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 > &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_b_data, std::vector< double > &local_Jac_data)
virtual void computeSecondaryVariable (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_prev, int const process_id)
virtual void preTimestep (std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, double const t, double const delta_t)
virtual void postTimestep (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, int const process_id)
void postNonLinearSolver (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, int const process_id)
virtual Eigen::Vector3d getFlux (MathLib::Point3d const &, double const, std::vector< double > const &) const
virtual Eigen::Vector3d getFlux (MathLib::Point3d const &, double const, std::vector< std::vector< double > > const &) const
 Fits to staggered scheme.
virtual std::optional< VectorSegmentgetVectorDeformationSegment () const
Public Member Functions inherited from NumLib::ExtrapolatableElement
virtual ~ExtrapolatableElement ()=default

Static Public Attributes

static int const KelvinVectorSize
static constexpr auto & N_u_op

Private 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 > getEpsilon0 () const override
virtual std::vector< double > const & getIntPtEpsilon0 (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

◆ 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
ShapeMatrixPolicyType< ShapeFunctionPressure, DisplacementDim > ShapeMatricesTypePressure
BMatrixPolicyType< ShapeFunctionDisplacement, DisplacementDim > BMatricesType

Definition at line 264 of file ThermoHydroMechanicsFEM.h.

◆ ShapeMatricesTypeDisplacement

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
using ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::ShapeMatricesTypeDisplacement
Initial value:
EigenFixedShapeMatrixPolicy< ShapeFunction, GlobalDim > ShapeMatrixPolicyType

Definition at line 56 of file ThermoHydroMechanicsFEM.h.

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

45 _element(e),
47{
48 unsigned const n_integration_points =
49 _integration_method.getNumberOfPoints();
50
54
55 auto const shape_matrices_u =
60
61 auto const shape_matrices_p =
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 =
92 _integration_method.getWeightedPoint(ip).getWeight() *
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
99 ip_data.dNdx = shape_matrices_p[ip].dNdx;
100
102 }
103}
#define OGS_FATAL(...)
Definition Error.h:26
ShapeMatrixPolicyType< ShapeFunctionDisplacement, DisplacementDim > ShapeMatricesTypeDisplacement
ThermoHydroMechanicsProcessData< DisplacementDim > & _process_data
std::vector< IntegrationPointDataForOutput< DisplacementDim >, Eigen::aligned_allocator< IntegrationPointDataForOutput< DisplacementDim > > > _ip_data_output
SecondaryData< typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType > _secondary_data
ShapeMatrixPolicyType< ShapeFunctionPressure, DisplacementDim > ShapeMatricesTypePressure
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)

References _element, _integration_method, _ip_data, _ip_data_output, _is_axially_symmetric, _process_data, _secondary_data, MeshLib::Element::getID(), NumLib::initShapeMatrices(), 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 > & local_rhs_data,
std::vector< double > & local_Jac_data )
overridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

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

569{
570 assert(local_x.size() ==
572
573 auto const x =
576 local_x_prev.size());
577
578 auto const [T, p, u] = localDOF(local_x);
579 auto const [T_prev, p_prev, u_prev] = localDOF(local_x_prev);
580
587
592
595
598
601
604
607
610
613
616
619 Kup;
621
622 auto const& medium = _process_data.media_map.getMedium(_element.getID());
623 bool const has_frozen_liquid_phase = medium->hasPhase("FrozenLiquid");
624
625 unsigned const n_integration_points =
626 _integration_method.getNumberOfPoints();
627
629 double average_velocity_norm = 0.0;
631
632 for (unsigned ip = 0; ip < n_integration_points; ip++)
633 {
634 auto const& N_u = _ip_data[ip].N_u;
636 std::nullopt, _element.getID(),
640 _element, N_u))};
641
642 auto const crv = updateConstitutiveRelations(
644
645 auto const& w = _ip_data[ip].integration_weight;
646
647 auto const& dNdx_u = _ip_data[ip].dNdx_u;
648
649 auto const& N = _ip_data[ip].N;
650 auto const& dNdx = _ip_data[ip].dNdx;
651
652 auto const T_int_pt = N.dot(T);
653
654 auto const x_coord =
657 _element, N_u);
658 auto const B =
663
664 auto const& b = _process_data.specific_body_force;
665 auto const velocity = _ip_data_output[ip].velocity;
666
667 //
668 // displacement equation, displacement part
669 //
670
672 {
676 .noalias() += B.transpose() * crv.J_uu_fr * B * w;
677
679 .noalias() -= B.transpose() * crv.r_u_fr * w;
680
684 .noalias() -= B.transpose() * crv.J_uT_fr * N * w;
685 }
686
690 .noalias() += B.transpose() * crv.C * B * w;
691
695 .noalias() -= B.transpose() * crv.C *
696 crv.solid_linear_thermal_expansion_coefficient * N *
697 w;
698
700 .noalias() -= (B.transpose() * _ip_data[ip].sigma_eff -
701 N_u_op(N_u).transpose() * crv.rho * b) *
702 w;
703
704 //
705 // displacement equation, pressure part (K_up)
706 //
707 Kup.noalias() +=
708 B.transpose() * crv.alpha_biot * Invariants::identity2 * N * w;
709
710 //
711 // pressure equation, pressure part (K_pp and M_pp).
712 //
713 laplace_p.noalias() += dNdx.transpose() * crv.K_over_mu * dNdx * w;
714
715 storage_p.noalias() +=
716 N.transpose() *
717 (_ip_data[ip].porosity * crv.fluid_compressibility +
718 (crv.alpha_biot - _ip_data[ip].porosity) * crv.beta_SR) *
719 N * w;
720
721 laplace_T.noalias() +=
722 dNdx.transpose() * crv.K_pT_thermal_osmosis * dNdx * w;
723 //
724 // RHS, pressure part
725 //
726 double const fluid_density = _ip_data_output[ip].fluid_density;
727 local_rhs.template segment<pressure_size>(pressure_index).noalias() +=
728 dNdx.transpose() * fluid_density * crv.K_over_mu * b * w;
729 //
730 // pressure equation, temperature part (M_pT)
731 //
732 storage_T.noalias() += N.transpose() * crv.beta * N * w;
733
734 //
735 // pressure equation, displacement part.
736 //
737 // Reusing Kup.transpose().
738
739 //
740 // temperature equation, temperature part.
741 //
742 KTT.noalias() +=
743 dNdx.transpose() * crv.effective_thermal_conductivity * dNdx * w;
744
745 ip_flux_vector.emplace_back(velocity * fluid_density * crv.c_f);
747
749 {
753 .noalias() -= N.transpose() * crv.J_TT_fr * N * w;
754 }
755
756 MTT.noalias() +=
757 N.transpose() * crv.effective_volumetric_heat_capacity * N * w;
758
759 //
760 // temperature equation, pressure part
761 //
762 KTp.noalias() +=
763 dNdx.transpose() * T_int_pt * crv.K_pT_thermal_osmosis * dNdx * w;
764
765 // linearized darcy
766 dKTT_dp.noalias() -= fluid_density * crv.c_f * N.transpose() *
767 (dNdx * T).transpose() * crv.K_over_mu * dNdx * w;
768
769 /* TODO (Joerg) Temperature changes due to thermal dilatation of the
770 * fluid, which are usually discarded as being very small.
771 * Zhou et al. (10.1016/S0020-7683(98)00089-4) states that:
772 * "Biot (1956) neglected this term and it is included here for
773 * completeness"
774 * Keeping the code here in the case these are needed for the named
775 * effects in the future.
776 if (fluid_compressibility != 0)
777 {
778 auto const C_el = _ip_data[ip].computeElasticTangentStiffness(
779 t, x_position, dt, static_cast<double>(T_int_pt));
780 auto const solid_skeleton_compressibility =
781 1 / solid_material.getBulkModulus(t, x_position, &C_el);
782 double const fluid_volumetric_thermal_expansion_coefficient =
783 MaterialPropertyLib::getLiquidThermalExpansivity(
784 liquid_phase, vars, fluid_density, x_position, t, dt);
785
786 KTT.noalias() +=
787 dNdx.transpose() *
788 (-T_int_pt * fluid_volumetric_thermal_expansion_coefficient *
789 K_pT_thermal_osmosis / fluid_compressibility) *
790 dNdx * w;
791
792 local_rhs.template segment<temperature_size>(temperature_index)
793 .noalias() +=
794 dNdx.transpose() *
795 (-T_int_pt * fluid_volumetric_thermal_expansion_coefficient /
796 fluid_compressibility) *
797 fluid_density * K_over_mu * b * w;
798 MTu part for rhs and Jacobian:
799 (-T_int_pt *
800 Invariants::trace(solid_linear_thermal_expansion_coefficient) /
801 solid_skeleton_compressibility) *
802 N.transpose() * identity2.transpose() * B * w;
803 KTp part for rhs and Jacobian:
804 dNdx.transpose() *
805 (T_int_pt * fluid_volumetric_thermal_expansion_coefficient *
806 K_over_mu / fluid_compressibility) *
807 dNdx * w;
808 }
809 */
810 }
811
814 average_velocity_norm / static_cast<double>(n_integration_points), KTT);
815
816 // temperature equation, temperature part
820 .noalias() += KTT + MTT / dt;
821
822 // temperature equation, pressure part
826 .noalias() += KTp + dKTT_dp;
827
828 // displacement equation, pressure part
832 .noalias() -= Kup;
833
834 // pressure equation, temperature part.
838 .noalias() += -storage_T / dt + laplace_T;
839
840 // pressure equation, pressure part.
844 .noalias() += laplace_p + storage_p / dt;
845
846 // pressure equation, displacement part.
850 .noalias() += Kup.transpose() / dt;
851
852 // pressure equation (f_p)
853 local_rhs.template segment<pressure_size>(pressure_index).noalias() -=
854 laplace_p * p + laplace_T * T + storage_p * (p - p_prev) / dt -
855 storage_T * (T - T_prev) / dt + Kup.transpose() * (u - u_prev) / dt;
856
857 // displacement equation (f_u)
859 .noalias() += Kup * p;
860
861 // temperature equation (f_T)
863 KTT * T + MTT * (T - T_prev) / dt;
864
866 KTp * p;
867}
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
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)
MatrixType< ShapeFunction::NPOINTS, ShapeFunction::NPOINTS > NodalMatrixType
static Eigen::Matrix< double, KelvinVectorSize, 1 > const identity2
Kelvin mapping of 2nd order identity tensor.

References _element, _integration_method, _ip_data, _ip_data_output, _is_axially_symmetric, _process_data, NumLib::detail::assembleAdvectionMatrix(), assembleWithJacobian(), ProcessLib::LinearBMatrix::computeBMatrix(), MathLib::createZeroedMatrix(), MathLib::createZeroedVector(), displacement_index, displacement_size, MathLib::KelvinVector::Invariants< KelvinVectorSize >::identity2, NumLib::interpolateCoordinates(), NumLib::interpolateXCoordinate(), localDOF(), N_u_op, pressure_index, pressure_size, temperature_index, temperature_size, and updateConstitutiveRelations().

Referenced by assembleWithJacobian().

◆ 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 928 of file ThermoHydroMechanicsFEM-impl.h.

932{
933 auto const p = local_x.template segment<pressure_size>(pressure_index);
934 auto const T =
936
937 unsigned const n_integration_points =
938 _integration_method.getNumberOfPoints();
939
940 double phi_fr_avg = 0;
941 double fluid_density_avg = 0;
942 double viscosity_avg = 0;
943
946
947 for (unsigned ip = 0; ip < n_integration_points; ip++)
948 {
949 auto& ip_data = _ip_data[ip];
950
951 phi_fr_avg += _ip_data[ip].phi_fr;
952 fluid_density_avg += _ip_data_output[ip].fluid_density;
953 viscosity_avg += _ip_data_output[ip].viscosity;
954 sigma_avg += ip_data.sigma_eff;
955 }
956
961
962 (*_process_data.element_phi_fr)[_element.getID()] = phi_fr_avg;
963 (*_process_data.element_fluid_density)[_element.getID()] =
965 (*_process_data.element_viscosity)[_element.getID()] = viscosity_avg;
966
967 Eigen::Map<KV>(&(*_process_data.element_stresses)[_element.getID() *
970
974 *_process_data.pressure_interpolated);
975
979 *_process_data.temperature_interpolated);
980}
Eigen::Matrix< double, 4, 1 > kelvinVectorToSymmetricTensor(Eigen::Matrix< double, 4, 1, Eigen::ColMajor, 4, 1 > const &v)

References _element, _integration_method, _ip_data, _ip_data_output, _is_axially_symmetric, _process_data, NumLib::interpolateToHigherOrderNodes(), MathLib::KelvinVector::kelvinVectorToSymmetricTensor(), pressure_index, and temperature_index.

◆ 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 342 of file ThermoHydroMechanicsFEM.h.

343 {
344 constexpr int kelvin_vector_size =
346
349 { return getIntPtEpsilon(0, {}, {}, values); });
350 }
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.
std::vector< double > transposeInPlace(StoreValuesFunction const &store_values_function)

References getIntPtEpsilon(), MathLib::KelvinVector::kelvin_vector_dimensions(), and ProcessLib::transposeInPlace().

◆ getEpsilon0()

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

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

Definition at line 303 of file ThermoHydroMechanicsFEM.h.

304 {
305 constexpr int kelvin_vector_size =
307
310 { return getIntPtEpsilon0(0, {}, {}, values); });
311 }
virtual std::vector< double > const & getIntPtEpsilon0(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override

References getIntPtEpsilon0(), MathLib::KelvinVector::kelvin_vector_dimensions(), and ProcessLib::transposeInPlace().

◆ 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 322 of file ThermoHydroMechanicsFEM.h.

323 {
324 constexpr int kelvin_vector_size =
326
329 { return getIntPtEpsilonM(0, {}, {}, values); });
330 }
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 getIntPtEpsilonM(), MathLib::KelvinVector::kelvin_vector_dimensions(), and ProcessLib::transposeInPlace().

◆ 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 872 of file ThermoHydroMechanicsFEM-impl.h.

878{
879 unsigned const n_integration_points =
880 _integration_method.getNumberOfPoints();
881
882 cache.clear();
886
887 for (unsigned ip = 0; ip < n_integration_points; ip++)
888 {
889 cache_matrix.col(ip).noalias() = _ip_data_output[ip].velocity;
890 }
891
892 return cache;
893}

References _integration_method, _ip_data_output, and 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

◆ getIntPtEpsilon0()

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
virtual std::vector< double > const & ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::getIntPtEpsilon0 ( 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

◆ getIntPtIceVolume()

◆ 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()

◆ 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

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

Definition at line 377 of file ThermoHydroMechanicsFEM.h.

378 {
379 return _process_data.material_ids == nullptr
380 ? 0
381 : (*_process_data.material_ids)[_element.getID()];
382 }

References _element, and _process_data.

◆ 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

◆ 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

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

Definition at line 397 of file ThermoHydroMechanicsFEM.h.

398 {
399 return *_ip_data[integration_point].material_state_variables;
400 }

References _ip_data.

◆ getNumberOfIntegrationPoints()

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

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

Definition at line 372 of file ThermoHydroMechanicsFEM.h.

373 {
374 return _integration_method.getNumberOfPoints();
375 }

References _integration_method.

◆ 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 221 of file ThermoHydroMechanicsFEM.h.

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

References _secondary_data.

◆ 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 239 of file ThermoHydroMechanicsFEM.h.

240 {
241 constexpr int kelvin_vector_size =
243
246 { return getIntPtSigma(0, {}, {}, values); });
247 }
std::vector< double > const & getIntPtSigma(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override

References getIntPtSigma(), MathLib::KelvinVector::kelvin_vector_dimensions(), and ProcessLib::transposeInPlace().

◆ 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 115 of file ThermoHydroMechanicsFEM.h.

116 {
117 unsigned const n_integration_points =
118 _integration_method.getNumberOfPoints();
119
120 for (unsigned ip = 0; ip < n_integration_points; ip++)
121 {
122 auto& ip_data = _ip_data[ip];
123
125 std::nullopt, _element.getID(),
130
132 if (_process_data.initial_stress.value)
133 {
134 ip_data.sigma_eff =
136 DisplacementDim>((*_process_data.initial_stress.value)(
138 independent
139 */
140 ,
141 x_position));
142 }
143
144 double const t = 0; // TODO (naumov) pass t from top
145 ip_data.solid_material.initializeInternalStateVariables(
146 t, x_position, *ip_data.material_state_variables);
147
148 ip_data.pushBackState();
149 }
150 }

References _element, _integration_method, _ip_data, _process_data, NumLib::interpolateCoordinates(), and MathLib::KelvinVector::symmetricTensorToKelvinVector().

◆ localDOF()

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
template<typename SolutionVector>
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 169 of file ThermoHydroMechanicsFEM.h.

173 {
174 unsigned const n_integration_points =
175 _integration_method.getNumberOfPoints();
176
177 auto const [T_prev, p_prev, u_prev] = localDOF(local_x_prev);
178
179 for (unsigned ip = 0; ip < n_integration_points; ip++)
180 {
181 auto& ip_data = _ip_data[ip];
182 auto const& N_u = ip_data.N_u;
183 auto const& dNdx_u = ip_data.dNdx_u;
184
186 std::nullopt, _element.getID(),
191
194
195 auto const x_coord =
198 _element, N_u);
199 auto const B = LinearBMatrix::computeBMatrix<
203
205
207 eps_prev = B * u_prev;
208
209 _ip_data[ip].eps0 =
210 _ip_data[ip].eps0_prev +
211 (1 - _ip_data[ip].phi_fr_prev / _ip_data[ip].porosity) *
213 _ip_data[ip].pushBackState();
214 }
215 }

References _element, _integration_method, _ip_data, _ip_data_output, _is_axially_symmetric, ProcessLib::LinearBMatrix::computeBMatrix(), NumLib::interpolateCoordinates(), NumLib::interpolateXCoordinate(), localDOF(), and 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

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 156 of file ThermoHydroMechanicsFEM.h.

158 {
159 unsigned const n_integration_points =
160 _integration_method.getNumberOfPoints();
161
162 for (unsigned ip = 0; ip < n_integration_points; ip++)
163 {
164 _ip_data_output[ip].velocity.setConstant(
166 }
167 }

References _integration_method, and _ip_data_output.

◆ 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 // TODO: For staggered scheme, overload
185 // LocalAssemblerInterface::setInitialConditions to enable local_x contains
186 // the primary variables from all coupled processes.
187 auto const [T, p, u] = localDOF(local_x);
188
189 double const dt = 0.0;
190
192 auto const& medium = _process_data.media_map.getMedium(_element.getID());
193 auto* const frozen_liquid_phase = medium->hasPhase("FrozenLiquid")
194 ? &medium->phase("FrozenLiquid")
195 : nullptr;
196
197 int const n_integration_points = _integration_method.getNumberOfPoints();
198 for (int ip = 0; ip < n_integration_points; ip++)
199 {
200 auto const& N = _ip_data[ip].N;
201 auto const& N_u = _ip_data[ip].N_u;
203 std::nullopt, _element.getID(),
207 _element, N_u))};
208
209 auto& sigma_eff = _ip_data[ip].sigma_eff;
210 if (_process_data.initial_stress.isTotalStress())
211 {
212 auto const alpha_b =
214 .template value<double>(vars, x_position, t, dt);
215
216 sigma_eff.noalias() += alpha_b * N.dot(p) * Invariants::identity2;
217 }
218 _ip_data[ip].sigma_eff_prev.noalias() = sigma_eff;
219
220 vars.temperature = N.dot(T);
222 {
223 _ip_data[ip].phi_fr =
225 .template value<double>(vars, x_position, t, dt);
226 }
227 }
228}

References _element, _integration_method, _ip_data, _process_data, MaterialPropertyLib::biot_coefficient, MathLib::KelvinVector::Invariants< KelvinVectorSize >::identity2, NumLib::interpolateCoordinates(), localDOF(), setInitialConditionsConcrete(), MaterialPropertyLib::VariableArray::temperature, and MaterialPropertyLib::volume_fraction.

Referenced by setInitialConditionsConcrete().

◆ 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
136 }
137 if (name == "epsilon_m")
138 {
141 }
142 if (name == "epsilon")
143 {
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(
157 [&name](auto const& iv) { return iv.name == name; });
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)
std::size_t setIntegrationPointKelvinVectorData(double const *values, IntegrationPointDataVector &ip_data_vector, MemberType IpData::*const member)

References _element, _integration_method, _ip_data, _process_data, DBUG(), ProcessLib::ThermoHydroMechanics::IntegrationPointData< BMatricesType, ShapeMatricesTypeDisplacement, ShapeMatricesTypePressure, DisplacementDim, ShapeFunctionDisplacement::NPOINTS >::eps, ProcessLib::ThermoHydroMechanics::IntegrationPointData< BMatricesType, ShapeMatricesTypeDisplacement, ShapeMatricesTypePressure, DisplacementDim, ShapeFunctionDisplacement::NPOINTS >::eps_m, ProcessLib::ThermoHydroMechanics::IntegrationPointData< BMatricesType, ShapeMatricesTypeDisplacement, ShapeMatricesTypePressure, DisplacementDim, ShapeFunctionDisplacement::NPOINTS >::material_state_variables, OGS_FATAL, ProcessLib::setIntegrationPointDataMaterialStateVariables(), ProcessLib::setIntegrationPointKelvinVectorData(), and ProcessLib::ThermoHydroMechanics::IntegrationPointData< BMatricesType, ShapeMatricesTypeDisplacement, ShapeMatricesTypePressure, DisplacementDim, ShapeFunctionDisplacement::NPOINTS >::sigma_eff.

◆ 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 233 of file ThermoHydroMechanicsFEM-impl.h.

240{
241 assert(local_x.size() ==
243
244 auto const [T, p, u] = localDOF(local_x);
245 auto const [T_prev, p_prev, u_prev] = localDOF(local_x_prev);
246
247 auto const& solid_material =
249 _process_data.solid_materials, _process_data.material_ids,
250 _element.getID());
251
252 auto const& medium = _process_data.media_map.getMedium(_element.getID());
253 auto const& liquid_phase = medium->phase("AqueousLiquid");
254 auto const& solid_phase = medium->phase("Solid");
255 auto* const frozen_liquid_phase = medium->hasPhase("FrozenLiquid")
256 ? &medium->phase("FrozenLiquid")
257 : nullptr;
259
260 auto const& N_u = ip_data.N_u;
261 auto const& dNdx_u = ip_data.dNdx_u;
262
263 auto const& N = ip_data.N;
264 auto const& dNdx = ip_data.dNdx;
265
266 auto const T_int_pt = N.dot(T);
267 auto const T_prev_int_pt = N.dot(T_prev);
268 double const dT_int_pt = T_int_pt - T_prev_int_pt;
269
270 auto const x_coord =
273 N_u);
274 auto const B =
279
281
282 auto& eps = ip_data.eps;
283 eps.noalias() = B * u;
285 B * u_prev;
286
287 vars.temperature = T_int_pt;
288 double const p_int_pt = N.dot(p);
289 vars.liquid_phase_pressure = p_int_pt;
290
291 vars.liquid_saturation = 1.0;
292
293 auto const solid_density =
295 .template value<double>(vars, x_position, t, dt);
296
297 auto const porosity =
299 .template value<double>(vars, x_position, t, dt);
300 vars.porosity = porosity;
301 ip_data.porosity = porosity;
302
303 crv.alpha_biot =
305 .template value<double>(vars, x_position, t, dt);
306 auto const& alpha = crv.alpha_biot;
307
308 auto const C_el = ip_data.computeElasticTangentStiffness(
309 t, x_position, dt, static_cast<double>(T_int_pt));
311 1 / solid_material.getBulkModulus(t, x_position, &C_el);
312
313 crv.beta_SR = (1 - alpha) * solid_skeleton_compressibility;
314
315 // Set mechanical variables for the intrinsic permeability model
316 // For stress dependent permeability.
317 {
318 auto const& identity2 = Invariants::identity2;
319 auto const sigma_total =
320 (ip_data.sigma_eff - alpha * p_int_pt * identity2).eval();
321 vars.total_stress.emplace<SymmetricTensor>(
323 }
324 // For strain dependent permeability
325 vars.volumetric_strain = Invariants::trace(ip_data.eps);
326 vars.equivalent_plastic_strain =
327 ip_data.material_state_variables->getEquivalentPlasticStrain();
328
329 auto const intrinsic_permeability =
332 .value(vars, x_position, t, dt));
333
334 auto const fluid_density =
336 .template value<double>(vars, x_position, t, dt);
337 ip_data_output.fluid_density = fluid_density;
338 vars.density = fluid_density;
339
340 auto const drho_dp =
342 .template dValue<double>(
344 x_position, t, dt);
345
346 crv.fluid_compressibility = 1 / fluid_density * drho_dp;
347
351
352 // Use the viscosity model to compute the viscosity
353 ip_data_output.viscosity =
355 .template value<double>(vars, x_position, t, dt);
356 crv.K_over_mu = intrinsic_permeability / ip_data_output.viscosity;
357
358 auto const& b = _process_data.specific_body_force;
359
360 // Consider also anisotropic thermal expansion.
361 crv.solid_linear_thermal_expansion_coefficient =
364 .property(
366 .value(vars, x_position, t, dt));
367
370 crv.solid_linear_thermal_expansion_coefficient * dT_int_pt;
371
372 crv.K_pT_thermal_osmosis =
373 (solid_phase.hasProperty(
379 .value(vars, x_position, t, dt))
381
382 GlobalDimVectorType const velocity = -crv.K_over_mu * dNdx * p -
383 crv.K_pT_thermal_osmosis * dNdx * T +
384 crv.K_over_mu * fluid_density * b;
385 ip_data_output.velocity = velocity;
386
387 //
388 // displacement equation, displacement part
389 //
390 auto& eps_m = ip_data.eps_m;
391 auto& eps_m_prev = ip_data.eps_m_prev;
392 eps_m.noalias() = eps_m_prev + eps - eps_prev - dthermal_strain;
393 vars.mechanical_strain
395 eps_m);
396
397 crv.C = ip_data.updateConstitutiveRelation(vars, t, x_position, dt,
399
401
402 crv.beta =
404 (alpha - porosity) *
405 Invariants::trace(crv.solid_linear_thermal_expansion_coefficient);
406
407 //
408 // pressure equation, displacement part.
409 //
410 // Reusing Kup.transpose().
411
412 //
413 // temperature equation, temperature part.
414 //
415 crv.c_f =
418 .template value<double>(vars, x_position, t, dt);
419 crv.effective_thermal_conductivity =
421 medium
422 ->property(
424 .value(vars, x_position, t, dt));
425
426 // Thermal conductivity is moved outside and zero matrix is passed instead
427 // due to multiplication with fluid's density times specific heat capacity.
428 crv.effective_thermal_conductivity.noalias() +=
429 fluid_density * crv.c_f *
431 _process_data.stabilizer, _element.getID(),
433 velocity, 0. /* phi */, 0. /* dispersivity_transversal */,
434 0. /*dispersivity_longitudinal*/);
435
436 double const c_s =
439 .template value<double>(vars, x_position, t, dt);
440
441 // Also modified by freezing terms.
442 crv.effective_volumetric_heat_capacity =
443 porosity * fluid_density * crv.c_f +
444 (1.0 - porosity) * solid_density * c_s;
445
447 {
449 double const phi_fr =
451 .template value<double>(vars, x_position, t, dt);
452 ip_data.phi_fr = phi_fr;
453
454 auto const frozen_liquid_value =
456 {
457 return (*frozen_liquid_phase)[p].template value<double>(
458 vars, x_position, t, dt);
459 };
460
461 double const rho_fr =
463
464 double const c_fr = frozen_liquid_value(
466
467 double const l_fr = frozen_liquid_value(
469
470 double const dphi_fr_dT =
472 .template dValue<double>(
474 x_position, t, dt);
475
476 double const phi_fr_prev = [&]()
477 {
479 vars_prev.temperature = T_prev_int_pt;
481 .template value<double>(vars_prev, x_position, t, dt);
482 }();
483 ip_data.phi_fr_prev = phi_fr_prev;
484
485 // alpha_T^I
490 ->property(
492 .value(vars, x_position, t, dt));
493
497
498 // alpha_{phi_I} -- linear expansion coeff. due to water-to-ice
499 // transition (phase change), and related phase_change_strain term
506 .value(vars, x_position, t, dt));
507
511
512 // eps0 ia a 'history variable' -- a solid matrix strain accrued
513 // prior to the onset of ice forming
514 auto& eps0 = ip_data.eps0;
515 auto const& eps0_prev = ip_data.eps0_prev;
516
517 // definition of eps_m_ice
518 auto& eps_m_ice = ip_data.eps_m_ice;
519 auto const& eps_m_ice_prev = ip_data.eps_m_ice_prev;
520
521 eps_m_ice.noalias() = eps_m_ice_prev + eps - eps_prev -
524
525 vars_ice.mechanical_strain
527 eps_m_ice);
528 auto const C_IR = ip_data.updateConstitutiveRelationIce(
529 *_process_data.ice_constitutive_relation, vars_ice, t, x_position,
531 crv.effective_volumetric_heat_capacity +=
532 -phi_fr * fluid_density * crv.c_f + phi_fr * rho_fr * c_fr -
534
535 // part of dMTT_dT derivative for freezing
536 double const d2phi_fr_dT2 =
538 .template d2Value<double>(
541 dt);
542
543 crv.J_uu_fr = phi_fr * C_IR;
544
545 auto const& sigma_eff_ice = ip_data.sigma_eff_ice;
546 crv.r_u_fr = phi_fr * sigma_eff_ice;
547
549
550 crv.J_TT_fr = ((rho_fr * c_fr - fluid_density * crv.c_f) * dphi_fr_dT +
552 dT_int_pt / dt;
553 }
554 return crv;
555}
typename ShapeMatricesTypePressure::GlobalDimVectorType GlobalDimVectorType
MathLib::KelvinVector::KelvinVectorType< GlobalDim > formKelvinVector(MaterialPropertyLib::PropertyDataType const &values)
A function to form a Kelvin vector from strain or stress alike property like thermal expansivity for ...
Eigen::Matrix< double, GlobalDim, GlobalDim > formEigenTensor(MaterialPropertyLib::PropertyDataType const &values)
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)
static double trace(Eigen::Matrix< double, KelvinVectorSize, 1 > const &v)
Trace of the corresponding tensor.

References _element, _is_axially_symmetric, _process_data, MaterialPropertyLib::biot_coefficient, ProcessLib::LinearBMatrix::computeBMatrix(), ProcessLib::ThermoHydroMechanics::IntegrationPointData< BMatricesType, ShapeMatrixTypeDisplacement, ShapeMatricesTypePressure, DisplacementDim, NPoints >::computeElasticTangentStiffness(), NumLib::computeHydrodynamicDispersion(), MaterialPropertyLib::density, MaterialPropertyLib::VariableArray::density, displacement_size, 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::formEigenTensor(), MaterialPropertyLib::formKelvinVector(), MaterialPropertyLib::getLiquidThermalExpansivity(), MathLib::KelvinVector::Invariants< KelvinVectorSize >::identity2, NumLib::interpolateXCoordinate(), MathLib::KelvinVector::kelvinVectorToSymmetricTensor(), MaterialPropertyLib::liquid_phase_pressure, MaterialPropertyLib::VariableArray::liquid_phase_pressure, MaterialPropertyLib::VariableArray::liquid_saturation, localDOF(), 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, pressure_size, 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, temperature_size, MaterialPropertyLib::thermal_conductivity, MaterialPropertyLib::thermal_expansivity, MaterialPropertyLib::thermal_osmosis_coefficient, MaterialPropertyLib::VariableArray::total_stress, MathLib::KelvinVector::Invariants< KelvinVectorSize >::trace(), 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 assembleWithJacobian(), and postTimestepConcrete().

Member Data Documentation

◆ _element

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
MeshLib::Element const& ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::_element
private

◆ _integration_method

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
NumLib::GenericIntegrationMethod const& ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::_integration_method
private

◆ _ip_data

◆ _ip_data_output

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

◆ _is_axially_symmetric

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

◆ _process_data

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
ThermoHydroMechanicsProcessData<DisplacementDim>& ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::_process_data
private

◆ _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 432 of file ThermoHydroMechanicsFEM.h.

Referenced by assembleWithJacobian().

◆ 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 433 of file ThermoHydroMechanicsFEM.h.

Referenced by assembleWithJacobian(), and updateConstitutiveRelations().

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

Referenced by assembleWithJacobian().

◆ pressure_index

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

◆ 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 431 of file ThermoHydroMechanicsFEM.h.

Referenced by assembleWithJacobian(), and updateConstitutiveRelations().

◆ temperature_index

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

◆ 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 429 of file ThermoHydroMechanicsFEM.h.

Referenced by assembleWithJacobian(), and updateConstitutiveRelations().


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