OGS
ProcessLib::RichardsComponentTransport::LocalAssemblerData< ShapeFunction, GlobalDim > Class Template Reference

Detailed Description

template<typename ShapeFunction, int GlobalDim>
class ProcessLib::RichardsComponentTransport::LocalAssemblerData< ShapeFunction, GlobalDim >

Definition at line 70 of file RichardsComponentTransportFEM.h.

#include <RichardsComponentTransportFEM.h>

Inheritance diagram for ProcessLib::RichardsComponentTransport::LocalAssemblerData< ShapeFunction, GlobalDim >:
[legend]
Collaboration diagram for ProcessLib::RichardsComponentTransport::LocalAssemblerData< ShapeFunction, GlobalDim >:
[legend]

Public Member Functions

 LocalAssemblerData (MeshLib::Element const &element, std::size_t const local_matrix_size, NumLib::GenericIntegrationMethod const &integration_method, bool is_axially_symmetric, RichardsComponentTransportProcessData const &process_data, ProcessVariable const &transport_process_variable)
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
std::vector< double > const & getIntPtDarcyVelocity (const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const 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 t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) 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 LocalMatrixType
using LocalVectorType
using NodalVectorType = typename ShapeMatricesType::NodalVectorType
using NodalRowVectorType = typename ShapeMatricesType::NodalRowVectorType
using GlobalDimVectorType = typename ShapeMatricesType::GlobalDimVectorType
using GlobalDimNodalMatrixType
using NodalMatrixType = typename ShapeMatricesType::NodalMatrixType
using GlobalDimMatrixType = typename ShapeMatricesType::GlobalDimMatrixType

Private Attributes

MeshLib::Element const & _element
RichardsComponentTransportProcessData const & _process_data
NumLib::GenericIntegrationMethod const & _integration_method
ProcessVariable const & _transport_process_variable
std::vector< IntegrationPointData< NodalRowVectorType, GlobalDimNodalMatrixType, NodalMatrixType >, Eigen::aligned_allocator< IntegrationPointData< NodalRowVectorType, GlobalDimNodalMatrixType, NodalMatrixType > > > _ip_data

Static Private Attributes

static const int concentration_index = 0
static const int concentration_size = ShapeFunction::NPOINTS
static const int pressure_index = ShapeFunction::NPOINTS
static const int pressure_size = ShapeFunction::NPOINTS

Member Typedef Documentation

◆ GlobalDimMatrixType

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::RichardsComponentTransport::LocalAssemblerData< ShapeFunction, GlobalDim >::GlobalDimMatrixType = typename ShapeMatricesType::GlobalDimMatrixType
private

Definition at line 90 of file RichardsComponentTransportFEM.h.

◆ GlobalDimNodalMatrixType

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

Definition at line 87 of file RichardsComponentTransportFEM.h.

◆ GlobalDimVectorType

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::RichardsComponentTransport::LocalAssemblerData< ShapeFunction, GlobalDim >::GlobalDimVectorType = typename ShapeMatricesType::GlobalDimVectorType
private

Definition at line 86 of file RichardsComponentTransportFEM.h.

◆ LocalMatrixType

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::RichardsComponentTransport::LocalAssemblerData< ShapeFunction, GlobalDim >::LocalMatrixType
private
Initial value:
typename ShapeMatricesType::template MatrixType<
NUM_NODAL_DOF * ShapeFunction::NPOINTS,
NUM_NODAL_DOF * ShapeFunction::NPOINTS>
const unsigned NUM_NODAL_DOF

Definition at line 76 of file RichardsComponentTransportFEM.h.

◆ LocalVectorType

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::RichardsComponentTransport::LocalAssemblerData< ShapeFunction, GlobalDim >::LocalVectorType
private
Initial value:
typename ShapeMatricesType::template VectorType<NUM_NODAL_DOF *
ShapeFunction::NPOINTS>

Definition at line 79 of file RichardsComponentTransportFEM.h.

◆ NodalMatrixType

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::RichardsComponentTransport::LocalAssemblerData< ShapeFunction, GlobalDim >::NodalMatrixType = typename ShapeMatricesType::NodalMatrixType
private

Definition at line 89 of file RichardsComponentTransportFEM.h.

◆ NodalRowVectorType

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::RichardsComponentTransport::LocalAssemblerData< ShapeFunction, GlobalDim >::NodalRowVectorType = typename ShapeMatricesType::NodalRowVectorType
private

Definition at line 84 of file RichardsComponentTransportFEM.h.

◆ NodalVectorType

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::RichardsComponentTransport::LocalAssemblerData< ShapeFunction, GlobalDim >::NodalVectorType = typename ShapeMatricesType::NodalVectorType
private

Definition at line 83 of file RichardsComponentTransportFEM.h.

◆ ShapeMatrices

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::RichardsComponentTransport::LocalAssemblerData< ShapeFunction, GlobalDim >::ShapeMatrices = typename ShapeMatricesType::ShapeMatrices
private

Definition at line 74 of file RichardsComponentTransportFEM.h.

◆ ShapeMatricesType

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::RichardsComponentTransport::LocalAssemblerData< ShapeFunction, GlobalDim >::ShapeMatricesType = ShapeMatrixPolicyType<ShapeFunction, GlobalDim>
private

Definition at line 73 of file RichardsComponentTransportFEM.h.

Constructor & Destructor Documentation

◆ LocalAssemblerData()

template<typename ShapeFunction, int GlobalDim>
ProcessLib::RichardsComponentTransport::LocalAssemblerData< ShapeFunction, GlobalDim >::LocalAssemblerData ( MeshLib::Element const & element,
std::size_t const local_matrix_size,
NumLib::GenericIntegrationMethod const & integration_method,
bool is_axially_symmetric,
RichardsComponentTransportProcessData const & process_data,
ProcessVariable const & transport_process_variable )

Definition at line 15 of file RichardsComponentTransportFEM-impl.h.

26{
27 // This assertion is valid only if all nodal d.o.f. use the same shape
28 // matrices.
31
32 unsigned const n_integration_points =
33 _integration_method.getNumberOfPoints();
35
36 auto const shape_matrices =
39
40 for (unsigned ip = 0; ip < n_integration_points; ip++)
41 {
42 auto const& sm = shape_matrices[ip];
43 double const integration_factor = sm.integralMeasure * sm.detJ;
44 _ip_data.emplace_back(
45 sm.N, sm.dNdx,
46 _integration_method.getWeightedPoint(ip).getWeight() *
48 sm.N.transpose() * sm.N * integration_factor *
49 _integration_method.getWeightedPoint(ip).getWeight());
50 }
51}
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 _element, _integration_method, _ip_data, _process_data, _transport_process_variable, NumLib::initShapeMatrices(), and ProcessLib::RichardsComponentTransport::NUM_NODAL_DOF.

Member Function Documentation

◆ assemble()

template<typename ShapeFunction, int GlobalDim>
void ProcessLib::RichardsComponentTransport::LocalAssemblerData< 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 54 of file RichardsComponentTransportFEM-impl.h.

59{
60 auto const local_matrix_size = local_x.size();
61 // This assertion is valid only if all nodal d.o.f. use the same shape
62 // matrices.
64
71
72 unsigned const n_integration_points =
73 _integration_method.getNumberOfPoints();
74
77
78 auto const& b = _process_data.specific_body_force;
79
81
84
85 // Get material properties
86 auto const& medium = *_process_data.media_map.getMedium(_element.getID());
87 auto const& phase = medium.phase("AqueousLiquid");
88 auto const& component =
89 phase.component(_transport_process_variable.getName());
90
100
102 {
103 auto const& ip_data = _ip_data[ip];
104 auto const& N = ip_data.N;
105 auto const& dNdx = ip_data.dNdx;
106 auto const& w = ip_data.integration_weight;
107
109 std::nullopt, _element.getID(),
112 _element, N))};
113
114 double C_int_pt = 0.0;
115 double p_int_pt = 0.0;
116 // Order matters: First C, then p!
118
119 vars.capillary_pressure = -p_int_pt;
121 .template value<double>(vars, pos, t, dt);
122
123 double const dSw_dpc =
125 .template dValue<double>(
127 pos, t, dt);
128
129 vars.concentration = C_int_pt;
130 vars.liquid_phase_pressure = p_int_pt;
131 // setting pG to 1 atm
132 // TODO : rewrite equations s.t. p_L = pG-p_cap
133 vars.gas_phase_pressure = 1.0e5;
134
135 // \todo the argument to getValue() has to be changed for non
136 // constant storage model
137 auto const specific_storage =
139 .template value<double>(vars, pos, t, dt);
140 // \todo the first argument has to be changed for non constant
141 // porosity model
142 auto const porosity =
144 .template value<double>(vars, pos, t, dt);
145
146 auto const retardation_factor =
148 .template value<double>(vars, pos, t, dt);
149
152 .template value<double>(vars, pos, t, dt);
155 .template value<double>(vars, pos, t, dt);
156
157 // Use the fluid density model to compute the density
159 .template value<double>(vars, pos, t, dt);
160 vars.density = density;
161
162 auto const decay_rate =
164 .template value<double>(vars, pos, t, dt);
165 auto const pore_diffusion_coefficient =
168 .value(vars, pos, t, dt));
169
172 vars, pos, t, dt));
173 vars.liquid_saturation = Sw;
174 auto const k_rel =
176 .template value<double>(vars, pos, t, dt);
177 // Use the viscosity model to compute the viscosity
179 .template value<double>(vars, pos, t, dt);
180 auto const K_times_k_rel_over_mu = K * (k_rel / mu);
181
183 _process_data.has_gravity
188
189 double const velocity_magnitude = velocity.norm();
191 velocity_magnitude != 0.0
197 velocity_magnitude * velocity * velocity.transpose())
201
202 // matrix assembly
203 KCC.noalias() +=
204 (dNdx.transpose() * hydrodynamic_dispersion * dNdx +
205 N.transpose() * velocity.transpose() * dNdx +
206 N.transpose() * decay_rate * porosity * retardation_factor * N) *
207 w;
208 MCC.noalias() += w * N.transpose() * porosity * retardation_factor * N;
209 Kpp.noalias() += w * dNdx.transpose() * K_times_k_rel_over_mu * dNdx;
210 // \TODO Extend to pressure dependent density.
211 double const drhow_dp(0.0);
212 Mpp.noalias() += (specific_storage * Sw + porosity * Sw * drhow_dp -
213 porosity * dSw_dpc) *
214 ip_data.mass_operator;
215
216 if (_process_data.has_gravity)
217 {
218 Bp += w * density * dNdx.transpose() * K_times_k_rel_over_mu * b;
219 }
220 /* with Oberbeck-Boussing assumption density difference only exists
221 * in buoyancy effects */
222 }
223}
typename ShapeMatricesType::GlobalDimVectorType GlobalDimVectorType
ShapeMatrixPolicyType< ShapeFunction, GlobalDim > ShapeMatricesType
typename ShapeMatricesType::GlobalDimMatrixType GlobalDimMatrixType
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, _process_data, _transport_process_variable, MaterialPropertyLib::capillary_pressure, MaterialPropertyLib::VariableArray::capillary_pressure, MaterialPropertyLib::VariableArray::concentration, concentration_index, MathLib::createZeroedMatrix(), MathLib::createZeroedVector(), MaterialPropertyLib::decay_rate, MaterialPropertyLib::density, MaterialPropertyLib::VariableArray::density, MaterialPropertyLib::formEigenTensor(), MaterialPropertyLib::VariableArray::gas_phase_pressure, NumLib::interpolateCoordinates(), MaterialPropertyLib::VariableArray::liquid_phase_pressure, MaterialPropertyLib::VariableArray::liquid_saturation, MaterialPropertyLib::longitudinal_dispersivity, ProcessLib::RichardsComponentTransport::NUM_NODAL_DOF, MaterialPropertyLib::permeability, MaterialPropertyLib::pore_diffusion, MaterialPropertyLib::porosity, pressure_index, pressure_size, MaterialPropertyLib::relative_permeability, MaterialPropertyLib::retardation_factor, MaterialPropertyLib::saturation, NumLib::detail::shapeFunctionInterpolate(), MaterialPropertyLib::storage, MaterialPropertyLib::transversal_dispersivity, and MaterialPropertyLib::viscosity.

◆ getIntPtDarcyVelocity()

template<typename ShapeFunction, int GlobalDim>
std::vector< double > const & ProcessLib::RichardsComponentTransport::LocalAssemblerData< ShapeFunction, GlobalDim >::getIntPtDarcyVelocity ( const double t,
std::vector< GlobalVector * > const & x,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & dof_table,
std::vector< double > & cache ) const
overridevirtual

Implements ProcessLib::RichardsComponentTransport::RichardsComponentTransportLocalAssemblerInterface.

Definition at line 227 of file RichardsComponentTransportFEM-impl.h.

232{
233 unsigned const n_integration_points =
234 _integration_method.getNumberOfPoints();
235
236 constexpr int process_id = 0; // monolithic scheme
237 auto const indices =
239 assert(!indices.empty());
240 auto const local_x = x[process_id]->get(indices);
241
242 cache.clear();
246
248
249 // Get material properties
250 auto const& medium = *_process_data.media_map.getMedium(_element.getID());
251 auto const& phase = medium.phase("AqueousLiquid");
252
255
256 for (unsigned ip = 0; ip < n_integration_points; ++ip)
257 {
258 auto const& ip_data = _ip_data[ip];
259 auto const& N = ip_data.N;
260 auto const& dNdx = ip_data.dNdx;
261
263 std::nullopt, _element.getID(),
266 _element, N))};
267
271 vars, pos, t, dt));
273 .template value<double>(vars, pos, t, dt);
274
275 double C_int_pt = 0.0;
276 double p_int_pt = 0.0;
278
279 // saturation
280 vars.capillary_pressure = -p_int_pt;
282 .template value<double>(vars, pos, t, dt);
283
284 vars.liquid_saturation = Sw;
285 auto const k_rel =
287 .template value<double>(vars, pos, t, dt);
288
289 cache_mat.col(ip).noalias() = -dNdx * p_nodal_values;
290 if (_process_data.has_gravity)
291 {
292 vars.concentration = C_int_pt;
293 vars.liquid_phase_pressure = p_int_pt;
295 .template value<double>(vars, pos, t, dt);
296 auto const b = _process_data.specific_body_force;
297 // here it is assumed that the vector b is directed 'downwards'
298 cache_mat.col(ip).noalias() += rho_w * b;
299 }
300 cache_mat.col(ip).noalias() = k_rel / mu * (K * cache_mat.col(ip));
301 }
302 return cache;
303}
std::vector< GlobalIndexType > getIndices(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)

References _element, _integration_method, _ip_data, _process_data, MathLib::createZeroedMatrix(), MaterialPropertyLib::density, MaterialPropertyLib::formEigenTensor(), NumLib::getIndices(), NumLib::interpolateCoordinates(), MaterialPropertyLib::permeability, MaterialPropertyLib::relative_permeability, MaterialPropertyLib::saturation, NumLib::detail::shapeFunctionInterpolate(), and MaterialPropertyLib::viscosity.

◆ getIntPtSaturation()

template<typename ShapeFunction, int GlobalDim>
std::vector< double > const & ProcessLib::RichardsComponentTransport::LocalAssemblerData< ShapeFunction, GlobalDim >::getIntPtSaturation ( const double t,
std::vector< GlobalVector * > const & x,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & dof_table,
std::vector< double > & cache ) const
overridevirtual

Implements ProcessLib::RichardsComponentTransport::RichardsComponentTransportLocalAssemblerInterface.

Definition at line 318 of file RichardsComponentTransportFEM-impl.h.

323{
325
326 auto const& medium = *_process_data.media_map.getMedium(_element.getID());
327
328 unsigned const n_integration_points =
329 _integration_method.getNumberOfPoints();
330
331 constexpr int process_id = 0; // monolithic scheme
332 auto const indices =
334 assert(!indices.empty());
335 auto const local_x = x[process_id]->get(indices);
336
337 cache.clear();
341
342 for (unsigned ip = 0; ip < n_integration_points; ++ip)
343 {
344 auto const& ip_data = _ip_data[ip];
345 auto const& N = ip_data.N;
346
348 std::nullopt, _element.getID(),
351 _element, N))};
352
353 double C_int_pt = 0.0;
354 double p_int_pt = 0.0;
356
357 // saturation
358 vars.capillary_pressure = -p_int_pt;
361 .template value<double>(vars, pos, t, dt);
362 cache_vec[ip] = Sw;
363 }
364
365 return cache;
366}

References _element, _integration_method, _ip_data, _process_data, MaterialPropertyLib::VariableArray::capillary_pressure, MathLib::createZeroedVector(), NumLib::getIndices(), NumLib::interpolateCoordinates(), MaterialPropertyLib::saturation, and NumLib::detail::shapeFunctionInterpolate().

◆ getShapeMatrix()

template<typename ShapeFunction, int GlobalDim>
Eigen::Map< const Eigen::RowVectorXd > ProcessLib::RichardsComponentTransport::LocalAssemblerData< ShapeFunction, GlobalDim >::getShapeMatrix ( const unsigned integration_point) const
overridevirtual

Provides the shape matrix at the given integration point.

Implements NumLib::ExtrapolatableElement.

Definition at line 307 of file RichardsComponentTransportFEM-impl.h.

309{
310 auto const& N = _ip_data[integration_point].N;
311
312 // assumes N is stored contiguously in memory
313 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
314}

References _ip_data.

Member Data Documentation

◆ _element

◆ _integration_method

◆ _ip_data

◆ _process_data

◆ _transport_process_variable

template<typename ShapeFunction, int GlobalDim>
ProcessVariable const& ProcessLib::RichardsComponentTransport::LocalAssemblerData< ShapeFunction, GlobalDim >::_transport_process_variable
private

Definition at line 128 of file RichardsComponentTransportFEM.h.

Referenced by LocalAssemblerData(), and assemble().

◆ concentration_index

template<typename ShapeFunction, int GlobalDim>
const int ProcessLib::RichardsComponentTransport::LocalAssemblerData< ShapeFunction, GlobalDim >::concentration_index = 0
staticprivate

Definition at line 137 of file RichardsComponentTransportFEM.h.

Referenced by assemble().

◆ concentration_size

template<typename ShapeFunction, int GlobalDim>
const int ProcessLib::RichardsComponentTransport::LocalAssemblerData< ShapeFunction, GlobalDim >::concentration_size = ShapeFunction::NPOINTS
staticprivate

Definition at line 138 of file RichardsComponentTransportFEM.h.

◆ pressure_index

template<typename ShapeFunction, int GlobalDim>
const int ProcessLib::RichardsComponentTransport::LocalAssemblerData< ShapeFunction, GlobalDim >::pressure_index = ShapeFunction::NPOINTS
staticprivate

Definition at line 139 of file RichardsComponentTransportFEM.h.

Referenced by assemble().

◆ pressure_size

template<typename ShapeFunction, int GlobalDim>
const int ProcessLib::RichardsComponentTransport::LocalAssemblerData< ShapeFunction, GlobalDim >::pressure_size = ShapeFunction::NPOINTS
staticprivate

Definition at line 140 of file RichardsComponentTransportFEM.h.

Referenced by assemble().


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