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 45 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
int getNumberOfVectorElementsForDeformation () 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.
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 56 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 59 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 64 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 260 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 49 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 66 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 29 of file ThermoHydroMechanicsFEM-impl.h.

38 _element(e),
40{
41 unsigned const n_integration_points =
42 _integration_method.getNumberOfPoints();
43
47
48 auto const shape_matrices_u =
53
54 auto const shape_matrices_p =
58
59 auto const& solid_material =
61 _process_data.solid_materials, _process_data.material_ids,
62 e.getID());
63
64 // Consistency check: if frozen liquid phase is given, then the constitutive
65 // relation for ice must also be given, and vice versa.
66 auto const& medium = _process_data.media_map.getMedium(_element.getID());
67 if (medium->hasPhase("FrozenLiquid") !=
68 (_process_data.ice_constitutive_relation != nullptr))
69 {
71 "Frozen liquid phase is {:s} and the solid material constitutive "
72 "relation for ice is {:s}. But both must be given (or both "
73 "omitted).",
74 medium->hasPhase("FrozenLiquid") ? "specified" : "not specified",
75 _process_data.ice_constitutive_relation != nullptr
76 ? "specified"
77 : "not specified");
78 }
79 for (unsigned ip = 0; ip < n_integration_points; ip++)
80 {
81 _ip_data.emplace_back(solid_material);
82 auto& ip_data = _ip_data[ip];
83 auto const& sm_u = shape_matrices_u[ip];
84 ip_data.integration_weight =
85 _integration_method.getWeightedPoint(ip).getWeight() *
86 sm_u.integralMeasure * sm_u.detJ;
87
88 ip_data.N_u = sm_u.N;
89 ip_data.dNdx_u = sm_u.dNdx;
90
92 ip_data.dNdx = shape_matrices_p[ip].dNdx;
93
95 }
96}
#define OGS_FATAL(...)
Definition Error.h:19
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 90 of file ThermoHydroMechanicsFEM.h.

96 {
98 "ThermoHydroMechanicsLocalAssembler: assembly without Jacobian is "
99 "not implemented.");
100 }

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

674{
675 assert(local_x.size() ==
677
678 auto const x =
681 local_x_prev.size());
682
683 auto const [T, p, u] = localDOF(local_x);
684 auto const [T_prev, p_prev, u_prev] = localDOF(local_x_prev);
685
692
697
700
703
706
709
712
715
718
721
724
727 Kup;
729
730 auto const& medium = _process_data.media_map.getMedium(_element.getID());
731 bool const has_frozen_liquid_phase = medium->hasPhase("FrozenLiquid");
732
733 unsigned const n_integration_points =
734 _integration_method.getNumberOfPoints();
735
737 double average_velocity_norm = 0.0;
739
740 for (unsigned ip = 0; ip < n_integration_points; ip++)
741 {
742 auto& ip_data = _ip_data[ip];
743 auto const& N_u = ip_data.N_u;
745 std::nullopt, _element.getID(),
749 _element, N_u))};
750
751 auto const crv = updateConstitutiveRelations(
753
754 auto const& w = ip_data.integration_weight;
755
756 auto const& dNdx_u = ip_data.dNdx_u;
757
758 auto const& N = ip_data.N;
759 auto const& dNdx = ip_data.dNdx;
760
761 auto const T_int_pt = N.dot(T);
762
763 auto const x_coord =
764 x_position.getCoordinates().value()[0]; // r for axisymmetry
765 auto const B =
770
771 auto const& b = _process_data.specific_body_force;
772 auto const velocity = _ip_data_output[ip].velocity;
773
774 //
775 // displacement equation, displacement part
776 //
777
779 {
783 .noalias() += B.transpose() * crv.J_uu_fr * B * w;
784
786 .noalias() -= B.transpose() * crv.r_u_fr * w;
787
791 .noalias() -= B.transpose() * crv.J_uT_fr * N * w;
792 }
793
797 .noalias() += B.transpose() * crv.C * B * w;
798
802 .noalias() -= B.transpose() * crv.C *
803 crv.solid_linear_thermal_expansion_coefficient * N *
804 w;
805
807 .noalias() -= (B.transpose() * ip_data.sigma_eff -
808 N_u_op(N_u).transpose() * crv.rho * b) *
809 w;
810
811 //
812 // displacement equation, pressure part (K_up)
813 //
814 Kup.noalias() +=
815 B.transpose() * Invariants::identity2 * N * (crv.alpha_biot * w);
816
817 double const fluid_density = _ip_data_output[ip].fluid_density;
819 {
820 Kup.noalias() +=
821 B.transpose() * Invariants::identity2 * N *
822 (crv.alpha_biot * ip_data.phi_fr / ip_data.porosity *
823 (_ip_data_output[ip].rho_fr / fluid_density - 1)) *
824 w;
825 }
826
827 //
828 // pressure equation, pressure part (K_pp and M_pp).
829 //
830 laplace_p.noalias() +=
831 dNdx.transpose() * crv.K_over_mu * dNdx * (crv.k_rel * w);
835 .noalias() += dNdx.transpose() * crv.K_over_mu * (dNdx * p) * N *
836 (crv.dk_rel_dT * w);
837
838 storage_p.noalias() +=
839 N.transpose() * N *
840 ((ip_data.porosity * crv.fluid_compressibility +
841 (crv.alpha_biot - ip_data.porosity) * crv.beta_SR) *
842 w);
843
845 {
846 storage_p.noalias() += N.transpose() * crv.storage_p_fr * N * w;
847
851 .noalias() += N.transpose() * crv.J_pT_fr * N * w;
852 }
853
854 laplace_T.noalias() +=
855 dNdx.transpose() * crv.K_pT_thermal_osmosis * dNdx * w;
856 //
857 // RHS, pressure part
858 //
859 local_rhs.template segment<pressure_size>(pressure_index).noalias() +=
860 dNdx.transpose() * crv.K_over_mu * b *
861 (fluid_density * crv.k_rel * w);
865 .noalias() -=
866 dNdx.transpose() * crv.K_over_mu * b * N *
867 ((fluid_density * crv.dk_rel_dT + crv.drho_LR_dT * crv.k_rel) * w);
868 //
869 // pressure equation, temperature part (M_pT)
870 //
871 storage_T.noalias() += N.transpose() * crv.beta * N * w;
872
874 {
875 storage_T.noalias() += N.transpose() * crv.storage_T_fr * N * w;
876 }
877
878 //
879 // pressure equation, displacement part.
880 //
881 // Reusing Kup.transpose().
882
883 //
884 // temperature equation, temperature part.
885 //
886 KTT.noalias() +=
887 dNdx.transpose() * crv.effective_thermal_conductivity * dNdx * w;
888 dKTT_dT_T.noalias() +=
889 dNdx.transpose() * crv.dlambda_eff_dT * dNdx * T * N * w;
890
891 ip_flux_vector.emplace_back(velocity * fluid_density * crv.c_f);
892 // Without any flux correction the flux derivative is as follows. The
893 // contribution to KTT is different if any stabilization scheme is used,
894 // but this is ignored for the moment.
896 crv.dvelocity_dT * fluid_density * crv.c_f +
897 velocity * crv.drho_LR_dT * crv.c_f;
898 dKTT_dT_T.noalias() +=
899 N.transpose() * dip_flux_vector_dT.transpose() * dNdx * T * N * w;
901
902 MTT.noalias() +=
903 N.transpose() * crv.effective_volumetric_heat_capacity * N * w;
907 .noalias() += N.transpose() * crv.J_TT * N * w;
908
909 //
910 // temperature equation, pressure part
911 //
912 KTp.noalias() +=
913 dNdx.transpose() * crv.K_pT_thermal_osmosis * dNdx * (T_int_pt * w);
914
915 // linearized darcy
916 dKTT_dp_T.noalias() -= N.transpose() * (dNdx * T).transpose() *
917 crv.K_over_mu * dNdx *
918 (fluid_density * crv.c_f * crv.k_rel * w);
919
920 /* TODO (Joerg) Temperature changes due to thermal dilatation of the
921 * fluid, which are usually discarded as being very small.
922 * Zhou et al. (10.1016/S0020-7683(98)00089-4) states that:
923 * "Biot (1956) neglected this term and it is included here for
924 * completeness"
925 * Keeping the code here in the case these are needed for the named
926 * effects in the future.
927 if (fluid_compressibility != 0)
928 {
929 auto const C_el = ip_data.computeElasticTangentStiffness(
930 t, x_position, dt, static_cast<double>(T_int_pt));
931 auto const solid_skeleton_compressibility =
932 1 / solid_material.getBulkModulus(t, x_position, &C_el);
933 double const fluid_volumetric_thermal_expansion_coefficient =
934 MaterialPropertyLib::getLiquidThermalExpansivity(
935 liquid_phase, vars, fluid_density, x_position, t, dt);
936
937 KTT.noalias() +=
938 dNdx.transpose() *
939 (-T_int_pt * fluid_volumetric_thermal_expansion_coefficient *
940 K_pT_thermal_osmosis / fluid_compressibility) *
941 dNdx * w;
942
943 local_rhs.template segment<temperature_size>(temperature_index)
944 .noalias() +=
945 dNdx.transpose() *
946 (-T_int_pt * fluid_volumetric_thermal_expansion_coefficient /
947 fluid_compressibility) *
948 fluid_density * crv.k_rel * K_over_mu * b * w;
949 MTu part for rhs and Jacobian:
950 (-T_int_pt *
951 Invariants::trace(solid_linear_thermal_expansion_coefficient) /
952 solid_skeleton_compressibility) *
953 N.transpose() * identity2.transpose() * B * w;
954 KTp part for rhs and Jacobian:
955 dNdx.transpose() *
956 (T_int_pt * fluid_volumetric_thermal_expansion_coefficient *
957 crv.k_rel * K_over_mu / fluid_compressibility) *
958 dNdx * w;
959 }
960 */
961 }
962
965 average_velocity_norm / static_cast<double>(n_integration_points), KTT);
966
967 // temperature equation, temperature part
971 .noalias() += KTT + dKTT_dT_T + MTT / dt;
972
973 // temperature equation, pressure part
977 .noalias() += KTp + dKTT_dp_T;
978
979 // displacement equation, pressure part
983 .noalias() -= Kup;
984
985 // pressure equation, temperature part.
989 .noalias() += -storage_T / dt + laplace_T;
990
991 // pressure equation, pressure part.
995 .noalias() += laplace_p + storage_p / dt;
996
997 // pressure equation, displacement part.
1001 .noalias() += Kup.transpose() / dt;
1002
1003 // pressure equation (f_p)
1004 local_rhs.template segment<pressure_size>(pressure_index).noalias() -=
1005 laplace_p * p + laplace_T * T + storage_p * (p - p_prev) / dt -
1006 storage_T * (T - T_prev) / dt + Kup.transpose() * (u - u_prev) / dt;
1007
1008 // displacement equation (f_u)
1010 .noalias() += Kup * p;
1011
1012 // temperature equation (f_T)
1014 KTT * T + MTT * (T - T_prev) / dt;
1015
1017 KTp * p;
1018}
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
typename ShapeMatricesTypePressure::GlobalDimVectorType GlobalDimVectorType
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, ParameterLib::SpatialPosition::getCoordinates(), MathLib::KelvinVector::Invariants< KelvinVectorSize >::identity2, NumLib::interpolateCoordinates(), 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 1079 of file ThermoHydroMechanicsFEM-impl.h.

1083{
1084 auto const p = local_x.template segment<pressure_size>(pressure_index);
1085 auto const T =
1087
1088 unsigned const n_integration_points =
1089 _integration_method.getNumberOfPoints();
1090
1091 double phi_fr_avg = 0;
1092 double fluid_density_avg = 0;
1093 double viscosity_avg = 0;
1094
1096 KV sigma_avg = KV::Zero();
1097
1098 for (unsigned ip = 0; ip < n_integration_points; ip++)
1099 {
1100 auto& ip_data = _ip_data[ip];
1101
1102 phi_fr_avg += ip_data.phi_fr;
1103 fluid_density_avg += _ip_data_output[ip].fluid_density;
1104 viscosity_avg += _ip_data_output[ip].viscosity;
1105 sigma_avg += ip_data.sigma_eff;
1106 }
1107
1112
1113 (*_process_data.element_phi_fr)[_element.getID()] = phi_fr_avg;
1114 (*_process_data.element_fluid_density)[_element.getID()] =
1116 (*_process_data.element_viscosity)[_element.getID()] = viscosity_avg;
1117
1118 Eigen::Map<KV>(&(*_process_data.element_stresses)[_element.getID() *
1121
1125 *_process_data.pressure_interpolated);
1126
1130 *_process_data.temperature_interpolated);
1131}
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 338 of file ThermoHydroMechanicsFEM.h.

339 {
340 constexpr int kelvin_vector_size =
342
345 { return getIntPtEpsilon(0, {}, {}, values); });
346 }
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 299 of file ThermoHydroMechanicsFEM.h.

300 {
301 constexpr int kelvin_vector_size =
303
306 { return getIntPtEpsilon0(0, {}, {}, values); });
307 }
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 318 of file ThermoHydroMechanicsFEM.h.

319 {
320 constexpr int kelvin_vector_size =
322
325 { return getIntPtEpsilonM(0, {}, {}, values); });
326 }
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 1023 of file ThermoHydroMechanicsFEM-impl.h.

1029{
1030 unsigned const n_integration_points =
1031 _integration_method.getNumberOfPoints();
1032
1033 cache.clear();
1037
1038 for (unsigned ip = 0; ip < n_integration_points; ip++)
1039 {
1040 cache_matrix.col(ip).noalias() = _ip_data_output[ip].velocity;
1041 }
1042
1043 return cache;
1044}

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

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

374 {
375 return _process_data.material_ids == nullptr
376 ? 0
377 : (*_process_data.material_ids)[_element.getID()];
378 }

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

394 {
395 return *_ip_data[integration_point].material_state_variables;
396 }

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

369 {
370 return _integration_method.getNumberOfPoints();
371 }

References _integration_method.

◆ getNumberOfVectorElementsForDeformation()

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

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 252 of file ThermoHydroMechanicsFEM.h.

253 {
254 return displacement_size;
255 }

References displacement_size.

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

214 {
215 auto const& N_u = _secondary_data.N_u[integration_point];
216
217 // assumes N is stored contiguously in memory
218 return Eigen::Map<const Eigen::RowVectorXd>(N_u.data(), N_u.size());
219 }

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

231 {
232 constexpr int kelvin_vector_size =
234
237 { return getIntPtSigma(0, {}, {}, values); });
238 }
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 108 of file ThermoHydroMechanicsFEM.h.

109 {
110 unsigned const n_integration_points =
111 _integration_method.getNumberOfPoints();
112
113 for (unsigned ip = 0; ip < n_integration_points; ip++)
114 {
115 auto& ip_data = _ip_data[ip];
116
118 std::nullopt, _element.getID(),
123
125 if (_process_data.initial_stress.value)
126 {
127 ip_data.sigma_eff =
129 DisplacementDim>((*_process_data.initial_stress.value)(
131 independent
132 */
133 ,
134 x_position));
135 }
136
137 double const t = 0; // TODO (naumov) pass t from top
138 ip_data.solid_material.initializeInternalStateVariables(
139 t, x_position, *ip_data.material_state_variables);
140
141 ip_data.pushBackState();
142 }
143 }

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

166 {
167 unsigned const n_integration_points =
168 _integration_method.getNumberOfPoints();
169
170 auto const [T_prev, p_prev, u_prev] = localDOF(local_x_prev);
171
172 for (unsigned ip = 0; ip < n_integration_points; ip++)
173 {
174 auto& ip_data = _ip_data[ip];
175 auto const& N_u = ip_data.N_u;
176 auto const& dNdx_u = ip_data.dNdx_u;
177
179 std::nullopt, _element.getID(),
184
187
188 auto const x_coord =
189 x_position.getCoordinates().value()[0]; // r for axisymmetry
190 auto const B = LinearBMatrix::computeBMatrix<
194
196
198 eps_prev = B * u_prev;
199
200 _ip_data[ip].eps0 =
201 _ip_data[ip].eps0_prev +
202 (1 - _ip_data[ip].phi_fr_prev / _ip_data[ip].porosity) *
204 _ip_data[ip].pushBackState();
205 }
206 }

References _element, _integration_method, _ip_data, _ip_data_output, _is_axially_symmetric, ProcessLib::LinearBMatrix::computeBMatrix(), ParameterLib::SpatialPosition::getCoordinates(), NumLib::interpolateCoordinates(), 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 149 of file ThermoHydroMechanicsFEM.h.

151 {
152 unsigned const n_integration_points =
153 _integration_method.getNumberOfPoints();
154
155 for (unsigned ip = 0; ip < n_integration_points; ip++)
156 {
157 _ip_data_output[ip].velocity.setConstant(
159 }
160 }

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

176{
177 // TODO: For staggered scheme, overload
178 // LocalAssemblerInterface::setInitialConditions to enable local_x contains
179 // the primary variables from all coupled processes.
180 auto const [T, p, u] = localDOF(local_x);
181
182 double const dt = 0.0;
183
185 auto const& medium = _process_data.media_map.getMedium(_element.getID());
186 auto* const frozen_liquid_phase = medium->hasPhase("FrozenLiquid")
187 ? &medium->phase("FrozenLiquid")
188 : nullptr;
189
190 int const n_integration_points = _integration_method.getNumberOfPoints();
191 for (int ip = 0; ip < n_integration_points; ip++)
192 {
193 auto& ip_data = _ip_data[ip];
194 auto const& N = ip_data.N;
195 auto const& N_u = ip_data.N_u;
197 std::nullopt, _element.getID(),
201 _element, N_u))};
202
203 auto& sigma_eff = ip_data.sigma_eff;
204 if (_process_data.initial_stress.isTotalStress())
205 {
206 auto const alpha_b =
208 .template value<double>(vars, x_position, t, dt);
209
210 sigma_eff.noalias() += alpha_b * N.dot(p) * Invariants::identity2;
211 }
212 ip_data.sigma_eff_prev.noalias() = sigma_eff;
213
214 vars.temperature = N.dot(T);
216 {
217 ip_data.phi_fr =
219 .template value<double>(vars, x_position, t, dt);
220 }
221 }
222}

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

105{
106 if (integration_order !=
107 static_cast<int>(_integration_method.getIntegrationOrder()))
108 {
109 OGS_FATAL(
110 "Setting integration point initial conditions; The integration "
111 "order of the local assembler for element {:d} is different from "
112 "the integration order in the initial condition.",
113 _element.getID());
114 }
115
116 if (name == "sigma")
117 {
118 if (_process_data.initial_stress.value)
119 {
120 OGS_FATAL(
121 "Setting initial conditions for stress from integration "
122 "point data and from a parameter '{:s}' is not possible "
123 "simultaneously.",
124 _process_data.initial_stress.value->name);
125 }
126
129 }
130 if (name == "epsilon_m")
131 {
134 }
135 if (name == "epsilon")
136 {
139 }
140 if (name.starts_with("material_state_variable_"))
141 {
142 name.remove_prefix(24);
143
144 // Using first ip data for solid material. TODO (naumov) move solid
145 // material into element, store only material state in IPs.
146 auto const& internal_variables =
147 _ip_data[0].solid_material.getInternalVariables();
148 if (auto const iv = std::find_if(
150 [&name](auto const& iv) { return iv.name == name; });
152 {
153 DBUG("Setting material state variable '{:s}'", name);
156 iv->reference);
157 }
158
159 int const element_id = _element.getID();
160 DBUG(
161 "The solid material of element {:d} (material ID {:d}) does not "
162 "have an internal state variable called {:s}.",
163 element_id, (*_process_data.material_ids)[element_id], name);
164 }
165
166 return 0;
167}
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:22
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 227 of file ThermoHydroMechanicsFEM-impl.h.

234{
235 assert(local_x.size() ==
237
238 auto const [T, p, u] = localDOF(local_x);
239 auto const [T_prev, p_prev, u_prev] = localDOF(local_x_prev);
240
241 auto const& solid_material =
243 _process_data.solid_materials, _process_data.material_ids,
244 _element.getID());
245
246 auto const& medium = _process_data.media_map.getMedium(_element.getID());
247 auto const& liquid_phase = medium->phase("AqueousLiquid");
248 auto const& solid_phase = medium->phase("Solid");
249 auto* const frozen_liquid_phase = medium->hasPhase("FrozenLiquid")
250 ? &medium->phase("FrozenLiquid")
251 : nullptr;
253
254 auto const& N_u = ip_data.N_u;
255 auto const& dNdx_u = ip_data.dNdx_u;
256
257 auto const& N = ip_data.N;
258 auto const& dNdx = ip_data.dNdx;
259
260 auto const T_int_pt = N.dot(T);
261 auto const T_prev_int_pt = N.dot(T_prev);
262 double const dT_int_pt = T_int_pt - T_prev_int_pt;
263
264 auto const x_coord =
265 x_position.getCoordinates().value()[0]; // r for axisymmetry
266 auto const B =
271
273
274 auto& eps = ip_data.eps;
275 eps.noalias() = B * u;
277 B * u_prev;
278
279 vars.temperature = T_int_pt;
280 double const p_int_pt = N.dot(p);
281 vars.liquid_phase_pressure = p_int_pt;
282 double const p_prev_int_pt = N.dot(p_prev);
283 double const dp_int_pt = p_int_pt - p_prev_int_pt;
284
285 vars.liquid_saturation = 1.0;
286
287 auto const solid_density =
289 .template value<double>(vars, x_position, t, dt);
290
291 auto const drho_SR_dT =
293 .template dValue<double>(vars,
295 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
348 crv.drho_LR_dT =
350 .template dValue<double>(vars,
352 x_position, t, dt);
353
357
358 // Use the viscosity model to compute the viscosity
359 ip_data_output.viscosity =
361 .template value<double>(vars, x_position, t, dt);
362 crv.K_over_mu = intrinsic_permeability / ip_data_output.viscosity;
363
364 crv.k_rel = 1.;
365 crv.dk_rel_dT = 0;
367 {
368 ip_data.phi_fr =
370 .template value<double>(vars, x_position, t, dt);
371
372 double const dphi_fr_dT =
374 .template dValue<double>(
376 x_position, t, dt);
377
378 // Set ice_volume_fraction variable for relative permeability
379 // calculation ice_volume_fraction = fraction of pore space occupied by
380 // ice
381 vars.ice_volume_fraction = ip_data.phi_fr / porosity;
382
383 crv.k_rel =
385 .property(
387 .template value<double>(vars, x_position, t, dt);
388
389 // dk_rel/dT = (dk_rel/dphi_ice) * (dphi_ice/dT)
390 // where dphi_ice/dT = dphi_fr_dT / porosity
391 auto const dk_rel_dphi_ice =
393 .property(
395 .template dValue<double>(
397 x_position, t, dt);
398 crv.dk_rel_dT = dk_rel_dphi_ice * dphi_fr_dT / porosity;
399 }
400
401 auto const& b = _process_data.specific_body_force;
402
403 // Consider also anisotropic thermal expansion.
404 crv.solid_linear_thermal_expansion_coefficient =
407 .property(
409 .value(vars, x_position, t, dt));
410
413 crv.solid_linear_thermal_expansion_coefficient * dT_int_pt;
414
415 crv.K_pT_thermal_osmosis =
416 (solid_phase.hasProperty(
422 .value(vars, x_position, t, dt))
424
426 -crv.k_rel * crv.K_over_mu * dNdx * p -
427 crv.K_pT_thermal_osmosis * dNdx * T +
428 (fluid_density * crv.k_rel) * crv.K_over_mu * b;
429 ip_data_output.velocity = velocity;
430 crv.dvelocity_dT =
431 -crv.dk_rel_dT * crv.K_over_mu * dNdx * p +
432 // TODO(naumov): - crv.K_pT_thermal_osmosis * dNdx * dT_dT +
433 (crv.dk_rel_dT * fluid_density + crv.k_rel * crv.drho_LR_dT) *
434 crv.K_over_mu * b;
435
436 //
437 // displacement equation, displacement part
438 //
439 auto& eps_m = ip_data.eps_m;
440 auto& eps_m_prev = ip_data.eps_m_prev;
441 eps_m.noalias() = eps_m_prev + eps - eps_prev - dthermal_strain;
442 vars.mechanical_strain
444 eps_m);
445
446 crv.C = ip_data.updateConstitutiveRelation(vars, t, x_position, dt,
448
450
451 crv.beta =
453 (alpha - porosity) *
454 Invariants::trace(crv.solid_linear_thermal_expansion_coefficient);
455
456 //
457 // pressure equation, displacement part.
458 //
459 // Reusing Kup.transpose().
460
461 //
462 // temperature equation, temperature part.
463 //
464 crv.c_f =
467 .template value<double>(vars, x_position, t, dt);
468 crv.effective_thermal_conductivity =
470 medium
471 ->property(
473 .value(vars, x_position, t, dt));
474
476 medium
479 x_position, t, dt));
480
481 // Thermal conductivity is moved outside and zero matrix is passed instead
482 // due to multiplication with fluid's density times specific heat capacity.
483 crv.effective_thermal_conductivity.noalias() +=
484 fluid_density * crv.c_f *
486 _process_data.stabilizer, _element.getID(),
488 velocity, 0. /* phi */, 0. /* dispersivity_transversal */,
489 0. /*dispersivity_longitudinal*/);
490
491 double const c_s =
494 .template value<double>(vars, x_position, t, dt);
495
496 // Also modified by freezing terms.
497 crv.effective_volumetric_heat_capacity =
498 porosity * fluid_density * crv.c_f +
499 (1.0 - porosity) * solid_density * c_s;
500 double dC_eff_dT = porosity * crv.drho_LR_dT * crv.c_f +
501 (1.0 - porosity) * drho_SR_dT * c_s;
502
504 {
506 double const phi_fr = ip_data.phi_fr;
507
508 auto const frozen_liquid_value =
510 {
511 return (*frozen_liquid_phase)[p].template value<double>(
512 vars, x_position, t, dt);
513 };
514
515 double const c_fr = frozen_liquid_value(
517
518 double const l_fr = frozen_liquid_value(
520
521 double const dphi_fr_dT =
523 .template dValue<double>(
525 x_position, t, dt);
526
527 double const d2phi_fr_dT2 =
529 .template d2Value<double>(
532 dt);
533
534 double const phi_fr_prev = [&]()
535 {
537 vars_prev.temperature = T_prev_int_pt;
539 .template value<double>(vars_prev, x_position, t, dt);
540 }();
541 ip_data.phi_fr_prev = phi_fr_prev;
542
543 double const rho_fr =
545 ip_data_output.rho_fr = rho_fr;
546
547 crv.rho += ip_data.phi_fr * rho_fr - ip_data.phi_fr * fluid_density;
548 crv.mass_exchange =
550 double const dmass_exchange_dT =
552 dphi_fr_dT * porosity * rho_fr * crv.drho_LR_dT /
554
555 // alpha_T^I
560 ->property(
562 .value(vars, x_position, t, dt));
563
567
568 crv.beta_T_SI =
569 porosity *
571 (alpha - porosity) *
573 crv.solid_linear_thermal_expansion_coefficient);
574
575 // alpha_{phi_I} -- linear expansion coeff. due to water-to-ice
576 // transition (phase change), and related phase_change_strain term
583 .value(vars, x_position, t, dt));
584
588
589 // eps0 ia a 'history variable' -- a solid matrix strain accrued
590 // prior to the onset of ice forming
591 auto& eps0 = ip_data.eps0;
592 auto const& eps0_prev = ip_data.eps0_prev;
593
594 // definition of eps_m_ice
595 auto& eps_m_ice = ip_data.eps_m_ice;
596 auto const& eps_m_ice_prev = ip_data.eps_m_ice_prev;
597
598 eps_m_ice.noalias() = eps_m_ice_prev + eps - eps_prev -
601
602 vars_ice.mechanical_strain
604 eps_m_ice);
605 auto const C_IR = ip_data.updateConstitutiveRelationIce(
606 *_process_data.ice_constitutive_relation, vars_ice, t, x_position,
608
609 auto const C_el_ice = ip_data.computeElasticTangentStiffnessIce(
610 *_process_data.ice_constitutive_relation, t, x_position, dt,
611 static_cast<double>(T_int_pt));
612 crv.beta_IR =
613 1. / _process_data.ice_constitutive_relation->getBulkModulus(
615
616 crv.effective_volumetric_heat_capacity +=
617 -phi_fr * fluid_density * crv.c_f + phi_fr * rho_fr * c_fr -
619
620 crv.J_uu_fr = phi_fr * C_IR;
621
622 auto const& sigma_eff_ice = ip_data.sigma_eff_ice;
623 crv.r_u_fr = phi_fr * sigma_eff_ice;
624
626
627 // part of dMTT_dT derivative for freezing
629 phi_fr * crv.drho_LR_dT * crv.c_f +
631 double const storage_p_fr_coeff =
632 (porosity * crv.beta_IR + (alpha - porosity) * crv.beta_SR) *
634 (porosity * crv.fluid_compressibility +
635 (alpha - porosity) * crv.beta_SR);
636 crv.storage_p_fr = phi_fr / porosity * storage_p_fr_coeff;
637
638 double const dstorage_p_fr_coeff_dT =
639 (porosity * crv.beta_IR + (alpha - porosity) * crv.beta_SR) *
640 crv.drho_LR_dT / (fluid_density * fluid_density);
641
642 crv.J_pT_fr = (dphi_fr_dT * storage_p_fr_coeff +
645
646 crv.storage_T_fr =
647 phi_fr / porosity *
648 (crv.beta_T_SI * rho_fr / fluid_density - crv.beta) -
649 crv.mass_exchange;
650 double const dstorage_T_fr_dT =
652 (crv.beta_T_SI * rho_fr / fluid_density - crv.beta) +
653 phi_fr / porosity * crv.beta_T_SI * rho_fr * crv.drho_LR_dT /
656 crv.J_pT_fr += dstorage_T_fr_dT * dT_int_pt / dt;
657 }
658 crv.J_TT = dC_eff_dT * dT_int_pt / dt;
659 return crv;
660}
constexpr Eigen::Matrix< double, GlobalDim, GlobalDim > formEigenTensor(MaterialPropertyLib::PropertyDataType const &values)
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 ...
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(), ProcessLib::ThermoHydroMechanics::IntegrationPointData< BMatricesType, ShapeMatrixTypeDisplacement, ShapeMatricesTypePressure, DisplacementDim, NPoints >::computeElasticTangentStiffnessIce(), 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(), ParameterLib::SpatialPosition::getCoordinates(), MaterialPropertyLib::getLiquidThermalExpansivity(), MaterialPropertyLib::ice_volume_fraction, MaterialPropertyLib::VariableArray::ice_volume_fraction, MathLib::KelvinVector::Invariants< KelvinVectorSize >::identity2, 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, MaterialPropertyLib::relative_permeability, ProcessLib::ThermoHydroMechanics::IntegrationPointDataForOutput< DisplacementDim >::rho_fr, 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 428 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 429 of file ThermoHydroMechanicsFEM.h.

Referenced by assembleWithJacobian(), getNumberOfVectorElementsForDeformation(), 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 62 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 68 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 427 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 425 of file ThermoHydroMechanicsFEM.h.

Referenced by assembleWithJacobian(), and updateConstitutiveRelations().


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