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