Loading [MathJax]/extensions/tex2jax.js
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_rhs_data,
120 std::vector<double>& local_Jac_data)
121{
122 auto const local_matrix_size = local_x.size();
123 assert(local_matrix_size == temperature_size + displacement_size);
124
125 auto T = Eigen::Map<typename ShapeMatricesType::template VectorType<
126 temperature_size> const>(local_x.data() + temperature_index,
127 temperature_size);
128
129 auto u = Eigen::Map<typename ShapeMatricesType::template VectorType<
130 displacement_size> const>(local_x.data() + displacement_index,
131 displacement_size);
132 bool const is_u_non_zero = u.norm() > 0.0;
133
134 auto T_prev = Eigen::Map<typename ShapeMatricesType::template VectorType<
135 temperature_size> const>(local_x_prev.data() + temperature_index,
136 temperature_size);
137
139 local_Jac_data, local_matrix_size, local_matrix_size);
140
141 auto local_rhs = MathLib::createZeroedVector<RhsVector>(local_rhs_data,
142 local_matrix_size);
143
144 typename ShapeMatricesType::template MatrixType<displacement_size,
145 temperature_size>
146 KuT;
147 KuT.setZero(displacement_size, temperature_size);
148
150 KTT.setZero(temperature_size, temperature_size);
151
153 DTT.setZero(temperature_size, temperature_size);
154
155 unsigned const n_integration_points =
156 _integration_method.getNumberOfPoints();
157
158 MPL::VariableArray variables;
159 MPL::VariableArray variables_prev;
161 x_position.setElementID(_element.getID());
162
163 auto const& medium = _process_data.media_map.getMedium(_element.getID());
164 auto const& solid_phase = medium->phase("Solid");
165
166 for (unsigned ip = 0; ip < n_integration_points; ip++)
167 {
168 auto const& w = _ip_data[ip].integration_weight;
169 auto const& N = _ip_data[ip].N;
170 auto const& dNdx = _ip_data[ip].dNdx;
171
172 ParameterLib::SpatialPosition const x_position{
173 std::nullopt, _element.getID(),
176 _element, N))};
177
178 auto const x_coord = x_position.getCoordinates().value()[0];
179 auto const& B =
180 LinearBMatrix::computeBMatrix<DisplacementDim,
181 ShapeFunction::NPOINTS,
183 dNdx, N, x_coord, _is_axially_symmetric);
184
185 auto& sigma = _ip_data[ip].sigma;
186 auto const& sigma_prev = _ip_data[ip].sigma_prev;
187
188 auto& eps = _ip_data[ip].eps;
189 auto const& eps_prev = _ip_data[ip].eps_prev;
190
191 auto& eps_m = _ip_data[ip].eps_m;
192 auto const& eps_m_prev = _ip_data[ip].eps_m_prev;
193
194 auto& state = _ip_data[ip].material_state_variables;
195
196 const double T_ip = N.dot(T); // T at integration point
197 double const T_prev_ip = N.dot(T_prev);
198
199 // Consider also anisotropic thermal expansion.
200 auto const solid_linear_thermal_expansivity_vector =
202 solid_phase
203 .property(
205 .value(variables, x_position, t, dt));
206
208 dthermal_strain =
209 solid_linear_thermal_expansivity_vector * (T_ip - T_prev_ip);
210
211 //
212 // displacement equation, displacement part
213 //
214 // For the restart computation, the displacement may not be
215 // reloaded but the initial strains are always available. For such case,
216 // the following computation is skipped.
217 if (is_u_non_zero)
218 {
219 eps.noalias() = B * u;
220 }
221
222 eps_m.noalias() = eps_m_prev + eps - eps_prev - dthermal_strain;
223
224 variables_prev.stress
226 sigma_prev);
227 variables_prev.mechanical_strain
229 eps_m_prev);
230 variables_prev.temperature = T_ip;
231 variables.mechanical_strain
233 eps_m);
234 variables.temperature = T_ip;
235
236 auto&& solution = _ip_data[ip].solid_material.integrateStress(
237 variables_prev, variables, t, x_position, dt, *state);
238
239 if (!solution)
240 {
241 OGS_FATAL("Computation of local constitutive relation failed.");
242 }
243
245 std::tie(sigma, state, C) = std::move(*solution);
246
247 local_Jac
248 .template block<displacement_size, displacement_size>(
249 displacement_index, displacement_index)
250 .noalias() += B.transpose() * C * B * w;
251
252 typename ShapeMatricesType::template MatrixType<DisplacementDim,
253 displacement_size>
254 N_u = ShapeMatricesType::template MatrixType<
255 DisplacementDim, displacement_size>::Zero(DisplacementDim,
256 displacement_size);
257
258 for (int i = 0; i < DisplacementDim; ++i)
259 {
260 N_u.template block<1, displacement_size / DisplacementDim>(
261 i, i * displacement_size / DisplacementDim)
262 .noalias() = N;
263 }
264
265 auto const rho_s =
266 solid_phase.property(MPL::PropertyType::density)
267 .template value<double>(variables, x_position, t, dt);
268
269 auto const& b = _process_data.specific_body_force;
270 local_rhs.template block<displacement_size, 1>(displacement_index, 0)
271 .noalias() -=
272 (B.transpose() * sigma - N_u.transpose() * rho_s * b) * w;
273
274 //
275 // displacement equation, temperature part
276 // The computation of KuT can be ignored.
277 auto const alpha_T_tensor =
279 solid_linear_thermal_expansivity_vector);
280 KuT.noalias() += B.transpose() * (C * alpha_T_tensor) * N * w;
281
282 if (_ip_data[ip].solid_material.getConstitutiveModel() ==
284 {
285 auto const s = Invariants::deviatoric_projection * sigma;
286 double const norm_s = Invariants::FrobeniusNorm(s);
287 const double creep_coefficient =
288 _ip_data[ip].solid_material.getTemperatureRelatedCoefficient(
289 t, dt, x_position, T_ip, norm_s);
290 KuT.noalias() += creep_coefficient * B.transpose() * s * N * w;
291 }
292
293 //
294 // temperature equation, temperature part;
295 //
296 auto const lambda =
297 solid_phase
298 .property(
300 .value(variables, x_position, t, dt);
301
302 GlobalDimMatrixType const thermal_conductivity =
304
305 KTT.noalias() += dNdx.transpose() * thermal_conductivity * dNdx * w;
306
307 auto const c =
308 solid_phase
309 .property(
311 .template value<double>(variables, x_position, t, dt);
312 DTT.noalias() += N.transpose() * rho_s * c * N * w;
313 }
314
315 // temperature equation, temperature part
316 local_Jac
317 .template block<temperature_size, temperature_size>(temperature_index,
318 temperature_index)
319 .noalias() += KTT + DTT / dt;
320
321 // displacement equation, temperature part
322 local_Jac
323 .template block<displacement_size, temperature_size>(displacement_index,
324 temperature_index)
325 .noalias() -= KuT;
326
327 local_rhs.template block<temperature_size, 1>(temperature_index, 0)
328 .noalias() -= KTT * T + DTT * (T - T_prev) / dt;
329}
330
331template <typename ShapeFunction, int DisplacementDim>
333 assembleWithJacobianForStaggeredScheme(const double t, double const dt,
334 Eigen::VectorXd const& local_x,
335 Eigen::VectorXd const& local_x_prev,
336 int const process_id,
337 std::vector<double>& local_b_data,
338 std::vector<double>& local_Jac_data)
339{
340 // For the equations with pressure
341 if (process_id == _process_data.heat_conduction_process_id)
342 {
343 assembleWithJacobianForHeatConductionEquations(
344 t, dt, local_x, local_x_prev, local_b_data, local_Jac_data);
345 return;
346 }
347
348 // For the equations with deformation
349 assembleWithJacobianForDeformationEquations(t, dt, local_x, local_x_prev,
350 local_b_data, local_Jac_data);
351}
352
353template <typename ShapeFunction, int DisplacementDim>
356 const double t, double const dt, Eigen::VectorXd const& local_x,
357 Eigen::VectorXd const& local_x_prev, std::vector<double>& local_b_data,
358 std::vector<double>& local_Jac_data)
359{
360 auto const local_T =
361 local_x.template segment<temperature_size>(temperature_index);
362
363 auto const local_T_prev =
364 local_x_prev.template segment<temperature_size>(temperature_index);
365
366 auto const local_u =
367 local_x.template segment<displacement_size>(displacement_index);
368 bool const is_u_non_zero = local_u.norm() > 0.0;
369
370 auto local_Jac = MathLib::createZeroedMatrix<
371 typename ShapeMatricesType::template MatrixType<displacement_size,
372 displacement_size>>(
373 local_Jac_data, displacement_size, displacement_size);
374
375 auto local_rhs = MathLib::createZeroedVector<
376 typename ShapeMatricesType::template VectorType<displacement_size>>(
377 local_b_data, displacement_size);
378
379 unsigned const n_integration_points =
380 _integration_method.getNumberOfPoints();
381
382 MPL::VariableArray variables;
383 MPL::VariableArray variables_prev;
384 auto const& medium = _process_data.media_map.getMedium(_element.getID());
385 auto const& solid_phase = medium->phase("Solid");
386
387 for (unsigned ip = 0; ip < n_integration_points; ip++)
388 {
389 auto const& w = _ip_data[ip].integration_weight;
390 auto const& N = _ip_data[ip].N;
391 auto const& dNdx = _ip_data[ip].dNdx;
392
393 ParameterLib::SpatialPosition const x_position{
394 std::nullopt, _element.getID(),
397 _element, N))};
398
399 auto const x_coord = x_position.getCoordinates().value()[0];
400
401 auto const& B =
402 LinearBMatrix::computeBMatrix<DisplacementDim,
403 ShapeFunction::NPOINTS,
405 dNdx, N, x_coord, _is_axially_symmetric);
406
407 auto& sigma = _ip_data[ip].sigma;
408 auto const& sigma_prev = _ip_data[ip].sigma_prev;
409
410 auto& eps = _ip_data[ip].eps;
411 auto const& eps_prev = _ip_data[ip].eps_prev;
412
413 auto& eps_m = _ip_data[ip].eps_m;
414 auto const& eps_m_prev = _ip_data[ip].eps_m_prev;
415
416 auto& state = _ip_data[ip].material_state_variables;
417
418 const double T_ip = N.dot(local_T); // T at integration point
419 variables.temperature = T_ip;
420 double const dT_ip = T_ip - N.dot(local_T_prev);
421
422 //
423 // displacement equation, displacement part
424 //
425 // For the restart computation, the displacement may not be
426 // reloaded but the initial strains are always available. For such case,
427 // the following computation is skipped.
428 if (is_u_non_zero)
429 {
430 eps.noalias() = B * local_u;
431 }
432
433 // Consider also anisotropic thermal expansion.
434 auto const solid_linear_thermal_expansivity_vector =
436 solid_phase
437 .property(
439 .value(variables, x_position, t, dt));
440
442 dthermal_strain = solid_linear_thermal_expansivity_vector * dT_ip;
443
444 eps_m.noalias() = eps_m_prev + eps - eps_prev - dthermal_strain;
445
446 variables_prev.stress
448 sigma_prev);
449 variables_prev.mechanical_strain
451 eps_m_prev);
452 variables_prev.temperature = T_ip;
453 variables.mechanical_strain
455 eps_m);
456
457 auto&& solution = _ip_data[ip].solid_material.integrateStress(
458 variables_prev, variables, t, x_position, dt, *state);
459
460 if (!solution)
461 {
462 OGS_FATAL("Computation of local constitutive relation failed.");
463 }
464
466 std::tie(sigma, state, C) = std::move(*solution);
467
468 local_Jac.noalias() += B.transpose() * C * B * w;
469
470 typename ShapeMatricesType::template MatrixType<DisplacementDim,
471 displacement_size>
472 N_u = ShapeMatricesType::template MatrixType<
473 DisplacementDim, displacement_size>::Zero(DisplacementDim,
474 displacement_size);
475
476 for (int i = 0; i < DisplacementDim; ++i)
477 {
478 N_u.template block<1, displacement_size / DisplacementDim>(
479 i, i * displacement_size / DisplacementDim)
480 .noalias() = N;
481 }
482
483 auto const rho_s =
484 solid_phase.property(MPL::PropertyType::density)
485 .template value<double>(variables, x_position, t, dt);
486
487 auto const& b = _process_data.specific_body_force;
488 local_rhs.noalias() -=
489 (B.transpose() * sigma - N_u.transpose() * rho_s * b) * w;
490 }
491}
492
493template <typename ShapeFunction, int DisplacementDim>
496 const double t, double const dt, Eigen::VectorXd const& local_x,
497 Eigen::VectorXd const& local_x_prev, std::vector<double>& local_b_data,
498 std::vector<double>& local_Jac_data)
499{
500 auto const local_T =
501 local_x.template segment<temperature_size>(temperature_index);
502
503 auto const local_T_prev =
504 local_x_prev.template segment<temperature_size>(temperature_index);
505
506 auto local_Jac = MathLib::createZeroedMatrix<
507 typename ShapeMatricesType::template MatrixType<temperature_size,
508 temperature_size>>(
509 local_Jac_data, temperature_size, temperature_size);
510
511 auto local_rhs = MathLib::createZeroedVector<
512 typename ShapeMatricesType::template VectorType<temperature_size>>(
513 local_b_data, temperature_size);
514
516 mass.setZero(temperature_size, temperature_size);
517
518 typename ShapeMatricesType::NodalMatrixType laplace;
519 laplace.setZero(temperature_size, temperature_size);
520
521 unsigned const n_integration_points =
522 _integration_method.getNumberOfPoints();
523
524 auto const& medium = _process_data.media_map.getMedium(_element.getID());
525 auto const& solid_phase = medium->phase("Solid");
526 MPL::VariableArray variables;
527
528 for (unsigned ip = 0; ip < n_integration_points; ip++)
529 {
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 ParameterLib::SpatialPosition const x_position{
535 std::nullopt, _element.getID(),
538 _element, N))};
539
540 const double T_ip = N.dot(local_T); // T at integration point
541 variables.temperature = T_ip;
542
543 auto const rho_s =
544 solid_phase.property(MPL::PropertyType::density)
545 .template value<double>(variables, x_position, t, dt);
546 auto const c_p =
547 solid_phase.property(MPL::PropertyType::specific_heat_capacity)
548 .template value<double>(variables, x_position, t, dt);
549
550 mass.noalias() += N.transpose() * rho_s * c_p * N * w;
551
552 auto const lambda =
553 solid_phase
554 .property(
556 .value(variables, x_position, t, dt);
557
558 GlobalDimMatrixType const thermal_conductivity =
560
561 laplace.noalias() += dNdx.transpose() * thermal_conductivity * dNdx * w;
562 }
563 local_Jac.noalias() += laplace + mass / dt;
564
565 local_rhs.noalias() -=
566 laplace * local_T + mass * (local_T - local_T_prev) / dt;
567}
568
569template <typename ShapeFunction, int DisplacementDim>
570std::size_t
577
578template <typename ShapeFunction, int DisplacementDim>
579std::vector<double>
581{
582 constexpr int kelvin_vector_size =
584
586 [this](std::vector<double>& values)
587 { return getIntPtSigma(0, {}, {}, values); });
588}
589
590template <typename ShapeFunction, int DisplacementDim>
591std::vector<double> const&
593 const double /*t*/,
594 std::vector<GlobalVector*> const& /*x*/,
595 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
596 std::vector<double>& cache) const
597{
599 _ip_data, &IpData::sigma, cache);
600}
601
602template <typename ShapeFunction, int DisplacementDim>
603std::size_t
610
611template <typename ShapeFunction, int DisplacementDim>
612std::vector<double> ThermoMechanicsLocalAssembler<
613 ShapeFunction, DisplacementDim>::getEpsilon() const
614{
615 constexpr int kelvin_vector_size =
617
619 [this](std::vector<double>& values)
620 { return getIntPtEpsilon(0, {}, {}, values); });
621}
622
623template <typename ShapeFunction, int DisplacementDim>
624std::vector<double> const&
626 const double /*t*/,
627 std::vector<GlobalVector*> const& /*x*/,
628 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
629 std::vector<double>& cache) const
630{
632 _ip_data, &IpData::eps, cache);
633}
634
635template <typename ShapeFunction, int DisplacementDim>
636std::vector<double> const&
639 const double /*t*/,
640 std::vector<GlobalVector*> const& /*x*/,
641 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
642 std::vector<double>& cache) const
643{
645 _ip_data, &IpData::eps_m, cache);
646}
647
648template <typename ShapeFunction, int DisplacementDim>
650 ShapeFunction, DisplacementDim>::setEpsilonMechanical(double const* values)
651{
653 values, _ip_data, &IpData::eps_m);
654}
655template <typename ShapeFunction, int DisplacementDim>
656std::vector<double> ThermoMechanicsLocalAssembler<
657 ShapeFunction, DisplacementDim>::getEpsilonMechanical() const
658{
659 constexpr int kelvin_vector_size =
661
663 [this](std::vector<double>& values)
664 { return getIntPtEpsilonMechanical(0, {}, {}, values); });
665}
666
667template <typename ShapeFunction, int DisplacementDim>
669 ShapeFunction, DisplacementDim>::getNumberOfIntegrationPoints() const
670{
671 return _integration_method.getNumberOfPoints();
672}
673
674template <typename ShapeFunction, int DisplacementDim>
675int ThermoMechanicsLocalAssembler<ShapeFunction,
676 DisplacementDim>::getMaterialID() const
677{
678 return _process_data.material_ids == nullptr
679 ? 0
680 : (*_process_data.material_ids)[_element.getID()];
681}
682
683template <typename ShapeFunction, int DisplacementDim>
685 DisplacementDim>::MaterialStateVariables const&
687 getMaterialStateVariablesAt(unsigned integration_point) const
688{
689 return *_ip_data[integration_point].material_state_variables;
690}
691} // namespace ThermoMechanics
692} // namespace ProcessLib
#define OGS_FATAL(...)
Definition Error.h:26
constexpr 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)
std::optional< MathLib::Point3d > const getCoordinates() const
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
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)
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_b_data, std::vector< double > &local_Jac_data) override
typename ShapeMatricesType::GlobalDimMatrixType GlobalDimMatrixType
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
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_rhs_data, std::vector< double > &local_Jac_data) override
auto & selectSolidConstitutiveRelation(SolidMaterialsMap const &constitutive_relations, MeshLib::PropertyVector< int > const *const material_ids, std::size_t const element_id)
MathLib::KelvinVector::KelvinVectorType< GlobalDim > formKelvinVector(MaterialPropertyLib::PropertyDataType const &values)
A function to form a Kelvin vector from strain or stress alike property like thermal expansivity for ...
Eigen::Matrix< double, GlobalDim, GlobalDim > formEigenTensor(MaterialPropertyLib::PropertyDataType const &values)
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)
std::array< double, 3 > interpolateCoordinates(MeshLib::Element const &e, typename ShapeMatricesType::ShapeMatrices::ShapeType const &N)
BMatrixType computeBMatrix(DNDX_Type const &dNdx, N_Type const &N, const double radius, const bool is_axially_symmetric)
Fills a B-matrix based on given shape function dN/dx values.
std::vector< double > transposeInPlace(StoreValuesFunction const &store_values_function)
std::vector< double > const & getIntegrationPointKelvinVectorData(IntegrationPointDataVector const &ip_data_vector, MemberType IpData::*const member, std::vector< double > &cache)
std::size_t setIntegrationPointKelvinVectorData(double const *values, IntegrationPointDataVector &ip_data_vector, MemberType IpData::*const member)
MatrixType< ShapeFunction::NPOINTS, ShapeFunction::NPOINTS > NodalMatrixType
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N