OGS
HTProcess.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
4#include "HTProcess.h"
5
6#include <cassert>
7
8#include "MonolithicHTFEM.h"
14#include "StaggeredHTFEM.h"
15
16namespace ProcessLib
17{
18namespace HT
19{
21 std::string name,
22 MeshLib::Mesh& mesh,
23 std::unique_ptr<ProcessLib::AbstractJacobianAssembler>&& jacobian_assembler,
24 std::vector<std::unique_ptr<ParameterLib::ParameterBase>> const& parameters,
25 unsigned const integration_order,
26 std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>&&
27 process_variables,
28 HTProcessData&& process_data,
29 SecondaryVariableCollection&& secondary_variables,
30 bool const use_monolithic_scheme,
31 std::unique_ptr<ProcessLib::SurfaceFluxData>&& surfaceflux)
32 : Process(std::move(name), mesh, std::move(jacobian_assembler), parameters,
33 integration_order, std::move(process_variables),
34 std::move(secondary_variables), use_monolithic_scheme),
35 _process_data(std::move(process_data)),
36 _surfaceflux(std::move(surfaceflux))
37{
38 this->_jacobian_assembler->checkPerturbationSize(2);
40 {
41 this->_jacobian_assembler->setNonDeformationComponentIDsNoSizeCheck(
42 {0, 1} /* two variables: pressure and temperature */);
43 }
44}
45
47 NumLib::LocalToGlobalIndexMap const& dof_table,
48 MeshLib::Mesh const& mesh,
49 unsigned const integration_order)
50{
51 int const mesh_space_dimension = _process_data.mesh_space_dimension;
52
54 {
56 mesh_space_dimension, mesh.getElements(), dof_table,
59 }
60 else
61 {
63 mesh_space_dimension, mesh.getElements(), dof_table,
66 }
67
68 _secondary_variables.addSecondaryVariable(
69 "darcy_velocity",
70 makeExtrapolator(mesh_space_dimension, getExtrapolator(),
73}
74
76 const double t, double const dt, std::vector<GlobalVector*> const& x,
77 std::vector<GlobalVector*> const& x_prev, int const process_id,
79{
80 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
82 {
83 DBUG("Assemble HTProcess.");
84 dof_tables.emplace_back(_local_to_global_index_map.get());
85 }
86 else
87 {
88 if (process_id == _process_data.heat_transport_process_id)
89 {
90 DBUG(
91 "Assemble the equations of heat transport process within "
92 "HTProcess.");
93 }
94 else
95 {
96 DBUG(
97 "Assemble the equations of single phase fully saturated "
98 "fluid flow process within HTProcess.");
99 }
100 dof_tables.emplace_back(_local_to_global_index_map.get());
101 dof_tables.emplace_back(_local_to_global_index_map.get());
102
103 // For numerical Jacobian assembler
104 // (only one variable per process in staggered scheme);
105 this->_jacobian_assembler->setNonDeformationComponentIDsNoSizeCheck(
106 {process_id});
107 }
108
109 // Call global assembler for each local assembly item.
112 getActiveElementIDs(), dof_tables, t, dt, x, x_prev, process_id, &M, &K,
113 &b);
114}
115
117 const double t, double const dt, std::vector<GlobalVector*> const& x,
118 std::vector<GlobalVector*> const& x_prev, int const process_id,
119 GlobalVector& b, GlobalMatrix& Jac)
120{
121 DBUG("AssembleWithJacobian HTProcess.");
122
123 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
125 {
126 dof_tables.emplace_back(_local_to_global_index_map.get());
127 }
128 else
129 {
130 dof_tables.emplace_back(_local_to_global_index_map.get());
131 dof_tables.emplace_back(_local_to_global_index_map.get());
132 }
133
134 // Call global assembler for each local assembly item.
137 _local_assemblers, getActiveElementIDs(), dof_tables, t, dt, x, x_prev,
138 process_id, &b, &Jac);
139}
140
141std::tuple<NumLib::LocalToGlobalIndexMap*, bool>
143{
145 {
146 // For single-variable-single-component processes reuse the existing DOF
147 // table.
148 const bool manage_storage = false;
149 return std::make_tuple(_local_to_global_index_map.get(),
150 manage_storage);
151 }
152
153 // Otherwise construct a new DOF table.
154 std::vector<MeshLib::MeshSubset> all_mesh_subsets_single_component{
156
157 const bool manage_storage = true;
158 return std::make_tuple(new NumLib::LocalToGlobalIndexMap(
159 std::move(all_mesh_subsets_single_component),
160 // by location order is needed for output
162 manage_storage);
163}
164
165Eigen::Vector3d HTProcess::getFlux(std::size_t element_id,
166 MathLib::Point3d const& p,
167 double const t,
168 std::vector<GlobalVector*> const& x) const
169{
170 // fetch local_x from primary variable
171 std::vector<GlobalIndexType> indices_cache;
172 auto const r_c_indices = NumLib::getRowColumnIndices(
173 element_id, *_local_to_global_index_map, indices_cache);
174 std::vector<std::vector<GlobalIndexType>> indices_of_all_coupled_processes{
175 x.size(), r_c_indices.rows};
176 auto const local_x =
177 getCoupledLocalSolutions(x, indices_of_all_coupled_processes);
178
179 return _local_assemblers[element_id]->getFlux(p, t, local_x);
180}
181
182// this is almost a copy of the implementation in the GroundwaterFlow
184 std::vector<GlobalVector*> const& x,
185 std::vector<GlobalVector*> const& /*x_prev*/,
186 const double t,
187 const double /*delta_t*/,
188 int const process_id)
189{
190 // For the monolithic scheme, process_id is always zero.
191 if (_use_monolithic_scheme && process_id != 0)
192 {
193 OGS_FATAL(
194 "The condition of process_id = 0 must be satisfied for monolithic "
195 "HTProcess, which is a single process.");
196 }
198 process_id != _process_data.hydraulic_process_id)
199 {
200 DBUG("This is the thermal part of the staggered HTProcess.");
201 return;
202 }
203 if (!_surfaceflux) // computing the surfaceflux is optional
204 {
205 return;
206 }
207
208 _surfaceflux->integrate(x, t, *this, process_id, _integration_order, _mesh,
210}
211} // namespace HT
212} // namespace ProcessLib
#define OGS_FATAL(...)
Definition Error.h:19
MathLib::EigenMatrix GlobalMatrix
MathLib::EigenVector GlobalVector
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
virtual std::vector< double > const & getIntPtDarcyVelocity(const double, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &) const =0
void initializeConcreteProcess(NumLib::LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh, unsigned const integration_order) override
Process specific initialization called by initialize().
Definition HTProcess.cpp:46
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
std::tuple< NumLib::LocalToGlobalIndexMap *, bool > getDOFTableForExtrapolatorData() const override
HTProcess(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, HTProcessData &&process_data, SecondaryVariableCollection &&secondary_variables, bool const use_monolithic_scheme, std::unique_ptr< ProcessLib::SurfaceFluxData > &&surfaceflux)
Definition HTProcess.cpp:20
std::unique_ptr< ProcessLib::SurfaceFluxData > _surfaceflux
Definition HTProcess.h:106
Eigen::Vector3d getFlux(std::size_t element_id, MathLib::Point3d const &p, double const t, std::vector< GlobalVector * > const &x) const override
HTProcessData _process_data
Definition HTProcess.h:102
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
Definition HTProcess.cpp:75
void postTimestepConcreteProcess(std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, const double t, const double delta_t, int const process_id) override
std::vector< std::unique_ptr< HTLocalAssemblerInterface > > _local_assemblers
Definition HTProcess.h:104
std::string const name
Definition Process.h:361
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::size_t > const & getActiveElementIDs() const
Definition Process.h:160
SecondaryVariableCollection _secondary_variables
Definition Process.h:369
VectorMatrixAssembler _global_assembler
Definition Process.h:376
unsigned const _integration_order
Definition Process.h:383
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
NumLib::Extrapolator & getExtrapolator() const
Definition Process.h:201
const bool _use_monolithic_scheme
Definition Process.h:378
Handles configuration of several secondary variables from the project file.
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)
NumLib::LocalToGlobalIndexMap::RowColumnIndices getRowColumnIndices(std::size_t const id, NumLib::LocalToGlobalIndexMap const &dof_table, std::vector< GlobalIndexType > &indices)
@ BY_LOCATION
Ordering data by spatial location.
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)
std::vector< double > getCoupledLocalSolutions(std::vector< GlobalVector * > const &global_solutions, std::vector< std::vector< GlobalIndexType > > const &indices)
SecondaryVariableFunctions makeExtrapolator(const unsigned num_components, NumLib::Extrapolator &extrapolator, LocalAssemblerCollection const &local_assemblers, typename NumLib::ExtrapolatableLocalAssemblerCollection< LocalAssemblerCollection >::IntegrationPointValuesMethod integration_point_values_method)
static void executeSelectedMemberDereferenced(Object &object, Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)