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 127 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.
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.
virtual std::optional< VectorSegmentgetVectorDeformationSegment () const
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 138 of file ThermoMechanicsFEM.h.

◆ GlobalDimMatrixType

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

Definition at line 137 of file ThermoMechanicsFEM.h.

◆ Invariants

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

Definition at line 142 of file ThermoMechanicsFEM.h.

◆ IpData

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

◆ NodalForceVectorType

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

Definition at line 144 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 145 of file ThermoMechanicsFEM.h.

◆ ShapeMatrices

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

Definition at line 136 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 131 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 27 of file ThermoMechanicsFEM-impl.h.

36 _element(e),
38{
39 unsigned const n_integration_points =
40 _integration_method.getNumberOfPoints();
41
44
45 auto const shape_matrices =
49
51 _process_data.solid_materials, _process_data.material_ids, e.getID());
52
53 for (unsigned ip = 0; ip < n_integration_points; ip++)
54 {
55 _ip_data.emplace_back(solid_material);
56 auto& ip_data = _ip_data[ip];
57 ip_data.integration_weight =
58 _integration_method.getWeightedPoint(ip).getWeight() *
59 shape_matrices[ip].integralMeasure * shape_matrices[ip].detJ;
60
61 static const int kelvin_vector_size =
63 // Initialize current time step values
64 ip_data.sigma.setZero(kelvin_vector_size);
65 ip_data.eps.setZero(kelvin_vector_size);
66
67 // Previous time step values are not initialized and are set later.
68 ip_data.sigma_prev.resize(kelvin_vector_size);
69 ip_data.eps_prev.resize(kelvin_vector_size);
70
71 ip_data.eps_m.setZero(kelvin_vector_size);
72 ip_data.eps_m_prev.setZero(kelvin_vector_size);
74 x_position.setElementID(_element.getID());
76 ip_data.dNdx = shape_matrices[ip].dNdx;
77
79 }
80}
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 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, 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 149 of file ThermoMechanicsFEM-impl.h.

155{
156 auto const local_matrix_size = local_x.size();
158
162
166 bool const is_u_non_zero = u.norm() > 0.0;
167
171
174
177
180 KuT;
182
185
188
189 unsigned const n_integration_points =
190 _integration_method.getNumberOfPoints();
191
195 x_position.setElementID(_element.getID());
196
197 auto const& medium = _process_data.media_map.getMedium(_element.getID());
198 auto const& solid_phase = medium->phase("Solid");
199
200 for (unsigned ip = 0; ip < n_integration_points; ip++)
201 {
202 auto const& w = _ip_data[ip].integration_weight;
203 auto const& N = _ip_data[ip].N;
204 auto const& dNdx = _ip_data[ip].dNdx;
205
207 std::nullopt, _element.getID(),
210 _element, N))};
211
212 auto const x_coord = x_position.getCoordinates().value()[0];
213 auto const& B =
218
219 auto& sigma = _ip_data[ip].sigma;
220 auto const& sigma_prev = _ip_data[ip].sigma_prev;
221
222 auto& eps = _ip_data[ip].eps;
223 auto const& eps_prev = _ip_data[ip].eps_prev;
224
225 auto& eps_m = _ip_data[ip].eps_m;
226 auto const& eps_m_prev = _ip_data[ip].eps_m_prev;
227
228 auto& state = _ip_data[ip].material_state_variables;
229
230 const double T_ip = N.dot(T); // T at integration point
231 double const T_prev_ip = N.dot(T_prev);
232
233 // Consider also anisotropic thermal expansion.
237 .property(
240
244
245 //
246 // displacement equation, displacement part
247 //
248 // For the restart computation, the displacement may not be
249 // reloaded but the initial strains are always available. For such case,
250 // the following computation is skipped.
251 if (is_u_non_zero)
252 {
253 eps.noalias() = B * u;
254 }
255
256 eps_m.noalias() = eps_m_prev + eps - eps_prev - dthermal_strain;
257
258 variables_prev.stress
260 sigma_prev);
261 variables_prev.mechanical_strain
263 eps_m_prev);
264 variables_prev.temperature = T_ip;
265 variables.mechanical_strain
267 eps_m);
268 variables.temperature = T_ip;
269
270 auto&& solution = _ip_data[ip].solid_material.integrateStress(
272
273 if (!solution)
274 {
275 OGS_FATAL("Computation of local constitutive relation failed.");
276 }
277
280
284 .noalias() += B.transpose() * C * B * w;
285
291
292 for (int i = 0; i < DisplacementDim; ++i)
293 {
296 .noalias() = N;
297 }
298
299 auto const rho_s =
301 .template value<double>(variables, x_position, t, dt);
302
303 auto const& b = _process_data.specific_body_force;
305 .noalias() -=
306 (B.transpose() * sigma - N_u.transpose() * rho_s * b) * w;
307
308 //
309 // displacement equation, temperature part
310 // The computation of KuT can be ignored.
311 auto const alpha_T_tensor =
314 KuT.noalias() += B.transpose() * (C * alpha_T_tensor) * N * w;
315
316 if (_ip_data[ip].solid_material.getConstitutiveModel() ==
318 {
320 double const norm_s = Invariants::FrobeniusNorm(s);
321 const double creep_coefficient =
322 _ip_data[ip].solid_material.getTemperatureRelatedCoefficient(
324 KuT.noalias() += creep_coefficient * B.transpose() * s * N * w;
325 }
326
327 //
328 // temperature equation, temperature part;
329 //
330 auto const lambda =
332 .property(
334 .value(variables, x_position, t, dt);
335
338
339 KTT.noalias() += dNdx.transpose() * thermal_conductivity * dNdx * w;
340
341 auto const c =
343 .property(
345 .template value<double>(variables, x_position, t, dt);
346 DTT.noalias() += N.transpose() * rho_s * c * N * w;
347 }
348
349 // temperature equation, temperature part
353 .noalias() += KTT + DTT / dt;
354
355 // displacement equation, temperature part
359 .noalias() -= KuT;
360
362 .noalias() -= KTT * T + DTT * (T - T_prev) / dt;
363}
MatrixType< _kelvin_vector_size, _number_of_dof > BMatrixType
typename ShapeMatricesType::GlobalDimMatrixType GlobalDimMatrixType
MathLib::KelvinVector::KelvinVectorType< GlobalDim > formKelvinVector(MaterialPropertyLib::PropertyDataType const &values)
A function to form a Kelvin vector from strain or stress alike property like thermal expansivity for ...
Eigen::Matrix< double, GlobalDim, GlobalDim > formEigenTensor(MaterialPropertyLib::PropertyDataType const &values)
Eigen::Matrix< double, 4, 1 > kelvinVectorToSymmetricTensor(Eigen::Matrix< double, 4, 1, Eigen::ColMajor, 4, 1 > const &v)
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::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 388 of file ThermoMechanicsFEM-impl.h.

393{
394 auto const local_T =
396
397 auto const local_T_prev =
399
400 auto const local_u =
402 bool const is_u_non_zero = local_u.norm() > 0.0;
403
408
412
413 unsigned const n_integration_points =
414 _integration_method.getNumberOfPoints();
415
418 auto const& medium = _process_data.media_map.getMedium(_element.getID());
419 auto const& solid_phase = medium->phase("Solid");
420
421 for (unsigned ip = 0; ip < n_integration_points; ip++)
422 {
423 auto const& w = _ip_data[ip].integration_weight;
424 auto const& N = _ip_data[ip].N;
425 auto const& dNdx = _ip_data[ip].dNdx;
426
428 std::nullopt, _element.getID(),
431 _element, N))};
432
433 auto const x_coord = x_position.getCoordinates().value()[0];
434
435 auto const& B =
440
441 auto& sigma = _ip_data[ip].sigma;
442 auto const& sigma_prev = _ip_data[ip].sigma_prev;
443
444 auto& eps = _ip_data[ip].eps;
445 auto const& eps_prev = _ip_data[ip].eps_prev;
446
447 auto& eps_m = _ip_data[ip].eps_m;
448 auto const& eps_m_prev = _ip_data[ip].eps_m_prev;
449
450 auto& state = _ip_data[ip].material_state_variables;
451
452 const double T_ip = N.dot(local_T); // T at integration point
453 variables.temperature = T_ip;
454 double const dT_ip = T_ip - N.dot(local_T_prev);
455
456 //
457 // displacement equation, displacement part
458 //
459 // For the restart computation, the displacement may not be
460 // reloaded but the initial strains are always available. For such case,
461 // the following computation is skipped.
462 if (is_u_non_zero)
463 {
464 eps.noalias() = B * local_u;
465 }
466
467 // Consider also anisotropic thermal expansion.
471 .property(
474
477
478 eps_m.noalias() = eps_m_prev + eps - eps_prev - dthermal_strain;
479
480 variables_prev.stress
482 sigma_prev);
483 variables_prev.mechanical_strain
485 eps_m_prev);
486 variables_prev.temperature = T_ip;
487 variables.mechanical_strain
489 eps_m);
490
491 auto&& solution = _ip_data[ip].solid_material.integrateStress(
493
494 if (!solution)
495 {
496 OGS_FATAL("Computation of local constitutive relation failed.");
497 }
498
501
502 local_Jac.noalias() += B.transpose() * C * B * w;
503
509
510 for (int i = 0; i < DisplacementDim; ++i)
511 {
514 .noalias() = N;
515 }
516
517 auto const rho_s =
519 .template value<double>(variables, x_position, t, dt);
520
521 auto const& b = _process_data.specific_body_force;
522 local_rhs.noalias() -=
523 (B.transpose() * sigma - N_u.transpose() * rho_s * b) * w;
524 }
525}

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::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 528 of file ThermoMechanicsFEM-impl.h.

533{
534 auto const local_T =
536
537 auto const local_T_prev =
539
544
548
551
554
555 unsigned const n_integration_points =
556 _integration_method.getNumberOfPoints();
557
558 auto const& medium = _process_data.media_map.getMedium(_element.getID());
559 auto const& solid_phase = medium->phase("Solid");
561
562 for (unsigned ip = 0; ip < n_integration_points; ip++)
563 {
564 auto const& w = _ip_data[ip].integration_weight;
565 auto const& N = _ip_data[ip].N;
566 auto const& dNdx = _ip_data[ip].dNdx;
567
569 std::nullopt, _element.getID(),
572 _element, N))};
573
574 const double T_ip = N.dot(local_T); // T at integration point
575 variables.temperature = T_ip;
576
577 auto const rho_s =
579 .template value<double>(variables, x_position, t, dt);
580 auto const c_p =
582 .template value<double>(variables, x_position, t, dt);
583
584 mass.noalias() += N.transpose() * rho_s * c_p * N * w;
585
586 auto const lambda =
588 .property(
590 .value(variables, x_position, t, dt);
591
594
595 laplace.noalias() += dNdx.transpose() * thermal_conductivity * dNdx * w;
596 }
597 local_Jac.noalias() += laplace + mass / dt;
598
599 local_rhs.noalias() -=
601}

References _element, _integration_method, _ip_data, _process_data, MathLib::createZeroedMatrix(), MathLib::createZeroedVector(), MaterialPropertyLib::density, MaterialPropertyLib::formEigenTensor(), NumLib::interpolateCoordinates(), 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 366 of file ThermoMechanicsFEM-impl.h.

373{
374 // For the equations with pressure
375 if (process_id == _process_data.heat_conduction_process_id)
376 {
379 return;
380 }
381
382 // For the equations with deformation
385}
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 647 of file ThermoMechanicsFEM-impl.h.

648{
649 constexpr int kelvin_vector_size =
651
654 { return getIntPtEpsilon(0, {}, {}, values); });
655}
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 691 of file ThermoMechanicsFEM-impl.h.

692{
693 constexpr int kelvin_vector_size =
695
698 { return getIntPtEpsilonMechanical(0, {}, {}, values); });
699}
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 710 of file ThermoMechanicsFEM-impl.h.

711{
712 return _process_data.material_ids == nullptr
713 ? 0
714 : (*_process_data.material_ids)[_element.getID()];
715}

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 720 of file ThermoMechanicsFEM-impl.h.

722{
723 return *_ip_data[integration_point].material_state_variables;
724}

References _ip_data.

◆ getNumberOfIntegrationPoints()

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

◆ 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 246 of file ThermoMechanicsFEM.h.

248 {
249 auto const& N = _secondary_data.N[integration_point];
250
251 // assumes N is stored contiguously in memory
252 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
253 }

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 614 of file ThermoMechanicsFEM-impl.h.

615{
616 constexpr int kelvin_vector_size =
618
621 { return getIntPtSigma(0, {}, {}, values); });
622}
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 194 of file ThermoMechanicsFEM.h.

195 {
196 unsigned const n_integration_points =
197 _integration_method.getNumberOfPoints();
198 for (unsigned ip = 0; ip < n_integration_points; ip++)
199 {
200 auto& ip_data = _ip_data[ip];
201
203 std::nullopt, _element.getID(),
207 _element, ip_data.N))};
208
210 if (_process_data.initial_stress != nullptr)
211 {
212 ip_data.sigma =
214 DisplacementDim>((*_process_data.initial_stress)(
216 double>::quiet_NaN() /* time independent */,
217 x_position));
218 }
219
220 double const t = 0; // TODO (naumov) pass t from top
221 ip_data.solid_material.initializeInternalStateVariables(
222 t, x_position, *ip_data.material_state_variables);
223
224 ip_data.pushBackState();
225 }
226 }

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 232 of file ThermoMechanicsFEM.h.

236 {
237 unsigned const n_integration_points =
238 _integration_method.getNumberOfPoints();
239
240 for (unsigned ip = 0; ip < n_integration_points; ip++)
241 {
242 _ip_data[ip].pushBackState();
243 }
244 }

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 638 of file ThermoMechanicsFEM-impl.h.

640{
643}
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

Reimplemented from ProcessLib::LocalAssemblerInterface.

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

87{
88 unsigned const n_integration_points =
89 this->_integration_method.getNumberOfPoints();
90 auto const local_u =
92
93 for (unsigned ip = 0; ip < n_integration_points; ip++)
94 {
95 auto const& N = _ip_data[ip].N;
96 auto const& dNdx = _ip_data[ip].dNdx;
97
99 std::nullopt, _element.getID(),
102 _element, N))};
103
104 auto const x_coord = x_position.getCoordinates().value()[0];
105
106 auto const& B =
111
112 _ip_data[ip].eps.noalias() = B * local_u;
113 }
114}

References _element, _integration_method, _ip_data, _is_axially_symmetric, ProcessLib::LinearBMatrix::computeBMatrix(), displacement_index, ParameterLib::SpatialPosition::getCoordinates(), and NumLib::interpolateCoordinates().

◆ 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 117 of file ThermoMechanicsFEM-impl.h.

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

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 356 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 362 of file ThermoMechanicsFEM.h.

Referenced by assembleWithJacobian(), and assembleWithJacobianForDeformationEquations().

◆ KelvinVectorSize

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

Definition at line 140 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: