OGS
ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim > Class Template Referencefinal

Detailed Description

template<int DisplacementDim>
class ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >

Definition at line 22 of file ThermoMechanicsProcess.h.

#include <ThermoMechanicsProcess.h>

Inheritance diagram for ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >:
[legend]
Collaboration diagram for ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >:
[legend]

Public Member Functions

 ThermoMechanicsProcess (std::string name, MeshLib::Mesh &mesh, std::unique_ptr< ProcessLib::AbstractJacobianAssembler > &&jacobian_assembler, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const &parameters, unsigned const integration_order, std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > &&process_variables, ThermoMechanicsProcessData< DisplacementDim > &&process_data, SecondaryVariableCollection &&secondary_variables, bool const use_monolithic_scheme)
 
MathLib::MatrixSpecifications getMatrixSpecifications (const int process_id) const override
 
ODESystem interface
bool isLinear () const override
 
- Public Member Functions inherited from ProcessLib::Process
 Process (std::string name_, MeshLib::Mesh &mesh, std::unique_ptr< AbstractJacobianAssembler > &&jacobian_assembler, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const &parameters, unsigned const integration_order, std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > &&process_variables, SecondaryVariableCollection &&secondary_variables, const bool use_monolithic_scheme=true)
 
void preTimestep (std::vector< GlobalVector * > const &x, const double t, const double delta_t, const int process_id)
 Preprocessing before starting assembly for new timestep.
 
void postTimestep (std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, const double t, const double delta_t, int const process_id)
 Postprocessing after a complete timestep.
 
void postNonLinearSolver (std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, const double t, double const dt, int const process_id)
 
void preIteration (const unsigned iter, GlobalVector const &x) final
 
void computeSecondaryVariable (double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_prev, int const process_id)
 compute secondary variables for the coupled equations or for output.
 
NumLib::IterationResult postIteration (GlobalVector const &x) final
 
