OGS
ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPLocalAssembler< ShapeFunction, GlobalDim > Class Template Reference

Detailed Description

template<typename ShapeFunction, int GlobalDim>
class ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPLocalAssembler< ShapeFunction, GlobalDim >

Definition at line 63 of file TwoPhaseFlowWithPPLocalAssembler.h.

#include <TwoPhaseFlowWithPPLocalAssembler.h>

Inheritance diagram for ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPLocalAssembler< ShapeFunction, GlobalDim >:
[legend]
Collaboration diagram for ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPLocalAssembler< ShapeFunction, GlobalDim >:
[legend]

Public Member Functions

 TwoPhaseFlowWithPPLocalAssembler (MeshLib::Element const &element, std::size_t const, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, TwoPhaseFlowWithPPProcessData const &process_data)
void setInitialConditionsConcrete (Eigen::VectorXd const local_x, double const t, int const process_id) override
void assemble (double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_x_prev, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data) override
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix (const unsigned integration_point) const override
 Provides the shape matrix at the given integration point.
std::vector< double > const & getIntPtSaturation (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &) const override
std::vector< double > const & getIntPtWetPressure (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &) const override
Public Member Functions inherited from ProcessLib::LocalAssemblerInterface
virtual ~LocalAssemblerInterface ()=default
virtual void setInitialConditions (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, double const t, int const process_id)
virtual void initialize (std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
virtual void preAssemble (double const, double const, std::vector< double > const &)
virtual void assembleForStaggeredScheme (double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, int const process_id, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data)
virtual void assembleWithJacobian (double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_x_prev, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data)
virtual void assembleWithJacobianForStaggeredScheme (double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, int const process_id, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data)
virtual void computeSecondaryVariable (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_prev, int const process_id)
virtual void preTimestep (std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, double const t, double const delta_t)
virtual void postTimestep (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, int const process_id)
void postNonLinearSolver (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, int const process_id)
virtual Eigen::Vector3d getFlux (MathLib::Point3d const &, double const, std::vector< double > const &) const
virtual Eigen::Vector3d getFlux (MathLib::Point3d const &, double const, std::vector< std::vector< double > > const &) const
 Fits to staggered scheme.
virtual int getNumberOfVectorElementsForDeformation () const
Public Member Functions inherited from NumLib::ExtrapolatableElement
virtual ~ExtrapolatableElement ()=default

Private Types

using ShapeMatricesType = ShapeMatrixPolicyType<ShapeFunction, GlobalDim>
using ShapeMatrices = typename ShapeMatricesType::ShapeMatrices
using LocalAssemblerTraits
using NodalRowVectorType = typename ShapeMatricesType::NodalRowVectorType
using GlobalDimNodalMatrixType
using NodalMatrixType = typename ShapeMatricesType::NodalMatrixType
using NodalVectorType = typename ShapeMatricesType::NodalVectorType
using GlobalDimMatrixType = typename ShapeMatricesType::GlobalDimMatrixType
using GlobalDimVectorType = typename ShapeMatricesType::GlobalDimVectorType
using LocalMatrixType = typename LocalAssemblerTraits::LocalMatrix
using LocalVectorType = typename LocalAssemblerTraits::LocalVector

Private Attributes

MeshLib::Element const & _element
NumLib::GenericIntegrationMethod const & _integration_method
TwoPhaseFlowWithPPProcessData const & _process_data
std::vector< IntegrationPointData< NodalRowVectorType, GlobalDimNodalMatrixType, NodalMatrixType >, Eigen::aligned_allocator< IntegrationPointData< NodalRowVectorType, GlobalDimNodalMatrixType, NodalMatrixType > > > _ip_data
std::vector< double > _saturation
std::vector< double > _pressure_wet

Static Private Attributes

static const int nonwet_pressure_coeff_index = 0
static const int cap_pressure_coeff_index = 1
static const int nonwet_pressure_matrix_index = 0
static const int cap_pressure_matrix_index = ShapeFunction::NPOINTS
static const int nonwet_pressure_size = ShapeFunction::NPOINTS
static const int cap_pressure_size = ShapeFunction::NPOINTS

Member Typedef Documentation

◆ GlobalDimMatrixType

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPLocalAssembler< ShapeFunction, GlobalDim >::GlobalDimMatrixType = typename ShapeMatricesType::GlobalDimMatrixType
private

Definition at line 77 of file TwoPhaseFlowWithPPLocalAssembler.h.

◆ GlobalDimNodalMatrixType

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPLocalAssembler< ShapeFunction, GlobalDim >::GlobalDimNodalMatrixType
private
Initial value:
MatrixType< GlobalDim, ShapeFunction::NPOINTS > GlobalDimNodalMatrixType

Definition at line 73 of file TwoPhaseFlowWithPPLocalAssembler.h.

◆ GlobalDimVectorType

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPLocalAssembler< ShapeFunction, GlobalDim >::GlobalDimVectorType = typename ShapeMatricesType::GlobalDimVectorType
private

Definition at line 78 of file TwoPhaseFlowWithPPLocalAssembler.h.

◆ LocalAssemblerTraits

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPLocalAssembler< ShapeFunction, GlobalDim >::LocalAssemblerTraits
private
Initial value:
ShapeMatricesType, ShapeFunction::NPOINTS, NUM_NODAL_DOF, GlobalDim>
ShapeMatrixPolicyType< ShapeFunction, GlobalDim > ShapeMatricesType
detail::LocalAssemblerTraitsFixed< ShpPol, NNodes, NodalDOF, Dim > LocalAssemblerTraits

Definition at line 69 of file TwoPhaseFlowWithPPLocalAssembler.h.

◆ LocalMatrixType

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPLocalAssembler< ShapeFunction, GlobalDim >::LocalMatrixType = typename LocalAssemblerTraits::LocalMatrix
private

Definition at line 79 of file TwoPhaseFlowWithPPLocalAssembler.h.

◆ LocalVectorType

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPLocalAssembler< ShapeFunction, GlobalDim >::LocalVectorType = typename LocalAssemblerTraits::LocalVector
private

Definition at line 80 of file TwoPhaseFlowWithPPLocalAssembler.h.

◆ NodalMatrixType

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPLocalAssembler< ShapeFunction, GlobalDim >::NodalMatrixType = typename ShapeMatricesType::NodalMatrixType
private

Definition at line 75 of file TwoPhaseFlowWithPPLocalAssembler.h.

◆ NodalRowVectorType

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPLocalAssembler< ShapeFunction, GlobalDim >::NodalRowVectorType = typename ShapeMatricesType::NodalRowVectorType
private

Definition at line 71 of file TwoPhaseFlowWithPPLocalAssembler.h.

◆ NodalVectorType

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPLocalAssembler< ShapeFunction, GlobalDim >::NodalVectorType = typename ShapeMatricesType::NodalVectorType
private

Definition at line 76 of file TwoPhaseFlowWithPPLocalAssembler.h.

◆ ShapeMatrices

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPLocalAssembler< ShapeFunction, GlobalDim >::ShapeMatrices = typename ShapeMatricesType::ShapeMatrices
private

Definition at line 67 of file TwoPhaseFlowWithPPLocalAssembler.h.

◆ ShapeMatricesType

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPLocalAssembler< ShapeFunction, GlobalDim >::ShapeMatricesType = ShapeMatrixPolicyType<ShapeFunction, GlobalDim>
private

Definition at line 66 of file TwoPhaseFlowWithPPLocalAssembler.h.

Constructor & Destructor Documentation

◆ TwoPhaseFlowWithPPLocalAssembler()

template<typename ShapeFunction, int GlobalDim>
ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPLocalAssembler< ShapeFunction, GlobalDim >::TwoPhaseFlowWithPPLocalAssembler ( MeshLib::Element const & element,
std::size_t const ,
NumLib::GenericIntegrationMethod const & integration_method,
bool const is_axially_symmetric,
TwoPhaseFlowWithPPProcessData const & process_data )
inline

Definition at line 83 of file TwoPhaseFlowWithPPLocalAssembler.h.

93 std::vector<double>(_integration_method.getNumberOfPoints())),
95 std::vector<double>(_integration_method.getNumberOfPoints()))
96 {
97 unsigned const n_integration_points =
98 _integration_method.getNumberOfPoints();
100 auto const shape_matrices =
104 for (unsigned ip = 0; ip < n_integration_points; ip++)
105 {
106 auto const& sm = shape_matrices[ip];
107 _ip_data.emplace_back(
108 sm.N, sm.dNdx,
109 sm.integralMeasure * sm.detJ *
110 _integration_method.getWeightedPoint(ip).getWeight());
111 }
112 }
std::vector< IntegrationPointData< NodalRowVectorType, GlobalDimNodalMatrixType, NodalMatrixType >, Eigen::aligned_allocator< IntegrationPointData< NodalRowVectorType, GlobalDimNodalMatrixType, NodalMatrixType > > > _ip_data

References _element, _integration_method, _ip_data, _pressure_wet, _process_data, _saturation, and NumLib::initShapeMatrices().

Member Function Documentation

◆ assemble()

template<typename ShapeFunction, int GlobalDim>
void ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPLocalAssembler< ShapeFunction, GlobalDim >::assemble ( double const t,
double const dt,
std::vector< double > const & local_x,
std::vector< double > const & local_x_prev,
std::vector< double > & local_M_data,
std::vector< double > & local_K_data,
std::vector< double > & local_b_data )
overridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 64 of file TwoPhaseFlowWithPPLocalAssembler-impl.h.

69{
70 auto const local_matrix_size = local_x.size();
71
73
80
81 auto Mgp =
86
89
92
93 auto Kgp =
96
99
102
105
106 auto Bl =
108
109 auto const& medium = *_process_data.media_map.getMedium(_element.getID());
110 auto const& liquid_phase = medium.phase("AqueousLiquid");
111 auto const& gas_phase = medium.phase("Gas");
112
113 unsigned const n_integration_points =
114 _integration_method.getNumberOfPoints();
115
117 pos.setElementID(_element.getID());
118
119 for (unsigned ip = 0; ip < n_integration_points; ip++)
120 {
121 auto const& w = _ip_data[ip].integration_weight;
122 auto const& N = _ip_data[ip].N;
123 auto const& dNdx = _ip_data[ip].dNdx;
124
125 NodalMatrixType M = N.transpose() * N * w;
126
127 double pc_int_pt = 0.;
128 double pn_int_pt = 0.;
130 pc_int_pt);
131
134
135 variables.gas_phase_pressure = pn_int_pt;
136 variables.capillary_pressure = pc_int_pt;
137 variables.temperature = _process_data.temperature(t, pos)[0];
138 variables.molar_mass =
140 .template value<double>(variables, pos, t, dt);
141
142 auto const rho_nonwet =
144 .template value<double>(variables, pos, t, dt);
145
146 auto const rho_wet = liquid_phase.property(MPL::PropertyType::density)
147 .template value<double>(variables, pos, t, dt);
148 auto& Sw = _saturation[ip];
150 .template value<double>(variables, pos, t, dt);
151
152 variables.liquid_saturation = Sw;
153 auto const dSw_dpc =
155 .template dValue<double>(
157
158 auto const porosity =
160 .template value<double>(variables, pos, t, dt);
161
162 auto const drhononwet_dpn =
164 .template dValue<double>(
166
167 auto const k_rel_wet =
169 .template value<double>(variables, pos, t, dt);
170 auto const k_rel_nonwet =
171 medium
172 .property(
174 .template value<double>(variables, pos, t, dt);
175
176 auto const mu_nonwet =
178 .template value<double>(variables, pos, t, dt);
179
180 auto const lambda_nonwet = k_rel_nonwet / mu_nonwet;
181
183 .template value<double>(variables, pos, t, dt);
184
185 auto const lambda_wet = k_rel_wet / mu_wet;
186
189 .template value<double>(variables, pos, t, dt));
190
191 Mgp.noalias() += porosity * (1 - Sw) * drhononwet_dpn * M;
192 Mgpc.noalias() += -porosity * rho_nonwet * dSw_dpc * M;
193 Mlpc.noalias() += porosity * dSw_dpc * rho_wet * M;
194
195 laplace_operator.noalias() = dNdx.transpose() * K * dNdx * w;
196
198 Klp.noalias() += rho_wet * lambda_wet * laplace_operator;
199 Klpc.noalias() += -rho_wet * lambda_wet * laplace_operator;
200
201 if (_process_data.has_gravity)
202 {
203 auto const& b = _process_data.specific_body_force;
204
205 NodalVectorType gravity_operator = dNdx.transpose() * K * b * w;
206 Bg.noalias() +=
208 Bl.noalias() += rho_wet * rho_wet * lambda_wet * gravity_operator;
209 } // end of has gravity
210 }
211 if (_process_data.has_mass_lumping)
212 {
213 for (unsigned row = 0; row < Mgpc.cols(); row++)
214 {
215 for (unsigned column = 0; column < Mgpc.cols(); column++)
216 {
217 if (row != column)
218 {
219 Mgpc(row, row) += Mgpc(row, column);
220 Mgpc(row, column) = 0.0;
221 Mgp(row, row) += Mgp(row, column);
222 Mgp(row, column) = 0.0;
223 Mlpc(row, row) += Mlpc(row, column);
224 Mlpc(row, column) = 0.0;
225 }
226 }
227 }
228 } // end of mass-lumping
229}
constexpr Eigen::Matrix< double, GlobalDim, GlobalDim > formEigenTensor(MaterialPropertyLib::PropertyDataType const &values)
Eigen::Map< Vector > createZeroedVector(std::vector< double > &data, Eigen::VectorXd::Index size)
Eigen::Map< Matrix > createZeroedMatrix(std::vector< double > &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
void shapeFunctionInterpolate(const NodalValues &, const ShapeMatrix &)

References _element, _integration_method, _ip_data, _pressure_wet, _process_data, _saturation, cap_pressure_matrix_index, MaterialPropertyLib::capillary_pressure, MaterialPropertyLib::VariableArray::capillary_pressure, MathLib::createZeroedMatrix(), MathLib::createZeroedVector(), MaterialPropertyLib::density, MaterialPropertyLib::formEigenTensor(), MaterialPropertyLib::gas_phase_pressure, MaterialPropertyLib::VariableArray::gas_phase_pressure, MaterialPropertyLib::VariableArray::liquid_saturation, MaterialPropertyLib::molar_mass, MaterialPropertyLib::VariableArray::molar_mass, nonwet_pressure_matrix_index, ProcessLib::TwoPhaseFlowWithPP::NUM_NODAL_DOF, MaterialPropertyLib::permeability, MaterialPropertyLib::porosity, MaterialPropertyLib::relative_permeability, MaterialPropertyLib::relative_permeability_nonwetting_phase, MaterialPropertyLib::saturation, ParameterLib::SpatialPosition::setElementID(), NumLib::detail::shapeFunctionInterpolate(), MaterialPropertyLib::VariableArray::temperature, and MaterialPropertyLib::viscosity.

◆ getIntPtSaturation()

template<typename ShapeFunction, int GlobalDim>
std::vector< double > const & ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPLocalAssembler< ShapeFunction, GlobalDim >::getIntPtSaturation ( const double ,
std::vector< GlobalVector * > const & ,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & ,
std::vector< double > &  ) const
inlineoverridevirtual

◆ getIntPtWetPressure()

template<typename ShapeFunction, int GlobalDim>
std::vector< double > const & ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPLocalAssembler< ShapeFunction, GlobalDim >::getIntPtWetPressure ( const double ,
std::vector< GlobalVector * > const & ,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & ,
std::vector< double > &  ) const
inlineoverridevirtual

◆ getShapeMatrix()

template<typename ShapeFunction, int GlobalDim>
Eigen::Map< const Eigen::RowVectorXd > ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPLocalAssembler< ShapeFunction, GlobalDim >::getShapeMatrix ( const unsigned integration_point) const
inlineoverridevirtual

Provides the shape matrix at the given integration point.

Implements NumLib::ExtrapolatableElement.

Definition at line 125 of file TwoPhaseFlowWithPPLocalAssembler.h.

127 {
128 auto const& N = _ip_data[integration_point].N;
129
130 // assumes N is stored contiguously in memory
131 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
132 }

References _ip_data.

◆ setInitialConditionsConcrete()

template<typename ShapeFunction, int GlobalDim>
void ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPLocalAssembler< ShapeFunction, GlobalDim >::setInitialConditionsConcrete ( Eigen::VectorXd const local_x,
double const t,
int const process_id )
overridevirtual

Member Data Documentation

◆ _element

◆ _integration_method

◆ _ip_data

◆ _pressure_wet

template<typename ShapeFunction, int GlobalDim>
std::vector<double> ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPLocalAssembler< ShapeFunction, GlobalDim >::_pressure_wet
private

◆ _process_data

◆ _saturation

template<typename ShapeFunction, int GlobalDim>
std::vector<double> ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPLocalAssembler< ShapeFunction, GlobalDim >::_saturation
private

◆ cap_pressure_coeff_index

template<typename ShapeFunction, int GlobalDim>
const int ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPLocalAssembler< ShapeFunction, GlobalDim >::cap_pressure_coeff_index = 1
staticprivate

Definition at line 174 of file TwoPhaseFlowWithPPLocalAssembler.h.

◆ cap_pressure_matrix_index

template<typename ShapeFunction, int GlobalDim>
const int ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPLocalAssembler< ShapeFunction, GlobalDim >::cap_pressure_matrix_index = ShapeFunction::NPOINTS
staticprivate

Definition at line 177 of file TwoPhaseFlowWithPPLocalAssembler.h.

Referenced by assemble(), and setInitialConditionsConcrete().

◆ cap_pressure_size

template<typename ShapeFunction, int GlobalDim>
const int ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPLocalAssembler< ShapeFunction, GlobalDim >::cap_pressure_size = ShapeFunction::NPOINTS
staticprivate

Definition at line 180 of file TwoPhaseFlowWithPPLocalAssembler.h.

Referenced by setInitialConditionsConcrete().

◆ nonwet_pressure_coeff_index

template<typename ShapeFunction, int GlobalDim>
const int ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPLocalAssembler< ShapeFunction, GlobalDim >::nonwet_pressure_coeff_index = 0
staticprivate

Definition at line 173 of file TwoPhaseFlowWithPPLocalAssembler.h.

◆ nonwet_pressure_matrix_index

template<typename ShapeFunction, int GlobalDim>
const int ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPLocalAssembler< ShapeFunction, GlobalDim >::nonwet_pressure_matrix_index = 0
staticprivate

Definition at line 176 of file TwoPhaseFlowWithPPLocalAssembler.h.

Referenced by assemble().

◆ nonwet_pressure_size

template<typename ShapeFunction, int GlobalDim>
const int ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPLocalAssembler< ShapeFunction, GlobalDim >::nonwet_pressure_size = ShapeFunction::NPOINTS
staticprivate

Definition at line 179 of file TwoPhaseFlowWithPPLocalAssembler.h.


The documentation for this class was generated from the following files: