OGS
ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim > Class Template Reference

Detailed Description

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
class ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >

Definition at line 52 of file ThermoHydroMechanicsFEM.h.

#include <ThermoHydroMechanicsFEM.h>

Inheritance diagram for ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >:
[legend]
Collaboration diagram for ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >:
[legend]

Public Types

using ShapeMatricesTypeDisplacement
 
using ShapeMatricesTypePressure
 
using GlobalDimMatrixType
 
using GlobalDimVectorType
 
using Invariants = MathLib::KelvinVector::Invariants<KelvinVectorSize>
 
using SymmetricTensor = Eigen::Matrix<double, KelvinVectorSize, 1>
 

Public Member Functions

 ThermoHydroMechanicsLocalAssembler (ThermoHydroMechanicsLocalAssembler const &)=delete
 
 ThermoHydroMechanicsLocalAssembler (ThermoHydroMechanicsLocalAssembler &&)=delete
 
 ThermoHydroMechanicsLocalAssembler (MeshLib::Element const &e, std::size_t const, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, ThermoHydroMechanicsProcessData< DisplacementDim > &process_data)
 
std::size_t setIPDataInitialConditions (std::string_view const name, double const *values, int const integration_order) override
 Returns number of read integration points.
 
void assemble (double const, double const, std::vector< double > const &, std::vector< double > const &, std::vector< double > &, std::vector< double > &, std::vector< double > &) 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 setInitialConditionsConcrete (Eigen::VectorXd const local_x, double const t, int const process_id) override
 
void preTimestepConcrete (std::vector< double > const &, double const, double const) override
 
void postTimestepConcrete (Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, double const t, double const dt, int const) 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.
 
