OGS
HydroMechanics/HydroMechanicsProcess.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
8#include "HydroMechanicsFEM.h"
14#include "ProcessLib/Process.h"
17
18namespace ProcessLib
19{
20namespace HydroMechanics
21{
22template <int DisplacementDim>
24 std::string name,
25 MeshLib::Mesh& mesh,
26 std::unique_ptr<ProcessLib::AbstractJacobianAssembler>&& jacobian_assembler,
27 std::vector<std::unique_ptr<ParameterLib::ParameterBase>> const& parameters,
28 unsigned const integration_order,
29 std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>&&
30 process_variables,
32 SecondaryVariableCollection&& secondary_variables,
33 bool const use_monolithic_scheme,
34 bool const is_linear)
35 : Process(std::move(name), mesh, std::move(jacobian_assembler), parameters,
36 integration_order, std::move(process_variables),
37 std::move(secondary_variables), use_monolithic_scheme),
38 AssemblyMixin<HydroMechanicsProcess<DisplacementDim>>{
39 *_jacobian_assembler, is_linear, use_monolithic_scheme},
40 process_data_(std::move(process_data))
41{
42 // For numerical Jacobian
43 if (this->_jacobian_assembler->isPerturbationEnabled() &&
44 _use_monolithic_scheme)
45 {
47 "Numerical Jacobian is not supported for the "
48 "HydroMechanicsProcess using the monolithic scheme yet.");
49 }
50
51 _integration_point_writer.emplace_back(
52 std::make_unique<MeshLib::IntegrationPointWriter>(
53 "sigma_ip",
54 static_cast<int>(mesh.getDimension() == 2 ? 4 : 6) /*n components*/,
55 integration_order, local_assemblers_, &LocalAssemblerIF::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 &LocalAssemblerIF::getEpsilon));
63
64 if (!_use_monolithic_scheme)
65 {
66 _integration_point_writer.emplace_back(
67 std::make_unique<MeshLib::IntegrationPointWriter>(
68 "strain_rate_variable_ip", 1, integration_order,
69 local_assemblers_, &LocalAssemblerIF::getStrainRateVariable));
70 }
71}
72
73template <int DisplacementDim>
78
79template <int DisplacementDim>
82 const int process_id) const
83{
84 // For the monolithic scheme or the M process (deformation) in the staggered
85 // scheme.
86 if (process_id == process_data_.mechanics_related_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 H process (pressure).
95 return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
96 &l.getGhostIndices(), &sparsity_pattern_with_linear_element_};
97}
98
99template <int DisplacementDim>
101{
102 // Create single component dof in every of the mesh's nodes.
103 _mesh_subset_all_nodes = std::make_unique<MeshLib::MeshSubset>(
104 _mesh, _mesh.getNodes(), process_data_.use_taylor_hood_elements);
105
106 // Create single component dof in the mesh's base nodes.
107 base_nodes_ = MeshLib::getBaseNodes(_mesh.getElements());
108 mesh_subset_base_nodes_ = std::make_unique<MeshLib::MeshSubset>(
109 _mesh, base_nodes_, process_data_.use_taylor_hood_elements);
110
111 // TODO move the two data members somewhere else.
112 // for extrapolation of secondary variables of stress or strain
113 std::vector<MeshLib::MeshSubset> all_mesh_subsets_single_component{
116 std::make_unique<NumLib::LocalToGlobalIndexMap>(
117 std::move(all_mesh_subsets_single_component),
118 // by location order is needed for output
120
121 if (process_data_.isMonolithicSchemeUsed())
122 {
123 // For pressure, which is the first
124 std::vector<MeshLib::MeshSubset> all_mesh_subsets{
126
127 // For displacement.
128 const int monolithic_process_id = 0;
129 std::generate_n(std::back_inserter(all_mesh_subsets),
130 getProcessVariables(monolithic_process_id)[1]
131 .get()
132 .getNumberOfGlobalComponents(),
133 [&]() { return *_mesh_subset_all_nodes; });
134
135 std::vector<int> const vec_n_components{1, DisplacementDim};
137 std::make_unique<NumLib::LocalToGlobalIndexMap>(
138 std::move(all_mesh_subsets), vec_n_components,
141 }
142 else
143 {
144 // For displacement equation.
145 const int process_id = 1;
146 std::vector<MeshLib::MeshSubset> all_mesh_subsets;
147 std::generate_n(std::back_inserter(all_mesh_subsets),
148 getProcessVariables(process_id)[0]
149 .get()
150 .getNumberOfGlobalComponents(),
151 [&]() { return *_mesh_subset_all_nodes; });
152
153 std::vector<int> const vec_n_components{DisplacementDim};
155 std::make_unique<NumLib::LocalToGlobalIndexMap>(
156 std::move(all_mesh_subsets), vec_n_components,
158
159 // For pressure equation.
160 // Collect the mesh subsets with base nodes in a vector.
161 std::vector<MeshLib::MeshSubset> all_mesh_subsets_base_nodes{
164 std::make_unique<NumLib::LocalToGlobalIndexMap>(
165 std::move(all_mesh_subsets_base_nodes),
166 // by location order is needed for output
168
171
174 }
175}
176
177template <int DisplacementDim>
179 NumLib::LocalToGlobalIndexMap const& dof_table,
180 MeshLib::Mesh const& mesh,
181 unsigned const integration_order)
182{
185 mesh.getElements(), dof_table, local_assemblers_,
186 NumLib::IntegrationOrder{integration_order}, mesh.isAxiallySymmetric(),
188
189 auto add_secondary_variable = [&](std::string const& name,
190 int const num_components,
191 auto get_ip_values_function)
192 {
193 _secondary_variables.addSecondaryVariable(
194 name,
195 makeExtrapolator(num_components, getExtrapolator(),
197 std::move(get_ip_values_function)));
198 };
199
200 add_secondary_variable("sigma",
202 DisplacementDim>::RowsAtCompileTime,
204
205 add_secondary_variable("epsilon",
207 DisplacementDim>::RowsAtCompileTime,
209
210 add_secondary_variable("velocity", DisplacementDim,
212
213 //
214 // enable output of internal variables defined by material models
215 //
217 LocalAssemblerIF>(process_data_.solid_materials,
218 add_secondary_variable);
219
220 process_data_.pressure_interpolated =
222 const_cast<MeshLib::Mesh&>(mesh), "pressure_interpolated",
224
225 process_data_.principal_stress_vector[0] =
227 const_cast<MeshLib::Mesh&>(mesh), "principal_stress_vector_1",
229
230 process_data_.principal_stress_vector[1] =
232 const_cast<MeshLib::Mesh&>(mesh), "principal_stress_vector_2",
234
235 process_data_.principal_stress_vector[2] =
237 const_cast<MeshLib::Mesh&>(mesh), "principal_stress_vector_3",
239
240 process_data_.principal_stress_values =
242 const_cast<MeshLib::Mesh&>(mesh), "principal_stress_values",
244
246 const_cast<MeshLib::Mesh&>(mesh), "permeability",
249
252
253 // Initialize local assemblers after all variables have been set.
257}
258
259template <int DisplacementDim>
261 std::map<int, std::shared_ptr<MaterialPropertyLib::Medium>> const& media)
262{
263 if (process_data_.isMonolithicSchemeUsed())
264 {
265 const int process_id_of_hydromechanics = 0;
267 *_local_to_global_index_map, process_id_of_hydromechanics, media);
268 return;
269 }
270
271 // Staggered scheme:
272 // for the equations of pressure
273 const int hydraulic_process_id = 0;
275 *local_to_global_index_map_with_base_nodes_, hydraulic_process_id,
276 media);
277
278 // for the equations of deformation.
279 const int mechanical_process_id = 1;
281 *_local_to_global_index_map, mechanical_process_id, media);
282}
283
284template <int DisplacementDim>
286 const double t, double const dt, std::vector<GlobalVector*> const& x,
287 std::vector<GlobalVector*> const& x_prev, int const process_id,
289{
290 DBUG("Assemble the equations for HydroMechanics");
291
292 // Note: This assembly function is for the Picard nonlinear solver. Since
293 // only the Newton-Raphson method is employed to simulate coupled HM
294 // processes in this class, this function is actually not used so far.
295
296 if (this->_jacobian_assembler->isPerturbationEnabled() &&
298 {
299 if (process_id == process_data_.hydraulic_process_id)
300 {
301 this->_jacobian_assembler->setNonDeformationComponentIDs(
302 {process_id} /* pressure variable id */);
303 }
304 else
305 {
306 OGS_FATAL(
307 "Numerical Jacobian is only supported for the "
308 "liquid fluid process in the staggered HydroMechanicsProcess.");
309 }
310 }
311
313 t, dt, x, x_prev, process_id, M, K, b);
314}
315
316template <int DisplacementDim>
319 const double t, double const dt, std::vector<GlobalVector*> const& x,
320 std::vector<GlobalVector*> const& x_prev, int const process_id,
321 GlobalVector& b, GlobalMatrix& Jac)
322{
323 // For the monolithic scheme
324 bool const use_monolithic_scheme = process_data_.isMonolithicSchemeUsed();
325 if (use_monolithic_scheme)
326 {
327 DBUG(
328 "Assemble the Jacobian of HydroMechanics for the monolithic "
329 "scheme.");
330 }
331 else
332 {
333 // For the staggered scheme
334 if (process_id == process_data_.hydraulic_process_id)
335 {
336 DBUG(
337 "Assemble the Jacobian equations of liquid fluid process in "
338 "HydroMechanics for the staggered scheme.");
339 }
340 else
341 {
342 DBUG(
343 "Assemble the Jacobian equations of mechanical process in "
344 "HydroMechanics for the staggered scheme.");
345 }
346 }
347
349 t, dt, x, x_prev, process_id, b, Jac);
350}
351
352template <int DisplacementDim>
354 std::vector<GlobalVector*> const& x, double const t, double const dt,
355 const int process_id)
356{
357 DBUG("PreTimestep HydroMechanicsProcess.");
358
359 if (hasMechanicalProcess(process_id))
360 {
364 t, dt);
365
368 }
369}
370
371template <int DisplacementDim>
372std::vector<std::vector<std::string>>
374 std::vector<std::reference_wrapper<MeshLib::Mesh>> const& meshes)
375{
376 INFO("HydroMechanicsProcess initializeSubmeshOutput().");
377 std::vector<std::vector<std::string>> per_process_residuum_names;
378 if (_process_variables.size() == 1) // monolithic
379 {
380 per_process_residuum_names = {{"MassFlowRate", "NodalForces"}};
381 }
382 else // staggered
383 {
384 per_process_residuum_names = {{"MassFlowRate"}, {"NodalForces"}};
385 }
386
388 initializeAssemblyOnSubmeshes(meshes, per_process_residuum_names);
389
390 return per_process_residuum_names;
391}
392
393template <int DisplacementDim>
395 std::vector<GlobalVector*> const& x,
396 std::vector<GlobalVector*> const& x_prev, double const t, double const dt,
397 const int process_id)
398{
399 if (process_id != process_data_.hydraulic_process_id)
400 {
401 return;
402 }
403
404 DBUG("PostTimestep HydroMechanicsProcess.");
405
408 getActiveElementIDs(), getDOFTables(x.size()), x, x_prev, t, dt,
409 process_id);
410}
411
412template <int DisplacementDim>
414 std::vector<GlobalVector*> const& x,
415 std::vector<GlobalVector*> const& x_prev, const double t, double const dt,
416 const int process_id)
417{
418 DBUG("PostNonLinearSolver HydroMechanicsProcess.");
419
420 // Calculate strain, stress or other internal variables of mechanics.
423 getActiveElementIDs(), getDOFTables(x.size()), x, x_prev, t, dt,
424 process_id);
425}
426
427template <int DisplacementDim>
429 setInitialConditionsConcreteProcess(std::vector<GlobalVector*>& x,
430 double const t,
431 int const process_id)
432{
433 // So far, this function only sets the initial stress using the input data.
434 if (process_id != process_data_.mechanics_related_process_id)
435 {
436 return;
437 }
438
439 DBUG("Set initial conditions of HydroMechanicsProcess.");
440
443 getActiveElementIDs(), getDOFTables(x.size()), x, t, process_id);
444}
445
446template <int DisplacementDim>
448 double const t, double const dt, std::vector<GlobalVector*> const& x,
449 GlobalVector const& x_prev, const int process_id)
450{
451 if (process_id != process_data_.hydraulic_process_id)
452 {
453 return;
454 }
455
456 DBUG("Compute the secondary variables for HydroMechanicsProcess.");
457
460 getActiveElementIDs(), getDOFTables(x.size()), t, dt, x, x_prev,
461 process_id);
462}
463
464template <int DisplacementDim>
465std::tuple<NumLib::LocalToGlobalIndexMap*, bool>
467{
468 const bool manage_storage = false;
469 return std::make_tuple(local_to_global_index_map_single_component_.get(),
470 manage_storage);
471}
472
473template <int DisplacementDim>
476{
477 if (hasMechanicalProcess(process_id))
478 {
480 }
481
482 // For the equation of pressure
484}
485
486template class HydroMechanicsProcess<2>;
487template class HydroMechanicsProcess<3>;
488
489} // namespace HydroMechanics
490} // 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
Properties & getProperties()
Definition Mesh.h:125
std::unique_ptr< MeshLib::MeshSubset const > mesh_subset_base_nodes_
void assembleConcreteProcess(const double t, double const, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b) override
std::unique_ptr< NumLib::LocalToGlobalIndexMap > local_to_global_index_map_with_base_nodes_
std::vector< std::vector< std::string > > initializeAssemblyOnSubmeshes(std::vector< std::reference_wrapper< MeshLib::Mesh > > const &meshes) override
HydroMechanicsProcess(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, HydroMechanicsProcessData< DisplacementDim > &&process_data, SecondaryVariableCollection &&secondary_variables, bool const use_monolithic_scheme, bool const is_linear)
void computeSecondaryVariableConcrete(double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_prev, const int process_id) override
void preTimestepConcreteProcess(std::vector< GlobalVector * > const &x, double const t, double const dt, const int process_id) override
void initializeBoundaryConditions(std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media) override
NumLib::LocalToGlobalIndexMap const & getDOFTable(const int process_id) const override
void assembleWithJacobianConcreteProcess(const double t, double const, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalVector &b, GlobalMatrix &Jac) override
std::unique_ptr< NumLib::LocalToGlobalIndexMap > local_to_global_index_map_single_component_
void postNonLinearSolverConcreteProcess(std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, const double t, double const dt, int const process_id) override
void setInitialConditionsConcreteProcess(std::vector< GlobalVector * > &x, double const t, int const process_id) override
MathLib::MatrixSpecifications getMatrixSpecifications(const int process_id) const override
std::vector< std::unique_ptr< LocalAssemblerIF > > local_assemblers_
void initializeConcreteProcess(NumLib::LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh, unsigned const integration_order) override
Process specific initialization called by initialize().
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
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 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
void assembleWithJacobian(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) final
Definition Process.cpp:279
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.
constexpr int kelvin_vector_dimensions(int const displacement_dim)
Kelvin vector dimensions for given displacement dimension.
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 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 & getIntPtSigma(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 & getIntPtEpsilon(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0