OGS
ThermoHydroMechanicsProcess.cpp
Go to the documentation of this file.
1
12
13#include <cassert>
14
19#include "ProcessLib/Process.h"
24
25namespace ProcessLib
26{
27namespace ThermoHydroMechanics
28{
29template <int DisplacementDim>
31 std::string name,
32 MeshLib::Mesh& mesh,
33 std::unique_ptr<ProcessLib::AbstractJacobianAssembler>&& jacobian_assembler,
34 std::vector<std::unique_ptr<ParameterLib::ParameterBase>> const& parameters,
35 unsigned const integration_order,
36 std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>&&
37 process_variables,
39 SecondaryVariableCollection&& secondary_variables,
40 bool const use_monolithic_scheme)
41 : Process(std::move(name), mesh, std::move(jacobian_assembler), parameters,
42 integration_order, std::move(process_variables),
43 std::move(secondary_variables), use_monolithic_scheme),
44 _process_data(std::move(process_data))
45{
47 mesh, "NodalForces", MeshLib::MeshItemType::Node, DisplacementDim);
48
50 mesh, "MassFlowRate", MeshLib::MeshItemType::Node, 1);
51
53 mesh, "HeatFlowRate", MeshLib::MeshItemType::Node, 1);
54
55 _integration_point_writer.emplace_back(
56 std::make_unique<MeshLib::IntegrationPointWriter>(
57 "sigma_ip",
58 static_cast<int>(mesh.getDimension() == 2 ? 4 : 6) /*n components*/,
59 integration_order, _local_assemblers,
61
62 _integration_point_writer.emplace_back(
63 std::make_unique<MeshLib::IntegrationPointWriter>(
64 "epsilon_m_ip",
65 static_cast<int>(mesh.getDimension() == 2 ? 4 : 6) /*n components*/,
66 integration_order, _local_assemblers,
68
69 _integration_point_writer.emplace_back(
70 std::make_unique<MeshLib::IntegrationPointWriter>(
71 "epsilon_ip",
72 static_cast<int>(mesh.getDimension() == 2 ? 4 : 6) /*n components*/,
73 integration_order, _local_assemblers,
75}
76
77template <int DisplacementDim>
79{
80 return false;
81}
82
83template <int DisplacementDim>
86 const int process_id) const
87{
88 // For the monolithic scheme or the M process (deformation) in the staggered
89 // scheme.
90 if (_use_monolithic_scheme || process_id == 2)
91 {
92 auto const& l = *_local_to_global_index_map;
93 return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
94 &l.getGhostIndices(), &this->_sparsity_pattern};
95 }
96
97 // For staggered scheme and T or H process (pressure).
98 auto const& l = *_local_to_global_index_map_with_base_nodes;
99 return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
100 &l.getGhostIndices(), &_sparsity_pattern_with_linear_element};
101}
102
103template <int DisplacementDim>
105{
106 // Create single component dof in every of the mesh's nodes.
107 _mesh_subset_all_nodes =
108 std::make_unique<MeshLib::MeshSubset>(_mesh, _mesh.getNodes());
109 // Create single component dof in the mesh's base nodes.
110 _base_nodes = MeshLib::getBaseNodes(_mesh.getElements());
111 _mesh_subset_base_nodes =
112 std::make_unique<MeshLib::MeshSubset>(_mesh, _base_nodes);
113
114 // TODO move the two data members somewhere else.
115 // for extrapolation of secondary variables of stress or strain
116 std::vector<MeshLib::MeshSubset> all_mesh_subsets_single_component{
117 *_mesh_subset_all_nodes};
118 _local_to_global_index_map_single_component =
119 std::make_unique<NumLib::LocalToGlobalIndexMap>(
120 std::move(all_mesh_subsets_single_component),
121 // by location order is needed for output
123
124 if (_use_monolithic_scheme)
125 {
126 // For temperature, which is the first
127 std::vector<MeshLib::MeshSubset> all_mesh_subsets{
128 *_mesh_subset_base_nodes};
129
130 // For pressure, which is the second
131 all_mesh_subsets.push_back(*_mesh_subset_base_nodes);
132
133 // For displacement.
134 const int monolithic_process_id = 0;
135 std::generate_n(std::back_inserter(all_mesh_subsets),
136 getProcessVariables(monolithic_process_id)[2]
137 .get()
138 .getNumberOfGlobalComponents(),
139 [&]() { return *_mesh_subset_all_nodes; });
140
141 std::vector<int> const vec_n_components{1, 1, DisplacementDim};
142 _local_to_global_index_map =
143 std::make_unique<NumLib::LocalToGlobalIndexMap>(
144 std::move(all_mesh_subsets), vec_n_components,
146 assert(_local_to_global_index_map);
147 }
148 else
149 {
150 // For displacement equation.
151 const int process_id = 2;
152 std::vector<MeshLib::MeshSubset> all_mesh_subsets;
153 std::generate_n(std::back_inserter(all_mesh_subsets),
154 getProcessVariables(process_id)[0]
155 .get()
156 .getNumberOfGlobalComponents(),
157 [&]() { return *_mesh_subset_all_nodes; });
158
159 std::vector<int> const vec_n_components{DisplacementDim};
160 _local_to_global_index_map =
161 std::make_unique<NumLib::LocalToGlobalIndexMap>(
162 std::move(all_mesh_subsets), vec_n_components,
164
165 // For pressure equation or temperature equation.
166 // Collect the mesh subsets with base nodes in a vector.
167 std::vector<MeshLib::MeshSubset> all_mesh_subsets_base_nodes{
168 *_mesh_subset_base_nodes};
169 _local_to_global_index_map_with_base_nodes =
170 std::make_unique<NumLib::LocalToGlobalIndexMap>(
171 std::move(all_mesh_subsets_base_nodes),
172 // by location order is needed for output
174
175 _sparsity_pattern_with_linear_element = NumLib::computeSparsityPattern(
176 *_local_to_global_index_map_with_base_nodes, _mesh);
177
178 assert(_local_to_global_index_map);
179 assert(_local_to_global_index_map_with_base_nodes);
180 }
181}
182
183template <int DisplacementDim>
185 NumLib::LocalToGlobalIndexMap const& dof_table,
186 MeshLib::Mesh const& mesh,
187 unsigned const integration_order)
188{
191 mesh.getElements(), dof_table, _local_assemblers,
192 NumLib::IntegrationOrder{integration_order}, mesh.isAxiallySymmetric(),
193 _process_data);
194
195 auto add_secondary_variable = [&](std::string const& name,
196 int const num_components,
197 auto get_ip_values_function)
198 {
199 _secondary_variables.addSecondaryVariable(
200 name,
201 makeExtrapolator(num_components, getExtrapolator(),
202 _local_assemblers,
203 std::move(get_ip_values_function)));
204 };
205
206 add_secondary_variable(
207 "sigma",
209 DisplacementDim>::RowsAtCompileTime,
211
212 add_secondary_variable(
213 "sigma_ice",
215 DisplacementDim>::RowsAtCompileTime,
217
218 add_secondary_variable(
219 "epsilon_m",
221 DisplacementDim>::RowsAtCompileTime,
223
224 add_secondary_variable(
225 "epsilon",
227 DisplacementDim>::RowsAtCompileTime,
229
230 add_secondary_variable(
231 "ice_volume_fraction", 1,
233
234 add_secondary_variable(
235 "velocity", mesh.getDimension(),
237
238 add_secondary_variable(
239 "fluid_density", 1,
241
242 add_secondary_variable(
243 "viscosity", 1,
245
246 _process_data.element_fluid_density =
248 const_cast<MeshLib::Mesh&>(mesh), "fluid_density_avg",
250
251 _process_data.element_viscosity = MeshLib::getOrCreateMeshProperty<double>(
252 const_cast<MeshLib::Mesh&>(mesh), "viscosity_avg",
254
255 _process_data.element_stresses = MeshLib::getOrCreateMeshProperty<double>(
256 const_cast<MeshLib::Mesh&>(mesh), "stress_avg",
259 DisplacementDim>::RowsAtCompileTime);
260
261 _process_data.pressure_interpolated =
263 const_cast<MeshLib::Mesh&>(mesh), "pressure_interpolated",
265
266 _process_data.temperature_interpolated =
268 const_cast<MeshLib::Mesh&>(mesh), "temperature_interpolated",
270
271 //
272 // enable output of internal variables defined by material models
273 //
275 LocalAssemblerInterface<DisplacementDim>>(_process_data.solid_materials,
276 add_secondary_variable);
277
280 _process_data.solid_materials, _local_assemblers,
281 _integration_point_writer, integration_order);
282
283 setIPDataInitialConditions(_integration_point_writer, mesh.getProperties(),
284 _local_assemblers);
285
286 // Initialize local assemblers after all variables have been set.
289 _local_assemblers, *_local_to_global_index_map);
290}
291
292template <int DisplacementDim>
294 std::map<int, std::shared_ptr<MaterialPropertyLib::Medium>> const& media)
295{
296 if (_use_monolithic_scheme)
297 {
298 const int process_id_of_thermohydromechancs = 0;
299 initializeProcessBoundaryConditionsAndSourceTerms(
300 *_local_to_global_index_map, process_id_of_thermohydromechancs,
301 media);
302 return;
303 }
304
305 // Staggered scheme:
306 // for the equations of heat transport
307 const int thermal_process_id = 0;
308 initializeProcessBoundaryConditionsAndSourceTerms(
309 *_local_to_global_index_map_with_base_nodes, thermal_process_id, media);
310
311 // for the equations of mass balance
312 const int hydraulic_process_id = 1;
313 initializeProcessBoundaryConditionsAndSourceTerms(
314 *_local_to_global_index_map_with_base_nodes, hydraulic_process_id,
315 media);
316
317 // for the equations of deformation.
318 const int mechanical_process_id = 2;
319 initializeProcessBoundaryConditionsAndSourceTerms(
320 *_local_to_global_index_map, mechanical_process_id, media);
321}
322
323template <int DisplacementDim>
325 setInitialConditionsConcreteProcess(std::vector<GlobalVector*>& x,
326 double const t,
327 int const process_id)
328{
329 DBUG("SetInitialConditions ThermoHydroMechanicsProcess.");
330
333 _local_assemblers, getDOFTables(x.size()), x, t, process_id);
334}
335
336template <int DisplacementDim>
338 const double t, double const dt, std::vector<GlobalVector*> const& x,
339 std::vector<GlobalVector*> const& x_prev, int const process_id,
341{
342 DBUG("Assemble the equations for ThermoHydroMechanics");
343
344 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_table = {
345 _local_to_global_index_map.get()};
346
348 _global_assembler, &VectorMatrixAssembler::assemble, _local_assemblers,
349 getActiveElementIDs(), dof_table, t, dt, x, x_prev, process_id, &M, &K,
350 &b);
351}
352
353template <int DisplacementDim>
356 const double t, double const dt, std::vector<GlobalVector*> const& x,
357 std::vector<GlobalVector*> const& x_prev, int const process_id,
358 GlobalVector& b, GlobalMatrix& Jac)
359{
360 // For the monolithic scheme
361 if (_use_monolithic_scheme)
362 {
363 DBUG(
364 "Assemble the Jacobian of ThermoHydroMechanics for the monolithic "
365 "scheme.");
366 }
367 else
368 {
369 // For the staggered scheme
370 if (process_id == 0)
371 {
372 DBUG(
373 "Assemble the Jacobian equations of heat transport process in "
374 "ThermoHydroMechanics for the staggered scheme.");
375 }
376 else if (process_id == 1)
377 {
378 DBUG(
379 "Assemble the Jacobian equations of liquid fluid process in "
380 "ThermoHydroMechanics for the staggered scheme.");
381 }
382 else
383 {
384 DBUG(
385 "Assemble the Jacobian equations of mechanical process in "
386 "ThermoHydroMechanics for the staggered scheme.");
387 }
388 }
389
390 auto const dof_tables = getDOFTables(x.size());
391
394 _local_assemblers, getActiveElementIDs(), dof_tables, t, dt, x, x_prev,
395 process_id, &b, &Jac);
396
397 auto copyRhs = [&](int const variable_id, auto& output_vector)
398 {
399 if (_use_monolithic_scheme)
400 {
401 transformVariableFromGlobalVector(b, variable_id, *dof_tables[0],
402 output_vector,
403 std::negate<double>());
404 }
405 else
406 {
407 transformVariableFromGlobalVector(b, 0, *dof_tables[process_id],
408 output_vector,
409 std::negate<double>());
410 }
411 };
412 if (_use_monolithic_scheme || process_id == 0)
413 {
414 copyRhs(0, *_heat_flux);
415 }
416 if (_use_monolithic_scheme || process_id == 1)
417 {
418 copyRhs(1, *_hydraulic_flow);
419 }
420 if (_use_monolithic_scheme || process_id == 2)
421 {
422 copyRhs(2, *_nodal_forces);
423 }
424}
425
426template <int DisplacementDim>
428 std::vector<GlobalVector*> const& x, double const t, double const dt,
429 const int process_id)
430{
431 DBUG("PreTimestep ThermoHydroMechanicsProcess.");
432
433 if (hasMechanicalProcess(process_id))
434 {
437 _local_assemblers, *_local_to_global_index_map, *x[process_id], t,
438 dt);
439 }
440}
441
442template <int DisplacementDim>
444 std::vector<GlobalVector*> const& x,
445 std::vector<GlobalVector*> const& x_prev, double const t, double const dt,
446 const int process_id)
447{
448 if (process_id != 0)
449 {
450 return;
451 }
452
453 DBUG("PostTimestep ThermoHydroMechanicsProcess.");
454
457 _local_assemblers, getActiveElementIDs(), getDOFTables(x.size()), x,
458 x_prev, t, dt, process_id);
459}
460
461template <int DisplacementDim>
463 computeSecondaryVariableConcrete(double const t, double const dt,
464 std::vector<GlobalVector*> const& x,
465 GlobalVector const& x_prev,
466 const int process_id)
467{
468 if (process_id != 0)
469 {
470 return;
471 }
472
473 DBUG("Compute the secondary variables for ThermoHydroMechanicsProcess.");
474
477 _local_assemblers, getActiveElementIDs(), getDOFTables(x.size()), t, dt,
478 x, x_prev, process_id);
479}
480
481template <int DisplacementDim>
482std::tuple<NumLib::LocalToGlobalIndexMap*, bool> ThermoHydroMechanicsProcess<
483 DisplacementDim>::getDOFTableForExtrapolatorData() const
484{
485 const bool manage_storage = false;
486 return std::make_tuple(_local_to_global_index_map_single_component.get(),
487 manage_storage);
488}
489
490template <int DisplacementDim>
493 const int process_id) const
494{
495 if (hasMechanicalProcess(process_id))
496 {
497 return *_local_to_global_index_map;
498 }
499
500 // For the equation of pressure
501 return *_local_to_global_index_map_with_base_nodes;
502}
503
504template class ThermoHydroMechanicsProcess<2>;
505template class ThermoHydroMechanicsProcess<3>;
506
507} // namespace ThermoHydroMechanics
508} // 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
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 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
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 preTimestepConcreteProcess(std::vector< GlobalVector * > const &x, double const t, double const dt, const int process_id) override
ThermoHydroMechanicsProcess(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, ThermoHydroMechanicsProcessData< DisplacementDim > &&process_data, SecondaryVariableCollection &&secondary_variables, bool const use_monolithic_scheme)
void initializeBoundaryConditions(std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media) override
std::vector< std::unique_ptr< LocalAssemblerInterface< DisplacementDim > > > _local_assemblers
void computeSecondaryVariableConcrete(double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_prev, const int process_id) 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 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
void setInitialConditionsConcreteProcess(std::vector< GlobalVector * > &x, double const t, int const process_id) override
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)
std::vector< Node * > getBaseNodes(std::vector< Element * > const &elements)
Definition Utils.h:26
@ 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 solidMaterialInternalVariablesToIntegrationPointWriter(std::map< int, std::shared_ptr< SolidMaterial > > const &solid_materials, std::vector< std::unique_ptr< LocalAssemblerInterface > > const &local_assemblers, std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > &integration_point_writer, int const integration_order)
void solidMaterialInternalToSecondaryVariables(std::map< int, std::shared_ptr< SolidMaterial > > const &solid_materials, AddSecondaryVariableCallback const &add_secondary_variable)
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)
void createLocalAssemblersHM(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)
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)