OGS
ThermoHydroMechanicsFEM-impl.h
Go to the documentation of this file.
1
11#pragma once
12
13#include <typeinfo>
14
28
29namespace ProcessLib
30{
31namespace ThermoHydroMechanics
32{
33template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
34 int DisplacementDim>
35ThermoHydroMechanicsLocalAssembler<ShapeFunctionDisplacement,
36 ShapeFunctionPressure, DisplacementDim>::
37 ThermoHydroMechanicsLocalAssembler(
38 MeshLib::Element const& e,
39 std::size_t const /*local_matrix_size*/,
40 NumLib::GenericIntegrationMethod const& integration_method,
41 bool const is_axially_symmetric,
43 : _process_data(process_data),
44 _integration_method(integration_method),
45 _element(e),
46 _is_axially_symmetric(is_axially_symmetric)
47{
48 unsigned const n_integration_points =
50
51 _ip_data.reserve(n_integration_points);
52 _ip_data_output.resize(n_integration_points);
53 _secondary_data.N_u.resize(n_integration_points);
54
55 auto const shape_matrices_u =
56 NumLib::initShapeMatrices<ShapeFunctionDisplacement,
58 DisplacementDim>(e, is_axially_symmetric,
60
61 auto const shape_matrices_p =
62 NumLib::initShapeMatrices<ShapeFunctionPressure,
63 ShapeMatricesTypePressure, DisplacementDim>(
64 e, is_axially_symmetric, _integration_method);
65
66 auto const& solid_material =
68 _process_data.solid_materials, _process_data.material_ids,
69 e.getID());
70
71 // Consistency check: if frozen liquid phase is given, then the constitutive
72 // relation for ice must also be given, and vice versa.
73 auto const& medium = _process_data.media_map.getMedium(_element.getID());
74 if (medium->hasPhase("FrozenLiquid") !=
75 (_process_data.ice_constitutive_relation != nullptr))
76 {
78 "Frozen liquid phase is {:s} and the solid material constitutive "
79 "relation for ice is {:s}. But both must be given (or both "
80 "omitted).",
81 medium->hasPhase("FrozenLiquid") ? "specified" : "not specified",
82 _process_data.ice_constitutive_relation != nullptr
83 ? "specified"
84 : "not specified");
85 }
86 for (unsigned ip = 0; ip < n_integration_points; ip++)
87 {
88 _ip_data.emplace_back(solid_material);
89 auto& ip_data = _ip_data[ip];
90 auto const& sm_u = shape_matrices_u[ip];
91 ip_data.integration_weight =
93 sm_u.integralMeasure * sm_u.detJ;
94
95 ip_data.N_u = sm_u.N;
96 ip_data.dNdx_u = sm_u.dNdx;
97
98 ip_data.N = shape_matrices_p[ip].N;
99 ip_data.dNdx = shape_matrices_p[ip].dNdx;
100
101 _secondary_data.N_u[ip] = shape_matrices_u[ip].N;
102 }
103}
104
105template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
106 int DisplacementDim>
108 ShapeFunctionDisplacement, ShapeFunctionPressure,
109 DisplacementDim>::setIPDataInitialConditions(std::string_view name,
110 double const* values,
111 int const integration_order)
112{
113 if (integration_order !=
114 static_cast<int>(_integration_method.getIntegrationOrder()))
115 {
116 OGS_FATAL(
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.",
120 _element.getID());
121 }
122
123 if (name == "sigma")
124 {
125 if (_process_data.initial_stress.value)
126 {
127 OGS_FATAL(
128 "Setting initial conditions for stress from integration "
129 "point data and from a parameter '{:s}' is not possible "
130 "simultaneously.",
131 _process_data.initial_stress.value->name);
132 }
133
135 values, _ip_data, &IpData::sigma_eff);
136 }
137 if (name == "epsilon_m")
138 {
140 values, _ip_data, &IpData::eps_m);
141 }
142 if (name == "epsilon")
143 {
145 values, _ip_data, &IpData::eps);
146 }
147 if (name.starts_with("material_state_variable_"))
148 {
149 name.remove_prefix(24);
150
151 // Using first ip data for solid material. TODO (naumov) move solid
152 // material into element, store only material state in IPs.
153 auto const& internal_variables =
154 _ip_data[0].solid_material.getInternalVariables();
155 if (auto const iv = std::find_if(
156 begin(internal_variables), end(internal_variables),
157 [&name](auto const& iv) { return iv.name == name; });
158 iv != end(internal_variables))
159 {
160 DBUG("Setting material state variable '{:s}'", name);
162 values, _ip_data, &IpData::material_state_variables,
163 iv->reference);
164 }
165
166 int const element_id = _element.getID();
167 DBUG(
168 "The solid material of element {:d} (material ID {:d}) does not "
169 "have an internal state variable called {:s}.",
170 element_id, (*_process_data.material_ids)[element_id], name);
171 }
172
173 return 0;
174}
175template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
176 int DisplacementDim>
178 ShapeFunctionDisplacement, ShapeFunctionPressure,
179 DisplacementDim>::setInitialConditionsConcrete(Eigen::VectorXd const
180 local_x,
181 double const t,
182 int const /*process_id*/)
183{
184 // TODO: For staggered scheme, overload
185 // LocalAssemblerInterface::setInitialConditions to enable local_x contains
186 // the primary variables from all coupled processes.
187 auto const [T, p, u] = localDOF(local_x);
188
189 double const dt = 0.0;
190
192 auto const& medium = _process_data.media_map.getMedium(_element.getID());
193 auto* const frozen_liquid_phase = medium->hasPhase("FrozenLiquid")
194 ? &medium->phase("FrozenLiquid")
195 : nullptr;
196
197 int const n_integration_points = _integration_method.getNumberOfPoints();
198 for (int ip = 0; ip < n_integration_points; ip++)
199 {
200 auto const& N = _ip_data[ip].N;
201 auto const& N_u = _ip_data[ip].N_u;
202 ParameterLib::SpatialPosition const x_position{
203 std::nullopt, _element.getID(),
205 NumLib::interpolateCoordinates<ShapeFunctionDisplacement,
207 _element, N_u))};
208
209 auto& sigma_eff = _ip_data[ip].sigma_eff;
210 if (_process_data.initial_stress.isTotalStress())
211 {
212 auto const alpha_b =
213 medium->property(MPL::PropertyType::biot_coefficient)
214 .template value<double>(vars, x_position, t, dt);
215
216 sigma_eff.noalias() += alpha_b * N.dot(p) * Invariants::identity2;
217 }
218 _ip_data[ip].sigma_eff_prev.noalias() = sigma_eff;
219
220 vars.temperature = N.dot(T);
221 if (frozen_liquid_phase)
222 {
223 _ip_data[ip].phi_fr =
225 .template value<double>(vars, x_position, t, dt);
226 }
227 }
228}
229
230template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
231 int DisplacementDim>
233 ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim>::
234 updateConstitutiveRelations(
235 Eigen::Ref<Eigen::VectorXd const> const local_x,
236 Eigen::Ref<Eigen::VectorXd const> const local_x_prev,
237 ParameterLib::SpatialPosition const& x_position, double const t,
238 double const dt, IpData& ip_data,
240{
241 assert(local_x.size() ==
242 pressure_size + displacement_size + temperature_size);
243
244 auto const [T, p, u] = localDOF(local_x);
245 auto const [T_prev, p_prev, u_prev] = localDOF(local_x_prev);
246
247 auto const& solid_material =
249 _process_data.solid_materials, _process_data.material_ids,
250 _element.getID());
251
252 auto const& medium = _process_data.media_map.getMedium(_element.getID());
253 auto const& liquid_phase = medium->phase("AqueousLiquid");
254 auto const& solid_phase = medium->phase("Solid");
255 auto* const frozen_liquid_phase = medium->hasPhase("FrozenLiquid")
256 ? &medium->phase("FrozenLiquid")
257 : nullptr;
259
260 auto const& N_u = ip_data.N_u;
261 auto const& dNdx_u = ip_data.dNdx_u;
262
263 auto const& N = ip_data.N;
264 auto const& dNdx = ip_data.dNdx;
265
266 auto const T_int_pt = N.dot(T);
267 auto const T_prev_int_pt = N.dot(T_prev);
268 double const dT_int_pt = T_int_pt - T_prev_int_pt;
269
270 auto const x_coord =
271 NumLib::interpolateXCoordinate<ShapeFunctionDisplacement,
273 N_u);
274 auto const B =
275 LinearBMatrix::computeBMatrix<DisplacementDim,
276 ShapeFunctionDisplacement::NPOINTS,
278 dNdx_u, N_u, x_coord, _is_axially_symmetric);
279
281
282 auto& eps = ip_data.eps;
283 eps.noalias() = B * u;
285 B * u_prev;
286
287 vars.temperature = T_int_pt;
288 double const p_int_pt = N.dot(p);
289 vars.liquid_phase_pressure = p_int_pt;
290
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 porosity =
299 .template value<double>(vars, x_position, t, dt);
300 vars.porosity = porosity;
301 ip_data.porosity = porosity;
302
303 crv.alpha_biot =
305 .template value<double>(vars, x_position, t, dt);
306 auto const& alpha = crv.alpha_biot;
307
308 auto const C_el = ip_data.computeElasticTangentStiffness(
309 t, x_position, dt, static_cast<double>(T_int_pt));
310 auto const solid_skeleton_compressibility =
311 1 / solid_material.getBulkModulus(t, x_position, &C_el);
312
313 crv.beta_SR = (1 - alpha) * solid_skeleton_compressibility;
314
315 // Set mechanical variables for the intrinsic permeability model
316 // For stress dependent permeability.
317 {
318 auto const& identity2 = Invariants::identity2;
319 auto const sigma_total =
320 (ip_data.sigma_eff - alpha * p_int_pt * identity2).eval();
321 vars.total_stress.emplace<SymmetricTensor>(
323 }
324 // For strain dependent permeability
325 vars.volumetric_strain = Invariants::trace(ip_data.eps);
327 ip_data.material_state_variables->getEquivalentPlasticStrain();
328
329 auto const intrinsic_permeability =
332 .value(vars, x_position, t, dt));
333
334 auto const fluid_density =
335 liquid_phase.property(MaterialPropertyLib::PropertyType::density)
336 .template value<double>(vars, x_position, t, dt);
337 ip_data_output.fluid_density = fluid_density;
338 vars.density = fluid_density;
339
340 auto const drho_dp =
341 liquid_phase.property(MaterialPropertyLib::PropertyType::density)
342 .template dValue<double>(
344 x_position, t, dt);
345
346 crv.fluid_compressibility = 1 / fluid_density * drho_dp;
347
348 double const fluid_volumetric_thermal_expansion_coefficient =
350 liquid_phase, vars, fluid_density, x_position, t, dt);
351
352 // Use the viscosity model to compute the viscosity
353 ip_data_output.viscosity =
355 .template value<double>(vars, x_position, t, dt);
356 crv.K_over_mu = intrinsic_permeability / ip_data_output.viscosity;
357
358 auto const& b = _process_data.specific_body_force;
359
360 // Consider also anisotropic thermal expansion.
361 crv.solid_linear_thermal_expansion_coefficient =
363 solid_phase
364 .property(
366 .value(vars, x_position, t, dt));
367
369 dthermal_strain =
370 crv.solid_linear_thermal_expansion_coefficient * dT_int_pt;
371
372 crv.K_pT_thermal_osmosis =
373 (solid_phase.hasProperty(
376 solid_phase
378 thermal_osmosis_coefficient)
379 .value(vars, x_position, t, dt))
380 : Eigen::MatrixXd::Zero(DisplacementDim, DisplacementDim));
381
382 GlobalDimVectorType const velocity = -crv.K_over_mu * dNdx * p -
383 crv.K_pT_thermal_osmosis * dNdx * T +
384 crv.K_over_mu * fluid_density * b;
385 ip_data_output.velocity = velocity;
386
387 //
388 // displacement equation, displacement part
389 //
390 auto& eps_m = ip_data.eps_m;
391 auto& eps_m_prev = ip_data.eps_m_prev;
392 eps_m.noalias() = eps_m_prev + eps - eps_prev - dthermal_strain;
395 eps_m);
396
397 crv.C = ip_data.updateConstitutiveRelation(vars, t, x_position, dt,
398 T_prev_int_pt);
399
400 crv.rho = solid_density * (1 - porosity) + porosity * fluid_density;
401
402 crv.beta =
403 porosity * fluid_volumetric_thermal_expansion_coefficient +
404 (alpha - porosity) *
405 Invariants::trace(crv.solid_linear_thermal_expansion_coefficient);
406
407 //
408 // pressure equation, displacement part.
409 //
410 // Reusing Kup.transpose().
411
412 //
413 // temperature equation, temperature part.
414 //
415 crv.c_f =
416 liquid_phase
418 .template value<double>(vars, x_position, t, dt);
419 crv.effective_thermal_conductivity =
421 medium
422 ->property(
424 .value(vars, x_position, t, dt));
425
426 // Thermal conductivity is moved outside and zero matrix is passed instead
427 // due to multiplication with fluid's density times specific heat capacity.
428 crv.effective_thermal_conductivity.noalias() +=
429 fluid_density * crv.c_f *
431 _process_data.stabilizer, _element.getID(),
432 GlobalDimMatrixType::Zero(DisplacementDim, DisplacementDim),
433 velocity, 0. /* phi */, 0. /* dispersivity_transversal */,
434 0. /*dispersivity_longitudinal*/);
435
436 double const c_s =
437 solid_phase
439 .template value<double>(vars, x_position, t, dt);
440
441 // Also modified by freezing terms.
442 crv.effective_volumetric_heat_capacity =
443 porosity * fluid_density * crv.c_f +
444 (1.0 - porosity) * solid_density * c_s;
445
446 if (frozen_liquid_phase)
447 {
449 double const phi_fr =
451 .template value<double>(vars, x_position, t, dt);
452 ip_data.phi_fr = phi_fr;
453
454 auto const frozen_liquid_value =
456 {
457 return (*frozen_liquid_phase)[p].template value<double>(
458 vars, x_position, t, dt);
459 };
460
461 double const rho_fr =
463
464 double const c_fr = frozen_liquid_value(
466
467 double const l_fr = frozen_liquid_value(
469
470 double const dphi_fr_dT =
472 .template dValue<double>(
474 x_position, t, dt);
475
476 double const phi_fr_prev = [&]()
477 {
479 vars_prev.temperature = T_prev_int_pt;
481 .template value<double>(vars_prev, x_position, t, dt);
482 }();
483 ip_data.phi_fr_prev = phi_fr_prev;
484
485 // alpha_T^I
487 DisplacementDim> const ice_linear_thermal_expansion_coefficient =
489 frozen_liquid_phase
490 ->property(
492 .value(vars, x_position, t, dt));
493
495 dthermal_strain_ice =
496 ice_linear_thermal_expansion_coefficient * dT_int_pt;
497
498 // alpha_{phi_I} -- linear expansion coeff. due to water-to-ice
499 // transition (phase change), and related phase_change_strain term
501 phase_change_expansion_coefficient =
503 frozen_liquid_phase
505 phase_change_expansivity)
506 .value(vars, x_position, t, dt));
507
509 dphase_change_strain = phase_change_expansion_coefficient *
510 (phi_fr - phi_fr_prev) / porosity;
511
512 // eps0 ia a 'history variable' -- a solid matrix strain accrued
513 // prior to the onset of ice forming
514 auto& eps0 = ip_data.eps0;
515 auto const& eps0_prev = ip_data.eps0_prev;
516
517 // definition of eps_m_ice
518 auto& eps_m_ice = ip_data.eps_m_ice;
519 auto const& eps_m_ice_prev = ip_data.eps_m_ice_prev;
520
521 eps_m_ice.noalias() = eps_m_ice_prev + eps - eps_prev -
522 (eps0 - eps0_prev) - dthermal_strain_ice -
523 dphase_change_strain;
524
525 vars_ice.mechanical_strain
527 eps_m_ice);
528 auto const C_IR = ip_data.updateConstitutiveRelationIce(
529 *_process_data.ice_constitutive_relation, vars_ice, t, x_position,
530 dt, T_prev_int_pt);
531 crv.effective_volumetric_heat_capacity +=
532 -phi_fr * fluid_density * crv.c_f + phi_fr * rho_fr * c_fr -
533 l_fr * rho_fr * dphi_fr_dT;
534
535 // part of dMTT_dT derivative for freezing
536 double const d2phi_fr_dT2 =
538 .template d2Value<double>(
541 dt);
542
543 crv.J_uu_fr = phi_fr * C_IR;
544
545 auto const& sigma_eff_ice = ip_data.sigma_eff_ice;
546 crv.r_u_fr = phi_fr * sigma_eff_ice;
547
548 crv.J_uT_fr = phi_fr * C_IR * ice_linear_thermal_expansion_coefficient;
549
550 crv.J_TT_fr = ((rho_fr * c_fr - fluid_density * crv.c_f) * dphi_fr_dT +
551 l_fr * rho_fr * d2phi_fr_dT2) *
552 dT_int_pt / dt;
553 }
554 return crv;
555}
556
557// Assembles the local Jacobian matrix. So far, the linearisation of HT part is
558// not considered as that in HT process.
559template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
560 int DisplacementDim>
562 ShapeFunctionDisplacement, ShapeFunctionPressure,
563 DisplacementDim>::assembleWithJacobian(double const t, double const dt,
564 std::vector<double> const& local_x,
565 std::vector<double> const&
566 local_x_prev,
567 std::vector<double>& local_rhs_data,
568 std::vector<double>& local_Jac_data)
569{
570 assert(local_x.size() ==
571 pressure_size + displacement_size + temperature_size);
572
573 auto const x =
574 Eigen::Map<Eigen::VectorXd const>(local_x.data(), local_x.size());
575 auto const x_prev = Eigen::Map<Eigen::VectorXd const>(local_x_prev.data(),
576 local_x_prev.size());
577
578 auto const [T, p, u] = localDOF(local_x);
579 auto const [T_prev, p_prev, u_prev] = localDOF(local_x_prev);
580
581 auto local_Jac = MathLib::createZeroedMatrix<
582 typename ShapeMatricesTypeDisplacement::template MatrixType<
583 temperature_size + displacement_size + pressure_size,
584 temperature_size + displacement_size + pressure_size>>(
585 local_Jac_data, displacement_size + pressure_size + temperature_size,
586 displacement_size + pressure_size + temperature_size);
587
588 auto local_rhs = MathLib::createZeroedVector<
589 typename ShapeMatricesTypeDisplacement::template VectorType<
590 displacement_size + pressure_size + temperature_size>>(
591 local_rhs_data, displacement_size + pressure_size + temperature_size);
592
594 MTT.setZero(temperature_size, temperature_size);
595
597 KTT.setZero(temperature_size, temperature_size);
598
600 KTp.setZero(temperature_size, pressure_size);
601
603 dKTT_dp.setZero(temperature_size, pressure_size);
604
606 laplace_p.setZero(pressure_size, pressure_size);
607
609 laplace_T.setZero(pressure_size, temperature_size);
610
612 storage_p.setZero(pressure_size, pressure_size);
613
615 storage_T.setZero(pressure_size, temperature_size);
616
617 typename ShapeMatricesTypeDisplacement::template MatrixType<
618 displacement_size, pressure_size>
619 Kup;
620 Kup.setZero(displacement_size, pressure_size);
621
622 auto const& medium = _process_data.media_map.getMedium(_element.getID());
623 bool const has_frozen_liquid_phase = medium->hasPhase("FrozenLiquid");
624
625 unsigned const n_integration_points =
626 _integration_method.getNumberOfPoints();
627
628 std::vector<GlobalDimVectorType> ip_flux_vector;
629 double average_velocity_norm = 0.0;
630 ip_flux_vector.reserve(n_integration_points);
631
632 for (unsigned ip = 0; ip < n_integration_points; ip++)
633 {
634 auto const& N_u = _ip_data[ip].N_u;
635 ParameterLib::SpatialPosition const x_position{
636 std::nullopt, _element.getID(),
638 NumLib::interpolateCoordinates<ShapeFunctionDisplacement,
640 _element, N_u))};
641
642 auto const crv = updateConstitutiveRelations(
643 x, x_prev, x_position, t, dt, _ip_data[ip], _ip_data_output[ip]);
644
645 auto const& w = _ip_data[ip].integration_weight;
646
647 auto const& dNdx_u = _ip_data[ip].dNdx_u;
648
649 auto const& N = _ip_data[ip].N;
650 auto const& dNdx = _ip_data[ip].dNdx;
651
652 auto const T_int_pt = N.dot(T);
653
654 auto const x_coord =
655 NumLib::interpolateXCoordinate<ShapeFunctionDisplacement,
657 _element, N_u);
658 auto const B =
659 LinearBMatrix::computeBMatrix<DisplacementDim,
660 ShapeFunctionDisplacement::NPOINTS,
662 dNdx_u, N_u, x_coord, _is_axially_symmetric);
663
664 auto const& b = _process_data.specific_body_force;
665 auto const velocity = _ip_data_output[ip].velocity;
666
667 //
668 // displacement equation, displacement part
669 //
670
671 if (has_frozen_liquid_phase)
672 {
673 local_Jac
674 .template block<displacement_size, displacement_size>(
675 displacement_index, displacement_index)
676 .noalias() += B.transpose() * crv.J_uu_fr * B * w;
677
678 local_rhs.template segment<displacement_size>(displacement_index)
679 .noalias() -= B.transpose() * crv.r_u_fr * w;
680
681 local_Jac
682 .template block<displacement_size, temperature_size>(
683 displacement_index, temperature_index)
684 .noalias() -= B.transpose() * crv.J_uT_fr * N * w;
685 }
686
687 local_Jac
688 .template block<displacement_size, displacement_size>(
689 displacement_index, displacement_index)
690 .noalias() += B.transpose() * crv.C * B * w;
691
692 local_Jac
693 .template block<displacement_size, temperature_size>(
694 displacement_index, temperature_index)
695 .noalias() -= B.transpose() * crv.C *
696 crv.solid_linear_thermal_expansion_coefficient * N *
697 w;
698
699 local_rhs.template segment<displacement_size>(displacement_index)
700 .noalias() -= (B.transpose() * _ip_data[ip].sigma_eff -
701 N_u_op(N_u).transpose() * crv.rho * b) *
702 w;
703
704 //
705 // displacement equation, pressure part (K_up)
706 //
707 Kup.noalias() +=
708 B.transpose() * crv.alpha_biot * Invariants::identity2 * N * w;
709
710 //
711 // pressure equation, pressure part (K_pp and M_pp).
712 //
713 laplace_p.noalias() += dNdx.transpose() * crv.K_over_mu * dNdx * w;
714
715 storage_p.noalias() +=
716 N.transpose() *
717 (_ip_data[ip].porosity * crv.fluid_compressibility +
718 (crv.alpha_biot - _ip_data[ip].porosity) * crv.beta_SR) *
719 N * w;
720
721 laplace_T.noalias() +=
722 dNdx.transpose() * crv.K_pT_thermal_osmosis * dNdx * w;
723 //
724 // RHS, pressure part
725 //
726 double const fluid_density = _ip_data_output[ip].fluid_density;
727 local_rhs.template segment<pressure_size>(pressure_index).noalias() +=
728 dNdx.transpose() * fluid_density * crv.K_over_mu * b * w;
729 //
730 // pressure equation, temperature part (M_pT)
731 //
732 storage_T.noalias() += N.transpose() * crv.beta * N * w;
733
734 //
735 // pressure equation, displacement part.
736 //
737 // Reusing Kup.transpose().
738
739 //
740 // temperature equation, temperature part.
741 //
742 KTT.noalias() +=
743 dNdx.transpose() * crv.effective_thermal_conductivity * dNdx * w;
744
745 ip_flux_vector.emplace_back(velocity * fluid_density * crv.c_f);
746 average_velocity_norm += velocity.norm();
747
748 if (has_frozen_liquid_phase)
749 {
750 local_Jac
751 .template block<temperature_size, temperature_size>(
752 temperature_index, temperature_index)
753 .noalias() -= N.transpose() * crv.J_TT_fr * N * w;
754 }
755
756 MTT.noalias() +=
757 N.transpose() * crv.effective_volumetric_heat_capacity * N * w;
758
759 //
760 // temperature equation, pressure part
761 //
762 KTp.noalias() +=
763 dNdx.transpose() * T_int_pt * crv.K_pT_thermal_osmosis * dNdx * w;
764
765 // linearized darcy
766 dKTT_dp.noalias() -= fluid_density * crv.c_f * N.transpose() *
767 (dNdx * T).transpose() * crv.K_over_mu * dNdx * w;
768
769 /* TODO (Joerg) Temperature changes due to thermal dilatation of the
770 * fluid, which are usually discarded as being very small.
771 * Zhou et al. (10.1016/S0020-7683(98)00089-4) states that:
772 * "Biot (1956) neglected this term and it is included here for
773 * completeness"
774 * Keeping the code here in the case these are needed for the named
775 * effects in the future.
776 if (fluid_compressibility != 0)
777 {
778 auto const C_el = _ip_data[ip].computeElasticTangentStiffness(
779 t, x_position, dt, static_cast<double>(T_int_pt));
780 auto const solid_skeleton_compressibility =
781 1 / solid_material.getBulkModulus(t, x_position, &C_el);
782 double const fluid_volumetric_thermal_expansion_coefficient =
783 MaterialPropertyLib::getLiquidThermalExpansivity(
784 liquid_phase, vars, fluid_density, x_position, t, dt);
785
786 KTT.noalias() +=
787 dNdx.transpose() *
788 (-T_int_pt * fluid_volumetric_thermal_expansion_coefficient *
789 K_pT_thermal_osmosis / fluid_compressibility) *
790 dNdx * w;
791
792 local_rhs.template segment<temperature_size>(temperature_index)
793 .noalias() +=
794 dNdx.transpose() *
795 (-T_int_pt * fluid_volumetric_thermal_expansion_coefficient /
796 fluid_compressibility) *
797 fluid_density * K_over_mu * b * w;
798 MTu part for rhs and Jacobian:
799 (-T_int_pt *
800 Invariants::trace(solid_linear_thermal_expansion_coefficient) /
801 solid_skeleton_compressibility) *
802 N.transpose() * identity2.transpose() * B * w;
803 KTp part for rhs and Jacobian:
804 dNdx.transpose() *
805 (T_int_pt * fluid_volumetric_thermal_expansion_coefficient *
806 K_over_mu / fluid_compressibility) *
807 dNdx * w;
808 }
809 */
810 }
811
813 _process_data.stabilizer, _ip_data, ip_flux_vector,
814 average_velocity_norm / static_cast<double>(n_integration_points), KTT);
815
816 // temperature equation, temperature part
817 local_Jac
818 .template block<temperature_size, temperature_size>(temperature_index,
819 temperature_index)
820 .noalias() += KTT + MTT / dt;
821
822 // temperature equation, pressure part
823 local_Jac
824 .template block<temperature_size, pressure_size>(temperature_index,
825 pressure_index)
826 .noalias() += KTp + dKTT_dp;
827
828 // displacement equation, pressure part
829 local_Jac
830 .template block<displacement_size, pressure_size>(displacement_index,
831 pressure_index)
832 .noalias() -= Kup;
833
834 // pressure equation, temperature part.
835 local_Jac
836 .template block<pressure_size, temperature_size>(pressure_index,
837 temperature_index)
838 .noalias() += -storage_T / dt + laplace_T;
839
840 // pressure equation, pressure part.
841 local_Jac
842 .template block<pressure_size, pressure_size>(pressure_index,
843 pressure_index)
844 .noalias() += laplace_p + storage_p / dt;
845
846 // pressure equation, displacement part.
847 local_Jac
848 .template block<pressure_size, displacement_size>(pressure_index,
849 displacement_index)
850 .noalias() += Kup.transpose() / dt;
851
852 // pressure equation (f_p)
853 local_rhs.template segment<pressure_size>(pressure_index).noalias() -=
854 laplace_p * p + laplace_T * T + storage_p * (p - p_prev) / dt -
855 storage_T * (T - T_prev) / dt + Kup.transpose() * (u - u_prev) / dt;
856
857 // displacement equation (f_u)
858 local_rhs.template segment<displacement_size>(displacement_index)
859 .noalias() += Kup * p;
860
861 // temperature equation (f_T)
862 local_rhs.template segment<temperature_size>(temperature_index).noalias() -=
863 KTT * T + MTT * (T - T_prev) / dt;
864
865 local_rhs.template segment<temperature_size>(temperature_index).noalias() -=
866 KTp * p;
867}
868
869template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
870 int DisplacementDim>
871std::vector<double> const& ThermoHydroMechanicsLocalAssembler<
872 ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim>::
873 getIntPtDarcyVelocity(
874 const double /*t*/,
875 std::vector<GlobalVector*> const& /*x*/,
876 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
877 std::vector<double>& cache) const
878{
879 unsigned const n_integration_points =
880 _integration_method.getNumberOfPoints();
881
882 cache.clear();
883 auto cache_matrix = MathLib::createZeroedMatrix<Eigen::Matrix<
884 double, DisplacementDim, Eigen::Dynamic, Eigen::RowMajor>>(
885 cache, DisplacementDim, n_integration_points);
886
887 for (unsigned ip = 0; ip < n_integration_points; ip++)
888 {
889 cache_matrix.col(ip).noalias() = _ip_data_output[ip].velocity;
890 }
891
892 return cache;
893}
894
895template <typename ShapeFunctionDisplacement, typename ShapeFunction,
896 int DisplacementDim>
897std::vector<double> const& ThermoHydroMechanicsLocalAssembler<
898 ShapeFunctionDisplacement, ShapeFunction, DisplacementDim>::
899 getIntPtFluidDensity(
900 const double /*t*/,
901 std::vector<GlobalVector*> const& /*x*/,
902 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
903 std::vector<double>& cache) const
904{
906 _ip_data_output,
908}
909
910template <typename ShapeFunctionDisplacement, typename ShapeFunction,
911 int DisplacementDim>
912std::vector<double> const& ThermoHydroMechanicsLocalAssembler<
913 ShapeFunctionDisplacement, ShapeFunction, DisplacementDim>::
914 getIntPtViscosity(
915 const double /*t*/,
916 std::vector<GlobalVector*> const& /*x*/,
917 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
918 std::vector<double>& cache) const
919{
921 _ip_data_output,
923}
924
925template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
926 int DisplacementDim>
928 ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim>::
929 computeSecondaryVariableConcrete(double const /*t*/, double const /*dt*/,
930 Eigen::VectorXd const& local_x,
931 Eigen::VectorXd const& /*local_x_prev*/)
932{
933 auto const p = local_x.template segment<pressure_size>(pressure_index);
934 auto const T =
935 local_x.template segment<temperature_size>(temperature_index);
936
937 unsigned const n_integration_points =
938 _integration_method.getNumberOfPoints();
939
940 double phi_fr_avg = 0;
941 double fluid_density_avg = 0;
942 double viscosity_avg = 0;
943
945 KV sigma_avg = KV::Zero();
946
947 for (unsigned ip = 0; ip < n_integration_points; ip++)
948 {
949 auto& ip_data = _ip_data[ip];
950
951 phi_fr_avg += _ip_data[ip].phi_fr;
952 fluid_density_avg += _ip_data_output[ip].fluid_density;
953 viscosity_avg += _ip_data_output[ip].viscosity;
954 sigma_avg += ip_data.sigma_eff;
955 }
956
957 phi_fr_avg /= n_integration_points;
958 fluid_density_avg /= n_integration_points;
959 viscosity_avg /= n_integration_points;
960 sigma_avg /= n_integration_points;
961
962 (*_process_data.element_phi_fr)[_element.getID()] = phi_fr_avg;
963 (*_process_data.element_fluid_density)[_element.getID()] =
964 fluid_density_avg;
965 (*_process_data.element_viscosity)[_element.getID()] = viscosity_avg;
966
967 Eigen::Map<KV>(&(*_process_data.element_stresses)[_element.getID() *
968 KV::RowsAtCompileTime]) =
970
972 ShapeFunctionPressure, typename ShapeFunctionDisplacement::MeshElement,
973 DisplacementDim>(_element, _is_axially_symmetric, p,
974 *_process_data.pressure_interpolated);
975
977 ShapeFunctionPressure, typename ShapeFunctionDisplacement::MeshElement,
978 DisplacementDim>(_element, _is_axially_symmetric, T,
979 *_process_data.temperature_interpolated);
980}
981} // namespace ThermoHydroMechanics
982} // namespace ProcessLib
#define OGS_FATAL(...)
Definition Error.h:26
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
constexpr double getWeight() const
std::size_t getID() const
Returns the ID of the element.
Definition Element.h:91
MathLib::WeightedPoint const & getWeightedPoint(unsigned const igp) const
MatrixType< _kelvin_vector_size, _number_of_dof > BMatrixType
ShapeMatrixPolicyType< ShapeFunctionDisplacement, DisplacementDim > ShapeMatricesTypeDisplacement
typename ShapeMatricesTypePressure::GlobalDimVectorType GlobalDimVectorType
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)
MathLib::KelvinVector::KelvinVectorType< GlobalDim > formKelvinVector(MaterialPropertyLib::PropertyDataType const &values)
A function to form a Kelvin vector from strain or stress alike property like thermal expansivity for ...
Eigen::Matrix< double, GlobalDim, GlobalDim > formEigenTensor(MaterialPropertyLib::PropertyDataType const &values)
double getLiquidThermalExpansivity(Phase const &phase, VariableArray const &vars, const double density, ParameterLib::SpatialPosition const &pos, double const t, double const dt)
Eigen::Matrix< double, 4, 1 > kelvinVectorToSymmetricTensor(Eigen::Matrix< double, 4, 1, Eigen::ColMajor, 4, 1 > const &v)
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), 1, Eigen::ColMajor > KelvinVectorType
Eigen::Map< Vector > createZeroedVector(std::vector< double > &data, Eigen::VectorXd::Index size)
Eigen::Map< Matrix > createZeroedMatrix(std::vector< double > &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
void assembleAdvectionMatrix(IPData const &ip_data_vector, NumLib::ShapeMatrixCache const &shape_matrix_cache, std::vector< FluxVectorType > const &ip_flux_vector, Eigen::MatrixBase< Derived > &laplacian_matrix)
double interpolateXCoordinate(MeshLib::Element const &e, typename ShapeMatricesType::ShapeMatrices::ShapeType const &N)
void interpolateToHigherOrderNodes(MeshLib::Element const &element, bool const is_axially_symmetric, Eigen::MatrixBase< EigenMatrixType > const &node_values, MeshLib::PropertyVector< double > &interpolated_values_global_vector)
std::vector< typename ShapeMatricesType::ShapeMatrices, Eigen::aligned_allocator< typename ShapeMatricesType::ShapeMatrices > > initShapeMatrices(MeshLib::Element const &e, bool const is_axially_symmetric, IntegrationMethod const &integration_method)
Eigen::MatrixXd computeHydrodynamicDispersion(NumericalStabilization const &stabilizer, std::size_t const element_id, Eigen::MatrixXd const &pore_diffusion_coefficient, Eigen::VectorXd const &velocity, double const porosity, double const solute_dispersivity_transverse, double const solute_dispersivity_longitudinal)
std::array< double, 3 > interpolateCoordinates(MeshLib::Element const &e, typename ShapeMatricesType::ShapeMatrices::ShapeType const &N)
BMatrixType computeBMatrix(DNDX_Type const &dNdx, N_Type const &N, const double radius, const bool is_axially_symmetric)
Fills a B-matrix based on given shape function dN/dx values.
std::vector< double > const & getIntegrationPointScalarData(IntegrationPointDataVector const &ip_data_vector, MemberType IpData::*const member, std::vector< double > &cache)
std::size_t setIntegrationPointDataMaterialStateVariables(double const *values, IntegrationPointDataVector &ip_data_vector, MemberType member, std::function< std::span< double >(MaterialStateVariables &)> get_values_span)
std::size_t setIntegrationPointKelvinVectorData(double const *values, IntegrationPointDataVector &ip_data_vector, MemberType IpData::*const member)
void setIPDataInitialConditions(std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > const &_integration_point_writer, MeshLib::Properties const &mesh_properties, LocalAssemblersVector &local_assemblers)
MatrixType< ShapeFunction::NPOINTS, ShapeFunction::NPOINTS > NodalMatrixType
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
std::unique_ptr< typename MaterialLib::Solids::MechanicsBase< DisplacementDim >::MaterialStateVariables > material_state_variables
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)
ShapeMatricesTypePressure::NodalRowVectorType N
ShapeMatrixTypeDisplacement::GlobalDimNodalMatrixType dNdx_u
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N_u