33template <
typename ShapeFunctionDisplacement,
typename ShapeFunctionPressure,
36 ShapeFunctionPressure, DisplacementDim>::
37 ThermoHydroMechanicsLocalAssembler(
41 bool const is_axially_symmetric,
48 unsigned const n_integration_points =
51 _ip_data.reserve(n_integration_points);
55 auto const shape_matrices_u =
58 DisplacementDim>(e, is_axially_symmetric,
61 auto const shape_matrices_p =
66 auto const& solid_material =
74 if (medium->hasPhase(
"FrozenLiquid") !=
78 "Frozen liquid phase is {:s} and the solid material constitutive "
79 "relation for ice is {:s}. But both must be given (or both "
81 medium->hasPhase(
"FrozenLiquid") ?
"specified" :
"not specified",
86 for (
unsigned ip = 0; ip < n_integration_points; ip++)
88 _ip_data.emplace_back(solid_material);
90 auto const& sm_u = shape_matrices_u[ip];
91 ip_data.integration_weight =
93 sm_u.integralMeasure * sm_u.detJ;
96 ip_data.dNdx_u = sm_u.dNdx;
98 ip_data.N = shape_matrices_p[ip].N;
99 ip_data.dNdx = shape_matrices_p[ip].dNdx;
105template <
typename ShapeFunctionDisplacement,
typename ShapeFunctionPressure,
108 ShapeFunctionDisplacement, ShapeFunctionPressure,
110 double const* values,
111 int const integration_order)
113 if (integration_order !=
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.",
128 "Setting initial conditions for stress from integration "
129 "point data and from a parameter '{:s}' is not possible "
137 if (name ==
"epsilon_m")
142 if (name ==
"epsilon")
147 if (name.starts_with(
"material_state_variable_"))
149 name.remove_prefix(24);
153 auto const& internal_variables =
154 _ip_data[0].solid_material.getInternalVariables();
155 if (
auto const iv = std::find_if(
156 begin(internal_variables), end(internal_variables),
157 [&name](
auto const& iv) {
return iv.name == name; });
158 iv != end(internal_variables))
160 DBUG(
"Setting material state variable '{:s}'", name);
166 int const element_id =
_element.getID();
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);
175template <
typename ShapeFunctionDisplacement,
typename ShapeFunctionPressure,
178 ShapeFunctionDisplacement, ShapeFunctionPressure,
187 auto const [T, p, u] =
localDOF(local_x);
189 double const dt = 0.0;
193 auto*
const frozen_liquid_phase = medium->hasPhase(
"FrozenLiquid")
194 ? &medium->phase(
"FrozenLiquid")
198 for (
int ip = 0; ip < n_integration_points; ip++)
209 auto& sigma_eff =
_ip_data[ip].sigma_eff;
214 .template value<double>(vars, x_position, t, dt);
218 _ip_data[ip].sigma_eff_prev.noalias() = sigma_eff;
221 if (frozen_liquid_phase)
225 .template value<double>(vars, x_position, t, dt);
230template <
typename ShapeFunctionDisplacement,
typename ShapeFunctionPressure,
233 ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim>::
234 updateConstitutiveRelations(
235 Eigen::Ref<Eigen::VectorXd const>
const local_x,
236 Eigen::Ref<Eigen::VectorXd const>
const local_x_prev,
238 double const dt,
IpData& ip_data,
241 assert(local_x.size() ==
244 auto const [T, p, u] =
localDOF(local_x);
245 auto const [T_prev, p_prev, u_prev] =
localDOF(local_x_prev);
247 auto const& solid_material =
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")
260 auto const& N_u = ip_data.
N_u;
261 auto const& dNdx_u = ip_data.
dNdx_u;
263 auto const& N = ip_data.
N;
264 auto const& dNdx = ip_data.
dNdx;
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;
276 ShapeFunctionDisplacement::NPOINTS,
282 auto& eps = ip_data.
eps;
283 eps.noalias() = B * u;
288 double const p_int_pt = N.dot(p);
293 auto const solid_density =
295 .template value<double>(vars, x_position, t, dt);
297 auto const porosity =
299 .template value<double>(vars, x_position, t, dt);
305 .template value<double>(vars, x_position, t, dt);
306 auto const& alpha = crv.alpha_biot;
309 t, x_position, dt,
static_cast<double>(T_int_pt));
310 auto const solid_skeleton_compressibility =
311 1 / solid_material.getBulkModulus(t, x_position, &C_el);
313 crv.beta_SR = (1 - alpha) * solid_skeleton_compressibility;
319 auto const sigma_total =
320 (ip_data.
sigma_eff - alpha * p_int_pt * identity2).eval();
329 auto const intrinsic_permeability =
332 .value(vars, x_position, t, dt));
334 auto const fluid_density =
336 .template value<double>(vars, x_position, t, dt);
342 .template dValue<double>(
346 crv.fluid_compressibility = 1 / fluid_density * drho_dp;
348 double const fluid_volumetric_thermal_expansion_coefficient =
350 liquid_phase, vars, fluid_density, x_position, t, dt);
355 .template value<double>(vars, x_position, t, dt);
356 crv.K_over_mu = intrinsic_permeability / ip_data_output.
viscosity;
361 crv.solid_linear_thermal_expansion_coefficient =
366 .value(vars, x_position, t, dt));
370 crv.solid_linear_thermal_expansion_coefficient * dT_int_pt;
372 crv.K_pT_thermal_osmosis =
373 (solid_phase.hasProperty(
378 thermal_osmosis_coefficient)
379 .value(vars, x_position, t, dt))
380 : Eigen::MatrixXd::Zero(DisplacementDim, DisplacementDim));
383 crv.K_pT_thermal_osmosis * dNdx * T +
384 crv.K_over_mu * fluid_density * b;
390 auto& eps_m = ip_data.
eps_m;
392 eps_m.noalias() = eps_m_prev + eps - eps_prev - dthermal_strain;
400 crv.rho = solid_density * (1 - porosity) + porosity * fluid_density;
403 porosity * fluid_volumetric_thermal_expansion_coefficient +
418 .template value<double>(vars, x_position, t, dt);
419 crv.effective_thermal_conductivity =
424 .value(vars, x_position, t, dt));
428 crv.effective_thermal_conductivity.noalias() +=
429 fluid_density * crv.c_f *
432 GlobalDimMatrixType::Zero(DisplacementDim, DisplacementDim),
439 .template value<double>(vars, x_position, t, dt);
442 crv.effective_volumetric_heat_capacity =
443 porosity * fluid_density * crv.c_f +
444 (1.0 - porosity) * solid_density * c_s;
446 if (frozen_liquid_phase)
449 double const phi_fr =
451 .template value<double>(vars, x_position, t, dt);
454 auto const frozen_liquid_value =
457 return (*frozen_liquid_phase)[p].template value<double>(
458 vars, x_position, t, dt);
461 double const rho_fr =
464 double const c_fr = frozen_liquid_value(
467 double const l_fr = frozen_liquid_value(
470 double const dphi_fr_dT =
472 .template dValue<double>(
476 double const phi_fr_prev = [&]()
481 .template value<double>(vars_prev, x_position, t, dt);
487 DisplacementDim>
const ice_linear_thermal_expansion_coefficient =
492 .value(vars, x_position, t, dt));
495 dthermal_strain_ice =
496 ice_linear_thermal_expansion_coefficient * dT_int_pt;
501 phase_change_expansion_coefficient =
505 phase_change_expansivity)
506 .value(vars, x_position, t, dt));
509 dphase_change_strain = phase_change_expansion_coefficient *
510 (phi_fr - phi_fr_prev) / porosity;
514 auto& eps0 = ip_data.
eps0;
515 auto const& eps0_prev = ip_data.
eps0_prev;
521 eps_m_ice.noalias() = eps_m_ice_prev + eps - eps_prev -
522 (eps0 - eps0_prev) - dthermal_strain_ice -
523 dphase_change_strain;
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 -
533 l_fr * rho_fr * dphi_fr_dT;
536 double const d2phi_fr_dT2 =
538 .template d2Value<double>(
543 crv.J_uu_fr = phi_fr * C_IR;
546 crv.r_u_fr = phi_fr * sigma_eff_ice;
548 crv.J_uT_fr = phi_fr * C_IR * ice_linear_thermal_expansion_coefficient;
550 crv.J_TT_fr = ((rho_fr * c_fr - fluid_density * crv.c_f) * dphi_fr_dT +
551 l_fr * rho_fr * d2phi_fr_dT2) *
559template <
typename ShapeFunctionDisplacement,
typename ShapeFunctionPressure,
562 ShapeFunctionDisplacement, ShapeFunctionPressure,
564 std::vector<double>
const& local_x,
565 std::vector<double>
const&
567 std::vector<double>& local_rhs_data,
568 std::vector<double>& local_Jac_data)
570 assert(local_x.size() ==
574 Eigen::Map<Eigen::VectorXd const>(local_x.data(), local_x.size());
575 auto const x_prev = Eigen::Map<Eigen::VectorXd const>(local_x_prev.data(),
576 local_x_prev.size());
578 auto const [T, p, u] =
localDOF(local_x);
579 auto const [T_prev, p_prev, u_prev] =
localDOF(local_x_prev);
582 typename ShapeMatricesTypeDisplacement::template MatrixType<
589 typename ShapeMatricesTypeDisplacement::template VectorType<
617 typename ShapeMatricesTypeDisplacement::template MatrixType<
623 bool const has_frozen_liquid_phase = medium->hasPhase(
"FrozenLiquid");
625 unsigned const n_integration_points =
628 std::vector<GlobalDimVectorType> ip_flux_vector;
629 double average_velocity_norm = 0.0;
630 ip_flux_vector.reserve(n_integration_points);
632 for (
unsigned ip = 0; ip < n_integration_points; ip++)
645 auto const& w =
_ip_data[ip].integration_weight;
647 auto const& dNdx_u =
_ip_data[ip].dNdx_u;
650 auto const& dNdx =
_ip_data[ip].dNdx;
652 auto const T_int_pt = N.dot(T);
660 ShapeFunctionDisplacement::NPOINTS,
671 if (has_frozen_liquid_phase)
674 .template block<displacement_size, displacement_size>(
676 .noalias() += B.transpose() * crv.J_uu_fr * B * w;
679 .noalias() -= B.transpose() * crv.r_u_fr * w;
682 .template block<displacement_size, temperature_size>(
684 .noalias() -= B.transpose() * crv.J_uT_fr * N * w;
688 .template block<displacement_size, displacement_size>(
690 .noalias() += B.transpose() * crv.C * B * w;
693 .template block<displacement_size, temperature_size>(
695 .noalias() -= B.transpose() * crv.C *
696 crv.solid_linear_thermal_expansion_coefficient * N *
700 .noalias() -= (B.transpose() *
_ip_data[ip].sigma_eff -
701 N_u_op(N_u).transpose() * crv.rho * b) *
713 laplace_p.noalias() += dNdx.transpose() * crv.K_over_mu * dNdx * w;
715 storage_p.noalias() +=
717 (
_ip_data[ip].porosity * crv.fluid_compressibility +
718 (crv.alpha_biot -
_ip_data[ip].porosity) * crv.beta_SR) *
721 laplace_T.noalias() +=
722 dNdx.transpose() * crv.K_pT_thermal_osmosis * dNdx * w;
727 local_rhs.template segment<pressure_size>(
pressure_index).noalias() +=
728 dNdx.transpose() * fluid_density * crv.K_over_mu * b * w;
732 storage_T.noalias() += N.transpose() * crv.beta * N * w;
743 dNdx.transpose() * crv.effective_thermal_conductivity * dNdx * w;
745 ip_flux_vector.emplace_back(velocity * fluid_density * crv.c_f);
746 average_velocity_norm += velocity.norm();
748 if (has_frozen_liquid_phase)
751 .template block<temperature_size, temperature_size>(
753 .noalias() -= N.transpose() * crv.J_TT_fr * N * w;
757 N.transpose() * crv.effective_volumetric_heat_capacity * N * w;
763 dNdx.transpose() * T_int_pt * crv.K_pT_thermal_osmosis * dNdx * w;
766 dKTT_dp.noalias() -= fluid_density * crv.c_f * N.transpose() *
767 (dNdx * T).transpose() * crv.K_over_mu * dNdx * w;
814 average_velocity_norm /
static_cast<double>(n_integration_points), KTT);
820 .noalias() += KTT + MTT / dt;
826 .noalias() += KTp + dKTT_dp;
838 .noalias() += -storage_T / dt + laplace_T;
844 .noalias() += laplace_p + storage_p / dt;
850 .noalias() += Kup.transpose() / dt;
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;
859 .noalias() += Kup * p;
863 KTT * T + MTT * (T - T_prev) / dt;
869template <
typename ShapeFunctionDisplacement,
typename ShapeFunctionPressure,
872 ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim>::
873 getIntPtDarcyVelocity(
875 std::vector<GlobalVector*>
const& ,
876 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& ,
877 std::vector<double>& cache)
const
879 unsigned const n_integration_points =
884 double, DisplacementDim, Eigen::Dynamic, Eigen::RowMajor>>(
885 cache, DisplacementDim, n_integration_points);
887 for (
unsigned ip = 0; ip < n_integration_points; ip++)
895template <
typename ShapeFunctionDisplacement,
typename ShapeFunction,
899 getIntPtFluidDensity(
901 std::vector<GlobalVector*>
const& ,
902 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& ,
903 std::vector<double>& cache)
const
910template <
typename ShapeFunctionDisplacement,
typename ShapeFunction,
916 std::vector<GlobalVector*>
const& ,
917 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& ,
918 std::vector<double>& cache)
const
925template <
typename ShapeFunctionDisplacement,
typename ShapeFunctionPressure,
928 ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim>::
929 computeSecondaryVariableConcrete(
double const ,
double const ,
930 Eigen::VectorXd
const& local_x,
931 Eigen::VectorXd
const& )
933 auto const p = local_x.template segment<pressure_size>(
pressure_index);
937 unsigned const n_integration_points =
940 double phi_fr_avg = 0;
941 double fluid_density_avg = 0;
942 double viscosity_avg = 0;
945 KV sigma_avg = KV::Zero();
947 for (
unsigned ip = 0; ip < n_integration_points; ip++)
954 sigma_avg += ip_data.sigma_eff;
957 phi_fr_avg /= n_integration_points;
958 fluid_density_avg /= n_integration_points;
959 viscosity_avg /= n_integration_points;
960 sigma_avg /= n_integration_points;
968 KV::RowsAtCompileTime]) =
972 ShapeFunctionPressure,
typename ShapeFunctionDisplacement::MeshElement,
977 ShapeFunctionPressure,
typename ShapeFunctionDisplacement::MeshElement,
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
KelvinVector mechanical_strain
KelvinVector total_stress
double equivalent_plastic_strain
double liquid_phase_pressure
std::size_t getID() const
Returns the ID of the element.
MatrixType< _kelvin_vector_size, _number_of_dof > BMatrixType
void setInitialConditionsConcrete(Eigen::VectorXd const local_x, double const t, int const process_id) override
static const int pressure_size
static const int pressure_index
static const int displacement_size
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
ShapeMatrixPolicyType< ShapeFunctionDisplacement, DisplacementDim > ShapeMatricesTypeDisplacement
NumLib::GenericIntegrationMethod const & _integration_method
static const int temperature_size
typename ShapeMatricesTypePressure::GlobalDimVectorType GlobalDimVectorType
Eigen::Matrix< double, KelvinVectorSize, 1 > SymmetricTensor
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
IntegrationPointData< BMatricesType, ShapeMatricesTypeDisplacement, ShapeMatricesTypePressure, DisplacementDim, ShapeFunctionDisplacement::NPOINTS > IpData
ThermoHydroMechanicsProcessData< DisplacementDim > & _process_data
std::vector< IntegrationPointDataForOutput< DisplacementDim >, Eigen::aligned_allocator< IntegrationPointDataForOutput< DisplacementDim > > > _ip_data_output
bool const _is_axially_symmetric
static constexpr auto & N_u_op
SecondaryData< typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType > _secondary_data
MeshLib::Element const & _element
ShapeMatrixPolicyType< ShapeFunctionPressure, DisplacementDim > ShapeMatricesTypePressure
std::vector< IpData, Eigen::aligned_allocator< IpData > > _ip_data
static const int temperature_index
static constexpr auto localDOF(SolutionVector const &x)
static const int displacement_index
auto & selectSolidConstitutiveRelation(SolidMaterialsMap const &constitutive_relations, MeshLib::PropertyVector< int > const *const material_ids, std::size_t const element_id)
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)
@ thermal_osmosis_coefficient
Eigen::Matrix< double, 4, 1 > kelvinVectorToSymmetricTensor(Eigen::Matrix< double, 4, 1, Eigen::ColMajor, 4, 1 > const &v)
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), 1, Eigen::ColMajor > KelvinVectorType
Eigen::Map< Vector > createZeroedVector(std::vector< double > &data, Eigen::VectorXd::Index size)
Eigen::Map< Matrix > createZeroedMatrix(std::vector< double > &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
void assembleAdvectionMatrix(IPData const &ip_data_vector, NumLib::ShapeMatrixCache const &shape_matrix_cache, std::vector< FluxVectorType > const &ip_flux_vector, Eigen::MatrixBase< Derived > &laplacian_matrix)
double interpolateXCoordinate(MeshLib::Element const &e, typename ShapeMatricesType::ShapeMatrices::ShapeType const &N)
void interpolateToHigherOrderNodes(MeshLib::Element const &element, bool const is_axially_symmetric, Eigen::MatrixBase< EigenMatrixType > const &node_values, MeshLib::PropertyVector< double > &interpolated_values_global_vector)
std::vector< typename ShapeMatricesType::ShapeMatrices, Eigen::aligned_allocator< typename ShapeMatricesType::ShapeMatrices > > initShapeMatrices(MeshLib::Element const &e, bool const is_axially_symmetric, IntegrationMethod const &integration_method)
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)
std::array< double, 3 > interpolateCoordinates(MeshLib::Element const &e, typename ShapeMatricesType::ShapeMatrices::ShapeType const &N)
BMatrixType computeBMatrix(DNDX_Type const &dNdx, N_Type const &N, const double radius, const bool is_axially_symmetric)
Fills a B-matrix based on given shape function dN/dx values.
std::vector< double > const & getIntegrationPointScalarData(IntegrationPointDataVector const &ip_data_vector, MemberType IpData::*const member, std::vector< double > &cache)
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)
void setIPDataInitialConditions(std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > const &_integration_point_writer, MeshLib::Properties const &mesh_properties, LocalAssemblersVector &local_assemblers)
MatrixType< ShapeFunction::NPOINTS, ShapeFunction::NPOINTS > NodalMatrixType
static Eigen::Matrix< double, KelvinVectorSize, 1 > const identity2
Kelvin mapping of 2nd order identity tensor.
static double trace(Eigen::Matrix< double, KelvinVectorSize, 1 > const &v)
Trace of the corresponding tensor.
BMatricesType::KelvinVectorType eps_m
BMatricesType::KelvinVectorType eps_m_ice_prev
BMatricesType::KelvinVectorType eps_m_prev
BMatricesType::KelvinMatrixType updateConstitutiveRelation(MaterialPropertyLib::VariableArray const &variable_array, double const t, ParameterLib::SpatialPosition const &x_position, double const dt, double const temperature_prev)
MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > computeElasticTangentStiffness(double const t, ParameterLib::SpatialPosition const &x_position, double const dt, double const temperature)
ShapeMatricesTypePressure::GlobalDimNodalMatrixType dNdx
BMatricesType::KelvinVectorType eps0_prev
std::unique_ptr< typename MaterialLib::Solids::MechanicsBase< DisplacementDim >::MaterialStateVariables > material_state_variables
BMatricesType::KelvinVectorType eps0
BMatricesType::KelvinVectorType sigma_eff
ShapeMatrixTypeDisplacement::NodalRowVectorType N_u
BMatricesType::KelvinMatrixType updateConstitutiveRelationIce(MaterialLib::Solids::MechanicsBase< DisplacementDim > const &ice_constitutive_relation, MaterialPropertyLib::VariableArray const &variable_array, double const t, ParameterLib::SpatialPosition const &x_position, double const dt, double const temperature_prev)
BMatricesType::KelvinVectorType eps_m_ice
BMatricesType::KelvinVectorType eps
ShapeMatricesTypePressure::NodalRowVectorType N
ShapeMatrixTypeDisplacement::GlobalDimNodalMatrixType dNdx_u
BMatricesType::KelvinVectorType sigma_eff_ice