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

Detailed Description

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

Definition at line 30 of file MonolithicHTFEM.h.

#include <MonolithicHTFEM.h>

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

Public Member Functions

 MonolithicHTFEM (MeshLib::Element const &element, std::size_t const local_matrix_size, NumLib::GenericIntegrationMethod const &integration_method, bool is_axially_symmetric, HTProcessData const &process_data)
void assemble (double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_x_prev, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data) override
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
Public Member Functions inherited from 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)
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
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< 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 NodalMatrixType = typename ShapeMatricesType::NodalMatrixType
using NodalRowVectorType = typename ShapeMatricesType::NodalRowVectorType
using GlobalDimVectorType = typename ShapeMatricesType::GlobalDimVectorType
using GlobalDimMatrixType = typename ShapeMatricesType::GlobalDimMatrixType

Additional Inherited Members

Protected Member Functions inherited from ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >
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 inherited from ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >
MeshLib::Element const & _element
HTProcessData const & _process_data
NumLib::GenericIntegrationMethod const & _integration_method
std::vector< IntegrationPointData< GlobalDimNodalMatrixType > > _ip_data
Static Protected Attributes inherited from ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >
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

Member Typedef Documentation

◆ GlobalDimMatrixType

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

Definition at line 47 of file MonolithicHTFEM.h.

◆ GlobalDimVectorType

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

Definition at line 46 of file MonolithicHTFEM.h.

◆ LocalMatrixType

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::HT::MonolithicHTFEM< 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 35 of file MonolithicHTFEM.h.

◆ LocalVectorType

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

Definition at line 38 of file MonolithicHTFEM.h.

◆ NodalMatrixType

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::HT::MonolithicHTFEM< ShapeFunction, GlobalDim >::NodalMatrixType = typename ShapeMatricesType::NodalMatrixType
private

Definition at line 43 of file MonolithicHTFEM.h.

◆ NodalRowVectorType

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

Definition at line 44 of file MonolithicHTFEM.h.

◆ NodalVectorType

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

Definition at line 42 of file MonolithicHTFEM.h.

◆ ShapeMatrices

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

Definition at line 33 of file MonolithicHTFEM.h.

◆ ShapeMatricesType

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

Definition at line 32 of file MonolithicHTFEM.h.

Constructor & Destructor Documentation

◆ MonolithicHTFEM()

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

Definition at line 50 of file MonolithicHTFEM.h.

58 {
59 }
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)
Definition HTFEM.h:41

References ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::HTFEM(), and ProcessLib::HT::NUM_NODAL_DOF.

Member Function Documentation

◆ assemble()

template<typename ShapeFunction, int GlobalDim>
void ProcessLib::HT::MonolithicHTFEM< 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 )
inlineoverridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 61 of file MonolithicHTFEM.h.

67 {
68 auto const local_matrix_size = local_x.size();
69 // This assertion is valid only if all nodal d.o.f. use the same shape
70 // matrices.
72
79
89
90 auto const& process_data = this->_process_data;
91
94
95 auto const& medium =
96 *process_data.media_map.getMedium(this->_element.getID());
97 auto const& liquid_phase =
99 auto const& solid_phase =
101
102 auto const& b =
104 .projected_specific_body_force_vectors[this->_element.getID()];
105
107
108 unsigned const n_integration_points =
109 this->_integration_method.getNumberOfPoints();
110
112 double average_velocity_norm = 0.0;
114
115 auto const& Ns =
116 process_data.shape_matrix_cache
118
119 for (unsigned ip(0); ip < n_integration_points; ip++)
120 {
121 auto const& ip_data = this->_ip_data[ip];
122 auto const& dNdx = ip_data.dNdx;
123 auto const& N = Ns[ip];
124 auto const& w = ip_data.integration_weight;
125
127 std::nullopt, this->_element.getID(),
131 this->_element, N))};
132
133 double T_int_pt = 0.0;
134 double p_int_pt = 0.0;
135 // Order matters: First T, then P!
137
138 vars.temperature = T_int_pt;
139 vars.liquid_phase_pressure = p_int_pt;
140
141 vars.liquid_saturation = 1.0;
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
148 auto const porosity =
150 .template value<double>(vars, pos, t, dt);
151 vars.porosity = porosity;
152
153 auto const intrinsic_permeability =
155 medium
156 .property(
158 .value(vars, pos, t, dt));
159
163 .template value<double>(vars, pos, t, dt);
164
165 // Use the fluid density model to compute the density
166 auto const fluid_density =
169 .template value<double>(vars, pos, t, dt);
170
171 vars.density = fluid_density;
172 const double dfluid_density_dp =
175 .template dValue<double>(
176 vars,
178 pos, t, dt);
179
180 // Use the viscosity model to compute the viscosity
181 auto const viscosity =
184 .template value<double>(vars, pos, t, dt);
186
188 process_data.has_gravity
190 fluid_density * b))
192
193 // matrix assembly
197 pos, t, dt);
198
199 KTT.noalias() +=
201
202 ip_flux_vector.emplace_back(velocity * fluid_density *
205
206 Kpp.noalias() += w * dNdx.transpose() * K_over_mu * dNdx;
207 MTT.noalias() += w *
211 N.transpose() * N;
212 Mpp.noalias() += w *
215 N.transpose() * N;
216 if (process_data.has_gravity)
217 {
218 Bp += w * fluid_density * dNdx.transpose() * K_over_mu * b;
219 }
220
221 if (process_data.has_fluid_thermal_expansion)
222 {
223 double const solid_thermal_expansion =
224 process_data.solid_thermal_expansion(t, pos)[0];
225 double const dfluid_density_dT =
228 .template dValue<double>(
230 pos, t, dt);
231 double const Tdot_int_pt =
232 (T_int_pt -
235 .dot(N)) /
236 dt;
237 double const biot_constant =
238 process_data.biot_constant(t, pos)[0];
239 double const eff_thermal_expansion =
242 Bp.noalias() += eff_thermal_expansion * Tdot_int_pt * w * N;
243 }
244 }
245
247 process_data.stabilizer, this->_ip_data,
248 process_data.shape_matrix_cache, ip_flux_vector,
249 average_velocity_norm / static_cast<double>(n_integration_points),
250 KTT);
251 }
NumLib::GenericIntegrationMethod const & _integration_method
Definition HTFEM.h:165
static const int temperature_index
Definition HTFEM.h:319
static const int temperature_size
Definition HTFEM.h:320
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)
Definition HTFEM.h:168
static const int pressure_size
Definition HTFEM.h:318
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)
Definition HTFEM.h:193
MeshLib::Element const & _element
Definition HTFEM.h:162
HTProcessData const & _process_data
Definition HTFEM.h:163
std::vector< IntegrationPointData< GlobalDimNodalMatrixType > > _ip_data
Definition HTFEM.h:166
static const int pressure_index
Definition HTFEM.h:317
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 &)
void assembleAdvectionMatrix(IPData const &ip_data_vector, NumLib::ShapeMatrixCache const &shape_matrix_cache, std::vector< FluxVectorType > const &ip_flux_vector, Eigen::MatrixBase< Derived > &laplacian_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, MaterialPropertyLib::AqueousLiquid, NumLib::detail::assembleAdvectionMatrix(), MathLib::createZeroedMatrix(), MathLib::createZeroedVector(), MaterialPropertyLib::density, MaterialPropertyLib::VariableArray::density, MaterialPropertyLib::formEigenTensor(), ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::getHeatEnergyCoefficient(), ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::getThermalConductivityDispersivity(), NumLib::interpolateCoordinates(), MaterialPropertyLib::liquid_phase_pressure, MaterialPropertyLib::VariableArray::liquid_phase_pressure, MaterialPropertyLib::VariableArray::liquid_saturation, ProcessLib::HT::NUM_NODAL_DOF, MaterialPropertyLib::permeability, MaterialPropertyLib::porosity, MaterialPropertyLib::VariableArray::porosity, ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::pressure_index, ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::pressure_size, NumLib::detail::shapeFunctionInterpolate(), MaterialPropertyLib::Solid, MaterialPropertyLib::specific_heat_capacity, MaterialPropertyLib::storage, MaterialPropertyLib::temperature, MaterialPropertyLib::VariableArray::temperature, ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::temperature_index, ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::temperature_size, and MaterialPropertyLib::viscosity.

◆ getIntPtDarcyVelocity()

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

Implements ProcessLib::HT::HTLocalAssemblerInterface.

Definition at line 253 of file MonolithicHTFEM.h.

258 {
259 int const process_id = 0; // monolithic case.
260 auto const indices =
262 assert(!indices.empty());
263 auto const& local_x = x[process_id]->get(indices);
264
266 }
std::vector< double > const & getIntPtDarcyVelocityLocal(const double t, std::vector< double > const &local_x, std::vector< double > &cache) const
Definition HTFEM.h:234
std::vector< GlobalIndexType > getIndices(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)

References ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::_element, NumLib::getIndices(), and ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::getIntPtDarcyVelocityLocal().


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