OGS
ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim > Class Template Reference

Detailed Description

template<typename ShapeFunction, typename IntegrationMethod, int DisplacementDim>
class ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >

Local assembler of ThermoMechanics process.

In this assembler, the solid density can be expressed as an exponential function of temperature, if a temperature dependent solid density is assumed. The theory behind this is given below.

During the temperature variation, the solid mass balance is

\[ \rho_s^{n+1} V^{n+1} = \rho_s^n V^n, \]

with \(\rho_s\) the density, \(V\), the volume, and \(n\) the index of the time step. Under pure thermo-mechanics condition, the volume change along with the temperature change is given by

\[ V^{n+1} = V^{n} + \alpha_T dT V^{n} = (1+\alpha_T dT)V^{n}, \]

where \(\alpha_T\) is the volumetric thermal expansivity of the solid. This gives

\[ \rho_s^{n+1} + \alpha_T dT \rho_s^{n+1} = \rho_s^n. \]

Therefore, we obtain the differential expression of the temperature dependent solid density as

\[ \frac{d \rho_s}{d T} = -\alpha_T \rho_s. \]

The solution of the above ODE, i.e the density expression, is given by

\[ \rho_s = {\rho_0} \mathrm{exp} (- \alpha_T (T-T0)), \]

with reference density \(\rho_0\) at a reference temperature of \(T_0\).

An MPL property with the type of Exponential (see MaterialPropertyLib::Exponential) can be used for the temperature dependent solid property.

Definition at line 127 of file ThermoMechanicsFEM.h.

#include <ThermoMechanicsFEM.h>

Inheritance diagram for ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >:
[legend]
Collaboration diagram for ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >:
[legend]

Public Types

using ShapeMatricesType = ShapeMatrixPolicyType< ShapeFunction, DisplacementDim >
 
using ShapeMatrices = typename ShapeMatricesType::ShapeMatrices
 
using GlobalDimMatrixType = typename ShapeMatricesType::GlobalDimMatrixType
 
using BMatricesType = BMatrixPolicyType< ShapeFunction, DisplacementDim >
 
using Invariants = MathLib::KelvinVector::Invariants< KelvinVectorSize >
 
using NodalForceVectorType = typename BMatricesType::NodalForceVectorType
 
using RhsVector = typename ShapeMatricesType::template VectorType< ShapeFunction::NPOINTS+ShapeFunction::NPOINTS *DisplacementDim >
 
using JacobianMatrix = typename ShapeMatricesType::template MatrixType< ShapeFunction::NPOINTS+ShapeFunction::NPOINTS *DisplacementDim, ShapeFunction::NPOINTS+ShapeFunction::NPOINTS *DisplacementDim >
 
using IpData = IntegrationPointData< BMatricesType, ShapeMatricesType, DisplacementDim >
 

Public Member Functions

 ThermoMechanicsLocalAssembler (ThermoMechanicsLocalAssembler const &)=delete
 
 ThermoMechanicsLocalAssembler (ThermoMechanicsLocalAssembler &&)=delete
 
 ThermoMechanicsLocalAssembler (MeshLib::Element const &e, std::size_t const, bool const is_axially_symmetric, unsigned const integration_order, ThermoMechanicsProcessData< DisplacementDim > &process_data)
 
std::size_t setIPDataInitialConditions (std::string const &name, double const *values, int const integration_order) override
 Returns number of read integration points. More...
 
void assemble (double const, double const, std::vector< double > const &, std::vector< double > const &, std::vector< double > &, std::vector< double > &, std::vector< double > &) override
 
void assembleWithJacobianForStaggeredScheme (double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_xdot, const double dxdot_dx, const double dx_dx, int const process_id, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data) override
 
void assembleWithJacobian (double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_xdot, const double dxdot_dx, const double dx_dx, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_rhs_data, std::vector< double > &local_Jac_data) override
 
void initializeConcrete () override
 
void postTimestepConcrete (Eigen::VectorXd const &, double const, double const) override
 
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix (const unsigned integration_point) const override
 Provides the shape matrix at the given integration point. More...
 
- Public Member Functions inherited from ProcessLib::LocalAssemblerInterface
virtual ~LocalAssemblerInterface ()=default
 
void setInitialConditions (std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, double const t, bool const use_monolithic_scheme, int const process_id)
 
virtual void initialize (std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
 
virtual void preAssemble (double const, double const, std::vector< double > const &)
 
virtual void assembleForStaggeredScheme (double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_xdot, int const process_id, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data)
 
virtual void computeSecondaryVariable (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_dot, int const process_id)
 
virtual void preTimestep (std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, double const t, double const delta_t)
 
virtual void postTimestep (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, double const t, double const dt)
 
void postNonLinearSolver (std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, GlobalVector const &xdot, double const t, double const dt, bool const use_monolithic_scheme, int const process_id)
 
virtual std::vector< double > interpolateNodalValuesToIntegrationPoints (std::vector< double > const &)
 
virtual Eigen::Vector3d getFlux (MathLib::Point3d const &, double const, std::vector< double > const &) const
 
virtual Eigen::Vector3d getFlux (MathLib::Point3d const &, double const, std::vector< std::vector< double >> const &) const
 Fits to staggered scheme. More...
 
- Public Member Functions inherited from NumLib::ExtrapolatableElement
virtual ~ExtrapolatableElement ()=default
 

Static Public Attributes

static int const KelvinVectorSize
 

Private Member Functions

void assembleWithJacobianForDeformationEquations (const double t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_xdot, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data)
 
void assembleWithJacobianForHeatConductionEquations (const double t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_xdot, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data)
 
std::size_t setSigma (double const *values)
 
std::vector< double > getSigma () const override
 
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
 
std::size_t setEpsilon (double const *values)
 
std::vector< double > getEpsilon () const 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
 
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
 
std::size_t setEpsilonMechanical (double const *values)
 
std::vector< double > getEpsilonMechanical () const override
 
unsigned getNumberOfIntegrationPoints () const override
 
MaterialLib::Solids::MechanicsBase< DisplacementDim >::MaterialStateVariables const & getMaterialStateVariablesAt (unsigned integration_point) const override
 

Private Attributes

ThermoMechanicsProcessData< DisplacementDim > & _process_data
 
std::vector< IpData, Eigen::aligned_allocator< IpData > > _ip_data
 
IntegrationMethod _integration_method
 
MeshLib::Element const & _element
 
SecondaryData< typename ShapeMatrices::ShapeType > _secondary_data
 
bool const _is_axially_symmetric
 

Static Private Attributes

static const int temperature_index = 0
 
static const int temperature_size = ShapeFunction::NPOINTS
 
static const int displacement_index = ShapeFunction::NPOINTS
 
static const int displacement_size
 

Member Typedef Documentation

◆ BMatricesType

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
using ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::BMatricesType = BMatrixPolicyType<ShapeFunction, DisplacementDim>

Definition at line 138 of file ThermoMechanicsFEM.h.

◆ GlobalDimMatrixType

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
using ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::GlobalDimMatrixType = typename ShapeMatricesType::GlobalDimMatrixType

Definition at line 137 of file ThermoMechanicsFEM.h.

◆ Invariants

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
using ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::Invariants = MathLib::KelvinVector::Invariants<KelvinVectorSize>

Definition at line 142 of file ThermoMechanicsFEM.h.

◆ IpData

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
using ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::IpData = IntegrationPointData<BMatricesType, ShapeMatricesType, DisplacementDim>

Definition at line 150 of file ThermoMechanicsFEM.h.

◆ JacobianMatrix

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
using ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::JacobianMatrix = typename ShapeMatricesType::template MatrixType< ShapeFunction::NPOINTS + ShapeFunction::NPOINTS * DisplacementDim, ShapeFunction::NPOINTS + ShapeFunction::NPOINTS * DisplacementDim>

Definition at line 147 of file ThermoMechanicsFEM.h.

◆ NodalForceVectorType

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
using ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::NodalForceVectorType = typename BMatricesType::NodalForceVectorType

Definition at line 144 of file ThermoMechanicsFEM.h.

◆ RhsVector

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
using ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::RhsVector = typename ShapeMatricesType::template VectorType< ShapeFunction::NPOINTS + ShapeFunction::NPOINTS * DisplacementDim>

Definition at line 145 of file ThermoMechanicsFEM.h.

◆ ShapeMatrices

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
using ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::ShapeMatrices = typename ShapeMatricesType::ShapeMatrices

Definition at line 136 of file ThermoMechanicsFEM.h.

◆ ShapeMatricesType

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
using ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::ShapeMatricesType = ShapeMatrixPolicyType<ShapeFunction, DisplacementDim>

Definition at line 131 of file ThermoMechanicsFEM.h.

Constructor & Destructor Documentation

◆ ThermoMechanicsLocalAssembler() [1/3]

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::ThermoMechanicsLocalAssembler ( ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim > const &  )
delete

◆ ThermoMechanicsLocalAssembler() [2/3]

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::ThermoMechanicsLocalAssembler ( ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim > &&  )
delete

◆ ThermoMechanicsLocalAssembler() [3/3]

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::ThermoMechanicsLocalAssembler ( MeshLib::Element const &  e,
std::size_t const  ,
bool const  is_axially_symmetric,
unsigned const  integration_order,
ThermoMechanicsProcessData< DisplacementDim > &  process_data 
)

Definition at line 29 of file ThermoMechanicsFEM-impl.h.

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 }
virtual std::size_t getID() const final
Returns the ID of the element.
Definition: Element.h:82
void setElementID(std::size_t element_id)
SecondaryData< typename ShapeMatrices::ShapeType > _secondary_data
ShapeMatrixPolicyType< ShapeFunction, DisplacementDim > ShapeMatricesType
ThermoMechanicsProcessData< DisplacementDim > & _process_data
std::vector< IpData, Eigen::aligned_allocator< IpData > > _ip_data
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)
constexpr int kelvin_vector_dimensions(int const displacement_dim)
Kelvin vector dimensions for given displacement dimension.
Definition: KelvinVector.h:23
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::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N
Definition: SecondaryData.h:28

References ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::_element, ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::_integration_method, ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::_ip_data, ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::_process_data, ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::_secondary_data, MeshLib::Element::getID(), NumLib::initShapeMatrices(), MathLib::KelvinVector::kelvin_vector_dimensions(), ProcessLib::HeatTransportBHE::SecondaryData< ShapeMatrixType >::N, MaterialLib::Solids::selectSolidConstitutiveRelation(), and ParameterLib::SpatialPosition::setElementID().

Member Function Documentation

◆ assemble()

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
void ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::assemble ( double const  ,
double const  ,
std::vector< double > const &  ,
std::vector< double > const &  ,
std::vector< double > &  ,
std::vector< double > &  ,
std::vector< double > &   
)
inlineoverridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 170 of file ThermoMechanicsFEM.h.

176  {
177  OGS_FATAL(
178  "ThermoMechanicsLocalAssembler: assembly without jacobian is not "
179  "implemented.");
180  }
#define OGS_FATAL(...)
Definition: Error.h:26

References OGS_FATAL.

◆ assembleWithJacobian()

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
void ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::assembleWithJacobian ( double const  t,
double const  dt,
std::vector< double > const &  local_x,
std::vector< double > const &  local_xdot,
const double  dxdot_dx,
const double  dx_dx,
std::vector< double > &  local_M_data,
std::vector< double > &  local_K_data,
std::vector< double > &  local_rhs_data,
std::vector< double > &  local_Jac_data 
)
overridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 121 of file ThermoMechanicsFEM-impl.h.

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,
137 
138  auto u = Eigen::Map<typename ShapeMatricesType::template VectorType<
139  displacement_size> const>(local_x.data() + displacement_index,
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,
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,
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>(
258  .noalias() += B.transpose() * C * B * w;
259 
260  typename ShapeMatricesType::template MatrixType<DisplacementDim,
262  N_u = ShapeMatricesType::template MatrixType<
263  DisplacementDim, displacement_size>::Zero(DisplacementDim,
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,
327  .noalias() += KTT + DTT / dt;
328 
329  // displacement equation, temperature part
330  local_Jac
331  .template block<displacement_size, temperature_size>(displacement_index,
333  .noalias() -= KuT;
334 
335  local_rhs.template block<temperature_size, 1>(temperature_index, 0)
336  .noalias() -= KTT * T + DTT * T_dot;
337 }
void setIntegrationPoint(unsigned integration_point)
MatrixType< _kelvin_vector_size, _number_of_dof > BMatrixType
Definition: BMatrixPolicy.h:55
typename ShapeMatricesType::GlobalDimMatrixType GlobalDimMatrixType
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
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), kelvin_vector_dimensions(DisplacementDim), Eigen::RowMajor > KelvinMatrixType
Definition: KelvinVector.h:56
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
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
Definition: KelvinVector.h:66

References MaterialPropertyLib::c, ProcessLib::LinearBMatrix::computeBMatrix(), MaterialLib::Solids::CreepBGRa, MaterialPropertyLib::density, MathLib::KelvinVector::kelvinVectorToSymmetricTensor(), OGS_FATAL, ParameterLib::SpatialPosition::setElementID(), ParameterLib::SpatialPosition::setIntegrationPoint(), MaterialPropertyLib::specific_heat_capacity, MaterialPropertyLib::temperature, MaterialPropertyLib::thermal_conductivity, and MaterialPropertyLib::thermal_expansivity.

◆ assembleWithJacobianForDeformationEquations()

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
void ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::assembleWithJacobianForDeformationEquations ( const double  t,
double const  dt,
Eigen::VectorXd const &  local_x,
Eigen::VectorXd const &  local_xdot,
std::vector< double > &  local_b_data,
std::vector< double > &  local_Jac_data 
)
private

Assemble local matrices and vectors arise from the linearized discretized weak form of the residual of the momentum balance equation,

\[ \nabla (\sigma - \mathbf{D} \alpha_T (T-T_0) \mathrm{I}) = f \]

where \( \sigma\) is the effective stress tensor, \(\mathbf{D}\) is the tangential operator, \(T\) is the temperature, \(T_0\) is the initial temperature, \(\alpha_T\) is the linear thermal expansion, \(\mathrm{I}\) is the identity tensor, and \(f\) is the body force.

Parameters
tTime
dtTime increment
local_xNodal values of \(x\) of all processes of an element.
local_xdotNodal values of \(\dot{x}\) of all processes of an element.
local_b_dataRight hand side vector of an element.
local_Jac_dataElement Jacobian matrix for the Newton-Raphson method.

Definition at line 367 of file ThermoMechanicsFEM-impl.h.

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,
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,
486  N_u = ShapeMatricesType::template MatrixType<
487  DisplacementDim, displacement_size>::Zero(DisplacementDim,
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 }
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)

References ProcessLib::LinearBMatrix::computeBMatrix(), MathLib::createZeroedMatrix(), MathLib::createZeroedVector(), MaterialPropertyLib::density, OGS_FATAL, ParameterLib::SpatialPosition::setElementID(), ParameterLib::SpatialPosition::setIntegrationPoint(), MaterialPropertyLib::temperature, and MaterialPropertyLib::thermal_expansivity.

◆ assembleWithJacobianForHeatConductionEquations()

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
void ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::assembleWithJacobianForHeatConductionEquations ( const double  t,
double const  dt,
Eigen::VectorXd const &  local_x,
Eigen::VectorXd const &  local_xdot,
std::vector< double > &  local_b_data,
std::vector< double > &  local_Jac_data 
)
private

Assemble local matrices and vectors arise from the linearized discretized weak form of the residual of the energy balance equation,

\[ \rho c_p \cdot{T} - \nabla (\mathbf{K} (\nabla T) = Q_T \]

where \( rho\) is the solid density, \( c_p\) is the specific heat capacity, \( \mathbf{K} \) is the thermal conductivity, and \( Q_T\) is the source/sink term.

Parameters
tTime
dtTime increment
local_xNodal values of \(x\) of all processes of an element.
local_xdotNodal values of \(\dot{x}\) of all processes of an element.
local_b_dataRight hand side vector of an element.
local_Jac_dataElement Jacobian matrix for the Newton-Raphson method.

Definition at line 510 of file ThermoMechanicsFEM-impl.h.

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,
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 }

References MathLib::createZeroedMatrix(), MathLib::createZeroedVector(), MaterialPropertyLib::density, ParameterLib::SpatialPosition::setElementID(), ParameterLib::SpatialPosition::setIntegrationPoint(), MaterialPropertyLib::specific_heat_capacity, MaterialPropertyLib::temperature, and MaterialPropertyLib::thermal_conductivity.

◆ assembleWithJacobianForStaggeredScheme()

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
void ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::assembleWithJacobianForStaggeredScheme ( double const  t,
double const  dt,
Eigen::VectorXd const &  local_x,
Eigen::VectorXd const &  local_xdot,
const double  dxdot_dx,
const double  dx_dx,
int const  process_id,
std::vector< double > &  local_M_data,
std::vector< double > &  local_K_data,
std::vector< double > &  local_b_data,
std::vector< double > &  local_Jac_data 
)
overridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 342 of file ThermoMechanicsFEM-impl.h.

350 {
351  // For the equations with pressure
352  if (process_id == _process_data.heat_conduction_process_id)
353  {
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 }
void assembleWithJacobianForHeatConductionEquations(const double t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_xdot, 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_xdot, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data)

◆ getEpsilon()

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
std::vector< double > ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::getEpsilon
overrideprivatevirtual

Implements ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssemblerInterface< DisplacementDim >.

Definition at line 632 of file ThermoMechanicsFEM-impl.h.

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 }
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

References MathLib::KelvinVector::kelvin_vector_dimensions().

◆ getEpsilonMechanical()

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
std::vector< double > ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::getEpsilonMechanical
overrideprivatevirtual

Implements ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssemblerInterface< DisplacementDim >.

Definition at line 683 of file ThermoMechanicsFEM-impl.h.

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 }
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

References MathLib::KelvinVector::kelvin_vector_dimensions().

◆ getIntPtEpsilon()

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
std::vector< double > const & ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::getIntPtEpsilon ( const double  t,
std::vector< GlobalVector * > const &  x,
std::vector< NumLib::LocalToGlobalIndexMap const * > const &  dof_table,
std::vector< double > &  cache 
) const
overrideprivatevirtual

Implements ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssemblerInterface< DisplacementDim >.

Definition at line 645 of file ThermoMechanicsFEM-impl.h.

651 {
652  return ProcessLib::getIntegrationPointKelvinVectorData<DisplacementDim>(
653  _ip_data, &IpData::eps, cache);
654 }

◆ getIntPtEpsilonMechanical()

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
std::vector< double > const & ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::getIntPtEpsilonMechanical ( const double  t,
std::vector< GlobalVector * > const &  x,
std::vector< NumLib::LocalToGlobalIndexMap const * > const &  dof_table,
std::vector< double > &  cache 
) const
overrideprivatevirtual

Implements ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssemblerInterface< DisplacementDim >.

Definition at line 659 of file ThermoMechanicsFEM-impl.h.

665 {
666  return ProcessLib::getIntegrationPointKelvinVectorData<DisplacementDim>(
667  _ip_data, &IpData::eps_m, cache);
668 }

