OGS
ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim > Class Template Reference

Detailed Description

template<typename ShapeFunction, int GlobalDim>
class ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >

Definition at line 34 of file HTFEM.h.

#include <HTFEM.h>

Inheritance diagram for ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >:
[legend]
Collaboration diagram for ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >:
[legend]

Public Member Functions

 HTFEM (MeshLib::Element const &element, std::size_t const local_matrix_size, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, HTProcessData const &process_data, const unsigned dof_per_node)
 
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix (const unsigned integration_point) const override
 Provides the shape matrix at the given integration point.
 
Eigen::Vector3d getFlux (MathLib::Point3d const &pnt_local_coords, double const t, std::vector< double > const &local_x) const override
 
- Public Member Functions inherited from ProcessLib::HT::HTLocalAssemblerInterface
 HTLocalAssemblerInterface ()=default
 
virtual std::vector< double > const & getIntPtDarcyVelocity (const double, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &) const =0
 
- 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 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)
 
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< std::vector< double > > const &) const
 Fits to staggered scheme.
 
- Public Member Functions inherited from NumLib::ExtrapolatableElement
virtual ~ExtrapolatableElement ()=default
 

Protected Member Functions

double getHeatEnergyCoefficient (MaterialPropertyLib::VariableArray const &vars, const double porosity, const double fluid_density, const double specific_heat_capacity_fluid, ParameterLib::SpatialPosition const &pos, double const t, double const dt)
 
GlobalDimMatrixType getThermalConductivityDispersivity (MaterialPropertyLib::VariableArray const &vars, const double fluid_density, const double specific_heat_capacity_fluid, const GlobalDimVectorType &velocity, ParameterLib::SpatialPosition const &pos, double const t, double const dt)
 
std::vector< double > const & getIntPtDarcyVelocityLocal (const double t, std::vector< double > const &local_x, std::vector< double > &cache) const
 

Protected Attributes

MeshLib::Element const & _element
 
HTProcessData const & _process_data
 
NumLib::GenericIntegrationMethod const & _integration_method
 
std::vector< IntegrationPointData< GlobalDimNodalMatrixType > > _ip_data
 

Static Protected Attributes

static const int pressure_index = ShapeFunction::NPOINTS
 
static const int pressure_size = ShapeFunction::NPOINTS
 
static const int temperature_index = 0
 
static const int temperature_size = ShapeFunction::NPOINTS
 

Private Types

using ShapeMatricesType = ShapeMatrixPolicyType<ShapeFunction, GlobalDim>
 
using ShapeMatrices = typename ShapeMatricesType::ShapeMatrices
 
using NodalVectorType = typename ShapeMatricesType::NodalVectorType
 
using NodalRowVectorType = typename ShapeMatricesType::NodalRowVectorType
 
using GlobalDimVectorType = typename ShapeMatricesType::GlobalDimVectorType
 
using GlobalDimNodalMatrixType
 
using GlobalDimMatrixType = typename ShapeMatricesType::GlobalDimMatrixType
 

Member Typedef Documentation

◆ GlobalDimMatrixType

template<typename ShapeFunction , int GlobalDim>
using ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::GlobalDimMatrixType = typename ShapeMatricesType::GlobalDimMatrixType
private

Definition at line 45 of file HTFEM.h.

◆ GlobalDimNodalMatrixType

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

Definition at line 43 of file HTFEM.h.

◆ GlobalDimVectorType

template<typename ShapeFunction , int GlobalDim>
using ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::GlobalDimVectorType = typename ShapeMatricesType::GlobalDimVectorType
private

Definition at line 42 of file HTFEM.h.

◆ NodalRowVectorType

template<typename ShapeFunction , int GlobalDim>
using ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::NodalRowVectorType = typename ShapeMatricesType::NodalRowVectorType
private

Definition at line 40 of file HTFEM.h.

◆ NodalVectorType

template<typename ShapeFunction , int GlobalDim>
using ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::NodalVectorType = typename ShapeMatricesType::NodalVectorType
private

Definition at line 39 of file HTFEM.h.

◆ ShapeMatrices

template<typename ShapeFunction , int GlobalDim>
using ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::ShapeMatrices = typename ShapeMatricesType::ShapeMatrices
private

Definition at line 37 of file HTFEM.h.

◆ ShapeMatricesType

template<typename ShapeFunction , int GlobalDim>
using ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::ShapeMatricesType = ShapeMatrixPolicyType<ShapeFunction, GlobalDim>
private

Definition at line 36 of file HTFEM.h.

Constructor & Destructor Documentation

◆ HTFEM()

template<typename ShapeFunction , int GlobalDim>
ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::HTFEM ( MeshLib::Element const & element,
std::size_t const local_matrix_size,
NumLib::GenericIntegrationMethod const & integration_method,
bool const is_axially_symmetric,
HTProcessData const & process_data,
const unsigned dof_per_node )
inline

Definition at line 48 of file HTFEM.h.

55 _element(element),
56 _process_data(process_data),
57 _integration_method(integration_method)
58 {
59 // This assertion is valid only if all nodal d.o.f. use the same shape
60 // matrices.
61 assert(local_matrix_size == ShapeFunction::NPOINTS * dof_per_node);
62 (void)local_matrix_size;
63 (void)dof_per_node;
64
65 unsigned const n_integration_points =
67 _ip_data.reserve(n_integration_points);
68
71
72 double const aperture_size = _process_data.aperture_size(0.0, pos)[0];
73
74 auto const shape_matrices =
76 GlobalDim>(element, is_axially_symmetric,
78
79 for (unsigned ip = 0; ip < n_integration_points; ip++)
80 {
81 _ip_data.emplace_back(
82 shape_matrices[ip].dNdx,
84 shape_matrices[ip].integralMeasure *
85 shape_matrices[ip].detJ * aperture_size);
86 }
87 }
double getWeight() const
std::size_t getID() const
Returns the ID of the element.
Definition Element.h:89
MathLib::WeightedPoint const & getWeightedPoint(unsigned const igp) const
void setElementID(std::size_t element_id)
NumLib::GenericIntegrationMethod const & _integration_method
Definition HTFEM.h:171
ShapeMatrixPolicyType< ShapeFunction, GlobalDim > ShapeMatricesType
Definition HTFEM.h:36
MeshLib::Element const & _element
Definition HTFEM.h:168
HTProcessData const & _process_data
Definition HTFEM.h:169
std::vector< IntegrationPointData< GlobalDimNodalMatrixType > > _ip_data
Definition HTFEM.h:172
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)
ParameterLib::Parameter< double > const & aperture_size

References ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::_element, ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::_integration_method, ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::_ip_data, ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::_process_data, ProcessLib::HT::HTProcessData::aperture_size, MeshLib::Element::getID(), NumLib::GenericIntegrationMethod::getNumberOfPoints(), MathLib::WeightedPoint::getWeight(), NumLib::GenericIntegrationMethod::getWeightedPoint(), NumLib::initShapeMatrices(), and ParameterLib::SpatialPosition::setElementID().

Member Function Documentation

◆ getFlux()

template<typename ShapeFunction , int GlobalDim>
Eigen::Vector3d ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::getFlux ( MathLib::Point3d const & pnt_local_coords,
double const t,
std::vector< double > const & local_x ) const
inlineoverridevirtual

Computes the flux in the point pnt_local_coords that is given in local coordinates using the values from local_x.

Implements ProcessLib::HT::HTLocalAssemblerInterface.

Definition at line 101 of file HTFEM.h.

104 {
105 // Eval shape matrices at given point
106 // Note: Axial symmetry is set to false here, because we only need dNdx
107 // here, which is not affected by axial symmetry.
108 auto const shape_matrices =
110 GlobalDim>(
111 _element, false /*is_axially_symmetric*/,
112 std::array{pnt_local_coords})[0];
113
115 pos.setElementID(this->_element.getID());
116
118
119 // local_x contains the local temperature and pressure values
120 double T_int_pt = 0.0;
121 double p_int_pt = 0.0;
122 NumLib::shapeFunctionInterpolate(local_x, shape_matrices.N, T_int_pt,
123 p_int_pt);
124
125 vars.temperature = T_int_pt;
126 vars.liquid_phase_pressure = p_int_pt;
127
128 auto const& medium =
130 auto const& liquid_phase = medium.phase("AqueousLiquid");
131
132 // TODO (naumov) Temporary value not used by current material models.
133 // Need extension of secondary variables interface.
134 double const dt = std::numeric_limits<double>::quiet_NaN();
135 // fetch permeability, viscosity, density
138 .value(vars, pos, t, dt));
139
140 auto const mu =
142 .template value<double>(vars, pos, t, dt);
143 GlobalDimMatrixType const K_over_mu = K / mu;
144
145 auto const p_nodal_values = Eigen::Map<const NodalVectorType>(
146 &local_x[local_x.size() / 2], ShapeFunction::NPOINTS);
148 -K_over_mu * shape_matrices.dNdx * p_nodal_values;
149
150 if (this->_process_data.has_gravity)
151 {
152 auto const rho_w =
153 liquid_phase
155 .template value<double>(vars, pos, t, dt);
156 auto const b =
158 [this->_element.getID()];
159 q += K_over_mu * rho_w * b;
160 }
161
162 Eigen::Vector3d flux;
163 flux.head<GlobalDim>() = q;
164 return flux;
165 }
Phase const & phase(std::size_t index) const
Definition Medium.cpp:33
typename ShapeMatricesType::GlobalDimMatrixType GlobalDimMatrixType
Definition HTFEM.h:45
typename ShapeMatricesType::GlobalDimVectorType GlobalDimVectorType
Definition HTFEM.h:42
Eigen::Matrix< double, GlobalDim, GlobalDim > formEigenTensor(MaterialPropertyLib::PropertyDataType const &values)
void shapeFunctionInterpolate(const NodalValues &, const ShapeMatrix &)
std::vector< typename ShapeMatricesType::ShapeMatrices, Eigen::aligned_allocator< typename ShapeMatricesType::ShapeMatrices > > computeShapeMatrices(MeshLib::Element const &e, bool const is_axially_symmetric, PointContainer const &points)
MaterialPropertyLib::MaterialSpatialDistributionMap media_map
std::vector< Eigen::VectorXd > const projected_specific_body_force_vectors
Projected specific body force vector: R * R^T * b.

References ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::_element, ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::_process_data, NumLib::computeShapeMatrices(), MaterialPropertyLib::density, MaterialPropertyLib::formEigenTensor(), MeshLib::Element::getID(), MaterialPropertyLib::MaterialSpatialDistributionMap::getMedium(), ProcessLib::HT::HTProcessData::has_gravity, MaterialPropertyLib::VariableArray::liquid_phase_pressure, ProcessLib::HT::HTProcessData::media_map, MaterialPropertyLib::permeability, MaterialPropertyLib::Medium::phase(), ProcessLib::HT::HTProcessData::projected_specific_body_force_vectors, ParameterLib::SpatialPosition::setElementID(), NumLib::detail::shapeFunctionInterpolate(), MaterialPropertyLib::VariableArray::temperature, and MaterialPropertyLib::viscosity.

◆ getHeatEnergyCoefficient()

template<typename ShapeFunction , int GlobalDim>
double ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::getHeatEnergyCoefficient ( MaterialPropertyLib::VariableArray const & vars,
const double porosity,
const double fluid_density,
const double specific_heat_capacity_fluid,
ParameterLib::SpatialPosition const & pos,
double const t,
double const dt )
inlineprotected

Definition at line 174 of file HTFEM.h.

