OGS
ProcessLib::TwoPhaseFlowWithPrho::TwoPhaseFlowWithPrhoLocalAssembler< ShapeFunction, GlobalDim > Class Template Reference

Detailed Description

template<typename ShapeFunction, int GlobalDim>
class ProcessLib::TwoPhaseFlowWithPrho::TwoPhaseFlowWithPrhoLocalAssembler< ShapeFunction, GlobalDim >

Definition at line 81 of file TwoPhaseFlowWithPrhoLocalAssembler.h.

#include <TwoPhaseFlowWithPrhoLocalAssembler.h>

Inheritance diagram for ProcessLib::TwoPhaseFlowWithPrho::TwoPhaseFlowWithPrhoLocalAssembler< ShapeFunction, GlobalDim >:
[legend]
Collaboration diagram for ProcessLib::TwoPhaseFlowWithPrho::TwoPhaseFlowWithPrhoLocalAssembler< ShapeFunction, GlobalDim >:
[legend]

Public Member Functions

 TwoPhaseFlowWithPrhoLocalAssembler (MeshLib::Element const &element, std::size_t const, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, TwoPhaseFlowWithPrhoProcessData const &process_data)
 
void assemble (double const t, double const, 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 & getIntPtNonWettingPressure (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_M_data, std::vector< double > &local_K_data, 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_M_data, std::vector< double > &local_K_data, 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.
 
- 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 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
 
std::vector< ShapeMatrices, Eigen::aligned_allocator< ShapeMatrices > > _shape_matrices
 
TwoPhaseFlowWithPrhoProcessData const & _process_data
 
std::vector< IntegrationPointData< NodalMatrixType >, Eigen::aligned_allocator< IntegrationPointData< NodalMatrixType > > > _ip_data
 
std::vector< double > _saturation
 
std::vector< double > _pressure_nonwetting
 used for secondary variable output
 

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::TwoPhaseFlowWithPrho::TwoPhaseFlowWithPrhoLocalAssembler< ShapeFunction, GlobalDim >::GlobalDimMatrixType = typename ShapeMatricesType::GlobalDimMatrixType
private

Definition at line 92 of file TwoPhaseFlowWithPrhoLocalAssembler.h.

◆ GlobalDimVectorType

template<typename ShapeFunction , int GlobalDim>
using ProcessLib::TwoPhaseFlowWithPrho::TwoPhaseFlowWithPrhoLocalAssembler< ShapeFunction, GlobalDim >::GlobalDimVectorType = typename ShapeMatricesType::GlobalDimVectorType
private

Definition at line 93 of file TwoPhaseFlowWithPrhoLocalAssembler.h.

◆ LocalAssemblerTraits

template<typename ShapeFunction , int GlobalDim>
using ProcessLib::TwoPhaseFlowWithPrho::TwoPhaseFlowWithPrhoLocalAssembler< ShapeFunction, GlobalDim >::LocalAssemblerTraits
private

◆ LocalMatrixType

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

Definition at line 94 of file TwoPhaseFlowWithPrhoLocalAssembler.h.

◆ LocalVectorType

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

Definition at line 95 of file TwoPhaseFlowWithPrhoLocalAssembler.h.

◆ NodalMatrixType

template<typename ShapeFunction , int GlobalDim>
using ProcessLib::TwoPhaseFlowWithPrho::TwoPhaseFlowWithPrhoLocalAssembler< ShapeFunction, GlobalDim >::NodalMatrixType = typename ShapeMatricesType::NodalMatrixType
private

Definition at line 90 of file TwoPhaseFlowWithPrhoLocalAssembler.h.

◆ NodalVectorType

template<typename ShapeFunction , int GlobalDim>
using ProcessLib::TwoPhaseFlowWithPrho::TwoPhaseFlowWithPrhoLocalAssembler< ShapeFunction, GlobalDim >::NodalVectorType = typename ShapeMatricesType::NodalVectorType
private

Definition at line 91 of file TwoPhaseFlowWithPrhoLocalAssembler.h.

◆ ShapeMatrices

template<typename ShapeFunction , int GlobalDim>
using ProcessLib::TwoPhaseFlowWithPrho::TwoPhaseFlowWithPrhoLocalAssembler< ShapeFunction, GlobalDim >::ShapeMatrices = typename ShapeMatricesType::ShapeMatrices
private

Definition at line 85 of file TwoPhaseFlowWithPrhoLocalAssembler.h.

◆ ShapeMatricesType

template<typename ShapeFunction , int GlobalDim>
using ProcessLib::TwoPhaseFlowWithPrho::TwoPhaseFlowWithPrhoLocalAssembler< ShapeFunction, GlobalDim >::ShapeMatricesType = ShapeMatrixPolicyType<ShapeFunction, GlobalDim>
private

Definition at line 84 of file TwoPhaseFlowWithPrhoLocalAssembler.h.

Constructor & Destructor Documentation

◆ TwoPhaseFlowWithPrhoLocalAssembler()

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

Definition at line 98 of file TwoPhaseFlowWithPrhoLocalAssembler.h.

104 : _element(element),
105 _integration_method(integration_method),
108 GlobalDim>(
109 element, is_axially_symmetric, _integration_method)),
110 _process_data(process_data),
112 std::vector<double>(_integration_method.getNumberOfPoints())),
114 std::vector<double>(_integration_method.getNumberOfPoints()))
115 {
116 unsigned const n_integration_points =
118 _ip_data.reserve(n_integration_points);
119 for (unsigned ip = 0; ip < n_integration_points; ip++)
120 {
121 _ip_data.emplace_back(*_process_data._material);
122 auto const& sm = _shape_matrices[ip];
123 _ip_data[ip].integration_weight =
124 sm.integralMeasure * sm.detJ *
126 _ip_data[ip].massOperator.setZero(ShapeFunction::NPOINTS,
127 ShapeFunction::NPOINTS);
128 _ip_data[ip].diffusionOperator.setZero(ShapeFunction::NPOINTS,
129 ShapeFunction::NPOINTS);
130 _ip_data[ip].massOperator.noalias() =
131 sm.N.transpose() * sm.N * _ip_data[ip].integration_weight;
132 _ip_data[ip].diffusionOperator.noalias() =
133 sm.dNdx.transpose() * sm.dNdx * _ip_data[ip].integration_weight;
134 }
135 }
double getWeight() const
MathLib::WeightedPoint const & getWeightedPoint(unsigned const igp) const
std::vector< IntegrationPointData< NodalMatrixType >, Eigen::aligned_allocator< IntegrationPointData< NodalMatrixType > > > _ip_data
std::vector< double > _pressure_nonwetting
used for secondary variable output
std::vector< ShapeMatrices, Eigen::aligned_allocator< ShapeMatrices > > _shape_matrices
std::vector< typename ShapeMatricesType::ShapeMatrices, Eigen::aligned_allocator< typename ShapeMatricesType::ShapeMatrices > > initShapeMatrices(MeshLib::Element const &e, bool const is_axially_symmetric, IntegrationMethod const &integration_method)
std::unique_ptr< TwoPhaseFlowWithPrhoMaterialProperties > _material

