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 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.
 
- 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
Initial value:
IntegrationPointData<BMatricesType, ShapeMatricesType, DisplacementDim>

Definition at line 150 of file ThermoMechanicsFEM.h.

◆ 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

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.

34 : _process_data(process_data),
35 _integration_method(integration_method),
36 _element(e),
37 _is_axially_symmetric(is_axially_symmetric)
38{
39 unsigned const n_integration_points =
41
42 _ip_data.reserve(n_integration_points);
43 _secondary_data.N.resize(n_integration_points);
44
45 auto const shape_matrices =
47 DisplacementDim>(e, is_axially_symmetric,
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 =
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());
75 ip_data.N = shape_matrices[ip].N;
76 ip_data.dNdx = shape_matrices[ip].dNdx;
77
78 _secondary_data.N[ip] = shape_matrices[ip].N;
79 }
80}
double getWeight() const
std::size_t getID() const
Returns the ID of the element.
Definition Element.h:89
MathLib::WeightedPoint const & getWeightedPoint(unsigned const igp) const
void setElementID(std::size_t element_id)
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.
std::vector< typename ShapeMatricesType::ShapeMatrices, Eigen::aligned_allocator< typename ShapeMatricesType::ShapeMatrices > > initShapeMatrices(MeshLib::Element const &e, bool const is_axially_symmetric, IntegrationMethod const &integration_method)
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N

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().

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

121{
122 auto const local_matrix_size = local_x.size();
123 assert(local_matrix_size == temperature_size + displacement_size);
124
125 auto T = Eigen::Map<typename ShapeMatricesType::template VectorType<
126 temperature_size> const>(local_x.data() + temperature_index,
128
129 auto u = Eigen::Map<typename ShapeMatricesType::template VectorType<
130 displacement_size> const>(local_x.data() + displacement_index,
132 bool const is_u_non_zero = u.norm() > 0.0;
133
134 auto T_prev = Eigen::Map<typename ShapeMatricesType::template VectorType<
135 temperature_size> const>(local_x_prev.data() + temperature_index,
137
139 local_Jac_data, local_matrix_size, local_matrix_size);
140
141 auto local_rhs = MathLib::createZeroedVector<RhsVector>(local_rhs_data,
142 local_matrix_size);
143
144 typename ShapeMatricesType::template MatrixType<displacement_size,
146 KuT;
148
151
154
155 unsigned const n_integration_points =
157
158 MPL::VariableArray variables;
159 MPL::VariableArray variables_prev;
161 x_position.setElementID(_element.getID());
162
163 auto const& medium = _process_data.media_map.getMedium(_element.getID());
164 auto const& solid_phase = medium->phase("Solid");
165
166 for (unsigned ip = 0; ip < n_integration_points; ip++)
167 {
168 x_position.setIntegrationPoint(ip);
169 auto const& w = _ip_data[ip].integration_weight;
170 auto const& N = _ip_data[ip].N;
171 auto const& dNdx = _ip_data[ip].dNdx;
172
173 auto const x_coord =
175 _element, N);
176 auto const& B =
177 LinearBMatrix::computeBMatrix<DisplacementDim,
178 ShapeFunction::NPOINTS,
180 dNdx, N, x_coord, _is_axially_symmetric);
181
182 auto& sigma = _ip_data[ip].sigma;
183 auto const& sigma_prev = _ip_data[ip].sigma_prev;
184
185 auto& eps = _ip_data[ip].eps;
186 auto const& eps_prev = _ip_data[ip].eps_prev;
187
188 auto& eps_m = _ip_data[ip].eps_m;
189 auto const& eps_m_prev = _ip_data[ip].eps_m_prev;
190
191 auto& state = _ip_data[ip].material_state_variables;
192
193 const double T_ip = N.dot(T); // T at integration point
194 double const T_prev_ip = N.dot(T_prev);
195
196 // Consider also anisotropic thermal expansion.
197 auto const solid_linear_thermal_expansivity_vector =
199 solid_phase
200 .property(
202 .value(variables, x_position, t, dt));
203
205 dthermal_strain =
206 solid_linear_thermal_expansivity_vector * (T_ip - T_prev_ip);
207
208 //
209 // displacement equation, displacement part
210 //
211 // For the restart computation, the displacement may not be
212 // reloaded but the initial strains are always available. For such case,
213 // the following computation is skipped.
214 if (is_u_non_zero)
215 {
216 eps.noalias() = B * u;
217 }
218
219 eps_m.noalias() = eps_m_prev + eps - eps_prev - dthermal_strain;
220
221 variables_prev.stress
223 sigma_prev);
224 variables_prev.mechanical_strain
226 eps_m_prev);
227 variables_prev.temperature = T_ip;
228 variables.mechanical_strain
230 eps_m);
231 variables.temperature = T_ip;
232
233 auto&& solution = _ip_data[ip].solid_material.integrateStress(
234 variables_prev, variables, t, x_position, dt, *state);
235
236 if (!solution)
237 {
238 OGS_FATAL("Computation of local constitutive relation failed.");
239 }
240
242 std::tie(sigma, state, C) = std::move(*solution);
243
244 local_Jac
245 .template block<displacement_size, displacement_size>(
247 .noalias() += B.transpose() * C * B * w;
248
249 typename ShapeMatricesType::template MatrixType<DisplacementDim,
251 N_u = ShapeMatricesType::template MatrixType<
252 DisplacementDim, displacement_size>::Zero(DisplacementDim,
254
255 for (int i = 0; i < DisplacementDim; ++i)
256 {
257 N_u.template block<1, displacement_size / DisplacementDim>(
258 i, i * displacement_size / DisplacementDim)
259 .noalias() = N;
260 }
261
262 auto const rho_s =
263 solid_phase.property(MPL::PropertyType::density)
264 .template value<double>(variables, x_position, t, dt);
265
266 auto const& b = _process_data.specific_body_force;
267 local_rhs.template block<displacement_size, 1>(displacement_index, 0)
268 .noalias() -=
269 (B.transpose() * sigma - N_u.transpose() * rho_s * b) * w;
270
271 //
272 // displacement equation, temperature part
273 // The computation of KuT can be ignored.
274 auto const alpha_T_tensor =
276 solid_linear_thermal_expansivity_vector);
277 KuT.noalias() += B.transpose() * (C * alpha_T_tensor) * N * w;
278
279 if (_ip_data[ip].solid_material.getConstitutiveModel() ==
281 {
282 auto const s = Invariants::deviatoric_projection * sigma;
283 double const norm_s = Invariants::FrobeniusNorm(s);
284 const double creep_coefficient =
285 _ip_data[ip].solid_material.getTemperatureRelatedCoefficient(
286 t, dt, x_position, T_ip, norm_s);
287 KuT.noalias() += creep_coefficient * B.transpose() * s * N * w;
288 }
289
290 //
291 // temperature equation, temperature part;
292 //
293 auto const lambda =
294 solid_phase
295 .property(
297 .value(variables, x_position, t, dt);
298
301
302 KTT.noalias() += dNdx.transpose() * thermal_conductivity * dNdx * w;
303
304 auto const c =
305 solid_phase
306 .property(
308 .template value<double>(variables, x_position, t, dt);
309 DTT.noalias() += N.transpose() * rho_s * c * N * w;
310 }
311
312 // temperature equation, temperature part
313 local_Jac
314 .template block<temperature_size, temperature_size>(temperature_index,
316 .noalias() += KTT + DTT / dt;
317
318 // displacement equation, temperature part
319 local_Jac
320 .template block<displacement_size, temperature_size>(displacement_index,
322 .noalias() -= KuT;
323
324 local_rhs.template block<temperature_size, 1>(temperature_index, 0)
325 .noalias() -= KTT * T + DTT * (T - T_prev) / dt;
326}
void setIntegrationPoint(unsigned integration_point)
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::Matrix< double, kelvin_vector_dimensions(DisplacementDim), 1, Eigen::ColMajor > KelvinVectorType
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), kelvin_vector_dimensions(DisplacementDim), Eigen::RowMajor > KelvinMatrixType
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)
double interpolateXCoordinate(MeshLib::Element const &e, typename ShapeMatricesType::ShapeMatrices::ShapeType const &N)
BMatrixType computeBMatrix(DNDX_Type const &dNdx, N_Type const &N, const double radius, const bool is_axially_symmetric)
Fills a B-matrix based on given shape function dN/dx values.
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 ProcessLib::LinearBMatrix::computeBMatrix(), MathLib::createZeroedMatrix(), MathLib::createZeroedVector(), MaterialLib::Solids::CreepBGRa, MaterialPropertyLib::formEigenTensor(), MaterialPropertyLib::formKelvinVector(), NumLib::interpolateXCoordinate(), MathLib::KelvinVector::kelvinVectorToSymmetricTensor(), MaterialPropertyLib::VariableArray::mechanical_strain, OGS_FATAL, ParameterLib::SpatialPosition::setElementID(), ParameterLib::SpatialPosition::setIntegrationPoint(), MaterialPropertyLib::specific_heat_capacity, MaterialPropertyLib::VariableArray::stress, MaterialPropertyLib::VariableArray::temperature, 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 351 of file ThermoMechanicsFEM-impl.h.

356{
357 auto const local_T =
358 local_x.template segment<temperature_size>(temperature_index);
359
360 auto const local_T_prev =
361 local_x_prev.template segment<temperature_size>(temperature_index);
362
363 auto const local_u =
364 local_x.template segment<displacement_size>(displacement_index);
365 bool const is_u_non_zero = local_u.norm() > 0.0;
366
367 auto local_Jac = MathLib::createZeroedMatrix<
368 typename ShapeMatricesType::template MatrixType<displacement_size,
370 local_Jac_data, displacement_size, displacement_size);
371
372 auto local_rhs = MathLib::createZeroedVector<
373 typename ShapeMatricesType::template VectorType<displacement_size>>(
374 local_b_data, displacement_size);
375
376 unsigned const n_integration_points =
378
379 MPL::VariableArray variables;
380 MPL::VariableArray variables_prev;
382 x_position.setElementID(_element.getID());
383 auto const& medium = _process_data.media_map.getMedium(_element.getID());
384 auto const& solid_phase = medium->phase("Solid");
385
386 for (unsigned ip = 0; ip < n_integration_points; ip++)
387 {
388 x_position.setIntegrationPoint(ip);
389 auto const& w = _ip_data[ip].integration_weight;
390 auto const& N = _ip_data[ip].N;
391 auto const& dNdx = _ip_data[ip].dNdx;
392
393 auto const x_coord =
395 _element, N);
396 auto const& B =
397 LinearBMatrix::computeBMatrix<DisplacementDim,
398 ShapeFunction::NPOINTS,
400 dNdx, N, x_coord, _is_axially_symmetric);
401
402 auto& sigma = _ip_data[ip].sigma;
403 auto const& sigma_prev = _ip_data[ip].sigma_prev;
404
405 auto& eps = _ip_data[ip].eps;
406 auto const& eps_prev = _ip_data[ip].eps_prev;
407
408 auto& eps_m = _ip_data[ip].eps_m;
409 auto const& eps_m_prev = _ip_data[ip].eps_m_prev;
410
411 auto& state = _ip_data[ip].material_state_variables;
412
413 const double T_ip = N.dot(local_T); // T at integration point
414 variables.temperature = T_ip;
415 double const dT_ip = T_ip - N.dot(local_T_prev);
416
417 //
418 // displacement equation, displacement part
419 //
420 // For the restart computation, the displacement may not be
421 // reloaded but the initial strains are always available. For such case,
422 // the following computation is skipped.
423 if (is_u_non_zero)
424 {
425 eps.noalias() = B * local_u;
426 }
427
428 // Consider also anisotropic thermal expansion.
429 auto const solid_linear_thermal_expansivity_vector =
431 solid_phase
432 .property(
434 .value(variables, x_position, t, dt));
435
437 dthermal_strain = solid_linear_thermal_expansivity_vector * dT_ip;
438
439 eps_m.noalias() = eps_m_prev + eps - eps_prev - dthermal_strain;
440
441 variables_prev.stress
443 sigma_prev);
444 variables_prev.mechanical_strain
446 eps_m_prev);
447 variables_prev.temperature = T_ip;
448 variables.mechanical_strain
450 eps_m);
451
452 auto&& solution = _ip_data[ip].solid_material.integrateStress(
453 variables_prev, variables, t, x_position, dt, *state);
454
455 if (!solution)
456 {
457 OGS_FATAL("Computation of local constitutive relation failed.");
458 }
459
461 std::tie(sigma, state, C) = std::move(*solution);
462
463 local_Jac.noalias() += B.transpose() * C * B * w;
464
465 typename ShapeMatricesType::template MatrixType<DisplacementDim,
467 N_u = ShapeMatricesType::template MatrixType<
468 DisplacementDim, displacement_size>::Zero(DisplacementDim,
470
471 for (int i = 0; i < DisplacementDim; ++i)
472 {
473 N_u.template block<1, displacement_size / DisplacementDim>(
474 i, i * displacement_size / DisplacementDim)
475 .noalias() = N;
476 }
477
478 auto const rho_s =
479 solid_phase.property(MPL::PropertyType::density)
480 .template value<double>(variables, x_position, t, dt);
481
482 auto const& b = _process_data.specific_body_force;
483 local_rhs.noalias() -=
484 (B.transpose() * sigma - N_u.transpose() * rho_s * b) * w;
485 }
486}

References ProcessLib::LinearBMatrix::computeBMatrix(), MathLib::createZeroedMatrix(), MathLib::createZeroedVector(), MaterialPropertyLib::formKelvinVector(), NumLib::interpolateXCoordinate(), MaterialPropertyLib::VariableArray::mechanical_strain, OGS_FATAL, ParameterLib::SpatialPosition::setElementID(), ParameterLib::SpatialPosition::setIntegrationPoint(), MaterialPropertyLib::VariableArray::stress, MaterialPropertyLib::VariableArray::temperature, and MaterialPropertyLib::thermal_expansivity.

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

494{
495 auto const local_T =
496 local_x.template segment<temperature_size>(temperature_index);
497
498 auto const local_T_prev =
499 local_x_prev.template segment<temperature_size>(temperature_index);
500
501 auto local_Jac = MathLib::createZeroedMatrix<
502 typename ShapeMatricesType::template MatrixType<temperature_size,
504 local_Jac_data, temperature_size, temperature_size);
505
506 auto local_rhs = MathLib::createZeroedVector<
507 typename ShapeMatricesType::template VectorType<temperature_size>>(
508 local_b_data, temperature_size);
509
512
513 typename ShapeMatricesType::NodalMatrixType laplace;
514 laplace.setZero(temperature_size, temperature_size);
515
516 unsigned const n_integration_points =
518
520 x_position.setElementID(_element.getID());
521 auto const& medium = _process_data.media_map.getMedium(_element.getID());
522 auto const& solid_phase = medium->phase("Solid");
523 MPL::VariableArray variables;
524
525 for (unsigned ip = 0; ip < n_integration_points; ip++)
526 {
527 x_position.setIntegrationPoint(ip);
528 auto const& w = _ip_data[ip].integration_weight;
529 auto const& N = _ip_data[ip].N;
530 auto const& dNdx = _ip_data[ip].dNdx;
531
532 const double T_ip = N.dot(local_T); // T at integration point
533 variables.temperature = T_ip;
534
535 auto const rho_s =
536 solid_phase.property(MPL::PropertyType::density)
537 .template value<double>(variables, x_position, t, dt);
538 auto const c_p =
539 solid_phase.property(MPL::PropertyType::specific_heat_capacity)
540 .template value<double>(variables, x_position, t, dt);
541
542 mass.noalias() += N.transpose() * rho_s * c_p * N * w;
543
544 auto const lambda =
545 solid_phase
546 .property(
548 .value(variables, x_position, t, dt);
549
552
553 laplace.noalias() += dNdx.transpose() * thermal_conductivity * dNdx * w;
554 }
555 local_Jac.noalias() += laplace + mass / dt;
556
557 local_rhs.noalias() -=
558 laplace * local_T + mass * (local_T - local_T_prev) / dt;
559}

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

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

