OGS
Process.cpp
Go to the documentation of this file.
1 
11 #include "Process.h"
12 
17 #include "ParameterLib/Parameter.h"
18 #include "ProcessVariable.h"
19 
20 namespace ProcessLib
21 {
23  std::string name_,
24  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,
30  SecondaryVariableCollection&& secondary_variables,
31  const bool use_monolithic_scheme)
32  : name(std::move(name_)),
33  _mesh(mesh),
34  _secondary_variables(std::move(secondary_variables)),
35  _global_assembler(std::move(jacobian_assembler)),
36  _use_monolithic_scheme(use_monolithic_scheme),
37  _coupled_solutions(nullptr),
38  _integration_order(integration_order),
39  _process_variables(std::move(process_variables)),
40  _boundary_conditions(
41  [&](const std::size_t number_of_process_variables)
42  -> std::vector<BoundaryConditionCollection>
43  {
44  std::vector<BoundaryConditionCollection> pcs_BCs;
45  pcs_BCs.reserve(number_of_process_variables);
46  for (std::size_t i = 0; i < number_of_process_variables; i++)
47  {
48  pcs_BCs.emplace_back(BoundaryConditionCollection(parameters));
49  }
50  return pcs_BCs;
51  }(_process_variables.size())),
52  _source_term_collections(
53  [&](const std::size_t number_of_processes)
54  -> std::vector<SourceTermCollection>
55  {
56  std::vector<SourceTermCollection> pcs_sts;
57  pcs_sts.reserve(number_of_processes);
58  for (std::size_t i = 0; i < number_of_processes; i++)
59  {
60  pcs_sts.emplace_back(SourceTermCollection(parameters));
61  }
62  return pcs_sts;
63  }(_process_variables.size()))
64 {
65 }
66 
68  const NumLib::LocalToGlobalIndexMap& dof_table, const int process_id)
69 {
70  auto const& per_process_variables = _process_variables[process_id];
71  auto& per_process_BCs = _boundary_conditions[process_id];
72 
73  per_process_BCs.addBCsForProcessVariables(per_process_variables, dof_table,
74  _integration_order, *this);
75 
76  auto& per_process_sts = _source_term_collections[process_id];
77  per_process_sts.addSourceTermsForProcessVariables(
78  per_process_variables, dof_table, _integration_order);
79 }
80 
82 {
83  // The number of processes is identical to the size of _process_variables,
84  // the vector contains variables for different processes. See the
85  // documentation of _process_variables.
86  const std::size_t number_of_processes = _process_variables.size();
87  for (std::size_t pcs_id = 0; pcs_id < number_of_processes; pcs_id++)
88  {
91  }
92 }
93 
95 {
96  DBUG("Initialize process.");
97 
98  DBUG("Construct dof mappings.");
100 
101  DBUG("Compute sparsity pattern");
103 
104  DBUG("Initialize the extrapolator");
106 
109 
110  DBUG("Initialize boundary conditions.");
112 }
113 
115  std::vector<GlobalVector*>& process_solutions,
116  std::vector<GlobalVector*> const& process_solutions_prev,
117  double const t,
118  int const process_id)
119 {
120  auto& x = *process_solutions[process_id];
121  auto& x_prev = *process_solutions_prev[process_id];
122 
123  // getDOFTableOfProcess can be overloaded by the specific process.
124  auto const& dof_table_of_process = getDOFTable(process_id);
125 
126  auto const& per_process_variables = _process_variables[process_id];
127  for (std::size_t variable_id = 0;
128  variable_id < per_process_variables.size();
129  variable_id++)
130  {
133 
134  auto const& pv = per_process_variables[variable_id];
135  DBUG("Set the initial condition of variable {:s} of process {:d}.",
136  pv.get().getName().data(), process_id);
137 
138  auto const& ic = pv.get().getInitialCondition();
139 
140  auto const num_comp = pv.get().getNumberOfGlobalComponents();
141 
142  for (int component_id = 0; component_id < num_comp; ++component_id)
143  {
144  auto const& mesh_subset =
145  dof_table_of_process.getMeshSubset(variable_id, component_id);
146  auto const mesh_id = mesh_subset.getMeshID();
147  for (auto const* node : mesh_subset.getNodes())
148  {
150  node->getID());
151 
152  pos.setNodeID(node->getID());
153  pos.setCoordinates(*node);
154  auto const& ic_value = ic(t, pos);
155 
156  auto global_index =
157  std::abs(dof_table_of_process.getGlobalIndex(l, variable_id,
158  component_id));
159 #ifdef USE_PETSC
160  // The global indices of the ghost entries of the global
161  // matrix or the global vectors need to be set as negative
162  // values for equation assembly, however the global indices
163  // start from zero. Therefore, any ghost entry with zero
164  // index is assigned an negative value of the vector size
165  // or the matrix dimension. To assign the initial value for
166  // the ghost entries, the negative indices of the ghost
167  // entries are restored to zero.
168  if (global_index == x.size())
169  global_index = 0;
170 #endif
171  x.set(global_index, ic_value[component_id]);
172  }
173  }
174  }
175 
177  MathLib::LinAlg::copy(x, x_prev); // pushState
178 
179  setInitialConditionsConcreteProcess(process_solutions, t, process_id);
180 }
181 
183  const int /*process_id*/) const
184 {
185  auto const& l = *_local_to_global_index_map;
186  return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
187  &l.getGhostIndices(), &_sparsity_pattern};
188 }
189 
191  const int process_id)
192 {
193  auto const& variables_per_process = getProcessVariables(process_id);
194  for (auto const& variable : variables_per_process)
195  {
196  variable.get().updateDeactivatedSubdomains(time);
197  }
198 }
199 
200 void Process::preAssemble(const double t, double const dt,
201  GlobalVector const& x)
202 {
203  preAssembleConcreteProcess(t, dt, x);
204 }
205 
206 void Process::assemble(const double t, double const dt,
207  std::vector<GlobalVector*> const& x,
208  std::vector<GlobalVector*> const& xdot,
209  int const process_id, GlobalMatrix& M, GlobalMatrix& K,
210  GlobalVector& b)
211 {
212  assert(x.size() == xdot.size());
213  for (std::size_t i = 0; i < x.size(); i++)
214  {
217  }
218 
219  assembleConcreteProcess(t, dt, x, xdot, process_id, M, K, b);
220 
221  // the last argument is for the jacobian, nullptr is for a unused jacobian
222  _boundary_conditions[process_id].applyNaturalBC(t, x, process_id, K, b,
223  nullptr);
224 
225  // the last argument is for the jacobian, nullptr is for a unused jacobian
226  _source_term_collections[process_id].integrate(t, *x[process_id], b,
227  nullptr);
228 }
229 
230 void Process::assembleWithJacobian(const double t, double const dt,
231  std::vector<GlobalVector*> const& x,
232  std::vector<GlobalVector*> const& xdot,
233  const double dxdot_dx, const double dx_dx,
234  int const process_id, GlobalMatrix& M,
235  GlobalMatrix& K, GlobalVector& b,
236  GlobalMatrix& Jac)
237 {
238  assert(x.size() == xdot.size());
239  for (std::size_t i = 0; i < x.size(); i++)
240  {
243  }
244 
245  assembleWithJacobianConcreteProcess(t, dt, x, xdot, dxdot_dx, dx_dx,
246  process_id, M, K, b, Jac);
247 
248  // TODO: apply BCs to Jacobian.
249  _boundary_conditions[process_id].applyNaturalBC(t, x, process_id, K, b,
250  &Jac);
251 
252  // the last argument is for the jacobian, nullptr is for a unused jacobian
253  _source_term_collections[process_id].integrate(t, *x[process_id], b, &Jac);
254 }
255 
257 {
259  {
261 
262  return;
263  }
264 
265  // For staggered scheme:
266  const int specified_process_id = 0;
268 }
269 
271 {
272  // Create single component dof in every of the mesh nodes.
274  std::make_unique<MeshLib::MeshSubset>(_mesh, _mesh.getNodes());
275 
276  // Vector of mesh subsets.
277  std::vector<MeshLib::MeshSubset> all_mesh_subsets;
278 
279  // Collect the mesh subsets in a vector for the components of each
280  // variables.
281  for (ProcessVariable const& pv : _process_variables[0])
282  {
283  std::generate_n(std::back_inserter(all_mesh_subsets),
284  pv.getNumberOfGlobalComponents(),
285  [&]() { return *_mesh_subset_all_nodes; });
286  }
287 
288  // Create a vector of the number of variable components
289  std::vector<int> vec_var_n_components;
290  transform(cbegin(_process_variables[0]), cend(_process_variables[0]),
291  back_inserter(vec_var_n_components),
292  [](ProcessVariable const& pv)
293  { return pv.getNumberOfGlobalComponents(); });
294 
296  std::make_unique<NumLib::LocalToGlobalIndexMap>(
297  std::move(all_mesh_subsets), vec_var_n_components,
299 
301 }
302 
304  const int specified_process_id)
305 {
306  // Create single component dof in every of the mesh nodes.
308  std::make_unique<MeshLib::MeshSubset>(_mesh, _mesh.getNodes());
309 
310  // Vector of mesh subsets.
311  std::vector<MeshLib::MeshSubset> all_mesh_subsets;
312 
313  // Vector of the number of variable components
314  std::vector<int> vec_var_n_components;
315  // Collect the mesh subsets in a vector for each variables' components.
316  std::generate_n(std::back_inserter(all_mesh_subsets),
317  _process_variables[specified_process_id][0]
318  .get()
319  .getNumberOfGlobalComponents(),
320  [&]() { return *_mesh_subset_all_nodes; });
321 
322  // Create a vector of the number of variable components.
323  vec_var_n_components.push_back(_process_variables[specified_process_id][0]
324  .get()
325  .getNumberOfGlobalComponents());
327  std::make_unique<NumLib::LocalToGlobalIndexMap>(
328  std::move(all_mesh_subsets), vec_var_n_components,
330 
332 }
333 
334 std::tuple<NumLib::LocalToGlobalIndexMap*, bool>
336 {
337  if (_local_to_global_index_map->getNumberOfGlobalComponents() == 1)
338  {
339  // For single-variable-single-component processes reuse the existing DOF
340  // table.
341  const bool manage_storage = false;
342  return std::make_tuple(_local_to_global_index_map.get(),
343  manage_storage);
344  }
345 
346  // Otherwise construct a new DOF table.
347  std::vector<MeshLib::MeshSubset> all_mesh_subsets_single_component;
348  all_mesh_subsets_single_component.emplace_back(*_mesh_subset_all_nodes);
349 
350  const bool manage_storage = true;
351 
352  return std::make_tuple(new NumLib::LocalToGlobalIndexMap(
353  std::move(all_mesh_subsets_single_component),
354  // by location order is needed for output
356  manage_storage);
357 }
358 
360 {
361  NumLib::LocalToGlobalIndexMap* dof_table_single_component;
362  bool manage_storage;
363 
364  std::tie(dof_table_single_component, manage_storage) =
366 
367  std::unique_ptr<NumLib::Extrapolator> extrapolator(
369  *dof_table_single_component));
370 
371  // TODO: Later on the DOF table can change during the simulation!
373  std::move(extrapolator), dof_table_single_component, manage_storage);
374 }
375 
377 {
380 }
381 
382 void Process::preTimestep(std::vector<GlobalVector*> const& x, const double t,
383  const double delta_t, const int process_id)
384 {
385  for (auto* const solution : x)
387  preTimestepConcreteProcess(x, t, delta_t, process_id);
388 
389  _boundary_conditions[process_id].preTimestep(t, x, process_id);
390 }
391 
392 void Process::postTimestep(std::vector<GlobalVector*> const& x, const double t,
393  const double delta_t, int const process_id)
394 {
395  for (auto* const solution : x)
397  postTimestepConcreteProcess(x, t, delta_t, process_id);
398 
399  _boundary_conditions[process_id].postTimestep(t, x, process_id);
400 }
401 
403  GlobalVector const& xdot, const double t,
404  double const dt, int const process_id)
405 {
408  postNonLinearSolverConcreteProcess(x, xdot, t, dt, process_id);
409 }
410 
412  double const dt,
413  std::vector<GlobalVector*> const& x,
414  GlobalVector const& x_dot,
415  int const process_id)
416 {
417  for (auto const* solution : x)
420 
421  computeSecondaryVariableConcrete(t, dt, x, x_dot, process_id);
422 }
423 
424 void Process::preIteration(const unsigned iter, const GlobalVector& x)
425 {
428 }
429 
431 {
434 }
435 
436 } // namespace ProcessLib
void DBUG(char const *fmt, Args const &... args)
Definition: Logging.h:27
Global vector based on Eigen vector.
Definition: EigenVector.h:26
std::vector< Node * > const & getNodes() const
Get the nodes-vector for the mesh.
Definition: Mesh.h:95
void setNodeID(std::size_t node_id)
void setCoordinates(MathLib::TemplatePoint< double, 3 > const &coordinates)
int getNumberOfGlobalComponents() const
Returns the number of components of the process variable.
virtual void initializeBoundaryConditions()
Definition: Process.cpp:81
virtual NumLib::IterationResult postIterationConcreteProcess(GlobalVector const &)
Definition: Process.h:275
virtual void preTimestepConcreteProcess(std::vector< GlobalVector * > const &, const double, const double, const int)
Definition: Process.h:239
std::vector< BoundaryConditionCollection > _boundary_conditions
Definition: Process.h:365
void constructDofTableOfSpecifiedProcessStaggeredScheme(const int specified_prosess_id)
Definition: Process.cpp:303
virtual void postTimestepConcreteProcess(std::vector< GlobalVector * > const &, const double, const double, int const)
Definition: Process.h:247
void postNonLinearSolver(GlobalVector const &x, GlobalVector const &xdot, const double t, double const dt, int const process_id)
Definition: Process.cpp:402
virtual void preIterationConcreteProcess(const unsigned, GlobalVector const &)
Definition: Process.h:261
virtual void assembleConcreteProcess(const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &xdot, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b)=0
void assembleWithJacobian(const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &xdot, const double dxdot_dx, const double dx_dx, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, GlobalMatrix &Jac) final
Definition: Process.cpp:230
void preAssemble(const double t, double const dt, GlobalVector const &x) final
Definition: Process.cpp:200
virtual void initializeConcreteProcess(NumLib::LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh, unsigned const integration_order)=0
Process specific initialization called by initialize().
void initializeProcessBoundaryConditionsAndSourceTerms(const NumLib::LocalToGlobalIndexMap &dof_table, const int process_id)
Definition: Process.cpp:67
void computeSecondaryVariable(double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_dot, int const process_id)
compute secondary variables for the coupled equations or for output.
Definition: Process.cpp:411
virtual void computeSecondaryVariableConcrete(double const, double const, std::vector< GlobalVector * > const &, GlobalVector const &, int const)
Definition: Process.h:266
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_all_nodes
Definition: Process.h:327
MeshLib::Mesh & _mesh
Definition: Process.h:326
std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > _process_variables
Definition: Process.h:360
virtual NumLib::LocalToGlobalIndexMap const & getDOFTable(const int) const
Definition: Process.h:137
virtual void postNonLinearSolverConcreteProcess(GlobalVector const &, GlobalVector const &, const double, double const, int const)
Definition: Process.h:255
virtual void setInitialConditionsConcreteProcess(std::vector< GlobalVector * > &, double const, int const)
Definition: Process.h:214
NumLib::IterationResult postIteration(GlobalVector const &x) final
Definition: Process.cpp:430
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:382
void initializeExtrapolator()
Definition: Process.cpp:359
std::vector< std::reference_wrapper< ProcessVariable > > const & getProcessVariables(const int process_id) const
Definition: Process.h:145
virtual void constructDofTable()
Definition: Process.cpp:256
void constructMonolithicProcessDofTable()
Definition: Process.cpp:270
virtual void preAssembleConcreteProcess(const double, double const, GlobalVector const &)
Definition: Process.h:221
virtual void assembleWithJacobianConcreteProcess(const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &xdot, const double dxdot_dx, const double dx_dx, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, GlobalMatrix &Jac)=0
unsigned const _integration_order
Definition: Process.h:344
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
Definition: Process.h:329
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:22
std::vector< SourceTermCollection > _source_term_collections
Definition: Process.h:371
void updateDeactivatedSubdomains(double const time, const int process_id)
Definition: Process.cpp:190
void computeSparsityPattern()
Definition: Process.cpp:376
MathLib::MatrixSpecifications getMatrixSpecifications(const int process_id) const override
Definition: Process.cpp:182
void setInitialConditions(std::vector< GlobalVector * > &process_solutions, std::vector< GlobalVector * > const &process_solutions_prev, double const t, int const process_id)
Definition: Process.cpp:114
void assemble(const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &xdot, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b) final
Definition: Process.cpp:206
void postTimestep(std::vector< GlobalVector * > const &x, const double t, const double delta_t, int const process_id)
Postprocessing after a complete timestep.
Definition: Process.cpp:392
ExtrapolatorData _extrapolator_data
Definition: Process.h:373
GlobalSparsityPattern _sparsity_pattern
Definition: Process.h:352
void preIteration(const unsigned iter, GlobalVector const &x) final
Definition: Process.cpp:424
const bool _use_monolithic_scheme
Definition: Process.h:335
virtual std::tuple< NumLib::LocalToGlobalIndexMap *, bool > getDOFTableForExtrapolatorData() const
Definition: Process.cpp:335
Handles configuration of several secondary variables from the project file.
IterationResult
Status flags telling the NonlinearSolver if an iteration succeeded.
void finalizeAssembly(PETScMatrix &A)
Definition: LinAlg.cpp:162
void copy(PETScVector const &x, PETScVector &y)
Definition: LinAlg.cpp:37
void setLocalAccessibleVector(PETScVector const &x)
Definition: LinAlg.cpp:27
@ 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.