OGS
ThermoMechanicsFEM-impl.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) OpenGeoSys Community (opengeosys.org)
2// SPDX-License-Identifier: BSD-3-Clause
3
4#pragma once
5
12
13namespace ProcessLib
14{
15namespace ThermoMechanics
16{
17template <typename ShapeFunction, int DisplacementDim>
20 MeshLib::Element const& e,
21 std::size_t const /*local_matrix_size*/,
22 NumLib::GenericIntegrationMethod const& integration_method,
23 bool const is_axially_symmetric,
25 : _process_data(process_data),
26 _integration_method(integration_method),
27 _element(e),
28 _is_axially_symmetric(is_axially_symmetric)
29{
30 unsigned const n_integration_points =
31 _integration_method.getNumberOfPoints();
32
33 _ip_data.reserve(n_integration_points);
34 _secondary_data.N.resize(n_integration_points);
35
36 auto const shape_matrices =
38 DisplacementDim>(e, is_axially_symmetric,
40
42 _process_data.solid_materials, _process_data.material_ids, e.getID());
43
44 for (unsigned ip = 0; ip < n_integration_points; ip++)
45 {
46 _ip_data.emplace_back(solid_material);
47 auto& ip_data = _ip_data[ip];
48 ip_data.integration_weight =
49 _integration_method.getWeightedPoint(ip).getWeight() *
50 shape_matrices[ip].integralMeasure * shape_matrices[ip].detJ;
51
52 static const int kelvin_vector_size =
54 // Initialize current time step values
55 ip_data.sigma.setZero(kelvin_vector_size);
56 ip_data.eps.setZero(kelvin_vector_size);
57
58 // Previous time step values are not initialized and are set later.
59 ip_data.sigma_prev.resize(kelvin_vector_size);
60 ip_data.eps_prev.resize(kelvin_vector_size);
61
62 ip_data.eps_m.setZero(kelvin_vector_size);
63 ip_data.eps_m_prev.setZero(kelvin_vector_size);
65 x_position.setElementID(_element.getID());
66 ip_data.N = shape_matrices[ip].N;
67 ip_data.dNdx = shape_matrices[ip].dNdx;
68
69 _secondary_data.N[ip] = shape_matrices[ip].N;
70 }
71}
72
73template <typename ShapeFunction, int DisplacementDim>
75 setInitialConditionsConcrete(Eigen::VectorXd const local_x,
76 double const /*t*/,
77 int const /*process_id*/)
78{
79 unsigned const n_integration_points =
80 this->_integration_method.getNumberOfPoints();
81 auto const local_u =
82 local_x.template segment<displacement_size>(displacement_index);
83
84 for (unsigned ip = 0; ip < n_integration_points; ip++)
85 {
86 auto const& N = _ip_data[ip].N;
87 auto const& dNdx = _ip_data[ip].dNdx;
88
89 ParameterLib::SpatialPosition const x_position{
90 std::nullopt, _element.getID(),
93 _element, N))};
94
95 auto const x_coord = x_position.getCoordinates().value()[0];
96
97 auto const& B =
98 LinearBMatrix::computeBMatrix<DisplacementDim,
99 ShapeFunction::NPOINTS,
101 dNdx, N, x_coord, _is_axially_symmetric);
102
103 _ip_data[ip].eps.noalias() = B * local_u;
104 }
105}
106
107template <typename ShapeFunction, int DisplacementDim>
109 setIPDataInitialConditions(std::string_view const name,
110 double const* values,
111 int const integration_order)
112{
113 if (integration_order !=
114 static_cast<int>(_integration_method.getIntegrationOrder()))
115 {
116 OGS_FATAL(
117 "Setting integration point initial conditions; The integration "
118 "order of the local assembler for element {:d} is different from "
119 "the integration order in the initial condition.",
120 _element.getID());
121 }
122
123 if (name == "sigma")
124 {
125 return setSigma(values);
126 }
127 if (name == "epsilon")
128 {
129 return setEpsilon(values);
130 }
131 if (name == "epsilon_m")
132 {
133 return setEpsilonMechanical(values);
134 }
135
136 return 0;
137}
138
139template <typename ShapeFunction, int DisplacementDim>
141 assembleWithJacobian(double const t, double const dt,
142 std::vector<double> const& local_x,
143 std::vector<double> const& local_x_prev,
144 std::vector<double>& local_rhs_data,
145 std::vector<double>& local_Jac_data)
146{
147 auto const local_matrix_size = local_x.size();
148 assert(local_matrix_size == temperature_size + displacement_size);
149
150 auto T = Eigen::Map<typename ShapeMatricesType::template VectorType<
151 temperature_size> const>(local_x.data() + temperature_index,
153
154 auto u = Eigen::Map<typename ShapeMatricesType::template VectorType<
155 displacement_size> const>(local_x.data() + displacement_index,
157 bool const is_u_non_zero = u.norm() > 0.0;
158
159 auto T_prev = Eigen::Map<typename ShapeMatricesType::template VectorType<
160 temperature_size> const>(local_x_prev.data() + temperature_index,
162
164 local_Jac_data, local_matrix_size, local_matrix_size);
165
166 auto local_rhs = MathLib::createZeroedVector<RhsVector>(local_rhs_data,
167 local_matrix_size);
168
169 typename ShapeMatricesType::template MatrixType<displacement_size,
171 KuT;
173
176
179
180 unsigned const n_integration_points =
181 _integration_method.getNumberOfPoints();
182
183 MPL::VariableArray variables;
184 MPL::VariableArray variables_prev;
186 x_position.setElementID(_element.getID());
187
188 auto const& medium = _process_data.media_map.getMedium(_element.getID());
189 auto const& solid_phase = medium->phase("Solid");
190
191 for (unsigned ip = 0; ip < n_integration_points; ip++)
192 {
193 auto const& w = _ip_data[ip].integration_weight;
194 auto const& N = _ip_data[ip].N;
195 auto const& dNdx = _ip_data[ip].dNdx;
196
197 ParameterLib::SpatialPosition const x_position{
198 std::nullopt, _element.getID(),
201 _element, N))};
202
203 auto const x_coord = x_position.getCoordinates().value()[0];
204 auto const& B =
205 LinearBMatrix::computeBMatrix<DisplacementDim,
206 ShapeFunction::NPOINTS,
208 dNdx, N, x_coord, _is_axially_symmetric);
209
210 auto& sigma = _ip_data[ip].sigma;
211 auto const& sigma_prev = _ip_data[ip].sigma_prev;
212
213 auto& eps = _ip_data[ip].eps;
214 auto const& eps_prev = _ip_data[ip].eps_prev;
215
216 auto& eps_m = _ip_data[ip].eps_m;
217 auto const& eps_m_prev = _ip_data[ip].eps_m_prev;
218
219 auto& state = _ip_data[ip].material_state_variables;
220
221 const double T_ip = N.dot(T); // T at integration point
222 double const T_prev_ip = N.dot(T_prev);
223
224 // Consider also anisotropic thermal expansion.
225 auto const solid_linear_thermal_expansivity_vector =
227 solid_phase
228 .property(
230 .value(variables, x_position, t, dt));
231
233 dthermal_strain =
234 solid_linear_thermal_expansivity_vector * (T_ip - T_prev_ip);
235
236 //
237 // displacement equation, displacement part
238 //
239 // For the restart computation, the displacement may not be
240 // reloaded but the initial strains are always available. For such case,
241 // the following computation is skipped.
242 if (is_u_non_zero)
243 {
244 eps.noalias() = B * u;
245 }
246
247 eps_m.noalias() = eps_m_prev + eps - eps_prev - dthermal_strain;
248
249 variables_prev.stress
251 sigma_prev);
252 variables_prev.mechanical_strain
254 eps_m_prev);
255 variables_prev.temperature = T_ip;
256 variables.mechanical_strain
258 eps_m);
259 variables.temperature = T_ip;
260
261 auto&& solution = _ip_data[ip].solid_material.integrateStress(
262 variables_prev, variables, t, x_position, dt, *state);
263
264 if (!solution)
265 {
266 OGS_FATAL("Computation of local constitutive relation failed.");
267 }
268
270 std::tie(sigma, state, C) = std::move(*solution);
271
272 local_Jac
273 .template block<displacement_size, displacement_size>(
275 .noalias() += B.transpose() * C * B * w;
276
277 typename ShapeMatricesType::template MatrixType<DisplacementDim,
279 N_u = ShapeMatricesType::template MatrixType<
280 DisplacementDim, displacement_size>::Zero(DisplacementDim,
282
283 for (int i = 0; i < DisplacementDim; ++i)
284 {
285 N_u.template block<1, displacement_size / DisplacementDim>(
286 i, i * displacement_size / DisplacementDim)
287 .noalias() = N;
288 }
289
290 auto const rho_s =
291 solid_phase.property(MPL::PropertyType::density)
292 .template value<double>(variables, x_position, t, dt);
293
294 auto const& b = _process_data.specific_body_force;
295 local_rhs.template block<displacement_size, 1>(displacement_index, 0)
296 .noalias() -=
297 (B.transpose() * sigma - N_u.transpose() * rho_s * b) * w;
298
299 //
300 // displacement equation, temperature part
301 // The computation of KuT can be ignored.
302 auto const alpha_T_tensor =
304 solid_linear_thermal_expansivity_vector);
305 KuT.noalias() += B.transpose() * (C * alpha_T_tensor) * N * w;
306
307 if (_ip_data[ip].solid_material.getConstitutiveModel() ==
309 {
310 auto const s = Invariants::deviatoric_projection * sigma;
311 double const norm_s = Invariants::FrobeniusNorm(s);
312 const double creep_coefficient =
313 _ip_data[ip].solid_material.getTemperatureRelatedCoefficient(
314 t, dt, x_position, T_ip, norm_s);
315 KuT.noalias() += creep_coefficient * B.transpose() * s * N * w;
316 }
317
318 //
319 // temperature equation, temperature part;
320 //
321 auto const lambda =
322 solid_phase
323 .property(
325 .value(variables, x_position, t, dt);
326
327 GlobalDimMatrixType const thermal_conductivity =
329
330 KTT.noalias() += dNdx.transpose() * thermal_conductivity * dNdx * w;
331
332 auto const c =
333 solid_phase
334 .property(
336 .template value<double>(variables, x_position, t, dt);
337 DTT.noalias() += N.transpose() * rho_s * c * N * w;
338 }
339
340 // temperature equation, temperature part
341 local_Jac
342 .template block<temperature_size, temperature_size>(temperature_index,
344 .noalias() += KTT + DTT / dt;
345
346 // displacement equation, temperature part
347 local_Jac
348 .template block<displacement_size, temperature_size>(displacement_index,
350 .noalias() -= KuT;
351
352 local_rhs.template block<temperature_size, 1>(temperature_index, 0)
353 .noalias() -= KTT * T + DTT * (T - T_prev) / dt;
354}
355
356template <typename ShapeFunction, int DisplacementDim>
358 assembleWithJacobianForStaggeredScheme(const double t, double const dt,
359 Eigen::VectorXd const& local_x,
360 Eigen::VectorXd const& local_x_prev,
361 int const process_id,
362 std::vector<double>& local_b_data,
363 std::vector<double>& local_Jac_data)
364{
365 // For the equations with pressure
366 if (process_id == _process_data.heat_conduction_process_id)
367 {
369 t, dt, local_x, local_x_prev, local_b_data, local_Jac_data);
370 return;
371 }
372
373 // For the equations with deformation
374 assembleWithJacobianForDeformationEquations(t, dt, local_x, local_x_prev,
375 local_b_data, local_Jac_data);
376}
377
378template <typename ShapeFunction, int DisplacementDim>
381 const double t, double const dt, Eigen::VectorXd const& local_x,
382 Eigen::VectorXd const& local_x_prev, std::vector<double>& local_b_data,
383 std::vector<double>& local_Jac_data)
384{
385 auto const local_T =
386 local_x.template segment<temperature_size>(temperature_index);
387
388 auto const local_T_prev =
389 local_x_prev.template segment<temperature_size>(temperature_index);
390
391 auto const local_u =
392 local_x.template segment<displacement_size>(displacement_index);
393 bool const is_u_non_zero = local_u.norm() > 0.0;
394
395 auto local_Jac = MathLib::createZeroedMatrix<
396 typename ShapeMatricesType::template MatrixType<displacement_size,
398 local_Jac_data, displacement_size, displacement_size);
399
400 auto local_rhs = MathLib::createZeroedVector<
401 typename ShapeMatricesType::template VectorType<displacement_size>>(
402 local_b_data, displacement_size);
403
404 unsigned const n_integration_points =
405 _integration_method.getNumberOfPoints();
406
407 MPL::VariableArray variables;
408 MPL::VariableArray variables_prev;
409 auto const& medium = _process_data.media_map.getMedium(_element.getID());
410 auto const& solid_phase = medium->phase("Solid");
411
412 for (unsigned ip = 0; ip < n_integration_points; ip++)
413 {
414 auto const& w = _ip_data[ip].integration_weight;
415 auto const& N = _ip_data[ip].N;
416 auto const& dNdx = _ip_data[ip].dNdx;
417
418 ParameterLib::SpatialPosition const x_position{
419 std::nullopt, _element.getID(),
422 _element, N))};
423
424 auto const x_coord = x_position.getCoordinates().value()[0];
425
426 auto const& B =
427 LinearBMatrix::computeBMatrix<DisplacementDim,
428 ShapeFunction::NPOINTS,
430 dNdx, N, x_coord, _is_axially_symmetric);
431
432 auto& sigma = _ip_data[ip].sigma;
433 auto const& sigma_prev = _ip_data[ip].sigma_prev;
434
435 auto& eps = _ip_data[ip].eps;
436 auto const& eps_prev = _ip_data[ip].eps_prev;
437
438 auto& eps_m = _ip_data[ip].eps_m;
439 auto const& eps_m_prev = _ip_data[ip].eps_m_prev;
440
441 auto& state = _ip_data[ip].material_state_variables;
442
443 const double T_ip = N.dot(local_T); // T at integration point
444 variables.temperature = T_ip;
445 double const dT_ip = T_ip - N.dot(local_T_prev);
446
447 //
448 // displacement equation, displacement part
449 //
450 // For the restart computation, the displacement may not be
451 // reloaded but the initial strains are always available. For such case,
452 // the following computation is skipped.
453 if (is_u_non_zero)
454 {
455 eps.noalias() = B * local_u;
456 }
457
458 // Consider also anisotropic thermal expansion.
459 auto const solid_linear_thermal_expansivity_vector =
461 solid_phase
462 .property(
464 .value(variables, x_position, t, dt));
465
467 dthermal_strain = solid_linear_thermal_expansivity_vector * dT_ip;
468
469 eps_m.noalias() = eps_m_prev + eps - eps_prev - dthermal_strain;
470
471 variables_prev.stress
473 sigma_prev);
474 variables_prev.mechanical_strain
476 eps_m_prev);
477 variables_prev.temperature = T_ip;
478 variables.mechanical_strain
480 eps_m);
481
482 auto&& solution = _ip_data[ip].solid_material.integrateStress(
483 variables_prev, variables, t, x_position, dt, *state);
484
485 if (!solution)
486 {
487 OGS_FATAL("Computation of local constitutive relation failed.");
488 }
489
491 std::tie(sigma, state, C) = std::move(*solution);
492
493 local_Jac.noalias() += B.transpose() * C * B * w;
494
495 typename ShapeMatricesType::template MatrixType<DisplacementDim,
497 N_u = ShapeMatricesType::template MatrixType<
498 DisplacementDim, displacement_size>::Zero(DisplacementDim,
500
501 for (int i = 0; i < DisplacementDim; ++i)
502 {
503 N_u.template block<1, displacement_size / DisplacementDim>(
504 i, i * displacement_size / DisplacementDim)
505 .noalias() = N;
506 }
507
508 auto const rho_s =
509 solid_phase.property(MPL::PropertyType::density)
510 .template value<double>(variables, x_position, t, dt);
511
512 auto const& b = _process_data.specific_body_force;
513 local_rhs.noalias() -=
514 (B.transpose() * sigma - N_u.transpose() * rho_s * b) * w;
515 }
516}
517
518template <typename ShapeFunction, int DisplacementDim>
521 const double t, double const dt, Eigen::VectorXd const& local_x,
522 Eigen::VectorXd const& local_x_prev, std::vector<double>& local_b_data,
523 std::vector<double>& local_Jac_data)
524{
525 auto const local_T =
526 local_x.template segment<temperature_size>(temperature_index);
527
528 auto const local_T_prev =
529 local_x_prev.template segment<temperature_size>(temperature_index);
530
531 auto local_Jac = MathLib::createZeroedMatrix<
532 typename ShapeMatricesType::template MatrixType<temperature_size,
534 local_Jac_data, temperature_size, temperature_size);
535
536 auto local_rhs = MathLib::createZeroedVector<
537 typename ShapeMatricesType::template VectorType<temperature_size>>(
538 local_b_data, temperature_size);
539
541 mass.setZero(temperature_size, temperature_size);
542
543 typename ShapeMatricesType::NodalMatrixType laplace;
544 laplace.setZero(temperature_size, temperature_size);
545
546 unsigned const n_integration_points =
547 _integration_method.getNumberOfPoints();
548
549 auto const& medium = _process_data.media_map.getMedium(_element.getID());
550 auto const& solid_phase = medium->phase("Solid");
551 MPL::VariableArray variables;
552
553 for (unsigned ip = 0; ip < n_integration_points; ip++)
554 {
555 auto const& w = _ip_data[ip].integration_weight;
556 auto const& N = _ip_data[ip].N;
557 auto const& dNdx = _ip_data[ip].dNdx;
558
559 ParameterLib::SpatialPosition const x_position{
560 std::nullopt, _element.getID(),
563 _element, N))};
564
565 const double T_ip = N.dot(local_T); // T at integration point
566 variables.temperature = T_ip;
567
568 auto const rho_s =
569 solid_phase.property(MPL::PropertyType::density)
570 .template value<double>(variables, x_position, t, dt);
571 auto const c_p =
573 .template value<double>(variables, x_position, t, dt);
574
575 mass.noalias() += N.transpose() * rho_s * c_p * N * w;
576
577 auto const lambda =
578 solid_phase
579 .property(
581 .value(variables, x_position, t, dt);
582
583 GlobalDimMatrixType const thermal_conductivity =
585
586 laplace.noalias() += dNdx.transpose() * thermal_conductivity * dNdx * w;
587 }
588 local_Jac.noalias() += laplace + mass / dt;
589
590 local_rhs.noalias() -=
591 laplace * local_T + mass * (local_T - local_T_prev) / dt;
592}
593
594template <typename ShapeFunction, int DisplacementDim>
595std::size_t
602
603template <typename ShapeFunction, int DisplacementDim>
604std::vector<double>
606{
607 constexpr int kelvin_vector_size =
609
611 [this](std::vector<double>& values)
612 { return getIntPtSigma(0, {}, {}, values); });
613}
614
615template <typename ShapeFunction, int DisplacementDim>
616std::vector<double> const&
618 const double /*t*/,
619 std::vector<GlobalVector*> const& /*x*/,
620 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
621 std::vector<double>& cache) const
622{
624 _ip_data, &IpData::sigma, cache);
625}
626
627template <typename ShapeFunction, int DisplacementDim>
628std::size_t
635
636template <typename ShapeFunction, int DisplacementDim>
637std::vector<double> ThermoMechanicsLocalAssembler<
638 ShapeFunction, DisplacementDim>::getEpsilon() const
639{
640 constexpr int kelvin_vector_size =
642
644 [this](std::vector<double>& values)
645 { return getIntPtEpsilon(0, {}, {}, values); });
646}
647
648template <typename ShapeFunction, int DisplacementDim>
649std::vector<double> const&
651 const double /*t*/,
652 std::vector<GlobalVector*> const& /*x*/,
653 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
654 std::vector<double>& cache) const
655{
657 _ip_data, &IpData::eps, cache);
658}
659
660template <typename ShapeFunction, int DisplacementDim>
661std::vector<double> const&
664 const double /*t*/,
665 std::vector<GlobalVector*> const& /*x*/,
666 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
667 std::vector<double>& cache) const
668{
670 _ip_data, &IpData::eps_m, cache);
671}
672
673template <typename ShapeFunction, int DisplacementDim>
675 ShapeFunction, DisplacementDim>::setEpsilonMechanical(double const* values)
676{
678 values, _ip_data, &IpData::eps_m);
679}
680template <typename ShapeFunction, int DisplacementDim>
681std::vector<double> ThermoMechanicsLocalAssembler<
682 ShapeFunction, DisplacementDim>::getEpsilonMechanical() const
683{
684 constexpr int kelvin_vector_size =
686
688 [this](std::vector<double>& values)
689 { return getIntPtEpsilonMechanical(0, {}, {}, values); });
690}
691
692template <typename ShapeFunction, int DisplacementDim>
695{
696 return _integration_method.getNumberOfPoints();
697}
698
699template <typename ShapeFunction, int DisplacementDim>
701 DisplacementDim>::getMaterialID() const
702{
703 return _process_data.material_ids == nullptr
704 ? 0
705 : (*_process_data.material_ids)[_element.getID()];
706}
707
708template <typename ShapeFunction, int DisplacementDim>
710 DisplacementDim>::MaterialStateVariables const&
712 getMaterialStateVariablesAt(unsigned integration_point) const
713{
714 return *_ip_data[integration_point].material_state_variables;
715}
716} // namespace ThermoMechanics
717} // namespace ProcessLib
#define OGS_FATAL(...)
Definition Error.h:19
std::size_t getID() const
Returns the ID of the element.
Definition Element.h:80
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
void setInitialConditionsConcrete(Eigen::VectorXd const local_x, double const t, int const process_id) override
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)
constexpr Eigen::Matrix< double, GlobalDim, GlobalDim > formEigenTensor(MaterialPropertyLib::PropertyDataType const &values)
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, 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
static double FrobeniusNorm(Eigen::Matrix< double, KelvinVectorSize, 1 > const &deviatoric_v)
Get the norm of the deviatoric stress.
static Eigen::Matrix< double, KelvinVectorSize, KelvinVectorSize > const deviatoric_projection