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 30 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 49 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 50 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 62 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 56 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 54 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 46 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 64 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 31 of file ThermoRichardsMechanicsFEM-impl.h.

41{
42 unsigned const n_integration_points =
43 this->integration_method_.getNumberOfPoints();
44
46
47 auto const shape_matrices_u =
52
53 auto const shape_matrices =
57
58 for (unsigned ip = 0; ip < n_integration_points; ip++)
59 {
60 auto& ip_data = ip_data_[ip];
61 auto const& sm_u = shape_matrices_u[ip];
62 ip_data_[ip].integration_weight =
63 integration_method.getWeightedPoint(ip).getWeight() *
64 sm_u.integralMeasure * sm_u.detJ;
65
66 ip_data.N_u = sm_u.N;
67 ip_data.dNdx_u = sm_u.dNdx;
68
69 // ip_data.N_p and ip_data.dNdx_p are used for both p and T variables
70 ip_data.N_p = shape_matrices[ip].N;
71 ip_data.dNdx_p = shape_matrices[ip].dNdx;
72 }
73}
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 273 of file ThermoRichardsMechanicsFEM-impl.h.

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

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

211{
212 auto& medium =
213 *this->process_data_.media_map.getMedium(this->element_.getID());
214
216 loc_mat.setZero();
218 loc_mat_current_ip.setZero(); // only to set the right matrix sizes
219
221
222 for (unsigned ip = 0; ip < this->integration_method_.getNumberOfPoints();
223 ++ip)
224 {
226 std::nullopt, this->element_.getID(),
230 this->element_, ip_data_[ip].N_u))};
231
233 t, dt, x_position, //
236 medium, //
238 this->current_states_[ip], this->prev_states_[ip],
239 this->material_states_[ip], this->output_data_[ip]);
241 }
242
244
247}
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 333 of file ThermoRichardsMechanicsFEM-impl.h.

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

384 {
387 }

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

413 {
415 }

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 393 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 398 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 527 of file ThermoRichardsMechanicsFEM-impl.h.

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

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

170{
176 this->process_data_.initial_stress.type ==
178 {
179 return;
180 }
181
183 this->process_data_.initial_stress.type == InitialStress::Type::Total)
184 {
185 return;
186 }
187
188 double const alpha_b =
190 .template value<double>(variables, x_position, t, 0.0 /*dt*/);
191
192 double const bishop =
194 .template value<double>(variables, x_position, t, 0.0 /*dt*/);
195
197 this->current_states_[ip], this->prev_states_[ip],
199}
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 243 of file ThermoRichardsMechanicsFEM.h.

244 {
245 return displacement_size;
246 }

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

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

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

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

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

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

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

83{
84 assert(local_x.size() ==
86
87 auto const p_L = local_x.template segment<pressure_size>(pressure_index);
88 auto const T =
90
91 constexpr double dt = std::numeric_limits<double>::quiet_NaN();
92 auto const& medium =
93 *this->process_data_.media_map.getMedium(this->element_.getID());
94
96
99 this->process_data_, this->solid_material_);
100
101 unsigned const n_integration_points =
102 this->integration_method_.getNumberOfPoints();
103 for (unsigned ip = 0; ip < n_integration_points; ip++)
104 {
105 // N is used for both T and p variables.
106 auto const& N = ip_data_[ip].N_p;
107
109 std::nullopt, this->element_.getID(),
113 this->element_, ip_data_[ip].N_u))};
114
115 double T_ip;
117
118 double p_cap_ip;
120
122
123 if (medium.hasProperty(MPL::PropertyType::saturation))
124 {
125 variables.capillary_pressure = p_cap_ip;
126 variables.liquid_phase_pressure = -p_cap_ip;
127 // setting pG to 1 atm
128 // TODO : rewrite equations s.t. p_L = pG-p_cap
129 variables.gas_phase_pressure = 1.0e5;
130 variables.temperature = T_ip;
131
132 double const S_L =
134 .template value<double>(variables, x_position, t, dt);
136 S_L;
137
138 variables.liquid_saturation = S_L;
139 }
140
142 {T_ip, 0, {}}, this->current_states_[ip],
143 this->prev_states_[ip]);
144
145 if (this->process_data_.initial_stress.value)
146 {
147 if (!medium.hasProperty(MPL::PropertyType::saturation))
148 {
149 OGS_FATAL(
150 "The medium has no saturation property required to compute "
151 "initial stress.");
152 }
154 -p_cap_ip);
155 }
156 }
157}
#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 37 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 60 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 66 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: