OGS
ThermoHydroMechanicsProcess.cpp
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) OpenGeoSys Community (opengeosys.org)
2// SPDX-License-Identifier: BSD-3-Clause
3
5
6#include <cassert>
7
12#include "ProcessLib/Process.h"
17
18namespace ProcessLib
19{
21{
22template <int DisplacementDim>
24 std::string name, MeshLib::Mesh& mesh,
25 std::unique_ptr<ProcessLib::AbstractJacobianAssembler>&& jacobian_assembler,
26 std::vector<std::unique_ptr<ParameterLib::ParameterBase>> const& parameters,
27 unsigned const integration_order,
28 std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>&&
29 process_variables,
31 SecondaryVariableCollection&& secondary_variables,
32 bool const use_monolithic_scheme, bool const is_linear)
33 : Process(std::move(name), mesh, std::move(jacobian_assembler), parameters,
34 integration_order, std::move(process_variables),
35 std::move(secondary_variables), use_monolithic_scheme),
37 *_jacobian_assembler, is_linear, use_monolithic_scheme},
38 _process_data(std::move(process_data))
39
40{
41 // For numerical Jacobian
42 if (this->_jacobian_assembler->isPerturbationEnabled())
43 {
45 "Numerical Jacobian is not supported for the "
46 "ThermoHydroMechanicsProcess yet.");
47 }
48
49 _integration_point_writer.emplace_back(
50 std::make_unique<MeshLib::IntegrationPointWriter>(
51 "sigma_ip",
52 static_cast<int>(mesh.getDimension() == 2 ? 4 : 6) /*n components*/,
53 integration_order, local_assemblers_,
55
56 _integration_point_writer.emplace_back(
57 std::make_unique<MeshLib::IntegrationPointWriter>(
58 "epsilon_m_ip",
59 static_cast<int>(mesh.getDimension() == 2 ? 4 : 6) /*n components*/,
60 integration_order, local_assemblers_,
62
63 _integration_point_writer.emplace_back(
64 std::make_unique<MeshLib::IntegrationPointWriter>(
65 "epsilon_ip",
66 static_cast<int>(mesh.getDimension() == 2 ? 4 : 6) /*n components*/,
67 integration_order, local_assemblers_,
69}
70
71template <int DisplacementDim>
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 || process_id == 2)
86 {
87 auto const& l = *_local_to_global_index_map;
88 return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
89 &l.getGhostIndices(), &this->_sparsity_pattern};
90 }
91
92 // For staggered scheme and T or H process (pressure).
94 return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
95 &l.getGhostIndices(), &_sparsity_pattern_with_linear_element};
96}
97
98template <int DisplacementDim>
100{
101 // Create single component dof in every of the mesh's nodes.
103 std::make_unique<MeshLib::MeshSubset>(_mesh, _mesh.getNodes());
104 // Create single component dof in the mesh's base nodes.
105 _base_nodes = MeshLib::getBaseNodes(_mesh.getElements());
107 std::make_unique<MeshLib::MeshSubset>(_mesh, _base_nodes);
108
109 // TODO move the two data members somewhere else.
110 // for extrapolation of secondary variables of stress or strain
111 std::vector<MeshLib::MeshSubset> all_mesh_subsets_single_component{
114 std::make_unique<NumLib::LocalToGlobalIndexMap>(
115 std::move(all_mesh_subsets_single_component),
116 // by location order is needed for output
118
120 {
121 // For temperature, which is the first
122 std::vector<MeshLib::MeshSubset> all_mesh_subsets{
124
125 // For pressure, which is the second
126 all_mesh_subsets.push_back(*_mesh_subset_base_nodes);
127
128 // For displacement.
129 const int monolithic_process_id = 0;
130 std::generate_n(std::back_inserter(all_mesh_subsets),
131 getProcessVariables(monolithic_process_id)[2]
132 .get()
133 .getNumberOfGlobalComponents(),
134 [&]() { return *_mesh_subset_all_nodes; });
135
136 std::vector<int> const vec_n_components{1, 1, DisplacementDim};
138 std::make_unique<NumLib::LocalToGlobalIndexMap>(
139 std::move(all_mesh_subsets), vec_n_components,
142 }
143 else
144 {
145 // For displacement equation.
146 const int process_id = 2;
147 std::vector<MeshLib::MeshSubset> all_mesh_subsets;
148 std::generate_n(std::back_inserter(all_mesh_subsets),
149 getProcessVariables(process_id)[0]
150 .get()
151 .getNumberOfGlobalComponents(),
152 [&]() { return *_mesh_subset_all_nodes; });
153
154 std::vector<int> const vec_n_components{DisplacementDim};
156 std::make_unique<NumLib::LocalToGlobalIndexMap>(
157 std::move(all_mesh_subsets), vec_n_components,
159
160 // For pressure equation or temperature equation.
161 // Collect the mesh subsets with base nodes in a vector.
162 std::vector<MeshLib::MeshSubset> all_mesh_subsets_base_nodes{
165 std::make_unique<NumLib::LocalToGlobalIndexMap>(
166 std::move(all_mesh_subsets_base_nodes),
167 // by location order is needed for output
169
172
175 }
176}
177
178template <int DisplacementDim>
180 NumLib::LocalToGlobalIndexMap const& dof_table,
181 MeshLib::Mesh const& mesh,
182 unsigned const integration_order)
183{
186 mesh.getElements(), dof_table, local_assemblers_,
187 NumLib::IntegrationOrder{integration_order}, mesh.isAxiallySymmetric(),
189
190 auto add_secondary_variable = [&](std::string const& name,
191 int const num_components,
192 auto get_ip_values_function)
193 {
194 _secondary_variables.addSecondaryVariable(
195 name,
196 makeExtrapolator(num_components, getExtrapolator(),
198 std::move(get_ip_values_function)));
199 };
200
201 add_secondary_variable(
202 "sigma",
204 DisplacementDim>::RowsAtCompileTime,
206
207 add_secondary_variable(
208 "sigma_ice",
210 DisplacementDim>::RowsAtCompileTime,
212
213 add_secondary_variable(
214 "epsilon_0",
216 DisplacementDim>::RowsAtCompileTime,
218
219 add_secondary_variable(
220 "epsilon_m",
222 DisplacementDim>::RowsAtCompileTime,
224
225 add_secondary_variable(
226 "epsilon",
228 DisplacementDim>::RowsAtCompileTime,
230
231 add_secondary_variable(
232 "ice_volume_fraction", 1,
234
235 add_secondary_variable(
236 "velocity", mesh.getDimension(),
238
239 add_secondary_variable(
240 "fluid_density", 1,
242
243 add_secondary_variable(
244 "viscosity", 1,
246
248 const_cast<MeshLib::Mesh&>(mesh), "ice_volume_fraction_avg",
250
251 _process_data.element_fluid_density =
253 const_cast<MeshLib::Mesh&>(mesh), "fluid_density_avg",
255
257 const_cast<MeshLib::Mesh&>(mesh), "viscosity_avg",
259
261 const_cast<MeshLib::Mesh&>(mesh), "stress_avg",
264 DisplacementDim>::RowsAtCompileTime);
265
266 _process_data.pressure_interpolated =
268 const_cast<MeshLib::Mesh&>(mesh), "pressure_interpolated",
270
271 _process_data.temperature_interpolated =
273 const_cast<MeshLib::Mesh&>(mesh), "temperature_interpolated",
275
276 //
277 // enable output of internal variables defined by material models
278 //
281 add_secondary_variable);
282
285 _process_data.solid_materials, local_assemblers_,
286 _integration_point_writer, integration_order);
287
290
291 // Initialize local assemblers after all variables have been set.
295}
296
297template <int DisplacementDim>
299 std::map<int, std::shared_ptr<MaterialPropertyLib::Medium>> const& media)
300{
302 {
303 const int process_id_of_thermohydromechancs = 0;
305 *_local_to_global_index_map, process_id_of_thermohydromechancs,
306 media);
307 return;
308 }
309
310 // Staggered scheme:
311 // for the equations of heat transport
312 const int thermal_process_id = 0;
314 *_local_to_global_index_map_with_base_nodes, thermal_process_id, media);
315
316 // for the equations of mass balance
317 const int hydraulic_process_id = 1;
319 *_local_to_global_index_map_with_base_nodes, hydraulic_process_id,
320 media);
321
322 // for the equations of deformation.
323 const int mechanical_process_id = 2;
325 *_local_to_global_index_map, mechanical_process_id, media);
326}
327
328template <int DisplacementDim>
330 setInitialConditionsConcreteProcess(std::vector<GlobalVector*>& x,
331 double const t,
332 int const process_id)
333{
334 DBUG("SetInitialConditions ThermoHydroMechanicsProcess.");
335
338 local_assemblers_, getDOFTables(x.size()), x, t, process_id);
339}
340
341template <int DisplacementDim>
343 const double t, double const dt, std::vector<GlobalVector*> const& x,
344 std::vector<GlobalVector*> const& x_prev, int const process_id,
346{
347 DBUG("Assemble the equations for ThermoHydroMechanics");
348
350 t, dt, x, x_prev, process_id, M, K, 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
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
391 assembleWithJacobian(t, dt, x, x_prev, process_id, b, Jac);
392}
393
394template <int DisplacementDim>
396 std::vector<GlobalVector*> const& x, double const t, double const dt,
397 const int process_id)
398{
399 DBUG("PreTimestep ThermoHydroMechanicsProcess.");
400
401 if (hasMechanicalProcess(process_id))
402 {
406 dt);
407
409 updateActiveElements();
410 }
411}
412
413template <int DisplacementDim>
415 std::vector<GlobalVector*> const& x,
416 std::vector<GlobalVector*> const& x_prev, double const t, double const dt,
417 const int process_id)
418{
419 if (process_id != 0)
420 {
421 return;
422 }
423
424 DBUG("PostTimestep ThermoHydroMechanicsProcess.");
425
429 x_prev, t, dt, process_id);
430}
431
432template <int DisplacementDim>
433std::vector<std::vector<std::string>>
435 std::vector<std::reference_wrapper<MeshLib::Mesh>> const& meshes)
436{
437 INFO("ThermoHydroMechanicsProcess process initializeSubmeshOutput().");
438 std::vector<std::vector<std::string>> per_process_residuum_names;
439 if (_process_variables.size() == 1) // monolithic
440 {
441 per_process_residuum_names = {
442 {"HeatFlowRate", "VolumetricFlowRate", "NodalForces"}};
443 }
444 else // staggered
445 {
446 per_process_residuum_names = {
447 {"HeatFlowRate"}, {"VolumetricFlowRate"}, {"NodalForces"}};
448 }
449
451 initializeAssemblyOnSubmeshes(meshes, per_process_residuum_names);
452
453 return per_process_residuum_names;
454}
455
456template <int DisplacementDim>
458 computeSecondaryVariableConcrete(double const t, double const dt,
459 std::vector<GlobalVector*> const& x,
460 GlobalVector const& x_prev,
461 const int process_id)
462{
463 if (process_id != 0)
464 {
465 return;
466 }
467
468 DBUG("Compute the secondary variables for ThermoHydroMechanicsProcess.");
469
473 x, x_prev, process_id);
474}
475
476template <int DisplacementDim>
477std::tuple<NumLib::LocalToGlobalIndexMap*, bool> ThermoHydroMechanicsProcess<
478 DisplacementDim>::getDOFTableForExtrapolatorData() const
479{
480 const bool manage_storage = false;
481 return std::make_tuple(_local_to_global_index_map_single_component.get(),
482 manage_storage);
483}
484
485template <int DisplacementDim>
488 const int process_id) const
489{
490 if (hasMechanicalProcess(process_id))
491 {
493 }
494
495 // For the equation of pressure
497}
498
499template class ThermoHydroMechanicsProcess<2>;
500template class ThermoHydroMechanicsProcess<3>;
501
502} // namespace ThermoHydroMechanics
503} // namespace ProcessLib
#define OGS_FATAL(...)
Definition Error.h:19
MathLib::EigenMatrix GlobalMatrix
MathLib::EigenVector GlobalVector
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:28
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:22
bool isAxiallySymmetric() const
Definition Mesh.h:128
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
Definition Mesh.h:100
unsigned getDimension() const
Returns the dimension of the mesh (determined by the maximum dimension over all elements).
Definition Mesh.h:79
Properties & getProperties()
Definition Mesh.h:125
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 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 initialize(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
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)
std::string const name
Definition Process.h:361
std::vector< NumLib::LocalToGlobalIndexMap const * > getDOFTables(int const number_of_processes) const
Definition Process.cpp:376
std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > _integration_point_writer
Definition Process.h:389
void initializeProcessBoundaryConditionsAndSourceTerms(const NumLib::LocalToGlobalIndexMap &dof_table, const int process_id, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
Definition Process.cpp:83
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
Definition Process.cpp:258
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:37
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_all_nodes
Definition Process.h:365
MeshLib::Mesh & _mesh
Definition Process.h:364
std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > _process_variables
Definition Process.h:399
std::vector< std::size_t > const & getActiveElementIDs() const
Definition Process.h:160
SecondaryVariableCollection _secondary_variables
Definition Process.h:369
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
Definition Process.h:367
std::unique_ptr< ProcessLib::AbstractJacobianAssembler > _jacobian_assembler
Definition Process.h:375
std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > const & getProcessVariables() const
Definition Process.h:149
NumLib::Extrapolator & getExtrapolator() const
Definition Process.h:201
GlobalSparsityPattern _sparsity_pattern
Definition Process.h:391
const bool _use_monolithic_scheme
Definition Process.h:378
Handles configuration of several secondary variables from the project file.
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map_with_base_nodes
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
THERMOHYDROMECHANICS_EXPORT NumLib::LocalToGlobalIndexMap const & getDOFTable(const int process_id) const override
THERMOHYDROMECHANICS_EXPORT 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, bool const is_linear)
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
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map_single_component
void initializeBoundaryConditions(std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media) override
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
std::tuple< NumLib::LocalToGlobalIndexMap *, bool > getDOFTableForExtrapolatorData() const 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
THERMOHYDROMECHANICS_EXPORT MathLib::MatrixSpecifications getMatrixSpecifications(const int process_id) const override
void setInitialConditionsConcreteProcess(std::vector< GlobalVector * > &x, double const t, int const process_id) override
std::vector< std::unique_ptr< LocalAssemblerInterface< DisplacementDim > > > local_assemblers_
ThermoHydroMechanicsProcessData< DisplacementDim > _process_data
std::vector< std::vector< std::string > > initializeAssemblyOnSubmeshes(std::vector< std::reference_wrapper< MeshLib::Mesh > > const &meshes) override
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)
@ 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 executeMemberOnDereferenced(Method method, Container const &container, Args &&... args)
virtual std::vector< double > const & getIntPtFluidDensity(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 & 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 & getIntPtEpsilon0(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 & getIntPtDarcyVelocity(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 & getIntPtSigmaIce(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 & getIntPtIceVolume(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 & getIntPtViscosity(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 & getIntPtEpsilonM(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