|  | OGS
    | 
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>
| 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< VectorSegment > | getVectorDeformationSegment () 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 | 
| using ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, DisplacementDim >::BMatricesType = BMatrixPolicyType<ShapeFunction, DisplacementDim> | 
Definition at line 138 of file ThermoMechanicsFEM.h.
| using ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, DisplacementDim >::GlobalDimMatrixType = typename ShapeMatricesType::GlobalDimMatrixType | 
Definition at line 137 of file ThermoMechanicsFEM.h.
| using ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, DisplacementDim >::Invariants = MathLib::KelvinVector::Invariants<KelvinVectorSize> | 
Definition at line 142 of file ThermoMechanicsFEM.h.
| using ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, DisplacementDim >::IpData | 
Definition at line 150 of file ThermoMechanicsFEM.h.
| using ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, DisplacementDim >::JacobianMatrix | 
Definition at line 147 of file ThermoMechanicsFEM.h.
| using ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, DisplacementDim >::NodalForceVectorType = typename BMatricesType::NodalForceVectorType | 
Definition at line 144 of file ThermoMechanicsFEM.h.
| using ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, DisplacementDim >::RhsVector | 
Definition at line 145 of file ThermoMechanicsFEM.h.
| using ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, DisplacementDim >::ShapeMatrices = typename ShapeMatricesType::ShapeMatrices | 
Definition at line 136 of file ThermoMechanicsFEM.h.
| using ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, DisplacementDim >::ShapeMatricesType | 
Definition at line 131 of file ThermoMechanicsFEM.h.
| 
 | delete | 
| 
 | delete | 
| 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.
References ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, DisplacementDim >::_element, ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, DisplacementDim >::_integration_method, ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, DisplacementDim >::_ip_data, ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, DisplacementDim >::_process_data, ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, DisplacementDim >::_secondary_data, MeshLib::Element::getID(), NumLib::GenericIntegrationMethod::getNumberOfPoints(), MathLib::WeightedPoint::getWeight(), NumLib::GenericIntegrationMethod::getWeightedPoint(), NumLib::initShapeMatrices(), MathLib::KelvinVector::kelvin_vector_dimensions(), ProcessLib::ThermoMechanics::SecondaryData< ShapeMatrixType >::N, MaterialLib::Solids::selectSolidConstitutiveRelation(), and ParameterLib::SpatialPosition::setElementID().
| 
 | inlineoverridevirtual | 
Reimplemented from ProcessLib::LocalAssemblerInterface.
Definition at line 170 of file ThermoMechanicsFEM.h.
References OGS_FATAL.
| 
 | overridevirtual | 
Reimplemented from ProcessLib::LocalAssemblerInterface.
Definition at line 149 of file ThermoMechanicsFEM-impl.h.
References ProcessLib::LinearBMatrix::computeBMatrix(), MathLib::createZeroedMatrix(), MathLib::createZeroedVector(), MaterialLib::Solids::CreepBGRa, MaterialPropertyLib::formEigenTensor(), MaterialPropertyLib::formKelvinVector(), 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, MaterialPropertyLib::thermal_conductivity, and MaterialPropertyLib::thermal_expansivity.
| 
 | 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.
| t | Time | 
| dt | Time increment | 
| local_x | Nodal values of \(x\) of all processes of an element. | 
| local_x_prev | Nodal values of \(x_{prev}\) of all processes of an element. | 
| local_b_data | Right hand side vector of an element. | 
| local_Jac_data | Element Jacobian matrix for the Newton-Raphson method. | 
Definition at line 388 of file ThermoMechanicsFEM-impl.h.
References ProcessLib::LinearBMatrix::computeBMatrix(), MathLib::createZeroedMatrix(), MathLib::createZeroedVector(), MaterialPropertyLib::formKelvinVector(), NumLib::interpolateCoordinates(), MaterialPropertyLib::VariableArray::mechanical_strain, OGS_FATAL, MaterialPropertyLib::VariableArray::stress, MaterialPropertyLib::VariableArray::temperature, and MaterialPropertyLib::thermal_expansivity.
| 
 | 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.
| t | Time | 
| dt | Time increment | 
| local_x | Nodal values of \(x\) of all processes of an element. | 
| local_x_prev | Nodal values of \(x_{prev}\) of all processes of an element. | 
| local_b_data | Right hand side vector of an element. | 
| local_Jac_data | Element Jacobian matrix for the Newton-Raphson method. | 
Definition at line 528 of file ThermoMechanicsFEM-impl.h.
References MathLib::createZeroedMatrix(), MathLib::createZeroedVector(), MaterialPropertyLib::formEigenTensor(), NumLib::interpolateCoordinates(), MaterialPropertyLib::VariableArray::temperature, and MaterialPropertyLib::thermal_conductivity.
| 
 | overridevirtual | 
Reimplemented from ProcessLib::LocalAssemblerInterface.
Definition at line 366 of file ThermoMechanicsFEM-impl.h.
| 
 | overrideprivatevirtual | 
Implements ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssemblerInterface< DisplacementDim >.
Definition at line 647 of file ThermoMechanicsFEM-impl.h.
References MathLib::KelvinVector::kelvin_vector_dimensions(), and ProcessLib::transposeInPlace().
| 
 | overrideprivatevirtual | 
Implements ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssemblerInterface< DisplacementDim >.
Definition at line 691 of file ThermoMechanicsFEM-impl.h.
References MathLib::KelvinVector::kelvin_vector_dimensions(), and ProcessLib::transposeInPlace().
| 
 | overrideprivatevirtual | 
Implements ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssemblerInterface< DisplacementDim >.
Definition at line 659 of file ThermoMechanicsFEM-impl.h.
References ProcessLib::getIntegrationPointKelvinVectorData().
| 
 | overrideprivatevirtual | 
Implements ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssemblerInterface< DisplacementDim >.
Definition at line 671 of file ThermoMechanicsFEM-impl.h.
References ProcessLib::getIntegrationPointKelvinVectorData().
| 
 | overrideprivatevirtual | 
Implements ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssemblerInterface< DisplacementDim >.
Definition at line 626 of file ThermoMechanicsFEM-impl.h.
References ProcessLib::getIntegrationPointKelvinVectorData().
| 
 | overrideprivatevirtual | 
Implements ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssemblerInterface< DisplacementDim >.
Definition at line 710 of file ThermoMechanicsFEM-impl.h.
| 
 | overrideprivatevirtual | 
Implements ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssemblerInterface< DisplacementDim >.
Definition at line 720 of file ThermoMechanicsFEM-impl.h.
| 
 | overrideprivatevirtual | 
Implements ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssemblerInterface< DisplacementDim >.
Definition at line 703 of file ThermoMechanicsFEM-impl.h.
| 
 | inlineoverridevirtual | 
Provides the shape matrix at the given integration point.
Implements NumLib::ExtrapolatableElement.
Definition at line 246 of file ThermoMechanicsFEM.h.
References ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, DisplacementDim >::_secondary_data, and ProcessLib::ThermoMechanics::SecondaryData< ShapeMatrixType >::N.
| 
 | overrideprivatevirtual | 
Implements ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssemblerInterface< DisplacementDim >.
Definition at line 614 of file ThermoMechanicsFEM-impl.h.
References MathLib::KelvinVector::kelvin_vector_dimensions(), and ProcessLib::transposeInPlace().
| 
 | inlineoverridevirtual | 
Set initial stress from parameter.
Reimplemented from ProcessLib::LocalAssemblerInterface.
Definition at line 194 of file ThermoMechanicsFEM.h.
References ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, DisplacementDim >::_element, ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, DisplacementDim >::_integration_method, ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, DisplacementDim >::_ip_data, ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, DisplacementDim >::_process_data, MeshLib::Element::getID(), NumLib::GenericIntegrationMethod::getNumberOfPoints(), NumLib::interpolateCoordinates(), and MathLib::KelvinVector::symmetricTensorToKelvinVector().
| 
 | inlineoverridevirtual | 
Reimplemented from ProcessLib::LocalAssemblerInterface.
Definition at line 232 of file ThermoMechanicsFEM.h.
References ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, DisplacementDim >::_integration_method, ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, DisplacementDim >::_ip_data, and NumLib::GenericIntegrationMethod::getNumberOfPoints().
| 
 | private | 
Definition at line 638 of file ThermoMechanicsFEM-impl.h.
References ProcessLib::setIntegrationPointKelvinVectorData().
| 
 | private | 
Definition at line 684 of file ThermoMechanicsFEM-impl.h.
References ProcessLib::setIntegrationPointKelvinVectorData().
| 
 | overridevirtual | 
Reimplemented from ProcessLib::LocalAssemblerInterface.
Definition at line 83 of file ThermoMechanicsFEM-impl.h.
References ProcessLib::LinearBMatrix::computeBMatrix(), and NumLib::interpolateCoordinates().
| 
 | overridevirtual | 
Returns number of read integration points.
Implements ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssemblerInterface< DisplacementDim >.
Definition at line 117 of file ThermoMechanicsFEM-impl.h.
References OGS_FATAL.
| 
 | private | 
Definition at line 605 of file ThermoMechanicsFEM-impl.h.
References ProcessLib::setIntegrationPointKelvinVectorData().
| 
 | private | 
Definition at line 355 of file ThermoMechanicsFEM.h.
Referenced by ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, DisplacementDim >::ThermoMechanicsLocalAssembler(), and ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, DisplacementDim >::initializeConcrete().
| 
 | private | 
Definition at line 354 of file ThermoMechanicsFEM.h.
Referenced by ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, DisplacementDim >::ThermoMechanicsLocalAssembler(), ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, DisplacementDim >::initializeConcrete(), and ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, DisplacementDim >::postTimestepConcrete().
| 
 | private | 
Definition at line 352 of file ThermoMechanicsFEM.h.
Referenced by ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, DisplacementDim >::ThermoMechanicsLocalAssembler(), ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, DisplacementDim >::initializeConcrete(), and ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, DisplacementDim >::postTimestepConcrete().
| 
 | private | 
Definition at line 357 of file ThermoMechanicsFEM.h.
| 
 | private | 
Definition at line 350 of file ThermoMechanicsFEM.h.
Referenced by ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, DisplacementDim >::ThermoMechanicsLocalAssembler(), and ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, DisplacementDim >::initializeConcrete().
| 
 | private | 
Definition at line 356 of file ThermoMechanicsFEM.h.
Referenced by ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, DisplacementDim >::ThermoMechanicsLocalAssembler(), and ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, DisplacementDim >::getShapeMatrix().
| 
 | staticprivate | 
Definition at line 361 of file ThermoMechanicsFEM.h.
| 
 | staticprivate | 
Definition at line 362 of file ThermoMechanicsFEM.h.
| 
 | static | 
Definition at line 140 of file ThermoMechanicsFEM.h.
| 
 | staticprivate | 
Definition at line 359 of file ThermoMechanicsFEM.h.
| 
 | staticprivate | 
Definition at line 360 of file ThermoMechanicsFEM.h.