OGS
Process.h
Go to the documentation of this file.
1
11#pragma once
12
13#include <string>
14#include <tuple>
15
27#include "ProcessVariable.h"
30#include "processlib_export.h"
31
32namespace MeshLib
33{
34class Mesh;
35}
36
37namespace ProcessLib
38{
39
41 : public NumLib::ODESystem< // TODO: later on use a simpler ODE system
42 NumLib::ODESystemTag::FirstOrderImplicitQuasilinear,
43 NumLib::NonlinearSolverTag::Newton>,
45{
46public:
48
49 Process(std::string name_,
50 MeshLib::Mesh& mesh,
51 std::unique_ptr<AbstractJacobianAssembler>&& jacobian_assembler,
52 std::vector<std::unique_ptr<ParameterLib::ParameterBase>> const&
53 parameters,
54 unsigned const integration_order,
55 std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>&&
56 process_variables,
57 SecondaryVariableCollection&& secondary_variables,
58 const bool use_monolithic_scheme = true);
59
61 void preTimestep(std::vector<GlobalVector*> const& x, const double t,
62 const double delta_t, const int process_id);
63
65 void postTimestep(std::vector<GlobalVector*> const& x,
66 std::vector<GlobalVector*> const& x_prev, const double t,
67 const double delta_t, int const process_id);
68
71 void postNonLinearSolver(std::vector<GlobalVector*> const& x,
72 std::vector<GlobalVector*> const& x_prev,
73 const double t, double const dt,
74 int const process_id);
75
76 void preIteration(const unsigned iter, GlobalVector const& x) final;
77
79 void computeSecondaryVariable(double const t,
80 double const dt,
81 std::vector<GlobalVector*> const& x,
82 GlobalVector const& x_prev,
83 int const process_id);
84
86
87 void initialize(
88 std::map<int, std::shared_ptr<MaterialPropertyLib::Medium>> const&
89 media);
90
92 std::vector<GlobalVector*>& process_solutions,
93 std::vector<GlobalVector*> const& process_solutions_prev,
94 double const t,
95 int const process_id);
96
98 const int process_id) const override;
99
100 void updateDeactivatedSubdomains(double const time, const int process_id);
101
102 virtual bool isMonolithicSchemeUsed() const
103 {
105 }
106
108 const double /*t*/,
109 std::vector<GlobalVector*> const& /*integration_point_values_vectors*/,
110 std::vector<GlobalVector*>& /*nodal_values_vectors*/)
111 {
112 }
113
114 void preAssemble(const double t, double const dt,
115 GlobalVector const& x) final;
116
117 void assemble(const double t, double const dt,
118 std::vector<GlobalVector*> const& x,
119 std::vector<GlobalVector*> const& x_prev,
120 int const process_id, GlobalMatrix& M, GlobalMatrix& K,
121 GlobalVector& b) final;
122
123 void assembleWithJacobian(const double t, double const dt,
124 std::vector<GlobalVector*> const& x,
125 std::vector<GlobalVector*> const& x_prev,
126 int const process_id, GlobalVector& b,
127 GlobalMatrix& Jac) final;
128
129 /* Computes data necessary for output.
130 *
131 * \pre Must be called before postTimestep() since the latter will overwrite
132 * previous states with current ones possibly affecting the data computed by
133 * this method.
134 */
135 void preOutput(const double t, double const dt,
136 std::vector<GlobalVector*> const& x,
137 std::vector<GlobalVector*> const& x_prev,
138 int const process_id);
139
140 std::vector<NumLib::IndexValueVector<GlobalIndexType>> const*
141 getKnownSolutions(double const t, GlobalVector const& x,
142 int const process_id) const final
143 {
144 return _boundary_conditions[process_id].getKnownSolutions(t, x);
145 }
146
148 const int /*process_id*/) const
149 {
151 }
152
153 MeshLib::Mesh& getMesh() const { return _mesh; }
154 std::vector<std::reference_wrapper<ProcessVariable>> const&
155 getProcessVariables(const int process_id) const
156 {
157 return _process_variables[process_id];
158 }
159
160 std::vector<std::size_t> const& getActiveElementIDs() const
161 {
163 }
164
169
170 std::vector<std::unique_ptr<MeshLib::IntegrationPointWriter>> const&
175
176 // Used as a call back for CalculateSurfaceFlux process.
177
178 virtual Eigen::Vector3d getFlux(
179 std::size_t /*element_id*/,
180 MathLib::Point3d const& /*p*/,
181 double const /*t*/,
182 std::vector<GlobalVector*> const& /*x*/) const
183 {
184 return Eigen::Vector3d{};
185 }
186
188 std::vector<GlobalVector*>& /*x*/,
189 std::vector<GlobalVector*> const& /*x_prev*/, double const /*t*/,
190 double const /*dt*/, NumLib::EquationSystem& /*ode_sys*/,
191 int const /*process_id*/)
192 {
193 }
194
195protected:
196 std::vector<NumLib::LocalToGlobalIndexMap const*> getDOFTables(
197 int const number_of_processes) const;
198
203
208
215 const NumLib::LocalToGlobalIndexMap& dof_table, const int process_id,
216 std::map<int, std::shared_ptr<MaterialPropertyLib::Medium>> const&
217 media);
218
219private:
222 NumLib::LocalToGlobalIndexMap const& dof_table,
223 MeshLib::Mesh const& mesh,
224 unsigned const integration_order) = 0;
225
228 virtual void initializeBoundaryConditions(
229 std::map<int, std::shared_ptr<MaterialPropertyLib::Medium>> const&
230 media);
231
233 std::vector<GlobalVector*>& /*x*/,
234 double const /*t*/,
235 int const /*process_id*/)
236 {
237 }
238
239 virtual void preAssembleConcreteProcess(const double /*t*/,
240 double const /*dt*/,
241 GlobalVector const& /*x*/)
242 {
243 }
244
246 const double t, double const dt, std::vector<GlobalVector*> const& x,
247 std::vector<GlobalVector*> const& x_prev, int const process_id,
248 GlobalMatrix& M, GlobalMatrix& K, GlobalVector& b) = 0;
249
251 const double t, double const dt, std::vector<GlobalVector*> const& x,
252 std::vector<GlobalVector*> const& x_prev, int const process_id,
253 GlobalVector& b, GlobalMatrix& Jac) = 0;
254
256 std::vector<GlobalVector*> const& /*x*/,
257 const double /*t*/,
258 const double /*dt*/,
259 const int /*process_id*/)
260 {
261 }
262
264 std::vector<GlobalVector*> const& /*x*/,
265 std::vector<GlobalVector*> const& /*x_prev*/,
266 const double /*t*/,
267 const double /*dt*/,
268 int const /*process_id*/)
269 {
270 }
271
273 std::vector<GlobalVector*> const& /*x*/,
274 std::vector<GlobalVector*> const& /*x_prev*/, const double /*t*/,
275 double const /*dt*/, int const /*process_id*/)
276 {
277 }
278
279 virtual void preIterationConcreteProcess(const unsigned /*iter*/,
280 GlobalVector const& /*x*/)
281 {
282 }
283
285 double const /*t*/,
286 double const /*dt*/,
287 std::vector<GlobalVector*> const& /*x*/,
288 GlobalVector const& /*x_prev*/,
289 int const /*process_id*/)
290 {
291 }
292
298
300 const double /*t*/, double const /*dt*/,
301 std::vector<GlobalVector*> const& /*x*/,
302 std::vector<GlobalVector*> const& /*x_prev*/, int const /*process_id*/)
303 {
304 }
305
306protected:
313 virtual void constructDofTable();
314
320
327 const int specified_process_id);
328
337 virtual std::tuple<NumLib::LocalToGlobalIndexMap*, bool>
339
342 std::vector<GlobalIndexType>
344
345private:
347
351
352public:
353 std::string const name;
354
355protected:
357 std::unique_ptr<MeshLib::MeshSubset const> _mesh_subset_all_nodes;
358
359 std::unique_ptr<NumLib::LocalToGlobalIndexMap> _local_to_global_index_map;
360
362
364
365 // TODO (CL) switch to the parallel vector matrix assembler here, once
366 // feature complete, or use the assembly mixin and remove these two members.
367 std::unique_ptr<ProcessLib::AbstractJacobianAssembler> _jacobian_assembler;
369
371
375 unsigned const _integration_order;
376
380 std::vector<std::unique_ptr<MeshLib::IntegrationPointWriter>>
382
384
385protected:
390 std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>
392
396 std::vector<BoundaryConditionCollection> _boundary_conditions;
397
398private:
402 std::vector<SourceTermCollection> _source_term_collections;
403
405
407 std::vector<std::size_t> _ids_of_active_elements;
408};
409
410} // namespace ProcessLib
MathLib::SparsityPattern< GlobalIndexType > GlobalSparsityPattern
Global vector based on Eigen vector.
Definition EigenVector.h:25
NumLib::LocalToGlobalIndexMap const & getDOFTable() const
NumLib::Extrapolator & getExtrapolator() const
std::string const name
Definition Process.h:353
std::vector< NumLib::IndexValueVector< GlobalIndexType > > const * getKnownSolutions(double const t, GlobalVector const &x, int const process_id) const final
Definition Process.h:141
std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > const & getIntegrationPointWriters() const
Definition Process.h:171
MeshLib::Mesh & getMesh() const
Definition Process.h:153
virtual NumLib::IterationResult postIterationConcreteProcess(GlobalVector const &)
Definition Process.h:293
virtual Eigen::Vector3d getFlux(std::size_t, MathLib::Point3d const &, double const, std::vector< GlobalVector * > const &) const
Definition Process.h:178
virtual void preTimestepConcreteProcess(std::vector< GlobalVector * > const &, const double, const double, const int)
Definition Process.h:255
std::vector< BoundaryConditionCollection > _boundary_conditions
Definition Process.h:396
std::vector< NumLib::LocalToGlobalIndexMap const * > getDOFTables(int const number_of_processes) const
Definition Process.cpp:385
std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > _integration_point_writer
Definition Process.h:381
SecondaryVariableCollection const & getSecondaryVariables() const
Definition Process.h:165
virtual void preIterationConcreteProcess(const unsigned, GlobalVector const &)
Definition Process.h:279
void constructDofTableOfSpecifiedProcessStaggeredScheme(const int specified_process_id)
Definition Process.cpp:354
static PROCESSLIB_EXPORT const std::string constant_one_parameter_name
Definition Process.h:47
void preAssemble(const double t, double const dt, GlobalVector const &x) final
Definition Process.cpp:260
virtual void initializeBoundaryConditions(std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
Definition Process.cpp:105
void initialize(std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
Definition Process.cpp:119
void initializeProcessBoundaryConditionsAndSourceTerms(const NumLib::LocalToGlobalIndexMap &dof_table, const int process_id, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
Definition Process.cpp:90
virtual void solveReactionEquation(std::vector< GlobalVector * > &, std::vector< GlobalVector * > const &, double const, double const, NumLib::EquationSystem &, int const)
Definition Process.h:187
virtual void initializeConcreteProcess(NumLib::LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh, unsigned const integration_order)=0
Process specific initialization called by initialize().
std::vector< std::size_t > _ids_of_active_elements
Union of active element ids per process variable.
Definition Process.h:407
virtual bool isMonolithicSchemeUsed() const
Definition Process.h:102
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:266
CellAverageData cell_average_data_
Definition Process.h:363
virtual void postTimestepConcreteProcess(std::vector< GlobalVector * > const &, std::vector< GlobalVector * > const &, const double, const double, int const)
Definition Process.h:263
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:44
virtual void computeSecondaryVariableConcrete(double const, double const, std::vector< GlobalVector * > const &, GlobalVector const &, int const)
Definition Process.h:284
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_all_nodes
Definition Process.h:357
void preOutput(const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id)
Definition Process.cpp:500
MeshLib::Mesh & _mesh
Definition Process.h:356
virtual 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)=0
std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > _process_variables
Definition Process.h:391
std::vector< std::size_t > const & getActiveElementIDs() const
Definition Process.h:160
NumLib::LocalToGlobalIndexMap const & getSingleComponentDOFTable() const
Definition Process.h:204
virtual NumLib::LocalToGlobalIndexMap const & getDOFTable(const int) const
Definition Process.h:147
virtual void setInitialConditionsConcreteProcess(std::vector< GlobalVector * > &, double const, int const)
Definition Process.h:232
NumLib::IterationResult postIteration(GlobalVector const &x) final
Definition Process.cpp:494
void postTimestep(std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, const double t, const double delta_t, int const process_id)
Postprocessing after a complete timestep.
Definition Process.cpp:453
SecondaryVariableCollection _secondary_variables
Definition Process.h:361
virtual 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)=0
void preTimestep(std::vector< GlobalVector * > const &x, const double t, const double delta_t, const int process_id)
Preprocessing before starting assembly for new timestep.
Definition Process.cpp:443
void computeSecondaryVariable(double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_prev, int const process_id)
compute secondary variables for the coupled equations or for output.
Definition Process.cpp:475
void initializeExtrapolator()
Definition Process.cpp:420
virtual void postNonLinearSolverConcreteProcess(std::vector< GlobalVector * > const &, std::vector< GlobalVector * > const &, const double, double const, int const)
Definition Process.h:272
std::vector< std::reference_wrapper< ProcessVariable > > const & getProcessVariables(const int process_id) const
Definition Process.h:155
virtual void constructDofTable()
Definition Process.cpp:307
void constructMonolithicProcessDofTable()
Definition Process.cpp:321
VectorMatrixAssembler _global_assembler
Definition Process.h:368
virtual void preAssembleConcreteProcess(const double, double const, GlobalVector const &)
Definition Process.h:239
void postNonLinearSolver(std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, const double t, double const dt, int const process_id)
Definition Process.cpp:465
std::vector< GlobalIndexType > getIndicesOfResiduumWithoutInitialCompensation() const override
Definition Process.cpp:512
unsigned const _integration_order
Definition Process.h:375
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
Definition Process.h:359
std::vector< SourceTermCollection > _source_term_collections
Definition Process.h:402
virtual void extrapolateIntegrationPointValuesToNodes(const double, std::vector< GlobalVector * > const &, std::vector< GlobalVector * > &)
Definition Process.h:107
std::unique_ptr< ProcessLib::AbstractJacobianAssembler > _jacobian_assembler
Definition Process.h:367
void updateDeactivatedSubdomains(double const time, const int process_id)
Definition Process.cpp:219
void computeSparsityPattern()
Definition Process.cpp:437
MathLib::MatrixSpecifications getMatrixSpecifications(const int process_id) const override
Definition Process.cpp:211
void setInitialConditions(std::vector< GlobalVector * > &process_solutions, std::vector< GlobalVector * > const &process_solutions_prev, double const t, int const process_id)
Definition Process.cpp:140
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:287
virtual void preOutputConcreteProcess(const double, double const, std::vector< GlobalVector * > const &, std::vector< GlobalVector * > const &, int const)
Definition Process.h:299
ExtrapolatorData _extrapolator_data
Definition Process.h:404
NumLib::Extrapolator & getExtrapolator() const
Definition Process.h:199
GlobalSparsityPattern _sparsity_pattern
Definition Process.h:383
void preIteration(const unsigned iter, GlobalVector const &x) final
Definition Process.cpp:488
const bool _use_monolithic_scheme
Definition Process.h:370
virtual std::tuple< NumLib::LocalToGlobalIndexMap *, bool > getDOFTableForExtrapolatorData() const
Definition Process.cpp:396
Handles configuration of several secondary variables from the project file.
IterationResult
Status flags telling the NonlinearSolver if an iteration succeeded.
#define PROCESSLIB_EXPORT