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