OGS
ThermoRichardsMechanicsFEM-impl.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) OpenGeoSys Community (opengeosys.org)
2// SPDX-License-Identifier: BSD-3-Clause
3
4#pragma once
5
6#include <spdlog/fmt/fmt.h>
7
8#include <cassert>
9#include <type_traits>
10
24
25namespace ProcessLib
26{
28{
29template <typename ShapeFunctionDisplacement, typename ShapeFunction,
30 int DisplacementDim, typename ConstitutiveTraits>
31ThermoRichardsMechanicsLocalAssembler<ShapeFunctionDisplacement, ShapeFunction,
32 DisplacementDim, ConstitutiveTraits>::
33 ThermoRichardsMechanicsLocalAssembler(
34 MeshLib::Element const& e,
35 std::size_t const /*local_matrix_size*/,
36 NumLib::GenericIntegrationMethod const& integration_method,
37 bool const is_axially_symmetric,
39 process_data)
40 : LocalAssemblerInterface<DisplacementDim, ConstitutiveTraits>(
41 e, integration_method, is_axially_symmetric, process_data)
42{
43 unsigned const n_integration_points =
44 this->integration_method_.getNumberOfPoints();
45
46 ip_data_.resize(n_integration_points);
47
48 auto const shape_matrices_u =
49 NumLib::initShapeMatrices<ShapeFunctionDisplacement,
51 DisplacementDim>(e, is_axially_symmetric,
52 integration_method);
53
54 auto const shape_matrices =
56 DisplacementDim>(e, is_axially_symmetric,
57 integration_method);
58
59 for (unsigned ip = 0; ip < n_integration_points; ip++)
60 {
61 auto& ip_data = ip_data_[ip];
62 auto const& sm_u = shape_matrices_u[ip];
63 ip_data_[ip].integration_weight =
64 integration_method.getWeightedPoint(ip).getWeight() *
65 sm_u.integralMeasure * sm_u.detJ;
66
67 ip_data.N_u = sm_u.N;
68 ip_data.dNdx_u = sm_u.dNdx;
69
70 // ip_data.N_p and ip_data.dNdx_p are used for both p and T variables
71 ip_data.N_p = shape_matrices[ip].N;
72 ip_data.dNdx_p = shape_matrices[ip].dNdx;
73 }
74}
75
76template <typename ShapeFunctionDisplacement, typename ShapeFunction,
77 int DisplacementDim, typename ConstitutiveTraits>
79 ShapeFunctionDisplacement, ShapeFunction, DisplacementDim,
80 ConstitutiveTraits>::setInitialConditionsConcrete(Eigen::VectorXd const
81 local_x,
82 double const t,
83 int const /*process_id*/)
84{
85 assert(local_x.size() ==
87
88 auto const p_L = local_x.template segment<pressure_size>(pressure_index);
89 auto const T =
90 local_x.template segment<temperature_size>(temperature_index);
91
92 constexpr double dt = std::numeric_limits<double>::quiet_NaN();
93 auto const& medium =
94 *this->process_data_.media_map.getMedium(this->element_.getID());
95
96 MediaData const media_data{medium};
97
98 typename ConstitutiveTraits::ConstitutiveSetting const constitutive_setting;
99 auto models = ConstitutiveTraits::createConstitutiveModels(
100 this->process_data_, this->solid_material_);
101
102 unsigned const n_integration_points =
103 this->integration_method_.getNumberOfPoints();
104 for (unsigned ip = 0; ip < n_integration_points; ip++)
105 {
106 // N is used for both T and p variables.
107 auto const& N = ip_data_[ip].N_p;
108
109 ParameterLib::SpatialPosition const x_position{
110 std::nullopt, this->element_.getID(),
112 NumLib::interpolateCoordinates<ShapeFunctionDisplacement,
114 this->element_, ip_data_[ip].N_u))};
115
116 double T_ip;
118
119 double p_cap_ip;
120 NumLib::shapeFunctionInterpolate(-p_L, N, p_cap_ip);
121
122 MPL::VariableArray variables;
123
124 if (medium.hasProperty(MPL::PropertyType::saturation))
125 {
126 variables.capillary_pressure = p_cap_ip;
127 variables.liquid_phase_pressure = -p_cap_ip;
128 // setting pG to 1 atm
129 // TODO : rewrite equations s.t. p_L = pG-p_cap
130 variables.gas_phase_pressure = 1.0e5;
131 variables.temperature = T_ip;
132
133 double const S_L =
134 medium.property(MPL::PropertyType::saturation)
135 .template value<double>(variables, x_position, t, dt);
136 std::get<PrevState<SaturationData>>(this->prev_states_[ip])->S_L =
137 S_L;
138
139 variables.liquid_saturation = S_L;
140 }
141
142 constitutive_setting.init(models, t, dt, x_position, media_data,
143 {T_ip, 0, {}}, this->current_states_[ip],
144 this->prev_states_[ip]);
145
146 if (this->process_data_.initial_stress.value)
147 {
148 if (!medium.hasProperty(MPL::PropertyType::saturation))
149 {
150 OGS_FATAL(
151 "The medium has no saturation property required to compute "
152 "initial stress.");
153 }
154 convertInitialStressType(ip, t, x_position, medium, variables,
155 -p_cap_ip);
156 }
157 }
158}
159
160template <typename ShapeFunctionDisplacement, typename ShapeFunction,
161 int DisplacementDim, typename ConstitutiveTraits>
162void ThermoRichardsMechanicsLocalAssembler<ShapeFunctionDisplacement,
163 ShapeFunction, DisplacementDim,
164 ConstitutiveTraits>::
165 convertInitialStressType(unsigned const ip,
166 double const t,
167 ParameterLib::SpatialPosition const x_position,
168 MaterialPropertyLib::Medium const& medium,
169 MPL::VariableArray const& variables,
170 double const p_at_ip)
171{
172 bool constexpr is_strain_temperature_constitutive =
174 DisplacementDim>,
175 ConstitutiveTraits>::value;
176 if (is_strain_temperature_constitutive &&
177 this->process_data_.initial_stress.type ==
179 {
180 return;
181 }
182
183 if (!is_strain_temperature_constitutive &&
184 this->process_data_.initial_stress.type == InitialStress::Type::Total)
185 {
186 return;
187 }
188
189 double const alpha_b =
191 .template value<double>(variables, x_position, t, 0.0 /*dt*/);
192
193 double const bishop =
195 .template value<double>(variables, x_position, t, 0.0 /*dt*/);
196
197 ConstitutiveTraits::ConstitutiveSetting::convertInitialStressType(
198 this->current_states_[ip], this->prev_states_[ip],
199 bishop * alpha_b * p_at_ip * Invariants::identity2);
200}
201
202template <typename ShapeFunctionDisplacement, typename ShapeFunction,
203 int DisplacementDim, typename ConstitutiveTraits>
204void ThermoRichardsMechanicsLocalAssembler<ShapeFunctionDisplacement,
205 ShapeFunction, DisplacementDim,
206 ConstitutiveTraits>::
207 assembleWithJacobian(double const t, double const dt,
208 std::vector<double> const& local_x,
209 std::vector<double> const& local_x_prev,
210 std::vector<double>& local_rhs_data,
211 std::vector<double>& local_Jac_data)
212{
213 auto& medium =
214 *this->process_data_.media_map.getMedium(this->element_.getID());
215
216 LocalMatrices loc_mat;
217 loc_mat.setZero();
218 LocalMatrices loc_mat_current_ip;
219 loc_mat_current_ip.setZero(); // only to set the right matrix sizes
220
221 typename ConstitutiveTraits::ConstitutiveSetting constitutive_setting;
222
223 for (unsigned ip = 0; ip < this->integration_method_.getNumberOfPoints();
224 ++ip)
225 {
226 ParameterLib::SpatialPosition const x_position{
227 std::nullopt, this->element_.getID(),
229 NumLib::interpolateCoordinates<ShapeFunctionDisplacement,
231 this->element_, ip_data_[ip].N_u))};
232
234 t, dt, x_position, //
235 local_x, local_x_prev, //
236 ip_data_[ip], constitutive_setting,
237 medium, //
238 loc_mat_current_ip, //
239 this->current_states_[ip], this->prev_states_[ip],
240 this->material_states_[ip], this->output_data_[ip]);
241 loc_mat += loc_mat_current_ip;
242 }
243
244 massLumping(loc_mat);
245
246 addToLocalMatrixData(dt, local_x, local_x_prev, loc_mat, local_rhs_data,
247 local_Jac_data);
248}
249
250template <typename ShapeFunctionDisplacement, typename ShapeFunction,
251 int DisplacementDim, typename ConstitutiveTraits>
252void ThermoRichardsMechanicsLocalAssembler<ShapeFunctionDisplacement,
253 ShapeFunction, DisplacementDim,
254 ConstitutiveTraits>::
255 massLumping(typename ThermoRichardsMechanicsLocalAssembler<
256 ShapeFunctionDisplacement, ShapeFunction, DisplacementDim,
257 ConstitutiveTraits>::LocalMatrices& loc_mat) const
258{
259 if (this->process_data_.apply_mass_lumping)
260 {
261 loc_mat.storage_p_a_p =
262 loc_mat.storage_p_a_p.colwise().sum().eval().asDiagonal();
263 loc_mat.storage_p_a_S =
264 loc_mat.storage_p_a_S.colwise().sum().eval().asDiagonal();
265 loc_mat.storage_p_a_S_Jpp =
266 loc_mat.storage_p_a_S_Jpp.colwise().sum().eval().asDiagonal();
267 }
268}
269
270template <typename ShapeFunctionDisplacement, typename ShapeFunction,
271 int DisplacementDim, typename ConstitutiveTraits>
272void ThermoRichardsMechanicsLocalAssembler<ShapeFunctionDisplacement,
273 ShapeFunction, DisplacementDim,
274 ConstitutiveTraits>::
275 addToLocalMatrixData(
276 double const dt,
277 std::vector<double> const& local_x,
278 std::vector<double> const& local_x_prev,
280 ShapeFunctionDisplacement, ShapeFunction, DisplacementDim,
281 ConstitutiveTraits>::LocalMatrices const& loc_mat,
282 std::vector<double>& local_rhs_data,
283 std::vector<double>& local_Jac_data) const
284{
285 constexpr auto local_matrix_dim =
287 assert(local_x.size() == local_matrix_dim);
288
289 auto local_Jac = MathLib::createZeroedMatrix<
290 typename ShapeMatricesTypeDisplacement::template MatrixType<
291 local_matrix_dim, local_matrix_dim>>(
292 local_Jac_data, local_matrix_dim, local_matrix_dim);
293
294 auto local_rhs =
295 MathLib::createZeroedVector<typename ShapeMatricesTypeDisplacement::
296 template VectorType<local_matrix_dim>>(
297 local_rhs_data, local_matrix_dim);
298
299 local_Jac.noalias() = loc_mat.Jac;
300 local_rhs.noalias() = -loc_mat.res;
301
302 //
303 // -- Jacobian
304 //
305 block_TT(local_Jac).noalias() += loc_mat.M_TT / dt + loc_mat.K_TT;
306 block_Tp(local_Jac).noalias() +=
307 loc_mat.K_Tp + loc_mat.dK_TT_dp + loc_mat.M_Tp / dt;
308
309 block_pT(local_Jac).noalias() += loc_mat.M_pT / dt + loc_mat.K_pT;
310 block_pp(local_Jac).noalias() +=
311 loc_mat.K_pp + loc_mat.storage_p_a_p / dt + loc_mat.storage_p_a_S_Jpp;
312 block_pu(local_Jac).noalias() = loc_mat.M_pu / dt;
313
314 //
315 // -- Residual
316 //
317 auto const [T, p_L, u] = localDOF(local_x);
318 auto const [T_prev, p_L_prev, u_prev] = localDOF(local_x_prev);
319
320 block_T(local_rhs).noalias() -= loc_mat.M_TT * (T - T_prev) / dt +
321 loc_mat.K_TT * T + loc_mat.K_Tp * p_L +
322 loc_mat.M_Tp * (p_L - p_L_prev) / dt;
323 block_p(local_rhs).noalias() -=
324 loc_mat.K_pp * p_L + loc_mat.K_pT * T +
325 (loc_mat.storage_p_a_p + loc_mat.storage_p_a_S) * (p_L - p_L_prev) /
326 dt +
327 loc_mat.M_pu * (u - u_prev) / dt + loc_mat.M_pT * (T - T_prev) / dt;
328}
329
330template <typename ShapeFunctionDisplacement, typename ShapeFunction,
331 int DisplacementDim, typename ConstitutiveTraits>
332void ThermoRichardsMechanicsLocalAssembler<ShapeFunctionDisplacement,
333 ShapeFunction, DisplacementDim,
334 ConstitutiveTraits>::
335 assembleWithJacobianSingleIP(
336 double const t, double const dt,
337 ParameterLib::SpatialPosition const& x_position,
338 std::vector<double> const& local_x,
339 std::vector<double> const& local_x_prev,
341 ShapeFunctionDisplacement, ShapeFunction, DisplacementDim,
342 ConstitutiveTraits>::IpData const& ip_data,
343 typename ConstitutiveTraits::ConstitutiveSetting& CS,
346 ShapeFunctionDisplacement, ShapeFunction, DisplacementDim,
347 ConstitutiveTraits>::LocalMatrices& out,
348 typename ConstitutiveTraits::StatefulData& current_state,
349 typename ConstitutiveTraits::StatefulDataPrev const& prev_state,
351 typename ConstitutiveTraits::OutputData& output_data) const
352{
353 auto const& N_u = ip_data.N_u;
354 auto const& dNdx_u = ip_data.dNdx_u;
355
356 // N and dNdx are used for both p and T variables
357 auto const& N = ip_data.N_p;
358 auto const& dNdx = ip_data.dNdx_p;
359
360 auto const B =
361 LinearBMatrix::computeBMatrix<DisplacementDim,
362 ShapeFunctionDisplacement::NPOINTS,
364 dNdx_u, N_u, (*x_position.getCoordinates())[0],
366
367 typename ConstitutiveTraits::ConstitutiveData CD;
368
369 auto const [T, p_L, u] = localDOF(local_x);
370 auto const [T_prev, p_L_prev, u_prev] = localDOF(local_x_prev);
371
372 {
373 auto models = ConstitutiveTraits::createConstitutiveModels(
374 this->process_data_, this->solid_material_);
375 typename ConstitutiveTraits::ConstitutiveTempData tmp;
376
377 double const T_ip = N * T;
378 double const T_prev_ip = N * T_prev;
379 GlobalDimVectorType const grad_T_ip = dNdx * T;
380
381 double const p_cap_ip = -N * p_L;
382 double const p_cap_prev_ip = -N * p_L_prev;
383 GlobalDimVectorType const grad_p_cap_ip = -dNdx * p_L;
384
385 KelvinVectorType eps = B * u;
386
387 CS.eval(models, t, dt, x_position, //
388 medium, //
389 {T_ip, T_prev_ip, grad_T_ip}, //
390 {p_cap_ip, p_cap_prev_ip, grad_p_cap_ip}, //
391 eps, current_state, prev_state, mat_state, tmp, output_data,
392 CD);
393 }
394
395 using NodalMatrix = typename ShapeMatricesType::NodalMatrixType;
396 NodalMatrix const NTN = N.transpose() * N;
397 NodalMatrix const dNTdN = dNdx.transpose() * dNdx;
398
399 // TODO is identity2.transpose() * B something like divergence?
400 auto const& identity2 = MathLib::KelvinVector::Invariants<
402 DisplacementDim)>::identity2;
403 typename ShapeMatricesTypeDisplacement::template MatrixType<
405 BTI2N = B.transpose() * identity2 * N;
406
407 /*
408 * Conventions:
409 *
410 * * use positive signs exclusively, any negative signs must be included in
411 * the coefficients coming from the constitutive setting
412 * * the used coefficients from the constitutive setting are named after the
413 * terms they appear in, e.g. K_TT_X_dNTdN means it is a scalar (X) that
414 * will be multiplied by dNdx.transpose() * dNdx. Placefolders for the
415 * coefficients are:
416 * * X -> scalar
417 * * V -> vector
418 * * K -> Kelvin vector
419 * * the Laplace terms have a special name, e.g., K_TT_Laplace
420 * * there shall be only one contribution to each of the LocalMatrices,
421 * assigned with = assignment; this point might be relaxed in the future
422 * * this method will overwrite the data in the passed LocalMatrices& out
423 * argument, not add to it
424 */
425
426 // residual, order T, p, u
427 block_p(out.res).noalias() =
428 dNdx.transpose() * std::get<EqPData<DisplacementDim>>(CD).rhs_p_dNT_V;
429 block_u(out.res).noalias() =
430 B.transpose() *
432 CD, current_state)
433 .sigma_total -
434 static_cast<int>(this->process_data_.apply_body_force_for_deformation) *
435 N_u_op(N_u).transpose() *
436 std::get<GravityData<DisplacementDim>>(CD).volumetric_body_force;
437
438 // Storage matrices
439 out.storage_p_a_p.noalias() =
440 std::get<EqPData<DisplacementDim>>(CD).storage_p_a_p_X_NTN * NTN;
441 out.storage_p_a_S.noalias() =
442 std::get<TRMStorageData>(CD).storage_p_a_S_X_NTN * NTN;
443 out.storage_p_a_S_Jpp.noalias() =
444 std::get<TRMStorageData>(CD).storage_p_a_S_Jpp_X_NTN * NTN;
445
446 // M matrices, order T, p, u
447 out.M_TT.noalias() =
448 std::get<EqTData<DisplacementDim>>(CD).M_TT_X_NTN * NTN;
449 out.M_Tp.noalias() =
450 std::get<TRMVaporDiffusionData<DisplacementDim>>(CD).M_Tp_X_NTN * NTN;
451
452 out.M_pT.noalias() =
453 std::get<EqPData<DisplacementDim>>(CD).M_pT_X_NTN * NTN;
454 out.M_pu.noalias() =
455 std::get<EqPData<DisplacementDim>>(CD).M_pu_X_BTI2N * BTI2N.transpose();
456
457 // K matrices, order T, p, u
458 out.K_TT.noalias() =
459 dNdx.transpose() *
460 std::get<TRMHeatStorageAndFluxData<DisplacementDim>>(CD)
461 .K_TT_Laplace *
462 dNdx +
463 N.transpose() *
464 (std::get<EqTData<DisplacementDim>>(CD).K_TT_NT_V_dN.transpose() *
465 dNdx) +
466 std::get<TRMVaporDiffusionData<DisplacementDim>>(CD).K_TT_X_dNTdN *
467 dNTdN;
468
469 out.dK_TT_dp.noalias() =
470 N.transpose() *
471 (std::get<TRMHeatStorageAndFluxData<DisplacementDim>>(CD)
472 .K_Tp_NT_V_dN.transpose() *
473 dNdx) +
474 std::get<TRMHeatStorageAndFluxData<DisplacementDim>>(CD).K_Tp_X_NTN *
475 NTN;
476 out.K_Tp.noalias() =
477 dNdx.transpose() *
478 std::get<ThermoOsmosisData<DisplacementDim>>(CD).K_Tp_Laplace *
479 dNdx +
481 dNTdN;
482
483 out.K_pp.noalias() =
484 dNdx.transpose() * std::get<EqPData<DisplacementDim>>(CD).K_pp_Laplace *
485 dNdx +
487 dNTdN;
488 out.K_pT.noalias() =
489 dNdx.transpose() *
490 std::get<ThermoOsmosisData<DisplacementDim>>(CD).K_pT_Laplace * dNdx;
491
492 // direct Jacobian contributions, order T, p, u
493 block_pT(out.Jac).noalias() =
494 std::get<TRMVaporDiffusionData<DisplacementDim>>(CD).J_pT_X_dNTdN *
495 dNTdN;
496 block_pp(out.Jac).noalias() =
497 std::get<TRMStorageData>(CD).J_pp_X_NTN * NTN +
498 std::get<EqPData<DisplacementDim>>(CD).J_pp_X_BTI2NT_u_dot_N *
499 BTI2N.transpose() * (u - u_prev) / dt *
500 N // TODO something with volumetric strain rate?
501 + dNdx.transpose() *
502 std::get<EqPData<DisplacementDim>>(CD).J_pp_dNT_V_N * N;
503
504 block_uT(out.Jac).noalias() =
505 B.transpose() *
506 std::get<SolidMechanicsDataStateless<DisplacementDim>>(CD).J_uT_BT_K_N *
507 N;
508 block_up(out.Jac).noalias() =
509 B.transpose() *
510 std::get<SolidMechanicsDataStateless<DisplacementDim>>(CD)
511 .J_up_BT_K_N *
512 N +
513 N_u_op(N_u).transpose() *
514 std::get<GravityData<DisplacementDim>>(CD).J_up_HT_V_N * N;
515 block_uu(out.Jac).noalias() =
516 B.transpose() *
517 std::get<SolidMechanicsDataStateless<DisplacementDim>>(CD)
518 .stiffness_tensor *
519 B;
520
521 out *= ip_data.integration_weight;
522}
523
524template <typename ShapeFunctionDisplacement, typename ShapeFunction,
525 int DisplacementDim, typename ConstitutiveTraits>
526void ThermoRichardsMechanicsLocalAssembler<ShapeFunctionDisplacement,
527 ShapeFunction, DisplacementDim,
528 ConstitutiveTraits>::
529 computeSecondaryVariableConcrete(double const t, double const dt,
530 Eigen::VectorXd const& local_x,
531 Eigen::VectorXd const& local_x_prev)
532{
533 auto const T = block_T(local_x);
534 auto const p_L = block_p(local_x);
535 auto const u = block_u(local_x);
536
537 auto const T_prev = block_T(local_x_prev);
538 auto const p_L_prev = block_p(local_x_prev);
539
540 auto const e_id = this->element_.getID();
541 auto const& process_data = this->process_data_;
542 auto& medium = *process_data.media_map.getMedium(e_id);
543
544 unsigned const n_integration_points =
545 this->integration_method_.getNumberOfPoints();
546
547 typename ConstitutiveTraits::ConstitutiveSetting constitutive_setting;
548
549 auto models = ConstitutiveTraits::createConstitutiveModels(
550 process_data, this->solid_material_);
551 typename ConstitutiveTraits::ConstitutiveTempData tmp;
552 typename ConstitutiveTraits::ConstitutiveData CD;
553
554 for (unsigned ip = 0; ip < n_integration_points; ip++)
555 {
556 auto& current_state = this->current_states_[ip];
557 auto& output_data = this->output_data_[ip];
558
559 auto const& ip_data = ip_data_[ip];
560
561 // N is used for both p and T variables
562 auto const& N = ip_data.N_p;
563 auto const& N_u = ip_data.N_u;
564 auto const& dNdx_u = ip_data.dNdx_u;
565 auto const& dNdx = ip_data.dNdx_p;
566
567 ParameterLib::SpatialPosition const x_position{
568 std::nullopt, this->element_.getID(),
570 NumLib::interpolateCoordinates<ShapeFunctionDisplacement,
572 this->element_, N_u))};
573
574 auto const x_coord =
575 x_position.getCoordinates().value()[0]; // r for axisymmetry
576 auto const B =
577 LinearBMatrix::computeBMatrix<DisplacementDim,
578 ShapeFunctionDisplacement::NPOINTS,
580 dNdx_u, N_u, x_coord, this->is_axially_symmetric_);
581
582 double const T_ip = N * T;
583 double const T_prev_ip = N * T_prev;
584 GlobalDimVectorType const grad_T_ip = dNdx * T;
585
586 double const p_cap_ip = -N * p_L;
587 double const p_cap_prev_ip = -N * p_L_prev;
588 GlobalDimVectorType const grad_p_cap_ip = -dNdx * p_L;
589
590 KelvinVectorType eps = B * u;
591
592 constitutive_setting.eval(models, //
593 t, dt, x_position, //
594 medium, //
595 {T_ip, T_prev_ip, grad_T_ip}, //
596 {p_cap_ip, p_cap_prev_ip, grad_p_cap_ip}, //
597 eps, current_state, this->prev_states_[ip],
598 this->material_states_[ip], tmp, output_data,
599 CD);
600 }
601
603 ShapeFunction, typename ShapeFunctionDisplacement::MeshElement,
604 DisplacementDim>(this->element_, this->is_axially_symmetric_, p_L,
605 *process_data.pressure_interpolated);
607 ShapeFunction, typename ShapeFunctionDisplacement::MeshElement,
608 DisplacementDim>(this->element_, this->is_axially_symmetric_, T,
609 *process_data.temperature_interpolated);
610}
611} // namespace ThermoRichardsMechanics
612} // namespace ProcessLib
#define OGS_FATAL(...)
Definition Error.h:19
Property const & property(PropertyType const &p) const
Definition Medium.cpp:49
constexpr double getWeight() const
MathLib::WeightedPoint const & getWeightedPoint(unsigned const igp) const
std::optional< MathLib::Point3d > const getCoordinates() const
MatrixType< _kelvin_vector_size, _number_of_dof > BMatrixType
void addToLocalMatrixData(double const dt, std::vector< double > const &local_x, std::vector< double > const &local_x_prev, LocalMatrices const &loc_mat, std::vector< double > &local_rhs_data, std::vector< double > &local_Jac_data) const
void assembleWithJacobianSingleIP(double const t, double const dt, ParameterLib::SpatialPosition const &x_position, std::vector< double > const &local_x, std::vector< double > const &local_x_prev, IpData const &ip_data, typename ConstitutiveTraits::ConstitutiveSetting &CS, MaterialPropertyLib::Medium &medium, LocalMatrices &out, typename ConstitutiveTraits::StatefulData &current_state, typename ConstitutiveTraits::StatefulDataPrev const &prev_state, MaterialStateData< DisplacementDim > &mat_state, typename ConstitutiveTraits::OutputData &output_data) const
ShapeMatrixPolicyType< ShapeFunctionDisplacement, DisplacementDim > ShapeMatricesTypeDisplacement
void convertInitialStressType(unsigned const ip, double const t, ParameterLib::SpatialPosition const x_position, MaterialPropertyLib::Medium const &medium, MPL::VariableArray const &variables, double const p_at_ip)
ThermoRichardsMechanicsLocalAssembler(ThermoRichardsMechanicsLocalAssembler const &)=delete
void setInitialConditionsConcrete(Eigen::VectorXd const local_x, double const t, int const process_id) override
ShapeMatrixPolicyType< ShapeFunction, DisplacementDim > ShapeMatricesType
IntegrationPointData< ShapeMatricesTypeDisplacement, ShapeMatricesType, DisplacementDim, ShapeFunctionDisplacement::NPOINTS > IpData
constexpr int kelvin_vector_dimensions(int const displacement_dim)
Kelvin vector dimensions for given displacement dimension.
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 shapeFunctionInterpolate(const NodalValues &, const ShapeMatrix &)
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)
auto & get(Tuples &... ts)
Definition Get.h:53
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.
MatrixType< ShapeFunction::NPOINTS, ShapeFunction::NPOINTS > NodalMatrixType
static Eigen::Matrix< double, KelvinVectorSize, 1 > const identity2
Kelvin mapping of 2nd order identity tensor.
ThermoRichardsMechanicsProcessData< DisplacementDim, ConstitutiveTraits > & process_data_
LocalAssemblerInterface(MeshLib::Element const &e, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, ThermoRichardsMechanicsProcessData< DisplacementDim, ConstitutiveTraits > &process_data)