OGS
ThermoMechanicsFEM-impl.h
Go to the documentation of this file.
1
13#pragma once
14
21
22namespace ProcessLib
23{
24namespace ThermoMechanics
25{
26template <typename ShapeFunction, int DisplacementDim>
29 MeshLib::Element const& e,
30 std::size_t const /*local_matrix_size*/,
31 NumLib::GenericIntegrationMethod const& integration_method,
32 bool const is_axially_symmetric,
34 : _process_data(process_data),
35 _integration_method(integration_method),
36 _element(e),
37 _is_axially_symmetric(is_axially_symmetric)
38{
39 unsigned const n_integration_points =
41
42 _ip_data.reserve(n_integration_points);
43 _secondary_data.N.resize(n_integration_points);
44
45 auto const shape_matrices =
47 DisplacementDim>(e, is_axially_symmetric,
49
51 _process_data.solid_materials, _process_data.material_ids, e.getID());
52
53 for (unsigned ip = 0; ip < n_integration_points; ip++)
54 {
55 _ip_data.emplace_back(solid_material);
56 auto& ip_data = _ip_data[ip];
57 ip_data.integration_weight =
59 shape_matrices[ip].integralMeasure * shape_matrices[ip].detJ;
60
61 static const int kelvin_vector_size =
63 // Initialize current time step values
64 ip_data.sigma.setZero(kelvin_vector_size);
65 ip_data.eps.setZero(kelvin_vector_size);
66
67 // Previous time step values are not initialized and are set later.
68 ip_data.sigma_prev.resize(kelvin_vector_size);
69 ip_data.eps_prev.resize(kelvin_vector_size);
70
71 ip_data.eps_m.setZero(kelvin_vector_size);
72 ip_data.eps_m_prev.setZero(kelvin_vector_size);
74 x_position.setElementID(_element.getID());
75 ip_data.N = shape_matrices[ip].N;
76 ip_data.dNdx = shape_matrices[ip].dNdx;
77
78 _secondary_data.N[ip] = shape_matrices[ip].N;
79 }
80}
81
82template <typename ShapeFunction, int DisplacementDim>
84 setIPDataInitialConditions(std::string_view const name,
85 double const* values,
86 int const integration_order)
87{
88 if (integration_order !=
89 static_cast<int>(_integration_method.getIntegrationOrder()))
90 {
92 "Setting integration point initial conditions; The integration "
93 "order of the local assembler for element {:d} is different from "
94 "the integration order in the initial condition.",
95 _element.getID());
96 }
97
98 if (name == "sigma")
99 {
100 return setSigma(values);
101 }
102 if (name == "epsilon")
103 {
104 return setEpsilon(values);
105 }
106 if (name == "epsilon_m")
107 {
108 return setEpsilonMechanical(values);
109 }
110
111 return 0;
112}
113
114template <typename ShapeFunction, int DisplacementDim>
116 assembleWithJacobian(double const t, double const dt,
117 std::vector<double> const& local_x,
118 std::vector<double> const& local_x_prev,
119 std::vector<double>& /*local_M_data*/,
120 std::vector<double>& /*local_K_data*/,
121 std::vector<double>& local_rhs_data,
122 std::vector<double>& local_Jac_data)
123{
124 auto const local_matrix_size = local_x.size();
125 assert(local_matrix_size == temperature_size + displacement_size);
126
127 auto T = Eigen::Map<typename ShapeMatricesType::template VectorType<
128 temperature_size> const>(local_x.data() + temperature_index,
129 temperature_size);
130
131 auto u = Eigen::Map<typename ShapeMatricesType::template VectorType<
132 displacement_size> const>(local_x.data() + displacement_index,
133 displacement_size);
134 bool const is_u_non_zero = u.norm() > 0.0;
135
136 auto T_prev = Eigen::Map<typename ShapeMatricesType::template VectorType<
137 temperature_size> const>(local_x_prev.data() + temperature_index,
138 temperature_size);
139
140 auto local_Jac = MathLib::createZeroedMatrix<JacobianMatrix>(
141 local_Jac_data, local_matrix_size, local_matrix_size);
142
143 auto local_rhs = MathLib::createZeroedVector<RhsVector>(local_rhs_data,
144 local_matrix_size);
145
146 typename ShapeMatricesType::template MatrixType<displacement_size,
147 temperature_size>
148 KuT;
149 KuT.setZero(displacement_size, temperature_size);
150
152 KTT.setZero(temperature_size, temperature_size);
153
155 DTT.setZero(temperature_size, temperature_size);
156
157 unsigned const n_integration_points =
158 _integration_method.getNumberOfPoints();
159
160 MPL::VariableArray variables;
161 MPL::VariableArray variables_prev;
163 x_position.setElementID(_element.getID());
164
165 auto const& medium = _process_data.media_map.getMedium(_element.getID());
166 auto const& solid_phase = medium->phase("Solid");
167
168 for (unsigned ip = 0; ip < n_integration_points; ip++)
169 {
170 x_position.setIntegrationPoint(ip);
171 auto const& w = _ip_data[ip].integration_weight;
172 auto const& N = _ip_data[ip].N;
173 auto const& dNdx = _ip_data[ip].dNdx;
174
175 auto const x_coord =
176 NumLib::interpolateXCoordinate<ShapeFunction, ShapeMatricesType>(
177 _element, N);
178 auto const& B =
179 LinearBMatrix::computeBMatrix<DisplacementDim,
180 ShapeFunction::NPOINTS,
182 dNdx, N, x_coord, _is_axially_symmetric);
183
184 auto& sigma = _ip_data[ip].sigma;
185 auto const& sigma_prev = _ip_data[ip].sigma_prev;
186
187 auto& eps = _ip_data[ip].eps;
188 auto const& eps_prev = _ip_data[ip].eps_prev;
189
190 auto& eps_m = _ip_data[ip].eps_m;
191 auto const& eps_m_prev = _ip_data[ip].eps_m_prev;
192
193 auto& state = _ip_data[ip].material_state_variables;
194
195 const double T_ip = N.dot(T); // T at integration point
196 double const T_prev_ip = N.dot(T_prev);
197
198 // Consider also anisotropic thermal expansion.
199 auto const solid_linear_thermal_expansivity_vector =
200 MPL::formKelvinVector<DisplacementDim>(
201 solid_phase
202 .property(
204 .value(variables, x_position, t, dt));
205
207 dthermal_strain =
208 solid_linear_thermal_expansivity_vector * (T_ip - T_prev_ip);
209
210 //
211 // displacement equation, displacement part
212 //
213 // For the restart computation, the displacement may not be
214 // reloaded but the initial strains are always available. For such case,
215 // the following computation is skipped.
216 if (is_u_non_zero)
217 {
218 eps.noalias() = B * u;
219 }
220
221 eps_m.noalias() = eps_m_prev + eps - eps_prev - dthermal_strain;
222
223 variables_prev.stress
225 sigma_prev);
226 variables_prev.mechanical_strain
228 eps_m_prev);
229 variables_prev.temperature = T_ip;
230 variables.mechanical_strain
232 eps_m);
233 variables.temperature = T_ip;
234
235 auto&& solution = _ip_data[ip].solid_material.integrateStress(
236 variables_prev, variables, t, x_position, dt, *state);
237
238 if (!solution)
239 {
240 OGS_FATAL("Computation of local constitutive relation failed.");
241 }
242
244 std::tie(sigma, state, C) = std::move(*solution);
245
246 local_Jac
247 .template block<displacement_size, displacement_size>(
248 displacement_index, displacement_index)
249 .noalias() += B.transpose() * C * B * w;
250
251 typename ShapeMatricesType::template MatrixType<DisplacementDim,
252 displacement_size>
253 N_u = ShapeMatricesType::template MatrixType<
254 DisplacementDim, displacement_size>::Zero(DisplacementDim,
255 displacement_size);
256
257 for (int i = 0; i < DisplacementDim; ++i)
258 {
259 N_u.template block<1, displacement_size / DisplacementDim>(
260 i, i * displacement_size / DisplacementDim)
261 .noalias() = N;
262 }
263
264 auto const rho_s =
265 solid_phase.property(MPL::PropertyType::density)
266 .template value<double>(variables, x_position, t, dt);
267
268 auto const& b = _process_data.specific_body_force;
269 local_rhs.template block<displacement_size, 1>(displacement_index, 0)
270 .noalias() -=
271 (B.transpose() * sigma - N_u.transpose() * rho_s * b) * w;
272
273 //
274 // displacement equation, temperature part
275 // The computation of KuT can be ignored.
276 auto const alpha_T_tensor =
278 solid_linear_thermal_expansivity_vector);
279 KuT.noalias() += B.transpose() * (C * alpha_T_tensor) * N * w;
280
281 if (_ip_data[ip].solid_material.getConstitutiveModel() ==
283 {
284 auto const s = Invariants::deviatoric_projection * sigma;
285 double const norm_s = Invariants::FrobeniusNorm(s);
286 const double creep_coefficient =
287 _ip_data[ip].solid_material.getTemperatureRelatedCoefficient(
288 t, dt, x_position, T_ip, norm_s);
289 KuT.noalias() += creep_coefficient * B.transpose() * s * N * w;
290 }
291
292 //
293 // temperature equation, temperature part;
294 //
295 auto const lambda =
296 solid_phase
297 .property(
299 .value(variables, x_position, t, dt);
300
301 GlobalDimMatrixType const thermal_conductivity =
302 MaterialPropertyLib::formEigenTensor<DisplacementDim>(lambda);
303
304 KTT.noalias() += dNdx.transpose() * thermal_conductivity * dNdx * w;
305
306 auto const c =
307 solid_phase
308 .property(
310 .template value<double>(variables, x_position, t, dt);
311 DTT.noalias() += N.transpose() * rho_s * c * N * w;
312 }
313
314 // temperature equation, temperature part
315 local_Jac
316 .template block<temperature_size, temperature_size>(temperature_index,
317 temperature_index)
318 .noalias() += KTT + DTT / dt;
319
320 // displacement equation, temperature part
321 local_Jac
322 .template block<displacement_size, temperature_size>(displacement_index,
323 temperature_index)
324 .noalias() -= KuT;
325
326 local_rhs.template block<temperature_size, 1>(temperature_index, 0)
327 .noalias() -= KTT * T + DTT * (T - T_prev) / dt;
328}
329
330template <typename ShapeFunction, int DisplacementDim>
333 const double t, double const dt, Eigen::VectorXd const& local_x,
334 Eigen::VectorXd const& local_x_prev, int const process_id,
335 std::vector<double>& /*local_M_data*/,
336 std::vector<double>& /*local_K_data*/,
337 std::vector<double>& local_b_data, std::vector<double>& local_Jac_data)
338{
339 // For the equations with pressure
340 if (process_id == _process_data.heat_conduction_process_id)
341 {
342 assembleWithJacobianForHeatConductionEquations(
343 t, dt, local_x, local_x_prev, local_b_data, local_Jac_data);
344 return;
345 }
346
347 // For the equations with deformation
348 assembleWithJacobianForDeformationEquations(t, dt, local_x, local_x_prev,
349 local_b_data, local_Jac_data);
350}
351
352template <typename ShapeFunction, int DisplacementDim>
355 const double t, double const dt, Eigen::VectorXd const& local_x,
356 Eigen::VectorXd const& local_x_prev, std::vector<double>& local_b_data,
357 std::vector<double>& local_Jac_data)
358{
359 auto const local_T =
360 local_x.template segment<temperature_size>(temperature_index);
361
362 auto const local_T_prev =
363 local_x_prev.template segment<temperature_size>(temperature_index);
364
365 auto const local_u =
366 local_x.template segment<displacement_size>(displacement_index);
367 bool const is_u_non_zero = local_u.norm() > 0.0;
368
369 auto local_Jac = MathLib::createZeroedMatrix<
370 typename ShapeMatricesType::template MatrixType<displacement_size,
371 displacement_size>>(
372 local_Jac_data, displacement_size, displacement_size);
373
374 auto local_rhs = MathLib::createZeroedVector<
375 typename ShapeMatricesType::template VectorType<displacement_size>>(
376 local_b_data, displacement_size);
377
378 unsigned const n_integration_points =
379 _integration_method.getNumberOfPoints();
380
381 MPL::VariableArray variables;
382 MPL::VariableArray variables_prev;
384 x_position.setElementID(_element.getID());
385 auto const& medium = _process_data.media_map.getMedium(_element.getID());
386 auto const& solid_phase = medium->phase("Solid");
387
388 for (unsigned ip = 0; ip < n_integration_points; ip++)
389 {
390 x_position.setIntegrationPoint(ip);
391 auto const& w = _ip_data[ip].integration_weight;
392 auto const& N = _ip_data[ip].N;
393 auto const& dNdx = _ip_data[ip].dNdx;
394
395 auto const x_coord =
396 NumLib::interpolateXCoordinate<ShapeFunction, ShapeMatricesType>(
397 _element, N);
398 auto const& B =
399 LinearBMatrix::computeBMatrix<DisplacementDim,
400 ShapeFunction::NPOINTS,
402 dNdx, N, x_coord, _is_axially_symmetric);
403
404 auto& sigma = _ip_data[ip].sigma;
405 auto const& sigma_prev = _ip_data[ip].sigma_prev;
406
407 auto& eps = _ip_data[ip].eps;
408 auto const& eps_prev = _ip_data[ip].eps_prev;
409
410 auto& eps_m = _ip_data[ip].eps_m;
411 auto const& eps_m_prev = _ip_data[ip].eps_m_prev;
412
413 auto& state = _ip_data[ip].material_state_variables;
414
415 const double T_ip = N.dot(local_T); // T at integration point
416 variables.temperature = T_ip;
417 double const dT_ip = T_ip - N.dot(local_T_prev);
418
419 //
420 // displacement equation, displacement part
421 //
422 // For the restart computation, the displacement may not be
423 // reloaded but the initial strains are always available. For such case,
424 // the following computation is skipped.
425 if (is_u_non_zero)
426 {
427 eps.noalias() = B * local_u;
428 }
429
430 // Consider also anisotropic thermal expansion.
431 auto const solid_linear_thermal_expansivity_vector =
432 MPL::formKelvinVector<DisplacementDim>(
433 solid_phase
434 .property(
436 .value(variables, x_position, t, dt));
437
439 dthermal_strain = solid_linear_thermal_expansivity_vector * dT_ip;
440
441 eps_m.noalias() = eps_m_prev + eps - eps_prev - dthermal_strain;
442
443 variables_prev.stress
445 sigma_prev);
446 variables_prev.mechanical_strain
448 eps_m_prev);
449 variables_prev.temperature = T_ip;
450 variables.mechanical_strain
452 eps_m);
453
454 auto&& solution = _ip_data[ip].solid_material.integrateStress(
455 variables_prev, variables, t, x_position, dt, *state);
456
457 if (!solution)
458 {
459 OGS_FATAL("Computation of local constitutive relation failed.");
460 }
461
463 std::tie(sigma, state, C) = std::move(*solution);
464
465 local_Jac.noalias() += B.transpose() * C * B * w;
466
467 typename ShapeMatricesType::template MatrixType<DisplacementDim,
468 displacement_size>
469 N_u = ShapeMatricesType::template MatrixType<
470 DisplacementDim, displacement_size>::Zero(DisplacementDim,
471 displacement_size);
472
473 for (int i = 0; i < DisplacementDim; ++i)
474 {
475 N_u.template block<1, displacement_size / DisplacementDim>(
476 i, i * displacement_size / DisplacementDim)
477 .noalias() = N;
478 }
479
480 auto const rho_s =
481 solid_phase.property(MPL::PropertyType::density)
482 .template value<double>(variables, x_position, t, dt);
483
484 auto const& b = _process_data.specific_body_force;
485 local_rhs.noalias() -=
486 (B.transpose() * sigma - N_u.transpose() * rho_s * b) * w;
487 }
488}
489
490template <typename ShapeFunction, int DisplacementDim>
493 const double t, double const dt, Eigen::VectorXd const& local_x,
494 Eigen::VectorXd const& local_x_prev, std::vector<double>& local_b_data,
495 std::vector<double>& local_Jac_data)
496{
497 auto const local_T =
498 local_x.template segment<temperature_size>(temperature_index);
499
500 auto const local_T_prev =
501 local_x_prev.template segment<temperature_size>(temperature_index);
502
503 auto local_Jac = MathLib::createZeroedMatrix<
504 typename ShapeMatricesType::template MatrixType<temperature_size,
505 temperature_size>>(
506 local_Jac_data, temperature_size, temperature_size);
507
508 auto local_rhs = MathLib::createZeroedVector<
509 typename ShapeMatricesType::template VectorType<temperature_size>>(
510 local_b_data, temperature_size);
511
513 mass.setZero(temperature_size, temperature_size);
514
515 typename ShapeMatricesType::NodalMatrixType laplace;
516 laplace.setZero(temperature_size, temperature_size);
517
518 unsigned const n_integration_points =
519 _integration_method.getNumberOfPoints();
520
522 x_position.setElementID(_element.getID());
523 auto const& medium = _process_data.media_map.getMedium(_element.getID());
524 auto const& solid_phase = medium->phase("Solid");
525 MPL::VariableArray variables;
526
527 for (unsigned ip = 0; ip < n_integration_points; ip++)
528 {
529 x_position.setIntegrationPoint(ip);
530 auto const& w = _ip_data[ip].integration_weight;
531 auto const& N = _ip_data[ip].N;
532 auto const& dNdx = _ip_data[ip].dNdx;
533
534 const double T_ip = N.dot(local_T); // T at integration point
535 variables.temperature = T_ip;
536
537 auto const rho_s =
538 solid_phase.property(MPL::PropertyType::density)
539 .template value<double>(variables, x_position, t, dt);
540 auto const c_p =
541 solid_phase.property(MPL::PropertyType::specific_heat_capacity)
542 .template value<double>(variables, x_position, t, dt);
543
544 mass.noalias() += N.transpose() * rho_s * c_p * N * w;
545
546 auto const lambda =
547 solid_phase
548 .property(
550 .value(variables, x_position, t, dt);
551
552 GlobalDimMatrixType const thermal_conductivity =
553 MaterialPropertyLib::formEigenTensor<DisplacementDim>(lambda);
554
555 laplace.noalias() += dNdx.transpose() * thermal_conductivity * dNdx * w;
556 }
557 local_Jac.noalias() += laplace + mass / dt;
558
559 local_rhs.noalias() -=
560 laplace * local_T + mass * (local_T - local_T_prev) / dt;
561}
562
563template <typename ShapeFunction, int DisplacementDim>
564std::size_t
566 double const* values)
567{
568 return ProcessLib::setIntegrationPointKelvinVectorData<DisplacementDim>(
569 values, _ip_data, &IpData::sigma);
570}
571
572template <typename ShapeFunction, int DisplacementDim>
573std::vector<double>
575{
576 constexpr int kelvin_vector_size =
578
579 return transposeInPlace<kelvin_vector_size>(
580 [this](std::vector<double>& values)
581 { return getIntPtSigma(0, {}, {}, values); });
582}
583
584template <typename ShapeFunction, int DisplacementDim>
585std::vector<double> const&
587 const double /*t*/,
588 std::vector<GlobalVector*> const& /*x*/,
589 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
590 std::vector<double>& cache) const
591{
592 return ProcessLib::getIntegrationPointKelvinVectorData<DisplacementDim>(
593 _ip_data, &IpData::sigma, cache);
594}
595
596template <typename ShapeFunction, int DisplacementDim>
597std::size_t
599 double const* values)
600{
601 return ProcessLib::setIntegrationPointKelvinVectorData<DisplacementDim>(
602 values, _ip_data, &IpData::eps);
603}
604
605template <typename ShapeFunction, int DisplacementDim>
606std::vector<double> ThermoMechanicsLocalAssembler<
607 ShapeFunction, DisplacementDim>::getEpsilon() const
608{
609 constexpr int kelvin_vector_size =
611
612 return transposeInPlace<kelvin_vector_size>(
613 [this](std::vector<double>& values)
614 { return getIntPtEpsilon(0, {}, {}, values); });
615}
616
617template <typename ShapeFunction, int DisplacementDim>
618std::vector<double> const&
620 const double /*t*/,
621 std::vector<GlobalVector*> const& /*x*/,
622 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
623 std::vector<double>& cache) const
624{
625 return ProcessLib::getIntegrationPointKelvinVectorData<DisplacementDim>(
626 _ip_data, &IpData::eps, cache);
627}
628
629template <typename ShapeFunction, int DisplacementDim>
630std::vector<double> const&
633 const double /*t*/,
634 std::vector<GlobalVector*> const& /*x*/,
635 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
636 std::vector<double>& cache) const
637{
638 return ProcessLib::getIntegrationPointKelvinVectorData<DisplacementDim>(
639 _ip_data, &IpData::eps_m, cache);
640}
641
642template <typename ShapeFunction, int DisplacementDim>
644 ShapeFunction, DisplacementDim>::setEpsilonMechanical(double const* values)
645{
646 return ProcessLib::setIntegrationPointKelvinVectorData<DisplacementDim>(
647 values, _ip_data, &IpData::eps_m);
648}
649template <typename ShapeFunction, int DisplacementDim>
650std::vector<double> ThermoMechanicsLocalAssembler<
651 ShapeFunction, DisplacementDim>::getEpsilonMechanical() const
652{
653 constexpr int kelvin_vector_size =
655
656 return transposeInPlace<kelvin_vector_size>(
657 [this](std::vector<double>& values)
658 { return getIntPtEpsilonMechanical(0, {}, {}, values); });
659}
660
661template <typename ShapeFunction, int DisplacementDim>
663 ShapeFunction, DisplacementDim>::getNumberOfIntegrationPoints() const
664{
665 return _integration_method.getNumberOfPoints();
666}
667
668template <typename ShapeFunction, int DisplacementDim>
669int ThermoMechanicsLocalAssembler<ShapeFunction,
670 DisplacementDim>::getMaterialID() const
671{
672 return _process_data.material_ids == nullptr
673 ? 0
674 : (*_process_data.material_ids)[_element.getID()];
675}
676
677template <typename ShapeFunction, int DisplacementDim>
679 DisplacementDim>::MaterialStateVariables const&
681 getMaterialStateVariablesAt(unsigned integration_point) const
682{
683 return *_ip_data[integration_point].material_state_variables;
684}
685} // namespace ThermoMechanics
686} // namespace ProcessLib
#define OGS_FATAL(...)
Definition Error.h:26
std::variant< std::monostate, Eigen::Matrix< double, 4, 1 >, Eigen::Matrix< double, 6, 1 > > mechanical_strain
std::variant< std::monostate, Eigen::Matrix< double, 4, 1 >, Eigen::Matrix< double, 6, 1 > > stress
double getWeight() const
std::size_t getID() const
Returns the ID of the element.
Definition Element.h:89
MathLib::WeightedPoint const & getWeightedPoint(unsigned const igp) const
void setElementID(std::size_t element_id)
void setIntegrationPoint(unsigned integration_point)
MatrixType< _kelvin_vector_size, _number_of_dof > BMatrixType
std::vector< double > const & getIntPtEpsilonMechanical(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
SecondaryData< typename ShapeMatrices::ShapeType > _secondary_data
std::size_t setIPDataInitialConditions(std::string_view const name, double const *values, int const integration_order) override
Returns number of read integration points.
std::vector< double > const & getIntPtSigma(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
ThermoMechanicsLocalAssembler(ThermoMechanicsLocalAssembler const &)=delete
MaterialLib::Solids::MechanicsBase< DisplacementDim >::MaterialStateVariables const & getMaterialStateVariablesAt(unsigned integration_point) const override
void assembleWithJacobianForStaggeredScheme(double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, int const process_id, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data) override
ShapeMatrixPolicyType< ShapeFunction, DisplacementDim > ShapeMatricesType
ThermoMechanicsProcessData< DisplacementDim > & _process_data
NumLib::GenericIntegrationMethod const & _integration_method
void assembleWithJacobianForHeatConductionEquations(const double t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data)
void assembleWithJacobianForDeformationEquations(const double t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data)
typename ShapeMatricesType::GlobalDimMatrixType GlobalDimMatrixType
void assembleWithJacobian(double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_x_prev, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_rhs_data, std::vector< double > &local_Jac_data) override
std::vector< IpData, Eigen::aligned_allocator< IpData > > _ip_data
std::vector< double > const & getIntPtEpsilon(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
auto & selectSolidConstitutiveRelation(SolidMaterialsMap const &constitutive_relations, MeshLib::PropertyVector< int > const *const material_ids, std::size_t const element_id)
Eigen::Matrix< double, 4, 1 > kelvinVectorToSymmetricTensor(Eigen::Matrix< double, 4, 1, Eigen::ColMajor, 4, 1 > const &v)
constexpr int kelvin_vector_dimensions(int const displacement_dim)
Kelvin vector dimensions for given displacement dimension.
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), 1, Eigen::ColMajor > KelvinVectorType
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), kelvin_vector_dimensions(DisplacementDim), Eigen::RowMajor > KelvinMatrixType
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)
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)
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
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N