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