References ProcessLib::TwoPhaseFlowWithPrho::TwoPhaseFlowWithPrhoLocalAssembler< ShapeFunction, GlobalDim >::_integration_method, ProcessLib::TwoPhaseFlowWithPrho::TwoPhaseFlowWithPrhoLocalAssembler< ShapeFunction, GlobalDim >::_ip_data, ProcessLib::TwoPhaseFlowWithPrho::TwoPhaseFlowWithPrhoProcessData::_material, ProcessLib::TwoPhaseFlowWithPrho::TwoPhaseFlowWithPrhoLocalAssembler< ShapeFunction, GlobalDim >::_process_data, ProcessLib::TwoPhaseFlowWithPrho::TwoPhaseFlowWithPrhoLocalAssembler< ShapeFunction, GlobalDim >::_shape_matrices, NumLib::GenericIntegrationMethod::getNumberOfPoints(), MathLib::WeightedPoint::getWeight(), and NumLib::GenericIntegrationMethod::getWeightedPoint().

Member Function Documentation

◆ assemble()

template<typename ShapeFunction , int GlobalDim>
void ProcessLib::TwoPhaseFlowWithPrho::TwoPhaseFlowWithPrhoLocalAssembler< 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

Here only consider one component in gas phase

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 26 of file TwoPhaseFlowWithPrhoLocalAssembler-impl.h.

31{
32 auto const local_matrix_size = local_x.size();
33
34 assert(local_matrix_size == ShapeFunction::NPOINTS * NUM_NODAL_DOF);
35
36 auto local_M = MathLib::createZeroedMatrix<LocalMatrixType>(
37 local_M_data, local_matrix_size, local_matrix_size);
38 auto local_K = MathLib::createZeroedMatrix<LocalMatrixType>(
39 local_K_data, local_matrix_size, local_matrix_size);
40 auto local_b = MathLib::createZeroedVector<LocalVectorType>(
41 local_b_data, local_matrix_size);
42
43 auto Mgp =
44 local_M.template block<nonwet_pressure_size, nonwet_pressure_size>(
46 auto Mgx = local_M.template block<nonwet_pressure_size, cap_pressure_size>(
48
49 auto Mlp = local_M.template block<cap_pressure_size, nonwet_pressure_size>(
51
52 auto Mlx = local_M.template block<cap_pressure_size, cap_pressure_size>(
54
55 NodalMatrixType laplace_operator =
56 NodalMatrixType::Zero(ShapeFunction::NPOINTS, ShapeFunction::NPOINTS);
57
58 auto Kgp =
59 local_K.template block<nonwet_pressure_size, nonwet_pressure_size>(
61
62 auto Kgx = local_K.template block<nonwet_pressure_size, cap_pressure_size>(
64
65 auto Klp = local_K.template block<cap_pressure_size, nonwet_pressure_size>(
67
68 auto Klx = local_K.template block<cap_pressure_size, cap_pressure_size>(
70
71 auto Bg = local_b.template segment<nonwet_pressure_size>(
73
74 auto Bl =
75 local_b.template segment<cap_pressure_size>(cap_pressure_matrix_index);
76
77 auto const& medium = *_process_data.media_map.getMedium(_element.getID());
78 auto const& liquid_phase = medium.phase("AqueousLiquid");
79 auto const& gas_phase = medium.phase("Gas");
80
81 unsigned const n_integration_points =
83
86 const int material_id =
87 _process_data._material->getMaterialID(pos.getElementID().value());
88
89 for (unsigned ip = 0; ip < n_integration_points; ip++)
90 {
91 auto const& sm = _shape_matrices[ip];
92
93 double pl_int_pt = 0.;
94 double totalrho_int_pt =
95 0.; // total mass density of the light component
96 NumLib::shapeFunctionInterpolate(local_x, sm.N, pl_int_pt,
97 totalrho_int_pt);
98
99 const double temperature = _process_data._temperature(t, pos)[0];
100
101 MPL::VariableArray variables;
102
103 variables.liquid_phase_pressure = pl_int_pt;
104 variables.temperature = temperature;
105 variables.molar_mass =
106 gas_phase.property(MPL::PropertyType::molar_mass)
107 .template value<double>(variables, pos, t, dt);
108
109 auto const permeability = MPL::formEigenTensor<GlobalDim>(
110 medium.property(MPL::PropertyType::permeability)
111 .value(variables, pos, t, dt));
112 auto const rho_gas = gas_phase.property(MPL::PropertyType::density)
113 .template value<double>(variables, pos, t, dt);
114 auto const rho_h2o = liquid_phase.property(MPL::PropertyType::density)
115 .template value<double>(variables, pos, t, dt);
116
117 double& Sw = _ip_data[ip].sw;
119 double const X_h2_nonwet = 1.0;
120 double& rho_h2_wet = _ip_data[ip].rho_m;
121 double& dSwdP = _ip_data[ip].dsw_dpg;
122 double& dSwdrho = _ip_data[ip].dsw_drho;
123 double& drhoh2wet = _ip_data[ip].drhom_dpg;
124 double& drhoh2wet_drho = _ip_data[ip].drhom_drho;
125 if (!_ip_data[ip].mat_property.computeConstitutiveRelation(
126 t,
127 pos,
128 material_id,
129 pl_int_pt,
130 totalrho_int_pt,
131 temperature,
132 Sw,
133 rho_h2_wet,
134 dSwdP,
135 dSwdrho,
136 drhoh2wet,
137 drhoh2wet_drho))
138 {
139 OGS_FATAL("Computation of local constitutive relation failed.");
140 }
141 double const pc = _process_data._material->getCapillaryPressure(
142 material_id, t, pos, pl_int_pt, temperature, Sw);
143 variables.capillary_pressure = pc;
144 variables.liquid_saturation = Sw;
145
146 double const rho_wet = rho_h2o + rho_h2_wet;
147 _saturation[ip] = Sw;
149 pl_int_pt + pc;
150
151 // Assemble M matrix
152 // nonwetting
153 double dPC_dSw =
154 _process_data._material->getCapillaryPressureDerivative(
155 material_id, t, pos, pl_int_pt, temperature, Sw);
156
157 double const porosity =
158 medium.property(MPL::PropertyType::porosity)
159 .template value<double>(variables, pos, t, dt);
160
161 Mgx.noalias() += porosity * _ip_data[ip].massOperator;
162
163 Mlp.noalias() += porosity * rho_h2o * dSwdP * _ip_data[ip].massOperator;
164
165 Mlx.noalias() +=
166 porosity * (1 + dSwdrho * rho_h2o) * _ip_data[ip].massOperator;
167 double const k_rel_gas =
168 _process_data._material->getNonwetRelativePermeability(
169 t, pos, _pressure_nonwetting[ip], temperature, Sw);
170 double const mu_gas =
171 gas_phase.property(MPL::PropertyType::viscosity)
172 .template value<double>(variables, pos, t, dt);
173 double const lambda_gas = k_rel_gas / mu_gas;
174 double const diffusion_coeff_component_h2 =
176
177 // wet
178 double const k_rel_wet =
179 _process_data._material->getWetRelativePermeability(
180 t, pos, pl_int_pt, temperature, Sw);
181 auto const mu_wet = liquid_phase.property(MPL::PropertyType::viscosity)
182 .template value<double>(variables, pos, t, dt);
183 double const lambda_wet = k_rel_wet / mu_wet;
184
185 laplace_operator.noalias() = sm.dNdx.transpose() * permeability *
186 sm.dNdx * _ip_data[ip].integration_weight;
187
188 Kgp.noalias() +=
189 (rho_gas * X_h2_nonwet * lambda_gas * (1 + dPC_dSw * dSwdP) +
190 rho_h2_wet * lambda_wet) *
191 laplace_operator +
192 (Sw * porosity * diffusion_coeff_component_h2 *
193 (rho_h2o / rho_wet) * drhoh2wet) *
194 _ip_data[ip].diffusionOperator;
195 Kgx.noalias() +=
196 (rho_gas * X_h2_nonwet * lambda_gas * dPC_dSw * dSwdrho) *
197 laplace_operator +
198 (Sw * porosity * diffusion_coeff_component_h2 *
199 (rho_h2o / rho_wet) * drhoh2wet_drho) *
200 _ip_data[ip].diffusionOperator;
201 Klp.noalias() += (rho_gas * lambda_gas * (1 + dPC_dSw * dSwdP) +
202 rho_wet * lambda_wet) *
203 laplace_operator;
204
205 Klx.noalias() +=
206 (rho_gas * lambda_gas * dPC_dSw * dSwdrho) * laplace_operator;
207
209 {
210 auto const& b = _process_data._specific_body_force;
211 Bg.noalias() += (rho_gas * rho_gas * lambda_gas +
212 rho_h2_wet * rho_wet * lambda_wet) *
213 sm.dNdx.transpose() * permeability * b *
214 _ip_data[ip].integration_weight;
215 Bl.noalias() += (rho_wet * lambda_wet * rho_wet +
216 rho_gas * rho_gas * lambda_gas) *
217 sm.dNdx.transpose() * permeability * b *
218 _ip_data[ip].integration_weight;
219
220 } // end of has gravity
221 }
223 {
224 for (unsigned row = 0; row < Mgp.cols(); row++)
225 {
226 for (unsigned column = 0; column < Mgp.cols(); column++)
227 {
228 if (row != column)
229 {
230 Mgx(row, row) += Mgx(row, column);
231 Mgx(row, column) = 0.0;
232 Mgp(row, row) += Mgp(row, column);
233 Mgp(row, column) = 0.0;
234 Mlx(row, row) += Mlx(row, column);
235 Mlx(row, column) = 0.0;
236 Mlp(row, row) += Mlp(row, column);
237 Mlp(row, column) = 0.0;
238 }
239 }
240 }
241 } // end of mass-lumping
242}
#define OGS_FATAL(...)
Definition Error.h:26
Phase const & phase(std::size_t index) const
Definition Medium.cpp:33
std::size_t getID() const
Returns the ID of the element.
Definition Element.h:89
std::optional< std::size_t > getElementID() const
void setElementID(std::size_t element_id)
void shapeFunctionInterpolate(const NodalValues &, const ShapeMatrix &)

References MaterialPropertyLib::VariableArray::capillary_pressure, MaterialPropertyLib::VariableArray::gas_phase_pressure, ParameterLib::SpatialPosition::getElementID(), MaterialPropertyLib::VariableArray::liquid_phase_pressure, MaterialPropertyLib::VariableArray::liquid_saturation, MaterialPropertyLib::VariableArray::molar_mass, ProcessLib::TwoPhaseFlowWithPrho::NUM_NODAL_DOF, OGS_FATAL, ParameterLib::SpatialPosition::setElementID(), NumLib::detail::shapeFunctionInterpolate(), and MaterialPropertyLib::VariableArray::temperature.

◆ getIntPtNonWettingPressure()

template<typename ShapeFunction , int GlobalDim>
std::vector< double > const & ProcessLib::TwoPhaseFlowWithPrho::TwoPhaseFlowWithPrhoLocalAssembler< ShapeFunction, GlobalDim >::getIntPtNonWettingPressure ( const double ,
std::vector< GlobalVector * > const & ,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & ,
std::vector< double > &  ) const
inlineoverridevirtual

◆ getIntPtSaturation()

template<typename ShapeFunction , int GlobalDim>
std::vector< double > const & ProcessLib::TwoPhaseFlowWithPrho::TwoPhaseFlowWithPrhoLocalAssembler< ShapeFunction, GlobalDim >::getIntPtSaturation ( 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::TwoPhaseFlowWithPrho::TwoPhaseFlowWithPrhoLocalAssembler< ShapeFunction, GlobalDim >::getShapeMatrix ( const unsigned integration_point) const
inlineoverridevirtual

Provides the shape matrix at the given integration point.

Implements NumLib::ExtrapolatableElement.

Definition at line 144 of file TwoPhaseFlowWithPrhoLocalAssembler.h.

146 {
147 auto const& N = _shape_matrices[integration_point].N;
148
149 // assumes N is stored contiguously in memory
150 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
151 }

References ProcessLib::TwoPhaseFlowWithPrho::TwoPhaseFlowWithPrhoLocalAssembler< ShapeFunction, GlobalDim >::_shape_matrices.

Member Data Documentation

◆ _element

template<typename ShapeFunction , int GlobalDim>
MeshLib::Element const& ProcessLib::TwoPhaseFlowWithPrho::TwoPhaseFlowWithPrhoLocalAssembler< ShapeFunction, GlobalDim >::_element
private

Definition at line 174 of file TwoPhaseFlowWithPrhoLocalAssembler.h.

◆ _integration_method

◆ _ip_data

◆ _pressure_nonwetting

template<typename ShapeFunction , int GlobalDim>
std::vector<double> ProcessLib::TwoPhaseFlowWithPrho::TwoPhaseFlowWithPrhoLocalAssembler< ShapeFunction, GlobalDim >::_pressure_nonwetting
private

◆ _process_data

◆ _saturation

template<typename ShapeFunction , int GlobalDim>
std::vector<double> ProcessLib::TwoPhaseFlowWithPrho::TwoPhaseFlowWithPrhoLocalAssembler< ShapeFunction, GlobalDim >::_saturation
private

◆ _shape_matrices

◆ cap_pressure_coeff_index

template<typename ShapeFunction , int GlobalDim>
const int ProcessLib::TwoPhaseFlowWithPrho::TwoPhaseFlowWithPrhoLocalAssembler< ShapeFunction, GlobalDim >::cap_pressure_coeff_index = 1
staticprivate

Definition at line 188 of file TwoPhaseFlowWithPrhoLocalAssembler.h.

◆ cap_pressure_matrix_index

template<typename ShapeFunction , int GlobalDim>
const int ProcessLib::TwoPhaseFlowWithPrho::TwoPhaseFlowWithPrhoLocalAssembler< ShapeFunction, GlobalDim >::cap_pressure_matrix_index = ShapeFunction::NPOINTS
staticprivate

Definition at line 191 of file TwoPhaseFlowWithPrhoLocalAssembler.h.

◆ cap_pressure_size

template<typename ShapeFunction , int GlobalDim>
const int ProcessLib::TwoPhaseFlowWithPrho::TwoPhaseFlowWithPrhoLocalAssembler< ShapeFunction, GlobalDim >::cap_pressure_size = ShapeFunction::NPOINTS
staticprivate

Definition at line 194 of file TwoPhaseFlowWithPrhoLocalAssembler.h.

◆ nonwet_pressure_coeff_index

template<typename ShapeFunction , int GlobalDim>
const int ProcessLib::TwoPhaseFlowWithPrho::TwoPhaseFlowWithPrhoLocalAssembler< ShapeFunction, GlobalDim >::nonwet_pressure_coeff_index = 0
staticprivate

Definition at line 187 of file TwoPhaseFlowWithPrhoLocalAssembler.h.

◆ nonwet_pressure_matrix_index

template<typename ShapeFunction , int GlobalDim>
const int ProcessLib::TwoPhaseFlowWithPrho::TwoPhaseFlowWithPrhoLocalAssembler< ShapeFunction, GlobalDim >::nonwet_pressure_matrix_index = 0
staticprivate

Definition at line 190 of file TwoPhaseFlowWithPrhoLocalAssembler.h.

◆ nonwet_pressure_size

template<typename ShapeFunction , int GlobalDim>
const int ProcessLib::TwoPhaseFlowWithPrho::TwoPhaseFlowWithPrhoLocalAssembler< ShapeFunction, GlobalDim >::nonwet_pressure_size = ShapeFunction::NPOINTS
staticprivate

Definition at line 193 of file TwoPhaseFlowWithPrhoLocalAssembler.h.


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