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

Detailed Description

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

Definition at line 37 of file ThermoRichardsMechanicsFEM.h.

#include <ThermoRichardsMechanicsFEM.h>

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

Classes

class  LocalMatrices
 

Public Types

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

Public Member Functions

 ThermoRichardsMechanicsLocalAssembler (ThermoRichardsMechanicsLocalAssembler const &)=delete
 
 ThermoRichardsMechanicsLocalAssembler (ThermoRichardsMechanicsLocalAssembler &&)=delete
 
 ThermoRichardsMechanicsLocalAssembler (MeshLib::Element const &e, std::size_t const, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, ThermoRichardsMechanicsProcessData< DisplacementDim, ConstitutiveTraits > &process_data)
 
void setInitialConditionsConcrete (Eigen::VectorXd const local_x, double const t, int const process_id) override
 
void assembleWithJacobian (double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_x_prev, std::vector< double > &local_rhs_data, std::vector< double > &local_Jac_data) override
 
void initializeConcrete () override
 
void computeSecondaryVariableConcrete (double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev) override
 
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix (const unsigned integration_point) const override
 Provides the shape matrix at the given integration point.
 
- Public Member Functions inherited from ProcessLib::ThermoRichardsMechanics::LocalAssemblerInterface< DisplacementDim, ConstitutiveTraits >
 LocalAssemblerInterface (MeshLib::Element const &e, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, ThermoRichardsMechanicsProcessData< DisplacementDim, ConstitutiveTraits > &process_data)
 
std::size_t setIPDataInitialConditions (std::string_view name, double const *values, int const integration_order)
 
std::vector< double > getMaterialStateVariableInternalState (std::function< std::span< double >(typename MaterialLib::Solids::MechanicsBase< DisplacementDim >::MaterialStateVariables &)> const &get_values_span, int const &n_components) const
 
unsigned getNumberOfIntegrationPoints () const
 
int getMaterialID () const
 
MaterialLib::Solids::MechanicsBase< DisplacementDim >::MaterialStateVariables const & getMaterialStateVariablesAt (unsigned integration_point) const
 
void postTimestepConcrete (Eigen::VectorXd const &, Eigen::VectorXd const &, double const, double const, int const) override
 
- Public Member Functions inherited from ProcessLib::LocalAssemblerInterface
virtual ~LocalAssemblerInterface ()=default
 
virtual void setInitialConditions (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, double const t, int const process_id)
 
virtual void initialize (std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
 
virtual void preAssemble (double const, double const, std::vector< double > const &)
 
virtual void assemble (double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_x_prev, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data)
 
virtual void assembleForStaggeredScheme (double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, int const process_id, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data)
 
virtual void assembleWithJacobianForStaggeredScheme (double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, int const process_id, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data)
 
virtual void computeSecondaryVariable (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_prev, int const process_id)
 
virtual void preTimestep (std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, double const t, double const delta_t)
 
virtual void postTimestep (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, int const process_id)
 
void postNonLinearSolver (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, int const process_id)
 
virtual Eigen::Vector3d getFlux (MathLib::Point3d const &, double const, std::vector< double > const &) const
 
virtual Eigen::Vector3d getFlux (MathLib::Point3d const &, double const, std::vector< std::vector< double > > const &) const
 Fits to staggered scheme.
 
- 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 39 of file ThermoRichardsMechanicsFEM-impl.h.

47 : LocalAssemblerInterface<DisplacementDim, ConstitutiveTraits>(
48 e, integration_method, is_axially_symmetric, process_data)
49{
50 unsigned const n_integration_points =
52
53 ip_data_.resize(n_integration_points);
54
55 auto const shape_matrices_u =
56 NumLib::initShapeMatrices<ShapeFunctionDisplacement,
58 DisplacementDim>(e, is_axially_symmetric,
59 integration_method);
60
61 auto const shape_matrices =
63 DisplacementDim>(e, is_axially_symmetric,
64 integration_method);
65
66 for (unsigned ip = 0; ip < n_integration_points; ip++)
67 {
68 auto& ip_data = ip_data_[ip];
69 auto const& sm_u = shape_matrices_u[ip];
70 ip_data_[ip].integration_weight =
71 integration_method.getWeightedPoint(ip).getWeight() *
72 sm_u.integralMeasure * sm_u.detJ;
73
74 ip_data.N_u = sm_u.N;
75 ip_data.dNdx_u = sm_u.dNdx;
76
77 // ip_data.N_p and ip_data.dNdx_p are used for both p and T variables
78 ip_data.N_p = shape_matrices[ip].N;
79 ip_data.dNdx_p = shape_matrices[ip].dNdx;
80 }
81}
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 269 of file ThermoRichardsMechanicsFEM-impl.h.

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

Reimplemented from ProcessLib::LocalAssemblerInterface.

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

207{
208 auto& medium =
209 *this->process_data_.media_map.getMedium(this->element_.getID());
210
211 LocalMatrices loc_mat;
212 loc_mat.setZero();
213 LocalMatrices loc_mat_current_ip;
214 loc_mat_current_ip.setZero(); // only to set the right matrix sizes
215
216 typename ConstitutiveTraits::ConstitutiveSetting constitutive_setting;
217
218 for (unsigned ip = 0; ip < this->integration_method_.getNumberOfPoints();
219 ++ip)
220 {
221 ParameterLib::SpatialPosition const x_position{
222 std::nullopt, this->element_.getID(), ip,
224 NumLib::interpolateCoordinates<ShapeFunctionDisplacement,
226 this->element_, ip_data_[ip].N_u))};
227
229 t, dt, x_position, //
230 local_x, local_x_prev, //
231 ip_data_[ip], constitutive_setting,
232 medium, //
233 loc_mat_current_ip, //
234 this->current_states_[ip], this->prev_states_[ip],
235 this->material_states_[ip], this->output_data_[ip]);
236 loc_mat += loc_mat_current_ip;
237 }
238
239 massLumping(loc_mat);
240
241 addToLocalMatrixData(dt, local_x, local_x_prev, loc_mat, local_rhs_data,
242 local_Jac_data);
243}
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 329 of file ThermoRichardsMechanicsFEM-impl.h.

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

527{
528 auto const T = block_T(local_x);
529 auto const p_L = block_p(local_x);
530 auto const u = block_u(local_x);
531
532 auto const T_prev = block_T(local_x_prev);
533 auto const p_L_prev = block_p(local_x_prev);
534
535 auto const e_id = this->element_.getID();
536 auto const& process_data = this->process_data_;
537 auto& medium = *process_data.media_map.getMedium(e_id);
538
539 unsigned const n_integration_points =
541
542 typename ConstitutiveTraits::ConstitutiveSetting constitutive_setting;
543
544 auto models = ConstitutiveTraits::createConstitutiveModels(
545 process_data, this->solid_material_);
546 typename ConstitutiveTraits::ConstitutiveTempData tmp;
547 typename ConstitutiveTraits::ConstitutiveData CD;
548
549 for (unsigned ip = 0; ip < n_integration_points; ip++)
550 {
551 auto& current_state = this->current_states_[ip];
552 auto& output_data = this->output_data_[ip];
553
554 auto const& ip_data = ip_data_[ip];
555
556 // N is used for both p and T variables
557 auto const& N = ip_data.N_p;
558 auto const& N_u = ip_data.N_u;
559 auto const& dNdx_u = ip_data.dNdx_u;
560 auto const& dNdx = ip_data.dNdx_p;
561
562 ParameterLib::SpatialPosition const x_position{
563 std::nullopt, this->element_.getID(), ip,
565 NumLib::interpolateCoordinates<ShapeFunctionDisplacement,
567 this->element_, N_u))};
568 auto const x_coord =
569 NumLib::interpolateXCoordinate<ShapeFunctionDisplacement,
571 this->element_, N_u);
572 auto const B =
573 LinearBMatrix::computeBMatrix<DisplacementDim,
574 ShapeFunctionDisplacement::NPOINTS,
576 dNdx_u, N_u, x_coord, this->is_axially_symmetric_);
577
578 double const T_ip = N * T;
579 double const T_prev_ip = N * T_prev;
580 GlobalDimVectorType const grad_T_ip = dNdx * T;
581
582 double const p_cap_ip = -N * p_L;
583 double const p_cap_prev_ip = -N * p_L_prev;
584 GlobalDimVectorType const grad_p_cap_ip = -dNdx * p_L;
585
586 KelvinVectorType eps = B * u;
587
588 constitutive_setting.eval(models, //
589 t, dt, x_position, //
590 medium, //
591 {T_ip, T_prev_ip, grad_T_ip}, //
592 {p_cap_ip, p_cap_prev_ip, grad_p_cap_ip}, //
593 eps, current_state, this->prev_states_[ip],
594 this->material_states_[ip], tmp, output_data,
595 CD);
596 }
597
599 ShapeFunction, typename ShapeFunctionDisplacement::MeshElement,
600 DisplacementDim>(this->element_, this->is_axially_symmetric_, p_L,
601 *process_data.pressure_interpolated);
603 ShapeFunction, typename ShapeFunctionDisplacement::MeshElement,
604 DisplacementDim>(this->element_, this->is_axially_symmetric_, T,
605 *process_data.temperature_interpolated);
606}
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(), and NumLib::interpolateXCoordinate().

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

166{
167 bool constexpr is_strain_temperature_constitutive =
168 std::is_same<ConstitutiveStress_StrainTemperature::ConstitutiveTraits<
169 DisplacementDim>,
170 ConstitutiveTraits>::value;
171 if (is_strain_temperature_constitutive &&
172 this->process_data_.initial_stress.type ==
174 {
175 return;
176 }
177
178 if (!is_strain_temperature_constitutive &&
179 this->process_data_.initial_stress.type == InitialStress::Type::Total)
180 {
181 return;
182 }
183
184 double const alpha_b =
185 medium.property(MPL::PropertyType::biot_coefficient)
186 .template value<double>(variables, x_position, t, 0.0 /*dt*/);
187
188 double const bishop =
189 medium.property(MPL::PropertyType::bishops_effective_stress)
190 .template value<double>(variables, x_position, t, 0.0 /*dt*/);
191
192 ConstitutiveTraits::ConstitutiveSetting::convertInitialStressType(
193 this->current_states_[ip], this->prev_states_[ip],
194 bishop * alpha_b * p_at_ip * Invariants::identity2);
195}
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 346 of file ThermoRichardsMechanicsFEM.h.

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

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

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

359 {
360 return NumLib::localDOF<
361 ShapeFunction, ShapeFunction,
363 }
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 249 of file ThermoRichardsMechanicsFEM-impl.h.

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

◆ setInitialConditionsConcrete()

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

Reimplemented from ProcessLib::LocalAssemblerInterface.

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

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