179 {
180 auto const& medium =
182 auto const& solid_phase = medium.phase("Solid");
183
184 auto const specific_heat_capacity_solid =
185 solid_phase
186 .property(
188 .template value<double>(vars, pos, t, dt);
189
190 auto const solid_density =
192 .template value<double>(vars, pos, t, dt);
193
194 return solid_density * specific_heat_capacity_solid * (1 - porosity) +
195 fluid_density * specific_heat_capacity_fluid * porosity;
196 }
Property const & property(PropertyType const &p) const
Definition Phase.cpp:53

References ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::_process_data, MaterialPropertyLib::density, MeshLib::Element::getID(), MaterialPropertyLib::MaterialSpatialDistributionMap::getMedium(), ProcessLib::HT::HTProcessData::media_map, MaterialPropertyLib::Medium::phase(), MaterialPropertyLib::Phase::property(), and MaterialPropertyLib::specific_heat_capacity.

Referenced by ProcessLib::HT::MonolithicHTFEM< ShapeFunction, GlobalDim >::assemble().

◆ getIntPtDarcyVelocityLocal()

template<typename ShapeFunction , int GlobalDim>
std::vector< double > const & ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::getIntPtDarcyVelocityLocal ( const double t,
std::vector< double > const & local_x,
std::vector< double > & cache ) const
inlineprotected

Definition at line 239 of file HTFEM.h.

242 {
243 std::vector<double> local_p{
244 local_x.data() + pressure_index,
245 local_x.data() + pressure_index + pressure_size};
246 std::vector<double> local_T{
247 local_x.data() + temperature_index,
248 local_x.data() + temperature_index + temperature_size};
249
250 auto const n_integration_points =
252
253 cache.clear();
254 auto cache_mat = MathLib::createZeroedMatrix<
255 Eigen::Matrix<double, GlobalDim, Eigen::Dynamic, Eigen::RowMajor>>(
256 cache, GlobalDim, n_integration_points);
257
260
262
263 auto const p_nodal_values = Eigen::Map<const NodalVectorType>(
264 &local_p[0], ShapeFunction::NPOINTS);
265
266 auto const& medium =
268 auto const& liquid_phase = medium.phase("AqueousLiquid");
269
270 auto const& Ns =
272 .NsHigherOrder<typename ShapeFunction::MeshElement>();
273
274 for (unsigned ip = 0; ip < n_integration_points; ++ip)
275 {
276 auto const& ip_data = _ip_data[ip];
277 auto const& dNdx = ip_data.dNdx;
278 auto const& N = Ns[ip];
279
280 pos.setIntegrationPoint(ip);
281
282 double T_int_pt = 0.0;
283 double p_int_pt = 0.0;
284 NumLib::shapeFunctionInterpolate(local_p, N, p_int_pt);
285 NumLib::shapeFunctionInterpolate(local_T, N, T_int_pt);
286
287 vars.temperature = T_int_pt;
288 vars.liquid_phase_pressure = p_int_pt;
289
290 // TODO (naumov) Temporary value not used by current material
291 // models. Need extension of secondary variables interface.
292 double const dt = std::numeric_limits<double>::quiet_NaN();
295 .value(vars, pos, t, dt));
296
297 auto const mu =
298 liquid_phase
300 .template value<double>(vars, pos, t, dt);
301 GlobalDimMatrixType const K_over_mu = K / mu;
302
303 cache_mat.col(ip).noalias() = -K_over_mu * dNdx * p_nodal_values;
304
306 {
307 auto const rho_w =
308 liquid_phase
310 .template value<double>(vars, pos, t, dt);
311 auto const b =
313 [_element.getID()];
314 // here it is assumed that the vector b is directed 'downwards'
315 cache_mat.col(ip).noalias() += K_over_mu * rho_w * b;
316 }
317 }
318
319 return cache;
320 }
auto const & NsHigherOrder() const
static const int temperature_index
Definition HTFEM.h:325
static const int temperature_size
Definition HTFEM.h:326
static const int pressure_size
Definition HTFEM.h:324
static const int pressure_index
Definition HTFEM.h:323
Eigen::Map< Matrix > createZeroedMatrix(std::vector< double > &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
NumLib::ShapeMatrixCache shape_matrix_cache
caches for each mesh element type the shape matrix

References ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::_element, ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::_integration_method, ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::_ip_data, ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::_process_data, MathLib::createZeroedMatrix(), MaterialPropertyLib::density, MaterialPropertyLib::formEigenTensor(), MeshLib::Element::getID(), MaterialPropertyLib::MaterialSpatialDistributionMap::getMedium(), NumLib::GenericIntegrationMethod::getNumberOfPoints(), ProcessLib::HT::HTProcessData::has_gravity, MaterialPropertyLib::VariableArray::liquid_phase_pressure, ProcessLib::HT::HTProcessData::media_map, NumLib::ShapeMatrixCache::NsHigherOrder(), MaterialPropertyLib::permeability, MaterialPropertyLib::Medium::phase(), ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::pressure_index, ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::pressure_size, ProcessLib::HT::HTProcessData::projected_specific_body_force_vectors, ParameterLib::SpatialPosition::setElementID(), ProcessLib::HT::HTProcessData::shape_matrix_cache, NumLib::detail::shapeFunctionInterpolate(), MaterialPropertyLib::VariableArray::temperature, ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::temperature_index, ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::temperature_size, and MaterialPropertyLib::viscosity.

Referenced by ProcessLib::HT::MonolithicHTFEM< ShapeFunction, GlobalDim >::getIntPtDarcyVelocity().

◆ getShapeMatrix()

template<typename ShapeFunction , int GlobalDim>
Eigen::Map< const Eigen::RowVectorXd > ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::getShapeMatrix ( const unsigned integration_point) const
inlineoverridevirtual

Provides the shape matrix at the given integration point.

Implements NumLib::ExtrapolatableElement.

Definition at line 89 of file HTFEM.h.

91 {
93 typename ShapeFunction::MeshElement>()[integration_point];
94
95 // assumes N is stored contiguously in memory
96 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
97 }

References ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::_process_data, NumLib::ShapeMatrixCache::NsHigherOrder(), and ProcessLib::HT::HTProcessData::shape_matrix_cache.

◆ getThermalConductivityDispersivity()

template<typename ShapeFunction , int GlobalDim>
GlobalDimMatrixType ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::getThermalConductivityDispersivity ( MaterialPropertyLib::VariableArray const & vars,
const double fluid_density,
const double specific_heat_capacity_fluid,
const GlobalDimVectorType & velocity,
ParameterLib::SpatialPosition const & pos,
double const t,
double const dt )
inlineprotected

Definition at line 198 of file HTFEM.h.

204 {
205 auto const& medium =
207
210 medium
211 .property(
213 .value(vars, pos, t, dt));
214
215 auto const thermal_dispersivity_transversal =
216 medium
218 thermal_transversal_dispersivity)
219 .template value<double>();
220
221 auto const thermal_dispersivity_longitudinal =
222 medium
224 thermal_longitudinal_dispersivity)
225 .template value<double>();
226
227 // Thermal conductivity is moved outside and zero matrix is passed
228 // instead due to multiplication with fluid's density times specific
229 // heat capacity.
230 return thermal_conductivity +
231 fluid_density * specific_heat_capacity_fluid *
234 GlobalDimMatrixType::Zero(GlobalDim, GlobalDim),
235 velocity, 0 /* phi */, thermal_dispersivity_transversal,
236 thermal_dispersivity_longitudinal);
237 }
Eigen::MatrixXd computeHydrodynamicDispersion(NumericalStabilization const &stabilizer, std::size_t const element_id, Eigen::MatrixXd const &pore_diffusion_coefficient, Eigen::VectorXd const &velocity, double const porosity, double const solute_dispersivity_transverse, double const solute_dispersivity_longitudinal)
double fluid_density(const double p, const double T, const double x)
NumLib::NumericalStabilization stabilizer

References ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::_element, ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::_process_data, NumLib::computeHydrodynamicDispersion(), MaterialPropertyLib::formEigenTensor(), MeshLib::Element::getID(), MaterialPropertyLib::MaterialSpatialDistributionMap::getMedium(), ProcessLib::HT::HTProcessData::media_map, ProcessLib::HT::HTProcessData::stabilizer, and MaterialPropertyLib::thermal_conductivity.

Referenced by ProcessLib::HT::MonolithicHTFEM< ShapeFunction, GlobalDim >::assemble().

Member Data Documentation

◆ _element

◆ _integration_method

◆ _ip_data

◆ _process_data

◆ pressure_index

template<typename ShapeFunction , int GlobalDim>
const int ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::pressure_index = ShapeFunction::NPOINTS
staticprotected

◆ pressure_size

template<typename ShapeFunction , int GlobalDim>
const int ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::pressure_size = ShapeFunction::NPOINTS
staticprotected

◆ temperature_index

template<typename ShapeFunction , int GlobalDim>
const int ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::temperature_index = 0
staticprotected

◆ temperature_size

template<typename ShapeFunction , int GlobalDim>
const int ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::temperature_size = ShapeFunction::NPOINTS
staticprotected

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