OGS
StokesFlowProcess.cpp
Go to the documentation of this file.
1
11#include "StokesFlowProcess.h"
12
13#include <cassert>
14
18
19namespace ProcessLib
20{
21namespace StokesFlow
22{
23template <int GlobalDim>
25 std::string name,
26 MeshLib::Mesh& mesh,
27 std::unique_ptr<ProcessLib::AbstractJacobianAssembler>&& jacobian_assembler,
28 std::vector<std::unique_ptr<ParameterLib::ParameterBase>> const& parameters,
29 unsigned const integration_order,
30 std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>&&
31 process_variables,
32 StokesFlowProcessData&& process_data,
33 SecondaryVariableCollection&& secondary_variables,
34 bool const use_monolithic_scheme)
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 _process_data(std::move(process_data))
39{
40}
41
42template <int GlobalDim>
44{
45 // Create single component dof in every of the mesh's nodes.
46 _mesh_subset_all_nodes =
47 std::make_unique<MeshLib::MeshSubset>(_mesh, _mesh.getNodes());
48 // Create single component dof in the mesh's base nodes.
49 _base_nodes = MeshLib::getBaseNodes(_mesh.getElements());
50 _mesh_subset_base_nodes =
51 std::make_unique<MeshLib::MeshSubset>(_mesh, _base_nodes);
52
53 // TODO move the two data members somewhere else.
54 // for extrapolation of secondary variables
55 std::vector<MeshLib::MeshSubset> all_mesh_subsets_single_component{
56 *_mesh_subset_all_nodes};
57 _local_to_global_index_map_single_component =
58 std::make_unique<NumLib::LocalToGlobalIndexMap>(
59 std::move(all_mesh_subsets_single_component),
60 // by location order is needed for output
62
63 assert(_use_monolithic_scheme);
64 {
65 // For vector variables, in this case liquid velocity.
66 int const process_id = 0;
67 std::vector<MeshLib::MeshSubset> all_mesh_subsets(
68 getProcessVariables(process_id)[0]
69 .get()
70 .getNumberOfGlobalComponents(),
71 *_mesh_subset_all_nodes);
72
73 // For scalar variables, in this case pressure.
74 all_mesh_subsets.push_back(*_mesh_subset_base_nodes);
75
76 std::vector<int> const vec_n_components{GlobalDim, 1};
77
78 _local_to_global_index_map =
79 std::make_unique<NumLib::LocalToGlobalIndexMap>(
80 std::move(all_mesh_subsets), vec_n_components,
82 assert(_local_to_global_index_map);
83 }
84}
85
86template <int GlobalDim>
88 NumLib::LocalToGlobalIndexMap const& dof_table,
89 MeshLib::Mesh const& mesh,
90 unsigned const integration_order)
91{
92 ProcessLib::createLocalAssemblersStokes<GlobalDim, LocalAssemblerData>(
93 mesh.getElements(), dof_table, _local_assemblers,
94 NumLib::IntegrationOrder{integration_order}, mesh.isAxiallySymmetric(),
95 _process_data);
96
97 _process_data.pressure_interpolated =
98 MeshLib::getOrCreateMeshProperty<double>(
99 const_cast<MeshLib::Mesh&>(mesh), "pressure_interpolated",
101}
102
103template <int GlobalDim>
105 std::map<int, std::shared_ptr<MaterialPropertyLib::Medium>> const& media)
106{
107 assert(_use_monolithic_scheme);
108 {
109 int const process_id = 0;
110 initializeProcessBoundaryConditionsAndSourceTerms(
111 *_local_to_global_index_map, process_id, media);
112 }
113}
114
115template <int GlobalDim>
118 const int /*process_id*/) const
119{
120 assert(_use_monolithic_scheme);
121 {
122 auto const& l = *_local_to_global_index_map;
123
124 return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
125 &l.getGhostIndices(), &this->_sparsity_pattern};
126 }
127}
128
129template <int GlobalDim>
131 const double t, double const dt, std::vector<GlobalVector*> const& x,
132 std::vector<GlobalVector*> const& x_prev, int const process_id,
134{
135 DBUG("Assemble StokesFlowProcess.");
136
137 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
138 assert(_use_monolithic_scheme);
139 {
140 dof_tables.push_back(_local_to_global_index_map.get());
141 }
142
143 ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
144 // Call global assembler for each local assembly item.
146 _global_assembler, &VectorMatrixAssembler::assemble, _local_assemblers,
147 pv.getActiveElementIDs(), dof_tables, t, dt, x, x_prev, process_id, M,
148 K, b);
149}
150
151template <int GlobalDim>
153 const double /*t*/, double const /*dt*/,
154 std::vector<GlobalVector*> const& /*x*/,
155 std::vector<GlobalVector*> const& /*x_prev*/, int const /*process_id*/,
156 GlobalMatrix& /*M*/, GlobalMatrix& /*K*/, GlobalVector& /*b*/,
157 GlobalMatrix& /*Jac*/)
158{
159 OGS_FATAL(
160 "Assembly of Jacobian matrix has not yet been implemented for "
161 "StokesFlowProcess.");
162}
163
164template <int GlobalDim>
166 double const t,
167 double const dt,
168 std::vector<GlobalVector*> const& x,
169 GlobalVector const& x_prev,
170 int const process_id)
171{
172 if (process_id != 0)
173 {
174 return;
175 }
176
177 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
178 dof_tables.reserve(x.size());
179 assert(_use_monolithic_scheme);
180 {
181 dof_tables.push_back(_local_to_global_index_map.get());
182 }
183
184 ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
187 _local_assemblers, pv.getActiveElementIDs(), dof_tables, t, dt, x,
188 x_prev, process_id);
189}
190
191template <int GlobalDim>
193 std::vector<GlobalVector*> const& x,
194 std::vector<GlobalVector*> const& x_prev,
195 const double t,
196 const double dt,
197 int const process_id)
198{
199 if (process_id != 0)
200 {
201 return;
202 }
203
204 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
205 dof_tables.reserve(x.size());
206 assert(_use_monolithic_scheme);
207 {
208 dof_tables.push_back(_local_to_global_index_map.get());
209 }
210
211 ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
214 pv.getActiveElementIDs(), dof_tables, x, x_prev, t, dt, process_id);
215}
216
217template <int GlobalDim>
219 const int /*process_id*/) const
220{
221 assert(_use_monolithic_scheme);
222 {
223 return *_local_to_global_index_map;
224 }
225}
226
227template class StokesFlowProcess<2>;
228} // namespace StokesFlow
229} // namespace ProcessLib
#define OGS_FATAL(...)
Definition Error.h:26
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
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)
std::vector< std::size_t > const & getActiveElementIDs() const
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, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, GlobalMatrix &Jac) override
void initializeConcreteProcess(NumLib::LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh, unsigned const integration_order) override
Process specific initialization called by initialize().
void computeSecondaryVariableConcrete(double const, double const, std::vector< GlobalVector * > const &x, GlobalVector const &, int const) 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
void initializeBoundaryConditions(std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media) override
NumLib::LocalToGlobalIndexMap const & getDOFTable(const int) const override
StokesFlowProcess(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, StokesFlowProcessData &&process_data, SecondaryVariableCollection &&secondary_variables, bool const use_monolithic_scheme)
MathLib::MatrixSpecifications getMatrixSpecifications(const int) const 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)
std::vector< Node * > getBaseNodes(std::vector< Element * > const &elements)
Definition Utils.h:26
@ BY_LOCATION
Ordering data by spatial location.
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)