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 77 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 std::optional< VectorSegmentgetVectorDeformationSegment () 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 97 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 94 of file RichardsComponentTransportFEM.h.

◆ GlobalDimVectorType

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

Definition at line 93 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 83 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 86 of file RichardsComponentTransportFEM.h.

◆ NodalMatrixType

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

Definition at line 96 of file RichardsComponentTransportFEM.h.

◆ NodalRowVectorType

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

Definition at line 91 of file RichardsComponentTransportFEM.h.

◆ NodalVectorType

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

Definition at line 90 of file RichardsComponentTransportFEM.h.

◆ ShapeMatrices

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

Definition at line 81 of file RichardsComponentTransportFEM.h.

◆ ShapeMatricesType

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

Definition at line 80 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 22 of file RichardsComponentTransportFEM-impl.h.

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

66{
67 auto const local_matrix_size = local_x.size();
68 // This assertion is valid only if all nodal d.o.f. use the same shape
69 // matrices.
71
78
79 unsigned const n_integration_points =
80 _integration_method.getNumberOfPoints();
81
84
85 auto const& b = _process_data.specific_body_force;
86
88
91
92 // Get material properties
93 auto const& medium = *_process_data.media_map.getMedium(_element.getID());
94 auto const& phase = medium.phase("AqueousLiquid");
95 auto const& component =
96 phase.component(_transport_process_variable.getName());
97
107
109 {
110 auto const& ip_data = _ip_data[ip];
111 auto const& N = ip_data.N;
112 auto const& dNdx = ip_data.dNdx;
113 auto const& w = ip_data.integration_weight;
114
116 std::nullopt, _element.getID(),
119 _element, N))};
120
121 double C_int_pt = 0.0;
122 double p_int_pt = 0.0;
123 // Order matters: First C, then p!
125
126 vars.capillary_pressure = -p_int_pt;
128 .template value<double>(vars, pos, t, dt);
129
130 double const dSw_dpc =
132 .template dValue<double>(
134 pos, t, dt);
135
136 vars.concentration = C_int_pt;
137 vars.liquid_phase_pressure = p_int_pt;
138 // setting pG to 1 atm
139 // TODO : rewrite equations s.t. p_L = pG-p_cap
140 vars.gas_phase_pressure = 1.0e5;
141
142 // \todo the argument to getValue() has to be changed for non
143 // constant storage model
144 auto const specific_storage =
146 .template value<double>(vars, pos, t, dt);
147 // \todo the first argument has to be changed for non constant
148 // porosity model
149 auto const porosity =
151 .template value<double>(vars, pos, t, dt);
152
153 auto const retardation_factor =
155 .template value<double>(vars, pos, t, dt);
156
159 .template value<double>(vars, pos, t, dt);
162 .template value<double>(vars, pos, t, dt);
163
164 // Use the fluid density model to compute the density
166 .template value<double>(vars, pos, t, dt);
167 vars.density = density;
168
169 auto const decay_rate =
171 .template value<double>(vars, pos, t, dt);
172 auto const pore_diffusion_coefficient =
175 .value(vars, pos, t, dt));
176
179 vars, pos, t, dt));
180 vars.liquid_saturation = Sw;
181 auto const k_rel =
183 .template value<double>(vars, pos, t, dt);
184 // Use the viscosity model to compute the viscosity
186 .template value<double>(vars, pos, t, dt);
187 auto const K_times_k_rel_over_mu = K * (k_rel / mu);
188
190 _process_data.has_gravity
195
196 double const velocity_magnitude = velocity.norm();
198 velocity_magnitude != 0.0
204 velocity_magnitude * velocity * velocity.transpose())
208
209 // matrix assembly
210 KCC.noalias() +=
211 (dNdx.transpose() * hydrodynamic_dispersion * dNdx +
212 N.transpose() * velocity.transpose() * dNdx +
213 N.transpose() * decay_rate * porosity * retardation_factor * N) *
214 w;
215 MCC.noalias() += w * N.transpose() * porosity * retardation_factor * N;
216 Kpp.noalias() += w * dNdx.transpose() * K_times_k_rel_over_mu * dNdx;
217 // \TODO Extend to pressure dependent density.
218 double const drhow_dp(0.0);
219 Mpp.noalias() += (specific_storage * Sw + porosity * Sw * drhow_dp -
220 porosity * dSw_dpc) *
221 ip_data.mass_operator;
222
223 if (_process_data.has_gravity)
224 {
225 Bp += w * density * dNdx.transpose() * K_times_k_rel_over_mu * b;
226 }
227 /* with Oberbeck-Boussing assumption density difference only exists
228 * in buoyancy effects */
229 }
230}
typename ShapeMatricesType::GlobalDimVectorType GlobalDimVectorType
ShapeMatrixPolicyType< ShapeFunction, GlobalDim > ShapeMatricesType
typename ShapeMatricesType::GlobalDimMatrixType GlobalDimMatrixType
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 234 of file RichardsComponentTransportFEM-impl.h.

239{
240 unsigned const n_integration_points =
241 _integration_method.getNumberOfPoints();
242
243 constexpr int process_id = 0; // monolithic scheme
244 auto const indices =
246 assert(!indices.empty());
247 auto const local_x = x[process_id]->get(indices);
248
249 cache.clear();
253
255
256 // Get material properties
257 auto const& medium = *_process_data.media_map.getMedium(_element.getID());
258 auto const& phase = medium.phase("AqueousLiquid");
259
262
263 for (unsigned ip = 0; ip < n_integration_points; ++ip)
264 {
265 auto const& ip_data = _ip_data[ip];
266 auto const& N = ip_data.N;
267 auto const& dNdx = ip_data.dNdx;
268
270 std::nullopt, _element.getID(),
273 _element, N))};
274
278 vars, pos, t, dt));
280 .template value<double>(vars, pos, t, dt);
281
282 double C_int_pt = 0.0;
283 double p_int_pt = 0.0;
285
286 // saturation
287 vars.capillary_pressure = -p_int_pt;
289 .template value<double>(vars, pos, t, dt);
290
291 vars.liquid_saturation = Sw;
292 auto const k_rel =
294 .template value<double>(vars, pos, t, dt);
295
296 cache_mat.col(ip).noalias() = -dNdx * p_nodal_values;
297 if (_process_data.has_gravity)
298 {
299 vars.concentration = C_int_pt;
300 vars.liquid_phase_pressure = p_int_pt;
302 .template value<double>(vars, pos, t, dt);
303 auto const b = _process_data.specific_body_force;
304 // here it is assumed that the vector b is directed 'downwards'
305 cache_mat.col(ip).noalias() += rho_w * b;
306 }
307 cache_mat.col(ip).noalias() = k_rel / mu * (K * cache_mat.col(ip));
308 }
309 return cache;
310}
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 325 of file RichardsComponentTransportFEM-impl.h.

330{
332
333 auto const& medium = *_process_data.media_map.getMedium(_element.getID());
334
335 unsigned const n_integration_points =
336 _integration_method.getNumberOfPoints();
337
338 constexpr int process_id = 0; // monolithic scheme
339 auto const indices =
341 assert(!indices.empty());
342 auto const local_x = x[process_id]->get(indices);
343
344 cache.clear();
348
349 for (unsigned ip = 0; ip < n_integration_points; ++ip)
350 {
351 auto const& ip_data = _ip_data[ip];
352 auto const& N = ip_data.N;
353
355 std::nullopt, _element.getID(),
358 _element, N))};
359
360 double C_int_pt = 0.0;
361 double p_int_pt = 0.0;
363
364 // saturation
365 vars.capillary_pressure = -p_int_pt;
368 .template value<double>(vars, pos, t, dt);
369 cache_vec[ip] = Sw;
370 }
371
372 return cache;
373}

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

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

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 135 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 144 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 145 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 146 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 147 of file RichardsComponentTransportFEM.h.

Referenced by assemble().


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