31namespace ThermoHydroMechanics
33template <
typename ShapeFunctionDisplacement,
typename ShapeFunctionPressure,
35ThermoHydroMechanicsLocalAssembler<ShapeFunctionDisplacement,
36 ShapeFunctionPressure, DisplacementDim>::
37 ThermoHydroMechanicsLocalAssembler(
41 bool const is_axially_symmetric,
43 : _process_data(process_data),
44 _integration_method(integration_method),
46 _is_axially_symmetric(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 !=
114 static_cast<int>(_integration_method.getIntegrationOrder()))
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.",
123 if (name ==
"sigma_ip")
125 if (_process_data.initial_stress !=
nullptr)
128 "Setting initial conditions for stress from integration "
129 "point data and from a parameter '{:s}' is not possible "
131 _process_data.initial_stress->name);
134 return ProcessLib::setIntegrationPointKelvinVectorData<DisplacementDim>(
135 values, _ip_data, &IpData::sigma_eff);
137 if (name ==
"epsilon_ip")
139 return ProcessLib::setIntegrationPointKelvinVectorData<DisplacementDim>(
140 values, _ip_data, &IpData::eps);
146template <
typename ShapeFunctionDisplacement,
typename ShapeFunctionPressure,
149 ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim>::
150 updateConstitutiveRelations(
151 Eigen::Ref<Eigen::VectorXd const>
const local_x,
152 Eigen::Ref<Eigen::VectorXd const>
const local_x_prev,
154 double const dt,
IpData& ip_data,
157 assert(local_x.size() ==
158 pressure_size + displacement_size + temperature_size);
160 auto const [T, p, u] = localDOF(local_x);
161 auto const [T_prev, p_prev, u_prev] = localDOF(local_x_prev);
163 auto const& solid_material =
165 _process_data.solid_materials, _process_data.material_ids,
168 auto const& medium = _process_data.media_map.getMedium(_element.getID());
169 auto const& liquid_phase = medium->phase(
"AqueousLiquid");
170 auto const& solid_phase = medium->phase(
"Solid");
171 auto*
const frozen_liquid_phase = medium->hasPhase(
"FrozenLiquid")
172 ? &medium->phase(
"FrozenLiquid")
176 auto const& N_u = ip_data.
N_u;
177 auto const& dNdx_u = ip_data.
dNdx_u;
179 auto const& N = ip_data.
N;
180 auto const& dNdx = ip_data.
dNdx;
182 auto const T_int_pt = N.dot(T);
183 auto const T_prev_int_pt = N.dot(T_prev);
184 double const dT_int_pt = T_int_pt - T_prev_int_pt;
192 ShapeFunctionDisplacement::NPOINTS,
194 dNdx_u, N_u, x_coord, _is_axially_symmetric);
198 auto& eps = ip_data.
eps;
199 eps.noalias() = B * u;
202 double const p_int_pt = N.dot(p);
207 auto const solid_density =
209 .template value<double>(vars, x_position, t, dt);
211 auto const porosity =
213 .template value<double>(vars, x_position, t, dt);
219 .template value<double>(vars, x_position, t, dt);
220 auto const& alpha = crv.alpha_biot;
223 t, x_position, dt,
static_cast<double>(T_int_pt));
224 auto const solid_skeleton_compressibility =
225 1 / solid_material.getBulkModulus(t, x_position, &C_el);
227 crv.beta_SR = (1 - alpha) * solid_skeleton_compressibility;
232 auto const& identity2 = Invariants::identity2;
233 auto const sigma_total =
234 (ip_data.
sigma_eff - alpha * p_int_pt * identity2).eval();
243 auto const intrinsic_permeability =
244 MaterialPropertyLib::formEigenTensor<DisplacementDim>(
246 .value(vars, x_position, t, dt));
248 auto const fluid_density =
250 .template value<double>(vars, x_position, t, dt);
256 .template dValue<double>(
260 crv.fluid_compressibility = 1 / fluid_density * drho_dp;
262 double const fluid_volumetric_thermal_expansion_coefficient =
264 liquid_phase, vars, fluid_density, x_position, t, dt);
269 .template value<double>(vars, x_position, t, dt);
270 crv.K_over_mu = intrinsic_permeability / ip_data_output.
viscosity;
272 auto const& b = _process_data.specific_body_force;
275 crv.solid_linear_thermal_expansion_coefficient =
276 MPL::formKelvinVector<DisplacementDim>(
280 .value(vars, x_position, t, dt));
284 crv.solid_linear_thermal_expansion_coefficient * dT_int_pt;
286 crv.K_pT_thermal_osmosis =
287 (solid_phase.hasProperty(
289 ? MaterialPropertyLib::formEigenTensor<DisplacementDim>(
292 thermal_osmosis_coefficient)
293 .value(vars, x_position, t, dt))
294 : Eigen::MatrixXd::Zero(DisplacementDim, DisplacementDim));
297 crv.K_pT_thermal_osmosis * dNdx * T +
298 crv.K_over_mu * fluid_density * b;
305 auto& eps_m = ip_data.
eps_m;
307 eps_m.noalias() = eps_m_prev + eps - eps_prev - dthermal_strain;
315 crv.rho = solid_density * (1 - porosity) + porosity * fluid_density;
318 porosity * fluid_volumetric_thermal_expansion_coefficient +
320 Invariants::trace(crv.solid_linear_thermal_expansion_coefficient);
333 .template value<double>(vars, x_position, t, dt);
334 crv.effective_thermal_conductivity =
335 MaterialPropertyLib::formEigenTensor<DisplacementDim>(
339 .value(vars, x_position, t, dt));
343 crv.effective_thermal_conductivity.noalias() +=
344 fluid_density * crv.c_f *
346 _process_data.stabilizer, _element.getID(),
347 GlobalDimMatrixType::Zero(DisplacementDim, DisplacementDim),
354 .template value<double>(vars, x_position, t, dt);
357 crv.effective_volumetric_heat_capacity =
358 porosity * fluid_density * crv.c_f +
359 (1.0 - porosity) * solid_density * c_s;
361 if (frozen_liquid_phase)
364 double const phi_fr =
366 .template value<double>(vars, x_position, t, dt);
369 auto const frozen_liquid_value =
372 return (*frozen_liquid_phase)[p].template value<double>(
373 vars, x_position, t, dt);
376 double const rho_fr =
379 double const c_fr = frozen_liquid_value(
382 double const l_fr = frozen_liquid_value(
385 double const dphi_fr_dT =
387 .template dValue<double>(
391 double const phi_fr_prev = [&]()
396 .template value<double>(vars_prev, x_position, t, dt);
402 DisplacementDim>
const ice_linear_thermal_expansion_coefficient =
403 MPL::formKelvinVector<DisplacementDim>(
407 .value(vars, x_position, t, dt));
410 dthermal_strain_ice =
411 ice_linear_thermal_expansion_coefficient * dT_int_pt;
416 phase_change_expansion_coefficient =
417 MPL::formKelvinVector<DisplacementDim>(
420 phase_change_expansivity)
421 .value(vars, x_position, t, dt));
424 dphase_change_strain = phase_change_expansion_coefficient *
425 (phi_fr - phi_fr_prev) / porosity;
429 auto& eps0 = ip_data.
eps0;
430 auto const& eps0_prev = ip_data.
eps0_prev;
436 eps_m_ice.noalias() = eps_m_ice_prev + eps - eps_prev -
437 (eps0 - eps0_prev) - dthermal_strain_ice -
438 dphase_change_strain;
444 *_process_data.ice_constitutive_relation, vars_ice, t, x_position,
446 crv.effective_volumetric_heat_capacity +=
447 -phi_fr * fluid_density * crv.c_f + phi_fr * rho_fr * c_fr -
448 l_fr * rho_fr * dphi_fr_dT;
451 double const d2phi_fr_dT2 =
453 .template d2Value<double>(
458 crv.J_uu_fr = phi_fr * C_IR;
461 crv.r_u_fr = phi_fr * sigma_eff_ice;
463 crv.J_uT_fr = phi_fr * C_IR * ice_linear_thermal_expansion_coefficient;
465 crv.J_TT_fr = ((rho_fr * c_fr - fluid_density * crv.c_f) * dphi_fr_dT +
466 l_fr * rho_fr * d2phi_fr_dT2) *
474template <
typename ShapeFunctionDisplacement,
typename ShapeFunctionPressure,
477 ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim>::
478 assembleWithJacobian(
double const t,
double const dt,
479 std::vector<double>
const& local_x,
480 std::vector<double>
const& local_x_prev,
481 std::vector<double>& ,
482 std::vector<double>& ,
483 std::vector<double>& local_rhs_data,
484 std::vector<double>& local_Jac_data)
486 assert(local_x.size() ==
487 pressure_size + displacement_size + temperature_size);
490 Eigen::Map<Eigen::VectorXd const>(local_x.data(), local_x.size());
491 auto const x_prev = Eigen::Map<Eigen::VectorXd const>(local_x_prev.data(),
492 local_x_prev.size());
494 auto const [T, p, u] = localDOF(local_x);
495 auto const [T_prev, p_prev, u_prev] = localDOF(local_x_prev);
498 typename ShapeMatricesTypeDisplacement::template MatrixType<
499 temperature_size + displacement_size + pressure_size,
500 temperature_size + displacement_size + pressure_size>>(
501 local_Jac_data, displacement_size + pressure_size + temperature_size,
502 displacement_size + pressure_size + temperature_size);
505 typename ShapeMatricesTypeDisplacement::template VectorType<
506 displacement_size + pressure_size + temperature_size>>(
507 local_rhs_data, displacement_size + pressure_size + temperature_size);
510 MTT.setZero(temperature_size, temperature_size);
513 KTT.setZero(temperature_size, temperature_size);
516 KTp.setZero(temperature_size, pressure_size);
519 dKTT_dp.setZero(temperature_size, pressure_size);
522 laplace_p.setZero(pressure_size, pressure_size);
525 laplace_T.setZero(pressure_size, temperature_size);
528 storage_p.setZero(pressure_size, pressure_size);
531 storage_T.setZero(pressure_size, temperature_size);
533 typename ShapeMatricesTypeDisplacement::template MatrixType<
534 displacement_size, pressure_size>
536 Kup.setZero(displacement_size, pressure_size);
538 auto const& medium = _process_data.media_map.getMedium(_element.getID());
539 bool const has_frozen_liquid_phase = medium->hasPhase(
"FrozenLiquid");
541 unsigned const n_integration_points =
542 _integration_method.getNumberOfPoints();
544 std::vector<GlobalDimVectorType> ip_flux_vector;
545 double average_velocity_norm = 0.0;
546 ip_flux_vector.reserve(n_integration_points);
548 for (
unsigned ip = 0; ip < n_integration_points; ip++)
550 auto const& N_u = _ip_data[ip].N_u;
552 std::nullopt, _element.getID(), ip,
558 auto const crv = updateConstitutiveRelations(
559 x, x_prev, x_position, t, dt, _ip_data[ip], _ip_data_output[ip]);
561 auto const& w = _ip_data[ip].integration_weight;
563 auto const& dNdx_u = _ip_data[ip].dNdx_u;
565 auto const& N = _ip_data[ip].N;
566 auto const& dNdx = _ip_data[ip].dNdx;
568 auto const T_int_pt = N.dot(T);
576 ShapeFunctionDisplacement::NPOINTS,
578 dNdx_u, N_u, x_coord, _is_axially_symmetric);
580 auto const& b = _process_data.specific_body_force;
581 auto const velocity = _ip_data_output[ip].velocity;
587 if (has_frozen_liquid_phase)
590 .template block<displacement_size, displacement_size>(
591 displacement_index, displacement_index)
592 .noalias() += B.transpose() * crv.J_uu_fr * B * w;
594 local_rhs.template segment<displacement_size>(displacement_index)
595 .noalias() -= B.transpose() * crv.r_u_fr * w;
598 .template block<displacement_size, temperature_size>(
599 displacement_index, temperature_index)
600 .noalias() -= B.transpose() * crv.J_uT_fr * N * w;
604 .template block<displacement_size, displacement_size>(
605 displacement_index, displacement_index)
606 .noalias() += B.transpose() * crv.C * B * w;
609 .template block<displacement_size, temperature_size>(
610 displacement_index, temperature_index)
611 .noalias() -= B.transpose() * crv.C *
612 crv.solid_linear_thermal_expansion_coefficient * N *
615 local_rhs.template segment<displacement_size>(displacement_index)
616 .noalias() -= (B.transpose() * _ip_data[ip].sigma_eff -
617 N_u_op(N_u).transpose() * crv.rho * b) *
624 B.transpose() * crv.alpha_biot * Invariants::identity2 * N * w;
629 laplace_p.noalias() += dNdx.transpose() * crv.K_over_mu * dNdx * w;
631 storage_p.noalias() +=
633 (_ip_data[ip].porosity * crv.fluid_compressibility +
634 (crv.alpha_biot - _ip_data[ip].porosity) * crv.beta_SR) *
637 laplace_T.noalias() +=
638 dNdx.transpose() * crv.K_pT_thermal_osmosis * dNdx * w;
642 double const fluid_density = _ip_data_output[ip].fluid_density;
643 local_rhs.template segment<pressure_size>(pressure_index).noalias() +=
644 dNdx.transpose() * fluid_density * crv.K_over_mu * b * w;
648 storage_T.noalias() += N.transpose() * crv.beta * N * w;
659 dNdx.transpose() * crv.effective_thermal_conductivity * dNdx * w;
661 ip_flux_vector.emplace_back(velocity * fluid_density * crv.c_f);
662 average_velocity_norm += velocity.norm();
664 if (has_frozen_liquid_phase)
667 .template block<temperature_size, temperature_size>(
668 temperature_index, temperature_index)
669 .noalias() -= N.transpose() * crv.J_TT_fr * N * w;
673 N.transpose() * crv.effective_volumetric_heat_capacity * N * w;
679 dNdx.transpose() * T_int_pt * crv.K_pT_thermal_osmosis * dNdx * w;
682 dKTT_dp.noalias() -= fluid_density * crv.c_f * N.transpose() *
683 (dNdx * T).transpose() * crv.K_over_mu * dNdx * w;
729 _process_data.stabilizer, _ip_data, ip_flux_vector,
730 average_velocity_norm /
static_cast<double>(n_integration_points), KTT);
734 .template block<temperature_size, temperature_size>(temperature_index,
736 .noalias() += KTT + MTT / dt;
740 .template block<temperature_size, pressure_size>(temperature_index,
742 .noalias() += KTp + dKTT_dp;
746 .template block<displacement_size, pressure_size>(displacement_index,
752 .template block<pressure_size, temperature_size>(pressure_index,
754 .noalias() -= storage_T / dt - laplace_T;
758 .template block<pressure_size, pressure_size>(pressure_index,
760 .noalias() += laplace_p + storage_p / dt;
764 .template block<pressure_size, displacement_size>(pressure_index,
766 .noalias() += Kup.transpose() / dt;
769 local_rhs.template segment<pressure_size>(pressure_index).noalias() -=
770 laplace_p * p + laplace_T * T + storage_p * (p - p_prev) / dt -
771 storage_T * (T - T_prev) / dt + Kup.transpose() * (u - u_prev) / dt;
774 local_rhs.template segment<displacement_size>(displacement_index)
775 .noalias() += Kup * p;
778 local_rhs.template segment<temperature_size>(temperature_index).noalias() -=
779 KTT * T + MTT * (T - T_prev) / dt;
781 local_rhs.template segment<temperature_size>(temperature_index).noalias() -=
785template <
typename ShapeFunctionDisplacement,
typename ShapeFunctionPressure,
788 ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim>::
789 getIntPtDarcyVelocity(
791 std::vector<GlobalVector*>
const& ,
792 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& ,
793 std::vector<double>& cache)
const
795 unsigned const n_integration_points =
796 _integration_method.getNumberOfPoints();
800 double, DisplacementDim, Eigen::Dynamic, Eigen::RowMajor>>(
801 cache, DisplacementDim, n_integration_points);
803 for (
unsigned ip = 0; ip < n_integration_points; ip++)
805 cache_matrix.col(ip).noalias() = _ip_data_output[ip].velocity;
811template <
typename ShapeFunctionDisplacement,
typename ShapeFunction,
814 ShapeFunctionDisplacement, ShapeFunction, DisplacementDim>::
815 getIntPtFluidDensity(
817 std::vector<GlobalVector*>
const& ,
818 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& ,
819 std::vector<double>& cache)
const
826template <
typename ShapeFunctionDisplacement,
typename ShapeFunction,
829 ShapeFunctionDisplacement, ShapeFunction, DisplacementDim>::
832 std::vector<GlobalVector*>
const& ,
833 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& ,
834 std::vector<double>& cache)
const
841template <
typename ShapeFunctionDisplacement,
typename ShapeFunctionPressure,
844 ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim>::
845 computeSecondaryVariableConcrete(
double const t,
double const dt,
846 Eigen::VectorXd
const& local_x,
847 Eigen::VectorXd
const& local_x_prev)
849 auto const p = local_x.template segment<pressure_size>(pressure_index);
851 local_x.template segment<temperature_size>(temperature_index);
853 unsigned const n_integration_points =
854 _integration_method.getNumberOfPoints();
856 double fluid_density_avg = 0;
857 double viscosity_avg = 0;
860 KV sigma_avg = KV::Zero();
862 for (
unsigned ip = 0; ip < n_integration_points; ip++)
864 auto& ip_data = _ip_data[ip];
865 auto const& N_u = ip_data.N_u;
868 std::nullopt, _element.getID(), ip,
874 updateConstitutiveRelations(local_x, local_x_prev, x_position, t, dt,
875 _ip_data[ip], _ip_data_output[ip]);
877 fluid_density_avg += _ip_data_output[ip].fluid_density;
878 viscosity_avg += _ip_data_output[ip].viscosity;
879 sigma_avg += ip_data.sigma_eff;
882 fluid_density_avg /= n_integration_points;
883 viscosity_avg /= n_integration_points;
884 sigma_avg /= n_integration_points;
886 (*_process_data.element_fluid_density)[_element.getID()] =
888 (*_process_data.element_viscosity)[_element.getID()] = viscosity_avg;
890 Eigen::Map<KV>(&(*_process_data.element_stresses)[_element.getID() *
891 KV::RowsAtCompileTime]) =
895 ShapeFunctionPressure,
typename ShapeFunctionDisplacement::MeshElement,
896 DisplacementDim>(_element, _is_axially_symmetric, p,
897 *_process_data.pressure_interpolated);
900 ShapeFunctionPressure,
typename ShapeFunctionDisplacement::MeshElement,
901 DisplacementDim>(_element, _is_axially_symmetric, T,
902 *_process_data.temperature_interpolated);
std::variant< std::monostate, Eigen::Matrix< double, 4, 1 >, Eigen::Matrix< double, 6, 1 > > mechanical_strain
double equivalent_plastic_strain
std::variant< std::monostate, Eigen::Matrix< double, 4, 1 >, Eigen::Matrix< double, 6, 1 > > total_stress
virtual std::size_t getID() const final
Returns the ID of the element.
MathLib::WeightedPoint const & getWeightedPoint(unsigned const igp) const
unsigned getNumberOfPoints() const
MatrixType< _kelvin_vector_size, _number_of_dof > BMatrixType
NumLib::GenericIntegrationMethod const & _integration_method
Eigen::Matrix< double, KelvinVectorSize, 1 > SymmetricTensor
ThermoHydroMechanicsProcessData< DisplacementDim > & _process_data
std::vector< IntegrationPointDataForOutput< DisplacementDim >, Eigen::aligned_allocator< IntegrationPointDataForOutput< DisplacementDim > > > _ip_data_output
ShapeMatrixPolicyType< ShapeFunctionPressure, DisplacementDim > ShapeMatricesTypePressure
SecondaryData< typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType > _secondary_data
MeshLib::Element const & _element
std::vector< IpData, Eigen::aligned_allocator< IpData > > _ip_data
typename ShapeMatricesTypePressure::GlobalDimVectorType GlobalDimVectorType
ShapeMatrixPolicyType< ShapeFunctionDisplacement, DisplacementDim > ShapeMatricesTypeDisplacement
auto & selectSolidConstitutiveRelation(SolidMaterialsMap const &constitutive_relations, MeshLib::PropertyVector< int > const *const material_ids, std::size_t const element_id)
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, kelvin_vector_dimensions(DisplacementDim), 1, Eigen::ColMajor > KelvinVectorType
Eigen::Matrix< double, 4, 1 > kelvinVectorToSymmetricTensor(Eigen::Matrix< double, 4, 1, Eigen::ColMajor, 4, 1 > const &v)
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, 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)
void setIPDataInitialConditions(std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > const &_integration_point_writer, MeshLib::Properties const &mesh_properties, LocalAssemblersVector &local_assemblers, bool const remove_name_suffix=false)
MatrixType< ShapeFunction::NPOINTS, ShapeFunction::NPOINTS > NodalMatrixType
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N_u
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
BMatricesType::KelvinVectorType eps_prev