6#include <spdlog/fmt/fmt.h>
29template <
typename ShapeFunctionDisplacement,
typename ShapeFunction,
30 int DisplacementDim,
typename ConstitutiveTraits>
32 DisplacementDim, ConstitutiveTraits>::
33 ThermoRichardsMechanicsLocalAssembler(
37 bool const is_axially_symmetric,
41 e, integration_method, is_axially_symmetric, process_data)
43 unsigned const n_integration_points =
46 ip_data_.resize(n_integration_points);
48 auto const shape_matrices_u =
51 DisplacementDim>(e, is_axially_symmetric,
54 auto const shape_matrices =
56 DisplacementDim>(e, is_axially_symmetric,
59 for (
unsigned ip = 0; ip < n_integration_points; ip++)
62 auto const& sm_u = shape_matrices_u[ip];
65 sm_u.integralMeasure * sm_u.detJ;
68 ip_data.dNdx_u = sm_u.dNdx;
71 ip_data.N_p = shape_matrices[ip].N;
72 ip_data.dNdx_p = shape_matrices[ip].dNdx;
76template <
typename ShapeFunctionDisplacement,
typename ShapeFunction,
77 int DisplacementDim,
typename ConstitutiveTraits>
85 assert(local_x.size() ==
88 auto const p_L = local_x.template segment<pressure_size>(
pressure_index);
92 constexpr double dt = std::numeric_limits<double>::quiet_NaN();
98 typename ConstitutiveTraits::ConstitutiveSetting
const constitutive_setting;
99 auto models = ConstitutiveTraits::createConstitutiveModels(
102 unsigned const n_integration_points =
104 for (
unsigned ip = 0; ip < n_integration_points; ip++)
110 std::nullopt, this->
element_.getID(),
135 .template value<double>(variables, x_position, t, dt);
136 std::get<PrevState<SaturationData>>(this->
prev_states_[ip])->S_L =
142 constitutive_setting.init(models, t, dt, x_position, media_data,
151 "The medium has no saturation property required to compute "
160template <
typename ShapeFunctionDisplacement,
typename ShapeFunction,
161 int DisplacementDim,
typename ConstitutiveTraits>
164 ConstitutiveTraits>::
165 convertInitialStressType(
unsigned const ip,
170 double const p_at_ip)
172 bool constexpr is_strain_temperature_constitutive =
175 ConstitutiveTraits>::value;
176 if (is_strain_temperature_constitutive &&
183 if (!is_strain_temperature_constitutive &&
189 double const alpha_b =
191 .template value<double>(variables, x_position, t, 0.0 );
193 double const bishop =
195 .template value<double>(variables, x_position, t, 0.0 );
197 ConstitutiveTraits::ConstitutiveSetting::convertInitialStressType(
202template <
typename ShapeFunctionDisplacement,
typename ShapeFunction,
203 int DisplacementDim,
typename ConstitutiveTraits>
206 ConstitutiveTraits>::
207 assembleWithJacobian(
double const t,
double const dt,
208 std::vector<double>
const& local_x,
209 std::vector<double>
const& local_x_prev,
210 std::vector<double>& local_rhs_data,
211 std::vector<double>& local_Jac_data)
221 typename ConstitutiveTraits::ConstitutiveSetting constitutive_setting;
227 std::nullopt, this->
element_.getID(),
235 local_x, local_x_prev,
241 loc_mat += loc_mat_current_ip;
250template <
typename ShapeFunctionDisplacement,
typename ShapeFunction,
251 int DisplacementDim,
typename ConstitutiveTraits>
254 ConstitutiveTraits>::
261 loc_mat.storage_p_a_p =
262 loc_mat.storage_p_a_p.colwise().sum().eval().asDiagonal();
263 loc_mat.storage_p_a_S =
264 loc_mat.storage_p_a_S.colwise().sum().eval().asDiagonal();
265 loc_mat.storage_p_a_S_Jpp =
266 loc_mat.storage_p_a_S_Jpp.colwise().sum().eval().asDiagonal();
270template <
typename ShapeFunctionDisplacement,
typename ShapeFunction,
271 int DisplacementDim,
typename ConstitutiveTraits>
274 ConstitutiveTraits>::
275 addToLocalMatrixData(
277 std::vector<double>
const& local_x,
278 std::vector<double>
const& local_x_prev,
282 std::vector<double>& local_rhs_data,
283 std::vector<double>& local_Jac_data)
const
285 constexpr auto local_matrix_dim =
287 assert(local_x.size() == local_matrix_dim);
290 typename ShapeMatricesTypeDisplacement::template MatrixType<
291 local_matrix_dim, local_matrix_dim>>(
292 local_Jac_data, local_matrix_dim, local_matrix_dim);
296 template VectorType<local_matrix_dim>>(
297 local_rhs_data, local_matrix_dim);
299 local_Jac.noalias() = loc_mat.Jac;
300 local_rhs.noalias() = -loc_mat.res;
305 block_TT(local_Jac).noalias() += loc_mat.M_TT / dt + loc_mat.K_TT;
307 loc_mat.K_Tp + loc_mat.dK_TT_dp + loc_mat.M_Tp / dt;
309 block_pT(local_Jac).noalias() += loc_mat.M_pT / dt + loc_mat.K_pT;
311 loc_mat.K_pp + loc_mat.storage_p_a_p / dt + loc_mat.storage_p_a_S_Jpp;
312 block_pu(local_Jac).noalias() = loc_mat.M_pu / dt;
317 auto const [T, p_L, u] =
localDOF(local_x);
318 auto const [T_prev, p_L_prev, u_prev] =
localDOF(local_x_prev);
320 block_T(local_rhs).noalias() -= loc_mat.M_TT * (T - T_prev) / dt +
321 loc_mat.K_TT * T + loc_mat.K_Tp * p_L +
322 loc_mat.M_Tp * (p_L - p_L_prev) / dt;
323 block_p(local_rhs).noalias() -=
324 loc_mat.K_pp * p_L + loc_mat.K_pT * T +
325 (loc_mat.storage_p_a_p + loc_mat.storage_p_a_S) * (p_L - p_L_prev) /
327 loc_mat.M_pu * (u - u_prev) / dt + loc_mat.M_pT * (T - T_prev) / dt;
330template <
typename ShapeFunctionDisplacement,
typename ShapeFunction,
331 int DisplacementDim,
typename ConstitutiveTraits>
334 ConstitutiveTraits>::
335 assembleWithJacobianSingleIP(
336 double const t,
double const dt,
338 std::vector<double>
const& local_x,
339 std::vector<double>
const& local_x_prev,
342 ConstitutiveTraits>
::IpData const& ip_data,
343 typename ConstitutiveTraits::ConstitutiveSetting& CS,
348 typename ConstitutiveTraits::StatefulData& current_state,
349 typename ConstitutiveTraits::StatefulDataPrev
const& prev_state,
351 typename ConstitutiveTraits::OutputData& output_data)
const
353 auto const& N_u = ip_data.N_u;
354 auto const& dNdx_u = ip_data.dNdx_u;
357 auto const& N = ip_data.N_p;
358 auto const& dNdx = ip_data.dNdx_p;
362 ShapeFunctionDisplacement::NPOINTS,
367 typename ConstitutiveTraits::ConstitutiveData CD;
369 auto const [T, p_L, u] =
localDOF(local_x);
370 auto const [T_prev, p_L_prev, u_prev] =
localDOF(local_x_prev);
373 auto models = ConstitutiveTraits::createConstitutiveModels(
375 typename ConstitutiveTraits::ConstitutiveTempData tmp;
377 double const T_ip = N * T;
378 double const T_prev_ip = N * T_prev;
381 double const p_cap_ip = -N * p_L;
382 double const p_cap_prev_ip = -N * p_L_prev;
387 CS.eval(models, t, dt, x_position,
389 {T_ip, T_prev_ip, grad_T_ip},
390 {p_cap_ip, p_cap_prev_ip, grad_p_cap_ip},
391 eps, current_state, prev_state, mat_state, tmp, output_data,
396 NodalMatrix
const NTN = N.transpose() * N;
397 NodalMatrix
const dNTdN = dNdx.transpose() * dNdx;
402 DisplacementDim)>::identity2;
403 typename ShapeMatricesTypeDisplacement::template MatrixType<
405 BTI2N = B.transpose() * identity2 * N;
428 dNdx.transpose() * std::get<EqPData<DisplacementDim>>(CD).rhs_p_dNT_V;
434 static_cast<int>(this->
process_data_.apply_body_force_for_deformation) *
436 std::get<GravityData<DisplacementDim>>(CD).volumetric_body_force;
439 out.storage_p_a_p.noalias() =
440 std::get<EqPData<DisplacementDim>>(CD).storage_p_a_p_X_NTN * NTN;
441 out.storage_p_a_S.noalias() =
442 std::get<TRMStorageData>(CD).storage_p_a_S_X_NTN * NTN;
443 out.storage_p_a_S_Jpp.noalias() =
444 std::get<TRMStorageData>(CD).storage_p_a_S_Jpp_X_NTN * NTN;
448 std::get<EqTData<DisplacementDim>>(CD).M_TT_X_NTN * NTN;
450 std::get<TRMVaporDiffusionData<DisplacementDim>>(CD).M_Tp_X_NTN * NTN;
453 std::get<EqPData<DisplacementDim>>(CD).M_pT_X_NTN * NTN;
455 std::get<EqPData<DisplacementDim>>(CD).M_pu_X_BTI2N * BTI2N.transpose();
460 std::get<TRMHeatStorageAndFluxData<DisplacementDim>>(CD)
464 (std::get<EqTData<DisplacementDim>>(CD).K_TT_NT_V_dN.transpose() *
466 std::get<TRMVaporDiffusionData<DisplacementDim>>(CD).K_TT_X_dNTdN *
469 out.dK_TT_dp.noalias() =
471 (std::get<TRMHeatStorageAndFluxData<DisplacementDim>>(CD)
472 .K_Tp_NT_V_dN.transpose() *
474 std::get<TRMHeatStorageAndFluxData<DisplacementDim>>(CD).K_Tp_X_NTN *
478 std::get<ThermoOsmosisData<DisplacementDim>>(CD).K_Tp_Laplace *
484 dNdx.transpose() * std::get<EqPData<DisplacementDim>>(CD).K_pp_Laplace *
490 std::get<ThermoOsmosisData<DisplacementDim>>(CD).K_pT_Laplace * dNdx;
494 std::get<TRMVaporDiffusionData<DisplacementDim>>(CD).J_pT_X_dNTdN *
497 std::get<TRMStorageData>(CD).J_pp_X_NTN * NTN +
498 std::get<EqPData<DisplacementDim>>(CD).J_pp_X_BTI2NT_u_dot_N *
499 BTI2N.transpose() * (u - u_prev) / dt *
502 std::get<EqPData<DisplacementDim>>(CD).J_pp_dNT_V_N * N;
506 std::get<SolidMechanicsDataStateless<DisplacementDim>>(CD).J_uT_BT_K_N *
510 std::get<SolidMechanicsDataStateless<DisplacementDim>>(CD)
514 std::get<GravityData<DisplacementDim>>(CD).J_up_HT_V_N * N;
517 std::get<SolidMechanicsDataStateless<DisplacementDim>>(CD)
521 out *= ip_data.integration_weight;
524template <
typename ShapeFunctionDisplacement,
typename ShapeFunction,
525 int DisplacementDim,
typename ConstitutiveTraits>
528 ConstitutiveTraits>::
529 computeSecondaryVariableConcrete(
double const t,
double const dt,
530 Eigen::VectorXd
const& local_x,
531 Eigen::VectorXd
const& local_x_prev)
533 auto const T =
block_T(local_x);
534 auto const p_L =
block_p(local_x);
535 auto const u =
block_u(local_x);
537 auto const T_prev =
block_T(local_x_prev);
538 auto const p_L_prev =
block_p(local_x_prev);
540 auto const e_id = this->
element_.getID();
542 auto& medium = *process_data.media_map.getMedium(e_id);
544 unsigned const n_integration_points =
547 typename ConstitutiveTraits::ConstitutiveSetting constitutive_setting;
549 auto models = ConstitutiveTraits::createConstitutiveModels(
551 typename ConstitutiveTraits::ConstitutiveTempData tmp;
552 typename ConstitutiveTraits::ConstitutiveData CD;
554 for (
unsigned ip = 0; ip < n_integration_points; ip++)
562 auto const& N = ip_data.N_p;
563 auto const& N_u = ip_data.N_u;
564 auto const& dNdx_u = ip_data.dNdx_u;
565 auto const& dNdx = ip_data.dNdx_p;
568 std::nullopt, this->
element_.getID(),
578 ShapeFunctionDisplacement::NPOINTS,
582 double const T_ip = N * T;
583 double const T_prev_ip = N * T_prev;
586 double const p_cap_ip = -N * p_L;
587 double const p_cap_prev_ip = -N * p_L_prev;
592 constitutive_setting.eval(models,
595 {T_ip, T_prev_ip, grad_T_ip},
596 {p_cap_ip, p_cap_prev_ip, grad_p_cap_ip},
603 ShapeFunction,
typename ShapeFunctionDisplacement::MeshElement,
605 *process_data.pressure_interpolated);
607 ShapeFunction,
typename ShapeFunctionDisplacement::MeshElement,
608 DisplacementDim>(this->element_, this->is_axially_symmetric_, T,
609 *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 &)
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_