Loading [MathJax]/extensions/MathMenu.js
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 37 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
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.
virtual std::optional< VectorSegmentgetVectorDeformationSegment () const
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 56 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 57 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 69 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 63 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 61 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 53 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 71 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 39 of file ThermoRichardsMechanicsFEM-impl.h.

49{
50 unsigned const n_integration_points =
51 this->integration_method_.getNumberOfPoints();
52
54
55 auto const shape_matrices_u =
60
61 auto const shape_matrices =
65
66 for (unsigned ip = 0; ip < n_integration_points; ip++)
67 {
68 auto& ip_data = ip_data_[ip];
69 auto const& sm_u = shape_matrices_u[ip];
70 ip_data_[ip].integration_weight =
71 integration_method.getWeightedPoint(ip).getWeight() *
72 sm_u.integralMeasure * sm_u.detJ;
73
74 ip_data.N_u = sm_u.N;
75 ip_data.dNdx_u = sm_u.dNdx;
76
77 // ip_data.N_p and ip_data.dNdx_p are used for both p and T variables
78 ip_data.N_p = shape_matrices[ip].N;
79 ip_data.dNdx_p = shape_matrices[ip].dNdx;
80 }
81}
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 281 of file ThermoRichardsMechanicsFEM-impl.h.

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

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 213 of file ThermoRichardsMechanicsFEM-impl.h.

219{
220 auto& medium =
221 *this->process_data_.media_map.getMedium(this->element_.getID());
222
224 loc_mat.setZero();
226 loc_mat_current_ip.setZero(); // only to set the right matrix sizes
227
229
230 for (unsigned ip = 0; ip < this->integration_method_.getNumberOfPoints();
231 ++ip)
232 {
234 std::nullopt, this->element_.getID(),
238 this->element_, ip_data_[ip].N_u))};
239
241 t, dt, x_position, //
244 medium, //
246 this->current_states_[ip], this->prev_states_[ip],
247 this->material_states_[ip], this->output_data_[ip]);
249 }
250
252
255}
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
std::vector< typename ConstitutiveTraits::OutputData > output_data_
ThermoRichardsMechanicsProcessData< DisplacementDim, ConstitutiveTraits > & process_data_
std::vector< typename ConstitutiveTraits::StatefulDataPrev > prev_states_
std::vector< MaterialStateData< DisplacementDim > > material_states_
std::vector< typename ConstitutiveTraits::StatefulData > current_states_

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 341 of file ThermoRichardsMechanicsFEM-impl.h.

359{
360 auto const& N_u = ip_data.N_u;
361 auto const& dNdx_u = ip_data.dNdx_u;
362
363 // N and dNdx are used for both p and T variables
364 auto const& N = ip_data.N_p;
365 auto const& dNdx = ip_data.dNdx_p;
366
367 auto const B =
371 dNdx_u, N_u, (*x_position.getCoordinates())[0],
373
375
376 auto const [T, p_L, u] = localDOF(local_x);
377 auto const [T_prev, p_L_prev, u_prev] = localDOF(local_x_prev);
378
379 {
381 this->process_data_, this->solid_material_);
383
384 double const T_ip = N * T;
385 double const T_prev_ip = N * T_prev;
387
388 double const p_cap_ip = -N * p_L;
389 double const p_cap_prev_ip = -N * p_L_prev;
391
393
394 CS.eval(models, t, dt, x_position, //
395 medium, //
396 {T_ip, T_prev_ip, grad_T_ip}, //
399 CD);
400 }
401
403 NodalMatrix const NTN = N.transpose() * N;
404 NodalMatrix const dNTdN = dNdx.transpose() * dNdx;
405
406 // TODO is identity2.transpose() * B something like divergence?
412 BTI2N = B.transpose() * identity2 * N;
413
414 /*
415 * Conventions:
416 *
417 * * use positive signs exclusively, any negative signs must be included in
418 * the coefficients coming from the constitutive setting
419 * * the used coefficients from the constitutive setting are named after the
420 * terms they appear in, e.g. K_TT_X_dNTdN means it is a scalar (X) that
421 * will be multiplied by dNdx.transpose() * dNdx. Placefolders for the
422 * coefficients are:
423 * * X -> scalar
424 * * V -> vector
425 * * K -> Kelvin vector
426 * * the Laplace terms have a special name, e.g., K_TT_Laplace
427 * * there shall be only one contribution to each of the LocalMatrices,
428 * assigned with = assignment; this point might be relaxed in the future
429 * * this method will overwrite the data in the passed LocalMatrices& out
430 * argument, not add to it
431 */
432
433 // residual, order T, p, u
434 block_p(out.res).noalias() =
436 block_u(out.res).noalias() =
437 B.transpose() *
440 .sigma_total -
441 static_cast<int>(this->process_data_.apply_body_force_for_deformation) *
442 N_u_op(N_u).transpose() *
444
445 // Storage matrices
446 out.storage_p_a_p.noalias() =
448 out.storage_p_a_S.noalias() =
449 std::get<TRMStorageData>(CD).storage_p_a_S_X_NTN * NTN;
450 out.storage_p_a_S_Jpp.noalias() =
451 std::get<TRMStorageData>(CD).storage_p_a_S_Jpp_X_NTN * NTN;
452
453 // M matrices, order T, p, u
454 out.M_TT.noalias() =
456 out.M_Tp.noalias() =
458
459 out.M_pT.noalias() =
461 out.M_pu.noalias() =
463
464 // K matrices, order T, p, u
465 out.K_TT.noalias() =
466 dNdx.transpose() *
468 .K_TT_Laplace *
469 dNdx +
470 N.transpose() *
471 (std::get<EqTData<DisplacementDim>>(CD).K_TT_NT_V_dN.transpose() *
472 dNdx) +
474 dNTdN;
475
476 out.dK_TT_dp.noalias() =
477 N.transpose() *
479 .K_Tp_NT_V_dN.transpose() *
480 dNdx) +
482 NTN;
483 out.K_Tp.noalias() =
484 dNdx.transpose() *
486 dNdx +
488 dNTdN;
489
490 out.K_pp.noalias() =
492 dNdx +
494 dNTdN;
495 out.K_pT.noalias() =
496 dNdx.transpose() *
498
499 // direct Jacobian contributions, order T, p, u
500 block_pT(out.Jac).noalias() =
502 dNTdN;
503 block_pp(out.Jac).noalias() =
504 std::get<TRMStorageData>(CD).J_pp_X_NTN * NTN +
506 BTI2N.transpose() * (u - u_prev) / dt *
507 N // TODO something with volumetric strain rate?
508 + dNdx.transpose() *
510
511 block_uT(out.Jac).noalias() =
512 B.transpose() *
514 N;
515 block_up(out.Jac).noalias() =
516 B.transpose() *
518 .J_up_BT_K_N *
519 N +
520 N_u_op(N_u).transpose() *
522 block_uu(out.Jac).noalias() =
523 B.transpose() *
526 B;
527
528 out *= ip_data.integration_weight;
529}
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
ConstitutiveTraits::SolidConstitutiveRelation const & solid_material_

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 385 of file ThermoRichardsMechanicsFEM.h.

386 {
389 }

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 414 of file ThermoRichardsMechanicsFEM.h.

415 {
417 }

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 395 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 400 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 535 of file ThermoRichardsMechanicsFEM-impl.h.

539{
540 auto const T = block_T(local_x);
541 auto const p_L = block_p(local_x);
542 auto const u = block_u(local_x);
543
544 auto const T_prev = block_T(local_x_prev);
545 auto const p_L_prev = block_p(local_x_prev);
546
547 auto const e_id = this->element_.getID();
548 auto const& process_data = this->process_data_;
549 auto& medium = *process_data.media_map.getMedium(e_id);
550
551 unsigned const n_integration_points =
552 this->integration_method_.getNumberOfPoints();
553
555
560
561 for (unsigned ip = 0; ip < n_integration_points; ip++)
562 {
563 auto& current_state = this->current_states_[ip];
564 auto& output_data = this->output_data_[ip];
565
566 auto const& ip_data = ip_data_[ip];
567
568 // N is used for both p and T variables
569 auto const& N = ip_data.N_p;
570 auto const& N_u = ip_data.N_u;
571 auto const& dNdx_u = ip_data.dNdx_u;
572 auto const& dNdx = ip_data.dNdx_p;
573
575 std::nullopt, this->element_.getID(),
579 this->element_, N_u))};
580 auto const x_coord =
583 this->element_, N_u);
584 auto const B =
589
590 double const T_ip = N * T;
591 double const T_prev_ip = N * T_prev;
593
594 double const p_cap_ip = -N * p_L;
595 double const p_cap_prev_ip = -N * p_L_prev;
597
599
601 t, dt, x_position, //
602 medium, //
603 {T_ip, T_prev_ip, grad_T_ip}, //
607 CD);
608 }
609
613 *process_data.pressure_interpolated);
617 *process_data.temperature_interpolated);
618}

References block_p(), block_T(), block_u(), ProcessLib::LinearBMatrix::computeBMatrix(), ProcessLib::ThermoRichardsMechanics::LocalAssemblerInterface< DisplacementDim, ConstitutiveTraits >::current_states_, ProcessLib::ThermoRichardsMechanics::LocalAssemblerInterface< DisplacementDim, ConstitutiveTraits >::element_, ProcessLib::ThermoRichardsMechanics::LocalAssemblerInterface< DisplacementDim, ConstitutiveTraits >::integration_method_, NumLib::interpolateCoordinates(), NumLib::interpolateToHigherOrderNodes(), NumLib::interpolateXCoordinate(), 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 171 of file ThermoRichardsMechanicsFEM-impl.h.

178{
184 this->process_data_.initial_stress.type ==
186 {
187 return;
188 }
189
191 this->process_data_.initial_stress.type == InitialStress::Type::Total)
192 {
193 return;
194 }
195
196 double const alpha_b =
198 .template value<double>(variables, x_position, t, 0.0 /*dt*/);
199
200 double const bishop =
202 .template value<double>(variables, x_position, t, 0.0 /*dt*/);
203
205 this->current_states_[ip], this->prev_states_[ip],
207}
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().

◆ 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 346 of file ThermoRichardsMechanicsFEM.h.

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

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 273 of file ThermoRichardsMechanicsFEM.h.

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

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 261 of file ThermoRichardsMechanicsFEM-impl.h.

265{
266 if (this->process_data_.apply_mass_lumping)
267 {
268 loc_mat.storage_p_a_p =
269 loc_mat.storage_p_a_p.colwise().sum().eval().asDiagonal();
270 loc_mat.storage_p_a_S =
271 loc_mat.storage_p_a_S.colwise().sum().eval().asDiagonal();
272 loc_mat.storage_p_a_S_Jpp =
273 loc_mat.storage_p_a_S_Jpp.colwise().sum().eval().asDiagonal();
274 }
275}

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 87 of file ThermoRichardsMechanicsFEM-impl.h.

91{
92 assert(local_x.size() ==
94
95 auto const p_L = local_x.template segment<pressure_size>(pressure_index);
96 auto const T =
98
99 constexpr double dt = std::numeric_limits<double>::quiet_NaN();
100 auto const& medium =
101 *this->process_data_.media_map.getMedium(this->element_.getID());
102
104
107 this->process_data_, this->solid_material_);
108
109 unsigned const n_integration_points =
110 this->integration_method_.getNumberOfPoints();
111 for (unsigned ip = 0; ip < n_integration_points; ip++)
112 {
113 // N is used for both T and p variables.
114 auto const& N = ip_data_[ip].N_p;
115
117 std::nullopt, this->element_.getID(),
121 this->element_, ip_data_[ip].N_u))};
122
123 double T_ip;
125
126 double p_cap_ip;
128
130
131 if (medium.hasProperty(MPL::PropertyType::saturation))
132 {
133 variables.capillary_pressure = p_cap_ip;
134 variables.liquid_phase_pressure = -p_cap_ip;
135 // setting pG to 1 atm
136 // TODO : rewrite equations s.t. p_L = pG-p_cap
137 variables.gas_phase_pressure = 1.0e5;
138 variables.temperature = T_ip;
139
140 double const S_L =
142 .template value<double>(variables, x_position, t, dt);
144 S_L;
145
146 variables.liquid_saturation = S_L;
147 }
148
150 {T_ip, 0, {}}, this->current_states_[ip],
151 this->prev_states_[ip]);
152
153 if (this->process_data_.initial_stress.value)
154 {
155 if (!medium.hasProperty(MPL::PropertyType::saturation))
156 {
157 OGS_FATAL(
158 "The medium has no saturation property required to compute "
159 "initial stress.");
160 }
162 -p_cap_ip);
163 }
164 }
165}
#define OGS_FATAL(...)
Definition Error.h:26
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 44 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 67 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 73 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: