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

Detailed Description

template<typename ShapeFunction, int DisplacementDim>
class ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, 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 120 of file ThermoMechanicsFEM.h.

#include <ThermoMechanicsFEM.h>

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

Public Types

using ShapeMatricesType
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
using JacobianMatrix
using IpData

Public Member Functions

 ThermoMechanicsLocalAssembler (ThermoMechanicsLocalAssembler const &)=delete
 ThermoMechanicsLocalAssembler (ThermoMechanicsLocalAssembler &&)=delete
 ThermoMechanicsLocalAssembler (MeshLib::Element const &e, std::size_t const, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, ThermoMechanicsProcessData< DisplacementDim > &process_data)
std::size_t setIPDataInitialConditions (std::string_view const name, double const *values, int const integration_order) override
 Returns number of read integration points.
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_x_prev, int const process_id, 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_x_prev, std::vector< double > &local_rhs_data, std::vector< double > &local_Jac_data) override
void initializeConcrete () override
void setInitialConditionsConcrete (Eigen::VectorXd const local_x, double const t, int const process_id) override
void postTimestepConcrete (Eigen::VectorXd const &, Eigen::VectorXd const &, double const, double const, int const) override
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix (const unsigned integration_point) const override
 Provides the shape matrix at the given integration point.
int getNumberOfVectorElementsForDeformation () const override
Public Member Functions inherited from ProcessLib::LocalAssemblerInterface
virtual ~LocalAssemblerInterface ()=default
virtual void setInitialConditions (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, double const t, 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_x_prev, 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_prev, 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, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, int const process_id)
void postNonLinearSolver (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, int const process_id)
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.
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_x_prev, 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_x_prev, 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
int getMaterialID () 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
NumLib::GenericIntegrationMethod const & _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, int DisplacementDim>
using ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, DisplacementDim >::BMatricesType = BMatrixPolicyType<ShapeFunction, DisplacementDim>

Definition at line 131 of file ThermoMechanicsFEM.h.

◆ GlobalDimMatrixType

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

Definition at line 130 of file ThermoMechanicsFEM.h.

◆ Invariants

Definition at line 135 of file ThermoMechanicsFEM.h.

◆ IpData

◆ JacobianMatrix

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

Definition at line 140 of file ThermoMechanicsFEM.h.

◆ NodalForceVectorType

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

Definition at line 137 of file ThermoMechanicsFEM.h.

◆ RhsVector

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

Definition at line 138 of file ThermoMechanicsFEM.h.

◆ ShapeMatrices

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

Definition at line 129 of file ThermoMechanicsFEM.h.

◆ ShapeMatricesType

template<typename ShapeFunction, int DisplacementDim>
using ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, DisplacementDim >::ShapeMatricesType
Initial value:
EigenFixedShapeMatrixPolicy< ShapeFunction, GlobalDim > ShapeMatrixPolicyType

Definition at line 124 of file ThermoMechanicsFEM.h.

Constructor & Destructor Documentation

◆ ThermoMechanicsLocalAssembler() [1/3]

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

◆ ThermoMechanicsLocalAssembler() [2/3]

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

◆ ThermoMechanicsLocalAssembler() [3/3]

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

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

28 _element(e),
30{
31 unsigned const n_integration_points =
32 _integration_method.getNumberOfPoints();
33
36
37 auto const shape_matrices =
41
43 _process_data.solid_materials, _process_data.material_ids, e.getID());
44
45 for (unsigned ip = 0; ip < n_integration_points; ip++)
46 {
47 _ip_data.emplace_back(solid_material);
48 auto& ip_data = _ip_data[ip];
49 ip_data.integration_weight =
50 _integration_method.getWeightedPoint(ip).getWeight() *
51 shape_matrices[ip].integralMeasure * shape_matrices[ip].detJ;
52
53 static const int kelvin_vector_size =
55 // Initialize current time step values
56 ip_data.sigma.setZero(kelvin_vector_size);
57 ip_data.eps.setZero(kelvin_vector_size);
58
59 // Previous time step values are not initialized and are set later.
60 ip_data.sigma_prev.resize(kelvin_vector_size);
61 ip_data.eps_prev.resize(kelvin_vector_size);
62
63 ip_data.eps_m.setZero(kelvin_vector_size);
64 ip_data.eps_m_prev.setZero(kelvin_vector_size);
66 x_position.setElementID(_element.getID());
68 ip_data.dNdx = shape_matrices[ip].dNdx;
69
71 }
72}
SecondaryData< typename ShapeMatrices::ShapeType > _secondary_data
ShapeMatrixPolicyType< ShapeFunction, DisplacementDim > ShapeMatricesType
ThermoMechanicsProcessData< DisplacementDim > & _process_data
NumLib::GenericIntegrationMethod const & _integration_method
std::vector< IpData, Eigen::aligned_allocator< IpData > > _ip_data
auto & selectSolidConstitutiveRelation(SolidMaterialsMap 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.

References _element, _integration_method, _ip_data, _is_axially_symmetric, _process_data, _secondary_data, MeshLib::Element::getID(), NumLib::initShapeMatrices(), MathLib::KelvinVector::kelvin_vector_dimensions(), MaterialLib::Solids::selectSolidConstitutiveRelation(), and ParameterLib::SpatialPosition::setElementID().

Member Function Documentation

◆ assemble()

template<typename ShapeFunction, int DisplacementDim>
void ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, 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 163 of file ThermoMechanicsFEM.h.

169 {
170 OGS_FATAL(
171 "ThermoMechanicsLocalAssembler: assembly without jacobian is not "
172 "implemented.");
173 }
#define OGS_FATAL(...)
Definition Error.h:19

References OGS_FATAL.

◆ assembleWithJacobian()

template<typename ShapeFunction, int DisplacementDim>
void ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, DisplacementDim >::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 )
overridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

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

147{
148 auto const local_matrix_size = local_x.size();
150
154
158 bool const is_u_non_zero = u.norm() > 0.0;
159
163
166
169
172 KuT;
174
177
180
181 unsigned const n_integration_points =
182 _integration_method.getNumberOfPoints();
183
187 x_position.setElementID(_element.getID());
188
189 auto const& medium = _process_data.media_map.getMedium(_element.getID());
190 auto const& solid_phase =
192
193 for (unsigned ip = 0; ip < n_integration_points; ip++)
194 {
195 auto const& w = _ip_data[ip].integration_weight;
196 auto const& N = _ip_data[ip].N;
197 auto const& dNdx = _ip_data[ip].dNdx;
198
200 std::nullopt, _element.getID(),
203 _element, N))};
204
205 auto const x_coord = x_position.getCoordinates().value()[0];
206 auto const& B =
211
212 auto& sigma = _ip_data[ip].sigma;
213 auto const& sigma_prev = _ip_data[ip].sigma_prev;
214
215 auto& eps = _ip_data[ip].eps;
216 auto const& eps_prev = _ip_data[ip].eps_prev;
217
218 auto& eps_m = _ip_data[ip].eps_m;
219 auto const& eps_m_prev = _ip_data[ip].eps_m_prev;
220
221 auto& state = _ip_data[ip].material_state_variables;
222
223 const double T_ip = N.dot(T); // T at integration point
224 double const T_prev_ip = N.dot(T_prev);
225
226 // Consider also anisotropic thermal expansion.
230 .property(
233
237
238 //
239 // displacement equation, displacement part
240 //
241 // For the restart computation, the displacement may not be
242 // reloaded but the initial strains are always available. For such case,
243 // the following computation is skipped.
244 if (is_u_non_zero)
245 {
246 eps.noalias() = B * u;
247 }
248
249 eps_m.noalias() = eps_m_prev + eps - eps_prev - dthermal_strain;
250
251 variables_prev.stress
253 sigma_prev);
254 variables_prev.mechanical_strain
256 eps_m_prev);
257 variables_prev.temperature = T_ip;
258 variables.mechanical_strain
260 eps_m);
261 variables.temperature = T_ip;
262
263 auto&& solution = _ip_data[ip].solid_material.integrateStress(
265
266 if (!solution)
267 {
268 OGS_FATAL("Computation of local constitutive relation failed.");
269 }
270
273
277 .noalias() += B.transpose() * C * B * w;
278
284
285 for (int i = 0; i < DisplacementDim; ++i)
286 {
289 .noalias() = N;
290 }
291
292 auto const rho_s =
294 .template value<double>(variables, x_position, t, dt);
295
296 auto const& b = _process_data.specific_body_force;
298 .noalias() -=
299 (B.transpose() * sigma - N_u.transpose() * rho_s * b) * w;
300
301 //
302 // displacement equation, temperature part
303 // The computation of KuT can be ignored.
304 auto const alpha_T_tensor =
307 KuT.noalias() += B.transpose() * (C * alpha_T_tensor) * N * w;
308
309 if (_ip_data[ip].solid_material.getConstitutiveModel() ==
311 {
313 double const norm_s = Invariants::FrobeniusNorm(s);
314 const double creep_coefficient =
315 _ip_data[ip].solid_material.getTemperatureRelatedCoefficient(
317 KuT.noalias() += creep_coefficient * B.transpose() * s * N * w;
318 }
319
320 //
321 // temperature equation, temperature part;
322 //
323 auto const lambda =
325 .property(
327 .value(variables, x_position, t, dt);
328
331
332 KTT.noalias() += dNdx.transpose() * thermal_conductivity * dNdx * w;
333
334 auto const c =
336 .property(
338 .template value<double>(variables, x_position, t, dt);
339 DTT.noalias() += N.transpose() * rho_s * c * N * w;
340 }
341
342 // temperature equation, temperature part
346 .noalias() += KTT + DTT / dt;
347
348 // displacement equation, temperature part
352 .noalias() -= KuT;
353
355 .noalias() -= KTT * T + DTT * (T - T_prev) / dt;
356}
MatrixType< _kelvin_vector_size, _number_of_dof > BMatrixType
typename ShapeMatricesType::GlobalDimMatrixType GlobalDimMatrixType
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)
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)
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

References _element, _integration_method, _ip_data, _is_axially_symmetric, _process_data, ProcessLib::LinearBMatrix::computeBMatrix(), MathLib::createZeroedMatrix(), MathLib::createZeroedVector(), MaterialLib::Solids::CreepBGRa, MaterialPropertyLib::density, MathLib::KelvinVector::Invariants< KelvinVectorSize >::deviatoric_projection, displacement_index, displacement_size, MaterialPropertyLib::formEigenTensor(), MaterialPropertyLib::formKelvinVector(), MathLib::KelvinVector::Invariants< KelvinVectorSize >::FrobeniusNorm(), ParameterLib::SpatialPosition::getCoordinates(), NumLib::interpolateCoordinates(), MathLib::KelvinVector::kelvinVectorToSymmetricTensor(), MaterialPropertyLib::VariableArray::mechanical_strain, OGS_FATAL, ParameterLib::SpatialPosition::setElementID(), MaterialPropertyLib::Solid, MaterialPropertyLib::specific_heat_capacity, MaterialPropertyLib::VariableArray::stress, MaterialPropertyLib::VariableArray::temperature, temperature_index, temperature_size, MaterialPropertyLib::thermal_conductivity, and MaterialPropertyLib::thermal_expansivity.

◆ assembleWithJacobianForDeformationEquations()

template<typename ShapeFunction, int DisplacementDim>
void ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, DisplacementDim >::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 )
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_x_prevNodal values of \(x_{prev}\) 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 381 of file ThermoMechanicsFEM-impl.h.

386{
387 auto const local_T =
389
390 auto const local_T_prev =
392
393 auto const local_u =
395 bool const is_u_non_zero = local_u.norm() > 0.0;
396
401
405
406 unsigned const n_integration_points =
407 _integration_method.getNumberOfPoints();
408
411 auto const& medium = _process_data.media_map.getMedium(_element.getID());
412 auto const& solid_phase =
414
415 for (unsigned ip = 0; ip < n_integration_points; ip++)
416 {
417 auto const& w = _ip_data[ip].integration_weight;
418 auto const& N = _ip_data[ip].N;
419 auto const& dNdx = _ip_data[ip].dNdx;
420
422 std::nullopt, _element.getID(),
425 _element, N))};
426
427 auto const x_coord = x_position.getCoordinates().value()[0];
428
429 auto const& B =
434
435 auto& sigma = _ip_data[ip].sigma;
436 auto const& sigma_prev = _ip_data[ip].sigma_prev;
437
438 auto& eps = _ip_data[ip].eps;
439 auto const& eps_prev = _ip_data[ip].eps_prev;
440
441 auto& eps_m = _ip_data[ip].eps_m;
442 auto const& eps_m_prev = _ip_data[ip].eps_m_prev;
443
444 auto& state = _ip_data[ip].material_state_variables;
445
446 const double T_ip = N.dot(local_T); // T at integration point
447 variables.temperature = T_ip;
448 double const dT_ip = T_ip - N.dot(local_T_prev);
449
450 //
451 // displacement equation, displacement part
452 //
453 // For the restart computation, the displacement may not be
454 // reloaded but the initial strains are always available. For such case,
455 // the following computation is skipped.
456 if (is_u_non_zero)
457 {
458 eps.noalias() = B * local_u;
459 }
460
461 // Consider also anisotropic thermal expansion.
465 .property(
468
471
472 eps_m.noalias() = eps_m_prev + eps - eps_prev - dthermal_strain;
473
474 variables_prev.stress
476 sigma_prev);
477 variables_prev.mechanical_strain
479 eps_m_prev);
480 variables_prev.temperature = T_ip;
481 variables.mechanical_strain
483 eps_m);
484
485 auto&& solution = _ip_data[ip].solid_material.integrateStress(
487
488 if (!solution)
489 {
490 OGS_FATAL("Computation of local constitutive relation failed.");
491 }
492
495
496 local_Jac.noalias() += B.transpose() * C * B * w;
497
503
504 for (int i = 0; i < DisplacementDim; ++i)
505 {
508 .noalias() = N;
509 }
510
511 auto const rho_s =
513 .template value<double>(variables, x_position, t, dt);
514
515 auto const& b = _process_data.specific_body_force;
516 local_rhs.noalias() -=
517 (B.transpose() * sigma - N_u.transpose() * rho_s * b) * w;
518 }
519}

References _element, _integration_method, _ip_data, _is_axially_symmetric, _process_data, ProcessLib::LinearBMatrix::computeBMatrix(), MathLib::createZeroedMatrix(), MathLib::createZeroedVector(), MaterialPropertyLib::density, displacement_index, displacement_size, MaterialPropertyLib::formKelvinVector(), ParameterLib::SpatialPosition::getCoordinates(), NumLib::interpolateCoordinates(), MaterialPropertyLib::VariableArray::mechanical_strain, OGS_FATAL, MaterialPropertyLib::Solid, MaterialPropertyLib::VariableArray::stress, MaterialPropertyLib::VariableArray::temperature, temperature_index, and MaterialPropertyLib::thermal_expansivity.

Referenced by assembleWithJacobianForStaggeredScheme().

◆ assembleWithJacobianForHeatConductionEquations()

template<typename ShapeFunction, int DisplacementDim>
void ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, DisplacementDim >::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 )
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_x_prevNodal values of \(x_{prev}\) 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 522 of file ThermoMechanicsFEM-impl.h.

527{
528 auto const local_T =
530
531 auto const local_T_prev =
533
538
542
545
548
549 unsigned const n_integration_points =
550 _integration_method.getNumberOfPoints();
551
552 auto const& medium = _process_data.media_map.getMedium(_element.getID());
553 auto const& solid_phase =
556
557 for (unsigned ip = 0; ip < n_integration_points; ip++)
558 {
559 auto const& w = _ip_data[ip].integration_weight;
560 auto const& N = _ip_data[ip].N;
561 auto const& dNdx = _ip_data[ip].dNdx;
562
564 std::nullopt, _element.getID(),
567 _element, N))};
568
569 const double T_ip = N.dot(local_T); // T at integration point
570 variables.temperature = T_ip;
571
572 auto const rho_s =
574 .template value<double>(variables, x_position, t, dt);
575 auto const c_p =
577 .template value<double>(variables, x_position, t, dt);
578
579 mass.noalias() += N.transpose() * rho_s * c_p * N * w;
580
581 auto const lambda =
583 .property(
585 .value(variables, x_position, t, dt);
586
589
590 laplace.noalias() += dNdx.transpose() * thermal_conductivity * dNdx * w;
591 }
592 local_Jac.noalias() += laplace + mass / dt;
593
594 local_rhs.noalias() -=
596}

References _element, _integration_method, _ip_data, _process_data, MathLib::createZeroedMatrix(), MathLib::createZeroedVector(), MaterialPropertyLib::density, MaterialPropertyLib::formEigenTensor(), NumLib::interpolateCoordinates(), MaterialPropertyLib::Solid, MaterialPropertyLib::specific_heat_capacity, MaterialPropertyLib::VariableArray::temperature, temperature_index, temperature_size, and MaterialPropertyLib::thermal_conductivity.

Referenced by assembleWithJacobianForStaggeredScheme().

◆ assembleWithJacobianForStaggeredScheme()

template<typename ShapeFunction, int DisplacementDim>
void ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, DisplacementDim >::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 )
overridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

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

366{
367 // For the equations with pressure
368 if (process_id == _process_data.heat_conduction_process_id)
369 {
372 return;
373 }
374
375 // For the equations with deformation
378}
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)

References _process_data, assembleWithJacobianForDeformationEquations(), and assembleWithJacobianForHeatConductionEquations().

◆ getEpsilon()

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

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

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

643{
644 constexpr int kelvin_vector_size =
646
649 { return getIntPtEpsilon(0, {}, {}, values); });
650}
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 > transposeInPlace(StoreValuesFunction const &store_values_function)

References getEpsilon(), getIntPtEpsilon(), MathLib::KelvinVector::kelvin_vector_dimensions(), and ProcessLib::transposeInPlace().

Referenced by getEpsilon().

◆ getEpsilonMechanical()

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

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

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

687{
688 constexpr int kelvin_vector_size =
690
693 { return getIntPtEpsilonMechanical(0, {}, {}, values); });
694}
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 getEpsilonMechanical(), getIntPtEpsilonMechanical(), MathLib::KelvinVector::kelvin_vector_dimensions(), and ProcessLib::transposeInPlace().

Referenced by getEpsilonMechanical().

◆ getIntPtEpsilon()

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

◆ getIntPtEpsilonMechanical()

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

◆ getIntPtSigma()

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

◆ getMaterialID()

template<typename ShapeFunction, int DisplacementDim>
int ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, DisplacementDim >::getMaterialID ( ) const
overrideprivatevirtual

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

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

706{
707 return _process_data.material_ids == nullptr
708 ? 0
709 : (*_process_data.material_ids)[_element.getID()];
710}

References _element, _process_data, and getMaterialID().

Referenced by getMaterialID().

◆ getMaterialStateVariablesAt()

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

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

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

717{
718 return *_ip_data[integration_point].material_state_variables;
719}

References _ip_data.

◆ getNumberOfIntegrationPoints()

template<typename ShapeFunction, int DisplacementDim>
unsigned ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, DisplacementDim >::getNumberOfIntegrationPoints ( ) const
overrideprivatevirtual

◆ getNumberOfVectorElementsForDeformation()

template<typename ShapeFunction, int DisplacementDim>
int ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, DisplacementDim >::getNumberOfVectorElementsForDeformation ( ) const
inlineoverridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 248 of file ThermoMechanicsFEM.h.

249 {
250 return displacement_size;
251 }

References displacement_size.

◆ getShapeMatrix()

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

Provides the shape matrix at the given integration point.

Implements NumLib::ExtrapolatableElement.

Definition at line 239 of file ThermoMechanicsFEM.h.

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

References _secondary_data.

◆ getSigma()

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

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

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

610{
611 constexpr int kelvin_vector_size =
613
616 { return getIntPtSigma(0, {}, {}, values); });
617}
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 getIntPtSigma(), MathLib::KelvinVector::kelvin_vector_dimensions(), and ProcessLib::transposeInPlace().

◆ initializeConcrete()

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

Set initial stress from parameter.

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 187 of file ThermoMechanicsFEM.h.

188 {
189 unsigned const n_integration_points =
190 _integration_method.getNumberOfPoints();
191 for (unsigned ip = 0; ip < n_integration_points; ip++)
192 {
193 auto& ip_data = _ip_data[ip];
194
196 std::nullopt, _element.getID(),
200 _element, ip_data.N))};
201
203 if (_process_data.initial_stress != nullptr)
204 {
205 ip_data.sigma =
207 DisplacementDim>((*_process_data.initial_stress)(
209 double>::quiet_NaN() /* time independent */,
210 x_position));
211 }
212
213 double const t = 0; // TODO (naumov) pass t from top
214 ip_data.solid_material.initializeInternalStateVariables(
215 t, x_position, *ip_data.material_state_variables);
216
217 ip_data.pushBackState();
218 }
219 }

References _element, _integration_method, _ip_data, _process_data, NumLib::interpolateCoordinates(), and MathLib::KelvinVector::symmetricTensorToKelvinVector().

◆ postTimestepConcrete()

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

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 225 of file ThermoMechanicsFEM.h.

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

References _integration_method, and _ip_data.

◆ setEpsilon()

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

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

635{
638}
std::size_t setIntegrationPointKelvinVectorData(double const *values, IntegrationPointDataVector &ip_data_vector, MemberType IpData::*const member)

References _ip_data, ProcessLib::ThermoMechanics::IntegrationPointData< BMatricesType, ShapeMatricesType, DisplacementDim >::eps, and ProcessLib::setIntegrationPointKelvinVectorData().

Referenced by setIPDataInitialConditions().

◆ setEpsilonMechanical()

◆ setInitialConditionsConcrete()

template<typename ShapeFunction, int DisplacementDim>
void ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, DisplacementDim >::setInitialConditionsConcrete ( Eigen::VectorXd const local_x,
double const t,
int const process_id )
overridevirtual

◆ setIPDataInitialConditions()

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

Returns number of read integration points.

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

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

113{
114 if (integration_order !=
115 static_cast<int>(_integration_method.getIntegrationOrder()))
116 {
117 OGS_FATAL(
118 "Setting integration point initial conditions; The integration "
119 "order of the local assembler for element {:d} is different from "
120 "the integration order in the initial condition.",
121 _element.getID());
122 }
123
124 if (name == "sigma")
125 {
126 return setSigma(values);
127 }
128 if (name == "epsilon")
129 {
130 return setEpsilon(values);
131 }
132 if (name == "epsilon_m")
133 {
135 }
136
137 return 0;
138}

References _element, _integration_method, OGS_FATAL, setEpsilon(), setEpsilonMechanical(), and setSigma().

◆ setSigma()

Member Data Documentation

◆ _element

◆ _integration_method

◆ _ip_data

◆ _is_axially_symmetric

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

◆ _process_data

◆ _secondary_data

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

Definition at line 354 of file ThermoMechanicsFEM.h.

Referenced by ThermoMechanicsLocalAssembler(), and getShapeMatrix().

◆ displacement_index

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

◆ displacement_size

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

Definition at line 360 of file ThermoMechanicsFEM.h.

Referenced by assembleWithJacobian(), assembleWithJacobianForDeformationEquations(), and getNumberOfVectorElementsForDeformation().

◆ KelvinVectorSize

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

Definition at line 133 of file ThermoMechanicsFEM.h.

◆ temperature_index

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

◆ temperature_size

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

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