std::vector< double > const & getIntPtDarcyVelocity (const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
 
std::vector< double > getSigma () const override
 
std::vector< double > const & getIntPtFluidDensity (const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
 
std::vector< double > const & getIntPtViscosity (const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) 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 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 Types

using BMatricesType
 
using IpData
 

Private Member Functions

ConstitutiveRelationsValues< DisplacementDim > updateConstitutiveRelations (Eigen::Ref< Eigen::VectorXd const > const local_x, Eigen::Ref< Eigen::VectorXd const > const local_x_prev, ParameterLib::SpatialPosition const &x_position, double const t, double const dt, IpData &ip_data, IntegrationPointDataForOutput< DisplacementDim > &ip_data_output) const
 
std::size_t setSigma (double const *values)
 
std::vector< double > const & getIntPtSigma (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
 
std::vector< double > const & getIntPtSigmaIce (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
 
std::vector< double > getEpsilonM () const override
 
virtual std::vector< double > const & getIntPtEpsilonM (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
 
std::vector< double > getEpsilon () const override
 
virtual std::vector< double > const & getIntPtEpsilon (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
 
std::vector< double > const & getIntPtIceVolume (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
 
unsigned getNumberOfIntegrationPoints () const override
 
int getMaterialID () const override
 
std::vector< double > getMaterialStateVariableInternalState (std::function< std::span< double >(typename MaterialLib::Solids::MechanicsBase< DisplacementDim >::MaterialStateVariables &)> const &get_values_span, int const &n_components) const override
 
MaterialLib::Solids::MechanicsBase< DisplacementDim >::MaterialStateVariables const & getMaterialStateVariablesAt (unsigned integration_point) const override
 

Static Private Member Functions

template<typename SolutionVector >
static constexpr auto localDOF (SolutionVector const &x)
 

Private Attributes

ThermoHydroMechanicsProcessData< DisplacementDim > & _process_data
 
std::vector< IpData, Eigen::aligned_allocator< IpData > > _ip_data
 
std::vector< IntegrationPointDataForOutput< DisplacementDim >, Eigen::aligned_allocator< IntegrationPointDataForOutput< DisplacementDim > > > _ip_data_output
 
NumLib::GenericIntegrationMethod const & _integration_method
 
MeshLib::Element const & _element
 
bool const _is_axially_symmetric
 
SecondaryData< typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType > _secondary_data
 

Static Private Attributes

static const int temperature_index = 0
 
static const int temperature_size = ShapeFunctionPressure::NPOINTS
 
static const int pressure_index = ShapeFunctionPressure::NPOINTS
 
static const int pressure_size = ShapeFunctionPressure::NPOINTS
 
static const int displacement_index = ShapeFunctionPressure::NPOINTS * 2
 
static const int displacement_size
 

Member Typedef Documentation

◆ BMatricesType

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
using ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::BMatricesType
private
Initial value:
BMatrixPolicyType<ShapeFunctionDisplacement, DisplacementDim>

Definition at line 264 of file ThermoHydroMechanicsFEM.h.

◆ GlobalDimMatrixType

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
using ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::GlobalDimMatrixType
Initial value:

Definition at line 63 of file ThermoHydroMechanicsFEM.h.

◆ GlobalDimVectorType

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
using ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::GlobalDimVectorType
Initial value:

Definition at line 66 of file ThermoHydroMechanicsFEM.h.

◆ Invariants

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
using ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::Invariants = MathLib::KelvinVector::Invariants<KelvinVectorSize>

Definition at line 71 of file ThermoHydroMechanicsFEM.h.

◆ IpData

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
using ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::IpData
private
Initial value:
ShapeMatricesTypePressure, DisplacementDim,
ShapeFunctionDisplacement::NPOINTS>
ShapeMatrixPolicyType< ShapeFunctionDisplacement, DisplacementDim > ShapeMatricesTypeDisplacement
BMatrixPolicyType< ShapeFunctionDisplacement, DisplacementDim > BMatricesType
ShapeMatrixPolicyType< ShapeFunctionPressure, DisplacementDim > ShapeMatricesTypePressure

Definition at line 266 of file ThermoHydroMechanicsFEM.h.

◆ ShapeMatricesTypeDisplacement

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
using ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::ShapeMatricesTypeDisplacement

◆ ShapeMatricesTypePressure

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
using ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::ShapeMatricesTypePressure

◆ SymmetricTensor

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
using ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::SymmetricTensor = Eigen::Matrix<double, KelvinVectorSize, 1>

Definition at line 73 of file ThermoHydroMechanicsFEM.h.

Constructor & Destructor Documentation

◆ ThermoHydroMechanicsLocalAssembler() [1/3]

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::ThermoHydroMechanicsLocalAssembler ( ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim > const & )
delete

◆ ThermoHydroMechanicsLocalAssembler() [2/3]

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::ThermoHydroMechanicsLocalAssembler ( ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim > && )
delete

◆ ThermoHydroMechanicsLocalAssembler() [3/3]

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::ThermoHydroMechanicsLocalAssembler ( MeshLib::Element const & e,
std::size_t const ,
NumLib::GenericIntegrationMethod const & integration_method,
bool const is_axially_symmetric,
ThermoHydroMechanicsProcessData< DisplacementDim > & process_data )

Definition at line 36 of file ThermoHydroMechanicsFEM-impl.h.

43 : _process_data(process_data),
44 _integration_method(integration_method),
45 _element(e),
46 _is_axially_symmetric(is_axially_symmetric)
47{
48 unsigned const n_integration_points =
50
51 _ip_data.reserve(n_integration_points);
52 _ip_data_output.resize(n_integration_points);
53 _secondary_data.N_u.resize(n_integration_points);
54
55 auto const shape_matrices_u =
56 NumLib::initShapeMatrices<ShapeFunctionDisplacement,
58 DisplacementDim>(e, is_axially_symmetric,
60
61 auto const shape_matrices_p =
62 NumLib::initShapeMatrices<ShapeFunctionPressure,
63 ShapeMatricesTypePressure, DisplacementDim>(
64 e, is_axially_symmetric, _integration_method);
65
66 auto const& solid_material =
68 _process_data.solid_materials, _process_data.material_ids,
69 e.getID());
70
71 // Consistency check: if frozen liquid phase is given, then the constitutive
72 // relation for ice must also be given, and vice versa.
73 auto const& medium = _process_data.media_map.getMedium(_element.getID());
74 if (medium->hasPhase("FrozenLiquid") !=
75 (_process_data.ice_constitutive_relation != nullptr))
76 {
78 "Frozen liquid phase is {:s} and the solid material constitutive "
79 "relation for ice is {:s}. But both must be given (or both "
80 "omitted).",
81 medium->hasPhase("FrozenLiquid") ? "specified" : "not specified",
82 _process_data.ice_constitutive_relation != nullptr
83 ? "specified"
84 : "not specified");
85 }
86 for (unsigned ip = 0; ip < n_integration_points; ip++)
87 {
88 _ip_data.emplace_back(solid_material);
89 auto& ip_data = _ip_data[ip];
90 auto const& sm_u = shape_matrices_u[ip];
91 ip_data.integration_weight =
93 sm_u.integralMeasure * sm_u.detJ;
94
95 ip_data.N_u = sm_u.N;
96 ip_data.dNdx_u = sm_u.dNdx;
97
98 ip_data.N = shape_matrices_p[ip].N;
99 ip_data.dNdx = shape_matrices_p[ip].dNdx;
100
101 _secondary_data.N_u[ip] = shape_matrices_u[ip].N;
102 }
103}
#define OGS_FATAL(...)
Definition Error.h:26
double getWeight() const
std::size_t getID() const
Returns the ID of the element.
Definition Element.h:89
MathLib::WeightedPoint const & getWeightedPoint(unsigned const igp) const
ThermoHydroMechanicsProcessData< DisplacementDim > & _process_data
std::vector< IntegrationPointDataForOutput< DisplacementDim >, Eigen::aligned_allocator< IntegrationPointDataForOutput< DisplacementDim > > > _ip_data_output
SecondaryData< typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType > _secondary_data
std::vector< IpData, Eigen::aligned_allocator< IpData > > _ip_data
auto & selectSolidConstitutiveRelation(SolidMaterialsMap const &constitutive_relations, MeshLib::PropertyVector< int > const *const material_ids, std::size_t const element_id)
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)
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N_u

References ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::_element, ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::_integration_method, ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::_ip_data, ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::_ip_data_output, ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::_process_data, ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::_secondary_data, MeshLib::Element::getID(), NumLib::GenericIntegrationMethod::getNumberOfPoints(), MathLib::WeightedPoint::getWeight(), NumLib::GenericIntegrationMethod::getWeightedPoint(), NumLib::initShapeMatrices(), ProcessLib::HydroMechanics::SecondaryData< ShapeMatrixType >::N_u, OGS_FATAL, and MaterialLib::Solids::selectSolidConstitutiveRelation().

Member Function Documentation

◆ assemble()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
void ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::assemble ( double const ,
double const ,
std::vector< double > const & ,
std::vector< double > const & ,
std::vector< double > & ,
std::vector< double > & ,
std::vector< double > &  )
inlineoverridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 97 of file ThermoHydroMechanicsFEM.h.

103 {
104 OGS_FATAL(
105 "ThermoHydroMechanicsLocalAssembler: assembly without Jacobian is "
106 "not implemented.");
107 }

References OGS_FATAL.

◆ assembleWithJacobian()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
void ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::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 552 of file ThermoHydroMechanicsFEM-impl.h.

560{
561 assert(local_x.size() ==
563
564 auto const x =
565 Eigen::Map<Eigen::VectorXd const>(local_x.data(), local_x.size());
566 auto const x_prev = Eigen::Map<Eigen::VectorXd const>(local_x_prev.data(),
567 local_x_prev.size());
568
569 auto const [T, p, u] = localDOF(local_x);
570 auto const [T_prev, p_prev, u_prev] = localDOF(local_x_prev);
571
572 auto local_Jac = MathLib::createZeroedMatrix<
573 typename ShapeMatricesTypeDisplacement::template MatrixType<
578
579 auto local_rhs = MathLib::createZeroedVector<
580 typename ShapeMatricesTypeDisplacement::template VectorType<
583
586
589
591 KTp.setZero(temperature_size, pressure_size);
592
594 dKTT_dp.setZero(temperature_size, pressure_size);
595
597 laplace_p.setZero(pressure_size, pressure_size);
598
600 laplace_T.setZero(pressure_size, temperature_size);
601
603 storage_p.setZero(pressure_size, pressure_size);
604
606 storage_T.setZero(pressure_size, temperature_size);
607
608 typename ShapeMatricesTypeDisplacement::template MatrixType<
610 Kup;
611 Kup.setZero(displacement_size, pressure_size);
612
613 auto const& medium = _process_data.media_map.getMedium(_element.getID());
614 bool const has_frozen_liquid_phase = medium->hasPhase("FrozenLiquid");
615
616 unsigned const n_integration_points =
618
619 std::vector<GlobalDimVectorType> ip_flux_vector;
620 double average_velocity_norm = 0.0;
621 ip_flux_vector.reserve(n_integration_points);
622
623 for (unsigned ip = 0; ip < n_integration_points; ip++)
624 {
625 auto const& N_u = _ip_data[ip].N_u;
626 ParameterLib::SpatialPosition const x_position{
627 std::nullopt, _element.getID(), ip,
629 NumLib::interpolateCoordinates<ShapeFunctionDisplacement,
631 _element, N_u))};
632
633 auto const crv = updateConstitutiveRelations(
634 x, x_prev, x_position, t, dt, _ip_data[ip], _ip_data_output[ip]);
635
636 auto const& w = _ip_data[ip].integration_weight;
637
638 auto const& dNdx_u = _ip_data[ip].dNdx_u;
639
640 auto const& N = _ip_data[ip].N;
641 auto const& dNdx = _ip_data[ip].dNdx;
642
643 auto const T_int_pt = N.dot(T);
644
645 auto const x_coord =
646 NumLib::interpolateXCoordinate<ShapeFunctionDisplacement,
648 _element, N_u);
649 auto const B =
650 LinearBMatrix::computeBMatrix<DisplacementDim,
651 ShapeFunctionDisplacement::NPOINTS,
653 dNdx_u, N_u, x_coord, _is_axially_symmetric);
654
655 auto const& b = _process_data.specific_body_force;
656 auto const velocity = _ip_data_output[ip].velocity;
657
658 //
659 // displacement equation, displacement part
660 //
661
662 if (has_frozen_liquid_phase)
663 {
664 local_Jac
665 .template block<displacement_size, displacement_size>(
667 .noalias() += B.transpose() * crv.J_uu_fr * B * w;
668
669 local_rhs.template segment<displacement_size>(displacement_index)
670 .noalias() -= B.transpose() * crv.r_u_fr * w;
671
672 local_Jac
673 .template block<displacement_size, temperature_size>(
675 .noalias() -= B.transpose() * crv.J_uT_fr * N * w;
676 }
677
678 local_Jac
679 .template block<displacement_size, displacement_size>(
681 .noalias() += B.transpose() * crv.C * B * w;
682
683 local_Jac
684 .template block<displacement_size, temperature_size>(
686 .noalias() -= B.transpose() * crv.C *
687 crv.solid_linear_thermal_expansion_coefficient * N *
688 w;
689
690 local_rhs.template segment<displacement_size>(displacement_index)
691 .noalias() -= (B.transpose() * _ip_data[ip].sigma_eff -
692 N_u_op(N_u).transpose() * crv.rho * b) *
693 w;
694
695 //
696 // displacement equation, pressure part (K_up)
697 //
698 Kup.noalias() +=
699 B.transpose() * crv.alpha_biot * Invariants::identity2 * N * w;
700
701 //
702 // pressure equation, pressure part (K_pp and M_pp).
703 //
704 laplace_p.noalias() += dNdx.transpose() * crv.K_over_mu * dNdx * w;
705
706 storage_p.noalias() +=
707 N.transpose() *
708 (_ip_data[ip].porosity * crv.fluid_compressibility +
709 (crv.alpha_biot - _ip_data[ip].porosity) * crv.beta_SR) *
710 N * w;
711
712 laplace_T.noalias() +=
713 dNdx.transpose() * crv.K_pT_thermal_osmosis * dNdx * w;
714 //
715 // RHS, pressure part
716 //
717 double const fluid_density = _ip_data_output[ip].fluid_density;
718 local_rhs.template segment<pressure_size>(pressure_index).noalias() +=
719 dNdx.transpose() * fluid_density * crv.K_over_mu * b * w;
720 //
721 // pressure equation, temperature part (M_pT)
722 //
723 storage_T.noalias() += N.transpose() * crv.beta * N * w;
724
725 //
726 // pressure equation, displacement part.
727 //
728 // Reusing Kup.transpose().
729
730 //
731 // temperature equation, temperature part.
732 //
733 KTT.noalias() +=
734 dNdx.transpose() * crv.effective_thermal_conductivity * dNdx * w;
735
736 ip_flux_vector.emplace_back(velocity * fluid_density * crv.c_f);
737 average_velocity_norm += velocity.norm();
738
739 if (has_frozen_liquid_phase)
740 {
741 local_Jac
742 .template block<temperature_size, temperature_size>(
744 .noalias() -= N.transpose() * crv.J_TT_fr * N * w;
745 }
746
747 MTT.noalias() +=
748 N.transpose() * crv.effective_volumetric_heat_capacity * N * w;
749
750 //
751 // temperature equation, pressure part
752 //
753 KTp.noalias() +=
754 dNdx.transpose() * T_int_pt * crv.K_pT_thermal_osmosis * dNdx * w;
755
756 // linearized darcy
757 dKTT_dp.noalias() -= fluid_density * crv.c_f * N.transpose() *
758 (dNdx * T).transpose() * crv.K_over_mu * dNdx * w;
759
760 /* TODO (Joerg) Temperature changes due to thermal dilatation of the
761 * fluid, which are usually discarded as being very small.
762 * Zhou et al. (10.1016/S0020-7683(98)00089-4) states that:
763 * "Biot (1956) neglected this term and it is included here for
764 * completeness"
765 * Keeping the code here in the case these are needed for the named
766 * effects in the future.
767 if (fluid_compressibility != 0)
768 {
769 auto const C_el = _ip_data[ip].computeElasticTangentStiffness(
770 t, x_position, dt, static_cast<double>(T_int_pt));
771 auto const solid_skeleton_compressibility =
772 1 / solid_material.getBulkModulus(t, x_position, &C_el);
773 double const fluid_volumetric_thermal_expansion_coefficient =
774 MaterialPropertyLib::getLiquidThermalExpansivity(
775 liquid_phase, vars, fluid_density, x_position, t, dt);
776
777 KTT.noalias() +=
778 dNdx.transpose() *
779 (-T_int_pt * fluid_volumetric_thermal_expansion_coefficient *
780 K_pT_thermal_osmosis / fluid_compressibility) *
781 dNdx * w;
782
783 local_rhs.template segment<temperature_size>(temperature_index)
784 .noalias() +=
785 dNdx.transpose() *
786 (-T_int_pt * fluid_volumetric_thermal_expansion_coefficient /
787 fluid_compressibility) *
788 fluid_density * K_over_mu * b * w;
789 MTu part for rhs and Jacobian:
790 (-T_int_pt *
791 Invariants::trace(solid_linear_thermal_expansion_coefficient) /
792 solid_skeleton_compressibility) *
793 N.transpose() * identity2.transpose() * B * w;
794 KTp part for rhs and Jacobian:
795 dNdx.transpose() *
796 (T_int_pt * fluid_volumetric_thermal_expansion_coefficient *
797 K_over_mu / fluid_compressibility) *
798 dNdx * w;
799 }
800 */
801 }
802
804 _process_data.stabilizer, _ip_data, ip_flux_vector,
805 average_velocity_norm / static_cast<double>(n_integration_points), KTT);
806
807 // temperature equation, temperature part
808 local_Jac
809 .template block<temperature_size, temperature_size>(temperature_index,
811 .noalias() += KTT + MTT / dt;
812
813 // temperature equation, pressure part
814 local_Jac
815 .template block<temperature_size, pressure_size>(temperature_index,
817 .noalias() += KTp + dKTT_dp;
818
819 // displacement equation, pressure part
820 local_Jac
821 .template block<displacement_size, pressure_size>(displacement_index,
823 .noalias() -= Kup;
824
825 // pressure equation, temperature part.
826 local_Jac
827 .template block<pressure_size, temperature_size>(pressure_index,
829 .noalias() -= storage_T / dt - laplace_T;
830
831 // pressure equation, pressure part.
832 local_Jac
833 .template block<pressure_size, pressure_size>(pressure_index,
835 .noalias() += laplace_p + storage_p / dt;
836
837 // pressure equation, displacement part.
838 local_Jac
839 .template block<pressure_size, displacement_size>(pressure_index,
841 .noalias() += Kup.transpose() / dt;
842
843 // pressure equation (f_p)
844 local_rhs.template segment<pressure_size>(pressure_index).noalias() -=
845 laplace_p * p + laplace_T * T + storage_p * (p - p_prev) / dt -
846 storage_T * (T - T_prev) / dt + Kup.transpose() * (u - u_prev) / dt;
847
848 // displacement equation (f_u)
849 local_rhs.template segment<displacement_size>(displacement_index)
850 .noalias() += Kup * p;
851
852 // temperature equation (f_T)
853 local_rhs.template segment<temperature_size>(temperature_index).noalias() -=
854 KTT * T + MTT * (T - T_prev) / dt;
855
856 local_rhs.template segment<temperature_size>(temperature_index).noalias() -=
857 KTp * p;
858}
MatrixType< _kelvin_vector_size, _number_of_dof > BMatrixType
ConstitutiveRelationsValues< DisplacementDim > updateConstitutiveRelations(Eigen::Ref< Eigen::VectorXd const > const local_x, Eigen::Ref< Eigen::VectorXd const > const local_x_prev, ParameterLib::SpatialPosition const &x_position, double const t, double const dt, IpData &ip_data, IntegrationPointDataForOutput< DisplacementDim > &ip_data_output) const
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)
void assembleAdvectionMatrix(IPData const &ip_data_vector, NumLib::ShapeMatrixCache const &shape_matrix_cache, std::vector< FluxVectorType > const &ip_flux_vector, Eigen::MatrixBase< Derived > &laplacian_matrix)
double interpolateXCoordinate(MeshLib::Element const &e, typename ShapeMatricesType::ShapeMatrices::ShapeType const &N)
std::array< double, 3 > interpolateCoordinates(MeshLib::Element const &e, typename ShapeMatricesType::ShapeMatrices::ShapeType const &N)
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.
double fluid_density(const double p, const double T, const double x)
MatrixType< ShapeFunction::NPOINTS, ShapeFunction::NPOINTS > NodalMatrixType
static Eigen::Matrix< double, KelvinVectorSize, 1 > const identity2
Kelvin mapping of 2nd order identity tensor.

References NumLib::detail::assembleAdvectionMatrix(), ProcessLib::LinearBMatrix::computeBMatrix(), MathLib::createZeroedMatrix(), MathLib::createZeroedVector(), NumLib::interpolateCoordinates(), and NumLib::interpolateXCoordinate().

◆ computeSecondaryVariableConcrete()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
void ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::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 919 of file ThermoHydroMechanicsFEM-impl.h.

923{
924 auto const p = local_x.template segment<pressure_size>(pressure_index);
925 auto const T =
926 local_x.template segment<temperature_size>(temperature_index);
927
928 unsigned const n_integration_points =
930
931 double fluid_density_avg = 0;
932 double viscosity_avg = 0;
933
935 KV sigma_avg = KV::Zero();
936
937 for (unsigned ip = 0; ip < n_integration_points; ip++)
938 {
939 auto& ip_data = _ip_data[ip];
940
941 fluid_density_avg += _ip_data_output[ip].fluid_density;
942 viscosity_avg += _ip_data_output[ip].viscosity;
943 sigma_avg += ip_data.sigma_eff;
944 }
945
946 fluid_density_avg /= n_integration_points;
947 viscosity_avg /= n_integration_points;
948 sigma_avg /= n_integration_points;
949
950 (*_process_data.element_fluid_density)[_element.getID()] =
951 fluid_density_avg;
952 (*_process_data.element_viscosity)[_element.getID()] = viscosity_avg;
953
954 Eigen::Map<KV>(&(*_process_data.element_stresses)[_element.getID() *
955 KV::RowsAtCompileTime]) =
957
959 ShapeFunctionPressure, typename ShapeFunctionDisplacement::MeshElement,
960 DisplacementDim>(_element, _is_axially_symmetric, p,
961 *_process_data.pressure_interpolated);
962
964 ShapeFunctionPressure, typename ShapeFunctionDisplacement::MeshElement,
965 DisplacementDim>(_element, _is_axially_symmetric, T,
966 *_process_data.temperature_interpolated);
967}
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
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 NumLib::interpolateToHigherOrderNodes(), and MathLib::KelvinVector::kelvinVectorToSymmetricTensor().

◆ getEpsilon()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
std::vector< double > ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::getEpsilon ( ) const
inlineoverrideprivatevirtual

Implements ProcessLib::ThermoHydroMechanics::LocalAssemblerInterface< DisplacementDim >.

Definition at line 325 of file ThermoHydroMechanicsFEM.h.

326 {
327 constexpr int kelvin_vector_size =
329
330 return transposeInPlace<kelvin_vector_size>(
331 [this](std::vector<double>& values)
332 { return getIntPtEpsilon(0, {}, {}, values); });
333 }
virtual std::vector< double > const & getIntPtEpsilon(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
constexpr int kelvin_vector_dimensions(int const displacement_dim)
Kelvin vector dimensions for given displacement dimension.

References ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::getIntPtEpsilon(), and MathLib::KelvinVector::kelvin_vector_dimensions().

◆ getEpsilonM()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
std::vector< double > ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::getEpsilonM ( ) const
inlineoverrideprivatevirtual

Implements ProcessLib::ThermoHydroMechanics::LocalAssemblerInterface< DisplacementDim >.

Definition at line 305 of file ThermoHydroMechanicsFEM.h.

306 {
307 constexpr int kelvin_vector_size =
309
310 return transposeInPlace<kelvin_vector_size>(
311 [this](std::vector<double>& values)
312 { return getIntPtEpsilonM(0, {}, {}, values); });
313 }
virtual std::vector< double > const & getIntPtEpsilonM(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override

References ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::getIntPtEpsilonM(), and MathLib::KelvinVector::kelvin_vector_dimensions().

◆ getIntPtDarcyVelocity()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
std::vector< double > const & ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::getIntPtDarcyVelocity ( const double t,
std::vector< GlobalVector * > const & x,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & dof_table,
std::vector< double > & cache ) const
overridevirtual

Implements ProcessLib::ThermoHydroMechanics::LocalAssemblerInterface< DisplacementDim >.

Definition at line 863 of file ThermoHydroMechanicsFEM-impl.h.

869{
870 unsigned const n_integration_points =
872
873 cache.clear();
874 auto cache_matrix = MathLib::createZeroedMatrix<Eigen::Matrix<
875 double, DisplacementDim, Eigen::Dynamic, Eigen::RowMajor>>(
876 cache, DisplacementDim, n_integration_points);
877
878 for (unsigned ip = 0; ip < n_integration_points; ip++)
879 {
880 cache_matrix.col(ip).noalias() = _ip_data_output[ip].velocity;
881 }
882
883 return cache;
884}

References MathLib::createZeroedMatrix().

◆ getIntPtEpsilon()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
virtual std::vector< double > const & ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::getIntPtEpsilon ( const double ,
std::vector< GlobalVector * > const & ,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & ,
std::vector< double > & cache ) const
inlineoverrideprivatevirtual

◆ getIntPtEpsilonM()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
virtual std::vector< double > const & ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::getIntPtEpsilonM ( const double ,
std::vector< GlobalVector * > const & ,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & ,
std::vector< double > & cache ) const
inlineoverrideprivatevirtual

◆ getIntPtFluidDensity()

template<typename ShapeFunctionDisplacement , typename ShapeFunction , int DisplacementDim>
std::vector< double > const & ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim >::getIntPtFluidDensity ( const double t,
std::vector< GlobalVector * > const & x,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & dof_table,
std::vector< double > & cache ) const
overridevirtual

Implements ProcessLib::ThermoHydroMechanics::LocalAssemblerInterface< DisplacementDim >.

Definition at line 889 of file ThermoHydroMechanicsFEM-impl.h.

895{
899}
std::vector< double > const & getIntegrationPointScalarData(IntegrationPointDataVector const &ip_data_vector, MemberType IpData::*const member, std::vector< double > &cache)

References ProcessLib::getIntegrationPointScalarData().

◆ getIntPtIceVolume()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
std::vector< double > const & ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::getIntPtIceVolume ( const double ,
std::vector< GlobalVector * > const & ,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & ,
std::vector< double > & cache ) const
inlineoverrideprivatevirtual

◆ getIntPtSigma()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
std::vector< double > const & ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::getIntPtSigma ( const double ,
std::vector< GlobalVector * > const & ,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & ,
std::vector< double > & cache ) const
inlineoverrideprivatevirtual

◆ getIntPtSigmaIce()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
std::vector< double > const & ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::getIntPtSigmaIce ( const double ,
std::vector< GlobalVector * > const & ,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & ,
std::vector< double > & cache ) const
inlineoverrideprivatevirtual

◆ getIntPtViscosity()

template<typename ShapeFunctionDisplacement , typename ShapeFunction , int DisplacementDim>
std::vector< double > const & ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim >::getIntPtViscosity ( const double t,
std::vector< GlobalVector * > const & x,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & dof_table,
std::vector< double > & cache ) const
overridevirtual

◆ getMaterialID()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
int ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::getMaterialID ( ) const
inlineoverrideprivatevirtual

◆ getMaterialStateVariableInternalState()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
std::vector< double > ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::getMaterialStateVariableInternalState ( std::function< std::span< double >(typename MaterialLib::Solids::MechanicsBase< DisplacementDim >::MaterialStateVariables &)> const & get_values_span,
int const & n_components ) const
inlineoverrideprivatevirtual

Implements ProcessLib::ThermoHydroMechanics::LocalAssemblerInterface< DisplacementDim >.

Definition at line 367 of file ThermoHydroMechanicsFEM.h.

372 {
375 n_components);
376 }
std::vector< double > getIntegrationPointDataMaterialStateVariables(IntegrationPointDataVector const &ip_data_vector, MemberType member, std::function< std::span< double >(MaterialStateVariables &)> get_values_span, int const n_components)
std::unique_ptr< typename MaterialLib::Solids::MechanicsBase< DisplacementDim >::MaterialStateVariables > material_state_variables

References ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::_ip_data, ProcessLib::getIntegrationPointDataMaterialStateVariables(), and ProcessLib::ThermoHydroMechanics::IntegrationPointData< BMatricesType, ShapeMatrixTypeDisplacement, ShapeMatricesTypePressure, DisplacementDim, NPoints >::material_state_variables.

◆ getMaterialStateVariablesAt()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
MaterialLib::Solids::MechanicsBase< DisplacementDim >::MaterialStateVariables const & ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::getMaterialStateVariablesAt ( unsigned integration_point) const
inlineoverrideprivatevirtual

◆ getNumberOfIntegrationPoints()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
unsigned ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::getNumberOfIntegrationPoints ( ) const
inlineoverrideprivatevirtual

◆ getShapeMatrix()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
Eigen::Map< const Eigen::RowVectorXd > ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::getShapeMatrix ( const unsigned integration_point) const
inlineoverridevirtual

Provides the shape matrix at the given integration point.

Implements NumLib::ExtrapolatableElement.

Definition at line 223 of file ThermoHydroMechanicsFEM.h.

225 {
226 auto const& N_u = _secondary_data.N_u[integration_point];
227
228 // assumes N is stored contiguously in memory
229 return Eigen::Map<const Eigen::RowVectorXd>(N_u.data(), N_u.size());
230 }

References ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::_secondary_data, and ProcessLib::HydroMechanics::SecondaryData< ShapeMatrixType >::N_u.

◆ getSigma()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
std::vector< double > ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::getSigma ( ) const
inlineoverridevirtual

Implements ProcessLib::ThermoHydroMechanics::LocalAssemblerInterface< DisplacementDim >.

Definition at line 241 of file ThermoHydroMechanicsFEM.h.

242 {
243 constexpr int kelvin_vector_size =
245
246 return transposeInPlace<kelvin_vector_size>(
247 [this](std::vector<double>& values)
248 { return getIntPtSigma(0, {}, {}, values); });
249 }
std::vector< double > const & getIntPtSigma(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override

References ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::getIntPtSigma(), and MathLib::KelvinVector::kelvin_vector_dimensions().

◆ initializeConcrete()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
void ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::initializeConcrete ( )
inlineoverridevirtual

Set initial stress from parameter.

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 117 of file ThermoHydroMechanicsFEM.h.

118 {
119 unsigned const n_integration_points =
121
122 for (unsigned ip = 0; ip < n_integration_points; ip++)
123 {
124 auto& ip_data = _ip_data[ip];
125
126 ParameterLib::SpatialPosition const x_position{
127 std::nullopt, _element.getID(), ip,
130 ShapeFunctionDisplacement,
132
134 if (_process_data.initial_stress.value)
135 {
136 ip_data.sigma_eff =
138 DisplacementDim>((*_process_data.initial_stress.value)(
139 std::numeric_limits<double>::quiet_NaN() /* time
140 independent
141 */
142 ,
143 x_position));
144 }
145
146 double const t = 0; // TODO (naumov) pass t from top
147 ip_data.solid_material.initializeInternalStateVariables(
148 t, x_position, *ip_data.material_state_variables);
149
150 ip_data.pushBackState();
151 }
152 }
Eigen::Matrix< double, Eigen::MatrixBase< Derived >::RowsAtCompileTime, 1 > symmetricTensorToKelvinVector(Eigen::MatrixBase< Derived > const &v)

References ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::_element, ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::_integration_method, ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::_ip_data, ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::_process_data, MeshLib::Element::getID(), NumLib::GenericIntegrationMethod::getNumberOfPoints(), NumLib::interpolateCoordinates(), and MathLib::KelvinVector::symmetricTensorToKelvinVector().

◆ localDOF()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
template<typename SolutionVector >
static constexpr auto ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::localDOF ( SolutionVector const & x)
inlinestaticconstexprprivate

◆ postTimestepConcrete()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
void ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::postTimestepConcrete ( Eigen::VectorXd const & local_x,
Eigen::VectorXd const & local_x_prev,
double const t,
double const dt,
int const  )
inlineoverridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 171 of file ThermoHydroMechanicsFEM.h.

175 {
176 unsigned const n_integration_points =
178
179 auto const [T_prev, p_prev, u_prev] = localDOF(local_x_prev);
180
181 for (unsigned ip = 0; ip < n_integration_points; ip++)
182 {
183 auto& ip_data = _ip_data[ip];
184 auto const& N_u = ip_data.N_u;
185 auto const& dNdx_u = ip_data.dNdx_u;
186
187 ParameterLib::SpatialPosition const x_position{
188 std::nullopt, _element.getID(), ip,
191 ShapeFunctionDisplacement,
193
194 updateConstitutiveRelations(local_x, local_x_prev, x_position, t,
195 dt, _ip_data[ip], _ip_data_output[ip]);
196
197 auto const x_coord =
198 NumLib::interpolateXCoordinate<ShapeFunctionDisplacement,
200 _element, N_u);
201 auto const B = LinearBMatrix::computeBMatrix<
202 DisplacementDim, ShapeFunctionDisplacement::NPOINTS,
203 typename BMatricesType::BMatrixType>(dNdx_u, N_u, x_coord,
205
206 ConstitutiveRelationsValues<DisplacementDim> crv;
207
209 eps_prev = B * u_prev;
210
211 _ip_data[ip].eps0 =
212 _ip_data[ip].eps0_prev +
213 (1 - _ip_data[ip].phi_fr_prev / _ip_data[ip].porosity) *
214 (eps_prev - _ip_data[ip].eps0_prev);
215 _ip_data[ip].pushBackState();
216 }
217 }

References ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::_element, ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::_integration_method, ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::_ip_data, ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::_ip_data_output, ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::_is_axially_symmetric, ProcessLib::LinearBMatrix::computeBMatrix(), MeshLib::Element::getID(), NumLib::GenericIntegrationMethod::getNumberOfPoints(), NumLib::interpolateCoordinates(), NumLib::interpolateXCoordinate(), ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::localDOF(), and ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::updateConstitutiveRelations().

◆ preTimestepConcrete()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
void ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::preTimestepConcrete ( std::vector< double > const & ,
double const ,
double const  )
inlineoverridevirtual

◆ setInitialConditionsConcrete()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
void ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::setInitialConditionsConcrete ( Eigen::VectorXd const local_x,
double const t,
int const process_id )
overridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 179 of file ThermoHydroMechanicsFEM-impl.h.

183{
184 if (!_process_data.initial_stress.isTotalStress())
185 {
186 return;
187 }
188
189 // TODO: For staggered scheme, overload
190 // LocalAssemblerInterface::setInitialConditions to enable local_x contains
191 // the primary variables from all coupled processes.
192 auto const p = local_x.template segment<pressure_size>(pressure_index);
193
194 double const dt = 0.0;
195
197 auto const& medium = _process_data.media_map.getMedium(_element.getID());
198
199 int const n_integration_points = _integration_method.getNumberOfPoints();
200 for (int ip = 0; ip < n_integration_points; ip++)
201 {
202 auto const& N = _ip_data[ip].N;
203 auto const& N_u = _ip_data[ip].N_u;
204 ParameterLib::SpatialPosition const x_position{
205 std::nullopt, _element.getID(), ip,
207 NumLib::interpolateCoordinates<ShapeFunctionDisplacement,
209 _element, N_u))};
210 auto const alpha_b =
211 medium->property(MPL::PropertyType::biot_coefficient)
212 .template value<double>(vars, x_position, t, dt);
213
214 auto& sigma_eff = _ip_data[ip].sigma_eff;
215 sigma_eff.noalias() += alpha_b * N.dot(p) * Invariants::identity2;
216 _ip_data[ip].sigma_eff_prev.noalias() = sigma_eff;
217 }
218}

References NumLib::interpolateCoordinates().

◆ setIPDataInitialConditions()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
std::size_t ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::setIPDataInitialConditions ( std::string_view const name,
double const * values,
int const integration_order )
overridevirtual

Returns number of read integration points.

Implements ProcessLib::ThermoHydroMechanics::LocalAssemblerInterface< DisplacementDim >.

Definition at line 109 of file ThermoHydroMechanicsFEM-impl.h.

112{
113 if (integration_order !=
114 static_cast<int>(_integration_method.getIntegrationOrder()))
115 {
116 OGS_FATAL(
117 "Setting integration point initial conditions; The integration "
118 "order of the local assembler for element {:d} is different from "
119 "the integration order in the initial condition.",
120 _element.getID());
121 }
122
123 if (name == "sigma")
124 {
125 if (_process_data.initial_stress.value)
126 {
127 OGS_FATAL(
128 "Setting initial conditions for stress from integration "
129 "point data and from a parameter '{:s}' is not possible "
130 "simultaneously.",
131 _process_data.initial_stress.value->name);
132 }
133
134 return ProcessLib::setIntegrationPointKelvinVectorData<DisplacementDim>(
135 values, _ip_data, &IpData::sigma_eff);
136 }
137 if (name == "epsilon_m")
138 {
139 return ProcessLib::setIntegrationPointKelvinVectorData<DisplacementDim>(
140 values, _ip_data, &IpData::eps_m);
141 }
142 if (name == "epsilon")
143 {
144 return ProcessLib::setIntegrationPointKelvinVectorData<DisplacementDim>(
145 values, _ip_data, &IpData::eps);
146 }
147 if (name.starts_with("material_state_variable_"))
148 {
149 name.remove_prefix(24);
150
151 // Using first ip data for solid material. TODO (naumov) move solid
152 // material into element, store only material state in IPs.
153 auto const& internal_variables =
154 _ip_data[0].solid_material.getInternalVariables();
155 if (auto const iv = std::find_if(
156 begin(internal_variables), end(internal_variables),
157 [&name](auto const& iv) { return iv.name == name; });
158 iv != end(internal_variables))
159 {
160 DBUG("Setting material state variable '{:s}'", name);
163 iv->reference);
164 }
165
166 int const element_id = _element.getID();
167 DBUG(
168 "The solid material of element {:d} (material ID {:d}) does not "
169 "have an internal state variable called {:s}.",
170 element_id, (*_process_data.material_ids)[element_id], name);
171 }
172
173 return 0;
174}
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
std::size_t setIntegrationPointDataMaterialStateVariables(double const *values, IntegrationPointDataVector &ip_data_vector, MemberType member, std::function< std::span< double >(MaterialStateVariables &)> get_values_span)

References DBUG(), OGS_FATAL, and ProcessLib::setIntegrationPointDataMaterialStateVariables().

◆ setSigma()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
std::size_t ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::setSigma ( double const * values)
inlineprivate

◆ updateConstitutiveRelations()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
ConstitutiveRelationsValues< DisplacementDim > ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::updateConstitutiveRelations ( Eigen::Ref< Eigen::VectorXd const > const local_x,
Eigen::Ref< Eigen::VectorXd const > const local_x_prev,
ParameterLib::SpatialPosition const & x_position,
double const t,
double const dt,
IpData & ip_data,
IntegrationPointDataForOutput< DisplacementDim > & ip_data_output ) const
private

Definition at line 223 of file ThermoHydroMechanicsFEM-impl.h.

230{
231 assert(local_x.size() ==
233
234 auto const [T, p, u] = localDOF(local_x);
235 auto const [T_prev, p_prev, u_prev] = localDOF(local_x_prev);
236
237 auto const& solid_material =
239 _process_data.solid_materials, _process_data.material_ids,
240 _element.getID());
241
242 auto const& medium = _process_data.media_map.getMedium(_element.getID());
243 auto const& liquid_phase = medium->phase("AqueousLiquid");
244 auto const& solid_phase = medium->phase("Solid");
245 auto* const frozen_liquid_phase = medium->hasPhase("FrozenLiquid")
246 ? &medium->phase("FrozenLiquid")
247 : nullptr;
249
250 auto const& N_u = ip_data.N_u;
251 auto const& dNdx_u = ip_data.dNdx_u;
252
253 auto const& N = ip_data.N;
254 auto const& dNdx = ip_data.dNdx;
255
256 auto const T_int_pt = N.dot(T);
257 auto const T_prev_int_pt = N.dot(T_prev);
258 double const dT_int_pt = T_int_pt - T_prev_int_pt;
259
260 auto const x_coord =
261 NumLib::interpolateXCoordinate<ShapeFunctionDisplacement,
263 N_u);
264 auto const B =
265 LinearBMatrix::computeBMatrix<DisplacementDim,
266 ShapeFunctionDisplacement::NPOINTS,
268 dNdx_u, N_u, x_coord, _is_axially_symmetric);
269
270 ConstitutiveRelationsValues<DisplacementDim> crv;
271
272 auto& eps = ip_data.eps;
273 eps.noalias() = B * u;
275 B * u_prev;
276
277 vars.temperature = T_int_pt;
278 double const p_int_pt = N.dot(p);
279 vars.liquid_phase_pressure = p_int_pt;
280
281 vars.liquid_saturation = 1.0;
282
283 auto const solid_density =
285 .template value<double>(vars, x_position, t, dt);
286
287 auto const porosity =
289 .template value<double>(vars, x_position, t, dt);
290 vars.porosity = porosity;
291 ip_data.porosity = porosity;
292
293 crv.alpha_biot =
295 .template value<double>(vars, x_position, t, dt);
296 auto const& alpha = crv.alpha_biot;
297
298 auto const C_el = ip_data.computeElasticTangentStiffness(
299 t, x_position, dt, static_cast<double>(T_int_pt));
300 auto const solid_skeleton_compressibility =
301 1 / solid_material.getBulkModulus(t, x_position, &C_el);
302
303 crv.beta_SR = (1 - alpha) * solid_skeleton_compressibility;
304
305 // Set mechanical variables for the intrinsic permeability model
306 // For stress dependent permeability.
307 {
308 auto const& identity2 = Invariants::identity2;
309 auto const sigma_total =
310 (ip_data.sigma_eff - alpha * p_int_pt * identity2).eval();
311 vars.total_stress.emplace<SymmetricTensor>(
313 }
314 // For strain dependent permeability
315 vars.volumetric_strain = Invariants::trace(ip_data.eps);
317 ip_data.material_state_variables->getEquivalentPlasticStrain();
318
319 auto const intrinsic_permeability =
320 MaterialPropertyLib::formEigenTensor<DisplacementDim>(
322 .value(vars, x_position, t, dt));
323
324 auto const fluid_density =
325 liquid_phase.property(MaterialPropertyLib::PropertyType::density)
326 .template value<double>(vars, x_position, t, dt);
327 ip_data_output.fluid_density = fluid_density;
328 vars.density = fluid_density;
329
330 auto const drho_dp =
331 liquid_phase.property(MaterialPropertyLib::PropertyType::density)
332 .template dValue<double>(
334 x_position, t, dt);
335
336 crv.fluid_compressibility = 1 / fluid_density * drho_dp;
337
338 double const fluid_volumetric_thermal_expansion_coefficient =
340 liquid_phase, vars, fluid_density, x_position, t, dt);
341
342 // Use the viscosity model to compute the viscosity
343 ip_data_output.viscosity =
345 .template value<double>(vars, x_position, t, dt);
346 crv.K_over_mu = intrinsic_permeability / ip_data_output.viscosity;
347
348 auto const& b = _process_data.specific_body_force;
349
350 // Consider also anisotropic thermal expansion.
351 crv.solid_linear_thermal_expansion_coefficient =
352 MPL::formKelvinVector<DisplacementDim>(
353 solid_phase
354 .property(
356 .value(vars, x_position, t, dt));
357
359 dthermal_strain =
360 crv.solid_linear_thermal_expansion_coefficient * dT_int_pt;
361
362 crv.K_pT_thermal_osmosis =
363 (solid_phase.hasProperty(
365 ? MaterialPropertyLib::formEigenTensor<DisplacementDim>(
366 solid_phase
369 .value(vars, x_position, t, dt))
370 : Eigen::MatrixXd::Zero(DisplacementDim, DisplacementDim));
371
372 GlobalDimVectorType const velocity = -crv.K_over_mu * dNdx * p -
373 crv.K_pT_thermal_osmosis * dNdx * T +
374 crv.K_over_mu * fluid_density * b;
375 ip_data_output.velocity = velocity;
376
377 //
378 // displacement equation, displacement part
379 //
380 auto& eps_m = ip_data.eps_m;
381 auto& eps_m_prev = ip_data.eps_m_prev;
382 eps_m.noalias() = eps_m_prev + eps - eps_prev - dthermal_strain;
385 eps_m);
386
387 crv.C = ip_data.updateConstitutiveRelation(vars, t, x_position, dt,
388 T_prev_int_pt);
389
390 crv.rho = solid_density * (1 - porosity) + porosity * fluid_density;
391
392 crv.beta =
393 porosity * fluid_volumetric_thermal_expansion_coefficient +
394 (alpha - porosity) *
395 Invariants::trace(crv.solid_linear_thermal_expansion_coefficient);
396
397 //
398 // pressure equation, displacement part.
399 //
400 // Reusing Kup.transpose().
401
402 //
403 // temperature equation, temperature part.
404 //
405 crv.c_f =
406 liquid_phase
408 .template value<double>(vars, x_position, t, dt);
409 crv.effective_thermal_conductivity =
410 MaterialPropertyLib::formEigenTensor<DisplacementDim>(
411 medium
412 ->property(
414 .value(vars, x_position, t, dt));
415
416 // Thermal conductivity is moved outside and zero matrix is passed instead
417 // due to multiplication with fluid's density times specific heat capacity.
418 crv.effective_thermal_conductivity.noalias() +=
419 fluid_density * crv.c_f *
421 _process_data.stabilizer, _element.getID(),
422 GlobalDimMatrixType::Zero(DisplacementDim, DisplacementDim),
423 velocity, 0. /* phi */, 0. /* dispersivity_transversal */,
424 0. /*dispersivity_longitudinal*/);
425
426 double const c_s =
427 solid_phase
429 .template value<double>(vars, x_position, t, dt);
430
431 // Also modified by freezing terms.
432 crv.effective_volumetric_heat_capacity =
433 porosity * fluid_density * crv.c_f +
434 (1.0 - porosity) * solid_density * c_s;
435
436 if (frozen_liquid_phase)
437 {
439 double const phi_fr =
441 .template value<double>(vars, x_position, t, dt);
442 ip_data.phi_fr = phi_fr;
443
444 auto const frozen_liquid_value =
446 {
447 return (*frozen_liquid_phase)[p].template value<double>(
448 vars, x_position, t, dt);
449 };
450
451 double const rho_fr =
453
454 double const c_fr = frozen_liquid_value(
456
457 double const l_fr = frozen_liquid_value(
459
460 double const dphi_fr_dT =
462 .template dValue<double>(
464 x_position, t, dt);
465
466 double const phi_fr_prev = [&]()
467 {
469 vars_prev.temperature = T_prev_int_pt;
471 .template value<double>(vars_prev, x_position, t, dt);
472 }();
473 ip_data.phi_fr_prev = phi_fr_prev;
474
475 // alpha_T^I
477 DisplacementDim> const ice_linear_thermal_expansion_coefficient =
478 MPL::formKelvinVector<DisplacementDim>(
479 frozen_liquid_phase
480 ->property(
482 .value(vars, x_position, t, dt));
483
485 dthermal_strain_ice =
486 ice_linear_thermal_expansion_coefficient * dT_int_pt;
487
488 // alpha_{phi_I} -- linear expansion coeff. due to water-to-ice
489 // transition (phase change), and related phase_change_strain term
491 phase_change_expansion_coefficient =
492 MPL::formKelvinVector<DisplacementDim>(
493 frozen_liquid_phase
496 .value(vars, x_position, t, dt));
497
499 dphase_change_strain = phase_change_expansion_coefficient *
500 (phi_fr - phi_fr_prev) / porosity;
501
502 // eps0 ia a 'history variable' -- a solid matrix strain accrued
503 // prior to the onset of ice forming
504 auto& eps0 = ip_data.eps0;
505 auto const& eps0_prev = ip_data.eps0_prev;
506
507 // definition of eps_m_ice
508 auto& eps_m_ice = ip_data.eps_m_ice;
509 auto const& eps_m_ice_prev = ip_data.eps_m_ice_prev;
510
511 eps_m_ice.noalias() = eps_m_ice_prev + eps - eps_prev -
512 (eps0 - eps0_prev) - dthermal_strain_ice -
513 dphase_change_strain;
514
515 vars_ice.mechanical_strain
517 eps_m_ice);
518 auto const C_IR = ip_data.updateConstitutiveRelationIce(
519 *_process_data.ice_constitutive_relation, vars_ice, t, x_position,
520 dt, T_prev_int_pt);
521 crv.effective_volumetric_heat_capacity +=
522 -phi_fr * fluid_density * crv.c_f + phi_fr * rho_fr * c_fr -
523 l_fr * rho_fr * dphi_fr_dT;
524
525 // part of dMTT_dT derivative for freezing
526 double const d2phi_fr_dT2 =
528 .template d2Value<double>(
531 dt);
532
533 crv.J_uu_fr = phi_fr * C_IR;
534
535 auto const& sigma_eff_ice = ip_data.sigma_eff_ice;
536 crv.r_u_fr = phi_fr * sigma_eff_ice;
537
538 crv.J_uT_fr = phi_fr * C_IR * ice_linear_thermal_expansion_coefficient;
539
540 crv.J_TT_fr = ((rho_fr * c_fr - fluid_density * crv.c_f) * dphi_fr_dT +
541 l_fr * rho_fr * d2phi_fr_dT2) *
542 dT_int_pt / dt;
543 }
544 return crv;
545}
std::variant< std::monostate, Eigen::Matrix< double, 4, 1 >, Eigen::Matrix< double, 6, 1 > > mechanical_strain
std::variant< std::monostate, Eigen::Matrix< double, 4, 1 >, Eigen::Matrix< double, 6, 1 > > total_stress
typename ShapeMatricesTypePressure::GlobalDimVectorType GlobalDimVectorType
double getLiquidThermalExpansivity(Phase const &phase, VariableArray const &vars, const double density, ParameterLib::SpatialPosition const &pos, double const t, double const dt)
Eigen::MatrixXd computeHydrodynamicDispersion(NumericalStabilization const &stabilizer, std::size_t const element_id, Eigen::MatrixXd const &pore_diffusion_coefficient, Eigen::VectorXd const &velocity, double const porosity, double const solute_dispersivity_transverse, double const solute_dispersivity_longitudinal)
auto eval(Function &f, Tuples &... ts) -> typename detail::GetFunctionReturnType< decltype(&Function::eval)>::type
Definition Apply.h:267
static double trace(Eigen::Matrix< double, KelvinVectorSize, 1 > const &v)
Trace of the corresponding tensor.

References MaterialPropertyLib::biot_coefficient, ProcessLib::LinearBMatrix::computeBMatrix(), ProcessLib::ThermoHydroMechanics::IntegrationPointData< BMatricesType, ShapeMatrixTypeDisplacement, ShapeMatricesTypePressure, DisplacementDim, NPoints >::computeElasticTangentStiffness(), NumLib::computeHydrodynamicDispersion(), MaterialPropertyLib::density, MaterialPropertyLib::VariableArray::density, ProcessLib::ThermoHydroMechanics::IntegrationPointData< BMatricesType, ShapeMatrixTypeDisplacement, ShapeMatricesTypePressure, DisplacementDim, NPoints >::dNdx, ProcessLib::ThermoHydroMechanics::IntegrationPointData< BMatricesType, ShapeMatrixTypeDisplacement, ShapeMatricesTypePressure, DisplacementDim, NPoints >::dNdx_u, ProcessLib::ThermoHydroMechanics::IntegrationPointData< BMatricesType, ShapeMatrixTypeDisplacement, ShapeMatricesTypePressure, DisplacementDim, NPoints >::eps, ProcessLib::ThermoHydroMechanics::IntegrationPointData< BMatricesType, ShapeMatrixTypeDisplacement, ShapeMatricesTypePressure, DisplacementDim, NPoints >::eps0, ProcessLib::ThermoHydroMechanics::IntegrationPointData< BMatricesType, ShapeMatrixTypeDisplacement, ShapeMatricesTypePressure, DisplacementDim, NPoints >::eps0_prev, ProcessLib::ThermoHydroMechanics::IntegrationPointData< BMatricesType, ShapeMatrixTypeDisplacement, ShapeMatricesTypePressure, DisplacementDim, NPoints >::eps_m, ProcessLib::ThermoHydroMechanics::IntegrationPointData< BMatricesType, ShapeMatrixTypeDisplacement, ShapeMatricesTypePressure, DisplacementDim, NPoints >::eps_m_ice, ProcessLib::ThermoHydroMechanics::IntegrationPointData< BMatricesType, ShapeMatrixTypeDisplacement, ShapeMatricesTypePressure, DisplacementDim, NPoints >::eps_m_ice_prev, ProcessLib::ThermoHydroMechanics::IntegrationPointData< BMatricesType, ShapeMatrixTypeDisplacement, ShapeMatricesTypePressure, DisplacementDim, NPoints >::eps_m_prev, MaterialPropertyLib::VariableArray::equivalent_plastic_strain, ProcessLib::ThermoHydroMechanics::IntegrationPointDataForOutput< DisplacementDim >::fluid_density, MaterialPropertyLib::getLiquidThermalExpansivity(), NumLib::interpolateXCoordinate(), MathLib::KelvinVector::kelvinVectorToSymmetricTensor(), MaterialPropertyLib::liquid_phase_pressure, MaterialPropertyLib::VariableArray::liquid_phase_pressure, MaterialPropertyLib::VariableArray::liquid_saturation, ProcessLib::ThermoHydroMechanics::IntegrationPointData< BMatricesType, ShapeMatrixTypeDisplacement, ShapeMatricesTypePressure, DisplacementDim, NPoints >::material_state_variables, MaterialPropertyLib::VariableArray::mechanical_strain, ProcessLib::ThermoHydroMechanics::IntegrationPointData< BMatricesType, ShapeMatrixTypeDisplacement, ShapeMatricesTypePressure, DisplacementDim, NPoints >::N, ProcessLib::ThermoHydroMechanics::IntegrationPointData< BMatricesType, ShapeMatrixTypeDisplacement, ShapeMatricesTypePressure, DisplacementDim, NPoints >::N_u, MaterialPropertyLib::permeability, ProcessLib::ThermoHydroMechanics::IntegrationPointData< BMatricesType, ShapeMatrixTypeDisplacement, ShapeMatricesTypePressure, DisplacementDim, NPoints >::phi_fr, ProcessLib::ThermoHydroMechanics::IntegrationPointData< BMatricesType, ShapeMatrixTypeDisplacement, ShapeMatricesTypePressure, DisplacementDim, NPoints >::phi_fr_prev, MaterialPropertyLib::porosity, MaterialPropertyLib::VariableArray::porosity, ProcessLib::ThermoHydroMechanics::IntegrationPointData< BMatricesType, ShapeMatrixTypeDisplacement, ShapeMatricesTypePressure, DisplacementDim, NPoints >::porosity, MaterialLib::Solids::selectSolidConstitutiveRelation(), ProcessLib::ThermoHydroMechanics::IntegrationPointData< BMatricesType, ShapeMatrixTypeDisplacement, ShapeMatricesTypePressure, DisplacementDim, NPoints >::sigma_eff, ProcessLib::ThermoHydroMechanics::IntegrationPointData< BMatricesType, ShapeMatrixTypeDisplacement, ShapeMatricesTypePressure, DisplacementDim, NPoints >::sigma_eff_ice, MaterialPropertyLib::specific_heat_capacity, MaterialPropertyLib::specific_latent_heat, MaterialPropertyLib::temperature, MaterialPropertyLib::VariableArray::temperature, MaterialPropertyLib::thermal_conductivity, MaterialPropertyLib::thermal_expansivity, MaterialPropertyLib::thermal_osmosis_coefficient, MaterialPropertyLib::VariableArray::total_stress, ProcessLib::ThermoHydroMechanics::IntegrationPointData< BMatricesType, ShapeMatrixTypeDisplacement, ShapeMatricesTypePressure, DisplacementDim, NPoints >::updateConstitutiveRelation(), ProcessLib::ThermoHydroMechanics::IntegrationPointData< BMatricesType, ShapeMatrixTypeDisplacement, ShapeMatricesTypePressure, DisplacementDim, NPoints >::updateConstitutiveRelationIce(), ProcessLib::ThermoHydroMechanics::IntegrationPointDataForOutput< DisplacementDim >::velocity, MaterialPropertyLib::viscosity, ProcessLib::ThermoHydroMechanics::IntegrationPointDataForOutput< DisplacementDim >::viscosity, MaterialPropertyLib::volume_fraction, and MaterialPropertyLib::VariableArray::volumetric_strain.

Referenced by ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::postTimestepConcrete().

Member Data Documentation

◆ _element

◆ _integration_method

◆ _ip_data

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
std::vector<IpData, Eigen::aligned_allocator<IpData> > ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::_ip_data
private

Definition at line 396 of file ThermoHydroMechanicsFEM.h.

Referenced by ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::ThermoHydroMechanicsLocalAssembler(), ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::getIntPtEpsilon(), ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::getIntPtEpsilonM(), ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::getIntPtIceVolume(), ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::getIntPtSigma(), ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::getIntPtSigmaIce(), ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::getMaterialStateVariableInternalState(), ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::getMaterialStateVariablesAt(), ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::initializeConcrete(), ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::postTimestepConcrete(), and ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::setSigma().

◆ _ip_data_output

◆ _is_axially_symmetric

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
bool const ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::_is_axially_symmetric
private

◆ _process_data

◆ _secondary_data

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
SecondaryData< typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType> ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::_secondary_data
private

◆ displacement_index

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
const int ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::displacement_index = ShapeFunctionPressure::NPOINTS * 2
staticprivate

Definition at line 415 of file ThermoHydroMechanicsFEM.h.

◆ displacement_size

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
const int ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::displacement_size
staticprivate
Initial value:
=
ShapeFunctionDisplacement::NPOINTS * DisplacementDim

Definition at line 416 of file ThermoHydroMechanicsFEM.h.

◆ KelvinVectorSize

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
int const ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::KelvinVectorSize
static
Initial value:

Definition at line 69 of file ThermoHydroMechanicsFEM.h.

◆ N_u_op

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
constexpr auto& ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::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 75 of file ThermoHydroMechanicsFEM.h.

◆ pressure_index

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
const int ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::pressure_index = ShapeFunctionPressure::NPOINTS
staticprivate

Definition at line 413 of file ThermoHydroMechanicsFEM.h.

◆ pressure_size

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
const int ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::pressure_size = ShapeFunctionPressure::NPOINTS
staticprivate

Definition at line 414 of file ThermoHydroMechanicsFEM.h.

◆ temperature_index

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
const int ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::temperature_index = 0
staticprivate

Definition at line 411 of file ThermoHydroMechanicsFEM.h.

◆ temperature_size

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
const int ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::temperature_size = ShapeFunctionPressure::NPOINTS
staticprivate

Definition at line 412 of file ThermoHydroMechanicsFEM.h.


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