336{
337 // For the equations with pressure
338 if (process_id == _process_data.heat_conduction_process_id)
339 {
341 t, dt, local_x, local_x_prev, local_b_data, local_Jac_data);
342 return;
343 }
344
345 // For the equations with deformation
346 assembleWithJacobianForDeformationEquations(t, dt, local_x, local_x_prev,
347 local_b_data, local_Jac_data);
348}
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)

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

606{
607 constexpr int kelvin_vector_size =
609
611 [this](std::vector<double>& values)
612 { return getIntPtEpsilon(0, {}, {}, values); });
613}
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 MathLib::KelvinVector::kelvin_vector_dimensions(), and ProcessLib::transposeInPlace().

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

650{
651 constexpr int kelvin_vector_size =
653
655 [this](std::vector<double>& values)
656 { return getIntPtEpsilonMechanical(0, {}, {}, values); });
657}
std::vector< double > const & getIntPtEpsilonMechanical(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override

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

◆ 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

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

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

622{
624 _ip_data, &IpData::eps, cache);
625}
std::vector< double > const & getIntegrationPointKelvinVectorData(IntegrationPointDataVector const &ip_data_vector, MemberType IpData::*const member, std::vector< double > &cache)

References ProcessLib::getIntegrationPointKelvinVectorData().

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

669{
670 return _process_data.material_ids == nullptr
671 ? 0
672 : (*_process_data.material_ids)[_element.getID()];
673}

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

680{
681 return *_ip_data[integration_point].material_state_variables;
682}

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

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

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

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

573{
574 constexpr int kelvin_vector_size =
576
578 [this](std::vector<double>& values)
579 { return getIntPtSigma(0, {}, {}, values); });
580}
std::vector< double > const & getIntPtSigma(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override

References MathLib::KelvinVector::kelvin_vector_dimensions(), 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 =
198 for (unsigned ip = 0; ip < n_integration_points; ip++)
199 {
200 auto& ip_data = _ip_data[ip];
201
202 ParameterLib::SpatialPosition const x_position{
203 std::nullopt, _element.getID(), ip,
205 NumLib::interpolateCoordinates<ShapeFunction,
207 _element, ip_data.N))};
208
210 if (_process_data.initial_stress != nullptr)
211 {
212 ip_data.sigma =
214 DisplacementDim>((*_process_data.initial_stress)(
215 std::numeric_limits<
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 }
Eigen::Matrix< double, Eigen::MatrixBase< Derived >::RowsAtCompileTime, 1 > symmetricTensorToKelvinVector(Eigen::MatrixBase< Derived > const &v)
std::array< double, 3 > interpolateCoordinates(MeshLib::Element const &e, typename ShapeMatricesType::ShapeMatrices::ShapeType const &N)

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().

◆ 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

◆ setEpsilon()

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

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

598{
600 values, _ip_data, &IpData::eps);
601}
std::size_t setIntegrationPointKelvinVectorData(double const *values, IntegrationPointDataVector &ip_data_vector, MemberType IpData::*const member)

References ProcessLib::setIntegrationPointKelvinVectorData().

◆ setEpsilonMechanical()

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

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

87{
88 if (integration_order !=
89 static_cast<int>(_integration_method.getIntegrationOrder()))
90 {
92 "Setting integration point initial conditions; The integration "
93 "order of the local assembler for element {:d} is different from "
94 "the integration order in the initial condition.",
95 _element.getID());
96 }
97
98 if (name == "sigma")
99 {
100 return setSigma(values);
101 }
102 if (name == "epsilon")
103 {
104 return setEpsilon(values);
105 }
106 if (name == "epsilon_m")
107 {
108 return setEpsilonMechanical(values);
109 }
110
111 return 0;
112}

References OGS_FATAL.

◆ setSigma()

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

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

Definition at line 353 of file ThermoMechanicsFEM.h.

◆ _process_data

◆ _secondary_data

◆ displacement_index

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

Definition at line 357 of file ThermoMechanicsFEM.h.

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

◆ 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

Definition at line 355 of file ThermoMechanicsFEM.h.

◆ temperature_size

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

Definition at line 356 of file ThermoMechanicsFEM.h.


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