OGS
ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits > Class Template Reference

Detailed Description

template<typename ShapeFunctionDisplacement, typename ShapeFunction, int DisplacementDim, typename ConstitutiveTraits>
class ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >

Definition at line 31 of file ThermoRichardsMechanicsFEM.h.

#include <ThermoRichardsMechanicsFEM.h>

Inheritance diagram for ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >:
[legend]
Collaboration diagram for ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >:
[legend]

Classes

class  LocalMatrices

Public Types

using ShapeMatricesTypeDisplacement
using ShapeMatricesType
using GlobalDimMatrixType = typename ShapeMatricesType::GlobalDimMatrixType
using GlobalDimVectorType = typename ShapeMatricesType::GlobalDimVectorType
using BMatricesType
using KelvinVectorType = typename BMatricesType::KelvinVectorType
using IpData
using Invariants = MathLib::KelvinVector::Invariants<KelvinVectorSize>
using SymmetricTensor = Eigen::Matrix<double, KelvinVectorSize, 1>

Public Member Functions

 ThermoRichardsMechanicsLocalAssembler (ThermoRichardsMechanicsLocalAssembler const &)=delete
 ThermoRichardsMechanicsLocalAssembler (ThermoRichardsMechanicsLocalAssembler &&)=delete
 ThermoRichardsMechanicsLocalAssembler (MeshLib::Element const &e, std::size_t const, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, ThermoRichardsMechanicsProcessData< DisplacementDim, ConstitutiveTraits > &process_data)
void setInitialConditionsConcrete (Eigen::VectorXd const local_x, double const t, int const process_id) override
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_rhs_data, std::vector< double > &local_Jac_data) override
int getNumberOfVectorElementsForDeformation () const override
void initializeConcrete () override
void computeSecondaryVariableConcrete (double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev) override
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix (const unsigned integration_point) const override
 Provides the shape matrix at the given integration point.
Public Member Functions inherited from ProcessLib::ThermoRichardsMechanics::LocalAssemblerInterface< DisplacementDim, ConstitutiveTraits >
 LocalAssemblerInterface (MeshLib::Element const &e, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, ThermoRichardsMechanicsProcessData< DisplacementDim, ConstitutiveTraits > &process_data)