void initialize (std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
 
void setInitialConditions (std::vector< GlobalVector * > &process_solutions, std::vector< GlobalVector * > const &process_solutions_prev, double const t, int const process_id)
 
MathLib::MatrixSpecifications getMatrixSpecifications (const int process_id) const override
 
void updateDeactivatedSubdomains (double const time, const int process_id)
 
virtual bool isMonolithicSchemeUsed () const
 
virtual void extrapolateIntegrationPointValuesToNodes (const double, std::vector< GlobalVector * > const &, std::vector< GlobalVector * > &)
 
void preAssemble (const double t, double const dt, GlobalVector const &x) final
 
void assemble (const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b) final
 
void assembleWithJacobian (const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, GlobalMatrix &Jac) final
 
void preOutput (const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id)
 
std::vector< NumLib::IndexValueVector< GlobalIndexType > > const * getKnownSolutions (double const t, GlobalVector const &x, int const process_id) const final
 
MeshLib::MeshgetMesh () const
 
std::vector< std::reference_wrapper< ProcessVariable > > const & getProcessVariables (const int process_id) const
 
std::vector< std::size_t > const & getActiveElementIDs () const
 
SecondaryVariableCollection const & getSecondaryVariables () const
 
std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > const & getIntegrationPointWriters () const
 
virtual Eigen::Vector3d getFlux (std::size_t, MathLib::Point3d const &, double const, std::vector< GlobalVector * > const &) const
 
virtual void solveReactionEquation (std::vector< GlobalVector * > &, std::vector< GlobalVector * > const &, double const, double const, NumLib::EquationSystem &, int const)
 
- Public Member Functions inherited from ProcessLib::SubmeshAssemblySupport
virtual std::vector< std::string > initializeAssemblyOnSubmeshes (std::vector< std::reference_wrapper< MeshLib::Mesh > > const &meshes)
 
virtual ~SubmeshAssemblySupport ()=default
 

Private Types

using LocalAssemblerInterface
 

Private Member Functions

void constructDofTable () override
 
void initializeConcreteProcess (NumLib::LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh, unsigned const integration_order) override
 Process specific initialization called by initialize().
 
void initializeBoundaryConditions (std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media) override
 
void assembleConcreteProcess (const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b) override
 
void assembleWithJacobianConcreteProcess (const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, GlobalMatrix &Jac) override
 
void preTimestepConcreteProcess (std::vector< GlobalVector * > const &x, double const t, double const dt, const int process_id) override
 
void postTimestepConcreteProcess (std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, const double t, const double dt, int const process_id) override
 
NumLib::LocalToGlobalIndexMap const & getDOFTable (const int process_id) const override
 

Private Attributes

ThermoMechanicsProcessData< DisplacementDim > _process_data
 
std::vector< std::unique_ptr< LocalAssemblerInterface > > _local_assemblers
 
std::unique_ptr< NumLib::LocalToGlobalIndexMap_local_to_global_index_map_single_component
 
GlobalSparsityPattern _sparsity_pattern_with_single_component
 
MeshLib::PropertyVector< double > * _nodal_forces = nullptr
 
MeshLib::PropertyVector< double > * _heat_flux = nullptr
 

Additional Inherited Members

- Public Attributes inherited from ProcessLib::Process
std::string const name
 
- Static Public Attributes inherited from ProcessLib::Process
static PROCESSLIB_EXPORT const std::string constant_one_parameter_name = "constant_one"
 
- Protected Member Functions inherited from ProcessLib::Process
std::vector< NumLib::LocalToGlobalIndexMap const * > getDOFTables (int const number_of_processes) const
 
NumLib::ExtrapolatorgetExtrapolator () const
 
NumLib::LocalToGlobalIndexMap const & getSingleComponentDOFTable () const
 
void initializeProcessBoundaryConditionsAndSourceTerms (const NumLib::LocalToGlobalIndexMap &dof_table, const int process_id, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
 
void constructMonolithicProcessDofTable ()
 
void constructDofTableOfSpecifiedProcessStaggeredScheme (const int specified_process_id)
 
virtual std::tuple< NumLib::LocalToGlobalIndexMap *, bool > getDOFTableForExtrapolatorData () const
 
std::vector< GlobalIndexTypegetIndicesOfResiduumWithoutInitialCompensation () const override
 
- Protected Attributes inherited from ProcessLib::Process
MeshLib::Mesh_mesh
 
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_all_nodes
 
std::unique_ptr< NumLib::LocalToGlobalIndexMap_local_to_global_index_map
 
SecondaryVariableCollection _secondary_variables
 
std::unique_ptr< ProcessLib::AbstractJacobianAssembler_jacobian_assembler
 
VectorMatrixAssembler _global_assembler
 
const bool _use_monolithic_scheme
 
unsigned const _integration_order
 
std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > _integration_point_writer
 
GlobalSparsityPattern _sparsity_pattern
 
std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > _process_variables
 
std::vector< BoundaryConditionCollection_boundary_conditions
 

Member Typedef Documentation

◆ LocalAssemblerInterface

template<int DisplacementDim>
using ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >::LocalAssemblerInterface
private
Initial value:
ThermoMechanicsLocalAssemblerInterface<DisplacementDim>

Definition at line 57 of file ThermoMechanicsProcess.h.

Constructor & Destructor Documentation

◆ ThermoMechanicsProcess()

template<int DisplacementDim>
ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >::ThermoMechanicsProcess ( std::string name,
MeshLib::Mesh & mesh,
std::unique_ptr< ProcessLib::AbstractJacobianAssembler > && jacobian_assembler,
std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const & parameters,
unsigned const integration_order,
std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > && process_variables,
ThermoMechanicsProcessData< DisplacementDim > && process_data,
SecondaryVariableCollection && secondary_variables,
bool const use_monolithic_scheme )

Definition at line 29 of file ThermoMechanicsProcess.cpp.

39 : Process(std::move(name), mesh, std::move(jacobian_assembler), parameters,
40 integration_order, std::move(process_variables),
41 std::move(secondary_variables), use_monolithic_scheme),
42 _process_data(std::move(process_data))
43{
44 _nodal_forces = MeshLib::getOrCreateMeshProperty<double>(
45 mesh, "NodalForces", MeshLib::MeshItemType::Node, DisplacementDim);
46
47 _heat_flux = MeshLib::getOrCreateMeshProperty<double>(
48 mesh, "HeatFlowRate", MeshLib::MeshItemType::Node, 1);
49
50 _integration_point_writer.emplace_back(
51 std::make_unique<MeshLib::IntegrationPointWriter>(
52 "sigma_ip",
53 static_cast<int>(mesh.getDimension() == 2 ? 4 : 6) /*n components*/,
54 integration_order, _local_assemblers,
55 &LocalAssemblerInterface::getSigma));
56
57 _integration_point_writer.emplace_back(
58 std::make_unique<MeshLib::IntegrationPointWriter>(
59 "epsilon_ip",
60 static_cast<int>(mesh.getDimension() == 2 ? 4 : 6) /*n components*/,
61 integration_order, _local_assemblers,
62 &LocalAssemblerInterface::getEpsilon));
63
64 _integration_point_writer.emplace_back(
65 std::make_unique<MeshLib::IntegrationPointWriter>(
66 "epsilon_m_ip",
67 static_cast<int>(mesh.getDimension() == 2 ? 4 : 6) /*n components*/,
68 integration_order, _local_assemblers,
69 &LocalAssemblerInterface::getEpsilonMechanical));
70}
unsigned getDimension() const
Returns the dimension of the mesh (determined by the maximum dimension over all elements).
Definition Mesh.h:88
std::string const name
Definition Process.h:354
std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > _integration_point_writer
Definition Process.h:380
Process(std::string name_, MeshLib::Mesh &mesh, std::unique_ptr< AbstractJacobianAssembler > &&jacobian_assembler, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const &parameters, unsigned const integration_order, std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > &&process_variables, SecondaryVariableCollection &&secondary_variables, const bool use_monolithic_scheme=true)
Definition Process.cpp:44
ThermoMechanicsLocalAssemblerInterface< DisplacementDim > LocalAssemblerInterface
std::vector< std::unique_ptr< LocalAssemblerInterface > > _local_assemblers
ThermoMechanicsProcessData< DisplacementDim > _process_data

References ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >::_heat_flux, ProcessLib::Process::_integration_point_writer, ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >::_local_assemblers, ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >::_nodal_forces, MeshLib::Mesh::getDimension(), ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssemblerInterface< DisplacementDim >::getEpsilon(), ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssemblerInterface< DisplacementDim >::getEpsilonMechanical(), ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssemblerInterface< DisplacementDim >::getSigma(), and MeshLib::Node.

Member Function Documentation

◆ assembleConcreteProcess()

template<int DisplacementDim>
void ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >::assembleConcreteProcess ( const double t,
double const dt,
std::vector< GlobalVector * > const & x,
std::vector< GlobalVector * > const & x_prev,
int const process_id,
GlobalMatrix & M,
GlobalMatrix & K,
GlobalVector & b )
overrideprivatevirtual

Implements ProcessLib::Process.

Definition at line 206 of file ThermoMechanicsProcess.cpp.

210{
211 DBUG("Assemble ThermoMechanicsProcess.");
212
213 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_table = {
215 ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
216
217 // Call global assembler for each local assembly item.
220 pv.getActiveElementIDs(), dof_table, t, dt, x, x_prev, process_id, M, K,
221 b);
222}
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
std::vector< std::size_t > const & getActiveElementIDs() const
std::vector< std::reference_wrapper< ProcessVariable > > const & getProcessVariables(const int process_id) const
Definition Process.h:155
VectorMatrixAssembler _global_assembler
Definition Process.h:367
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
Definition Process.h:360
void assemble(std::size_t const mesh_item_id, LocalAssemblerInterface &local_assembler, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b)
static void executeSelectedMemberDereferenced(Object &object, Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)

References ProcessLib::VectorMatrixAssembler::assemble(), DBUG(), NumLib::SerialExecutor::executeSelectedMemberDereferenced(), and ProcessLib::ProcessVariable::getActiveElementIDs().

◆ assembleWithJacobianConcreteProcess()

template<int DisplacementDim>
void ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >::assembleWithJacobianConcreteProcess ( const double t,
double const dt,
std::vector< GlobalVector * > const & x,
std::vector< GlobalVector * > const & x_prev,
int const process_id,
GlobalMatrix & M,
GlobalMatrix & K,
GlobalVector & b,
GlobalMatrix & Jac )
overrideprivatevirtual

Implements ProcessLib::Process.

Definition at line 225 of file ThermoMechanicsProcess.cpp.

230{
231 DBUG("AssembleJacobian ThermoMechanicsProcess.");
232
233 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
234 // For the monolithic scheme
236 {
237 DBUG(
238 "Assemble the Jacobian of ThermoMechanics for the monolithic"
239 " scheme.");
240 dof_tables.emplace_back(_local_to_global_index_map.get());
241 }
242 else
243 {
244 // For the staggered scheme
245 if (process_id == _process_data.heat_conduction_process_id)
246 {
247 DBUG(
248 "Assemble the Jacobian equations of heat conduction process in "
249 "ThermoMechanics for the staggered scheme.");
250 }
251 else
252 {
253 DBUG(
254 "Assemble the Jacobian equations of mechanical process in "
255 "ThermoMechanics for the staggered scheme.");
256 }
257
258 // For the flexible appearance order of processes in the coupling.
259 if (_process_data.heat_conduction_process_id ==
260 0) // First: the heat conduction process
261 {
262 dof_tables.emplace_back(
264 dof_tables.emplace_back(_local_to_global_index_map.get());
265 }
266 else // vice versa
267 {
268 dof_tables.emplace_back(_local_to_global_index_map.get());
269 dof_tables.emplace_back(
271 }
272 }
273
274 ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
275
278 _local_assemblers, pv.getActiveElementIDs(), dof_tables, t, dt, x,
279 x_prev, process_id, M, K, b, Jac);
280
281 // TODO (naumov): Refactor the copy rhs part. This is copy from HM.
282 auto copyRhs = [&](int const variable_id, auto& output_vector)
283 {
285 {
286 transformVariableFromGlobalVector(b, variable_id, *dof_tables[0],
287 output_vector,
288 std::negate<double>());
289 }
290 else
291 {
292 transformVariableFromGlobalVector(b, 0, *dof_tables[process_id],
293 output_vector,
294 std::negate<double>());
295 }
296 };
298 process_id == _process_data.heat_conduction_process_id)
299 {
300 copyRhs(0, *_heat_flux);
301 }
303 process_id == _process_data.mechanics_process_id)
304 {
305 copyRhs(1, *_nodal_forces);
306 }
307}
const bool _use_monolithic_scheme
Definition Process.h:369
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map_single_component
void assembleWithJacobian(std::size_t const mesh_item_id, LocalAssemblerInterface &local_assembler, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, GlobalMatrix &Jac)
void transformVariableFromGlobalVector(GlobalVector const &input_vector, int const variable_id, NumLib::LocalToGlobalIndexMap const &local_to_global_index_map, MeshLib::PropertyVector< double > &output_vector, Functor map_function)

References ProcessLib::VectorMatrixAssembler::assembleWithJacobian(), DBUG(), NumLib::SerialExecutor::executeSelectedMemberDereferenced(), and ProcessLib::ProcessVariable::getActiveElementIDs().

◆ constructDofTable()

template<int DisplacementDim>
void ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >::constructDofTable ( )
overrideprivatevirtual

This function is for general cases, in which all equations of the coupled processes have the same number of unknowns. For the general cases with the staggered scheme, all equations of the coupled processes share one DOF table hold by _local_to_global_index_map. Other cases can be considered by overloading this member function in the derived class.

Reimplemented from ProcessLib::Process.

Definition at line 102 of file ThermoMechanicsProcess.cpp.

103{
104 // Note: the heat conduction process and the mechanical process use the same
105 // order of shape functions.
106
108 {
110 return;
111 }
113 _process_data.mechanics_process_id);
114
115 // TODO move the two data members somewhere else.
116 // for extrapolation of secondary variables of stress or strain
117 std::vector<MeshLib::MeshSubset> all_mesh_subsets_single_component{
121 std::move(all_mesh_subsets_single_component),
122 // by location order is needed for output
124
126 {
130 }
131}
void constructDofTableOfSpecifiedProcessStaggeredScheme(const int specified_process_id)
Definition Process.cpp:355
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_all_nodes
Definition Process.h:358
MeshLib::Mesh & _mesh
Definition Process.h:357
void constructMonolithicProcessDofTable()
Definition Process.cpp:322
@ BY_LOCATION
Ordering data by spatial location.
GlobalSparsityPattern computeSparsityPattern(LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh)
Computes a sparsity pattern for the given inputs.

References NumLib::BY_LOCATION, and NumLib::computeSparsityPattern().

◆ getDOFTable()

template<int DisplacementDim>
NumLib::LocalToGlobalIndexMap const & ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >::getDOFTable ( const int process_id) const
overrideprivatevirtual

Reimplemented from ProcessLib::Process.

Definition at line 358 of file ThermoMechanicsProcess.cpp.

359{
360 if (_process_data.mechanics_process_id == process_id)
361 {
363 }
364
365 // For the equation of pressure
367}

◆ getMatrixSpecifications()

template<int DisplacementDim>
MathLib::MatrixSpecifications ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >::getMatrixSpecifications ( const int process_id) const
override

Get the size and the sparse pattern of the global matrix in order to create the global matrices and vectors for the system equations of this process.

Parameters
process_idProcess ID. If the monolithic scheme is applied, process_id = 0.
Returns
Matrix specifications including size and sparse pattern.

Definition at line 80 of file ThermoMechanicsProcess.cpp.

82{
83 // For the monolithic scheme or the M process (deformation) in the staggered
84 // scheme.
86 process_id == _process_data.mechanics_process_id)
87 {
88 auto const& l = *_local_to_global_index_map;
89 return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
90 &l.getGhostIndices(), &this->_sparsity_pattern};
91 }
92
93 // For staggered scheme and T process.
95 return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
96 &l.getGhostIndices(), &_sparsity_pattern_with_single_component};
97}
GlobalSparsityPattern _sparsity_pattern
Definition Process.h:382

◆ initializeBoundaryConditions()

template<int DisplacementDim>
void ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >::initializeBoundaryConditions ( std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const & media)
overrideprivatevirtual

Member function to initialize the boundary conditions for all coupled processes. It is called by initialize().

Reimplemented from ProcessLib::Process.

Definition at line 183 of file ThermoMechanicsProcess.cpp.

185{
187 {
188 const int process_id_of_thermomechanics = 0;
190 *_local_to_global_index_map, process_id_of_thermomechanics, media);
191 return;
192 }
193
194 // Staggered scheme:
195 // for the equations of heat conduction
198 _process_data.heat_conduction_process_id, media);
199
200 // for the equations of deformation.
202 *_local_to_global_index_map, _process_data.mechanics_process_id, media);
203}
void initializeProcessBoundaryConditionsAndSourceTerms(const NumLib::LocalToGlobalIndexMap &dof_table, const int process_id, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
Definition Process.cpp:89

◆ initializeConcreteProcess()

template<int DisplacementDim>
void ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >::initializeConcreteProcess ( NumLib::LocalToGlobalIndexMap const & dof_table,
MeshLib::Mesh const & mesh,
unsigned const integration_order )
overrideprivatevirtual

Process specific initialization called by initialize().

Implements ProcessLib::Process.

Definition at line 134 of file ThermoMechanicsProcess.cpp.

138{
140 DisplacementDim, ThermoMechanicsLocalAssembler>(
141 mesh.getElements(), dof_table, _local_assemblers,
142 NumLib::IntegrationOrder{integration_order}, mesh.isAxiallySymmetric(),
144
145 auto add_secondary_variable = [&](std::string const& name,
146 int const num_components,
147 auto get_ip_values_function)
148 {
150 name,
151 makeExtrapolator(num_components, getExtrapolator(),
153 std::move(get_ip_values_function)));
154 };
155
156 add_secondary_variable("sigma",
158 DisplacementDim>::RowsAtCompileTime,
160
161 add_secondary_variable("epsilon",
163 DisplacementDim>::RowsAtCompileTime,
165
166 //
167 // enable output of internal variables defined by material models
168 //
170 LocalAssemblerInterface>(_process_data.solid_materials,
171 add_secondary_variable);
172
175
176 // Initialize local assemblers after all variables have been set.
180}
virtual void initialize(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
SecondaryVariableCollection _secondary_variables
Definition Process.h:362
NumLib::Extrapolator & getExtrapolator() const
Definition Process.h:199
void addSecondaryVariable(std::string const &internal_name, SecondaryVariableFunctions &&fcts)
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), 1, Eigen::ColMajor > KelvinVectorType
void solidMaterialInternalToSecondaryVariables(std::map< int, std::unique_ptr< SolidMaterial > > const &solid_materials, AddSecondaryVariableCallback const &add_secondary_variable)
void createLocalAssemblers(std::vector< MeshLib::Element * > const &mesh_elements, NumLib::LocalToGlobalIndexMap const &dof_table, std::vector< std::unique_ptr< LocalAssemblerInterface > > &local_assemblers, ProviderOrOrder const &provider_or_order, ExtraCtorArgs &&... extra_ctor_args)
SecondaryVariableFunctions makeExtrapolator(const unsigned num_components, NumLib::Extrapolator &extrapolator, LocalAssemblerCollection const &local_assemblers, typename NumLib::ExtrapolatableLocalAssemblerCollection< LocalAssemblerCollection >::IntegrationPointValuesMethod integration_point_values_method)
void setIPDataInitialConditions(std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > const &_integration_point_writer, MeshLib::Properties const &mesh_properties, LocalAssemblersVector &local_assemblers)
static void executeMemberOnDereferenced(Method method, Container const &container, Args &&... args)
virtual 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 =0
virtual 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 =0

References ProcessLib::SmallDeformation::createLocalAssemblers(), NumLib::SerialExecutor::executeMemberOnDereferenced(), MeshLib::Mesh::getElements(), MeshLib::Mesh::getProperties(), ProcessLib::LocalAssemblerInterface::initialize(), MeshLib::Mesh::isAxiallySymmetric(), ProcessLib::makeExtrapolator(), ProcessLib::setIPDataInitialConditions(), and ProcessLib::Deformation::solidMaterialInternalToSecondaryVariables().

◆ isLinear()

template<int DisplacementDim>
bool ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >::isLinear ( ) const
override

Definition at line 73 of file ThermoMechanicsProcess.cpp.

74{
75 return false;
76}

◆ postTimestepConcreteProcess()

template<int DisplacementDim>
void ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >::postTimestepConcreteProcess ( std::vector< GlobalVector * > const & x,
std::vector< GlobalVector * > const & x_prev,
const double t,
const double dt,
int const process_id )
overrideprivatevirtual

Reimplemented from ProcessLib::Process.

Definition at line 331 of file ThermoMechanicsProcess.cpp.

335{
336 if (process_id != 0)
337 {
338 return;
339 }
340
341 DBUG("PostTimestep ThermoMechanicsProcess.");
342 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
343 auto const n_processes = x.size();
344 dof_tables.reserve(n_processes);
345 for (std::size_t process_id = 0; process_id < n_processes; ++process_id)
346 {
347 dof_tables.push_back(&getDOFTable(process_id));
348 }
349
350 ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
353 pv.getActiveElementIDs(), dof_tables, x, x_prev, t, dt, process_id);
354}
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)
NumLib::LocalToGlobalIndexMap const & getDOFTable(const int process_id) const override
static void executeSelectedMemberOnDereferenced(Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)

References DBUG(), NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), ProcessLib::ProcessVariable::getActiveElementIDs(), and ProcessLib::LocalAssemblerInterface::postTimestep().