◆ getIntPtSigma()

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
std::vector< double > const & ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::getIntPtSigma ( const double  t,
std::vector< GlobalVector * > const &  x,
std::vector< NumLib::LocalToGlobalIndexMap const * > const &  dof_table,
std::vector< double > &  cache 
) const
overrideprivatevirtual

Implements ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssemblerInterface< DisplacementDim >.

Definition at line 608 of file ThermoMechanicsFEM-impl.h.

614 {
615  return ProcessLib::getIntegrationPointKelvinVectorData<DisplacementDim>(
616  _ip_data, &IpData::sigma, cache);
617 }

◆ getMaterialStateVariablesAt()

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
MaterialLib::Solids::MechanicsBase< DisplacementDim >::MaterialStateVariables const & ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::getMaterialStateVariablesAt ( unsigned  integration_point) const
overrideprivatevirtual

Implements ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssemblerInterface< DisplacementDim >.

Definition at line 707 of file ThermoMechanicsFEM-impl.h.

709 {
710  return *_ip_data[integration_point].material_state_variables;
711 }

◆ getNumberOfIntegrationPoints()

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
unsigned ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::getNumberOfIntegrationPoints
overrideprivatevirtual

Implements ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssemblerInterface< DisplacementDim >.

Definition at line 697 of file ThermoMechanicsFEM-impl.h.

698 {
699  return _integration_method.getNumberOfPoints();
700 }

◆ getShapeMatrix()

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
Eigen::Map<const Eigen::RowVectorXd> ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::getShapeMatrix ( const unsigned  integration_point) const
inlineoverridevirtual

Provides the shape matrix at the given integration point.

Implements NumLib::ExtrapolatableElement.

Definition at line 241 of file ThermoMechanicsFEM.h.

243  {
244  auto const& N = _secondary_data.N[integration_point];
245 
246  // assumes N is stored contiguously in memory
247  return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
248  }

References ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::_secondary_data, and ProcessLib::HeatTransportBHE::SecondaryData< ShapeMatrixType >::N.

◆ getSigma()

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
std::vector< double > ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::getSigma
overrideprivatevirtual

Implements ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssemblerInterface< DisplacementDim >.

Definition at line 595 of file ThermoMechanicsFEM-impl.h.

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 }
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

References MathLib::KelvinVector::kelvin_vector_dimensions().

◆ initializeConcrete()

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
void ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::initializeConcrete ( )
inlineoverridevirtual

Set initial stress from parameter.

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 199 of file ThermoMechanicsFEM.h.