std::size_t setIPDataInitialConditions (std::string_view name, double const *values, int const integration_order)
std::vector< double > getMaterialStateVariableInternalState (std::function< std::span< double >(typename MaterialLib::Solids::MechanicsBase< DisplacementDim >::MaterialStateVariables &)> const &get_values_span, int const &n_components) const
unsigned getNumberOfIntegrationPoints () const
int getMaterialID () const
MaterialLib::Solids::MechanicsBase< DisplacementDim >::MaterialStateVariables const & getMaterialStateVariablesAt (unsigned integration_point) const
void postTimestepConcrete (Eigen::VectorXd const &, Eigen::VectorXd const &, double const, double const, int 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 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 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.
Public Member Functions inherited from NumLib::ExtrapolatableElement
virtual ~ExtrapolatableElement ()=default

Static Public Attributes

static int const KelvinVectorSize
static constexpr auto & N_u_op

Private Member Functions

void assembleWithJacobianSingleIP (double const t, double const dt, ParameterLib::SpatialPosition const &x_position, std::vector< double > const &local_x, std::vector< double > const &local_x_prev, IpData const &ip_data, typename ConstitutiveTraits::ConstitutiveSetting &CS, MaterialPropertyLib::Medium &medium, LocalMatrices &out, typename ConstitutiveTraits::StatefulData &current_state, typename ConstitutiveTraits::StatefulDataPrev const &prev_state, MaterialStateData< DisplacementDim > &mat_state, typename ConstitutiveTraits::OutputData &output_data) const
void addToLocalMatrixData (double const dt, std::vector< double > const &local_x, std::vector< double > const &local_x_prev, LocalMatrices const &loc_mat, std::vector< double > &local_rhs_data, std::vector< double > &local_Jac_data) const
void massLumping (LocalMatrices &loc_mat) const
void convertInitialStressType (unsigned const ip, double const t, ParameterLib::SpatialPosition const x_position, MaterialPropertyLib::Medium const &medium, MPL::VariableArray const &variables, double const p_at_ip)

Static Private Member Functions

static constexpr auto localDOF (std::vector< double > const &x)
static auto block_uu (auto &mat)
static auto block_up (auto &mat)
static auto block_uT (auto &mat)
static auto block_pu (auto &mat)
static auto block_pp (auto &mat)
static auto block_pT (auto &mat)
static auto block_Tp (auto &mat)
static auto block_TT (auto &mat)
static auto block_u (auto &vec)
static auto block_p (auto &vec)
static auto block_T (auto &vec)

Private Attributes

std::vector< IpDataip_data_

Static Private Attributes

static constexpr int temperature_index = 0
static constexpr int temperature_size = ShapeFunction::NPOINTS
static constexpr int pressure_index = temperature_size
static constexpr int pressure_size = ShapeFunction::NPOINTS
static constexpr int displacement_index = 2 * ShapeFunction::NPOINTS
static constexpr int displacement_size

Additional Inherited Members

Static Public Member Functions inherited from ProcessLib::ThermoRichardsMechanics::LocalAssemblerInterface< DisplacementDim, ConstitutiveTraits >
static auto getReflectionDataForOutput ()
Protected Attributes inherited from ProcessLib::ThermoRichardsMechanics::LocalAssemblerInterface< DisplacementDim, ConstitutiveTraits >
ThermoRichardsMechanicsProcessData< DisplacementDim, ConstitutiveTraits > & process_data_
std::vector< typename ConstitutiveTraits::StatefulData > current_states_
std::vector< typename ConstitutiveTraits::StatefulDataPrev > prev_states_
std::vector< MaterialStateData< DisplacementDim > > material_states_
NumLib::GenericIntegrationMethod const & integration_method_
MeshLib::Element const & element_
bool const is_axially_symmetric_
ConstitutiveTraits::SolidConstitutiveRelation const & solid_material_
std::vector< typename ConstitutiveTraits::OutputData > output_data_

Member Typedef Documentation

◆ BMatricesType

template<typename ShapeFunctionDisplacement, typename ShapeFunction, int DisplacementDim, typename ConstitutiveTraits>
using ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::BMatricesType

◆ GlobalDimMatrixType

template<typename ShapeFunctionDisplacement, typename ShapeFunction, int DisplacementDim, typename ConstitutiveTraits>
using ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::GlobalDimMatrixType = typename ShapeMatricesType::GlobalDimMatrixType

Definition at line 50 of file ThermoRichardsMechanicsFEM.h.

◆ GlobalDimVectorType

template<typename ShapeFunctionDisplacement, typename ShapeFunction, int DisplacementDim, typename ConstitutiveTraits>
using ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::GlobalDimVectorType = typename ShapeMatricesType::GlobalDimVectorType

Definition at line 51 of file ThermoRichardsMechanicsFEM.h.

◆ Invariants

template<typename ShapeFunctionDisplacement, typename ShapeFunction, int DisplacementDim, typename ConstitutiveTraits>
using ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::Invariants = MathLib::KelvinVector::Invariants<KelvinVectorSize>

Definition at line 63 of file ThermoRichardsMechanicsFEM.h.

◆ IpData

template<typename ShapeFunctionDisplacement, typename ShapeFunction, int DisplacementDim, typename ConstitutiveTraits>
using ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::IpData
Initial value:
ShapeMatricesType, DisplacementDim,
ShapeFunctionDisplacement::NPOINTS>
ShapeMatrixPolicyType< ShapeFunctionDisplacement, DisplacementDim > ShapeMatricesTypeDisplacement
ShapeMatrixPolicyType< ShapeFunction, DisplacementDim > ShapeMatricesType

Definition at line 57 of file ThermoRichardsMechanicsFEM.h.

◆ KelvinVectorType

template<typename ShapeFunctionDisplacement, typename ShapeFunction, int DisplacementDim, typename ConstitutiveTraits>
using ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::KelvinVectorType = typename BMatricesType::KelvinVectorType

Definition at line 55 of file ThermoRichardsMechanicsFEM.h.

◆ ShapeMatricesType

template<typename ShapeFunctionDisplacement, typename ShapeFunction, int DisplacementDim, typename ConstitutiveTraits>
using ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::ShapeMatricesType
Initial value:
EigenFixedShapeMatrixPolicy< ShapeFunction, GlobalDim > ShapeMatrixPolicyType

Definition at line 47 of file ThermoRichardsMechanicsFEM.h.

◆ ShapeMatricesTypeDisplacement

template<typename ShapeFunctionDisplacement, typename ShapeFunction, int DisplacementDim, typename ConstitutiveTraits>
using ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::ShapeMatricesTypeDisplacement

◆ SymmetricTensor

template<typename ShapeFunctionDisplacement, typename ShapeFunction, int DisplacementDim, typename ConstitutiveTraits>
using ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::SymmetricTensor = Eigen::Matrix<double, KelvinVectorSize, 1>

Definition at line 65 of file ThermoRichardsMechanicsFEM.h.

Constructor & Destructor Documentation

◆ ThermoRichardsMechanicsLocalAssembler() [1/3]

template<typename ShapeFunctionDisplacement, typename ShapeFunction, int DisplacementDim, typename ConstitutiveTraits>
ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::ThermoRichardsMechanicsLocalAssembler ( ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits > const & )
delete

◆ ThermoRichardsMechanicsLocalAssembler() [2/3]

template<typename ShapeFunctionDisplacement, typename ShapeFunction, int DisplacementDim, typename ConstitutiveTraits>
ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::ThermoRichardsMechanicsLocalAssembler ( ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits > && )
delete

◆ ThermoRichardsMechanicsLocalAssembler() [3/3]

template<typename ShapeFunctionDisplacement, typename ShapeFunction, int DisplacementDim, typename ConstitutiveTraits>
ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::ThermoRichardsMechanicsLocalAssembler ( MeshLib::Element const & e,
std::size_t const ,
NumLib::GenericIntegrationMethod const & integration_method,
bool const is_axially_symmetric,
ThermoRichardsMechanicsProcessData< DisplacementDim, ConstitutiveTraits > & process_data )

Definition at line 32 of file ThermoRichardsMechanicsFEM-impl.h.

42{
43 unsigned const n_integration_points =
44 this->integration_method_.getNumberOfPoints();
45
47
48 auto const shape_matrices_u =
53
54 auto const shape_matrices =
58
59 for (unsigned ip = 0; ip < n_integration_points; ip++)
60 {
61 auto& ip_data = ip_data_[ip];
62 auto const& sm_u = shape_matrices_u[ip];
63 ip_data_[ip].integration_weight =
64 integration_method.getWeightedPoint(ip).getWeight() *
65 sm_u.integralMeasure * sm_u.detJ;
66
67 ip_data.N_u = sm_u.N;
68 ip_data.dNdx_u = sm_u.dNdx;
69
70 // ip_data.N_p and ip_data.dNdx_p are used for both p and T variables
71 ip_data.N_p = shape_matrices[ip].N;
72 ip_data.dNdx_p = shape_matrices[ip].dNdx;
73 }
74}
ShapeMatrixPolicyType< ShapeFunctionDisplacement, DisplacementDim > ShapeMatricesTypeDisplacement
LocalAssemblerInterface(MeshLib::Element const &e, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, ThermoRichardsMechanicsProcessData< DisplacementDim, ConstitutiveTraits > &process_data)

References ProcessLib::ThermoRichardsMechanics::LocalAssemblerInterface< DisplacementDim, ConstitutiveTraits >::LocalAssemblerInterface(), MathLib::WeightedPoint::getWeight(), NumLib::GenericIntegrationMethod::getWeightedPoint(), NumLib::initShapeMatrices(), ProcessLib::ThermoRichardsMechanics::LocalAssemblerInterface< DisplacementDim, ConstitutiveTraits >::integration_method_, and ip_data_.

Member Function Documentation

◆ addToLocalMatrixData()

template<typename ShapeFunctionDisplacement, typename ShapeFunction, int DisplacementDim, typename ConstitutiveTraits>
void ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::addToLocalMatrixData ( double const dt,
std::vector< double > const & local_x,
std::vector< double > const & local_x_prev,
LocalMatrices const & loc_mat,
std::vector< double > & local_rhs_data,
std::vector< double > & local_Jac_data ) const
private

Definition at line 274 of file ThermoRichardsMechanicsFEM-impl.h.

284{
285 constexpr auto local_matrix_dim =
288
293
294 auto local_rhs =
298
299 local_Jac.noalias() = loc_mat.Jac;
300 local_rhs.noalias() = -loc_mat.res;
301
302 //
303 // -- Jacobian
304 //
305 block_TT(local_Jac).noalias() += loc_mat.M_TT / dt + loc_mat.K_TT;
306 block_Tp(local_Jac).noalias() +=
307 loc_mat.K_Tp + loc_mat.dK_TT_dp + loc_mat.M_Tp / dt;
308
309 block_pT(local_Jac).noalias() += loc_mat.M_pT / dt + loc_mat.K_pT;
310 block_pp(local_Jac).noalias() +=
311 loc_mat.K_pp + loc_mat.storage_p_a_p / dt + loc_mat.storage_p_a_S_Jpp;
312 block_pu(local_Jac).noalias() = loc_mat.M_pu / dt;
313
314 //
315 // -- Residual
316 //
317 auto const [T, p_L, u] = localDOF(local_x);
318 auto const [T_prev, p_L_prev, u_prev] = localDOF(local_x_prev);
319
320 block_T(local_rhs).noalias() -= loc_mat.M_TT * (T - T_prev) / dt +
321 loc_mat.K_TT * T + loc_mat.K_Tp * p_L +
322 loc_mat.M_Tp * (p_L - p_L_prev) / dt;
323 block_p(local_rhs).noalias() -=
324 loc_mat.K_pp * p_L + loc_mat.K_pT * T +
325 (loc_mat.storage_p_a_p + loc_mat.storage_p_a_S) * (p_L - p_L_prev) /
326 dt +
327 loc_mat.M_pu * (u - u_prev) / dt + loc_mat.M_pT * (T - T_prev) / dt;
328}

References ThermoRichardsMechanicsLocalAssembler(), block_p(), block_pp(), block_pT(), block_pu(), block_T(), block_Tp(), block_TT(), MathLib::createZeroedMatrix(), MathLib::createZeroedVector(), displacement_size, localDOF(), pressure_size, and temperature_size.

Referenced by assembleWithJacobian().

◆ assembleWithJacobian()

template<typename ShapeFunctionDisplacement, typename ShapeFunction, int DisplacementDim, typename ConstitutiveTraits>
void ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::assembleWithJacobian ( double const t,
double const dt,
std::vector< double > const & local_x,
std::vector< double > const & local_x_prev,
std::vector< double > & local_rhs_data,
std::vector< double > & local_Jac_data )
overridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 206 of file ThermoRichardsMechanicsFEM-impl.h.

212{
213 auto& medium =
214 *this->process_data_.media_map.getMedium(this->element_.getID());
215
217 loc_mat.setZero();
219 loc_mat_current_ip.setZero(); // only to set the right matrix sizes
220
222
223 for (unsigned ip = 0; ip < this->integration_method_.getNumberOfPoints();
224 ++ip)
225 {
227 std::nullopt, this->element_.getID(),
231 this->element_, ip_data_[ip].N_u))};
232
234 t, dt, x_position, //
237 medium, //
239 this->current_states_[ip], this->prev_states_[ip],
240 this->material_states_[ip], this->output_data_[ip]);
242 }
243
245
248}
void addToLocalMatrixData(double const dt, std::vector< double > const &local_x, std::vector< double > const &local_x_prev, LocalMatrices const &loc_mat, std::vector< double > &local_rhs_data, std::vector< double > &local_Jac_data) const
void assembleWithJacobianSingleIP(double const t, double const dt, ParameterLib::SpatialPosition const &x_position, std::vector< double > const &local_x, std::vector< double > const &local_x_prev, IpData const &ip_data, typename ConstitutiveTraits::ConstitutiveSetting &CS, MaterialPropertyLib::Medium &medium, LocalMatrices &out, typename ConstitutiveTraits::StatefulData &current_state, typename ConstitutiveTraits::StatefulDataPrev const &prev_state, MaterialStateData< DisplacementDim > &mat_state, typename ConstitutiveTraits::OutputData &output_data) const
ThermoRichardsMechanicsProcessData< DisplacementDim, ConstitutiveTraits > & process_data_

References addToLocalMatrixData(), assembleWithJacobianSingleIP(), ProcessLib::ThermoRichardsMechanics::LocalAssemblerInterface< DisplacementDim, ConstitutiveTraits >::current_states_, ProcessLib::ThermoRichardsMechanics::LocalAssemblerInterface< DisplacementDim, ConstitutiveTraits >::element_, ProcessLib::ThermoRichardsMechanics::LocalAssemblerInterface< DisplacementDim, ConstitutiveTraits >::integration_method_, NumLib::interpolateCoordinates(), ip_data_, massLumping(), ProcessLib::ThermoRichardsMechanics::LocalAssemblerInterface< DisplacementDim, ConstitutiveTraits >::material_states_, ProcessLib::ThermoRichardsMechanics::LocalAssemblerInterface< DisplacementDim, ConstitutiveTraits >::output_data_, ProcessLib::ThermoRichardsMechanics::LocalAssemblerInterface< DisplacementDim, ConstitutiveTraits >::prev_states_, ProcessLib::ThermoRichardsMechanics::LocalAssemblerInterface< DisplacementDim, ConstitutiveTraits >::process_data_, and ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::LocalMatrices::setZero().

◆ assembleWithJacobianSingleIP()

template<typename ShapeFunctionDisplacement, typename ShapeFunction, int DisplacementDim, typename ConstitutiveTraits>
void ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::assembleWithJacobianSingleIP ( double const t,
double const dt,
ParameterLib::SpatialPosition const & x_position,
std::vector< double > const & local_x,
std::vector< double > const & local_x_prev,
IpData const & ip_data,
typename ConstitutiveTraits::ConstitutiveSetting & CS,
MaterialPropertyLib::Medium & medium,
LocalMatrices & out,
typename ConstitutiveTraits::StatefulData & current_state,
typename ConstitutiveTraits::StatefulDataPrev const & prev_state,
MaterialStateData< DisplacementDim > & mat_state,
typename ConstitutiveTraits::OutputData & output_data ) const
private

Definition at line 334 of file ThermoRichardsMechanicsFEM-impl.h.

352{
353 auto const& N_u = ip_data.N_u;
354 auto const& dNdx_u = ip_data.dNdx_u;
355
356 // N and dNdx are used for both p and T variables
357 auto const& N = ip_data.N_p;
358 auto const& dNdx = ip_data.dNdx_p;
359
360 auto const B =
364 dNdx_u, N_u, (*x_position.getCoordinates())[0],
366
368
369 auto const [T, p_L, u] = localDOF(local_x);
370 auto const [T_prev, p_L_prev, u_prev] = localDOF(local_x_prev);
371
372 {
374 this->process_data_, this->solid_material_);
376
377 double const T_ip = N * T;
378 double const T_prev_ip = N * T_prev;
380
381 double const p_cap_ip = -N * p_L;
382 double const p_cap_prev_ip = -N * p_L_prev;
384
386
387 CS.eval(models, t, dt, x_position, //
388 medium, //
389 {T_ip, T_prev_ip, grad_T_ip}, //
392 CD);
393 }
394
396 NodalMatrix const NTN = N.transpose() * N;
397 NodalMatrix const dNTdN = dNdx.transpose() * dNdx;
398
399 // TODO is identity2.transpose() * B something like divergence?
405 BTI2N = B.transpose() * identity2 * N;
406
407 /*
408 * Conventions:
409 *
410 * * use positive signs exclusively, any negative signs must be included in
411 * the coefficients coming from the constitutive setting
412 * * the used coefficients from the constitutive setting are named after the
413 * terms they appear in, e.g. K_TT_X_dNTdN means it is a scalar (X) that
414 * will be multiplied by dNdx.transpose() * dNdx. Placefolders for the
415 * coefficients are:
416 * * X -> scalar
417 * * V -> vector
418 * * K -> Kelvin vector
419 * * the Laplace terms have a special name, e.g., K_TT_Laplace
420 * * there shall be only one contribution to each of the LocalMatrices,
421 * assigned with = assignment; this point might be relaxed in the future
422 * * this method will overwrite the data in the passed LocalMatrices& out
423 * argument, not add to it
424 */
425
426 // residual, order T, p, u
427 block_p(out.res).noalias() =
429 block_u(out.res).noalias() =
430 B.transpose() *
433 .sigma_total -
434 static_cast<int>(this->process_data_.apply_body_force_for_deformation) *
435 N_u_op(N_u).transpose() *
437
438 // Storage matrices
439 out.storage_p_a_p.noalias() =
441 out.storage_p_a_S.noalias() =
442 std::get<TRMStorageData>(CD).storage_p_a_S_X_NTN * NTN;
443 out.storage_p_a_S_Jpp.noalias() =
444 std::get<TRMStorageData>(CD).storage_p_a_S_Jpp_X_NTN * NTN;
445
446 // M matrices, order T, p, u
447 out.M_TT.noalias() =
449 out.M_Tp.noalias() =
451
452 out.M_pT.noalias() =
454 out.M_pu.noalias() =
456
457 // K matrices, order T, p, u
458 out.K_TT.noalias() =
459 dNdx.transpose() *
461 .K_TT_Laplace *
462 dNdx +
463 N.transpose() *
464 (std::get<EqTData<DisplacementDim>>(CD).K_TT_NT_V_dN.transpose() *
465 dNdx) +
467 dNTdN;
468
469 out.dK_TT_dp.noalias() =
470 N.transpose() *
472 .K_Tp_NT_V_dN.transpose() *
473 dNdx) +
475 NTN;
476 out.K_Tp.noalias() =
477 dNdx.transpose() *
479 dNdx +
481 dNTdN;
482
483 out.K_pp.noalias() =
485 dNdx +
487 dNTdN;
488 out.K_pT.noalias() =
489 dNdx.transpose() *
491
492 // direct Jacobian contributions, order T, p, u
493 block_pT(out.Jac).noalias() =
495 dNTdN;
496 block_pp(out.Jac).noalias() =
497 std::get<TRMStorageData>(CD).J_pp_X_NTN * NTN +
499 BTI2N.transpose() * (u - u_prev) / dt *
500 N // TODO something with volumetric strain rate?
501 + dNdx.transpose() *
503
504 block_uT(out.Jac).noalias() =
505 B.transpose() *
507 N;
508 block_up(out.Jac).noalias() =
509 B.transpose() *
511 .J_up_BT_K_N *
512 N +
513 N_u_op(N_u).transpose() *
515 block_uu(out.Jac).noalias() =
516 B.transpose() *
519 B;
520
521 out *= ip_data.integration_weight;
522}
MatrixType< _kelvin_vector_size, _number_of_dof > BMatrixType
constexpr int kelvin_vector_dimensions(int const displacement_dim)
Kelvin vector dimensions for given displacement dimension.
MatrixType< ShapeFunction::NPOINTS, ShapeFunction::NPOINTS > NodalMatrixType

References ThermoRichardsMechanicsLocalAssembler(), block_p(), block_pp(), block_pT(), block_u(), block_up(), block_uT(), block_uu(), ProcessLib::LinearBMatrix::computeBMatrix(), displacement_size, ProcessLib::Graph::get(), ParameterLib::SpatialPosition::getCoordinates(), ProcessLib::ThermoRichardsMechanics::LocalAssemblerInterface< DisplacementDim, ConstitutiveTraits >::is_axially_symmetric_, ProcessLib::ThermoRichardsMechanics::TRMVaporDiffusionData< DisplacementDim >::K_pp_X_dNTdN, ProcessLib::ThermoRichardsMechanics::TRMVaporDiffusionData< DisplacementDim >::K_Tp_X_dNTdN, MathLib::KelvinVector::kelvin_vector_dimensions(), localDOF(), N_u_op, pressure_size, ProcessLib::ThermoRichardsMechanics::LocalAssemblerInterface< DisplacementDim, ConstitutiveTraits >::process_data_, and ProcessLib::ThermoRichardsMechanics::LocalAssemblerInterface< DisplacementDim, ConstitutiveTraits >::solid_material_.

Referenced by assembleWithJacobian().

◆ block_p()

template<typename ShapeFunctionDisplacement, typename ShapeFunction, int DisplacementDim, typename ConstitutiveTraits>
auto ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::block_p ( auto & vec)
inlinestaticprivate

◆ block_pp()

template<typename ShapeFunctionDisplacement, typename ShapeFunction, int DisplacementDim, typename ConstitutiveTraits>
auto ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::block_pp ( auto & mat)
inlinestaticprivate

Definition at line 384 of file ThermoRichardsMechanicsFEM.h.

385 {
388 }

References pressure_index.

Referenced by addToLocalMatrixData(), and assembleWithJacobianSingleIP().

◆ block_pT()

template<typename ShapeFunctionDisplacement, typename ShapeFunction, int DisplacementDim, typename ConstitutiveTraits>
auto ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::block_pT ( auto & mat)
inlinestaticprivate

◆ block_pu()

template<typename ShapeFunctionDisplacement, typename ShapeFunction, int DisplacementDim, typename ConstitutiveTraits>
auto ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::block_pu ( auto & mat)
inlinestaticprivate

◆ block_T()

template<typename ShapeFunctionDisplacement, typename ShapeFunction, int DisplacementDim, typename ConstitutiveTraits>
auto ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::block_T ( auto & vec)
inlinestaticprivate

Definition at line 413 of file ThermoRichardsMechanicsFEM.h.

414 {
416 }

References temperature_index.

Referenced by addToLocalMatrixData(), and computeSecondaryVariableConcrete().

◆ block_Tp()

template<typename ShapeFunctionDisplacement, typename ShapeFunction, int DisplacementDim, typename ConstitutiveTraits>
auto ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::block_Tp ( auto & mat)
inlinestaticprivate

Definition at line 394 of file ThermoRichardsMechanicsFEM.h.

References pressure_index, and temperature_index.

Referenced by addToLocalMatrixData().

◆ block_TT()

template<typename ShapeFunctionDisplacement, typename ShapeFunction, int DisplacementDim, typename ConstitutiveTraits>
auto ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::block_TT ( auto & mat)
inlinestaticprivate

Definition at line 399 of file ThermoRichardsMechanicsFEM.h.

References temperature_index.

Referenced by addToLocalMatrixData().

◆ block_u()

template<typename ShapeFunctionDisplacement, typename ShapeFunction, int DisplacementDim, typename ConstitutiveTraits>
auto ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::block_u ( auto & vec)
inlinestaticprivate

◆ block_up()

template<typename ShapeFunctionDisplacement, typename ShapeFunction, int DisplacementDim, typename ConstitutiveTraits>
auto ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::block_up ( auto & mat)
inlinestaticprivate

◆ block_uT()

template<typename ShapeFunctionDisplacement, typename ShapeFunction, int DisplacementDim, typename ConstitutiveTraits>
auto ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::block_uT ( auto & mat)
inlinestaticprivate

◆ block_uu()

template<typename ShapeFunctionDisplacement, typename ShapeFunction, int DisplacementDim, typename ConstitutiveTraits>
auto ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::block_uu ( auto & mat)
inlinestaticprivate

◆ computeSecondaryVariableConcrete()

template<typename ShapeFunctionDisplacement, typename ShapeFunction, int DisplacementDim, typename ConstitutiveTraits>
void ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::computeSecondaryVariableConcrete ( double const t,
double const dt,
Eigen::VectorXd const & local_x,
Eigen::VectorXd const & local_x_prev )
overridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 528 of file ThermoRichardsMechanicsFEM-impl.h.

532{
533 auto const T = block_T(local_x);
534 auto const p_L = block_p(local_x);
535 auto const u = block_u(local_x);
536
537 auto const T_prev = block_T(local_x_prev);
538 auto const p_L_prev = block_p(local_x_prev);
539
540 auto const e_id = this->element_.getID();
541 auto const& process_data = this->process_data_;
542 auto& medium = *process_data.media_map.getMedium(e_id);
543
544 unsigned const n_integration_points =
545 this->integration_method_.getNumberOfPoints();
546
548
553
554 for (unsigned ip = 0; ip < n_integration_points; ip++)
555 {
556 auto& current_state = this->current_states_[ip];
557 auto& output_data = this->output_data_[ip];
558
559 auto const& ip_data = ip_data_[ip];
560
561 // N is used for both p and T variables
562 auto const& N = ip_data.N_p;
563 auto const& N_u = ip_data.N_u;
564 auto const& dNdx_u = ip_data.dNdx_u;
565 auto const& dNdx = ip_data.dNdx_p;
566
568 std::nullopt, this->element_.getID(),
572 this->element_, N_u))};
573
574 auto const x_coord =
575 x_position.getCoordinates().value()[0]; // r for axisymmetry
576 auto const B =
581
582 double const T_ip = N * T;
583 double const T_prev_ip = N * T_prev;
585
586 double const p_cap_ip = -N * p_L;
587 double const p_cap_prev_ip = -N * p_L_prev;
589
591
593 t, dt, x_position, //
594 medium, //
595 {T_ip, T_prev_ip, grad_T_ip}, //
599 CD);
600 }
601
605 *process_data.pressure_interpolated);
609 *process_data.temperature_interpolated);
610}

References block_p(), block_T(), block_u(), ProcessLib::LinearBMatrix::computeBMatrix(), ProcessLib::ThermoRichardsMechanics::LocalAssemblerInterface< DisplacementDim, ConstitutiveTraits >::current_states_, ProcessLib::ThermoRichardsMechanics::LocalAssemblerInterface< DisplacementDim, ConstitutiveTraits >::element_, ParameterLib::SpatialPosition::getCoordinates(), ProcessLib::ThermoRichardsMechanics::LocalAssemblerInterface< DisplacementDim, ConstitutiveTraits >::integration_method_, NumLib::interpolateCoordinates(), NumLib::interpolateToHigherOrderNodes(), ip_data_, ProcessLib::ThermoRichardsMechanics::LocalAssemblerInterface< DisplacementDim, ConstitutiveTraits >::is_axially_symmetric_, ProcessLib::ThermoRichardsMechanics::LocalAssemblerInterface< DisplacementDim, ConstitutiveTraits >::material_states_, ProcessLib::ThermoRichardsMechanics::LocalAssemblerInterface< DisplacementDim, ConstitutiveTraits >::output_data_, ProcessLib::ThermoRichardsMechanics::LocalAssemblerInterface< DisplacementDim, ConstitutiveTraits >::prev_states_, ProcessLib::ThermoRichardsMechanics::LocalAssemblerInterface< DisplacementDim, ConstitutiveTraits >::process_data_, and ProcessLib::ThermoRichardsMechanics::LocalAssemblerInterface< DisplacementDim, ConstitutiveTraits >::solid_material_.

◆ convertInitialStressType()

template<typename ShapeFunctionDisplacement, typename ShapeFunction, int DisplacementDim, typename ConstitutiveTraits>
void ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::convertInitialStressType ( unsigned const ip,
double const t,
ParameterLib::SpatialPosition const x_position,
MaterialPropertyLib::Medium const & medium,
MPL::VariableArray const & variables,
double const p_at_ip )
private

This function resets the initial stress type according to the input initial stress type, either total or effective. If subtype = StressSaturation_StrainPressureTemperature is bing used in the process setting, the initial effective stress is converted to total stress. Otherwise, the initial total stress is converted to effective stress.

Definition at line 164 of file ThermoRichardsMechanicsFEM-impl.h.

171{
177 this->process_data_.initial_stress.type ==
179 {
180 return;
181 }
182
184 this->process_data_.initial_stress.type == InitialStress::Type::Total)
185 {
186 return;
187 }
188
189 double const alpha_b =
191 .template value<double>(variables, x_position, t, 0.0 /*dt*/);
192
193 double const bishop =
195 .template value<double>(variables, x_position, t, 0.0 /*dt*/);
196
198 this->current_states_[ip], this->prev_states_[ip],
200}
static Eigen::Matrix< double, KelvinVectorSize, 1 > const identity2
Kelvin mapping of 2nd order identity tensor.

References MaterialPropertyLib::biot_coefficient, MaterialPropertyLib::bishops_effective_stress, ProcessLib::ThermoRichardsMechanics::LocalAssemblerInterface< DisplacementDim, ConstitutiveTraits >::current_states_, ProcessLib::InitialStress::Effective, MathLib::KelvinVector::Invariants< KelvinVectorSize >::identity2, ProcessLib::ThermoRichardsMechanics::LocalAssemblerInterface< DisplacementDim, ConstitutiveTraits >::prev_states_, ProcessLib::ThermoRichardsMechanics::LocalAssemblerInterface< DisplacementDim, ConstitutiveTraits >::process_data_, MaterialPropertyLib::Medium::property(), and ProcessLib::InitialStress::Total.

Referenced by setInitialConditionsConcrete().

◆ getNumberOfVectorElementsForDeformation()

template<typename ShapeFunctionDisplacement, typename ShapeFunction, int DisplacementDim, typename ConstitutiveTraits>
int ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::getNumberOfVectorElementsForDeformation ( ) const
inlineoverridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 244 of file ThermoRichardsMechanicsFEM.h.

245 {
246 return displacement_size;
247 }

References displacement_size.

◆ getShapeMatrix()

template<typename ShapeFunctionDisplacement, typename ShapeFunction, int DisplacementDim, typename ConstitutiveTraits>
Eigen::Map< const Eigen::RowVectorXd > ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::getShapeMatrix ( const unsigned integration_point) const
inlineoverridevirtual

Provides the shape matrix at the given integration point.

Implements NumLib::ExtrapolatableElement.

Definition at line 345 of file ThermoRichardsMechanicsFEM.h.

347 {
348 auto const& N_u = ip_data_[integration_point].N_u;
349
350 // assumes N is stored contiguously in memory
351 return Eigen::Map<const Eigen::RowVectorXd>(N_u.data(), N_u.size());
352 }

References ip_data_.

◆ initializeConcrete()

template<typename ShapeFunctionDisplacement, typename ShapeFunction, int DisplacementDim, typename ConstitutiveTraits>
void ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::initializeConcrete ( )
inlineoverridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 272 of file ThermoRichardsMechanicsFEM.h.

273 {
274 unsigned const n_integration_points =
275 this->integration_method_.getNumberOfPoints();
277 auto const& medium =
278 *this->process_data_.media_map.getMedium(this->element_.getID());
279
280 for (unsigned ip = 0; ip < n_integration_points; ip++)
281 {
283 std::nullopt, this->element_.getID(),
287 this->element_, this->ip_data_[ip].N_u))};
288 auto& current_state = this->current_states_[ip];
289
290 // Set initial stress from parameter.
291 if (this->process_data_.initial_stress.value)
292 {
297 (*this->process_data_.initial_stress.value)(
299 }
300
301 if (*this->process_data_.initialize_porosity_from_medium_property)
302 {
303 // Initial porosity. Could be read from integration point data
304 // or mesh.
306 medium.property(MPL::porosity)
309
311 {
316 }
317 else
318 {
321 }
322 }
323
324 double const t = 0; // TODO (naumov) pass t from top
325 this->solid_material_.initializeInternalStateVariables(
326 t, x_position,
328 }
329
330 for (unsigned ip = 0; ip < n_integration_points; ip++)
331 {
332 this->material_states_[ip].pushBackState();
333 }
334
335 for (unsigned ip = 0; ip < n_integration_points; ip++)
336 {
337 this->prev_states_[ip] = this->current_states_[ip];
338 }
339 }

References ProcessLib::ThermoRichardsMechanics::LocalAssemblerInterface< DisplacementDim, ConstitutiveTraits >::current_states_, ProcessLib::ThermoRichardsMechanics::LocalAssemblerInterface< DisplacementDim, ConstitutiveTraits >::element_, ProcessLib::ThermoRichardsMechanics::LocalAssemblerInterface< DisplacementDim, ConstitutiveTraits >::integration_method_, NumLib::interpolateCoordinates(), ip_data_, ProcessLib::ThermoRichardsMechanics::LocalAssemblerInterface< DisplacementDim, ConstitutiveTraits >::material_states_, MaterialPropertyLib::porosity, ProcessLib::ThermoRichardsMechanics::LocalAssemblerInterface< DisplacementDim, ConstitutiveTraits >::prev_states_, ProcessLib::ThermoRichardsMechanics::LocalAssemblerInterface< DisplacementDim, ConstitutiveTraits >::process_data_, ProcessLib::ThermoRichardsMechanics::LocalAssemblerInterface< DisplacementDim, ConstitutiveTraits >::solid_material_, MathLib::KelvinVector::symmetricTensorToKelvinVector(), and MaterialPropertyLib::transport_porosity.

◆ localDOF()

template<typename ShapeFunctionDisplacement, typename ShapeFunction, int DisplacementDim, typename ConstitutiveTraits>
constexpr auto ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::localDOF ( std::vector< double > const & x)
inlinestaticconstexprprivate

◆ massLumping()

template<typename ShapeFunctionDisplacement, typename ShapeFunction, int DisplacementDim, typename ConstitutiveTraits>
void ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::massLumping ( LocalMatrices & loc_mat) const
private

Definition at line 254 of file ThermoRichardsMechanicsFEM-impl.h.

258{
259 if (this->process_data_.apply_mass_lumping)
260 {
261 loc_mat.storage_p_a_p =
262 loc_mat.storage_p_a_p.colwise().sum().eval().asDiagonal();
263 loc_mat.storage_p_a_S =
264 loc_mat.storage_p_a_S.colwise().sum().eval().asDiagonal();
265 loc_mat.storage_p_a_S_Jpp =
266 loc_mat.storage_p_a_S_Jpp.colwise().sum().eval().asDiagonal();
267 }
268}

References ThermoRichardsMechanicsLocalAssembler(), and ProcessLib::ThermoRichardsMechanics::LocalAssemblerInterface< DisplacementDim, ConstitutiveTraits >::process_data_.

Referenced by assembleWithJacobian().

◆ setInitialConditionsConcrete()

template<typename ShapeFunctionDisplacement, typename ShapeFunction, int DisplacementDim, typename ConstitutiveTraits>
void ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::setInitialConditionsConcrete ( Eigen::VectorXd const local_x,
double const t,
int const process_id )
overridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 80 of file ThermoRichardsMechanicsFEM-impl.h.

84{
85 assert(local_x.size() ==
87
88 auto const p_L = local_x.template segment<pressure_size>(pressure_index);
89 auto const T =
91
92 constexpr double dt = std::numeric_limits<double>::quiet_NaN();
93 auto const& medium =
94 *this->process_data_.media_map.getMedium(this->element_.getID());
95
97
100 this->process_data_, this->solid_material_);
101
102 unsigned const n_integration_points =
103 this->integration_method_.getNumberOfPoints();
104 for (unsigned ip = 0; ip < n_integration_points; ip++)
105 {
106 // N is used for both T and p variables.
107 auto const& N = ip_data_[ip].N_p;
108
110 std::nullopt, this->element_.getID(),
114 this->element_, ip_data_[ip].N_u))};
115
116 double T_ip;
118
119 double p_cap_ip;
121
123
124 if (medium.hasProperty(MPL::PropertyType::saturation))
125 {
126 variables.capillary_pressure = p_cap_ip;
127 variables.liquid_phase_pressure = -p_cap_ip;
128 // setting pG to 1 atm
129 // TODO : rewrite equations s.t. p_L = pG-p_cap
130 variables.gas_phase_pressure = 1.0e5;
131 variables.temperature = T_ip;
132
133 double const S_L =
135 .template value<double>(variables, x_position, t, dt);
137 S_L;
138
139 variables.liquid_saturation = S_L;
140 }
141
143 {T_ip, 0, {}}, this->current_states_[ip],
144 this->prev_states_[ip]);
145
146 if (this->process_data_.initial_stress.value)
147 {
148 if (!medium.hasProperty(MPL::PropertyType::saturation))
149 {
150 OGS_FATAL(
151 "The medium has no saturation property required to compute "
152 "initial stress.");
153 }
155 -p_cap_ip);
156 }
157 }
158}
#define OGS_FATAL(...)
Definition Error.h:19
void convertInitialStressType(unsigned const ip, double const t, ParameterLib::SpatialPosition const x_position, MaterialPropertyLib::Medium const &medium, MPL::VariableArray const &variables, double const p_at_ip)
void shapeFunctionInterpolate(const NodalValues &, const ShapeMatrix &)

References MaterialPropertyLib::VariableArray::capillary_pressure, convertInitialStressType(), ProcessLib::ThermoRichardsMechanics::LocalAssemblerInterface< DisplacementDim, ConstitutiveTraits >::current_states_, displacement_size, ProcessLib::ThermoRichardsMechanics::LocalAssemblerInterface< DisplacementDim, ConstitutiveTraits >::element_, MaterialPropertyLib::VariableArray::gas_phase_pressure, ProcessLib::ThermoRichardsMechanics::LocalAssemblerInterface< DisplacementDim, ConstitutiveTraits >::integration_method_, NumLib::interpolateCoordinates(), ip_data_, MaterialPropertyLib::VariableArray::liquid_phase_pressure, MaterialPropertyLib::VariableArray::liquid_saturation, OGS_FATAL, pressure_index, pressure_size, ProcessLib::ThermoRichardsMechanics::LocalAssemblerInterface< DisplacementDim, ConstitutiveTraits >::prev_states_, ProcessLib::ThermoRichardsMechanics::LocalAssemblerInterface< DisplacementDim, ConstitutiveTraits >::process_data_, MaterialPropertyLib::saturation, setInitialConditionsConcrete(), NumLib::detail::shapeFunctionInterpolate(), ProcessLib::ThermoRichardsMechanics::LocalAssemblerInterface< DisplacementDim, ConstitutiveTraits >::solid_material_, MaterialPropertyLib::VariableArray::temperature, temperature_index, and temperature_size.

Referenced by setInitialConditionsConcrete().

Member Data Documentation

◆ displacement_index

template<typename ShapeFunctionDisplacement, typename ShapeFunction, int DisplacementDim, typename ConstitutiveTraits>
int ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::displacement_index = 2 * ShapeFunction::NPOINTS
staticconstexprprivate

Definition at line 38 of file ThermoRichardsMechanicsFEM.h.

Referenced by block_pu(), block_u(), block_up(), block_uT(), and block_uu().

◆ displacement_size

template<typename ShapeFunctionDisplacement, typename ShapeFunction, int DisplacementDim, typename ConstitutiveTraits>
int ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::displacement_size
staticconstexprprivate

◆ ip_data_

template<typename ShapeFunctionDisplacement, typename ShapeFunction, int DisplacementDim, typename ConstitutiveTraits>
std::vector<IpData> ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::ip_data_
private

◆ KelvinVectorSize

template<typename ShapeFunctionDisplacement, typename ShapeFunction, int DisplacementDim, typename ConstitutiveTraits>
int const ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::KelvinVectorSize
static
Initial value:

Definition at line 61 of file ThermoRichardsMechanicsFEM.h.

◆ N_u_op

template<typename ShapeFunctionDisplacement, typename ShapeFunction, int DisplacementDim, typename ConstitutiveTraits>
auto& ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::N_u_op
staticconstexpr
Initial value:
DisplacementDim,
constexpr Eigen::CwiseNullaryOp< EigenBlockMatrixViewFunctor< D, M >, typename EigenBlockMatrixViewFunctor< D, M >::Matrix > eigenBlockMatrixView(const Eigen::MatrixBase< M > &matrix)
RowVectorType< ShapeFunction::NPOINTS > NodalRowVectorType

Definition at line 67 of file ThermoRichardsMechanicsFEM.h.

Referenced by assembleWithJacobianSingleIP().

◆ pressure_index

template<typename ShapeFunctionDisplacement, typename ShapeFunction, int DisplacementDim, typename ConstitutiveTraits>
int ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::pressure_index = temperature_size
staticconstexprprivate

◆ pressure_size

template<typename ShapeFunctionDisplacement, typename ShapeFunction, int DisplacementDim, typename ConstitutiveTraits>
int ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::pressure_size = ShapeFunction::NPOINTS
staticconstexprprivate

◆ temperature_index

template<typename ShapeFunctionDisplacement, typename ShapeFunction, int DisplacementDim, typename ConstitutiveTraits>
int ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::temperature_index = 0
staticconstexprprivate

◆ temperature_size

template<typename ShapeFunctionDisplacement, typename ShapeFunction, int DisplacementDim, typename ConstitutiveTraits>
int ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::temperature_size = ShapeFunction::NPOINTS
staticconstexprprivate

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