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 > &, std::vector< double > &, 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_M_data, std::vector< double > &local_K_data, 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
Initial value:
BMatrixPolicyType<ShapeFunctionDisplacement, DisplacementDim>

Definition at line 59 of file ThermoRichardsMechanicsFEM.h.

◆ 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:
IntegrationPointData<ShapeMatricesTypeDisplacement,
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

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

48 : LocalAssemblerInterface<DisplacementDim, ConstitutiveTraits>(
49 e, integration_method, is_axially_symmetric, process_data)
50{
51 unsigned const n_integration_points =
53
54 ip_data_.resize(n_integration_points);
55
56 auto const shape_matrices_u =
57 NumLib::initShapeMatrices<ShapeFunctionDisplacement,
59 DisplacementDim>(e, is_axially_symmetric,
60 integration_method);
61
62 auto const shape_matrices =
64 DisplacementDim>(e, is_axially_symmetric,
65 integration_method);
66
67 for (unsigned ip = 0; ip < n_integration_points; ip++)
68 {
69 auto& ip_data = ip_data_[ip];
70 auto const& sm_u = shape_matrices_u[ip];
71 ip_data_[ip].integration_weight =
72 integration_method.getWeightedPoint(ip).getWeight() *
73 sm_u.integralMeasure * sm_u.detJ;
74
75 ip_data.N_u = sm_u.N;
76 ip_data.dNdx_u = sm_u.dNdx;
77
78 // ip_data.N_p and ip_data.dNdx_p are used for both p and T variables
79 ip_data.N_p = shape_matrices[ip].N;
80 ip_data.dNdx_p = shape_matrices[ip].dNdx;
81 }
82}
std::vector< typename ShapeMatricesType::ShapeMatrices, Eigen::aligned_allocator< typename ShapeMatricesType::ShapeMatrices > > initShapeMatrices(MeshLib::Element const &e, bool const is_axially_symmetric, IntegrationMethod const &integration_method)

References NumLib::GenericIntegrationMethod::getNumberOfPoints(), MathLib::WeightedPoint::getWeight(), NumLib::GenericIntegrationMethod::getWeightedPoint(), NumLib::initShapeMatrices(), ProcessLib::ThermoRichardsMechanics::LocalAssemblerInterface< DisplacementDim, ConstitutiveTraits >::integration_method_, and ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::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 272 of file ThermoRichardsMechanicsFEM-impl.h.

282{
283 constexpr auto local_matrix_dim =
285 assert(local_x.size() == local_matrix_dim);
286
287 auto local_Jac = MathLib::createZeroedMatrix<
288 typename ShapeMatricesTypeDisplacement::template MatrixType<
289 local_matrix_dim, local_matrix_dim>>(
290 local_Jac_data, local_matrix_dim, local_matrix_dim);
291
292 auto local_rhs =
293 MathLib::createZeroedVector<typename ShapeMatricesTypeDisplacement::
294 template VectorType<local_matrix_dim>>(
295 local_rhs_data, local_matrix_dim);
296
297 local_Jac.noalias() = loc_mat.Jac;
298 local_rhs.noalias() = -loc_mat.res;
299
300 //
301 // -- Jacobian
302 //
303 block_TT(local_Jac).noalias() += loc_mat.M_TT / dt + loc_mat.K_TT;
304 block_Tp(local_Jac).noalias() +=
305 loc_mat.K_Tp + loc_mat.dK_TT_dp + loc_mat.M_Tp / dt;
306
307 block_pT(local_Jac).noalias() += loc_mat.M_pT / dt + loc_mat.K_pT;
308 block_pp(local_Jac).noalias() +=
309 loc_mat.K_pp + loc_mat.storage_p_a_p / dt + loc_mat.storage_p_a_S_Jpp;
310 block_pu(local_Jac).noalias() = loc_mat.M_pu / dt;
311
312 //
313 // -- Residual
314 //
315 auto const [T, p_L, u] = localDOF(local_x);
316 auto const [T_prev, p_L_prev, u_prev] = localDOF(local_x_prev);
317
318 block_T(local_rhs).noalias() -= loc_mat.M_TT * (T - T_prev) / dt +
319 loc_mat.K_TT * T + loc_mat.K_Tp * p_L +
320 loc_mat.M_Tp * (p_L - p_L_prev) / dt;
321 block_p(local_rhs).noalias() -=
322 loc_mat.K_pp * p_L + loc_mat.K_pT * T +
323 (loc_mat.storage_p_a_p + loc_mat.storage_p_a_S) * (p_L - p_L_prev) /
324 dt +
325 loc_mat.M_pu * (u - u_prev) / dt + loc_mat.M_pT * (T - T_prev) / dt;
326}
Eigen::Map< Vector > createZeroedVector(std::vector< double > &data, Eigen::VectorXd::Index size)
Eigen::Map< Matrix > createZeroedMatrix(std::vector< double > &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)

References MathLib::createZeroedMatrix(), and MathLib::createZeroedVector().

◆ 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 > & ,
std::vector< double > & ,
std::vector< double > & local_rhs_data,
std::vector< double > & local_Jac_data )
overridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

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

210{
211 auto& medium =
212 *this->process_data_.media_map.getMedium(this->element_.getID());
213
214 LocalMatrices loc_mat;
215 loc_mat.setZero();
216 LocalMatrices loc_mat_current_ip;
217 loc_mat_current_ip.setZero(); // only to set the right matrix sizes
218
219 typename ConstitutiveTraits::ConstitutiveSetting constitutive_setting;
220
221 for (unsigned ip = 0; ip < this->integration_method_.getNumberOfPoints();
222 ++ip)
223 {
224 ParameterLib::SpatialPosition const x_position{
225 std::nullopt, this->element_.getID(), ip,
227 NumLib::interpolateCoordinates<ShapeFunctionDisplacement,
229 this->element_, ip_data_[ip].N_u))};
230
232 t, dt, x_position, //
233 local_x, local_x_prev, //
234 ip_data_[ip], constitutive_setting,
235 medium, //
236 loc_mat_current_ip, //
237 this->current_states_[ip], this->prev_states_[ip],
238 this->material_states_[ip], this->output_data_[ip]);
239 loc_mat += loc_mat_current_ip;
240 }
241
242 massLumping(loc_mat);
243
244 addToLocalMatrixData(dt, local_x, local_x_prev, loc_mat, local_rhs_data,
245 local_Jac_data);
246}
std::size_t getID() const
Returns the ID of the element.
Definition Element.h:89
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::array< double, 3 > interpolateCoordinates(MeshLib::Element const &e, typename ShapeMatricesType::ShapeMatrices::ShapeType const &N)
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 NumLib::interpolateCoordinates(), 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 332 of file ThermoRichardsMechanicsFEM-impl.h.

350{
351 auto const& N_u = ip_data.N_u;
352 auto const& dNdx_u = ip_data.dNdx_u;
353
354 // N and dNdx are used for both p and T variables
355 auto const& N = ip_data.N_p;
356 auto const& dNdx = ip_data.dNdx_p;
357
358 auto const B =
359 LinearBMatrix::computeBMatrix<DisplacementDim,
360 ShapeFunctionDisplacement::NPOINTS,
362 dNdx_u, N_u, (*x_position.getCoordinates())[0],
364
365 auto const [T, p_L, u] = localDOF(local_x);
366 auto const [T_prev, p_L_prev, u_prev] = localDOF(local_x_prev);
367
368 GlobalDimVectorType const grad_T_ip = dNdx * T;
369
370 typename ConstitutiveTraits::ConstitutiveModels models(
371 this->process_data_, this->solid_material_);
372 typename ConstitutiveTraits::ConstitutiveTempData tmp;
373 typename ConstitutiveTraits::ConstitutiveData CD;
374
375 {
376 double const T_ip = N * T;
377 double const T_prev_ip = N * T_prev;
378
379 double const p_cap_ip = -N * p_L;
380 double const p_cap_prev_ip = -N * p_L_prev;
381 GlobalDimVectorType const grad_p_cap_ip = -dNdx * p_L;
382
383 KelvinVectorType eps = B * u;
384
385 CS.eval(models, t, dt, x_position, //
386 medium, //
387 {T_ip, T_prev_ip, grad_T_ip}, //
388 {p_cap_ip, p_cap_prev_ip, grad_p_cap_ip}, //
389 eps, current_state, prev_state, mat_state, tmp, output_data,
390 CD);
391 }
392
393 using NodalMatrix = typename ShapeMatricesType::NodalMatrixType;
394 NodalMatrix const NTN = N.transpose() * N;
395 NodalMatrix const dNTdN = dNdx.transpose() * dNdx;
396
397 // TODO is identity2.transpose() * B something like divergence?
398 auto const& identity2 = MathLib::KelvinVector::Invariants<
400 DisplacementDim)>::identity2;
401 typename ShapeMatricesTypeDisplacement::template MatrixType<
403 BTI2N = B.transpose() * identity2 * N;
404
405 /*
406 * Conventions:
407 *
408 * * use positive signs exclusively, any negative signs must be included in
409 * the coefficients coming from the constitutive setting
410 * * the used coefficients from the constitutive setting are named after the
411 * terms they appear in, e.g. K_TT_X_dNTdN means it is a scalar (X) that
412 * will be multiplied by dNdx.transpose() * dNdx. Placefolders for the
413 * coefficients are:
414 * * X -> scalar
415 * * V -> vector
416 * * K -> Kelvin vector
417 * * the Laplace terms have a special name, e.g., K_TT_Laplace
418 * * there shall be only one contribution to each of the LocalMatrices,
419 * assigned with = assignment; this point might be relaxed in the future
420 * * this method will overwrite the data in the passed LocalMatrices& out
421 * argument, not add to it
422 */
423
424 // residual, order T, p, u
425 block_p(out.res).noalias() =
426 dNdx.transpose() * std::get<EqPData<DisplacementDim>>(CD).rhs_p_dNT_V;
427 block_u(out.res).noalias() =
428 B.transpose() *
429 ProcessLib::Graph::get<TotalStressData<DisplacementDim>>(
430 CD, current_state)
431 .sigma_total -
432 static_cast<int>(this->process_data_.apply_body_force_for_deformation) *
433 N_u_op(N_u).transpose() *
434 std::get<GravityData<DisplacementDim>>(CD).volumetric_body_force;
435
436 // Storage matrices
437 out.storage_p_a_p.noalias() =
438 std::get<EqPData<DisplacementDim>>(CD).storage_p_a_p_X_NTN * NTN;
439 out.storage_p_a_S.noalias() =
440 std::get<TRMStorageData>(CD).storage_p_a_S_X_NTN * NTN;
441 out.storage_p_a_S_Jpp.noalias() =
442 std::get<TRMStorageData>(CD).storage_p_a_S_Jpp_X_NTN * NTN;
443
444 // M matrices, order T, p, u
445 out.M_TT.noalias() =
446 std::get<EqTData<DisplacementDim>>(CD).M_TT_X_NTN * NTN;
447 out.M_Tp.noalias() =
448 std::get<TRMVaporDiffusionData<DisplacementDim>>(CD).M_Tp_X_NTN * NTN;
449
450 out.M_pT.noalias() =
451 std::get<EqPData<DisplacementDim>>(CD).M_pT_X_NTN * NTN;
452 out.M_pu.noalias() =
453 std::get<EqPData<DisplacementDim>>(CD).M_pu_X_BTI2N * BTI2N.transpose();
454
455 // K matrices, order T, p, u
456 out.K_TT.noalias() =
457 dNdx.transpose() *
458 std::get<TRMHeatStorageAndFluxData<DisplacementDim>>(CD)
459 .K_TT_Laplace *
460 dNdx +
461 N.transpose() *
462 (std::get<EqTData<DisplacementDim>>(CD).K_TT_NT_V_dN.transpose() *
463 dNdx) +
464 std::get<TRMVaporDiffusionData<DisplacementDim>>(CD).K_TT_X_dNTdN *
465 dNTdN;
466
467 out.dK_TT_dp.noalias() =
468 N.transpose() *
469 (std::get<TRMHeatStorageAndFluxData<DisplacementDim>>(CD)
470 .K_Tp_NT_V_dN.transpose() *
471 dNdx) +
472 std::get<TRMHeatStorageAndFluxData<DisplacementDim>>(CD).K_Tp_X_NTN *
473 NTN;
474 out.K_Tp.noalias() =
475 dNdx.transpose() *
476 std::get<ThermoOsmosisData<DisplacementDim>>(CD).K_Tp_Laplace *
477 dNdx +
478 std::get<TRMVaporDiffusionData<DisplacementDim>>(CD).K_Tp_X_dNTdN *
479 dNTdN;
480
481 out.K_pp.noalias() =
482 dNdx.transpose() * std::get<EqPData<DisplacementDim>>(CD).K_pp_Laplace *
483 dNdx +
484 std::get<TRMVaporDiffusionData<DisplacementDim>>(CD).K_pp_X_dNTdN *
485 dNTdN;
486 out.K_pT.noalias() =
487 dNdx.transpose() *
488 std::get<ThermoOsmosisData<DisplacementDim>>(CD).K_pT_Laplace * dNdx;
489
490 // direct Jacobian contributions, order T, p, u
491 block_pT(out.Jac).noalias() =
492 std::get<TRMVaporDiffusionData<DisplacementDim>>(CD).J_pT_X_dNTdN *
493 dNTdN;
494 block_pp(out.Jac).noalias() =
495 std::get<TRMStorageData>(CD).J_pp_X_NTN * NTN +
496 std::get<EqPData<DisplacementDim>>(CD).J_pp_X_BTI2NT_u_dot_N *
497 BTI2N.transpose() * (u - u_prev) / dt *
498 N // TODO something with volumetric strain rate?
499 + dNdx.transpose() *
500 std::get<EqPData<DisplacementDim>>(CD).J_pp_dNT_V_N * N;
501
502 block_uT(out.Jac).noalias() =
503 B.transpose() *
504 std::get<SolidMechanicsDataStateless<DisplacementDim>>(CD).J_uT_BT_K_N *
505 N;
506 block_up(out.Jac).noalias() =
507 B.transpose() *
508 std::get<SolidMechanicsDataStateless<DisplacementDim>>(CD)
509 .J_up_BT_K_N *
510 N +
511 N_u_op(N_u).transpose() *
512 std::get<GravityData<DisplacementDim>>(CD).J_up_HT_V_N * N;
513 block_uu(out.Jac).noalias() =
514 B.transpose() *
515 std::get<SolidMechanicsDataStateless<DisplacementDim>>(CD)
516 .stiffness_tensor *
517 B;
518
519 out *= ip_data.integration_weight;
520}
MatrixType< _kelvin_vector_size, _number_of_dof > BMatrixType
constexpr int kelvin_vector_dimensions(int const displacement_dim)
Kelvin vector dimensions for given displacement dimension.
BMatrixType computeBMatrix(DNDX_Type const &dNdx, N_Type const &N, const double radius, const bool is_axially_symmetric)
Fills a B-matrix based on given shape function dN/dx values.
MatrixType< ShapeFunction::NPOINTS, ShapeFunction::NPOINTS > NodalMatrixType
ConstitutiveTraits::SolidConstitutiveRelation const & solid_material_

References ProcessLib::LinearBMatrix::computeBMatrix(), ParameterLib::SpatialPosition::getCoordinates(), ProcessLib::ThermoRichardsMechanics::TRMVaporDiffusionData< DisplacementDim >::K_pp_X_dNTdN, ProcessLib::ThermoRichardsMechanics::TRMVaporDiffusionData< DisplacementDim >::K_Tp_X_dNTdN, ProcessLib::ThermoRichardsMechanics::TRMHeatStorageAndFluxData< DisplacementDim >::K_Tp_X_NTN, ProcessLib::ThermoRichardsMechanics::TRMVaporDiffusionData< DisplacementDim >::K_TT_X_dNTdN, and MathLib::KelvinVector::kelvin_vector_dimensions().

◆ block_p()

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

◆ block_pp()

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

◆ block_pT()

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

◆ block_pu()

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

◆ block_T()

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

◆ block_Tp()

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

◆ block_TT()

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

◆ block_u()

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

◆ block_up()

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

◆ block_uT()

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

◆ block_uu()

template<typename ShapeFunctionDisplacement , typename ShapeFunction , int DisplacementDim, typename ConstitutiveTraits >
static 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 526 of file ThermoRichardsMechanicsFEM-impl.h.

530{
531 auto const T = block_T(local_x);
532 auto const p_L = block_p(local_x);
533 auto const u = block_u(local_x);
534
535 auto const T_prev = block_T(local_x_prev);
536 auto const p_L_prev = block_p(local_x_prev);
537
538 auto const e_id = this->element_.getID();
539 auto const& process_data = this->process_data_;
540 auto& medium = *process_data.media_map.getMedium(e_id);
541
542 unsigned const n_integration_points =
544
545 double saturation_avg = 0;
546 double porosity_avg = 0;
547 double liquid_density_avg = 0;
548 double viscosity_avg = 0;
549
551 KV sigma_avg = KV::Zero();
552
553 typename ConstitutiveTraits::ConstitutiveSetting constitutive_setting;
554
555 typename ConstitutiveTraits::ConstitutiveModels models(
556 process_data, this->solid_material_);
557 typename ConstitutiveTraits::ConstitutiveTempData tmp;
558 typename ConstitutiveTraits::ConstitutiveData CD;
559
560 for (unsigned ip = 0; ip < n_integration_points; ip++)
561 {
562 auto& current_state = this->current_states_[ip];
563 auto& output_data = this->output_data_[ip];
564
565 auto const& ip_data = ip_data_[ip];
566
567 // N is used for both p and T variables
568 auto const& N = ip_data.N_p;
569 auto const& N_u = ip_data.N_u;
570 auto const& dNdx_u = ip_data.dNdx_u;
571 auto const& dNdx = ip_data.dNdx_p;
572
573 ParameterLib::SpatialPosition const x_position{
574 std::nullopt, this->element_.getID(), ip,
576 NumLib::interpolateCoordinates<ShapeFunctionDisplacement,
578 this->element_, N_u))};
579 auto const x_coord =
580 NumLib::interpolateXCoordinate<ShapeFunctionDisplacement,
582 this->element_, N_u);
583 auto const B =
584 LinearBMatrix::computeBMatrix<DisplacementDim,
585 ShapeFunctionDisplacement::NPOINTS,
587 dNdx_u, N_u, x_coord, this->is_axially_symmetric_);
588
589 double const T_ip = N * T;
590 double const T_prev_ip = N * T_prev;
591 GlobalDimVectorType const grad_T_ip = dNdx * T;
592
593 double const p_cap_ip = -N * p_L;
594 double const p_cap_prev_ip = -N * p_L_prev;
595 GlobalDimVectorType const grad_p_cap_ip = -dNdx * p_L;
596
597 KelvinVectorType eps = B * u;
598
599 constitutive_setting.eval(models, //
600 t, dt, x_position, //
601 medium, //
602 {T_ip, T_prev_ip, grad_T_ip}, //
603 {p_cap_ip, p_cap_prev_ip, grad_p_cap_ip}, //
604 eps, current_state, this->prev_states_[ip],
605 this->material_states_[ip], tmp, output_data,
606 CD);
607
608 saturation_avg += std::get<SaturationData>(current_state).S_L;
609 porosity_avg += std::get<PorosityData>(current_state).phi;
610
611 liquid_density_avg += std::get<LiquidDensityData>(output_data).rho_LR;
612 viscosity_avg += std::get<LiquidViscosityData>(output_data).viscosity;
613 sigma_avg += ConstitutiveTraits::ConstitutiveSetting::statefulStress(
614 current_state);
615 }
616 saturation_avg /= n_integration_points;
617 porosity_avg /= n_integration_points;
618 viscosity_avg /= n_integration_points;
619 liquid_density_avg /= n_integration_points;
620 sigma_avg /= n_integration_points;
621
622 (*process_data.element_saturation)[e_id] = saturation_avg;
623 (*process_data.element_porosity)[e_id] = porosity_avg;
624 (*process_data.element_liquid_density)[e_id] = liquid_density_avg;
625 (*process_data.element_viscosity)[e_id] = viscosity_avg;
626
627 Eigen::Map<KV>(
628 &(*process_data.element_stresses)[e_id * KV::RowsAtCompileTime]) =
630
632 ShapeFunction, typename ShapeFunctionDisplacement::MeshElement,
633 DisplacementDim>(this->element_, this->is_axially_symmetric_, p_L,
634 *process_data.pressure_interpolated);
636 ShapeFunction, typename ShapeFunctionDisplacement::MeshElement,
637 DisplacementDim>(this->element_, this->is_axially_symmetric_, T,
638 *process_data.temperature_interpolated);
639}
Eigen::Matrix< double, 4, 1 > kelvinVectorToSymmetricTensor(Eigen::Matrix< double, 4, 1, Eigen::ColMajor, 4, 1 > const &v)
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), 1, Eigen::ColMajor > KelvinVectorType
double interpolateXCoordinate(MeshLib::Element const &e, typename ShapeMatricesType::ShapeMatrices::ShapeType const &N)
void interpolateToHigherOrderNodes(MeshLib::Element const &element, bool const is_axially_symmetric, Eigen::MatrixBase< EigenMatrixType > const &node_values, MeshLib::PropertyVector< double > &interpolated_values_global_vector)

References ProcessLib::LinearBMatrix::computeBMatrix(), NumLib::interpolateCoordinates(), NumLib::interpolateToHigherOrderNodes(), NumLib::interpolateXCoordinate(), and MathLib::KelvinVector::kelvinVectorToSymmetricTensor().

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

167{
168 bool constexpr is_strain_temperature_constitutive =
169 std::is_same<ConstitutiveStress_StrainTemperature::ConstitutiveTraits<
170 DisplacementDim>,
171 ConstitutiveTraits>::value;
172 if (is_strain_temperature_constitutive &&
173 this->process_data_.initial_stress.type ==
175 {
176 return;
177 }
178
179 if (!is_strain_temperature_constitutive &&
180 this->process_data_.initial_stress.type == InitialStress::Type::Total)
181 {
182 return;
183 }
184
185 double const alpha_b =
186 medium.property(MPL::PropertyType::biot_coefficient)
187 .template value<double>(variables, x_position, t, 0.0 /*dt*/);
188
189 double const bishop =
190 medium.property(MPL::PropertyType::bishops_effective_stress)
191 .template value<double>(variables, x_position, t, 0.0 /*dt*/);
192
193 ConstitutiveTraits::ConstitutiveSetting::convertInitialStressType(
194 this->current_states_[ip], this->prev_states_[ip],
195 bishop * alpha_b * p_at_ip * Invariants::identity2);
196}
static Eigen::Matrix< double, KelvinVectorSize, 1 > const identity2
Kelvin mapping of 2nd order identity tensor.

References ProcessLib::InitialStress::Effective, MaterialPropertyLib::Medium::property(), and ProcessLib::InitialStress::Total.

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

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

References ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::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 275 of file ThermoRichardsMechanicsFEM.h.

276 {
277 unsigned const n_integration_points =
279 auto const time_independent = std::numeric_limits<double>::quiet_NaN();
280 auto const& medium =
281 *this->process_data_.media_map.getMedium(this->element_.getID());
282
283 for (unsigned ip = 0; ip < n_integration_points; ip++)
284 {
285 ParameterLib::SpatialPosition const x_position{
286 std::nullopt, this->element_.getID(), ip,
288 ShapeFunctionDisplacement,
290 this->element_, this->ip_data_[ip].N_u))};
291 auto& current_state = this->current_states_[ip];
292
293 // Set initial stress from parameter.
294 if (this->process_data_.initial_stress.value)
295 {
296 ConstitutiveTraits::ConstitutiveSetting::statefulStress(
297 current_state) =
299 DisplacementDim>(
300 (*this->process_data_.initial_stress.value)(
301 time_independent, x_position));
302 }
303
304 if (*this->process_data_.initialize_porosity_from_medium_property)
305 {
306 // Initial porosity. Could be read from integration point data
307 // or mesh.
308 std::get<PorosityData>(current_state).phi =
309 medium.property(MPL::porosity)
310 .template initialValue<double>(x_position,
311 time_independent);
312
313 if (medium.hasProperty(MPL::PropertyType::transport_porosity))
314 {
315 std::get<TransportPorosityData>(current_state).phi =
316 medium.property(MPL::transport_porosity)
317 .template initialValue<double>(x_position,
318 time_independent);
319 }
320 else
321 {
322 std::get<TransportPorosityData>(current_state).phi =
323 std::get<PorosityData>(current_state).phi;
324 }
325 }
326
327 double const t = 0; // TODO (naumov) pass t from top
328 this->solid_material_.initializeInternalStateVariables(
329 t, x_position,
330 *this->material_states_[ip].material_state_variables);
331 }
332
333 for (unsigned ip = 0; ip < n_integration_points; ip++)
334 {
335 this->material_states_[ip].pushBackState();
336 }
337
338 for (unsigned ip = 0; ip < n_integration_points; ip++)
339 {
340 this->prev_states_[ip] = this->current_states_[ip];
341 }
342 }
Eigen::Matrix< double, Eigen::MatrixBase< Derived >::RowsAtCompileTime, 1 > symmetricTensorToKelvinVector(Eigen::MatrixBase< Derived > const &v)

References ProcessLib::ThermoRichardsMechanics::LocalAssemblerInterface< DisplacementDim, ConstitutiveTraits >::current_states_, ProcessLib::ThermoRichardsMechanics::LocalAssemblerInterface< DisplacementDim, ConstitutiveTraits >::element_, MeshLib::Element::getID(), NumLib::GenericIntegrationMethod::getNumberOfPoints(), ProcessLib::ThermoRichardsMechanics::LocalAssemblerInterface< DisplacementDim, ConstitutiveTraits >::integration_method_, NumLib::interpolateCoordinates(), ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::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 >
static constexpr auto ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::localDOF ( std::vector< double > const & x)
inlinestaticconstexprprivate

Definition at line 360 of file ThermoRichardsMechanicsFEM.h.

361 {
362 return NumLib::localDOF<
363 ShapeFunction, ShapeFunction,
365 }
auto localDOF(ElementDOFVector const &x)
Definition LocalDOF.h:64

References NumLib::localDOF().

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

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

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

92{
93 assert(local_x.size() ==
95
96 auto const p_L = local_x.template segment<pressure_size>(pressure_index);
97 auto const T =
98 local_x.template segment<temperature_size>(temperature_index);
99
100 constexpr double dt = std::numeric_limits<double>::quiet_NaN();
101 auto const& medium =
102 *this->process_data_.media_map.getMedium(this->element_.getID());
103
104 MediaData const media_data{medium};
105
106 typename ConstitutiveTraits::ConstitutiveSetting const constitutive_setting;
107 typename ConstitutiveTraits::ConstitutiveModels models(
108 this->process_data_, this->solid_material_);
109
110 unsigned const n_integration_points =
112 for (unsigned ip = 0; ip < n_integration_points; ip++)
113 {
114 // N is used for both T and p variables.
115 auto const& N = ip_data_[ip].N_p;
116
117 ParameterLib::SpatialPosition const x_position{
118 std::nullopt, this->element_.getID(), ip,
120 NumLib::interpolateCoordinates<ShapeFunctionDisplacement,
122 this->element_, ip_data_[ip].N_u))};
123
124 double T_ip;
126
127 double p_cap_ip;
128 NumLib::shapeFunctionInterpolate(-p_L, N, p_cap_ip);
129
130 MPL::VariableArray variables;
131 variables.capillary_pressure = p_cap_ip;
132 variables.liquid_phase_pressure = -p_cap_ip;
133 // setting pG to 1 atm
134 // TODO : rewrite equations s.t. p_L = pG-p_cap
135 variables.gas_phase_pressure = 1.0e5;
136 variables.temperature = T_ip;
137
138 double const S_L =
139 medium.property(MPL::PropertyType::saturation)
140 .template value<double>(variables, x_position, t, dt);
141 std::get<PrevState<SaturationData>>(this->prev_states_[ip])->S_L = S_L;
142
143 constitutive_setting.init(models, t, dt, x_position, media_data,
144 {T_ip, 0, {}}, this->current_states_[ip],
145 this->prev_states_[ip]);
146
147 if (this->process_data_.initial_stress.value)
148 {
149 variables.liquid_saturation = S_L;
150 convertInitialStressType(ip, t, x_position, medium, variables,
151 -p_cap_ip);
152 }
153 }
154}
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, MaterialPropertyLib::VariableArray::gas_phase_pressure, NumLib::interpolateCoordinates(), MaterialPropertyLib::VariableArray::liquid_phase_pressure, MaterialPropertyLib::VariableArray::liquid_saturation, NumLib::detail::shapeFunctionInterpolate(), and MaterialPropertyLib::VariableArray::temperature.

Member Data Documentation

◆ displacement_index

◆ displacement_size

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

◆ ip_data_

◆ 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 >
constexpr 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.

◆ pressure_index

◆ pressure_size

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

◆ temperature_index

◆ temperature_size

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

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