Loading [MathJax]/extensions/tex2jax.js
OGS
TH2MFEM-impl.h
Go to the documentation of this file.
1
11#pragma once
12
13#include <Eigen/LU>
14
25#include "TH2MProcessData.h"
26
27namespace ProcessLib
28{
29namespace TH2M
30{
31namespace MPL = MaterialPropertyLib;
32
33template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
34 int DisplacementDim>
35TH2MLocalAssembler<ShapeFunctionDisplacement, ShapeFunctionPressure,
36 DisplacementDim>::
37 TH2MLocalAssembler(
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 : LocalAssemblerInterface<DisplacementDim>(
44 e, integration_method, is_axially_symmetric, process_data)
45{
46 unsigned const n_integration_points =
48
49 _ip_data.resize(n_integration_points);
50 _secondary_data.N_u.resize(n_integration_points);
51
52 auto const shape_matrices_u =
53 NumLib::initShapeMatrices<ShapeFunctionDisplacement,
55 DisplacementDim>(e, is_axially_symmetric,
57
58 auto const shape_matrices_p =
59 NumLib::initShapeMatrices<ShapeFunctionPressure,
60 ShapeMatricesTypePressure, DisplacementDim>(
61 e, is_axially_symmetric, this->integration_method_);
62
63 for (unsigned ip = 0; ip < n_integration_points; ip++)
64 {
65 auto& ip_data = _ip_data[ip];
66 auto const& sm_u = shape_matrices_u[ip];
67 ip_data.integration_weight =
69 sm_u.integralMeasure * sm_u.detJ;
70
71 ip_data.N_u = sm_u.N;
72 ip_data.dNdx_u = sm_u.dNdx;
73
74 ip_data.N_p = shape_matrices_p[ip].N;
75 ip_data.dNdx_p = shape_matrices_p[ip].dNdx;
76
77 _secondary_data.N_u[ip] = shape_matrices_u[ip].N;
78 }
79}
80
81template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
82 int DisplacementDim>
83std::tuple<
84 std::vector<ConstitutiveRelations::ConstitutiveData<DisplacementDim>>,
85 std::vector<ConstitutiveRelations::ConstitutiveTempData<DisplacementDim>>>
86TH2MLocalAssembler<ShapeFunctionDisplacement, ShapeFunctionPressure,
87 DisplacementDim>::
88 updateConstitutiveVariables(
89 Eigen::VectorXd const& local_x, Eigen::VectorXd const& local_x_prev,
90 double const t, double const dt,
92 models)
93{
94 [[maybe_unused]] auto const matrix_size =
95 gas_pressure_size + capillary_pressure_size + temperature_size +
96 displacement_size;
97
98 assert(local_x.size() == matrix_size);
99
100 auto const gas_pressure =
101 local_x.template segment<gas_pressure_size>(gas_pressure_index);
102 auto const gas_pressure_prev =
103 local_x_prev.template segment<gas_pressure_size>(gas_pressure_index);
104 auto const capillary_pressure =
105 local_x.template segment<capillary_pressure_size>(
106 capillary_pressure_index);
107 auto const capillary_pressure_prev =
108 local_x_prev.template segment<capillary_pressure_size>(
109 capillary_pressure_index);
110
111 auto const temperature =
112 local_x.template segment<temperature_size>(temperature_index);
113 auto const temperature_prev =
114 local_x_prev.template segment<temperature_size>(temperature_index);
115
116 auto const displacement =
117 local_x.template segment<displacement_size>(displacement_index);
118 auto const displacement_prev =
119 local_x_prev.template segment<displacement_size>(displacement_index);
120
121 auto const& medium =
122 *this->process_data_.media_map.getMedium(this->element_.getID());
123 ConstitutiveRelations::MediaData media_data{medium};
124
125 unsigned const n_integration_points =
126 this->integration_method_.getNumberOfPoints();
127
128 std::vector<ConstitutiveRelations::ConstitutiveData<DisplacementDim>>
129 ip_constitutive_data(n_integration_points);
130 std::vector<ConstitutiveRelations::ConstitutiveTempData<DisplacementDim>>
131 ip_constitutive_variables(n_integration_points);
132
133 for (unsigned ip = 0; ip < n_integration_points; ip++)
134 {
135 auto& ip_data = _ip_data[ip];
136 auto& ip_cv = ip_constitutive_variables[ip];
137 auto& ip_cd = ip_constitutive_data[ip];
138 auto& ip_out = this->output_data_[ip];
139 auto& current_state = this->current_states_[ip];
140 auto& prev_state = this->prev_states_[ip];
141
142 auto const& Np = ip_data.N_p;
143 auto const& NT = Np;
144 auto const& Nu = ip_data.N_u;
145 auto const& gradNu = ip_data.dNdx_u;
146 auto const& gradNp = ip_data.dNdx_p;
148 std::nullopt, this->element_.getID(),
150 NumLib::interpolateCoordinates<ShapeFunctionDisplacement,
152 this->element_, Nu))};
153 auto const x_coord =
154 NumLib::interpolateXCoordinate<ShapeFunctionDisplacement,
156 this->element_, Nu);
157
158 double const T = NT.dot(temperature);
159 double const T_prev = NT.dot(temperature_prev);
160 double const pG = Np.dot(gas_pressure);
161 double const pG_prev = Np.dot(gas_pressure_prev);
162 double const pCap = Np.dot(capillary_pressure);
163 double const pCap_prev = Np.dot(capillary_pressure_prev);
164 ConstitutiveRelations::TemperatureData const T_data{T, T_prev};
165 ConstitutiveRelations::GasPressureData const pGR_data{pG, pG_prev};
167 pCap_prev};
169 this->process_data_.reference_temperature(t, pos)[0]};
171 grad_p_GR{gradNp * gas_pressure};
173 DisplacementDim> const grad_p_cap{gradNp * capillary_pressure};
175 grad_T{gradNp * temperature};
176
177 // medium properties
178 models.elastic_tangent_stiffness_model.eval({pos, t, dt}, T_data,
179 ip_cv.C_el_data);
180
181 models.biot_model.eval({pos, t, dt}, media_data, ip_cv.biot_data);
182
183 auto const Bu =
184 LinearBMatrix::computeBMatrix<DisplacementDim,
185 ShapeFunctionDisplacement::NPOINTS,
187 gradNu, Nu, x_coord, this->is_axially_symmetric_);
188
189 ip_out.eps_data.eps.noalias() = Bu * displacement;
190 models.S_L_model.eval({pos, t, dt}, media_data, pCap_data,
191 current_state.S_L_data);
192
193 models.chi_S_L_model.eval({pos, t, dt}, media_data,
194 current_state.S_L_data,
195 current_state.chi_S_L);
196
197 models.chi_S_L_prev_model.eval({pos, t, dt}, media_data,
198 prev_state.S_L_data, prev_state.chi_S_L);
199
200 // solid phase compressibility
201 models.beta_p_SR_model.eval({pos, t, dt}, ip_cv.biot_data,
202 ip_cv.C_el_data, ip_cv.beta_p_SR);
203
204 // If there is swelling stress rate, compute swelling stress.
205 models.swelling_model.eval(
206 {pos, t, dt}, media_data, ip_cv.C_el_data, current_state.S_L_data,
207 prev_state.S_L_data, prev_state.swelling_data,
208 current_state.swelling_data, ip_cv.swelling_data);
209
210 // solid phase linear thermal expansion coefficient
211 models.s_therm_exp_model.eval({pos, t, dt}, media_data, T_data, T0,
212 ip_cv.s_therm_exp_data);
213
214 models.mechanical_strain_model.eval(
215 T_data, ip_cv.s_therm_exp_data, ip_out.eps_data,
216 Bu * displacement_prev, prev_state.mechanical_strain_data,
217 ip_cv.swelling_data, current_state.mechanical_strain_data);
218
219 models.s_mech_model.eval(
220 {pos, t, dt}, T_data, current_state.mechanical_strain_data,
221 prev_state.mechanical_strain_data, prev_state.eff_stress_data,
222 current_state.eff_stress_data, this->material_states_[ip],
223 ip_cd.s_mech_data, ip_cv.equivalent_plastic_strain_data);
224
225 models.total_stress_model.eval(current_state.eff_stress_data,
226 ip_cv.biot_data, current_state.chi_S_L,
227 pGR_data, pCap_data,
228 ip_cv.total_stress_data);
229
230 models.pure_liquid_density_model.eval({pos, t, dt}, media_data,
231 pGR_data, pCap_data, T_data,
232 current_state.rho_W_LR);
233
235 {pos, t, dt}, media_data, pGR_data, pCap_data, T_data,
236 current_state.rho_W_LR, ip_out.fluid_enthalpy_data,
237 ip_out.mass_mole_fractions_data, ip_out.fluid_density_data,
238 ip_out.vapour_pressure_data, current_state.constituent_density_data,
239 ip_cv.phase_transition_data);
240
241 models.viscosity_model.eval({pos, t, dt}, media_data, T_data,
242 ip_out.mass_mole_fractions_data,
243 ip_cv.viscosity_data);
244
245 models.porosity_model.eval(
246 {pos, t, dt}, media_data, current_state.S_L_data,
247 prev_state.S_L_data, pCap_data, pGR_data, current_state.chi_S_L,
248 prev_state.chi_S_L, ip_cv.beta_p_SR, ip_out.eps_data,
249 Bu * displacement_prev, prev_state.porosity_data,
250 current_state.porosity_data);
251
252 if (medium.hasProperty(MPL::PropertyType::transport_porosity))
253 {
254 models.transport_porosity_model.eval(
255 {pos, t, dt}, media_data, current_state.S_L_data,
256 prev_state.S_L_data, pCap_data, pGR_data, current_state.chi_S_L,
257 prev_state.chi_S_L, ip_cv.beta_p_SR,
258 current_state.mechanical_strain_data,
259 prev_state.mechanical_strain_data,
260 prev_state.transport_porosity_data, current_state.porosity_data,
261 current_state.transport_porosity_data);
262 }
263 else
264 {
265 current_state.transport_porosity_data.phi =
266 current_state.porosity_data.phi;
267 }
268
269 models.permeability_model.eval(
270 {pos, t, dt}, media_data, current_state.S_L_data, pGR_data,
271 pCap_data, T_data, current_state.transport_porosity_data,
272 ip_cv.total_stress_data, current_state.mechanical_strain_data,
273 ip_out.eps_data, ip_cv.equivalent_plastic_strain_data,
274 ip_out.permeability_data);
275
276 models.solid_density_model.eval(
277 {pos, t, dt}, media_data, T_data, current_state.eff_stress_data,
278 pCap_data, pGR_data, current_state.chi_S_L,
279 current_state.porosity_data, ip_out.solid_density_data);
280
281 models.solid_heat_capacity_model.eval({pos, t, dt}, media_data, T_data,
282 ip_cv.solid_heat_capacity_data);
283
284 models.thermal_conductivity_model.eval(
285 {pos, t, dt}, media_data, T_data, current_state.porosity_data,
286 current_state.S_L_data, ip_cv.thermal_conductivity_data);
287
288 models.advection_model.eval(current_state.constituent_density_data,
289 ip_out.permeability_data,
290 current_state.rho_W_LR,
291 ip_cv.viscosity_data,
292 ip_cv.advection_data);
293
294 models.gravity_model.eval(
295 ip_out.fluid_density_data,
296 current_state.porosity_data,
297 current_state.S_L_data,
298 ip_out.solid_density_data,
300 this->process_data_.specific_body_force},
301 ip_cv.volumetric_body_force);
302
303 models.diffusion_velocity_model.eval(grad_p_cap,
304 grad_p_GR,
305 ip_out.mass_mole_fractions_data,
306 ip_cv.phase_transition_data,
307 current_state.porosity_data,
308 current_state.S_L_data,
309 grad_T,
310 ip_out.diffusion_velocity_data);
311
312 models.solid_enthalpy_model.eval(ip_cv.solid_heat_capacity_data, T_data,
313 ip_out.solid_enthalpy_data);
314
315 models.internal_energy_model.eval(ip_out.fluid_density_data,
316 ip_cv.phase_transition_data,
317 current_state.porosity_data,
318 current_state.S_L_data,
319 ip_out.solid_density_data,
320 ip_out.solid_enthalpy_data,
321 current_state.internal_energy_data);
322
324 ip_out.fluid_density_data,
325 ip_out.fluid_enthalpy_data,
326 current_state.porosity_data,
327 current_state.S_L_data,
328 ip_out.solid_density_data,
329 ip_out.solid_enthalpy_data,
330 ip_cv.effective_volumetric_enthalpy_data);
331
332 models.fC_1_model.eval(ip_cv.advection_data, ip_out.fluid_density_data,
333 ip_cv.fC_1);
334
335 if (!this->process_data_.apply_mass_lumping)
336 {
337 models.fC_2a_model.eval(ip_cv.biot_data,
338 pCap_data,
339 current_state.constituent_density_data,
340 current_state.porosity_data,
341 current_state.S_L_data,
342 ip_cv.beta_p_SR,
343 ip_cv.fC_2a);
344 }
345 models.fC_3a_model.eval(dt,
346 current_state.constituent_density_data,
347 prev_state.constituent_density_data,
348 current_state.S_L_data,
349 ip_cv.fC_3a);
350
351 models.fC_4_LCpG_model.eval(ip_cv.advection_data,
352 ip_out.fluid_density_data,
353 ip_cv.phase_transition_data,
354 current_state.porosity_data,
355 current_state.S_L_data,
356 ip_cv.fC_4_LCpG);
357
358 models.fC_4_LCpC_model.eval(ip_cv.advection_data,
359 ip_out.fluid_density_data,
360 ip_cv.phase_transition_data,
361 current_state.porosity_data,
362 current_state.S_L_data,
363 ip_cv.fC_4_LCpC);
364
365 models.fC_4_LCT_model.eval(ip_out.fluid_density_data,
366 ip_cv.phase_transition_data,
367 current_state.porosity_data,
368 current_state.S_L_data,
369 ip_cv.fC_4_LCT);
370
371 models.fC_4_MCpG_model.eval(ip_cv.biot_data,
372 current_state.constituent_density_data,
373 current_state.porosity_data,
374 current_state.S_L_data,
375 ip_cv.beta_p_SR,
376 ip_cv.fC_4_MCpG);
377
378 models.fC_4_MCpC_model.eval(ip_cv.biot_data,
379 pCap_data,
380 current_state.constituent_density_data,
381 current_state.porosity_data,
382 prev_state.S_L_data,
383 current_state.S_L_data,
384 ip_cv.beta_p_SR,
385 ip_cv.fC_4_MCpC);
386
387 models.fC_4_MCT_model.eval(ip_cv.biot_data,
388 current_state.constituent_density_data,
389 current_state.porosity_data,
390 current_state.S_L_data,
391 ip_cv.s_therm_exp_data,
392 ip_cv.fC_4_MCT);
393
394 models.fC_4_MCu_model.eval(ip_cv.biot_data,
395 current_state.constituent_density_data,
396 current_state.S_L_data,
397 ip_cv.fC_4_MCu);
398
399 models.fW_1_model.eval(ip_cv.advection_data, ip_out.fluid_density_data,
400 ip_cv.fW_1);
401
402 if (!this->process_data_.apply_mass_lumping)
403 {
404 models.fW_2_model.eval(ip_cv.biot_data,
405 pCap_data,
406 current_state.constituent_density_data,
407 current_state.porosity_data,
408 current_state.rho_W_LR,
409 current_state.S_L_data,
410 ip_cv.beta_p_SR,
411 ip_cv.fW_2);
412 }
413 models.fW_3a_model.eval(dt,
414 current_state.constituent_density_data,
415 prev_state.constituent_density_data,
416 prev_state.rho_W_LR,
417 current_state.rho_W_LR,
418 current_state.S_L_data,
419 ip_cv.fW_3a);
420
421 models.fW_4_LWpG_model.eval(ip_cv.advection_data,
422 ip_out.fluid_density_data,
423 ip_cv.phase_transition_data,
424 current_state.porosity_data,
425 current_state.S_L_data,
426 ip_cv.fW_4_LWpG);
427
428 models.fW_4_LWpC_model.eval(ip_cv.advection_data,
429 ip_out.fluid_density_data,
430 ip_cv.phase_transition_data,
431 current_state.porosity_data,
432 current_state.S_L_data,
433 ip_cv.fW_4_LWpC);
434
435 models.fW_4_LWT_model.eval(ip_out.fluid_density_data,
436 ip_cv.phase_transition_data,
437 current_state.porosity_data,
438 current_state.S_L_data,
439 ip_cv.fW_4_LWT);
440
441 models.fW_4_MWpG_model.eval(ip_cv.biot_data,
442 current_state.constituent_density_data,
443 current_state.porosity_data,
444 current_state.rho_W_LR,
445 current_state.S_L_data,
446 ip_cv.beta_p_SR,
447 ip_cv.fW_4_MWpG);
448
449 models.fW_4_MWpC_model.eval(ip_cv.biot_data,
450 pCap_data,
451 current_state.constituent_density_data,
452 current_state.porosity_data,
453 prev_state.S_L_data,
454 current_state.rho_W_LR,
455 current_state.S_L_data,
456 ip_cv.beta_p_SR,
457 ip_cv.fW_4_MWpC);
458
459 models.fW_4_MWT_model.eval(ip_cv.biot_data,
460 current_state.constituent_density_data,
461 current_state.porosity_data,
462 current_state.rho_W_LR,
463 current_state.S_L_data,
464 ip_cv.s_therm_exp_data,
465 ip_cv.fW_4_MWT);
466
467 models.fW_4_MWu_model.eval(ip_cv.biot_data,
468 current_state.constituent_density_data,
469 current_state.rho_W_LR,
470 current_state.S_L_data,
471 ip_cv.fW_4_MWu);
472
473 models.fT_1_model.eval(dt,
474 current_state.internal_energy_data,
475 prev_state.internal_energy_data,
476 ip_cv.fT_1);
477
478 // ---------------------------------------------------------------------
479 // Derivatives for Jacobian
480 // ---------------------------------------------------------------------
481
482 models.darcy_velocity_model.eval(
483 grad_p_cap,
484 ip_out.fluid_density_data,
485 grad_p_GR,
486 ip_out.permeability_data,
488 this->process_data_.specific_body_force},
489 ip_cv.viscosity_data,
490 ip_out.darcy_velocity_data);
491
492 models.fT_2_model.eval(ip_out.darcy_velocity_data,
493 ip_out.fluid_density_data,
494 ip_out.fluid_enthalpy_data,
495 ip_cv.fT_2);
496
497 models.fT_3_model.eval(
498 current_state.constituent_density_data,
499 ip_out.darcy_velocity_data,
500 ip_out.diffusion_velocity_data,
501 ip_out.fluid_density_data,
502 ip_cv.phase_transition_data,
504 this->process_data_.specific_body_force},
505 ip_cv.fT_3);
506
507 models.fu_2_KupC_model.eval(ip_cv.biot_data, current_state.chi_S_L,
508 ip_cv.fu_2_KupC);
509 }
510
511 return {ip_constitutive_data, ip_constitutive_variables};
512}
513
514template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
515 int DisplacementDim>
516std::vector<ConstitutiveRelations::DerivativesData<DisplacementDim>>
517TH2MLocalAssembler<ShapeFunctionDisplacement, ShapeFunctionPressure,
518 DisplacementDim>::
519 updateConstitutiveVariablesDerivatives(
520 Eigen::VectorXd const& local_x, Eigen::VectorXd const& local_x_prev,
521 double const t, double const dt,
522 std::vector<
524 ip_constitutive_data,
525 std::vector<
527 ip_constitutive_variables,
529 models)
530{
531 [[maybe_unused]] auto const matrix_size =
532 gas_pressure_size + capillary_pressure_size + temperature_size +
533 displacement_size;
534
535 assert(local_x.size() == matrix_size);
536
537 auto const gas_pressure =
538 local_x.template segment<gas_pressure_size>(gas_pressure_index);
539 auto const gas_pressure_prev =
540 local_x_prev.template segment<gas_pressure_size>(gas_pressure_index);
541 auto const temperature =
542 local_x.template segment<temperature_size>(temperature_index);
543 auto const temperature_prev =
544 local_x_prev.template segment<temperature_size>(temperature_index);
545 auto const displacement_prev =
546 local_x_prev.template segment<displacement_size>(displacement_index);
547
548 auto const capillary_pressure =
549 local_x.template segment<capillary_pressure_size>(
550 capillary_pressure_index);
551 auto const capillary_pressure_prev =
552 local_x_prev.template segment<capillary_pressure_size>(
553 capillary_pressure_index);
554
555 auto const& medium =
556 *this->process_data_.media_map.getMedium(this->element_.getID());
557 ConstitutiveRelations::MediaData media_data{medium};
558
559 unsigned const n_integration_points =
560 this->integration_method_.getNumberOfPoints();
561
562 std::vector<ConstitutiveRelations::DerivativesData<DisplacementDim>>
563 ip_d_data(n_integration_points);
564
565 for (unsigned ip = 0; ip < n_integration_points; ip++)
566 {
567 auto const& ip_data = _ip_data[ip];
568 auto& ip_dd = ip_d_data[ip];
569 auto const& ip_cd = ip_constitutive_data[ip];
570 auto const& ip_cv = ip_constitutive_variables[ip];
571 auto const& ip_out = this->output_data_[ip];
572 auto const& current_state = this->current_states_[ip];
573 auto const& prev_state = this->prev_states_[ip];
574
575 auto const& Nu = ip_data.N_u;
576 auto const& Np = ip_data.N_p;
577 auto const& NT = Np;
578 auto const& gradNu = ip_data.dNdx_u;
579 auto const x_coord =
580 NumLib::interpolateXCoordinate<ShapeFunctionDisplacement,
582 this->element_, Nu);
583
585 std::nullopt, this->element_.getID(),
587 NumLib::interpolateCoordinates<ShapeFunctionDisplacement,
589 this->element_, Nu))};
590
591 double const T = NT.dot(temperature);
592 double const T_prev = NT.dot(temperature_prev);
593 double const pG = Np.dot(gas_pressure);
594 double const pG_prev = Np.dot(gas_pressure_prev);
595 double const pCap = Np.dot(capillary_pressure);
596 double const pCap_prev = Np.dot(capillary_pressure_prev);
597 ConstitutiveRelations::TemperatureData const T_data{T, T_prev};
598 ConstitutiveRelations::GasPressureData const pGR_data{pG, pG_prev};
600 pCap_prev};
601
602 auto const Bu =
603 LinearBMatrix::computeBMatrix<DisplacementDim,
604 ShapeFunctionDisplacement::NPOINTS,
606 gradNu, Nu, x_coord, this->is_axially_symmetric_);
607
608 models.S_L_model.dEval({pos, t, dt}, media_data, pCap_data,
609 ip_dd.dS_L_dp_cap);
610
611 models.advection_model.dEval(current_state.constituent_density_data,
612 ip_out.permeability_data,
613 ip_cv.viscosity_data,
614 ip_dd.dS_L_dp_cap,
615 ip_cv.phase_transition_data,
616 ip_dd.advection_d_data);
617
618 models.porosity_model.dEval(
619 {pos, t, dt}, media_data, current_state.S_L_data,
620 prev_state.S_L_data, pCap_data, pGR_data, current_state.chi_S_L,
621 prev_state.chi_S_L, ip_cv.beta_p_SR, ip_out.eps_data,
622 Bu * displacement_prev, prev_state.porosity_data,
623 ip_dd.porosity_d_data);
624
625 models.thermal_conductivity_model.dEval(
626 {pos, t, dt}, media_data, T_data, current_state.porosity_data,
627 ip_dd.porosity_d_data, current_state.S_L_data,
628 ip_dd.thermal_conductivity_d_data);
629
630 models.solid_density_model.dEval(
631 {pos, t, dt}, media_data, T_data, current_state.eff_stress_data,
632 pCap_data, pGR_data, current_state.chi_S_L,
633 current_state.porosity_data, ip_dd.solid_density_d_data);
634
636 ip_out.fluid_density_data,
637 ip_cv.phase_transition_data,
638 current_state.porosity_data,
639 ip_dd.porosity_d_data,
640 current_state.S_L_data,
641 ip_out.solid_density_data,
642 ip_dd.solid_density_d_data,
643 ip_out.solid_enthalpy_data,
644 ip_cv.solid_heat_capacity_data,
645 ip_dd.effective_volumetric_internal_energy_d_data);
646
648 ip_out.fluid_density_data,
649 ip_out.fluid_enthalpy_data,
650 ip_cv.phase_transition_data,
651 current_state.porosity_data,
652 ip_dd.porosity_d_data,
653 current_state.S_L_data,
654 ip_out.solid_density_data,
655 ip_dd.solid_density_d_data,
656 ip_out.solid_enthalpy_data,
657 ip_cv.solid_heat_capacity_data,
658 ip_dd.effective_volumetric_enthalpy_d_data);
659 if (!this->process_data_.apply_mass_lumping)
660 {
661 models.fC_2a_model.dEval(ip_cv.biot_data,
662 pCap_data,
663 current_state.constituent_density_data,
664 ip_cv.phase_transition_data,
665 current_state.porosity_data,
666 ip_dd.porosity_d_data,
667 current_state.S_L_data,
668 ip_dd.dS_L_dp_cap,
669 ip_cv.beta_p_SR,
670 ip_dd.dfC_2a);
671 }
672 models.fC_3a_model.dEval(dt,
673 current_state.constituent_density_data,
674 prev_state.constituent_density_data,
675 ip_cv.phase_transition_data,
676 current_state.S_L_data,
677 ip_dd.dS_L_dp_cap,
678 ip_dd.dfC_3a);
679
680 models.fC_4_LCpG_model.dEval(ip_out.permeability_data,
681 ip_cv.viscosity_data,
682 ip_cv.phase_transition_data,
683 ip_dd.advection_d_data,
684 ip_dd.dfC_4_LCpG);
685
686 models.fC_4_LCpC_model.dEval(current_state.constituent_density_data,
687 ip_out.permeability_data,
688 ip_cv.phase_transition_data,
689 ip_dd.dS_L_dp_cap,
690 ip_cv.viscosity_data,
691 ip_dd.dfC_4_LCpC);
692
693 models.fC_4_MCpG_model.dEval(ip_cv.biot_data,
694 current_state.constituent_density_data,
695 ip_cv.phase_transition_data,
696 current_state.porosity_data,
697 ip_dd.porosity_d_data,
698 current_state.S_L_data,
699 ip_cv.beta_p_SR,
700 ip_dd.dfC_4_MCpG);
701
702 models.fC_4_MCT_model.dEval(ip_cv.biot_data,
703 current_state.constituent_density_data,
704 ip_cv.phase_transition_data,
705 current_state.porosity_data,
706 ip_dd.porosity_d_data,
707 current_state.S_L_data,
708 ip_cv.s_therm_exp_data,
709 ip_dd.dfC_4_MCT);
710
711 models.fC_4_MCu_model.dEval(ip_cv.biot_data,
712 ip_cv.phase_transition_data,
713 current_state.S_L_data,
714 ip_dd.dfC_4_MCu);
715
716 if (!this->process_data_.apply_mass_lumping)
717 {
718 models.fW_2_model.dEval(ip_cv.biot_data,
719 pCap_data,
720 current_state.constituent_density_data,
721 ip_cv.phase_transition_data,
722 current_state.porosity_data,
723 ip_dd.porosity_d_data,
724 current_state.rho_W_LR,
725 current_state.S_L_data,
726 ip_dd.dS_L_dp_cap,
727 ip_cv.beta_p_SR,
728 ip_dd.dfW_2);
729 }
730
731 models.fW_3a_model.dEval(dt,
732 current_state.constituent_density_data,
733 ip_cv.phase_transition_data,
734 prev_state.constituent_density_data,
735 prev_state.rho_W_LR,
736 current_state.rho_W_LR,
737 current_state.S_L_data,
738 ip_dd.dS_L_dp_cap,
739 ip_dd.dfW_3a);
740
741 models.fW_4_LWpG_model.dEval(current_state.constituent_density_data,
742 ip_out.permeability_data,
743 ip_cv.phase_transition_data,
744 current_state.rho_W_LR,
745 ip_dd.dS_L_dp_cap,
746 ip_cv.viscosity_data,
747 ip_dd.dfW_4_LWpG);
748
749 models.fW_4_LWpC_model.dEval(ip_cv.advection_data,
750 ip_out.fluid_density_data,
751 ip_out.permeability_data,
752 ip_cv.phase_transition_data,
753 current_state.porosity_data,
754 current_state.rho_W_LR,
755 current_state.S_L_data,
756 ip_dd.dS_L_dp_cap,
757 ip_cv.viscosity_data,
758 ip_dd.dfW_4_LWpC);
759
760 models.fT_1_model.dEval(
761 dt, ip_dd.effective_volumetric_internal_energy_d_data, ip_dd.dfT_1);
762
763 models.fT_2_model.dEval(
764 ip_out.darcy_velocity_data,
765 ip_out.fluid_density_data,
766 ip_out.fluid_enthalpy_data,
767 ip_out.permeability_data,
768 ip_cv.phase_transition_data,
770 this->process_data_.specific_body_force},
771 ip_cv.viscosity_data,
772 ip_dd.dfT_2);
773
774 models.fu_1_KuT_model.dEval(ip_cd.s_mech_data, ip_cv.s_therm_exp_data,
775 ip_dd.dfu_1_KuT);
776
777 models.fu_2_KupC_model.dEval(ip_cv.biot_data,
778 current_state.chi_S_L,
779 pCap_data,
780 ip_dd.dS_L_dp_cap,
781 ip_dd.dfu_2_KupC);
782 }
783
784 return ip_d_data;
785}
786
787template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
788 int DisplacementDim>
789std::size_t TH2MLocalAssembler<
790 ShapeFunctionDisplacement, ShapeFunctionPressure,
791 DisplacementDim>::setIPDataInitialConditions(std::string_view name,
792 double const* values,
793 int const integration_order)
794{
795 if (integration_order !=
796 static_cast<int>(this->integration_method_.getIntegrationOrder()))
797 {
798 OGS_FATAL(
799 "Setting integration point initial conditions; The integration "
800 "order of the local assembler for element {:d} is different "
801 "from the integration order in the initial condition.",
802 this->element_.getID());
803 }
804
805 if (name == "sigma" && this->process_data_.initial_stress.value)
806 {
807 OGS_FATAL(
808 "Setting initial conditions for stress from integration "
809 "point data and from a parameter '{:s}' is not possible "
810 "simultaneously.",
811 this->process_data_.initial_stress.value->name);
812 }
813
814 if (name.starts_with("material_state_variable_"))
815 {
816 name.remove_prefix(24);
817 DBUG("Setting material state variable '{:s}'", name);
818
819 auto const& internal_variables =
820 this->solid_material_.getInternalVariables();
821 if (auto const iv = std::find_if(
822 begin(internal_variables), end(internal_variables),
823 [&name](auto const& iv) { return iv.name == name; });
824 iv != end(internal_variables))
825 {
826 DBUG("Setting material state variable '{:s}'", name);
828 values, this->material_states_,
830 DisplacementDim>::material_state_variables,
831 iv->reference);
832 }
833
834 WARN(
835 "Could not find variable {:s} in solid material model's "
836 "internal variables.",
837 name);
838 return 0;
839 }
840
841 // TODO this logic could be pulled out of the local assembler into the
842 // process. That might lead to a slightly better performance due to less
843 // string comparisons.
845 name, values, this->current_states_);
846}
847
848template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
849 int DisplacementDim>
850void TH2MLocalAssembler<ShapeFunctionDisplacement, ShapeFunctionPressure,
851 DisplacementDim>::
852 setInitialConditionsConcrete(Eigen::VectorXd const local_x,
853 double const t,
854 int const /*process_id*/)
855{
856 [[maybe_unused]] auto const matrix_size =
857 gas_pressure_size + capillary_pressure_size + temperature_size +
858 displacement_size;
859
860 assert(local_x.size() == matrix_size);
861
862 auto const capillary_pressure =
863 local_x.template segment<capillary_pressure_size>(
864 capillary_pressure_index);
865
866 auto const p_GR =
867 local_x.template segment<gas_pressure_size>(gas_pressure_index);
868
869 auto const temperature =
870 local_x.template segment<temperature_size>(temperature_index);
871
872 auto const displacement =
873 local_x.template segment<displacement_size>(displacement_index);
874
875 constexpr double dt = std::numeric_limits<double>::quiet_NaN();
876 auto const& medium =
877 *this->process_data_.media_map.getMedium(this->element_.getID());
878 auto const& solid_phase = medium.phase("Solid");
879
881 this->solid_material_, *this->process_data_.phase_transition_model_};
882
883 unsigned const n_integration_points =
884 this->integration_method_.getNumberOfPoints();
885
886 for (unsigned ip = 0; ip < n_integration_points; ip++)
887 {
889
890 auto& ip_data = _ip_data[ip];
891 auto& ip_out = this->output_data_[ip];
892 auto& prev_state = this->prev_states_[ip];
893 auto const& Np = ip_data.N_p;
894 auto const& NT = Np;
895 auto const& Nu = ip_data.N_u;
896 auto const& gradNu = ip_data.dNdx_u;
897 auto const x_coord =
898 NumLib::interpolateXCoordinate<ShapeFunctionDisplacement,
900 this->element_, Nu);
902 std::nullopt, this->element_.getID(),
904 NumLib::interpolateCoordinates<ShapeFunctionDisplacement,
906 this->element_, ip_data.N_u))};
907
908 double const pCap = Np.dot(capillary_pressure);
909 vars.capillary_pressure = pCap;
910
911 double const T = NT.dot(temperature);
913 T, T}; // T_prev = T in initialization.
914 vars.temperature = T;
915
916 auto const Bu =
917 LinearBMatrix::computeBMatrix<DisplacementDim,
918 ShapeFunctionDisplacement::NPOINTS,
920 gradNu, Nu, x_coord, this->is_axially_symmetric_);
921
922 auto& eps = ip_out.eps_data.eps;
923 eps.noalias() = Bu * displacement;
924
925 // Set volumetric strain rate for the general case without swelling.
926 vars.volumetric_strain = Invariants::trace(eps);
927
928 double const S_L =
929 medium.property(MPL::PropertyType::saturation)
930 .template value<double>(
931 vars, pos, t, std::numeric_limits<double>::quiet_NaN());
932 this->prev_states_[ip].S_L_data->S_L = S_L;
933
934 // TODO (naumov) Double computation of C_el might be avoided if
935 // updateConstitutiveVariables is called before. But it might interfere
936 // with eps_m initialization.
938 C_el_data;
939 models.elastic_tangent_stiffness_model.eval({pos, t, dt}, T_data,
940 C_el_data);
941 auto const& C_el = C_el_data.stiffness_tensor;
942
943 // Set eps_m_prev from potentially non-zero eps and sigma_sw from
944 // restart.
945 auto const& sigma_sw = this->current_states_[ip].swelling_data.sigma_sw;
946 prev_state.mechanical_strain_data->eps_m.noalias() =
947 solid_phase.hasProperty(MPL::PropertyType::swelling_stress_rate)
948 ? eps + C_el.inverse() * sigma_sw
949 : eps;
950
951 if (this->process_data_.initial_stress.isTotalStress())
952 {
953 auto const alpha_b =
954 medium.property(MPL::PropertyType::biot_coefficient)
955 .template value<double>(vars, pos, t, 0.0 /*dt*/);
956
957 vars.liquid_saturation = S_L;
958 double const bishop =
959 medium.property(MPL::PropertyType::bishops_effective_stress)
960 .template value<double>(vars, pos, t, 0.0 /*dt*/);
961
962 this->current_states_[ip].eff_stress_data.sigma_eff.noalias() +=
963 alpha_b * Np.dot(p_GR - bishop * capillary_pressure) *
964 Invariants::identity2;
965 this->prev_states_[ip].eff_stress_data =
966 this->current_states_[ip].eff_stress_data;
967 }
968 }
969
970 // local_x_prev equal to local_x s.t. the local_x_dot is zero.
971 updateConstitutiveVariables(local_x, local_x, t, 0, models);
972
973 for (unsigned ip = 0; ip < n_integration_points; ip++)
974 {
975 this->material_states_[ip].pushBackState();
976 this->prev_states_[ip] = this->current_states_[ip];
977 }
978}
979
980template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
981 int DisplacementDim>
983 ShapeFunctionDisplacement, ShapeFunctionPressure,
984 DisplacementDim>::assemble(double const t, double const dt,
985 std::vector<double> const& local_x,
986 std::vector<double> const& local_x_prev,
987 std::vector<double>& local_M_data,
988 std::vector<double>& local_K_data,
989 std::vector<double>& local_rhs_data)
990{
991 auto const matrix_size = gas_pressure_size + capillary_pressure_size +
992 temperature_size + displacement_size;
993 assert(local_x.size() == matrix_size);
994
995 auto const capillary_pressure =
996 Eigen::Map<VectorType<capillary_pressure_size> const>(
997 local_x.data() + capillary_pressure_index, capillary_pressure_size);
998
999 auto const capillary_pressure_prev =
1000 Eigen::Map<VectorType<capillary_pressure_size> const>(
1001 local_x_prev.data() + capillary_pressure_index,
1002 capillary_pressure_size);
1003
1004 // pointer to local_M_data vector
1005 auto local_M =
1007 local_M_data, matrix_size, matrix_size);
1008
1009 // pointer to local_K_data vector
1010 auto local_K =
1012 local_K_data, matrix_size, matrix_size);
1013
1014 // pointer to local_rhs_data vector
1016 local_rhs_data, matrix_size);
1017
1018 // component-formulation
1019 // W - liquid phase main component
1020 // C - gas phase main component
1021 // pointer-matrices to the mass matrix - C component equation
1022 auto MCpG = local_M.template block<C_size, gas_pressure_size>(
1023 C_index, gas_pressure_index);
1024 auto MCpC = local_M.template block<C_size, capillary_pressure_size>(
1025 C_index, capillary_pressure_index);
1026 auto MCT = local_M.template block<C_size, temperature_size>(
1027 C_index, temperature_index);
1028 auto MCu = local_M.template block<C_size, displacement_size>(
1029 C_index, displacement_index);
1030
1031 // pointer-matrices to the stiffness matrix - C component equation
1032 auto LCpG = local_K.template block<C_size, gas_pressure_size>(
1033 C_index, gas_pressure_index);
1034 auto LCpC = local_K.template block<C_size, capillary_pressure_size>(
1035 C_index, capillary_pressure_index);
1036 auto LCT = local_K.template block<C_size, temperature_size>(
1037 C_index, temperature_index);
1038
1039 // pointer-matrices to the mass matrix - W component equation
1040 auto MWpG = local_M.template block<W_size, gas_pressure_size>(
1041 W_index, gas_pressure_index);
1042 auto MWpC = local_M.template block<W_size, capillary_pressure_size>(
1043 W_index, capillary_pressure_index);
1044 auto MWT = local_M.template block<W_size, temperature_size>(
1045 W_index, temperature_index);
1046 auto MWu = local_M.template block<W_size, displacement_size>(
1047 W_index, displacement_index);
1048
1049 // pointer-matrices to the stiffness matrix - W component equation
1050 auto LWpG = local_K.template block<W_size, gas_pressure_size>(
1051 W_index, gas_pressure_index);
1052 auto LWpC = local_K.template block<W_size, capillary_pressure_size>(
1053 W_index, capillary_pressure_index);
1054 auto LWT = local_K.template block<W_size, temperature_size>(
1055 W_index, temperature_index);
1056
1057 // pointer-matrices to the mass matrix - temperature equation
1058 auto MTu = local_M.template block<temperature_size, displacement_size>(
1059 temperature_index, displacement_index);
1060
1061 // pointer-matrices to the stiffness matrix - temperature equation
1062 auto KTT = local_K.template block<temperature_size, temperature_size>(
1063 temperature_index, temperature_index);
1064
1065 // pointer-matrices to the stiffness matrix - displacement equation
1066 auto KUpG = local_K.template block<displacement_size, gas_pressure_size>(
1067 displacement_index, gas_pressure_index);
1068 auto KUpC =
1069 local_K.template block<displacement_size, capillary_pressure_size>(
1070 displacement_index, capillary_pressure_index);
1071
1072 auto KUU = local_K.template block<displacement_size, displacement_size>(
1073 displacement_index, displacement_index);
1074
1075 // pointer-vectors to the right hand side terms - C-component equation
1076 auto fC = local_f.template segment<C_size>(C_index);
1077 // pointer-vectors to the right hand side terms - W-component equation
1078 auto fW = local_f.template segment<W_size>(W_index);
1079 // pointer-vectors to the right hand side terms - temperature equation
1080 auto fT = local_f.template segment<temperature_size>(temperature_index);
1081 // pointer-vectors to the right hand side terms - displacement equation
1082 auto fU = local_f.template segment<displacement_size>(displacement_index);
1083
1084 unsigned const n_integration_points =
1085 this->integration_method_.getNumberOfPoints();
1086
1088 this->solid_material_, *this->process_data_.phase_transition_model_};
1089
1090 auto const [ip_constitutive_data, ip_constitutive_variables] =
1091 updateConstitutiveVariables(
1092 Eigen::Map<Eigen::VectorXd const>(local_x.data(), local_x.size()),
1093 Eigen::Map<Eigen::VectorXd const>(local_x_prev.data(),
1094 local_x_prev.size()),
1095 t, dt, models);
1096
1097 for (unsigned int_point = 0; int_point < n_integration_points; int_point++)
1098 {
1099 auto& ip = _ip_data[int_point];
1100 auto& ip_cv = ip_constitutive_variables[int_point];
1101 auto& ip_cd = ip_constitutive_data[int_point];
1102
1103 auto& current_state = this->current_states_[int_point];
1104 auto const& prev_state = this->prev_states_[int_point];
1105
1106 auto const& Np = ip.N_p;
1107 auto const& NT = Np;
1108 auto const& Nu = ip.N_u;
1110 std::nullopt, this->element_.getID(),
1112 NumLib::interpolateCoordinates<ShapeFunctionDisplacement,
1114 this->element_, Nu))};
1115
1116 auto const& NpT = Np.transpose().eval();
1117 auto const& NTT = NT.transpose().eval();
1118
1119 auto const& gradNp = ip.dNdx_p;
1120 auto const& gradNT = gradNp;
1121 auto const& gradNu = ip.dNdx_u;
1122
1123 auto const& gradNpT = gradNp.transpose().eval();
1124 auto const& gradNTT = gradNT.transpose().eval();
1125
1126 auto const& w = ip.integration_weight;
1127
1128 auto const x_coord =
1129 NumLib::interpolateXCoordinate<ShapeFunctionDisplacement,
1131 this->element_, Nu);
1132
1133 auto const Bu =
1134 LinearBMatrix::computeBMatrix<DisplacementDim,
1135 ShapeFunctionDisplacement::NPOINTS,
1137 gradNu, Nu, x_coord, this->is_axially_symmetric_);
1138
1139 auto const NTN = (Np.transpose() * Np).eval();
1140 auto const BTI2N = (Bu.transpose() * Invariants::identity2 * Np).eval();
1141
1142 double const pCap = Np.dot(capillary_pressure);
1143 double const pCap_prev = Np.dot(capillary_pressure_prev);
1144
1145 auto const s_L = current_state.S_L_data.S_L;
1146 auto const s_L_dot = (s_L - prev_state.S_L_data->S_L) / dt;
1147
1148 auto const& b = this->process_data_.specific_body_force;
1149
1150 // ---------------------------------------------------------------------
1151 // C-component equation
1152 // ---------------------------------------------------------------------
1153
1154 MCpG.noalias() += NTN * (ip_cv.fC_4_MCpG.m * w);
1155 MCpC.noalias() += NTN * (ip_cv.fC_4_MCpC.m * w);
1156
1157 if (this->process_data_.apply_mass_lumping)
1158 {
1159 if (pCap - pCap_prev != 0.) // avoid division by Zero
1160 {
1161 MCpC.noalias() +=
1162 NTN * (ip_cv.fC_4_MCpC.ml / (pCap - pCap_prev) * w);
1163 }
1164 }
1165
1166 MCT.noalias() += NTN * (ip_cv.fC_4_MCT.m * w);
1167 MCu.noalias() += BTI2N.transpose() * (ip_cv.fC_4_MCu.m * w);
1168
1169 LCpG.noalias() += gradNpT * ip_cv.fC_4_LCpG.L * gradNp * w;
1170
1171 LCpC.noalias() += gradNpT * ip_cv.fC_4_LCpC.L * gradNp * w;
1172
1173 LCT.noalias() += gradNpT * ip_cv.fC_4_LCT.L * gradNp * w;
1174
1175 fC.noalias() += gradNpT * ip_cv.fC_1.A * b * w;
1176
1177 if (!this->process_data_.apply_mass_lumping)
1178 {
1179 fC.noalias() -= NpT * (ip_cv.fC_2a.a * s_L_dot * w);
1180 }
1181 // fC_III
1182 fC.noalias() -=
1183 NpT * (current_state.porosity_data.phi * ip_cv.fC_3a.a * w);
1184
1185 // ---------------------------------------------------------------------
1186 // W-component equation
1187 // ---------------------------------------------------------------------
1188
1189 MWpG.noalias() += NTN * (ip_cv.fW_4_MWpG.m * w);
1190 MWpC.noalias() += NTN * (ip_cv.fW_4_MWpC.m * w);
1191
1192 if (this->process_data_.apply_mass_lumping)
1193 {
1194 if (pCap - pCap_prev != 0.) // avoid division by Zero
1195 {
1196 MWpC.noalias() +=
1197 NTN * (ip_cv.fW_4_MWpC.ml / (pCap - pCap_prev) * w);
1198 }
1199 }
1200
1201 MWT.noalias() += NTN * (ip_cv.fW_4_MWT.m * w);
1202
1203 MWu.noalias() += BTI2N.transpose() * (ip_cv.fW_4_MWu.m * w);
1204
1205 LWpG.noalias() += gradNpT * ip_cv.fW_4_LWpG.L * gradNp * w;
1206
1207 LWpC.noalias() += gradNpT * ip_cv.fW_4_LWpC.L * gradNp * w;
1208
1209 LWT.noalias() += gradNpT * ip_cv.fW_4_LWT.L * gradNp * w;
1210
1211 fW.noalias() += gradNpT * ip_cv.fW_1.A * b * w;
1212
1213 if (!this->process_data_.apply_mass_lumping)
1214 {
1215 fW.noalias() -= NpT * (ip_cv.fW_2.a * s_L_dot * w);
1216 }
1217
1218 fW.noalias() -=
1219 NpT * (current_state.porosity_data.phi * ip_cv.fW_3a.a * w);
1220
1221 // ---------------------------------------------------------------------
1222 // - temperature equation
1223 // ---------------------------------------------------------------------
1224
1225 MTu.noalias() +=
1226 BTI2N.transpose() *
1227 (ip_cv.effective_volumetric_enthalpy_data.rho_h_eff * w);
1228
1229 KTT.noalias() +=
1230 gradNTT * ip_cv.thermal_conductivity_data.lambda * gradNT * w;
1231
1232 fT.noalias() -= NTT * (ip_cv.fT_1.m * w);
1233
1234 fT.noalias() += gradNTT * ip_cv.fT_2.A * w;
1235
1236 fT.noalias() += gradNTT * ip_cv.fT_3.gradN * w;
1237
1238 fT.noalias() += NTT * (ip_cv.fT_3.N * w);
1239
1240 // ---------------------------------------------------------------------
1241 // - displacement equation
1242 // ---------------------------------------------------------------------
1243
1244 KUpG.noalias() -= BTI2N * (ip_cv.biot_data() * w);
1245
1246 KUpC.noalias() += BTI2N * (ip_cv.fu_2_KupC.m * w);
1247
1248 KUU.noalias() +=
1249 Bu.transpose() * ip_cd.s_mech_data.stiffness_tensor * Bu * w;
1250
1251 fU.noalias() -=
1252 (Bu.transpose() * current_state.eff_stress_data.sigma_eff -
1253 N_u_op(Nu).transpose() * ip_cv.volumetric_body_force()) *
1254 w;
1255
1256 if (this->process_data_.apply_mass_lumping)
1257 {
1258 MCpG = MCpG.colwise().sum().eval().asDiagonal();
1259 MCpC = MCpC.colwise().sum().eval().asDiagonal();
1260 MWpG = MWpG.colwise().sum().eval().asDiagonal();
1261 MWpC = MWpC.colwise().sum().eval().asDiagonal();
1262 }
1263 } // int_point-loop
1264}
1265
1266// Assembles the local Jacobian matrix. So far, the linearisation of HT part is
1267// not considered as that in HT process.
1268template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
1269 int DisplacementDim>
1270void TH2MLocalAssembler<ShapeFunctionDisplacement, ShapeFunctionPressure,
1271 DisplacementDim>::
1272 assembleWithJacobian(double const t, double const dt,
1273 std::vector<double> const& local_x,
1274 std::vector<double> const& local_x_prev,
1275 std::vector<double>& local_rhs_data,
1276 std::vector<double>& local_Jac_data)
1277{
1278 auto const matrix_size = gas_pressure_size + capillary_pressure_size +
1279 temperature_size + displacement_size;
1280 assert(local_x.size() == matrix_size);
1281
1282 auto const temperature = Eigen::Map<VectorType<temperature_size> const>(
1283 local_x.data() + temperature_index, temperature_size);
1284
1285 auto const gas_pressure = Eigen::Map<VectorType<gas_pressure_size> const>(
1286 local_x.data() + gas_pressure_index, gas_pressure_size);
1287
1288 auto const capillary_pressure =
1289 Eigen::Map<VectorType<capillary_pressure_size> const>(
1290 local_x.data() + capillary_pressure_index, capillary_pressure_size);
1291
1292 auto const displacement = Eigen::Map<VectorType<displacement_size> const>(
1293 local_x.data() + displacement_index, displacement_size);
1294
1295 auto const gas_pressure_prev =
1296 Eigen::Map<VectorType<gas_pressure_size> const>(
1297 local_x_prev.data() + gas_pressure_index, gas_pressure_size);
1298
1299 auto const capillary_pressure_prev =
1300 Eigen::Map<VectorType<capillary_pressure_size> const>(
1301 local_x_prev.data() + capillary_pressure_index,
1302 capillary_pressure_size);
1303
1304 auto const temperature_prev =
1305 Eigen::Map<VectorType<temperature_size> const>(
1306 local_x_prev.data() + temperature_index, temperature_size);
1307
1308 auto const displacement_prev =
1309 Eigen::Map<VectorType<displacement_size> const>(
1310 local_x_prev.data() + displacement_index, displacement_size);
1311
1312 auto local_Jac =
1314 local_Jac_data, matrix_size, matrix_size);
1315
1317 local_rhs_data, matrix_size);
1318
1319 // component-formulation
1320 // W - liquid phase main component
1321 // C - gas phase main component
1322
1323 // C component equation matrices
1325 MatrixType<C_size, gas_pressure_size>::Zero(C_size, gas_pressure_size);
1328 C_size, capillary_pressure_size);
1330 MatrixType<C_size, temperature_size>::Zero(C_size, temperature_size);
1332 MatrixType<C_size, displacement_size>::Zero(C_size, displacement_size);
1333
1335 MatrixType<C_size, gas_pressure_size>::Zero(C_size, gas_pressure_size);
1338 C_size, capillary_pressure_size);
1340 MatrixType<C_size, temperature_size>::Zero(C_size, temperature_size);
1341
1342 // mass matrix - W component equation
1344 MatrixType<W_size, gas_pressure_size>::Zero(W_size, gas_pressure_size);
1347 W_size, capillary_pressure_size);
1349 MatrixType<W_size, temperature_size>::Zero(W_size, temperature_size);
1351 MatrixType<W_size, displacement_size>::Zero(W_size, displacement_size);
1352
1353 // stiffness matrix - W component equation
1355 MatrixType<W_size, gas_pressure_size>::Zero(W_size, gas_pressure_size);
1358 W_size, capillary_pressure_size);
1360 MatrixType<W_size, temperature_size>::Zero(W_size, temperature_size);
1361
1362 // mass matrix - temperature equation
1365 temperature_size, displacement_size);
1366
1367 // stiffness matrix - temperature equation
1370 temperature_size);
1371
1372 // stiffness matrices - displacement equation coupling into pressures
1375 displacement_size, gas_pressure_size);
1378 displacement_size, capillary_pressure_size);
1379
1380 // pointer-vectors to the right hand side terms - C-component equation
1381 auto fC = local_f.template segment<C_size>(C_index);
1382 // pointer-vectors to the right hand side terms - W-component equation
1383 auto fW = local_f.template segment<W_size>(W_index);
1384 // pointer-vectors to the right hand side terms - temperature equation
1385 auto fT = local_f.template segment<temperature_size>(temperature_index);
1386 // pointer-vectors to the right hand side terms - displacement equation
1387 auto fU = local_f.template segment<displacement_size>(displacement_index);
1388
1389 unsigned const n_integration_points =
1390 this->integration_method_.getNumberOfPoints();
1391
1393 this->solid_material_, *this->process_data_.phase_transition_model_};
1394
1395 auto const [ip_constitutive_data, ip_constitutive_variables] =
1396 updateConstitutiveVariables(
1397 Eigen::Map<Eigen::VectorXd const>(local_x.data(), local_x.size()),
1398 Eigen::Map<Eigen::VectorXd const>(local_x_prev.data(),
1399 local_x_prev.size()),
1400 t, dt, models);
1401
1402 auto const ip_d_data = updateConstitutiveVariablesDerivatives(
1403 Eigen::Map<Eigen::VectorXd const>(local_x.data(), local_x.size()),
1404 Eigen::Map<Eigen::VectorXd const>(local_x_prev.data(),
1405 local_x_prev.size()),
1406 t, dt, ip_constitutive_data, ip_constitutive_variables, models);
1407
1408 for (unsigned int_point = 0; int_point < n_integration_points; int_point++)
1409 {
1410 auto& ip = _ip_data[int_point];
1411 auto& ip_cd = ip_constitutive_data[int_point];
1412 auto& ip_dd = ip_d_data[int_point];
1413 auto& ip_cv = ip_constitutive_variables[int_point];
1414 auto& current_state = this->current_states_[int_point];
1415 auto& prev_state = this->prev_states_[int_point];
1416
1417 auto const& Np = ip.N_p;
1418 auto const& NT = Np;
1419 auto const& Nu = ip.N_u;
1421 std::nullopt, this->element_.getID(),
1423 NumLib::interpolateCoordinates<ShapeFunctionDisplacement,
1425 this->element_, Nu))};
1426
1427 auto const& NpT = Np.transpose().eval();
1428 auto const& NTT = NT.transpose().eval();
1429
1430 auto const& gradNp = ip.dNdx_p;
1431 auto const& gradNT = gradNp;
1432 auto const& gradNu = ip.dNdx_u;
1433
1434 auto const& gradNpT = gradNp.transpose().eval();
1435 auto const& gradNTT = gradNT.transpose().eval();
1436
1437 auto const& w = ip.integration_weight;
1438
1439 auto const x_coord =
1440 NumLib::interpolateXCoordinate<ShapeFunctionDisplacement,
1442 this->element_, Nu);
1443
1444 auto const Bu =
1445 LinearBMatrix::computeBMatrix<DisplacementDim,
1446 ShapeFunctionDisplacement::NPOINTS,
1448 gradNu, Nu, x_coord, this->is_axially_symmetric_);
1449
1450 auto const NTN = (Np.transpose() * Np).eval();
1451 auto const BTI2N = (Bu.transpose() * Invariants::identity2 * Np).eval();
1452
1453 double const div_u_dot =
1454 Invariants::trace(Bu * (displacement - displacement_prev) / dt);
1455
1456 double const pGR = Np.dot(gas_pressure);
1457 double const pCap = Np.dot(capillary_pressure);
1458 double const T = NT.dot(temperature);
1459
1460 GlobalDimVectorType const gradpGR = gradNp * gas_pressure;
1461 GlobalDimVectorType const gradpCap = gradNp * capillary_pressure;
1462 GlobalDimVectorType const gradT = gradNT * temperature;
1463
1464 double const pGR_prev = Np.dot(gas_pressure_prev);
1465 double const pCap_prev = Np.dot(capillary_pressure_prev);
1466 double const T_prev = NT.dot(temperature_prev);
1467
1468 auto const& s_L = current_state.S_L_data.S_L;
1469 auto const s_L_dot = (s_L - prev_state.S_L_data->S_L) / dt;
1470
1471 auto const& b = this->process_data_.specific_body_force;
1472
1473 // ---------------------------------------------------------------------
1474 // C-component equation
1475 // ---------------------------------------------------------------------
1476
1477 MCpG.noalias() += NTN * (ip_cv.fC_4_MCpG.m * w);
1478 MCpC.noalias() += NTN * (ip_cv.fC_4_MCpC.m * w);
1479
1480 if (this->process_data_.apply_mass_lumping)
1481 {
1482 if (pCap - pCap_prev != 0.) // avoid division by Zero
1483 {
1484 MCpC.noalias() +=
1485 NTN * (ip_cv.fC_4_MCpC.ml / (pCap - pCap_prev) * w);
1486 }
1487 }
1488
1489 MCT.noalias() += NTN * (ip_cv.fC_4_MCT.m * w);
1490 // d (fC_4_MCT * T_dot)/d T
1491 local_Jac
1492 .template block<C_size, temperature_size>(C_index,
1493 temperature_index)
1494 .noalias() += NTN * (ip_dd.dfC_4_MCT.dT * (T - T_prev) / dt * w);
1495
1496 MCu.noalias() += BTI2N.transpose() * (ip_cv.fC_4_MCu.m * w);
1497 // d (fC_4_MCu * u_dot)/d T
1498 local_Jac
1499 .template block<C_size, temperature_size>(C_index,
1500 temperature_index)
1501 .noalias() += NTN * (ip_dd.dfC_4_MCu.dT * div_u_dot * w);
1502
1503 LCpG.noalias() += gradNpT * ip_cv.fC_4_LCpG.L * gradNp * w;
1504
1505 // d (fC_4_LCpG * grad p_GR)/d p_GR
1506 local_Jac.template block<C_size, C_size>(C_index, C_index).noalias() +=
1507 gradNpT * ip_dd.dfC_4_LCpG.dp_GR * gradpGR * Np * w;
1508
1509 // d (fC_4_LCpG * grad p_GR)/d p_cap
1510 local_Jac.template block<C_size, W_size>(C_index, W_index).noalias() +=
1511 gradNpT * ip_dd.dfC_4_LCpG.dp_cap * gradpGR * Np * w;
1512
1513 // d (fC_4_LCpG * grad p_GR)/d T
1514 local_Jac
1515 .template block<C_size, temperature_size>(C_index,
1516 temperature_index)
1517 .noalias() += gradNpT * ip_dd.dfC_4_LCpG.dT * gradpGR * NT * w;
1518
1519 // d (fC_4_MCpG * p_GR_dot)/d p_GR
1520 local_Jac.template block<C_size, C_size>(C_index, C_index).noalias() +=
1521 NTN * (ip_dd.dfC_4_MCpG.dp_GR * (pGR - pGR_prev) / dt * w);
1522
1523 // d (fC_4_MCpG * p_GR_dot)/d T
1524 local_Jac
1525 .template block<C_size, temperature_size>(C_index,
1526 temperature_index)
1527 .noalias() +=
1528 NTN * (ip_dd.dfC_4_MCpG.dT * (pGR - pGR_prev) / dt * w);
1529
1530 LCpC.noalias() -= gradNpT * ip_cv.fC_4_LCpC.L * gradNp * w;
1531
1532 /* TODO (naumov) This part is not tested by any of the current ctests.
1533 // d (fC_4_LCpC * grad p_cap)/d p_GR
1534 local_Jac.template block<C_size, C_size>(C_index, C_index).noalias() +=
1535 gradNpT * ip_dd.dfC_4_LCpC.dp_GR * gradpCap * Np * w;
1536 // d (fC_4_LCpC * grad p_cap)/d p_cap
1537 local_Jac.template block<C_size, W_size>(C_index, W_index).noalias() +=
1538 gradNpT * ip_dd.dfC_4_LCpC.dp_cap * gradpCap * Np * w;
1539
1540 local_Jac
1541 .template block<C_size, temperature_size>(C_index,
1542 temperature_index)
1543 .noalias() += gradNpT * ip_dd.dfC_4_LCpC.dT * gradpCap * Np * w;
1544 */
1545
1546 LCT.noalias() += gradNpT * ip_cv.fC_4_LCT.L * gradNp * w;
1547
1548 // fC_1
1549 fC.noalias() += gradNpT * ip_cv.fC_1.A * b * w;
1550
1551 if (!this->process_data_.apply_mass_lumping)
1552 {
1553 // fC_2 = \int a * s_L_dot
1554 fC.noalias() -= NpT * (ip_cv.fC_2a.a * s_L_dot * w);
1555
1556 local_Jac.template block<C_size, C_size>(C_index, C_index)
1557 .noalias() +=
1558 NTN * ((ip_dd.dfC_2a.dp_GR * s_L_dot
1559 /*- ip_cv.fC_2a.a * (ds_L_dp_GR = 0) / dt*/) *
1560 w);
1561
1562 local_Jac.template block<C_size, W_size>(C_index, W_index)
1563 .noalias() +=
1564 NTN * ((ip_dd.dfC_2a.dp_cap * s_L_dot +
1565 ip_cv.fC_2a.a * ip_dd.dS_L_dp_cap() / dt) *
1566 w);
1567
1568 local_Jac
1569 .template block<C_size, temperature_size>(C_index,
1570 temperature_index)
1571 .noalias() += NTN * (ip_dd.dfC_2a.dT * s_L_dot * w);
1572 }
1573 {
1574 // fC_3 = \int phi * a
1575 fC.noalias() -=
1576 NpT * (current_state.porosity_data.phi * ip_cv.fC_3a.a * w);
1577
1578 local_Jac.template block<C_size, C_size>(C_index, C_index)
1579 .noalias() += NTN * (current_state.porosity_data.phi *
1580 ip_dd.dfC_3a.dp_GR * w);
1581
1582 local_Jac.template block<C_size, W_size>(C_index, W_index)
1583 .noalias() += NTN * (current_state.porosity_data.phi *
1584 ip_dd.dfC_3a.dp_cap * w);
1585
1586 local_Jac
1587 .template block<C_size, temperature_size>(C_index,
1588 temperature_index)
1589 .noalias() +=
1590 NTN * ((ip_dd.porosity_d_data.dphi_dT * ip_cv.fC_3a.a +
1591 current_state.porosity_data.phi * ip_dd.dfC_3a.dT) *
1592 w);
1593 }
1594 // ---------------------------------------------------------------------
1595 // W-component equation
1596 // ---------------------------------------------------------------------
1597
1598 MWpG.noalias() += NTN * (ip_cv.fW_4_MWpG.m * w);
1599 MWpC.noalias() += NTN * (ip_cv.fW_4_MWpC.m * w);
1600
1601 if (this->process_data_.apply_mass_lumping)
1602 {
1603 if (pCap - pCap_prev != 0.) // avoid division by Zero
1604 {
1605 MWpC.noalias() +=
1606 NTN * (ip_cv.fW_4_MWpC.ml / (pCap - pCap_prev) * w);
1607 }
1608 }
1609
1610 MWT.noalias() += NTN * (ip_cv.fW_4_MWT.m * w);
1611
1612 MWu.noalias() += BTI2N.transpose() * (ip_cv.fW_4_MWu.m * w);
1613
1614 LWpG.noalias() += gradNpT * ip_cv.fW_4_LWpG.L * gradNp * w;
1615
1616 // fW_4 LWpG' parts; LWpG = \int grad (a + d) grad
1617 local_Jac.template block<W_size, C_size>(W_index, C_index).noalias() +=
1618 gradNpT * ip_dd.dfW_4_LWpG.dp_GR * gradpGR * Np * w;
1619
1620 local_Jac.template block<W_size, W_size>(W_index, W_index).noalias() +=
1621 gradNpT * ip_dd.dfW_4_LWpG.dp_cap * gradpGR * Np * w;
1622
1623 local_Jac
1624 .template block<W_size, temperature_size>(W_index,
1625 temperature_index)
1626 .noalias() += gradNpT * ip_dd.dfW_4_LWpG.dT * gradpGR * NT * w;
1627
1628 LWpC.noalias() += gradNpT * ip_cv.fW_4_LWpC.L * gradNp * w;
1629
1630 // fW_4 LWp_cap' parts; LWpC = \int grad (a + d) grad
1631 local_Jac.template block<W_size, C_size>(W_index, C_index).noalias() -=
1632 gradNpT * ip_dd.dfW_4_LWpC.dp_GR * gradpCap * Np * w;
1633
1634 local_Jac.template block<W_size, W_size>(W_index, W_index).noalias() -=
1635 gradNpT * ip_dd.dfW_4_LWpC.dp_cap * gradpCap * Np * w;
1636
1637 local_Jac
1638 .template block<W_size, temperature_size>(W_index,
1639 temperature_index)
1640 .noalias() -= gradNpT * ip_dd.dfW_4_LWpC.dT * gradpCap * NT * w;
1641
1642 LWT.noalias() += gradNpT * ip_cv.fW_4_LWT.L * gradNp * w;
1643
1644 // fW_1
1645 fW.noalias() += gradNpT * ip_cv.fW_1.A * b * w;
1646
1647 // fW_2 = \int a * s_L_dot
1648 if (!this->process_data_.apply_mass_lumping)
1649 {
1650 fW.noalias() -= NpT * (ip_cv.fW_2.a * s_L_dot * w);
1651
1652 local_Jac.template block<W_size, C_size>(W_index, C_index)
1653 .noalias() += NTN * (ip_dd.dfW_2.dp_GR * s_L_dot * w);
1654
1655 // sign negated because of dp_cap = -dp_LR
1656 // TODO (naumov) Had to change the sign to get equal Jacobian WW
1657 // blocks in A2 Test. Where is the error?
1658 local_Jac.template block<W_size, W_size>(W_index, W_index)
1659 .noalias() += NTN * ((ip_dd.dfW_2.dp_cap * s_L_dot +
1660 ip_cv.fW_2.a * ip_dd.dS_L_dp_cap() / dt) *
1661 w);
1662
1663 local_Jac
1664 .template block<W_size, temperature_size>(W_index,
1665 temperature_index)
1666 .noalias() += NTN * (ip_dd.dfW_2.dT * s_L_dot * w);
1667 }
1668
1669 // fW_3 = \int phi * a
1670 fW.noalias() -=
1671 NpT * (current_state.porosity_data.phi * ip_cv.fW_3a.a * w);
1672
1673 local_Jac.template block<W_size, C_size>(W_index, C_index).noalias() +=
1674 NTN * (current_state.porosity_data.phi * ip_dd.dfW_3a.dp_GR * w);
1675
1676 local_Jac.template block<W_size, W_size>(W_index, W_index).noalias() +=
1677 NTN * (current_state.porosity_data.phi * ip_dd.dfW_3a.dp_cap * w);
1678
1679 local_Jac
1680 .template block<W_size, temperature_size>(W_index,
1681 temperature_index)
1682 .noalias() +=
1683 NTN * ((ip_dd.porosity_d_data.dphi_dT * ip_cv.fW_3a.a +
1684 current_state.porosity_data.phi * ip_dd.dfW_3a.dT) *
1685 w);
1686
1687 // ---------------------------------------------------------------------
1688 // - temperature equation
1689 // ---------------------------------------------------------------------
1690
1691 MTu.noalias() +=
1692 BTI2N.transpose() *
1693 (ip_cv.effective_volumetric_enthalpy_data.rho_h_eff * w);
1694
1695 // dfT_4/dp_GR
1696 // d (MTu * u_dot)/dp_GR
1697 local_Jac
1698 .template block<temperature_size, C_size>(temperature_index,
1699 C_index)
1700 .noalias() +=
1701 NTN * (ip_dd.effective_volumetric_enthalpy_d_data.drho_h_eff_dp_GR *
1702 div_u_dot * w);
1703
1704 // dfT_4/dp_cap
1705 // d (MTu * u_dot)/dp_cap
1706 local_Jac
1707 .template block<temperature_size, W_size>(temperature_index,
1708 W_index)
1709 .noalias() -=
1710 NTN *
1711 (ip_dd.effective_volumetric_enthalpy_d_data.drho_h_eff_dp_cap *
1712 div_u_dot * w);
1713
1714 // dfT_4/dT
1715 // d (MTu * u_dot)/dT
1716 local_Jac
1717 .template block<temperature_size, temperature_size>(
1718 temperature_index, temperature_index)
1719 .noalias() +=
1720 NTN * (ip_dd.effective_volumetric_enthalpy_d_data.drho_h_eff_dT *
1721 div_u_dot * w);
1722
1723 KTT.noalias() +=
1724 gradNTT * ip_cv.thermal_conductivity_data.lambda * gradNT * w;
1725
1726 // d KTT/dp_GR * T
1727 // TODO (naumov) always zero if lambda_xR have no derivatives wrt. p_GR.
1728 // dlambda_dp_GR =
1729 // (dphi_G_dp_GR = 0) * lambdaGR + phi_G * dlambda_GR_dp_GR +
1730 // (dphi_L_dp_GR = 0) * lambdaLR + phi_L * dlambda_LR_dp_GR +
1731 // (dphi_S_dp_GR = 0) * lambdaSR + phi_S * dlambda_SR_dp_GR +
1732 // = 0
1733 //
1734 // Since dlambda/dp_GR is 0 the derivative is omitted:
1735 // local_Jac
1736 // .template block<temperature_size, C_size>(temperature_index,
1737 // C_index)
1738 // .noalias() += gradNTT * dlambda_dp_GR * gradT * Np * w;
1739
1740 // d KTT/dp_cap * T
1741 local_Jac
1742 .template block<temperature_size, W_size>(temperature_index,
1743 W_index)
1744 .noalias() += gradNTT *
1745 ip_dd.thermal_conductivity_d_data.dlambda_dp_cap *
1746 gradT * Np * w;
1747
1748 // d KTT/dT * T
1749 local_Jac
1750 .template block<temperature_size, temperature_size>(
1751 temperature_index, temperature_index)
1752 .noalias() += gradNTT *
1753 ip_dd.thermal_conductivity_d_data.dlambda_dT * gradT *
1754 NT * w;
1755
1756 // fT_1
1757 fT.noalias() -= NTT * (ip_cv.fT_1.m * w);
1758
1759 // dfT_1/dp_GR
1760 local_Jac
1761 .template block<temperature_size, C_size>(temperature_index,
1762 C_index)
1763 .noalias() += NTN * (ip_dd.dfT_1.dp_GR * w);
1764
1765 // dfT_1/dp_cap
1766 local_Jac
1767 .template block<temperature_size, W_size>(temperature_index,
1768 W_index)
1769 .noalias() += NTN * (ip_dd.dfT_1.dp_cap * w);
1770
1771 // dfT_1/dT
1772 // MTT
1773 local_Jac
1774 .template block<temperature_size, temperature_size>(
1775 temperature_index, temperature_index)
1776 .noalias() += NTN * (ip_dd.dfT_1.dT * w);
1777
1778 // fT_2
1779 fT.noalias() += gradNTT * ip_cv.fT_2.A * w;
1780
1781 // dfT_2/dp_GR
1782 local_Jac
1783 .template block<temperature_size, C_size>(temperature_index,
1784 C_index)
1785 .noalias() -=
1786 // dfT_2/dp_GR first part
1787 gradNTT * ip_dd.dfT_2.dp_GR_Npart * Np * w +
1788 // dfT_2/dp_GR second part
1789 gradNTT * ip_dd.dfT_2.dp_GR_gradNpart * gradNp * w;
1790
1791 // dfT_2/dp_cap
1792 local_Jac
1793 .template block<temperature_size, W_size>(temperature_index,
1794 W_index)
1795 .noalias() -=
1796 // first part of dfT_2/dp_cap
1797 gradNTT * (-ip_dd.dfT_2.dp_cap_Npart) * Np * w +
1798 // second part of dfT_2/dp_cap
1799 gradNTT * (-ip_dd.dfT_2.dp_cap_gradNpart) * gradNp * w;
1800
1801 // dfT_2/dT
1802 local_Jac
1803 .template block<temperature_size, temperature_size>(
1804 temperature_index, temperature_index)
1805 .noalias() -= gradNTT * ip_dd.dfT_2.dT * NT * w;
1806
1807 // fT_3
1808 fT.noalias() += NTT * (ip_cv.fT_3.N * w);
1809
1810 fT.noalias() += gradNTT * ip_cv.fT_3.gradN * w;
1811
1812 // ---------------------------------------------------------------------
1813 // - displacement equation
1814 // ---------------------------------------------------------------------
1815
1816 KUpG.noalias() -= BTI2N * (ip_cv.biot_data() * w);
1817
1818 // dfU_2/dp_GR = dKUpG/dp_GR * p_GR + KUpG. The former is zero, the
1819 // latter is handled below.
1820
1821 KUpC.noalias() += BTI2N * (ip_cv.fu_2_KupC.m * w);
1822
1823 // dfU_2/dp_cap = dKUpC/dp_cap * p_cap + KUpC. The former is handled
1824 // here, the latter below.
1825 local_Jac
1826 .template block<displacement_size, W_size>(displacement_index,
1827 W_index)
1828 .noalias() += BTI2N * (ip_dd.dfu_2_KupC.dp_cap * w);
1829
1830 local_Jac
1831 .template block<displacement_size, displacement_size>(
1832 displacement_index, displacement_index)
1833 .noalias() +=
1834 Bu.transpose() * ip_cd.s_mech_data.stiffness_tensor * Bu * w;
1835
1836 // fU_1
1837 fU.noalias() -=
1838 (Bu.transpose() * current_state.eff_stress_data.sigma_eff -
1839 N_u_op(Nu).transpose() * ip_cv.volumetric_body_force()) *
1840 w;
1841
1842 // KuT
1843 local_Jac
1844 .template block<displacement_size, temperature_size>(
1845 displacement_index, temperature_index)
1846 .noalias() -= Bu.transpose() * ip_dd.dfu_1_KuT.dT * NT * w;
1847
1848 /* TODO (naumov) Test with gravity needed to check this Jacobian part.
1849 local_Jac
1850 .template block<displacement_size, temperature_size>(
1851 displacement_index, temperature_index)
1852 .noalias() += N_u_op(Nu).transpose() * ip_cv.drho_dT * b *
1853 N_u_op(Nu).transpose() * w;
1854 */
1855
1856 if (this->process_data_.apply_mass_lumping)
1857 {
1858 MCpG = MCpG.colwise().sum().eval().asDiagonal();
1859 MCpC = MCpC.colwise().sum().eval().asDiagonal();
1860 MWpG = MWpG.colwise().sum().eval().asDiagonal();
1861 MWpC = MWpC.colwise().sum().eval().asDiagonal();
1862 }
1863 } // int_point-loop
1864
1865 // --- Gas ---
1866 // fC_4
1867 fC.noalias() -= LCpG * gas_pressure + LCpC * capillary_pressure +
1868 LCT * temperature +
1869 MCpG * (gas_pressure - gas_pressure_prev) / dt +
1870 MCpC * (capillary_pressure - capillary_pressure_prev) / dt +
1871 MCT * (temperature - temperature_prev) / dt +
1872 MCu * (displacement - displacement_prev) / dt;
1873
1874 local_Jac.template block<C_size, C_size>(C_index, C_index).noalias() +=
1875 LCpG + MCpG / dt;
1876 local_Jac.template block<C_size, W_size>(C_index, W_index).noalias() +=
1877 LCpC + MCpC / dt;
1878 local_Jac
1879 .template block<C_size, temperature_size>(C_index, temperature_index)
1880 .noalias() += LCT + MCT / dt;
1881 local_Jac
1882 .template block<C_size, displacement_size>(C_index, displacement_index)
1883 .noalias() += MCu / dt;
1884
1885 // --- Capillary pressure ---
1886 // fW_4
1887 fW.noalias() -= LWpG * gas_pressure + LWpC * capillary_pressure +
1888 LWT * temperature +
1889 MWpG * (gas_pressure - gas_pressure_prev) / dt +
1890 MWpC * (capillary_pressure - capillary_pressure_prev) / dt +
1891 MWT * (temperature - temperature_prev) / dt +
1892 MWu * (displacement - displacement_prev) / dt;
1893
1894 local_Jac.template block<W_size, W_size>(W_index, W_index).noalias() +=
1895 LWpC + MWpC / dt;
1896 local_Jac.template block<W_size, C_size>(W_index, C_index).noalias() +=
1897 LWpG + MWpG / dt;
1898 local_Jac
1899 .template block<W_size, temperature_size>(W_index, temperature_index)
1900 .noalias() += LWT + MWT / dt;
1901 local_Jac
1902 .template block<W_size, displacement_size>(W_index, displacement_index)
1903 .noalias() += MWu / dt;
1904
1905 // --- Temperature ---
1906 // fT_4
1907 fT.noalias() -=
1908 KTT * temperature + MTu * (displacement - displacement_prev) / dt;
1909
1910 local_Jac
1911 .template block<temperature_size, temperature_size>(temperature_index,
1912 temperature_index)
1913 .noalias() += KTT;
1914 local_Jac
1915 .template block<temperature_size, displacement_size>(temperature_index,
1916 displacement_index)
1917 .noalias() += MTu / dt;
1918
1919 // --- Displacement ---
1920 // fU_2
1921 fU.noalias() -= KUpG * gas_pressure + KUpC * capillary_pressure;
1922
1923 local_Jac
1924 .template block<displacement_size, C_size>(displacement_index, C_index)
1925 .noalias() += KUpG;
1926 local_Jac
1927 .template block<displacement_size, W_size>(displacement_index, W_index)
1928 .noalias() += KUpC;
1929}
1930
1931template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
1932 int DisplacementDim>
1933void TH2MLocalAssembler<ShapeFunctionDisplacement, ShapeFunctionPressure,
1934 DisplacementDim>::
1935 computeSecondaryVariableConcrete(double const t, double const dt,
1936 Eigen::VectorXd const& local_x,
1937 Eigen::VectorXd const& local_x_prev)
1938{
1939 auto const gas_pressure =
1940 local_x.template segment<gas_pressure_size>(gas_pressure_index);
1941 auto const capillary_pressure =
1942 local_x.template segment<capillary_pressure_size>(
1943 capillary_pressure_index);
1944 auto const liquid_pressure = (gas_pressure - capillary_pressure).eval();
1945
1947 ShapeFunctionPressure, typename ShapeFunctionDisplacement::MeshElement,
1948 DisplacementDim>(this->element_, this->is_axially_symmetric_,
1949 gas_pressure,
1950 *this->process_data_.gas_pressure_interpolated);
1951
1953 ShapeFunctionPressure, typename ShapeFunctionDisplacement::MeshElement,
1954 DisplacementDim>(this->element_, this->is_axially_symmetric_,
1955 capillary_pressure,
1956 *this->process_data_.capillary_pressure_interpolated);
1957
1959 ShapeFunctionPressure, typename ShapeFunctionDisplacement::MeshElement,
1960 DisplacementDim>(this->element_, this->is_axially_symmetric_,
1961 liquid_pressure,
1962 *this->process_data_.liquid_pressure_interpolated);
1963
1964 auto const temperature =
1965 local_x.template segment<temperature_size>(temperature_index);
1966
1968 ShapeFunctionPressure, typename ShapeFunctionDisplacement::MeshElement,
1969 DisplacementDim>(this->element_, this->is_axially_symmetric_,
1970 temperature,
1971 *this->process_data_.temperature_interpolated);
1972
1974 this->solid_material_, *this->process_data_.phase_transition_model_};
1975
1976 updateConstitutiveVariables(local_x, local_x_prev, t, dt, models);
1977}
1978
1979} // namespace TH2M
1980} // namespace ProcessLib
#define OGS_FATAL(...)
Definition Error.h:26
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
void WARN(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:40
constexpr double getWeight() const
MathLib::WeightedPoint const & getWeightedPoint(unsigned const igp) const
MatrixType< _kelvin_vector_size, _number_of_dof > BMatrixType
ShapeMatrixPolicyType< ShapeFunctionDisplacement, DisplacementDim > ShapeMatricesTypeDisplacement
Definition TH2MFEM.h:53
SecondaryData< typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType > _secondary_data
Definition TH2MFEM.h:258
typename ShapeMatricesTypePressure::template MatrixType< M, N > MatrixType
Definition TH2MFEM.h:64
ShapeMatrixPolicyType< ShapeFunctionPressure, DisplacementDim > ShapeMatricesTypePressure
Definition TH2MFEM.h:56
typename ShapeMatricesTypePressure::GlobalDimVectorType GlobalDimVectorType
Definition TH2MFEM.h:70
std::vector< IpData > _ip_data
Definition TH2MFEM.h:254
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)
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)
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::size_t reflectSetIPData(std::string_view const name, double const *values, std::vector< IPData > &ip_data_vector)
std::size_t setIntegrationPointDataMaterialStateVariables(double const *values, IntegrationPointDataVector &ip_data_vector, MemberType member, std::function< std::span< double >(MaterialStateVariables &)> get_values_span)
void setIPDataInitialConditions(std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > const &_integration_point_writer, MeshLib::Properties const &mesh_properties, LocalAssemblersVector &local_assemblers)
void eval(SpaceTimeData const &x_t, MediaData const &media_data, BiotData &out) const
Definition Biot.cpp:16
void eval(SpaceTimeData const &x_t, MediaData const &media_data, SaturationData const &S_L_data, BishopsData &out) const
Definition Bishops.cpp:36
void eval(SpaceTimeData const &x_t, MediaData const &media_data, PrevState< SaturationData > const &S_L_data, PrevState< BishopsData > &out) const
Definition Bishops.cpp:41
Data that is needed for the equation system assembly.
DiffusionVelocityModel< DisplacementDim > diffusion_velocity_model
ElasticTangentStiffnessModel< DisplacementDim > elastic_tangent_stiffness_model
ThermalConductivityModel< DisplacementDim > thermal_conductivity_model
TransportPorosityModel< DisplacementDim > transport_porosity_model
SolidCompressibilityModel< DisplacementDim, SolidConstitutiveRelation< DisplacementDim > > beta_p_SR_model
SolidThermalExpansionModel< DisplacementDim > s_therm_exp_model
MechanicalStrainModel< DisplacementDim > mechanical_strain_model
void dEval(FluidDensityData const &fluid_density_data, FluidEnthalpyData const &fluid_enthalpy_data, PhaseTransitionData const &phase_transition_data, PorosityData const &porosity_data, PorosityDerivativeData const &porosity_d_data, SaturationData const &S_L_data, SolidDensityData const &solid_density_data, SolidDensityDerivativeData const &solid_density_d_data, SolidEnthalpyData const &solid_enthalpy_data, SolidHeatCapacityData const &solid_heat_capacity_data, EffectiveVolumetricEnthalpyDerivatives &effective_volumetric_enthalpy_d_data) const
Definition Enthalpy.cpp:35
void eval(FluidDensityData const &fluid_density_data, FluidEnthalpyData const &fluid_enthalpy_data, PorosityData const &porosity_data, SaturationData const &S_L_data, SolidDensityData const &solid_density_data, SolidEnthalpyData const &solid_enthalpy_data, EffectiveVolumetricEnthalpy &effective_volumetric_enthalpy_data) const
Definition Enthalpy.cpp:16
void dEval(BiotData const &biot_data, CapillaryPressureData const pCap, ConstituentDensityData const &constituent_density_data, PhaseTransitionData const &phase_transition_data, PorosityData const &porosity_data, PorosityDerivativeData const &porosity_d_data, SaturationData const &S_L_data, SaturationDataDeriv const &dS_L_dp_cap, SolidCompressibilityData const &beta_p_SR, FC2aDerivativeData &dfC_2a) const
Definition CEquation.cpp:47
void eval(BiotData const biot_data, CapillaryPressureData const pCap, ConstituentDensityData const &constituent_density_data, PorosityData const &porosity_data, SaturationData const &S_L_data, SolidCompressibilityData const beta_p_SR, FC2aData &fC_2a) const
Definition CEquation.cpp:29
void eval(double const dt, ConstituentDensityData const &constituent_density_data, PrevState< ConstituentDensityData > const &constituent_density_data_prev, SaturationData const &S_L_data, FC3aData &fC_3a) const
void dEval(double const dt, ConstituentDensityData const &constituent_density_data, PrevState< ConstituentDensityData > const &constituent_density_data_prev, PhaseTransitionData const &phase_transition_data, SaturationData const &S_L_data, SaturationDataDeriv const &dS_L_dp_cap, FC3aDerivativeData &dfC_3a) const
void eval(BiotData const &biot_data, CapillaryPressureData const pCap, ConstituentDensityData const &constituent_density_data, PorosityData const &porosity_data, PrevState< SaturationData > const &S_L_data_prev, SaturationData const &S_L_data, SolidCompressibilityData const &beta_p_SR, FC4MCpCData &fC_4_MCpC) const
void eval(BiotData const &biot_data, ConstituentDensityData const &constituent_density_data, PorosityData const &porosity_data, SaturationData const &S_L_data, SolidCompressibilityData const &beta_p_SR, FC4MCpGData &fC_4_MCpG) const
void dEval(BiotData const &biot_data, ConstituentDensityData const &constituent_density_data, PhaseTransitionData const &phase_transition_data, PorosityData const &porosity_data, PorosityDerivativeData const &porosity_d_data, SaturationData const &S_L_data, SolidCompressibilityData const &beta_p_SR, FC4MCpGDerivativeData &dfC_4_MCpG) const
void eval(BiotData const &biot_data, ConstituentDensityData const &constituent_density_data, SaturationData const &S_L_data, FC4MCuData &fC_4_MCu) const
void dEval(BiotData const &biot_data, PhaseTransitionData const &phase_transition_data, SaturationData const &S_L_data, FC4MCuDerivativeData &dfC_4_MCu) const
void dEval(double const dt, EffectiveVolumetricInternalEnergyDerivatives const &effective_volumetric_internal_energy_d_data, FT1DerivativeData &dfT_1) const
Definition TEquation.cpp:33
void eval(double const dt, InternalEnergyData const &internal_energy_data, PrevState< InternalEnergyData > const &internal_energy_data_prev, FT1Data &fT_1) const
Definition TEquation.cpp:16
void dEval(BiotData const &biot_data, BishopsData const &chi_S_L, CapillaryPressureData const &p_cap, SaturationDataDeriv const &dS_L_dp_cap, FU2KUpCDerivativeData &dfu_2_KupC) const
Definition UEquation.cpp:36
void eval(BiotData const &biot_data, BishopsData const &chi_S_L, FU2KUpCData &fu_2_KupC) const
Definition UEquation.cpp:29
void dEval(BiotData const &biot_data, CapillaryPressureData const pCap, ConstituentDensityData const &constituent_density_data, PhaseTransitionData const &phase_transition_data, PorosityData const &porosity_data, PorosityDerivativeData const &porosity_d_data, PureLiquidDensityData const &rho_W_LR, SaturationData const &S_L_data, SaturationDataDeriv const &dS_L_dp_cap, SolidCompressibilityData const &beta_p_SR, FW2DerivativeData &dfW_2) const
Definition WEquation.cpp:48
void eval(BiotData const biot_data, CapillaryPressureData const pCap, ConstituentDensityData const &constituent_density_data, PorosityData const &porosity_data, PureLiquidDensityData const &rho_W_LR, SaturationData const &S_L_data, SolidCompressibilityData const beta_p_SR, FW2Data &fW_2) const
Definition WEquation.cpp:29
void dEval(double const dt, ConstituentDensityData const &constituent_density_data, PhaseTransitionData const &phase_transition_data, PrevState< ConstituentDensityData > const &constituent_density_data_prev, PrevState< PureLiquidDensityData > const &rho_W_LR_prev, PureLiquidDensityData const &rho_W_LR, SaturationData const &S_L_data, SaturationDataDeriv const &dS_L_dp_cap, FW3aDerivativeData &dfW_3a) const
void eval(double const dt, ConstituentDensityData const &constituent_density_data, PrevState< ConstituentDensityData > const &constituent_density_data_prev, PrevState< PureLiquidDensityData > const &rho_W_LR_prev, PureLiquidDensityData const &rho_W_LR, SaturationData const &S_L_data, FW3aData &fW_3a) const
void eval(BiotData const &biot_data, CapillaryPressureData const pCap, ConstituentDensityData const &constituent_density_data, PorosityData const &porosity_data, PrevState< SaturationData > const &S_L_data_prev, PureLiquidDensityData const &rho_W_LR, SaturationData const &S_L_data, SolidCompressibilityData const &beta_p_SR, FW4MWpCData &fW_4_MWpC) const
void eval(BiotData const &biot_data, ConstituentDensityData const &constituent_density_data, PorosityData const &porosity_data, PureLiquidDensityData const &rho_W_LR, SaturationData const &S_L_data, SolidCompressibilityData const &beta_p_SR, FW4MWpGData &fW_4_MWpG) const
void eval(BiotData const &biot_data, ConstituentDensityData const &constituent_density_data, PureLiquidDensityData const &rho_W_LR, SaturationData const &S_L_data, FW4MWuData &fW_4_MWu) const
void eval(FluidDensityData const &fluid_density_data, PhaseTransitionData const &phase_transition_data, PorosityData const &porosity_data, SaturationData const &S_L_data, SolidDensityData const &solid_density_data, SolidEnthalpyData const &solid_enthalpy_data, InternalEnergyData &internal_energy_data) const
void dEval(FluidDensityData const &fluid_density_data, PhaseTransitionData const &phase_transition_data, PorosityData const &porosity_data, PorosityDerivativeData const &porosity_d_data, SaturationData const &S_L_data, SolidDensityData const &solid_density_data, SolidDensityDerivativeData const &solid_density_d_data, SolidEnthalpyData const &solid_enthalpy_data, SolidHeatCapacityData const &solid_heat_capacity_data, EffectiveVolumetricInternalEnergyDerivatives &effective_volumetric_internal_energy_d_data) const
virtual void eval(SpaceTimeData const &x_t, MediaData const &media_data, GasPressureData const &p_GR, CapillaryPressureData const &p_cap, TemperatureData const &T_data, PureLiquidDensityData const &rho_W_LR, FluidEnthalpyData &fluid_enthalpy_data, MassMoleFractionsData &mass_mole_fractions_data, FluidDensityData &fluid_density_data, VapourPartialPressureData &vapour_pressure_data, ConstituentDensityData &constituent_density_data, PhaseTransitionData &cv) const =0
void eval(SpaceTimeData const &x_t, MediaData const &media_data, GasPressureData const &p_GR, CapillaryPressureData const &p_cap, TemperatureData const &T_data, PureLiquidDensityData &out) const
void dEval(SpaceTimeData const &x_t, MediaData const &media_data, CapillaryPressureData const &p_cap, SaturationDataDeriv &dS_L_data) const
void eval(SpaceTimeData const &x_t, MediaData const &media_data, CapillaryPressureData const &p_cap, SaturationData &S_L_data) const
void eval(SolidHeatCapacityData const &solid_heat_capacity_data, TemperatureData const &T_data, SolidEnthalpyData &solid_enthalpy_data) const
Definition Enthalpy.cpp:92
void eval(SpaceTimeData const &x_t, MediaData const &media_data, TemperatureData const &T_data, SolidHeatCapacityData &solid_heat_capacity) const
void eval(SpaceTimeData const &x_t, MediaData const &media_data, TemperatureData const &T_data, MassMoleFractionsData const &mass_mole_fractions_data, ViscosityData &viscosity_data) const
Definition Viscosity.cpp:16
NumLib::GenericIntegrationMethod const & integration_method_
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N_u
Definition TH2MFEM.h:45