◆ preTimestepConcreteProcess()

template<int DisplacementDim>
void ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >::preTimestepConcreteProcess ( std::vector< GlobalVector * > const & x,
double const t,
double const dt,
const int process_id )
overrideprivatevirtual

Reimplemented from ProcessLib::Process.

Definition at line 310 of file ThermoMechanicsProcess.cpp.

313{
314 DBUG("PreTimestep ThermoMechanicsProcess.");
315
316 ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
317
318 assert(process_id < 2);
319
320 if (process_id == _process_data.mechanics_process_id)
321 {
325 *x[process_id], t, dt);
326 return;
327 }
328}
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)

References DBUG(), NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), ProcessLib::ProcessVariable::getActiveElementIDs(), and ProcessLib::LocalAssemblerInterface::preTimestep().

Member Data Documentation

◆ _heat_flux

template<int DisplacementDim>
MeshLib::PropertyVector<double>* ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >::_heat_flux = nullptr
private

◆ _local_assemblers

template<int DisplacementDim>
std::vector<std::unique_ptr<LocalAssemblerInterface> > ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >::_local_assemblers
private

◆ _local_to_global_index_map_single_component

template<int DisplacementDim>
std::unique_ptr<NumLib::LocalToGlobalIndexMap> ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >::_local_to_global_index_map_single_component
private

Definition at line 101 of file ThermoMechanicsProcess.h.

◆ _nodal_forces

template<int DisplacementDim>
MeshLib::PropertyVector<double>* ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >::_nodal_forces = nullptr
private

◆ _process_data

template<int DisplacementDim>
ThermoMechanicsProcessData<DisplacementDim> ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >::_process_data
private

Definition at line 96 of file ThermoMechanicsProcess.h.

◆ _sparsity_pattern_with_single_component

template<int DisplacementDim>
GlobalSparsityPattern ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >::_sparsity_pattern_with_single_component
private

Sparsity pattern for the heat conduction equation, and it is initialized only if the staggered scheme is used.

Definition at line 105 of file ThermoMechanicsProcess.h.


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