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 auto const x_coord =
174 _element, N);
175 auto const& B =
176 LinearBMatrix::computeBMatrix<DisplacementDim,
177 ShapeFunction::NPOINTS,
179 dNdx, N, x_coord, _is_axially_symmetric);
180
181 auto& sigma = _ip_data[ip].sigma;
182 auto const& sigma_prev = _ip_data[ip].sigma_prev;
183
184 auto& eps = _ip_data[ip].eps;
185 auto const& eps_prev = _ip_data[ip].eps_prev;
186
187 auto& eps_m = _ip_data[ip].eps_m;
188 auto const& eps_m_prev = _ip_data[ip].eps_m_prev;
189
190 auto& state = _ip_data[ip].material_state_variables;
191
192 const double T_ip = N.dot(T); // T at integration point
193 double const T_prev_ip = N.dot(T_prev);
194
195 // Consider also anisotropic thermal expansion.
196 auto const solid_linear_thermal_expansivity_vector =
198 solid_phase
199 .property(
201 .value(variables, x_position, t, dt));
202
204 dthermal_strain =
205 solid_linear_thermal_expansivity_vector * (T_ip - T_prev_ip);
206
207 //
208 // displacement equation, displacement part
209 //
210 // For the restart computation, the displacement may not be
211 // reloaded but the initial strains are always available. For such case,
212 // the following computation is skipped.
213 if (is_u_non_zero)
214 {
215 eps.noalias() = B * u;
216 }
217
218 eps_m.noalias() = eps_m_prev + eps - eps_prev - dthermal_strain;
219
220 variables_prev.stress
222 sigma_prev);
223 variables_prev.mechanical_strain
225 eps_m_prev);
226 variables_prev.temperature = T_ip;
227 variables.mechanical_strain
229 eps_m);
230 variables.temperature = T_ip;
231
232 auto&& solution = _ip_data[ip].solid_material.integrateStress(
233 variables_prev, variables, t, x_position, dt, *state);
234
235 if (!solution)
236 {
237 OGS_FATAL("Computation of local constitutive relation failed.");
238 }
239
241 std::tie(sigma, state, C) = std::move(*solution);
242
243 local_Jac
244 .template block<displacement_size, displacement_size>(
245 displacement_index, displacement_index)
246 .noalias() += B.transpose() * C * B * w;
247
248 typename ShapeMatricesType::template MatrixType<DisplacementDim,
249 displacement_size>
250 N_u = ShapeMatricesType::template MatrixType<
251 DisplacementDim, displacement_size>::Zero(DisplacementDim,
252 displacement_size);
253
254 for (int i = 0; i < DisplacementDim; ++i)
255 {
256 N_u.template block<1, displacement_size / DisplacementDim>(
257 i, i * displacement_size / DisplacementDim)
258 .noalias() = N;
259 }
260
261 auto const rho_s =
262 solid_phase.property(MPL::PropertyType::density)
263 .template value<double>(variables, x_position, t, dt);
264
265 auto const& b = _process_data.specific_body_force;
266 local_rhs.template block<displacement_size, 1>(displacement_index, 0)
267 .noalias() -=
268 (B.transpose() * sigma - N_u.transpose() * rho_s * b) * w;
269
270 //
271 // displacement equation, temperature part
272 // The computation of KuT can be ignored.
273 auto const alpha_T_tensor =
275 solid_linear_thermal_expansivity_vector);
276 KuT.noalias() += B.transpose() * (C * alpha_T_tensor) * N * w;
277
278 if (_ip_data[ip].solid_material.getConstitutiveModel() ==
280 {
281 auto const s = Invariants::deviatoric_projection * sigma;
282 double const norm_s = Invariants::FrobeniusNorm(s);
283 const double creep_coefficient =
284 _ip_data[ip].solid_material.getTemperatureRelatedCoefficient(
285 t, dt, x_position, T_ip, norm_s);
286 KuT.noalias() += creep_coefficient * B.transpose() * s * N * w;
287 }
288
289 //
290 // temperature equation, temperature part;
291 //
292 auto const lambda =
293 solid_phase
294 .property(
296 .value(variables, x_position, t, dt);
297
298 GlobalDimMatrixType const thermal_conductivity =
300
301 KTT.noalias() += dNdx.transpose() * thermal_conductivity * dNdx * w;
302
303 auto const c =
304 solid_phase
305 .property(
307 .template value<double>(variables, x_position, t, dt);
308 DTT.noalias() += N.transpose() * rho_s * c * N * w;
309 }
310
311 // temperature equation, temperature part
312 local_Jac
313 .template block<temperature_size, temperature_size>(temperature_index,
314 temperature_index)
315 .noalias() += KTT + DTT / dt;
316
317 // displacement equation, temperature part
318 local_Jac
319 .template block<displacement_size, temperature_size>(displacement_index,
320 temperature_index)
321 .noalias() -= KuT;
322
323 local_rhs.template block<temperature_size, 1>(temperature_index, 0)
324 .noalias() -= KTT * T + DTT * (T - T_prev) / dt;
325}
326
327template <typename ShapeFunction, int DisplacementDim>
329 assembleWithJacobianForStaggeredScheme(const double t, double const dt,
330 Eigen::VectorXd const& local_x,
331 Eigen::VectorXd const& local_x_prev,
332 int const process_id,
333 std::vector<double>& local_b_data,
334 std::vector<double>& local_Jac_data)
335{
336 // For the equations with pressure
337 if (process_id == _process_data.heat_conduction_process_id)
338 {
339 assembleWithJacobianForHeatConductionEquations(
340 t, dt, local_x, local_x_prev, local_b_data, local_Jac_data);
341 return;
342 }
343
344 // For the equations with deformation
345 assembleWithJacobianForDeformationEquations(t, dt, local_x, local_x_prev,
346 local_b_data, local_Jac_data);
347}
348
349template <typename ShapeFunction, int DisplacementDim>
352 const double t, double const dt, Eigen::VectorXd const& local_x,
353 Eigen::VectorXd const& local_x_prev, std::vector<double>& local_b_data,
354 std::vector<double>& local_Jac_data)
355{
356 auto const local_T =
357 local_x.template segment<temperature_size>(temperature_index);
358
359 auto const local_T_prev =
360 local_x_prev.template segment<temperature_size>(temperature_index);
361
362 auto const local_u =
363 local_x.template segment<displacement_size>(displacement_index);
364 bool const is_u_non_zero = local_u.norm() > 0.0;
365
366 auto local_Jac = MathLib::createZeroedMatrix<
367 typename ShapeMatricesType::template MatrixType<displacement_size,
368 displacement_size>>(
369 local_Jac_data, displacement_size, displacement_size);
370
371 auto local_rhs = MathLib::createZeroedVector<
372 typename ShapeMatricesType::template VectorType<displacement_size>>(
373 local_b_data, displacement_size);
374
375 unsigned const n_integration_points =
376 _integration_method.getNumberOfPoints();
377
378 MPL::VariableArray variables;
379 MPL::VariableArray variables_prev;
381 x_position.setElementID(_element.getID());
382 auto const& medium = _process_data.media_map.getMedium(_element.getID());
383 auto const& solid_phase = medium->phase("Solid");
384
385 for (unsigned ip = 0; ip < n_integration_points; ip++)
386 {
387 auto const& w = _ip_data[ip].integration_weight;
388 auto const& N = _ip_data[ip].N;
389 auto const& dNdx = _ip_data[ip].dNdx;
390
391 auto const x_coord =
393 _element, N);
394 auto const& B =
395 LinearBMatrix::computeBMatrix<DisplacementDim,
396 ShapeFunction::NPOINTS,
398 dNdx, N, x_coord, _is_axially_symmetric);
399
400 auto& sigma = _ip_data[ip].sigma;
401 auto const& sigma_prev = _ip_data[ip].sigma_prev;
402
403 auto& eps = _ip_data[ip].eps;
404 auto const& eps_prev = _ip_data[ip].eps_prev;
405
406 auto& eps_m = _ip_data[ip].eps_m;
407 auto const& eps_m_prev = _ip_data[ip].eps_m_prev;
408
409 auto& state = _ip_data[ip].material_state_variables;
410
411 const double T_ip = N.dot(local_T); // T at integration point
412 variables.temperature = T_ip;
413 double const dT_ip = T_ip - N.dot(local_T_prev);
414
415 //
416 // displacement equation, displacement part
417 //
418 // For the restart computation, the displacement may not be
419 // reloaded but the initial strains are always available. For such case,
420 // the following computation is skipped.
421 if (is_u_non_zero)
422 {
423 eps.noalias() = B * local_u;
424 }
425
426 // Consider also anisotropic thermal expansion.
427 auto const solid_linear_thermal_expansivity_vector =
429 solid_phase
430 .property(
432 .value(variables, x_position, t, dt));
433
435 dthermal_strain = solid_linear_thermal_expansivity_vector * dT_ip;
436
437 eps_m.noalias() = eps_m_prev + eps - eps_prev - dthermal_strain;
438
439 variables_prev.stress
441 sigma_prev);
442 variables_prev.mechanical_strain
444 eps_m_prev);
445 variables_prev.temperature = T_ip;
446 variables.mechanical_strain
448 eps_m);
449
450 auto&& solution = _ip_data[ip].solid_material.integrateStress(
451 variables_prev, variables, t, x_position, dt, *state);
452
453 if (!solution)
454 {
455 OGS_FATAL("Computation of local constitutive relation failed.");
456 }
457
459 std::tie(sigma, state, C) = std::move(*solution);
460
461 local_Jac.noalias() += B.transpose() * C * B * w;
462
463 typename ShapeMatricesType::template MatrixType<DisplacementDim,
464 displacement_size>
465 N_u = ShapeMatricesType::template MatrixType<
466 DisplacementDim, displacement_size>::Zero(DisplacementDim,
467 displacement_size);
468
469 for (int i = 0; i < DisplacementDim; ++i)
470 {
471 N_u.template block<1, displacement_size / DisplacementDim>(
472 i, i * displacement_size / DisplacementDim)
473 .noalias() = N;
474 }
475
476 auto const rho_s =
477 solid_phase.property(MPL::PropertyType::density)
478 .template value<double>(variables, x_position, t, dt);
479
480 auto const& b = _process_data.specific_body_force;
481 local_rhs.noalias() -=
482 (B.transpose() * sigma - N_u.transpose() * rho_s * b) * w;
483 }
484}
485
486template <typename ShapeFunction, int DisplacementDim>
489 const double t, double const dt, Eigen::VectorXd const& local_x,
490 Eigen::VectorXd const& local_x_prev, std::vector<double>& local_b_data,
491 std::vector<double>& local_Jac_data)
492{
493 auto const local_T =
494 local_x.template segment<temperature_size>(temperature_index);
495
496 auto const local_T_prev =
497 local_x_prev.template segment<temperature_size>(temperature_index);
498
499 auto local_Jac = MathLib::createZeroedMatrix<
500 typename ShapeMatricesType::template MatrixType<temperature_size,
501 temperature_size>>(
502 local_Jac_data, temperature_size, temperature_size);
503
504 auto local_rhs = MathLib::createZeroedVector<
505 typename ShapeMatricesType::template VectorType<temperature_size>>(
506 local_b_data, temperature_size);
507
509 mass.setZero(temperature_size, temperature_size);
510
511 typename ShapeMatricesType::NodalMatrixType laplace;
512 laplace.setZero(temperature_size, temperature_size);
513
514 unsigned const n_integration_points =
515 _integration_method.getNumberOfPoints();
516
518 x_position.setElementID(_element.getID());
519 auto const& medium = _process_data.media_map.getMedium(_element.getID());
520 auto const& solid_phase = medium->phase("Solid");
521 MPL::VariableArray variables;
522
523 for (unsigned ip = 0; ip < n_integration_points; ip++)
524 {
525 auto const& w = _ip_data[ip].integration_weight;
526 auto const& N = _ip_data[ip].N;
527 auto const& dNdx = _ip_data[ip].dNdx;
528
529 const double T_ip = N.dot(local_T); // T at integration point
530 variables.temperature = T_ip;
531
532 auto const rho_s =
533 solid_phase.property(MPL::PropertyType::density)
534 .template value<double>(variables, x_position, t, dt);
535 auto const c_p =
536 solid_phase.property(MPL::PropertyType::specific_heat_capacity)
537 .template value<double>(variables, x_position, t, dt);
538
539 mass.noalias() += N.transpose() * rho_s * c_p * N * w;
540
541 auto const lambda =
542 solid_phase
543 .property(
545 .value(variables, x_position, t, dt);
546
547 GlobalDimMatrixType const thermal_conductivity =
549
550 laplace.noalias() += dNdx.transpose() * thermal_conductivity * dNdx * w;
551 }
552 local_Jac.noalias() += laplace + mass / dt;
553
554 local_rhs.noalias() -=
555 laplace * local_T + mass * (local_T - local_T_prev) / dt;
556}
557
558template <typename ShapeFunction, int DisplacementDim>
559std::size_t
566
567template <typename ShapeFunction, int DisplacementDim>
568std::vector<double>
570{
571 constexpr int kelvin_vector_size =
573
575 [this](std::vector<double>& values)
576 { return getIntPtSigma(0, {}, {}, values); });
577}
578
579template <typename ShapeFunction, int DisplacementDim>
580std::vector<double> const&
582 const double /*t*/,
583 std::vector<GlobalVector*> const& /*x*/,
584 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
585 std::vector<double>& cache) const
586{
588 _ip_data, &IpData::sigma, cache);
589}
590
591template <typename ShapeFunction, int DisplacementDim>
592std::size_t
599
600template <typename ShapeFunction, int DisplacementDim>
601std::vector<double> ThermoMechanicsLocalAssembler<
602 ShapeFunction, DisplacementDim>::getEpsilon() const
603{
604 constexpr int kelvin_vector_size =
606
608 [this](std::vector<double>& values)
609 { return getIntPtEpsilon(0, {}, {}, values); });
610}
611
612template <typename ShapeFunction, int DisplacementDim>
613std::vector<double> const&
615 const double /*t*/,
616 std::vector<GlobalVector*> const& /*x*/,
617 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
618 std::vector<double>& cache) const
619{
621 _ip_data, &IpData::eps, cache);
622}
623
624template <typename ShapeFunction, int DisplacementDim>
625std::vector<double> const&
628 const double /*t*/,
629 std::vector<GlobalVector*> const& /*x*/,
630 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
631 std::vector<double>& cache) const
632{
634 _ip_data, &IpData::eps_m, cache);
635}
636
637template <typename ShapeFunction, int DisplacementDim>
639 ShapeFunction, DisplacementDim>::setEpsilonMechanical(double const* values)
640{
642 values, _ip_data, &IpData::eps_m);
643}
644template <typename ShapeFunction, int DisplacementDim>
645std::vector<double> ThermoMechanicsLocalAssembler<
646 ShapeFunction, DisplacementDim>::getEpsilonMechanical() const
647{
648 constexpr int kelvin_vector_size =
650
652 [this](std::vector<double>& values)
653 { return getIntPtEpsilonMechanical(0, {}, {}, values); });
654}
655
656template <typename ShapeFunction, int DisplacementDim>
658 ShapeFunction, DisplacementDim>::getNumberOfIntegrationPoints() const
659{
660 return _integration_method.getNumberOfPoints();
661}
662
663template <typename ShapeFunction, int DisplacementDim>
664int ThermoMechanicsLocalAssembler<ShapeFunction,
665 DisplacementDim>::getMaterialID() const
666{
667 return _process_data.material_ids == nullptr
668 ? 0
669 : (*_process_data.material_ids)[_element.getID()];
670}
671
672template <typename ShapeFunction, int DisplacementDim>
674 DisplacementDim>::MaterialStateVariables const&
676 getMaterialStateVariablesAt(unsigned integration_point) const
677{
678 return *_ip_data[integration_point].material_state_variables;
679}
680} // namespace ThermoMechanics
681} // 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)
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)
double interpolateXCoordinate(MeshLib::Element const &e, typename ShapeMatricesType::ShapeMatrices::ShapeType const &N)
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.
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