14#include <spdlog/fmt/fmt.h>
36template <
typename ShapeFunctionDisplacement,
typename ShapeFunction,
37 int DisplacementDim,
typename ConstitutiveTraits>
39 DisplacementDim, ConstitutiveTraits>::
40 ThermoRichardsMechanicsLocalAssembler(
44 bool const is_axially_symmetric,
48 e, integration_method, is_axially_symmetric, process_data)
50 unsigned const n_integration_points =
53 ip_data_.resize(n_integration_points);
55 auto const shape_matrices_u =
58 DisplacementDim>(e, is_axially_symmetric,
61 auto const shape_matrices =
63 DisplacementDim>(e, is_axially_symmetric,
66 for (
unsigned ip = 0; ip < n_integration_points; ip++)
69 auto const& sm_u = shape_matrices_u[ip];
72 sm_u.integralMeasure * sm_u.detJ;
75 ip_data.dNdx_u = sm_u.dNdx;
78 ip_data.N_p = shape_matrices[ip].N;
79 ip_data.dNdx_p = shape_matrices[ip].dNdx;
83template <
typename ShapeFunctionDisplacement,
typename ShapeFunction,
84 int DisplacementDim,
typename ConstitutiveTraits>
92 assert(local_x.size() ==
95 auto const p_L = local_x.template segment<pressure_size>(
pressure_index);
99 constexpr double dt = std::numeric_limits<double>::quiet_NaN();
105 typename ConstitutiveTraits::ConstitutiveSetting
const constitutive_setting;
106 auto models = ConstitutiveTraits::createConstitutiveModels(
109 unsigned const n_integration_points =
111 for (
unsigned ip = 0; ip < n_integration_points; ip++)
117 std::nullopt, this->
element_.getID(),
142 .template value<double>(variables, x_position, t, dt);
143 std::get<PrevState<SaturationData>>(this->
prev_states_[ip])->S_L =
149 constitutive_setting.init(models, t, dt, x_position, media_data,
158 "The medium has no saturation property required to compute "
167template <
typename ShapeFunctionDisplacement,
typename ShapeFunction,
168 int DisplacementDim,
typename ConstitutiveTraits>
171 ConstitutiveTraits>::
172 convertInitialStressType(
unsigned const ip,
177 double const p_at_ip)
179 bool constexpr is_strain_temperature_constitutive =
182 ConstitutiveTraits>::value;
183 if (is_strain_temperature_constitutive &&
190 if (!is_strain_temperature_constitutive &&
196 double const alpha_b =
198 .template value<double>(variables, x_position, t, 0.0 );
200 double const bishop =
202 .template value<double>(variables, x_position, t, 0.0 );
204 ConstitutiveTraits::ConstitutiveSetting::convertInitialStressType(
209template <
typename ShapeFunctionDisplacement,
typename ShapeFunction,
210 int DisplacementDim,
typename ConstitutiveTraits>
213 ConstitutiveTraits>::
214 assembleWithJacobian(
double const t,
double const dt,
215 std::vector<double>
const& local_x,
216 std::vector<double>
const& local_x_prev,
217 std::vector<double>& local_rhs_data,
218 std::vector<double>& local_Jac_data)
228 typename ConstitutiveTraits::ConstitutiveSetting constitutive_setting;
234 std::nullopt, this->
element_.getID(),
242 local_x, local_x_prev,
248 loc_mat += loc_mat_current_ip;
257template <
typename ShapeFunctionDisplacement,
typename ShapeFunction,
258 int DisplacementDim,
typename ConstitutiveTraits>
261 ConstitutiveTraits>::
268 loc_mat.storage_p_a_p =
269 loc_mat.storage_p_a_p.colwise().sum().eval().asDiagonal();
270 loc_mat.storage_p_a_S =
271 loc_mat.storage_p_a_S.colwise().sum().eval().asDiagonal();
272 loc_mat.storage_p_a_S_Jpp =
273 loc_mat.storage_p_a_S_Jpp.colwise().sum().eval().asDiagonal();
277template <
typename ShapeFunctionDisplacement,
typename ShapeFunction,
278 int DisplacementDim,
typename ConstitutiveTraits>
281 ConstitutiveTraits>::
282 addToLocalMatrixData(
284 std::vector<double>
const& local_x,
285 std::vector<double>
const& local_x_prev,
289 std::vector<double>& local_rhs_data,
290 std::vector<double>& local_Jac_data)
const
292 constexpr auto local_matrix_dim =
294 assert(local_x.size() == local_matrix_dim);
297 typename ShapeMatricesTypeDisplacement::template MatrixType<
298 local_matrix_dim, local_matrix_dim>>(
299 local_Jac_data, local_matrix_dim, local_matrix_dim);
303 template VectorType<local_matrix_dim>>(
304 local_rhs_data, local_matrix_dim);
306 local_Jac.noalias() = loc_mat.Jac;
307 local_rhs.noalias() = -loc_mat.res;
312 block_TT(local_Jac).noalias() += loc_mat.M_TT / dt + loc_mat.K_TT;
314 loc_mat.K_Tp + loc_mat.dK_TT_dp + loc_mat.M_Tp / dt;
316 block_pT(local_Jac).noalias() += loc_mat.M_pT / dt + loc_mat.K_pT;
318 loc_mat.K_pp + loc_mat.storage_p_a_p / dt + loc_mat.storage_p_a_S_Jpp;
319 block_pu(local_Jac).noalias() = loc_mat.M_pu / dt;
324 auto const [T, p_L, u] =
localDOF(local_x);
325 auto const [T_prev, p_L_prev, u_prev] =
localDOF(local_x_prev);
327 block_T(local_rhs).noalias() -= loc_mat.M_TT * (T - T_prev) / dt +
328 loc_mat.K_TT * T + loc_mat.K_Tp * p_L +
329 loc_mat.M_Tp * (p_L - p_L_prev) / dt;
330 block_p(local_rhs).noalias() -=
331 loc_mat.K_pp * p_L + loc_mat.K_pT * T +
332 (loc_mat.storage_p_a_p + loc_mat.storage_p_a_S) * (p_L - p_L_prev) /
334 loc_mat.M_pu * (u - u_prev) / dt + loc_mat.M_pT * (T - T_prev) / dt;
337template <
typename ShapeFunctionDisplacement,
typename ShapeFunction,
338 int DisplacementDim,
typename ConstitutiveTraits>
341 ConstitutiveTraits>::
342 assembleWithJacobianSingleIP(
343 double const t,
double const dt,
345 std::vector<double>
const& local_x,
346 std::vector<double>
const& local_x_prev,
349 ConstitutiveTraits>
::IpData const& ip_data,
350 typename ConstitutiveTraits::ConstitutiveSetting& CS,
355 typename ConstitutiveTraits::StatefulData& current_state,
356 typename ConstitutiveTraits::StatefulDataPrev
const& prev_state,
358 typename ConstitutiveTraits::OutputData& output_data)
const
360 auto const& N_u = ip_data.N_u;
361 auto const& dNdx_u = ip_data.dNdx_u;
364 auto const& N = ip_data.N_p;
365 auto const& dNdx = ip_data.dNdx_p;
369 ShapeFunctionDisplacement::NPOINTS,
374 typename ConstitutiveTraits::ConstitutiveData CD;
376 auto const [T, p_L, u] =
localDOF(local_x);
377 auto const [T_prev, p_L_prev, u_prev] =
localDOF(local_x_prev);
380 auto models = ConstitutiveTraits::createConstitutiveModels(
382 typename ConstitutiveTraits::ConstitutiveTempData tmp;
384 double const T_ip = N * T;
385 double const T_prev_ip = N * T_prev;
388 double const p_cap_ip = -N * p_L;
389 double const p_cap_prev_ip = -N * p_L_prev;
394 CS.eval(models, t, dt, x_position,
396 {T_ip, T_prev_ip, grad_T_ip},
397 {p_cap_ip, p_cap_prev_ip, grad_p_cap_ip},
398 eps, current_state, prev_state, mat_state, tmp, output_data,
403 NodalMatrix
const NTN = N.transpose() * N;
404 NodalMatrix
const dNTdN = dNdx.transpose() * dNdx;
409 DisplacementDim)>::identity2;
410 typename ShapeMatricesTypeDisplacement::template MatrixType<
412 BTI2N = B.transpose() * identity2 * N;
435 dNdx.transpose() * std::get<EqPData<DisplacementDim>>(CD).rhs_p_dNT_V;
441 static_cast<int>(this->
process_data_.apply_body_force_for_deformation) *
443 std::get<GravityData<DisplacementDim>>(CD).volumetric_body_force;
446 out.storage_p_a_p.noalias() =
447 std::get<EqPData<DisplacementDim>>(CD).storage_p_a_p_X_NTN * NTN;
448 out.storage_p_a_S.noalias() =
449 std::get<TRMStorageData>(CD).storage_p_a_S_X_NTN * NTN;
450 out.storage_p_a_S_Jpp.noalias() =
451 std::get<TRMStorageData>(CD).storage_p_a_S_Jpp_X_NTN * NTN;
455 std::get<EqTData<DisplacementDim>>(CD).M_TT_X_NTN * NTN;
457 std::get<TRMVaporDiffusionData<DisplacementDim>>(CD).M_Tp_X_NTN * NTN;
460 std::get<EqPData<DisplacementDim>>(CD).M_pT_X_NTN * NTN;
462 std::get<EqPData<DisplacementDim>>(CD).M_pu_X_BTI2N * BTI2N.transpose();
467 std::get<TRMHeatStorageAndFluxData<DisplacementDim>>(CD)
471 (std::get<EqTData<DisplacementDim>>(CD).K_TT_NT_V_dN.transpose() *
473 std::get<TRMVaporDiffusionData<DisplacementDim>>(CD).K_TT_X_dNTdN *
476 out.dK_TT_dp.noalias() =
478 (std::get<TRMHeatStorageAndFluxData<DisplacementDim>>(CD)
479 .K_Tp_NT_V_dN.transpose() *
481 std::get<TRMHeatStorageAndFluxData<DisplacementDim>>(CD).K_Tp_X_NTN *
485 std::get<ThermoOsmosisData<DisplacementDim>>(CD).K_Tp_Laplace *
491 dNdx.transpose() * std::get<EqPData<DisplacementDim>>(CD).K_pp_Laplace *
497 std::get<ThermoOsmosisData<DisplacementDim>>(CD).K_pT_Laplace * dNdx;
501 std::get<TRMVaporDiffusionData<DisplacementDim>>(CD).J_pT_X_dNTdN *
504 std::get<TRMStorageData>(CD).J_pp_X_NTN * NTN +
505 std::get<EqPData<DisplacementDim>>(CD).J_pp_X_BTI2NT_u_dot_N *
506 BTI2N.transpose() * (u - u_prev) / dt *
509 std::get<EqPData<DisplacementDim>>(CD).J_pp_dNT_V_N * N;
513 std::get<SolidMechanicsDataStateless<DisplacementDim>>(CD).J_uT_BT_K_N *
517 std::get<SolidMechanicsDataStateless<DisplacementDim>>(CD)
521 std::get<GravityData<DisplacementDim>>(CD).J_up_HT_V_N * N;
524 std::get<SolidMechanicsDataStateless<DisplacementDim>>(CD)
528 out *= ip_data.integration_weight;
531template <
typename ShapeFunctionDisplacement,
typename ShapeFunction,
532 int DisplacementDim,
typename ConstitutiveTraits>
535 ConstitutiveTraits>::
536 computeSecondaryVariableConcrete(
double const t,
double const dt,
537 Eigen::VectorXd
const& local_x,
538 Eigen::VectorXd
const& local_x_prev)
540 auto const T =
block_T(local_x);
541 auto const p_L =
block_p(local_x);
542 auto const u =
block_u(local_x);
544 auto const T_prev =
block_T(local_x_prev);
545 auto const p_L_prev =
block_p(local_x_prev);
547 auto const e_id = this->
element_.getID();
549 auto& medium = *process_data.media_map.getMedium(e_id);
551 unsigned const n_integration_points =
554 typename ConstitutiveTraits::ConstitutiveSetting constitutive_setting;
556 auto models = ConstitutiveTraits::createConstitutiveModels(
558 typename ConstitutiveTraits::ConstitutiveTempData tmp;
559 typename ConstitutiveTraits::ConstitutiveData CD;
561 for (
unsigned ip = 0; ip < n_integration_points; ip++)
569 auto const& N = ip_data.N_p;
570 auto const& N_u = ip_data.N_u;
571 auto const& dNdx_u = ip_data.dNdx_u;
572 auto const& dNdx = ip_data.dNdx_p;
575 std::nullopt, this->
element_.getID(),
586 ShapeFunctionDisplacement::NPOINTS,
590 double const T_ip = N * T;
591 double const T_prev_ip = N * T_prev;
594 double const p_cap_ip = -N * p_L;
595 double const p_cap_prev_ip = -N * p_L_prev;
600 constitutive_setting.eval(models,
603 {T_ip, T_prev_ip, grad_T_ip},
604 {p_cap_ip, p_cap_prev_ip, grad_p_cap_ip},
611 ShapeFunction,
typename ShapeFunctionDisplacement::MeshElement,
613 *process_data.pressure_interpolated);
615 ShapeFunction,
typename ShapeFunctionDisplacement::MeshElement,
616 DisplacementDim>(this->element_, this->is_axially_symmetric_, T,
617 *process_data.temperature_interpolated);
Property const & property(PropertyType const &p) const
double gas_phase_pressure
double capillary_pressure
double liquid_phase_pressure
constexpr double getWeight() const
MathLib::WeightedPoint const & getWeightedPoint(unsigned const igp) const
std::optional< MathLib::Point3d > const getCoordinates() const
MatrixType< _kelvin_vector_size, _number_of_dof > BMatrixType
static auto block_TT(auto &mat)
static constexpr auto & N_u_op
void addToLocalMatrixData(double const dt, std::vector< double > const &local_x, std::vector< double > const &local_x_prev, LocalMatrices const &loc_mat, std::vector< double > &local_rhs_data, std::vector< double > &local_Jac_data) const
static constexpr int temperature_size
static constexpr int temperature_index
static constexpr int pressure_size
void assembleWithJacobianSingleIP(double const t, double const dt, ParameterLib::SpatialPosition const &x_position, std::vector< double > const &local_x, std::vector< double > const &local_x_prev, IpData const &ip_data, typename ConstitutiveTraits::ConstitutiveSetting &CS, MaterialPropertyLib::Medium &medium, LocalMatrices &out, typename ConstitutiveTraits::StatefulData ¤t_state, typename ConstitutiveTraits::StatefulDataPrev const &prev_state, MaterialStateData< DisplacementDim > &mat_state, typename ConstitutiveTraits::OutputData &output_data) const
ShapeMatrixPolicyType< ShapeFunctionDisplacement, DisplacementDim > ShapeMatricesTypeDisplacement
static auto block_p(auto &vec)
typename ShapeMatricesType::GlobalDimVectorType GlobalDimVectorType
static auto block_pu(auto &mat)
std::vector< IpData > ip_data_
void convertInitialStressType(unsigned const ip, double const t, ParameterLib::SpatialPosition const x_position, MaterialPropertyLib::Medium const &medium, MPL::VariableArray const &variables, double const p_at_ip)
typename BMatricesType::KelvinVectorType KelvinVectorType
static auto block_up(auto &mat)
static auto block_Tp(auto &mat)
ThermoRichardsMechanicsLocalAssembler(ThermoRichardsMechanicsLocalAssembler const &)=delete
void setInitialConditionsConcrete(Eigen::VectorXd const local_x, double const t, int const process_id) override
ShapeMatrixPolicyType< ShapeFunction, DisplacementDim > ShapeMatricesType
static auto block_pp(auto &mat)
IntegrationPointData< ShapeMatricesTypeDisplacement, ShapeMatricesType, DisplacementDim, ShapeFunctionDisplacement::NPOINTS > IpData
static auto block_pT(auto &mat)
static auto block_uT(auto &mat)
static auto block_uu(auto &mat)
static constexpr auto localDOF(std::vector< double > const &x)
static constexpr int displacement_size
static constexpr int pressure_index
static auto block_u(auto &vec)
static auto block_T(auto &vec)
void massLumping(LocalMatrices &loc_mat) const
@ bishops_effective_stress
constexpr int kelvin_vector_dimensions(int const displacement_dim)
Kelvin vector dimensions for given displacement dimension.
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 shapeFunctionInterpolate(const NodalValues &, const ShapeMatrix &)
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)
std::array< double, 3 > interpolateCoordinates(MeshLib::Element const &e, typename ShapeMatricesType::ShapeMatrices::ShapeType const &N)
auto & get(Tuples &... ts)
BMatrixType computeBMatrix(DNDX_Type const &dNdx, N_Type const &N, const double radius, const bool is_axially_symmetric)
Fills a B-matrix based on given shape function dN/dx values.
MatrixType< ShapeFunction::NPOINTS, ShapeFunction::NPOINTS > NodalMatrixType
static Eigen::Matrix< double, KelvinVectorSize, 1 > const identity2
Kelvin mapping of 2nd order identity tensor.
bool const is_axially_symmetric_
NumLib::GenericIntegrationMethod const & integration_method_
std::vector< typename ConstitutiveTraits::OutputData > output_data_
ThermoRichardsMechanicsProcessData< DisplacementDim, ConstitutiveTraits > & process_data_
LocalAssemblerInterface(MeshLib::Element const &e, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, ThermoRichardsMechanicsProcessData< DisplacementDim, ConstitutiveTraits > &process_data)
std::vector< typename ConstitutiveTraits::StatefulDataPrev > prev_states_
std::vector< MaterialStateData< DisplacementDim > > material_states_
ConstitutiveTraits::SolidConstitutiveRelation const & solid_material_
std::vector< typename ConstitutiveTraits::StatefulData > current_states_
MeshLib::Element const & element_