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 16 of file RichardsComponentTransportFEM-impl.h.

27{
28 // This assertion is valid only if all nodal d.o.f. use the same shape
29 // matrices.
32
33 unsigned const n_integration_points =
34 _integration_method.getNumberOfPoints();
36
37 auto const shape_matrices =
40
41 for (unsigned ip = 0; ip < n_integration_points; ip++)
42 {
43 auto const& sm = shape_matrices[ip];
44 double const integration_factor = sm.integralMeasure * sm.detJ;
45 _ip_data.emplace_back(
46 sm.N, sm.dNdx,
47 _integration_method.getWeightedPoint(ip).getWeight() *
49 sm.N.transpose() * sm.N * integration_factor *
50 _integration_method.getWeightedPoint(ip).getWeight());
51 }
52}
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 55 of file RichardsComponentTransportFEM-impl.h.

60{
61 auto const local_matrix_size = local_x.size();
62 // This assertion is valid only if all nodal d.o.f. use the same shape
63 // matrices.
65
72
73 unsigned const n_integration_points =
74 _integration_method.getNumberOfPoints();
75
78
79 auto const& b = _process_data.specific_body_force;
80
82
85
86 // Get material properties
87 auto const& medium = *_process_data.media_map.getMedium(_element.getID());
88 auto const& phase =
90 auto const& component =
91 phase.component(_transport_process_variable.getName());
92
102
104 {
105 auto const& ip_data = _ip_data[ip];
106 auto const& N = ip_data.N;
107 auto const& dNdx = ip_data.dNdx;
108 auto const& w = ip_data.integration_weight;
109
111 std::nullopt, _element.getID(),
114 _element, N))};
115
116 double C_int_pt = 0.0;
117 double p_int_pt = 0.0;
118 // Order matters: First C, then p!
120
121 vars.capillary_pressure = -p_int_pt;
123 .template value<double>(vars, pos, t, dt);
124
125 double const dSw_dpc =
127 .template dValue<double>(
129 pos, t, dt);
130
131 vars.concentration = C_int_pt;
132 vars.liquid_phase_pressure = p_int_pt;
133 // setting pG to 1 atm
134 // TODO : rewrite equations s.t. p_L = pG-p_cap
135 vars.gas_phase_pressure = 1.0e5;
136
137 // \todo the argument to getValue() has to be changed for non
138 // constant storage model
139 auto const specific_storage =
141 .template value<double>(vars, pos, t, dt);
142 // \todo the first argument has to be changed for non constant
143 // porosity model
144 auto const porosity =
146 .template value<double>(vars, pos, t, dt);
147
148 auto const retardation_factor =
150 .template value<double>(vars, pos, t, dt);
151
154 .template value<double>(vars, pos, t, dt);
157 .template value<double>(vars, pos, t, dt);
158
159 // Use the fluid density model to compute the density
161 .template value<double>(vars, pos, t, dt);
162 vars.density = density;
163
164 auto const decay_rate =
166 .template value<double>(vars, pos, t, dt);
167 auto const pore_diffusion_coefficient =
170 .value(vars, pos, t, dt));
171
174 vars, pos, t, dt));
175 vars.liquid_saturation = Sw;
176 auto const k_rel =
178 .template value<double>(vars, pos, t, dt);
179 // Use the viscosity model to compute the viscosity
181 .template value<double>(vars, pos, t, dt);
182 auto const K_times_k_rel_over_mu = K * (k_rel / mu);
183
185 _process_data.has_gravity
190
191 double const velocity_magnitude = velocity.norm();
193 velocity_magnitude != 0.0
199 velocity_magnitude * velocity * velocity.transpose())
203
204 // matrix assembly
205 KCC.noalias() +=
206 (dNdx.transpose() * hydrodynamic_dispersion * dNdx +
207 N.transpose() * velocity.transpose() * dNdx +
208 N.transpose() * decay_rate * porosity * retardation_factor * N) *
209 w;
210 MCC.noalias() += w * N.transpose() * porosity * retardation_factor * N;
211 Kpp.noalias() += w * dNdx.transpose() * K_times_k_rel_over_mu * dNdx;
212 // \TODO Extend to pressure dependent density.
213 double const drhow_dp(0.0);
214 Mpp.noalias() += (specific_storage * Sw + porosity * Sw * drhow_dp -
215 porosity * dSw_dpc) *
216 ip_data.mass_operator;
217
218 if (_process_data.has_gravity)
219 {
220 Bp += w * density * dNdx.transpose() * K_times_k_rel_over_mu * b;
221 }
222 /* with Oberbeck-Boussing assumption density difference only exists
223 * in buoyancy effects */
224 }
225}
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::AqueousLiquid, 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 229 of file RichardsComponentTransportFEM-impl.h.

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

References _element, _integration_method, _ip_data, _process_data, MaterialPropertyLib::AqueousLiquid, 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 321 of file RichardsComponentTransportFEM-impl.h.

326{
328
329 auto const& medium = *_process_data.media_map.getMedium(_element.getID());
330
331 unsigned const n_integration_points =
332 _integration_method.getNumberOfPoints();
333
334 constexpr int process_id = 0; // monolithic scheme
335 auto const indices =
337 assert(!indices.empty());
338 auto const local_x = x[process_id]->get(indices);
339
340 cache.clear();
344
345 for (unsigned ip = 0; ip < n_integration_points; ++ip)
346 {
347 auto const& ip_data = _ip_data[ip];
348 auto const& N = ip_data.N;
349
351 std::nullopt, _element.getID(),
354 _element, N))};
355
356 double C_int_pt = 0.0;
357 double p_int_pt = 0.0;
359
360 // saturation
361 vars.capillary_pressure = -p_int_pt;
364 .template value<double>(vars, pos, t, dt);
365 cache_vec[ip] = Sw;
366 }
367
368 return cache;
369}

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 310 of file RichardsComponentTransportFEM-impl.h.

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

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: