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 if (!_process_data.initial_stress.isTotalStress())
185 {
186 return;
187 }
188
189 // TODO: For staggered scheme, overload
190 // LocalAssemblerInterface::setInitialConditions to enable local_x contains
191 // the primary variables from all coupled processes.
192 auto const p = local_x.template segment<pressure_size>(pressure_index);
193
194 double const dt = 0.0;
195
197 auto const& medium = _process_data.media_map.getMedium(_element.getID());
198
199 int const n_integration_points = _integration_method.getNumberOfPoints();
200 for (int ip = 0; ip < n_integration_points; ip++)
201 {
202 auto const& N = _ip_data[ip].N;
203 auto const& N_u = _ip_data[ip].N_u;
204 ParameterLib::SpatialPosition const x_position{
205 std::nullopt, _element.getID(), ip,
207 NumLib::interpolateCoordinates<ShapeFunctionDisplacement,
209 _element, N_u))};
210 auto const alpha_b =
211 medium->property(MPL::PropertyType::biot_coefficient)
212 .template value<double>(vars, x_position, t, dt);
213
214 auto& sigma_eff = _ip_data[ip].sigma_eff;
215 sigma_eff.noalias() += alpha_b * N.dot(p) * Invariants::identity2;
216 _ip_data[ip].sigma_eff_prev.noalias() = sigma_eff;
217 }
218}
219
220template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
221 int DisplacementDim>
223 ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim>::
224 updateConstitutiveRelations(
225 Eigen::Ref<Eigen::VectorXd const> const local_x,
226 Eigen::Ref<Eigen::VectorXd const> const local_x_prev,
227 ParameterLib::SpatialPosition const& x_position, double const t,
228 double const dt, IpData& ip_data,
230{
231 assert(local_x.size() ==
232 pressure_size + displacement_size + temperature_size);
233
234 auto const [T, p, u] = localDOF(local_x);
235 auto const [T_prev, p_prev, u_prev] = localDOF(local_x_prev);
236
237 auto const& solid_material =
239 _process_data.solid_materials, _process_data.material_ids,
240 _element.getID());
241
242 auto const& medium = _process_data.media_map.getMedium(_element.getID());
243 auto const& liquid_phase = medium->phase("AqueousLiquid");
244 auto const& solid_phase = medium->phase("Solid");
245 auto* const frozen_liquid_phase = medium->hasPhase("FrozenLiquid")
246 ? &medium->phase("FrozenLiquid")
247 : nullptr;
249
250 auto const& N_u = ip_data.N_u;
251 auto const& dNdx_u = ip_data.dNdx_u;
252
253 auto const& N = ip_data.N;
254 auto const& dNdx = ip_data.dNdx;
255
256 auto const T_int_pt = N.dot(T);
257 auto const T_prev_int_pt = N.dot(T_prev);
258 double const dT_int_pt = T_int_pt - T_prev_int_pt;
259
260 auto const x_coord =
261 NumLib::interpolateXCoordinate<ShapeFunctionDisplacement,
263 N_u);
264 auto const B =
265 LinearBMatrix::computeBMatrix<DisplacementDim,
266 ShapeFunctionDisplacement::NPOINTS,
268 dNdx_u, N_u, x_coord, _is_axially_symmetric);
269
271
272 auto& eps = ip_data.eps;
273 eps.noalias() = B * u;
275 B * u_prev;
276
277 vars.temperature = T_int_pt;
278 double const p_int_pt = N.dot(p);
279 vars.liquid_phase_pressure = p_int_pt;
280
281 vars.liquid_saturation = 1.0;
282
283 auto const solid_density =
285 .template value<double>(vars, x_position, t, dt);
286
287 auto const porosity =
289 .template value<double>(vars, x_position, t, dt);
290 vars.porosity = porosity;
291 ip_data.porosity = porosity;
292
293 crv.alpha_biot =
295 .template value<double>(vars, x_position, t, dt);
296 auto const& alpha = crv.alpha_biot;
297
298 auto const C_el = ip_data.computeElasticTangentStiffness(
299 t, x_position, dt, static_cast<double>(T_int_pt));
300 auto const solid_skeleton_compressibility =
301 1 / solid_material.getBulkModulus(t, x_position, &C_el);
302
303 crv.beta_SR = (1 - alpha) * solid_skeleton_compressibility;
304
305 // Set mechanical variables for the intrinsic permeability model
306 // For stress dependent permeability.
307 {
308 auto const& identity2 = Invariants::identity2;
309 auto const sigma_total =
310 (ip_data.sigma_eff - alpha * p_int_pt * identity2).eval();
311 vars.total_stress.emplace<SymmetricTensor>(
313 }
314 // For strain dependent permeability
315 vars.volumetric_strain = Invariants::trace(ip_data.eps);
317 ip_data.material_state_variables->getEquivalentPlasticStrain();
318
319 auto const intrinsic_permeability =
322 .value(vars, x_position, t, dt));
323
324 auto const fluid_density =
325 liquid_phase.property(MaterialPropertyLib::PropertyType::density)
326 .template value<double>(vars, x_position, t, dt);
327 ip_data_output.fluid_density = fluid_density;
328 vars.density = fluid_density;
329
330 auto const drho_dp =
331 liquid_phase.property(MaterialPropertyLib::PropertyType::density)
332 .template dValue<double>(
334 x_position, t, dt);
335
336 crv.fluid_compressibility = 1 / fluid_density * drho_dp;
337
338 double const fluid_volumetric_thermal_expansion_coefficient =
340 liquid_phase, vars, fluid_density, x_position, t, dt);
341
342 // Use the viscosity model to compute the viscosity
343 ip_data_output.viscosity =
345 .template value<double>(vars, x_position, t, dt);
346 crv.K_over_mu = intrinsic_permeability / ip_data_output.viscosity;
347
348 auto const& b = _process_data.specific_body_force;
349
350 // Consider also anisotropic thermal expansion.
351 crv.solid_linear_thermal_expansion_coefficient =
353 solid_phase
354 .property(
356 .value(vars, x_position, t, dt));
357
359 dthermal_strain =
360 crv.solid_linear_thermal_expansion_coefficient * dT_int_pt;
361
362 crv.K_pT_thermal_osmosis =
363 (solid_phase.hasProperty(
366 solid_phase
368 thermal_osmosis_coefficient)
369 .value(vars, x_position, t, dt))
370 : Eigen::MatrixXd::Zero(DisplacementDim, DisplacementDim));
371
372 GlobalDimVectorType const velocity = -crv.K_over_mu * dNdx * p -
373 crv.K_pT_thermal_osmosis * dNdx * T +
374 crv.K_over_mu * fluid_density * b;
375 ip_data_output.velocity = velocity;
376
377 //
378 // displacement equation, displacement part
379 //
380 auto& eps_m = ip_data.eps_m;
381 auto& eps_m_prev = ip_data.eps_m_prev;
382 eps_m.noalias() = eps_m_prev + eps - eps_prev - dthermal_strain;
385 eps_m);
386
387 crv.C = ip_data.updateConstitutiveRelation(vars, t, x_position, dt,
388 T_prev_int_pt);
389
390 crv.rho = solid_density * (1 - porosity) + porosity * fluid_density;
391
392 crv.beta =
393 porosity * fluid_volumetric_thermal_expansion_coefficient +
394 (alpha - porosity) *
395 Invariants::trace(crv.solid_linear_thermal_expansion_coefficient);
396
397 //
398 // pressure equation, displacement part.
399 //
400 // Reusing Kup.transpose().
401
402 //
403 // temperature equation, temperature part.
404 //
405 crv.c_f =
406 liquid_phase
408 .template value<double>(vars, x_position, t, dt);
409 crv.effective_thermal_conductivity =
411 medium
412 ->property(
414 .value(vars, x_position, t, dt));
415
416 // Thermal conductivity is moved outside and zero matrix is passed instead
417 // due to multiplication with fluid's density times specific heat capacity.
418 crv.effective_thermal_conductivity.noalias() +=
419 fluid_density * crv.c_f *
421 _process_data.stabilizer, _element.getID(),
422 GlobalDimMatrixType::Zero(DisplacementDim, DisplacementDim),
423 velocity, 0. /* phi */, 0. /* dispersivity_transversal */,
424 0. /*dispersivity_longitudinal*/);
425
426 double const c_s =
427 solid_phase
429 .template value<double>(vars, x_position, t, dt);
430
431 // Also modified by freezing terms.
432 crv.effective_volumetric_heat_capacity =
433 porosity * fluid_density * crv.c_f +
434 (1.0 - porosity) * solid_density * c_s;
435
436 if (frozen_liquid_phase)
437 {
439 double const phi_fr =
441 .template value<double>(vars, x_position, t, dt);
442 ip_data.phi_fr = phi_fr;
443
444 auto const frozen_liquid_value =
446 {
447 return (*frozen_liquid_phase)[p].template value<double>(
448 vars, x_position, t, dt);
449 };
450
451 double const rho_fr =
453
454 double const c_fr = frozen_liquid_value(
456
457 double const l_fr = frozen_liquid_value(
459
460 double const dphi_fr_dT =
462 .template dValue<double>(
464 x_position, t, dt);
465
466 double const phi_fr_prev = [&]()
467 {
469 vars_prev.temperature = T_prev_int_pt;
471 .template value<double>(vars_prev, x_position, t, dt);
472 }();
473 ip_data.phi_fr_prev = phi_fr_prev;
474
475 // alpha_T^I
477 DisplacementDim> const ice_linear_thermal_expansion_coefficient =
479 frozen_liquid_phase
480 ->property(
482 .value(vars, x_position, t, dt));
483
485 dthermal_strain_ice =
486 ice_linear_thermal_expansion_coefficient * dT_int_pt;
487
488 // alpha_{phi_I} -- linear expansion coeff. due to water-to-ice
489 // transition (phase change), and related phase_change_strain term
491 phase_change_expansion_coefficient =
493 frozen_liquid_phase
495 phase_change_expansivity)
496 .value(vars, x_position, t, dt));
497
499 dphase_change_strain = phase_change_expansion_coefficient *
500 (phi_fr - phi_fr_prev) / porosity;
501
502 // eps0 ia a 'history variable' -- a solid matrix strain accrued
503 // prior to the onset of ice forming
504 auto& eps0 = ip_data.eps0;
505 auto const& eps0_prev = ip_data.eps0_prev;
506
507 // definition of eps_m_ice
508 auto& eps_m_ice = ip_data.eps_m_ice;
509 auto const& eps_m_ice_prev = ip_data.eps_m_ice_prev;
510
511 eps_m_ice.noalias() = eps_m_ice_prev + eps - eps_prev -
512 (eps0 - eps0_prev) - dthermal_strain_ice -
513 dphase_change_strain;
514
515 vars_ice.mechanical_strain
517 eps_m_ice);
518 auto const C_IR = ip_data.updateConstitutiveRelationIce(
519 *_process_data.ice_constitutive_relation, vars_ice, t, x_position,
520 dt, T_prev_int_pt);
521 crv.effective_volumetric_heat_capacity +=
522 -phi_fr * fluid_density * crv.c_f + phi_fr * rho_fr * c_fr -
523 l_fr * rho_fr * dphi_fr_dT;
524
525 // part of dMTT_dT derivative for freezing
526 double const d2phi_fr_dT2 =
528 .template d2Value<double>(
531 dt);
532
533 crv.J_uu_fr = phi_fr * C_IR;
534
535 auto const& sigma_eff_ice = ip_data.sigma_eff_ice;
536 crv.r_u_fr = phi_fr * sigma_eff_ice;
537
538 crv.J_uT_fr = phi_fr * C_IR * ice_linear_thermal_expansion_coefficient;
539
540 crv.J_TT_fr = ((rho_fr * c_fr - fluid_density * crv.c_f) * dphi_fr_dT +
541 l_fr * rho_fr * d2phi_fr_dT2) *
542 dT_int_pt / dt;
543 }
544 return crv;
545}
546
547// Assembles the local Jacobian matrix. So far, the linearisation of HT part is
548// not considered as that in HT process.
549template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
550 int DisplacementDim>
552 ShapeFunctionDisplacement, ShapeFunctionPressure,
553 DisplacementDim>::assembleWithJacobian(double const t, double const dt,
554 std::vector<double> const& local_x,
555 std::vector<double> const&
556 local_x_prev,
557 std::vector<double>& local_rhs_data,
558 std::vector<double>& local_Jac_data)
559{
560 assert(local_x.size() ==
561 pressure_size + displacement_size + temperature_size);
562
563 auto const x =
564 Eigen::Map<Eigen::VectorXd const>(local_x.data(), local_x.size());
565 auto const x_prev = Eigen::Map<Eigen::VectorXd const>(local_x_prev.data(),
566 local_x_prev.size());
567
568 auto const [T, p, u] = localDOF(local_x);
569 auto const [T_prev, p_prev, u_prev] = localDOF(local_x_prev);
570
571 auto local_Jac = MathLib::createZeroedMatrix<
572 typename ShapeMatricesTypeDisplacement::template MatrixType<
573 temperature_size + displacement_size + pressure_size,
574 temperature_size + displacement_size + pressure_size>>(
575 local_Jac_data, displacement_size + pressure_size + temperature_size,
576 displacement_size + pressure_size + temperature_size);
577
578 auto local_rhs = MathLib::createZeroedVector<
579 typename ShapeMatricesTypeDisplacement::template VectorType<
580 displacement_size + pressure_size + temperature_size>>(
581 local_rhs_data, displacement_size + pressure_size + temperature_size);
582
584 MTT.setZero(temperature_size, temperature_size);
585
587 KTT.setZero(temperature_size, temperature_size);
588
590 KTp.setZero(temperature_size, pressure_size);
591
593 dKTT_dp.setZero(temperature_size, pressure_size);
594
596 laplace_p.setZero(pressure_size, pressure_size);
597
599 laplace_T.setZero(pressure_size, temperature_size);
600
602 storage_p.setZero(pressure_size, pressure_size);
603
605 storage_T.setZero(pressure_size, temperature_size);
606
607 typename ShapeMatricesTypeDisplacement::template MatrixType<
608 displacement_size, pressure_size>
609 Kup;
610 Kup.setZero(displacement_size, pressure_size);
611
612 auto const& medium = _process_data.media_map.getMedium(_element.getID());
613 bool const has_frozen_liquid_phase = medium->hasPhase("FrozenLiquid");
614
615 unsigned const n_integration_points =
616 _integration_method.getNumberOfPoints();
617
618 std::vector<GlobalDimVectorType> ip_flux_vector;
619 double average_velocity_norm = 0.0;
620 ip_flux_vector.reserve(n_integration_points);
621
622 for (unsigned ip = 0; ip < n_integration_points; ip++)
623 {
624 auto const& N_u = _ip_data[ip].N_u;
625 ParameterLib::SpatialPosition const x_position{
626 std::nullopt, _element.getID(), ip,
628 NumLib::interpolateCoordinates<ShapeFunctionDisplacement,
630 _element, N_u))};
631
632 auto const crv = updateConstitutiveRelations(
633 x, x_prev, x_position, t, dt, _ip_data[ip], _ip_data_output[ip]);
634
635 auto const& w = _ip_data[ip].integration_weight;
636
637 auto const& dNdx_u = _ip_data[ip].dNdx_u;
638
639 auto const& N = _ip_data[ip].N;
640 auto const& dNdx = _ip_data[ip].dNdx;
641
642 auto const T_int_pt = N.dot(T);
643
644 auto const x_coord =
645 NumLib::interpolateXCoordinate<ShapeFunctionDisplacement,
647 _element, N_u);
648 auto const B =
649 LinearBMatrix::computeBMatrix<DisplacementDim,
650 ShapeFunctionDisplacement::NPOINTS,
652 dNdx_u, N_u, x_coord, _is_axially_symmetric);
653
654 auto const& b = _process_data.specific_body_force;
655 auto const velocity = _ip_data_output[ip].velocity;
656
657 //
658 // displacement equation, displacement part
659 //
660
661 if (has_frozen_liquid_phase)
662 {
663 local_Jac
664 .template block<displacement_size, displacement_size>(
665 displacement_index, displacement_index)
666 .noalias() += B.transpose() * crv.J_uu_fr * B * w;
667
668 local_rhs.template segment<displacement_size>(displacement_index)
669 .noalias() -= B.transpose() * crv.r_u_fr * w;
670
671 local_Jac
672 .template block<displacement_size, temperature_size>(
673 displacement_index, temperature_index)
674 .noalias() -= B.transpose() * crv.J_uT_fr * N * w;
675 }
676
677 local_Jac
678 .template block<displacement_size, displacement_size>(
679 displacement_index, displacement_index)
680 .noalias() += B.transpose() * crv.C * B * w;
681
682 local_Jac
683 .template block<displacement_size, temperature_size>(
684 displacement_index, temperature_index)
685 .noalias() -= B.transpose() * crv.C *
686 crv.solid_linear_thermal_expansion_coefficient * N *
687 w;
688
689 local_rhs.template segment<displacement_size>(displacement_index)
690 .noalias() -= (B.transpose() * _ip_data[ip].sigma_eff -
691 N_u_op(N_u).transpose() * crv.rho * b) *
692 w;
693
694 //
695 // displacement equation, pressure part (K_up)
696 //
697 Kup.noalias() +=
698 B.transpose() * crv.alpha_biot * Invariants::identity2 * N * w;
699
700 //
701 // pressure equation, pressure part (K_pp and M_pp).
702 //
703 laplace_p.noalias() += dNdx.transpose() * crv.K_over_mu * dNdx * w;
704
705 storage_p.noalias() +=
706 N.transpose() *
707 (_ip_data[ip].porosity * crv.fluid_compressibility +
708 (crv.alpha_biot - _ip_data[ip].porosity) * crv.beta_SR) *
709 N * w;
710
711 laplace_T.noalias() +=
712 dNdx.transpose() * crv.K_pT_thermal_osmosis * dNdx * w;
713 //
714 // RHS, pressure part
715 //
716 double const fluid_density = _ip_data_output[ip].fluid_density;
717 local_rhs.template segment<pressure_size>(pressure_index).noalias() +=
718 dNdx.transpose() * fluid_density * crv.K_over_mu * b * w;
719 //
720 // pressure equation, temperature part (M_pT)
721 //
722 storage_T.noalias() += N.transpose() * crv.beta * N * w;
723
724 //
725 // pressure equation, displacement part.
726 //
727 // Reusing Kup.transpose().
728
729 //
730 // temperature equation, temperature part.
731 //
732 KTT.noalias() +=
733 dNdx.transpose() * crv.effective_thermal_conductivity * dNdx * w;
734
735 ip_flux_vector.emplace_back(velocity * fluid_density * crv.c_f);
736 average_velocity_norm += velocity.norm();
737
738 if (has_frozen_liquid_phase)
739 {
740 local_Jac
741 .template block<temperature_size, temperature_size>(
742 temperature_index, temperature_index)
743 .noalias() -= N.transpose() * crv.J_TT_fr * N * w;
744 }
745
746 MTT.noalias() +=
747 N.transpose() * crv.effective_volumetric_heat_capacity * N * w;
748
749 //
750 // temperature equation, pressure part
751 //
752 KTp.noalias() +=
753 dNdx.transpose() * T_int_pt * crv.K_pT_thermal_osmosis * dNdx * w;
754
755 // linearized darcy
756 dKTT_dp.noalias() -= fluid_density * crv.c_f * N.transpose() *
757 (dNdx * T).transpose() * crv.K_over_mu * dNdx * w;
758
759 /* TODO (Joerg) Temperature changes due to thermal dilatation of the
760 * fluid, which are usually discarded as being very small.
761 * Zhou et al. (10.1016/S0020-7683(98)00089-4) states that:
762 * "Biot (1956) neglected this term and it is included here for
763 * completeness"
764 * Keeping the code here in the case these are needed for the named
765 * effects in the future.
766 if (fluid_compressibility != 0)
767 {
768 auto const C_el = _ip_data[ip].computeElasticTangentStiffness(
769 t, x_position, dt, static_cast<double>(T_int_pt));
770 auto const solid_skeleton_compressibility =
771 1 / solid_material.getBulkModulus(t, x_position, &C_el);
772 double const fluid_volumetric_thermal_expansion_coefficient =
773 MaterialPropertyLib::getLiquidThermalExpansivity(
774 liquid_phase, vars, fluid_density, x_position, t, dt);
775
776 KTT.noalias() +=
777 dNdx.transpose() *
778 (-T_int_pt * fluid_volumetric_thermal_expansion_coefficient *
779 K_pT_thermal_osmosis / fluid_compressibility) *
780 dNdx * w;
781
782 local_rhs.template segment<temperature_size>(temperature_index)
783 .noalias() +=
784 dNdx.transpose() *
785 (-T_int_pt * fluid_volumetric_thermal_expansion_coefficient /
786 fluid_compressibility) *
787 fluid_density * K_over_mu * b * w;
788 MTu part for rhs and Jacobian:
789 (-T_int_pt *
790 Invariants::trace(solid_linear_thermal_expansion_coefficient) /
791 solid_skeleton_compressibility) *
792 N.transpose() * identity2.transpose() * B * w;
793 KTp part for rhs and Jacobian:
794 dNdx.transpose() *
795 (T_int_pt * fluid_volumetric_thermal_expansion_coefficient *
796 K_over_mu / fluid_compressibility) *
797 dNdx * w;
798 }
799 */
800 }
801
803 _process_data.stabilizer, _ip_data, ip_flux_vector,
804 average_velocity_norm / static_cast<double>(n_integration_points), KTT);
805
806 // temperature equation, temperature part
807 local_Jac
808 .template block<temperature_size, temperature_size>(temperature_index,
809 temperature_index)
810 .noalias() += KTT + MTT / dt;
811
812 // temperature equation, pressure part
813 local_Jac
814 .template block<temperature_size, pressure_size>(temperature_index,
815 pressure_index)
816 .noalias() += KTp + dKTT_dp;
817
818 // displacement equation, pressure part
819 local_Jac
820 .template block<displacement_size, pressure_size>(displacement_index,
821 pressure_index)
822 .noalias() -= Kup;
823
824 // pressure equation, temperature part.
825 local_Jac
826 .template block<pressure_size, temperature_size>(pressure_index,
827 temperature_index)
828 .noalias() -= storage_T / dt - laplace_T;
829
830 // pressure equation, pressure part.
831 local_Jac
832 .template block<pressure_size, pressure_size>(pressure_index,
833 pressure_index)
834 .noalias() += laplace_p + storage_p / dt;
835
836 // pressure equation, displacement part.
837 local_Jac
838 .template block<pressure_size, displacement_size>(pressure_index,
839 displacement_index)
840 .noalias() += Kup.transpose() / dt;
841
842 // pressure equation (f_p)
843 local_rhs.template segment<pressure_size>(pressure_index).noalias() -=
844 laplace_p * p + laplace_T * T + storage_p * (p - p_prev) / dt -
845 storage_T * (T - T_prev) / dt + Kup.transpose() * (u - u_prev) / dt;
846
847 // displacement equation (f_u)
848 local_rhs.template segment<displacement_size>(displacement_index)
849 .noalias() += Kup * p;
850
851 // temperature equation (f_T)
852 local_rhs.template segment<temperature_size>(temperature_index).noalias() -=
853 KTT * T + MTT * (T - T_prev) / dt;
854
855 local_rhs.template segment<temperature_size>(temperature_index).noalias() -=
856 KTp * p;
857}
858
859template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
860 int DisplacementDim>
861std::vector<double> const& ThermoHydroMechanicsLocalAssembler<
862 ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim>::
863 getIntPtDarcyVelocity(
864 const double /*t*/,
865 std::vector<GlobalVector*> const& /*x*/,
866 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
867 std::vector<double>& cache) const
868{
869 unsigned const n_integration_points =
870 _integration_method.getNumberOfPoints();
871
872 cache.clear();
873 auto cache_matrix = MathLib::createZeroedMatrix<Eigen::Matrix<
874 double, DisplacementDim, Eigen::Dynamic, Eigen::RowMajor>>(
875 cache, DisplacementDim, n_integration_points);
876
877 for (unsigned ip = 0; ip < n_integration_points; ip++)
878 {
879 cache_matrix.col(ip).noalias() = _ip_data_output[ip].velocity;
880 }
881
882 return cache;
883}
884
885template <typename ShapeFunctionDisplacement, typename ShapeFunction,
886 int DisplacementDim>
887std::vector<double> const& ThermoHydroMechanicsLocalAssembler<
888 ShapeFunctionDisplacement, ShapeFunction, DisplacementDim>::
889 getIntPtFluidDensity(
890 const double /*t*/,
891 std::vector<GlobalVector*> const& /*x*/,
892 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
893 std::vector<double>& cache) const
894{
896 _ip_data_output,
898}
899
900template <typename ShapeFunctionDisplacement, typename ShapeFunction,
901 int DisplacementDim>
902std::vector<double> const& ThermoHydroMechanicsLocalAssembler<
903 ShapeFunctionDisplacement, ShapeFunction, DisplacementDim>::
904 getIntPtViscosity(
905 const double /*t*/,
906 std::vector<GlobalVector*> const& /*x*/,
907 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
908 std::vector<double>& cache) const
909{
911 _ip_data_output,
913}
914
915template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
916 int DisplacementDim>
918 ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim>::
919 computeSecondaryVariableConcrete(double const /*t*/, double const /*dt*/,
920 Eigen::VectorXd const& local_x,
921 Eigen::VectorXd const& /*local_x_prev*/)
922{
923 auto const p = local_x.template segment<pressure_size>(pressure_index);
924 auto const T =
925 local_x.template segment<temperature_size>(temperature_index);
926
927 unsigned const n_integration_points =
928 _integration_method.getNumberOfPoints();
929
930 double fluid_density_avg = 0;
931 double viscosity_avg = 0;
932
934 KV sigma_avg = KV::Zero();
935
936 for (unsigned ip = 0; ip < n_integration_points; ip++)
937 {
938 auto& ip_data = _ip_data[ip];
939
940 fluid_density_avg += _ip_data_output[ip].fluid_density;
941 viscosity_avg += _ip_data_output[ip].viscosity;
942 sigma_avg += ip_data.sigma_eff;
943 }
944
945 fluid_density_avg /= n_integration_points;
946 viscosity_avg /= n_integration_points;
947 sigma_avg /= n_integration_points;
948
949 (*_process_data.element_fluid_density)[_element.getID()] =
950 fluid_density_avg;
951 (*_process_data.element_viscosity)[_element.getID()] = viscosity_avg;
952
953 Eigen::Map<KV>(&(*_process_data.element_stresses)[_element.getID() *
954 KV::RowsAtCompileTime]) =
956
958 ShapeFunctionPressure, typename ShapeFunctionDisplacement::MeshElement,
959 DisplacementDim>(_element, _is_axially_symmetric, p,
960 *_process_data.pressure_interpolated);
961
963 ShapeFunctionPressure, typename ShapeFunctionDisplacement::MeshElement,
964 DisplacementDim>(_element, _is_axially_symmetric, T,
965 *_process_data.temperature_interpolated);
966}
967} // namespace ThermoHydroMechanics
968} // namespace ProcessLib
#define OGS_FATAL(...)
Definition Error.h:26
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
double getWeight() const
std::size_t getID() const
Returns the ID of the element.
Definition Element.h:89
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