OGS
ThermoMechanicsProcess.cpp
Go to the documentation of this file.
1
12
13#include <cassert>
14
22#include "ThermoMechanicsFEM.h"
23
24namespace ProcessLib
25{
26namespace ThermoMechanics
27{
28template <int DisplacementDim>
30 std::string name, MeshLib::Mesh& mesh,
31 std::unique_ptr<ProcessLib::AbstractJacobianAssembler>&& jacobian_assembler,
32 std::vector<std::unique_ptr<ParameterLib::ParameterBase>> const& parameters,
33 unsigned const integration_order,
34 std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>&&
35 process_variables,
37 SecondaryVariableCollection&& secondary_variables,
38 bool const use_monolithic_scheme)
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{
45 mesh, "NodalForces", MeshLib::MeshItemType::Node, DisplacementDim);
46
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,
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,
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,
70}
71
72template <int DisplacementDim>
74{
75 return false;
76}
77
78template <int DisplacementDim>
81 const int process_id) const
82{
83 // For the monolithic scheme or the M process (deformation) in the staggered
84 // scheme.
85 if (_use_monolithic_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.
94 auto const& l = *_local_to_global_index_map_single_component;
95 return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
96 &l.getGhostIndices(), &_sparsity_pattern_with_single_component};
97}
98
99// TODO [WW]: remove if (_use_monolithic_scheme) during the refactoring of the
100// coupling part.
101template <int DisplacementDim>
103{
104 // Note: the heat conduction process and the mechanical process use the same
105 // order of shape functions.
106
107 if (_use_monolithic_scheme)
108 {
109 constructMonolithicProcessDofTable();
110 return;
111 }
112 constructDofTableOfSpecifiedProcessStaggeredScheme(
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{
118 *_mesh_subset_all_nodes};
119 _local_to_global_index_map_single_component.reset(
121 std::move(all_mesh_subsets_single_component),
122 // by location order is needed for output
124
125 if (!_use_monolithic_scheme)
126 {
127 _sparsity_pattern_with_single_component =
129 *_local_to_global_index_map_single_component, _mesh);
130 }
131}
132
133template <int DisplacementDim>
135 NumLib::LocalToGlobalIndexMap const& dof_table,
136 MeshLib::Mesh const& mesh,
137 unsigned const integration_order)
138{
140 DisplacementDim, ThermoMechanicsLocalAssembler>(
141 mesh.getElements(), dof_table, _local_assemblers,
142 NumLib::IntegrationOrder{integration_order}, mesh.isAxiallySymmetric(),
143 _process_data);
144
145 auto add_secondary_variable = [&](std::string const& name,
146 int const num_components,
147 auto get_ip_values_function)
148 {
149 _secondary_variables.addSecondaryVariable(
150 name,
151 makeExtrapolator(num_components, getExtrapolator(),
152 _local_assemblers,
153 std::move(get_ip_values_function)));
154 };
155
156 add_secondary_variable("sigma",
158 DisplacementDim>::RowsAtCompileTime,
159 &LocalAssemblerInterface::getIntPtSigma);
160
161 add_secondary_variable("epsilon",
163 DisplacementDim>::RowsAtCompileTime,
164 &LocalAssemblerInterface::getIntPtEpsilon);
165
166 //
167 // enable output of internal variables defined by material models
168 //
170 LocalAssemblerInterface>(_process_data.solid_materials,
171 add_secondary_variable);
172
173 setIPDataInitialConditions(_integration_point_writer, mesh.getProperties(),
174 _local_assemblers);
175
176 // Initialize local assemblers after all variables have been set.
178 &LocalAssemblerInterface::initialize, _local_assemblers,
179 *_local_to_global_index_map);
180}
181
182template <int DisplacementDim>
184 std::map<int, std::shared_ptr<MaterialPropertyLib::Medium>> const& media)
185{
186 if (_use_monolithic_scheme)
187 {
188 const int process_id_of_thermomechanics = 0;
189 initializeProcessBoundaryConditionsAndSourceTerms(
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
196 initializeProcessBoundaryConditionsAndSourceTerms(
197 *_local_to_global_index_map_single_component,
198 _process_data.heat_conduction_process_id, media);
199
200 // for the equations of deformation.
201 initializeProcessBoundaryConditionsAndSourceTerms(
202 *_local_to_global_index_map, _process_data.mechanics_process_id, media);
203}
204
205template <int DisplacementDim>
207 const double t, double const dt, std::vector<GlobalVector*> const& x,
208 std::vector<GlobalVector*> const& x_prev, int const process_id,
210{
211 DBUG("Assemble ThermoMechanicsProcess.");
212
213 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_table = {
214 _local_to_global_index_map.get()};
215
216 // Call global assembler for each local assembly item.
218 _global_assembler, &VectorMatrixAssembler::assemble, _local_assemblers,
219 getActiveElementIDs(), dof_table, t, dt, x, x_prev, process_id, &M, &K,
220 &b);
221}
222
223template <int DisplacementDim>
226 const double t, double const dt, std::vector<GlobalVector*> const& x,
227 std::vector<GlobalVector*> const& x_prev, int const process_id,
228 GlobalVector& b, GlobalMatrix& Jac)
229{
230 DBUG("AssembleJacobian ThermoMechanicsProcess.");
231
232 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
233 // For the monolithic scheme
234 if (_use_monolithic_scheme)
235 {
236 DBUG(
237 "Assemble the Jacobian of ThermoMechanics for the monolithic"
238 " scheme.");
239 dof_tables.emplace_back(_local_to_global_index_map.get());
240 }
241 else
242 {
243 // For the staggered scheme
244 if (process_id == _process_data.heat_conduction_process_id)
245 {
246 DBUG(
247 "Assemble the Jacobian equations of heat conduction process in "
248 "ThermoMechanics for the staggered scheme.");
249 }
250 else
251 {
252 DBUG(
253 "Assemble the Jacobian equations of mechanical process in "
254 "ThermoMechanics for the staggered scheme.");
255 }
256
257 // For the flexible appearance order of processes in the coupling.
258 if (_process_data.heat_conduction_process_id ==
259 0) // First: the heat conduction process
260 {
261 dof_tables.emplace_back(
262 _local_to_global_index_map_single_component.get());
263 dof_tables.emplace_back(_local_to_global_index_map.get());
264 }
265 else // vice versa
266 {
267 dof_tables.emplace_back(_local_to_global_index_map.get());
268 dof_tables.emplace_back(
269 _local_to_global_index_map_single_component.get());
270 }
271 }
272
275 _local_assemblers, getActiveElementIDs(), dof_tables, t, dt, x, x_prev,
276 process_id, &b, &Jac);
277
278 // TODO (naumov): Refactor the copy rhs part. This is copy from HM.
279 auto copyRhs = [&](int const variable_id, auto& output_vector)
280 {
281 if (_use_monolithic_scheme)
282 {
283 transformVariableFromGlobalVector(b, variable_id, *dof_tables[0],
284 output_vector,
285 std::negate<double>());
286 }
287 else
288 {
289 transformVariableFromGlobalVector(b, 0, *dof_tables[process_id],
290 output_vector,
291 std::negate<double>());
292 }
293 };
294 if (_use_monolithic_scheme ||
295 process_id == _process_data.heat_conduction_process_id)
296 {
297 copyRhs(0, *_heat_flux);
298 }
299 if (_use_monolithic_scheme ||
300 process_id == _process_data.mechanics_process_id)
301 {
302 copyRhs(1, *_nodal_forces);
303 }
304}
305
306template <int DisplacementDim>
308 std::vector<GlobalVector*> const& x, double const t, double const dt,
309 const int process_id)
310{
311 DBUG("PreTimestep ThermoMechanicsProcess.");
312
313 assert(process_id < 2);
314
315 if (process_id == _process_data.mechanics_process_id)
316 {
318 &LocalAssemblerInterface::preTimestep, _local_assemblers,
319 getActiveElementIDs(), *_local_to_global_index_map, *x[process_id],
320 t, dt);
321 return;
322 }
323}
324
325template <int DisplacementDim>
327 std::vector<GlobalVector*> const& x,
328 std::vector<GlobalVector*> const& x_prev, double const t, double const dt,
329 int const process_id)
330{
331 if (process_id != 0)
332 {
333 return;
334 }
335
336 DBUG("PostTimestep ThermoMechanicsProcess.");
337 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
338 auto const n_processes = x.size();
339 dof_tables.reserve(n_processes);
340 for (std::size_t process_id = 0; process_id < n_processes; ++process_id)
341 {
342 dof_tables.push_back(&getDOFTable(process_id));
343 }
344
346 &LocalAssemblerInterface::postTimestep, _local_assemblers,
347 getActiveElementIDs(), dof_tables, x, x_prev, t, dt, process_id);
348}
349
350template <int DisplacementDim>
353{
354 if (_process_data.mechanics_process_id == process_id)
355 {
356 return *_local_to_global_index_map;
357 }
358
359 // For the equation of pressure
360 return *_local_to_global_index_map_single_component;
361}
362
363template class ThermoMechanicsProcess<2>;
364template class ThermoMechanicsProcess<3>;
365
366} // namespace ThermoMechanics
367} // namespace ProcessLib
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
Global vector based on Eigen vector.
Definition EigenVector.h:25
bool isAxiallySymmetric() const
Definition Mesh.h:137
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
Definition Mesh.h:109
unsigned getDimension() const
Returns the dimension of the mesh (determined by the maximum dimension over all elements).
Definition Mesh.h:88
Properties & getProperties()
Definition Mesh.h:134
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)
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 initialize(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > _integration_point_writer
Definition Process.h:390
Handles configuration of several secondary variables from the project file.
void initializeBoundaryConditions(std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media) override
void preTimestepConcreteProcess(std::vector< GlobalVector * > const &x, double const t, double const dt, const int process_id) override
void assembleWithJacobianConcreteProcess(const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalVector &b, GlobalMatrix &Jac) 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
MathLib::MatrixSpecifications getMatrixSpecifications(const int process_id) const override
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)
std::vector< std::unique_ptr< LocalAssemblerInterface > > _local_assemblers
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
NumLib::LocalToGlobalIndexMap const & getDOFTable(const int process_id) const override
void initializeConcreteProcess(NumLib::LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh, unsigned const integration_order) override
Process specific initialization called by initialize().
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)
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, GlobalVector *b, GlobalMatrix *Jac)
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), 1, Eigen::ColMajor > KelvinVectorType
PropertyVector< T > * getOrCreateMeshProperty(Mesh &mesh, std::string const &property_name, MeshItemType const item_type, int const number_of_components)
@ 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.
void solidMaterialInternalToSecondaryVariables(std::map< int, std::shared_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 executeSelectedMemberOnDereferenced(Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)
static void executeSelectedMemberDereferenced(Object &object, Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)
static void executeMemberOnDereferenced(Method method, Container const &container, Args &&... args)
virtual std::vector< double > getEpsilonMechanical() const =0