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 73 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 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_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 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 87 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 83 of file TwoPhaseFlowWithPPLocalAssembler.h.

◆ GlobalDimVectorType

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

Definition at line 88 of file TwoPhaseFlowWithPPLocalAssembler.h.

◆ LocalAssemblerTraits

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

◆ LocalMatrixType

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

Definition at line 89 of file TwoPhaseFlowWithPPLocalAssembler.h.

◆ LocalVectorType

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

Definition at line 90 of file TwoPhaseFlowWithPPLocalAssembler.h.

◆ NodalMatrixType

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

Definition at line 85 of file TwoPhaseFlowWithPPLocalAssembler.h.

◆ NodalRowVectorType

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

Definition at line 81 of file TwoPhaseFlowWithPPLocalAssembler.h.

◆ NodalVectorType

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

Definition at line 86 of file TwoPhaseFlowWithPPLocalAssembler.h.

◆ ShapeMatrices

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

Definition at line 77 of file TwoPhaseFlowWithPPLocalAssembler.h.

◆ ShapeMatricesType

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

Definition at line 76 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 93 of file TwoPhaseFlowWithPPLocalAssembler.h.

99 : _element(element),
100 _integration_method(integration_method),
101 _process_data(process_data),
103 std::vector<double>(_integration_method.getNumberOfPoints())),
105 std::vector<double>(_integration_method.getNumberOfPoints()))
106 {
107 unsigned const n_integration_points =
109 _ip_data.reserve(n_integration_points);
110 auto const shape_matrices =
112 GlobalDim>(element, is_axially_symmetric,
114 for (unsigned ip = 0; ip < n_integration_points; ip++)
115 {
116 auto const& sm = shape_matrices[ip];
117 _ip_data.emplace_back(
118 sm.N, sm.dNdx,
119 sm.integralMeasure * sm.detJ *
121 sm.N.transpose() * sm.N * sm.integralMeasure * sm.detJ *
123 }
124 }
double getWeight() const
MathLib::WeightedPoint const & getWeightedPoint(unsigned const igp) const
std::vector< IntegrationPointData< NodalRowVectorType, GlobalDimNodalMatrixType, NodalMatrixType >, Eigen::aligned_allocator< IntegrationPointData< NodalRowVectorType, GlobalDimNodalMatrixType, NodalMatrixType > > > _ip_data
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)

References ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPLocalAssembler< ShapeFunction, GlobalDim >::_integration_method, ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPLocalAssembler< ShapeFunction, GlobalDim >::_ip_data, NumLib::GenericIntegrationMethod::getNumberOfPoints(), MathLib::WeightedPoint::getWeight(), NumLib::GenericIntegrationMethod::getWeightedPoint(), 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 46 of file TwoPhaseFlowWithPPLocalAssembler-impl.h.

51{
52 auto const local_matrix_size = local_x.size();
53
54 assert(local_matrix_size == ShapeFunction::NPOINTS * NUM_NODAL_DOF);
55
56 auto local_M = MathLib::createZeroedMatrix<LocalMatrixType>(
57 local_M_data, local_matrix_size, local_matrix_size);
58 auto local_K = MathLib::createZeroedMatrix<LocalMatrixType>(
59 local_K_data, local_matrix_size, local_matrix_size);
60 auto local_b = MathLib::createZeroedVector<LocalVectorType>(
61 local_b_data, local_matrix_size);
62
63 auto Mgp =
64 local_M.template block<nonwet_pressure_size, nonwet_pressure_size>(
66 auto Mgpc = local_M.template block<nonwet_pressure_size, cap_pressure_size>(
68
69 auto Mlpc = local_M.template block<cap_pressure_size, cap_pressure_size>(
71
72 NodalMatrixType laplace_operator =
73 NodalMatrixType::Zero(ShapeFunction::NPOINTS, ShapeFunction::NPOINTS);
74
75 auto Kgp =
76 local_K.template block<nonwet_pressure_size, nonwet_pressure_size>(
78
79 auto Klp = local_K.template block<cap_pressure_size, nonwet_pressure_size>(
81
82 auto Klpc = local_K.template block<cap_pressure_size, cap_pressure_size>(
84
85 auto Bg = local_b.template segment<nonwet_pressure_size>(
87
88 auto Bl =
89 local_b.template segment<cap_pressure_size>(cap_pressure_matrix_index);
90
91 auto const& medium = *_process_data.media_map.getMedium(_element.getID());
92 auto const& liquid_phase = medium.phase("AqueousLiquid");
93 auto const& gas_phase = medium.phase("Gas");
94
95 unsigned const n_integration_points =
97
100
101 for (unsigned ip = 0; ip < n_integration_points; ip++)
102 {
103 pos.setIntegrationPoint(ip);
104 auto const& w = _ip_data[ip].integration_weight;
105 auto const& dNdx = _ip_data[ip].dNdx;
106 auto const& M = _ip_data[ip].massOperator;
107
108 double pc_int_pt = 0.;
109 double pn_int_pt = 0.;
110 NumLib::shapeFunctionInterpolate(local_x, _ip_data[ip].N, pn_int_pt,
111 pc_int_pt);
112
113 _pressure_wet[ip] = pn_int_pt - pc_int_pt;
114 MPL::VariableArray variables;
115
116 variables.gas_phase_pressure = pn_int_pt;
117 variables.capillary_pressure = pc_int_pt;
118 variables.temperature = _process_data.temperature(t, pos)[0];
119 variables.molar_mass =
120 gas_phase.property(MPL::PropertyType::molar_mass)
121 .template value<double>(variables, pos, t, dt);
122
123 auto const rho_nonwet =
124 gas_phase.property(MPL::PropertyType::density)
125 .template value<double>(variables, pos, t, dt);
126
127 auto const rho_wet = liquid_phase.property(MPL::PropertyType::density)
128 .template value<double>(variables, pos, t, dt);
129 auto& Sw = _saturation[ip];
130 Sw = medium.property(MPL::PropertyType::saturation)
131 .template value<double>(variables, pos, t, dt);
132
133 auto const dSw_dpc =
134 medium.property(MPL::PropertyType::saturation)
135 .template dValue<double>(
136 variables, MPL::Variable::capillary_pressure, pos, t, dt);
137
138 auto const porosity =
139 medium.property(MPL::PropertyType::porosity)
140 .template value<double>(variables, pos, t, dt);
141
142 auto const drhononwet_dpn =
143 gas_phase.property(MPL::PropertyType::density)
144 .template dValue<double>(
145 variables, MPL::Variable::gas_phase_pressure, pos, t, dt);
146
147 auto const k_rel_wet =
148 medium.property(MPL::PropertyType::relative_permeability)
149 .template value<double>(variables, pos, t, dt);
150 auto const k_rel_nonwet =
151 medium
152 .property(
153 MPL::PropertyType::relative_permeability_nonwetting_phase)
154 .template value<double>(variables, pos, t, dt);
155
156 auto const mu_nonwet =
157 gas_phase.property(MPL::PropertyType::viscosity)
158 .template value<double>(variables, pos, t, dt);
159
160 auto const lambda_nonwet = k_rel_nonwet / mu_nonwet;
161
162 auto const mu_wet = liquid_phase.property(MPL::PropertyType::viscosity)
163 .template value<double>(variables, pos, t, dt);
164
165 auto const lambda_wet = k_rel_wet / mu_wet;
166
167 auto const K = MPL::formEigenTensor<GlobalDim>(
168 medium.property(MPL::PropertyType::permeability)
169 .template value<double>(variables, pos, t, dt));
170
171 Mgp.noalias() += porosity * (1 - Sw) * drhononwet_dpn * M;
172 Mgpc.noalias() += -porosity * rho_nonwet * dSw_dpc * M;
173 Mlpc.noalias() += porosity * dSw_dpc * rho_wet * M;
174
175 laplace_operator.noalias() = dNdx.transpose() * K * dNdx * w;
176
177 Kgp.noalias() += rho_nonwet * lambda_nonwet * laplace_operator;
178 Klp.noalias() += rho_wet * lambda_wet * laplace_operator;
179 Klpc.noalias() += -rho_wet * lambda_wet * laplace_operator;
180
182 {
183 auto const& b = _process_data.specific_body_force;
184
185 NodalVectorType gravity_operator = dNdx.transpose() * K * b * w;
186 Bg.noalias() +=
187 rho_nonwet * rho_nonwet * lambda_nonwet * gravity_operator;
188 Bl.noalias() += rho_wet * rho_wet * lambda_wet * gravity_operator;
189 } // end of has gravity
190 }
192 {
193 for (unsigned row = 0; row < Mgpc.cols(); row++)
194 {
195 for (unsigned column = 0; column < Mgpc.cols(); column++)
196 {
197 if (row != column)
198 {
199 Mgpc(row, row) += Mgpc(row, column);
200 Mgpc(row, column) = 0.0;
201 Mgp(row, row) += Mgp(row, column);
202 Mgp(row, column) = 0.0;
203 Mlpc(row, row) += Mlpc(row, column);
204 Mlpc(row, column) = 0.0;
205 }
206 }
207 }
208 } // end of mass-lumping
209}
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
void setElementID(std::size_t element_id)
void setIntegrationPoint(unsigned integration_point)
void shapeFunctionInterpolate(const NodalValues &, const ShapeMatrix &)
MaterialPropertyLib::MaterialSpatialDistributionMap media_map

References MaterialPropertyLib::VariableArray::capillary_pressure, MaterialPropertyLib::VariableArray::gas_phase_pressure, MaterialPropertyLib::VariableArray::molar_mass, ProcessLib::TwoPhaseFlowWithPP::NUM_NODAL_DOF, ParameterLib::SpatialPosition::setElementID(), ParameterLib::SpatialPosition::setIntegrationPoint(), NumLib::detail::shapeFunctionInterpolate(), and MaterialPropertyLib::VariableArray::temperature.

◆ 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 133 of file TwoPhaseFlowWithPPLocalAssembler.h.

135 {
136 auto const& N = _ip_data[integration_point].N;
137
138 // assumes N is stored contiguously in memory
139 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
140 }

References ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPLocalAssembler< ShapeFunction, GlobalDim >::_ip_data.

Member Data Documentation

◆ _element

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

Definition at line 163 of file TwoPhaseFlowWithPPLocalAssembler.h.

◆ _integration_method

◆ _ip_data

◆ _pressure_wet

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

◆ _process_data

template<typename ShapeFunction , int GlobalDim>
TwoPhaseFlowWithPPProcessData const& ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPLocalAssembler< ShapeFunction, GlobalDim >::_process_data
private

Definition at line 167 of file TwoPhaseFlowWithPPLocalAssembler.h.

◆ _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 182 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 185 of file TwoPhaseFlowWithPPLocalAssembler.h.

◆ cap_pressure_size

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

Definition at line 188 of file TwoPhaseFlowWithPPLocalAssembler.h.

◆ 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 181 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 184 of file TwoPhaseFlowWithPPLocalAssembler.h.

◆ nonwet_pressure_size

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

Definition at line 187 of file TwoPhaseFlowWithPPLocalAssembler.h.


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