200  {
201  unsigned const n_integration_points =
202  _integration_method.getNumberOfPoints();
203  for (unsigned ip = 0; ip < n_integration_points; ip++)
204  {
205  auto& ip_data = _ip_data[ip];
206 
208  if (_process_data.initial_stress != nullptr)
209  {
210  ParameterLib::SpatialPosition const x_position{
211  std::nullopt, _element.getID(), ip,
213  NumLib::interpolateCoordinates<ShapeFunction,
215  _element, ip_data.N))};
216 
217  ip_data.sigma =
219  DisplacementDim>((*_process_data.initial_stress)(
220  std::numeric_limits<
221  double>::quiet_NaN() /* time independent */,
222  x_position));
223  }
224 
225  ip_data.pushBackState();
226  }
227  }
Eigen::Matrix< double, Eigen::MatrixBase< Derived >::RowsAtCompileTime, 1 > symmetricTensorToKelvinVector(Eigen::MatrixBase< Derived > const &v)
Definition: KelvinVector.h:170
MathLib::TemplatePoint< double, 3 > Point3d
std::array< double, 3 > interpolateCoordinates(MeshLib::Element const &e, typename ShapeMatricesType::ShapeMatrices::ShapeType const &N)

References ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::_element, ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::_integration_method, ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::_ip_data, ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::_process_data, MeshLib::Element::getID(), NumLib::interpolateCoordinates(), and MathLib::KelvinVector::symmetricTensorToKelvinVector().

◆ postTimestepConcrete()

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
void ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::postTimestepConcrete ( Eigen::VectorXd const &  ,
double const  ,
double const   
)
inlineoverridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 229 of file ThermoMechanicsFEM.h.

231  {
232  unsigned const n_integration_points =
233  _integration_method.getNumberOfPoints();
234 
235  for (unsigned ip = 0; ip < n_integration_points; ip++)
236  {
237  _ip_data[ip].pushBackState();
238  }
239  }

References ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::_integration_method, and ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::_ip_data.

◆ setEpsilon()

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
std::size_t ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::setEpsilon ( double const *  values)
private

Definition at line 623 of file ThermoMechanicsFEM-impl.h.

624 {
625  return ProcessLib::setIntegrationPointKelvinVectorData<DisplacementDim>(
626  values, _ip_data, &IpData::eps);
627 }

◆ setEpsilonMechanical()

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
std::size_t ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::setEpsilonMechanical ( double const *  values)
private

Definition at line 674 of file ThermoMechanicsFEM-impl.h.

675 {
676  return ProcessLib::setIntegrationPointKelvinVectorData<DisplacementDim>(
677  values, _ip_data, &IpData::eps_m);
678 }

◆ setIPDataInitialConditions()

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
std::size_t ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::setIPDataInitialConditions ( std::string const &  name,
double const *  values,
int const  integration_order 
)
overridevirtual

Returns number of read integration points.

Implements ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssemblerInterface< DisplacementDim >.

Definition at line 88 of file ThermoMechanicsFEM-impl.h.

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 }

References MaterialPropertyLib::name, and OGS_FATAL.

◆ setSigma()

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
std::size_t ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::setSigma ( double const *  values)
private

Definition at line 586 of file ThermoMechanicsFEM-impl.h.

587 {
588  return ProcessLib::setIntegrationPointKelvinVectorData<DisplacementDim>(
589  values, _ip_data, &IpData::sigma);
590 }

Member Data Documentation

◆ _element

◆ _integration_method

◆ _ip_data

◆ _is_axially_symmetric

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
bool const ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::_is_axially_symmetric
private

Definition at line 351 of file ThermoMechanicsFEM.h.

◆ _process_data

◆ _secondary_data

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
SecondaryData<typename ShapeMatrices::ShapeType> ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::_secondary_data
private

◆ displacement_index

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
const int ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::displacement_index = ShapeFunction::NPOINTS
staticprivate

Definition at line 355 of file ThermoMechanicsFEM.h.

◆ displacement_size

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
const int ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::displacement_size
staticprivate
Initial value:
=
ShapeFunction::NPOINTS * DisplacementDim

Definition at line 356 of file ThermoMechanicsFEM.h.

◆ KelvinVectorSize

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
int const ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::KelvinVectorSize
static
Initial value:

Definition at line 140 of file ThermoMechanicsFEM.h.

◆ temperature_index

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
const int ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::temperature_index = 0
staticprivate

Definition at line 353 of file ThermoMechanicsFEM.h.

◆ temperature_size

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
const int ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::temperature_size = ShapeFunction::NPOINTS
staticprivate

Definition at line 354 of file ThermoMechanicsFEM.h.


The documentation for this class was generated from the following files: