OGS
ThermoHydroMechanicsFEM-impl.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) OpenGeoSys Community (opengeosys.org)
2// SPDX-License-Identifier: BSD-3-Clause
3
4#pragma once
5
6#include <typeinfo>
7
22
23namespace ProcessLib
24{
26{
27template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
28 int DisplacementDim>
29ThermoHydroMechanicsLocalAssembler<ShapeFunctionDisplacement,
30 ShapeFunctionPressure, DisplacementDim>::
31 ThermoHydroMechanicsLocalAssembler(
32 MeshLib::Element const& e,
33 std::size_t const /*local_matrix_size*/,
34 NumLib::GenericIntegrationMethod const& integration_method,
35 bool const is_axially_symmetric,
37 : _process_data(process_data),
38 _integration_method(integration_method),
39 _element(e),
40 _is_axially_symmetric(is_axially_symmetric)
41{
42 unsigned const n_integration_points =
43 _integration_method.getNumberOfPoints();
44
45 _ip_data.reserve(n_integration_points);
46 _ip_data_output.resize(n_integration_points);
47 _secondary_data.N_u.resize(n_integration_points);
48
49 auto const shape_matrices_u =
50 NumLib::initShapeMatrices<ShapeFunctionDisplacement,
52 DisplacementDim>(e, is_axially_symmetric,
54
55 auto const shape_matrices_p =
56 NumLib::initShapeMatrices<ShapeFunctionPressure,
57 ShapeMatricesTypePressure, DisplacementDim>(
58 e, is_axially_symmetric, _integration_method);
59
60 auto const& solid_material =
62 _process_data.solid_materials, _process_data.material_ids,
63 e.getID());
64
65 // Consistency check: if frozen liquid phase is given, then the constitutive
66 // relation for ice must also be given, and vice versa.
67 auto const& medium = _process_data.media_map.getMedium(_element.getID());
69 (_process_data.ice_constitutive_relation != nullptr))
70 {
72 "Frozen liquid phase is {:s} and the solid material constitutive "
73 "relation for ice is {:s}. But both must be given (or both "
74 "omitted).",
76 ? "specified"
77 : "not specified",
78 _process_data.ice_constitutive_relation != nullptr
79 ? "specified"
80 : "not specified");
81 }
82 for (unsigned ip = 0; ip < n_integration_points; ip++)
83 {
84 _ip_data.emplace_back(solid_material);
85 auto& ip_data = _ip_data[ip];
86 auto const& sm_u = shape_matrices_u[ip];
87 ip_data.integration_weight =
88 _integration_method.getWeightedPoint(ip).getWeight() *
89 sm_u.integralMeasure * sm_u.detJ;
90
91 ip_data.N_u = sm_u.N;
92 ip_data.dNdx_u = sm_u.dNdx;
93
94 ip_data.N = shape_matrices_p[ip].N;
95 ip_data.dNdx = shape_matrices_p[ip].dNdx;
96
97 _secondary_data.N_u[ip] = shape_matrices_u[ip].N;
98 }
99}
100
101template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
102 int DisplacementDim>
104 ShapeFunctionDisplacement, ShapeFunctionPressure,
105 DisplacementDim>::setIPDataInitialConditions(std::string_view name,
106 double const* values,
107 int const integration_order)
108{
109 if (integration_order !=
110 static_cast<int>(_integration_method.getIntegrationOrder()))
111 {
112 OGS_FATAL(
113 "Setting integration point initial conditions; The integration "
114 "order of the local assembler for element {:d} is different from "
115 "the integration order in the initial condition.",
116 _element.getID());
117 }
118
119 if (name == "sigma")
120 {
121 if (_process_data.initial_stress.value)
122 {
123 OGS_FATAL(
124 "Setting initial conditions for stress from integration "
125 "point data and from a parameter '{:s}' is not possible "
126 "simultaneously.",
127 _process_data.initial_stress.value->name);
128 }
129
131 values, _ip_data, &IpData::sigma_eff);
132 }
133 if (name == "epsilon_m")
134 {
136 values, _ip_data, &IpData::eps_m);
137 }
138 if (name == "epsilon")
139 {
141 values, _ip_data, &IpData::eps);
142 }
143 if (name.starts_with("material_state_variable_"))
144 {
145 name.remove_prefix(24);
146
147 // Using first ip data for solid material. TODO (naumov) move solid
148 // material into element, store only material state in IPs.
149 auto const& internal_variables =
150 _ip_data[0].solid_material.getInternalVariables();
151 if (auto const iv = std::find_if(
152 begin(internal_variables), end(internal_variables),
153 [&name](auto const& iv) { return iv.name == name; });
154 iv != end(internal_variables))
155 {
156 DBUG("Setting material state variable '{:s}'", name);
159 iv->reference);
160 }
161
162 int const element_id = _element.getID();
163 DBUG(
164 "The solid material of element {:d} (material ID {:d}) does not "
165 "have an internal state variable called {:s}.",
166 element_id, (*_process_data.material_ids)[element_id], name);
167 }
168
169 return 0;
170}
171template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
172 int DisplacementDim>
174 ShapeFunctionDisplacement, ShapeFunctionPressure,
175 DisplacementDim>::setInitialConditionsConcrete(Eigen::VectorXd const
176 local_x,
177 double const t,
178 int const /*process_id*/)
179{
180 // TODO: For staggered scheme, overload
181 // LocalAssemblerInterface::setInitialConditions to enable local_x contains
182 // the primary variables from all coupled processes.
183 auto const [T, p, u] = localDOF(local_x);
184
185 double const dt = 0.0;
186
188 auto const& medium = _process_data.media_map.getMedium(_element.getID());
189 auto* const frozen_liquid_phase =
192 : nullptr;
193
194 int const n_integration_points = _integration_method.getNumberOfPoints();
195 for (int ip = 0; ip < n_integration_points; ip++)
196 {
197 auto& ip_data = _ip_data[ip];
198 auto const& N = ip_data.N;
199 auto const& N_u = ip_data.N_u;
200 ParameterLib::SpatialPosition const x_position{
201 std::nullopt, _element.getID(),
203 NumLib::interpolateCoordinates<ShapeFunctionDisplacement,
205 _element, N_u))};
206
207 auto& sigma_eff = ip_data.sigma_eff;
208 if (_process_data.initial_stress.isTotalStress())
209 {
210 auto const alpha_b =
212 .template value<double>(vars, x_position, t, dt);
213
214 sigma_eff.noalias() += alpha_b * N.dot(p) * Invariants::identity2;
215 }
216 ip_data.sigma_eff_prev.noalias() = sigma_eff;
217
218 vars.temperature = N.dot(T);
219 if (frozen_liquid_phase)
220 {
221 ip_data.phi_fr =
223 .template value<double>(vars, x_position, t, dt);
224 }
225 }
226}
227
228template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
229 int DisplacementDim>
231 ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim>::
232 updateConstitutiveRelations(
233 Eigen::Ref<Eigen::VectorXd const> const local_x,
234 Eigen::Ref<Eigen::VectorXd const> const local_x_prev,
235 ParameterLib::SpatialPosition const& x_position, double const t,
236 double const dt, IpData& ip_data,
238{
239 assert(local_x.size() ==
241
242 auto const [T, p, u] = localDOF(local_x);
243 auto const [T_prev, p_prev, u_prev] = localDOF(local_x_prev);
244
245 auto const& solid_material =
247 _process_data.solid_materials, _process_data.material_ids,
248 _element.getID());
249
250 auto const& medium = _process_data.media_map.getMedium(_element.getID());
251 auto const& liquid_phase =
253 auto const& solid_phase =
255 auto* const frozen_liquid_phase =
258 : nullptr;
260
261 auto const& N_u = ip_data.N_u;
262 auto const& dNdx_u = ip_data.dNdx_u;
263
264 auto const& N = ip_data.N;
265 auto const& dNdx = ip_data.dNdx;
266
267 auto const T_int_pt = N.dot(T);
268 auto const T_prev_int_pt = N.dot(T_prev);
269 double const dT_int_pt = T_int_pt - T_prev_int_pt;
270
271 auto const x_coord =
272 x_position.getCoordinates().value()[0]; // r for axisymmetry
273 auto const B =
274 LinearBMatrix::computeBMatrix<DisplacementDim,
275 ShapeFunctionDisplacement::NPOINTS,
277 dNdx_u, N_u, x_coord, _is_axially_symmetric);
278
280
281 auto& eps = ip_data.eps;
282 eps.noalias() = B * u;
284 B * u_prev;
285
286 vars.temperature = T_int_pt;
287 double const p_int_pt = N.dot(p);
288 vars.liquid_phase_pressure = p_int_pt;
289 double const p_prev_int_pt = N.dot(p_prev);
290 double const dp_int_pt = p_int_pt - p_prev_int_pt;
291 vars.liquid_saturation = 1.0;
292
293 auto const solid_density =
295 .template value<double>(vars, x_position, t, dt);
296
297 auto const drho_SR_dT =
299 .template dValue<double>(vars,
301 x_position, t, dt);
302
303 auto const porosity =
305 .template value<double>(vars, x_position, t, dt);
306 vars.porosity = porosity;
307 ip_data.porosity = porosity;
308
309 crv.alpha_biot =
311 .template value<double>(vars, x_position, t, dt);
312 auto const& alpha = crv.alpha_biot;
313
314 auto const C_el = ip_data.computeElasticTangentStiffness(
315 t, x_position, dt, static_cast<double>(T_int_pt));
316 auto const solid_skeleton_compressibility =
317 1 / solid_material.getBulkModulus(t, x_position, &C_el);
318
319 crv.beta_SR = (1 - alpha) * solid_skeleton_compressibility;
320
321 // Set mechanical variables for the intrinsic permeability model
322 // For stress dependent permeability.
323 {
324 auto const& identity2 = Invariants::identity2;
325 auto const sigma_total =
326 (ip_data.sigma_eff - alpha * p_int_pt * identity2).eval();
327 vars.total_stress.emplace<SymmetricTensor>(
329 }
330 // For strain dependent permeability
333 ip_data.material_state_variables->getEquivalentPlasticStrain();
334
335 auto const intrinsic_permeability =
338 .value(vars, x_position, t, dt));
339
340 auto const fluid_density =
341 liquid_phase.property(MaterialPropertyLib::PropertyType::density)
342 .template value<double>(vars, x_position, t, dt);
343 ip_data_output.fluid_density = fluid_density;
344 vars.density = fluid_density;
345
346 auto const drho_dp =
347 liquid_phase.property(MaterialPropertyLib::PropertyType::density)
348 .template dValue<double>(
350 x_position, t, dt);
351 crv.drho_LR_dp = drho_dp;
352
353 crv.fluid_compressibility = 1 / fluid_density * drho_dp;
354
355 crv.drho_LR_dT =
356 liquid_phase.property(MaterialPropertyLib::PropertyType::density)
357 .template dValue<double>(vars,
359 x_position, t, dt);
360
361 double const fluid_volumetric_thermal_expansion_coefficient =
363 liquid_phase, vars, fluid_density, x_position, t, dt);
364
365 // Use the viscosity model to compute the viscosity
366 ip_data_output.viscosity =
368 .template value<double>(vars, x_position, t, dt);
369 crv.K_over_mu = intrinsic_permeability / ip_data_output.viscosity;
370
371 crv.k_rel = 1.;
372 crv.dk_rel_dT = 0;
373 if (frozen_liquid_phase)
374 {
375 ip_data.phi_fr =
377 .template value<double>(vars, x_position, t, dt);
378
379 double const dphi_fr_dT =
381 .template dValue<double>(
383 x_position, t, dt);
384
385 // Set ice_volume_fraction variable for relative permeability
386 // calculation ice_volume_fraction = fraction of pore space occupied by
387 // ice
388 vars.ice_volume_fraction = ip_data.phi_fr / porosity;
389
390 crv.k_rel =
391 liquid_phase
392 .property(
394 .template value<double>(vars, x_position, t, dt);
395
396 // dk_rel/dT = (dk_rel/dphi_ice) * (dphi_ice/dT)
397 // where dphi_ice/dT = dphi_fr_dT / porosity
398 auto const dk_rel_dphi_ice =
399 liquid_phase
400 .property(
402 .template dValue<double>(
404 x_position, t, dt);
405 crv.dk_rel_dT = dk_rel_dphi_ice * dphi_fr_dT / porosity;
406 }
407
408 auto const& b = _process_data.specific_body_force;
409
410 // Consider also anisotropic thermal expansion.
411 crv.solid_linear_thermal_expansion_coefficient =
413 solid_phase
414 .property(
416 .value(vars, x_position, t, dt));
417
419 dthermal_strain =
420 crv.solid_linear_thermal_expansion_coefficient * dT_int_pt;
421
422 crv.K_pT_thermal_osmosis =
423 (solid_phase.hasProperty(
426 solid_phase
428 thermal_osmosis_coefficient)
429 .value(vars, x_position, t, dt))
430 : Eigen::MatrixXd::Zero(DisplacementDim, DisplacementDim));
431
432 GlobalDimVectorType const velocity =
433 -crv.k_rel * crv.K_over_mu * dNdx * p -
434 crv.K_pT_thermal_osmosis * dNdx * T +
435 (fluid_density * crv.k_rel) * crv.K_over_mu * b;
436 ip_data_output.velocity = velocity;
437 crv.dvelocity_dT =
438 -crv.dk_rel_dT * crv.K_over_mu * dNdx * p +
439 // TODO(naumov): - crv.K_pT_thermal_osmosis * dNdx * dT_dT +
440 (crv.dk_rel_dT * fluid_density + crv.k_rel * crv.drho_LR_dT) *
441 crv.K_over_mu * b;
442
443 //
444 // displacement equation, displacement part
445 //
446 auto& eps_m = ip_data.eps_m;
447 auto& eps_m_prev = ip_data.eps_m_prev;
448 eps_m.noalias() = eps_m_prev + eps - eps_prev - dthermal_strain;
449
450 crv.eps_v_dot = (vars.volumetric_strain - Invariants::trace(eps_prev)) / dt;
451
454 eps_m);
455
456 crv.C = ip_data.updateConstitutiveRelation(vars, t, x_position, dt,
457 T_prev_int_pt);
458
459 crv.rho = solid_density * (1 - porosity) + porosity * fluid_density;
460
461 crv.beta =
462 porosity * fluid_volumetric_thermal_expansion_coefficient +
463 (alpha - porosity) *
464 Invariants::trace(crv.solid_linear_thermal_expansion_coefficient);
465
466 //
467 // pressure equation, displacement part.
468 //
469 // Reusing Kup.transpose().
470
471 //
472 // temperature equation, temperature part.
473 //
474 crv.c_f =
475 liquid_phase
477 .template value<double>(vars, x_position, t, dt);
478 crv.effective_thermal_conductivity =
480 medium
481 ->property(
483 .value(vars, x_position, t, dt));
484
486 medium
489 x_position, t, dt));
490
491 // Thermal conductivity is moved outside and zero matrix is passed instead
492 // due to multiplication with fluid's density times specific heat capacity.
493 crv.effective_thermal_conductivity.noalias() +=
494 fluid_density * crv.c_f *
496 _process_data.stabilizer, _element.getID(),
497 GlobalDimMatrixType::Zero(DisplacementDim, DisplacementDim),
498 velocity, 0. /* phi */, 0. /* dispersivity_transversal */,
499 0. /*dispersivity_longitudinal*/);
500
501 double const c_s =
502 solid_phase
504 .template value<double>(vars, x_position, t, dt);
505
506 // Also modified by freezing terms.
507 crv.effective_volumetric_heat_capacity =
508 porosity * fluid_density * crv.c_f +
509 (1.0 - porosity) * solid_density * c_s;
510 double dC_eff_dT = porosity * crv.drho_LR_dT * crv.c_f +
511 (1.0 - porosity) * drho_SR_dT * c_s;
512
513 if (frozen_liquid_phase)
514 {
516 double const phi_fr = ip_data.phi_fr;
517
518 auto const frozen_liquid_value =
520 {
521 return (*frozen_liquid_phase)[p].template value<double>(
522 vars, x_position, t, dt);
523 };
524
525 double const c_fr = frozen_liquid_value(
527
528 double const l_fr = frozen_liquid_value(
530
531 double const dphi_fr_dT =
533 .template dValue<double>(
535 x_position, t, dt);
536
537 double const d2phi_fr_dT2 =
539 .template d2Value<double>(
542 dt);
543
544 double const phi_fr_prev = [&]()
545 {
547 vars_prev.temperature = T_prev_int_pt;
549 .template value<double>(vars_prev, x_position, t, dt);
550 }();
551 ip_data.phi_fr_prev = phi_fr_prev;
552
553 double const rho_fr =
555 ip_data_output.rho_fr = rho_fr;
556
557 crv.rho += ip_data.phi_fr * rho_fr - ip_data.phi_fr * fluid_density;
558 crv.mass_exchange =
559 -dphi_fr_dT * porosity * (1. - rho_fr / fluid_density);
560 double const dmass_exchange_dT =
561 -d2phi_fr_dT2 * porosity * (1. - rho_fr / fluid_density) +
562 dphi_fr_dT * porosity * rho_fr * crv.drho_LR_dT /
563 (fluid_density * fluid_density);
564
565 // alpha_T^I
567 DisplacementDim> const ice_linear_thermal_expansion_coefficient =
569 frozen_liquid_phase
570 ->property(
572 .value(vars, x_position, t, dt));
573
575 dthermal_strain_ice =
576 ice_linear_thermal_expansion_coefficient * dT_int_pt;
577
578 crv.beta_T_SI =
579 porosity *
580 Invariants::trace(ice_linear_thermal_expansion_coefficient) +
581 (alpha - porosity) *
583 crv.solid_linear_thermal_expansion_coefficient);
584
585 // alpha_{phi_I} -- linear expansion coeff. due to water-to-ice
586 // transition (phase change), and related phase_change_strain term
588 phase_change_expansion_coefficient =
590 frozen_liquid_phase
592 phase_change_expansivity)
593 .value(vars, x_position, t, dt));
594
596 dphase_change_strain = phase_change_expansion_coefficient *
597 (phi_fr - phi_fr_prev) / porosity;
598
599 // eps0 ia a 'history variable' -- a solid matrix strain accrued
600 // prior to the onset of ice forming
601 auto& eps0 = ip_data.eps0;
602 auto const& eps0_prev = ip_data.eps0_prev;
603
604 // definition of eps_m_ice
605 auto& eps_m_ice = ip_data.eps_m_ice;
606 auto const& eps_m_ice_prev = ip_data.eps_m_ice_prev;
607
608 eps_m_ice.noalias() = eps_m_ice_prev + eps - eps_prev -
609 (eps0 - eps0_prev) - dthermal_strain_ice -
610 dphase_change_strain;
611
612 vars_ice.mechanical_strain
614 eps_m_ice);
615 auto const C_IR = ip_data.updateConstitutiveRelationIce(
616 *_process_data.ice_constitutive_relation, vars_ice, t, x_position,
617 dt, T_prev_int_pt);
618
619 auto const C_el_ice = ip_data.computeElasticTangentStiffnessIce(
620 *_process_data.ice_constitutive_relation, t, x_position, dt,
621 static_cast<double>(T_int_pt));
622 crv.beta_IR =
623 1. / _process_data.ice_constitutive_relation->getBulkModulus(
624 t, x_position, &C_el_ice);
625
626 crv.effective_volumetric_heat_capacity +=
627 -phi_fr * fluid_density * crv.c_f + phi_fr * rho_fr * c_fr -
628 l_fr * rho_fr * dphi_fr_dT;
629
630 crv.J_uu_fr = phi_fr * C_IR;
631
632 auto const& sigma_eff_ice = ip_data.sigma_eff_ice;
633 crv.r_u_fr = phi_fr * sigma_eff_ice;
634
635 crv.J_uT_fr = phi_fr * C_IR * ice_linear_thermal_expansion_coefficient;
636
637 // part of dMTT_dT derivative for freezing
638 dC_eff_dT += -dphi_fr_dT * fluid_density * crv.c_f -
639 phi_fr * crv.drho_LR_dT * crv.c_f +
640 dphi_fr_dT * rho_fr * c_fr - l_fr * rho_fr * d2phi_fr_dT2;
641 double const storage_p_fr_coeff =
642 (porosity * crv.beta_IR + (alpha - porosity) * crv.beta_SR) *
643 rho_fr / fluid_density -
644 (porosity * crv.fluid_compressibility +
645 (alpha - porosity) * crv.beta_SR);
646 crv.storage_p_fr = phi_fr / porosity * storage_p_fr_coeff;
647
648 double const dstorage_p_fr_coeff_dT =
649 (porosity * crv.beta_IR + (alpha - porosity) * crv.beta_SR) *
650 crv.drho_LR_dT / (fluid_density * fluid_density);
651
652 crv.J_pT_fr = (dphi_fr_dT * storage_p_fr_coeff +
653 phi_fr * dstorage_p_fr_coeff_dT) /
654 porosity * dp_int_pt / dt;
655
656 crv.storage_T_fr =
657 phi_fr / porosity *
658 (crv.beta_T_SI * rho_fr / fluid_density - crv.beta) -
659 crv.mass_exchange;
660 double const dstorage_T_fr_dT =
661 dphi_fr_dT / porosity *
662 (crv.beta_T_SI * rho_fr / fluid_density - crv.beta) +
663 phi_fr / porosity * crv.beta_T_SI * rho_fr * crv.drho_LR_dT /
664 (fluid_density * fluid_density) -
665 dmass_exchange_dT;
666 crv.J_pT_fr += dstorage_T_fr_dT * dT_int_pt / dt;
667 }
668 crv.J_TT = dC_eff_dT * dT_int_pt / dt;
669 return crv;
670}
671
672// Assembles the local Jacobian matrix. So far, the linearisation of HT part is
673// not considered as that in HT process.
674template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
675 int DisplacementDim>
677 ShapeFunctionDisplacement, ShapeFunctionPressure,
678 DisplacementDim>::assembleWithJacobian(double const t, double const dt,
679 std::vector<double> const& local_x,
680 std::vector<double> const&
681 local_x_prev,
682 std::vector<double>& local_rhs_data,
683 std::vector<double>& local_Jac_data)
684{
685 assert(local_x.size() ==
687
688 auto const x =
689 Eigen::Map<Eigen::VectorXd const>(local_x.data(), local_x.size());
690 auto const x_prev = Eigen::Map<Eigen::VectorXd const>(local_x_prev.data(),
691 local_x_prev.size());
692
693 auto const [T, p, u] = localDOF(local_x);
694 auto const [T_prev, p_prev, u_prev] = localDOF(local_x_prev);
695
696 auto local_Jac = MathLib::createZeroedMatrix<
697 typename ShapeMatricesTypeDisplacement::template MatrixType<
702
703 auto local_rhs = MathLib::createZeroedVector<
704 typename ShapeMatricesTypeDisplacement::template VectorType<
707
710
713
715 KTp.setZero(temperature_size, pressure_size);
716
718 dKTT_dT_T.setZero(temperature_size, pressure_size);
719
721 dKTT_dp_T.setZero(temperature_size, pressure_size);
722
724 laplace_p.setZero(pressure_size, pressure_size);
725
727 laplace_T.setZero(pressure_size, temperature_size);
728
730 storage_p.setZero(pressure_size, pressure_size);
731
733 storage_T.setZero(pressure_size, temperature_size);
734
735 typename ShapeMatricesTypeDisplacement::template MatrixType<
737 Kup;
738 Kup.setZero(displacement_size, pressure_size);
739
740 typename ShapeMatricesTypeDisplacement::template MatrixType<
742 Kpu;
743 if (!_process_data.is_volume_balance_equation_type)
744 {
745 Kpu.setZero(pressure_size, displacement_size);
746 }
747
748 auto const& medium = _process_data.media_map.getMedium(_element.getID());
749 bool const has_frozen_liquid_phase =
751
752 unsigned const n_integration_points =
753 _integration_method.getNumberOfPoints();
754
755 std::vector<GlobalDimVectorType> ip_flux_vector;
756 double average_velocity_norm = 0.0;
757 ip_flux_vector.reserve(n_integration_points);
758
759 for (unsigned ip = 0; ip < n_integration_points; ip++)
760 {
761 auto& ip_data = _ip_data[ip];
762 auto const& N_u = ip_data.N_u;
763 ParameterLib::SpatialPosition const x_position{
764 std::nullopt, _element.getID(),
766 NumLib::interpolateCoordinates<ShapeFunctionDisplacement,
768 _element, N_u))};
769
770 auto const crv = updateConstitutiveRelations(
771 x, x_prev, x_position, t, dt, ip_data, _ip_data_output[ip]);
772
773 auto const& w = ip_data.integration_weight;
774
775 auto const& dNdx_u = ip_data.dNdx_u;
776
777 auto const& N = ip_data.N;
778 auto const& dNdx = ip_data.dNdx;
779
780 auto const T_int_pt = N.dot(T);
781
782 auto const x_coord =
783 x_position.getCoordinates().value()[0]; // r for axisymmetry
784 auto const B =
785 LinearBMatrix::computeBMatrix<DisplacementDim,
786 ShapeFunctionDisplacement::NPOINTS,
788 dNdx_u, N_u, x_coord, _is_axially_symmetric);
789
790 auto const& b = _process_data.specific_body_force;
791 auto const velocity = _ip_data_output[ip].velocity;
792
793 //
794 // displacement equation, displacement part
795 //
796
797 auto const C_eff = has_frozen_liquid_phase
798 ? (crv.C + crv.J_uu_fr).eval()
799 : crv.C.eval();
800 local_Jac
801 .template block<displacement_size, displacement_size>(
803 .noalias() += B.transpose() * C_eff * B * w;
804
805 auto const uT_coeff =
806 has_frozen_liquid_phase
807 ? (crv.J_uT_fr +
808 crv.C * crv.solid_linear_thermal_expansion_coefficient)
809 .eval()
810 : (crv.C * crv.solid_linear_thermal_expansion_coefficient)
811 .eval();
812
813 if (has_frozen_liquid_phase)
814 {
815 local_rhs.template segment<displacement_size>(displacement_index)
816 .noalias() -= B.transpose() * crv.r_u_fr * w;
817 }
818
819 local_Jac
820 .template block<displacement_size, temperature_size>(
822 .noalias() -= B.transpose() * uT_coeff * N * w;
823
824 local_rhs.template segment<displacement_size>(displacement_index)
825 .noalias() -= (B.transpose() * ip_data.sigma_eff -
826 N_u_op(N_u).transpose() * crv.rho * b) *
827 w;
828
829 //
830 // displacement equation, pressure part (K_up)
831 //
832 double const fluid_density = _ip_data_output[ip].fluid_density;
833 double const up_coeff =
834 has_frozen_liquid_phase
835 ? crv.alpha_biot * ip_data.phi_fr / ip_data.porosity *
836 (_ip_data_output[ip].rho_fr / fluid_density - 1) +
837 crv.alpha_biot
838 : crv.alpha_biot;
839 Kup.noalias() +=
840 B.transpose() * Invariants::identity2 * N * (up_coeff * w);
841
842 //
843 // pressure equation, pressure part (K_pp and M_pp).
844 //
845 double const scaling_factor =
846 _process_data.is_volume_balance_equation_type ? 1.0 : fluid_density;
847
848 laplace_p.noalias() += dNdx.transpose() * crv.K_over_mu * dNdx *
849 (crv.k_rel * scaling_factor * w);
850 local_Jac
851 .template block<pressure_size, temperature_size>(pressure_index,
853 .noalias() += dNdx.transpose() * crv.K_over_mu * (dNdx * p) * N *
854 (crv.dk_rel_dT * scaling_factor * w);
855
856 double const storage_p_coeff_no_fr =
857 ip_data.porosity * crv.fluid_compressibility +
858 (crv.alpha_biot - ip_data.porosity) * crv.beta_SR;
859 double const storage_p_coeff =
860 has_frozen_liquid_phase ? crv.storage_p_fr + storage_p_coeff_no_fr
861 : storage_p_coeff_no_fr;
862
863 storage_p.noalias() +=
864 N.transpose() * N * (storage_p_coeff * scaling_factor * w);
865
866 if (has_frozen_liquid_phase)
867 {
868 local_Jac
869 .template block<pressure_size, temperature_size>(
871 .noalias() +=
872 N.transpose() * crv.J_pT_fr * N * scaling_factor * w;
873 }
874
875 laplace_T.noalias() += dNdx.transpose() * crv.K_pT_thermal_osmosis *
876 dNdx * scaling_factor * w;
877 //
878 // RHS, pressure part
879 //
880 local_rhs.template segment<pressure_size>(pressure_index).noalias() -=
881 N * (up_coeff * crv.eps_v_dot * scaling_factor * w);
882
883 local_rhs.template segment<pressure_size>(pressure_index).noalias() +=
884 dNdx.transpose() * crv.K_over_mu * b *
885 (fluid_density * crv.k_rel * scaling_factor * w);
886 local_Jac
887 .template block<pressure_size, temperature_size>(pressure_index,
889 .noalias() -= dNdx.transpose() * crv.K_over_mu * b * N *
890 (fluid_density * crv.dk_rel_dT * scaling_factor * w);
891
892 //
893 // pressure equation, temperature part (M_pT)
894 //
895
896 double const storage_T_coeff =
897 has_frozen_liquid_phase ? crv.storage_T_fr + crv.beta : crv.beta;
898
899 storage_T.noalias() +=
900 N.transpose() * storage_T_coeff * N * scaling_factor * w;
901
902 //
903 // pressure equation, displacement part.
904 //
905 // reusing Kup.transpose() if the equation balance type is not volume.
906
907 if (!_process_data.is_volume_balance_equation_type)
908 {
909 Kpu.noalias() += N.transpose() * Invariants::identity2.transpose() *
910 B * (up_coeff * scaling_factor * w);
911
912 //
913 // The contribution to Jacobian from d()/ drho
914 // drho/dp, d()/ drho drho/dT:
915 //
916 double const storage_p_solid_coeff =
917 (crv.alpha_biot - ip_data.porosity) * crv.beta_SR;
918
919 double const p_dot = N.dot(p - p_prev) / dt;
920 double const T_dot = N.dot(T - T_prev) / dt;
921 double const drho_dp_coeff = storage_p_solid_coeff * p_dot +
922 storage_T_coeff * T_dot +
923 up_coeff * crv.eps_v_dot;
924
925 local_Jac
926 .template block<pressure_size, pressure_size>(pressure_index,
928 .noalias() +=
929 // TODO (WW) : Add ip_data.porosity * d2rho_LR_dp2 * w.
930 N.transpose() * N * (drho_dp_coeff * crv.drho_LR_dp * w);
931 local_Jac
932 .template block<pressure_size, temperature_size>(
934 .noalias() +=
935 // TODO (WW) : Add ip_data.porosity * d2rho_LR_dpdT * w.
936 N.transpose() * N * (drho_dp_coeff * crv.drho_LR_dT * w);
937
938 // The term from d (rho_L K(grad p - rho_L b)/dp:
939 // derivative of rhp_L * K_over_mu * k_rel (grad p + rho_l g) with
940 // respect to pressure and temperature.
941 auto const dlaplace_temporal_factor =
942 (-velocity - (fluid_density * crv.k_rel) * crv.K_over_mu * b);
943 local_Jac
944 .template block<pressure_size, pressure_size>(pressure_index,
946 .noalias() += dNdx.transpose() * dlaplace_temporal_factor * N *
947 (crv.drho_LR_dp * w);
948 local_Jac
949 .template block<pressure_size, temperature_size>(
951 .noalias() += dNdx.transpose() * dlaplace_temporal_factor * N *
952 (crv.drho_LR_dT * w);
953 }
954
955 //
956 // temperature equation, temperature part.
957 //
958 KTT.noalias() +=
959 dNdx.transpose() * crv.effective_thermal_conductivity * dNdx * w;
960 dKTT_dT_T.noalias() +=
961 dNdx.transpose() * crv.dlambda_eff_dT * dNdx * T * N * w;
962
963 ip_flux_vector.emplace_back(velocity * fluid_density * crv.c_f);
964 // Without any flux correction the flux derivative is as follows. The
965 // contribution to KTT is different if any stabilization scheme is used,
966 // but this is ignored for the moment.
967 GlobalDimVectorType const dip_flux_vector_dT =
968 crv.dvelocity_dT * fluid_density * crv.c_f +
969 velocity * crv.drho_LR_dT * crv.c_f;
970 dKTT_dT_T.noalias() +=
971 N.transpose() * dip_flux_vector_dT.transpose() * dNdx * T * N * w;
972 average_velocity_norm += velocity.norm();
973
974 MTT.noalias() +=
975 N.transpose() * crv.effective_volumetric_heat_capacity * N * w;
976 local_Jac
977 .template block<temperature_size, temperature_size>(
979 .noalias() += N.transpose() * crv.J_TT * N * w;
980
981 //
982 // temperature equation, pressure part
983 //
984 KTp.noalias() +=
985 dNdx.transpose() * crv.K_pT_thermal_osmosis * dNdx * (T_int_pt * w);
986
987 // linearized darcy
988 dKTT_dp_T.noalias() -= N.transpose() * (dNdx * T).transpose() *
989 crv.K_over_mu * dNdx *
990 (fluid_density * crv.c_f * crv.k_rel * w);
991
992 /* TODO (Joerg) Temperature changes due to thermal dilatation of the
993 * fluid, which are usually discarded as being very small.
994 * Zhou et al. (10.1016/S0020-7683(98)00089-4) states that:
995 * "Biot (1956) neglected this term and it is included here for
996 * completeness"
997 * Keeping the code here in the case these are needed for the named
998 * effects in the future.
999 if (fluid_compressibility != 0)
1000 {
1001 auto const C_el = ip_data.computeElasticTangentStiffness(
1002 t, x_position, dt, static_cast<double>(T_int_pt));
1003 auto const solid_skeleton_compressibility =
1004 1 / solid_material.getBulkModulus(t, x_position, &C_el);
1005 double const fluid_volumetric_thermal_expansion_coefficient =
1006 MaterialPropertyLib::getLiquidThermalExpansivity(
1007 liquid_phase, vars, fluid_density, x_position, t, dt);
1008
1009 KTT.noalias() +=
1010 dNdx.transpose() *
1011 (-T_int_pt * fluid_volumetric_thermal_expansion_coefficient *
1012 K_pT_thermal_osmosis / fluid_compressibility) *
1013 dNdx * w;
1014
1015 local_rhs.template segment<temperature_size>(temperature_index)
1016 .noalias() +=
1017 dNdx.transpose() *
1018 (-T_int_pt * fluid_volumetric_thermal_expansion_coefficient /
1019 fluid_compressibility) *
1020 fluid_density * crv.k_rel * K_over_mu * b * w;
1021 MTu part for rhs and Jacobian:
1022 (-T_int_pt *
1023 Invariants::trace(solid_linear_thermal_expansion_coefficient) /
1024 solid_skeleton_compressibility) *
1025 N.transpose() * identity2.transpose() * B * w;
1026 KTp part for rhs and Jacobian:
1027 dNdx.transpose() *
1028 (T_int_pt * fluid_volumetric_thermal_expansion_coefficient *
1029 crv.k_rel * K_over_mu / fluid_compressibility) *
1030 dNdx * w;
1031 }
1032 */
1033 }
1034
1036 _process_data.stabilizer, _ip_data, ip_flux_vector,
1037 average_velocity_norm / static_cast<double>(n_integration_points), KTT);
1038
1039 // temperature equation, temperature part
1040 local_Jac
1041 .template block<temperature_size, temperature_size>(temperature_index,
1043 .noalias() += KTT + dKTT_dT_T + MTT / dt;
1044
1045 // temperature equation, pressure part
1046 local_Jac
1047 .template block<temperature_size, pressure_size>(temperature_index,
1049 .noalias() += KTp + dKTT_dp_T;
1050
1051 // displacement equation, pressure part
1052 local_Jac
1053 .template block<displacement_size, pressure_size>(displacement_index,
1055 .noalias() -= Kup;
1056
1057 // pressure equation, temperature part.
1058 local_Jac
1059 .template block<pressure_size, temperature_size>(pressure_index,
1061 .noalias() += -storage_T / dt + laplace_T;
1062
1063 // pressure equation, pressure part.
1064 local_Jac
1065 .template block<pressure_size, pressure_size>(pressure_index,
1067 .noalias() += laplace_p + storage_p / dt;
1068
1069 // pressure equation, displacement part.
1070 if (_process_data.is_volume_balance_equation_type)
1071 {
1072 local_Jac
1073 .template block<pressure_size, displacement_size>(
1075 .noalias() += Kup.transpose() / dt;
1076 }
1077 else
1078 {
1079 local_Jac
1080 .template block<pressure_size, displacement_size>(
1082 .noalias() += Kpu / dt;
1083 }
1084
1085 // pressure equation (f_p)
1086 local_rhs.template segment<pressure_size>(pressure_index).noalias() -=
1087 laplace_p * p + laplace_T * T + storage_p * (p - p_prev) / dt -
1088 storage_T * (T - T_prev) / dt;
1089
1090 // displacement equation (f_u)
1091 local_rhs.template segment<displacement_size>(displacement_index)
1092 .noalias() += Kup * p;
1093
1094 // temperature equation (f_T)
1095 local_rhs.template segment<temperature_size>(temperature_index).noalias() -=
1096 KTT * T + MTT * (T - T_prev) / dt;
1097
1098 local_rhs.template segment<temperature_size>(temperature_index).noalias() -=
1099 KTp * p;
1100}
1101
1102template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
1103 int DisplacementDim>
1104std::vector<double> const& ThermoHydroMechanicsLocalAssembler<
1105 ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim>::
1106 getIntPtDarcyVelocity(
1107 const double /*t*/,
1108 std::vector<GlobalVector*> const& /*x*/,
1109 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
1110 std::vector<double>& cache) const
1111{
1112 unsigned const n_integration_points =
1113 _integration_method.getNumberOfPoints();
1114
1115 cache.clear();
1116 auto cache_matrix = MathLib::createZeroedMatrix<Eigen::Matrix<
1117 double, DisplacementDim, Eigen::Dynamic, Eigen::RowMajor>>(
1118 cache, DisplacementDim, n_integration_points);
1119
1120 for (unsigned ip = 0; ip < n_integration_points; ip++)
1121 {
1122 cache_matrix.col(ip).noalias() = _ip_data_output[ip].velocity;
1123 }
1124
1125 return cache;
1126}
1127
1128template <typename ShapeFunctionDisplacement, typename ShapeFunction,
1129 int DisplacementDim>
1130std::vector<double> const& ThermoHydroMechanicsLocalAssembler<
1131 ShapeFunctionDisplacement, ShapeFunction, DisplacementDim>::
1132 getIntPtFluidDensity(
1133 const double /*t*/,
1134 std::vector<GlobalVector*> const& /*x*/,
1135 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
1136 std::vector<double>& cache) const
1137{
1141}
1142
1143template <typename ShapeFunctionDisplacement, typename ShapeFunction,
1144 int DisplacementDim>
1145std::vector<double> const& ThermoHydroMechanicsLocalAssembler<
1146 ShapeFunctionDisplacement, ShapeFunction, DisplacementDim>::
1147 getIntPtViscosity(
1148 const double /*t*/,
1149 std::vector<GlobalVector*> const& /*x*/,
1150 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
1151 std::vector<double>& cache) const
1152{
1156}
1157
1158template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
1159 int DisplacementDim>
1161 ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim>::
1162 computeSecondaryVariableConcrete(double const /*t*/, double const /*dt*/,
1163 Eigen::VectorXd const& local_x,
1164 Eigen::VectorXd const& /*local_x_prev*/)
1165{
1166 auto const p = local_x.template segment<pressure_size>(pressure_index);
1167 auto const T =
1168 local_x.template segment<temperature_size>(temperature_index);
1169
1170 unsigned const n_integration_points =
1171 _integration_method.getNumberOfPoints();
1172
1173 double phi_fr_avg = 0;
1174 double fluid_density_avg = 0;
1175 double viscosity_avg = 0;
1176
1178 KV sigma_avg = KV::Zero();
1179 KV sigma_ice_avg = KV::Zero();
1180
1181 for (unsigned ip = 0; ip < n_integration_points; ip++)
1182 {
1183 auto& ip_data = _ip_data[ip];
1184
1185 phi_fr_avg += ip_data.phi_fr;
1186 fluid_density_avg += _ip_data_output[ip].fluid_density;
1187 viscosity_avg += _ip_data_output[ip].viscosity;
1188 sigma_avg += ip_data.sigma_eff;
1189 sigma_ice_avg += ip_data.sigma_eff_ice;
1190 }
1191
1192 phi_fr_avg /= n_integration_points;
1193 fluid_density_avg /= n_integration_points;
1194 viscosity_avg /= n_integration_points;
1195 sigma_avg /= n_integration_points;
1196 sigma_ice_avg /= n_integration_points;
1197
1198 (*_process_data.element_phi_fr)[_element.getID()] = phi_fr_avg;
1199 (*_process_data.element_fluid_density)[_element.getID()] =
1200 fluid_density_avg;
1201 (*_process_data.element_viscosity)[_element.getID()] = viscosity_avg;
1202
1203 Eigen::Map<KV>(&(*_process_data.element_stresses)[_element.getID() *
1204 KV::RowsAtCompileTime]) =
1206
1207 Eigen::Map<KV>(&(
1209 .element_ice_stresses)[_element.getID() * KV::RowsAtCompileTime]) =
1211
1213 ShapeFunctionPressure, typename ShapeFunctionDisplacement::MeshElement,
1214 DisplacementDim>(_element, _is_axially_symmetric, p,
1215 *_process_data.pressure_interpolated);
1216
1218 ShapeFunctionPressure, typename ShapeFunctionDisplacement::MeshElement,
1219 DisplacementDim>(_element, _is_axially_symmetric, T,
1220 *_process_data.temperature_interpolated);
1221}
1222} // namespace ThermoHydroMechanics
1223} // namespace ProcessLib
#define OGS_FATAL(...)
Definition Error.h:19
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:22
std::size_t getID() const
Returns the ID of the element.
Definition Element.h:80
std::optional< MathLib::Point3d > const getCoordinates() const
MatrixType< _kelvin_vector_size, _number_of_dof > BMatrixType
void setInitialConditionsConcrete(Eigen::VectorXd const local_x, double const t, int const process_id) override
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
typename ShapeMatricesTypePressure::GlobalDimVectorType GlobalDimVectorType
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
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)
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::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)
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::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)
MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > computeElasticTangentStiffnessIce(MaterialLib::Solids::MechanicsBase< DisplacementDim > const &ice_constitutive_relation, double const t, ParameterLib::SpatialPosition const &x_position, double const dt, double const temperature)
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)