OGS
Process.cpp
Go to the documentation of this file.
1
11#include "Process.h"
12
13#include <range/v3/algorithm/any_of.hpp>
14#include <range/v3/algorithm/set_algorithm.hpp>
15#include <range/v3/view/drop.hpp>
16#include <range/v3/view/transform.hpp>
17
22#include "ProcessVariable.h"
23
24namespace
25{
26void setLocalAccessibleVectors(std::vector<GlobalVector*> const& x,
27 std::vector<GlobalVector*> const& x_prev)
28{
29 for (auto* const solution : x)
30 {
32 }
33 for (auto* const solution : x_prev)
34 {
36 }
37}
38} // namespace
39
40namespace ProcessLib
41{
42const std::string Process::constant_one_parameter_name = "constant_one";
43
45 std::string name_,
46 MeshLib::Mesh& mesh,
47 std::unique_ptr<ProcessLib::AbstractJacobianAssembler>&& jacobian_assembler,
48 std::vector<std::unique_ptr<ParameterLib::ParameterBase>> const& parameters,
49 unsigned const integration_order,
50 std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>&&
51 process_variables,
52 SecondaryVariableCollection&& secondary_variables,
53 const bool use_monolithic_scheme)
54 : name(std::move(name_)),
55 _mesh(mesh),
56 _secondary_variables(std::move(secondary_variables)),
57 _jacobian_assembler(std::move(jacobian_assembler)),
58 _global_assembler(*_jacobian_assembler),
59 _use_monolithic_scheme(use_monolithic_scheme),
60 _integration_order(integration_order),
61 _process_variables(std::move(process_variables)),
62 _boundary_conditions(
63 [&](const std::size_t number_of_process_variables)
64 -> std::vector<BoundaryConditionCollection>
65 {
66 std::vector<BoundaryConditionCollection> pcs_BCs;
67 pcs_BCs.reserve(number_of_process_variables);
68 for (std::size_t i = 0; i < number_of_process_variables; i++)
69 {
70 pcs_BCs.emplace_back(BoundaryConditionCollection(parameters));
71 }
72 return pcs_BCs;
73 }(_process_variables.size())),
74 _source_term_collections(
75 [&](const std::size_t number_of_processes)
76 -> std::vector<SourceTermCollection>
77 {
78 std::vector<SourceTermCollection> pcs_sts;
79 pcs_sts.reserve(number_of_processes);
80 for (std::size_t i = 0; i < number_of_processes; i++)
81 {
82 pcs_sts.emplace_back(SourceTermCollection(parameters));
83 }
84 return pcs_sts;
85 }(_process_variables.size()))
86{
87}
88
90 const NumLib::LocalToGlobalIndexMap& dof_table, const int process_id,
91 std::map<int, std::shared_ptr<MaterialPropertyLib::Medium>> const& media)
92{
93 auto const& per_process_variables = _process_variables[process_id];
94 auto& per_process_BCs = _boundary_conditions[process_id];
95
96 per_process_BCs.addBCsForProcessVariables(per_process_variables, dof_table,
97 _integration_order, *this, media);
98
99 auto& per_process_sts = _source_term_collections[process_id];
100 per_process_sts.addSourceTermsForProcessVariables(
101 per_process_variables, dof_table, _integration_order);
102}
103
105 std::map<int, std::shared_ptr<MaterialPropertyLib::Medium>> const& media)
106{
107 // The number of processes is identical to the size of _process_variables,
108 // the vector contains variables for different processes. See the
109 // documentation of _process_variables.
110 const std::size_t number_of_processes = _process_variables.size();
111 for (std::size_t pcs_id = 0; pcs_id < number_of_processes; pcs_id++)
112 {
114 *_local_to_global_index_map, pcs_id, media);
115 }
116}
117
119 std::map<int, std::shared_ptr<MaterialPropertyLib::Medium>> const& media)
120{
121 DBUG("Initialize process.");
122
123 DBUG("Construct dof mappings.");
125
126 DBUG("Compute sparsity pattern");
128
129 DBUG("Initialize the extrapolator");
131
134
135 DBUG("Initialize boundary conditions.");
137}
138
140 std::vector<GlobalVector*>& process_solutions,
141 std::vector<GlobalVector*> const& process_solutions_prev,
142 double const t,
143 int const process_id)
144{
145 auto& x = *process_solutions[process_id];
146 auto& x_prev = *process_solutions_prev[process_id];
147
148 // getDOFTableOfProcess can be overloaded by the specific process.
149 auto const& dof_table_of_process = getDOFTable(process_id);
150
151 auto const& per_process_variables = _process_variables[process_id];
152 for (std::size_t variable_id = 0;
153 variable_id < per_process_variables.size();
154 variable_id++)
155 {
158
159 auto const& pv = per_process_variables[variable_id];
160 DBUG("Set the initial condition of variable {:s} of process {:d}.",
161 pv.get().getName().data(), process_id);
162
163 auto const& ic = pv.get().getInitialCondition();
164
165 auto const num_comp = pv.get().getNumberOfGlobalComponents();
166
167 for (int component_id = 0; component_id < num_comp; ++component_id)
168 {
169 auto const& mesh_subset =
170 dof_table_of_process.getMeshSubset(variable_id, component_id);
171 auto const mesh_id = mesh_subset.getMeshID();
172 for (auto const* node : mesh_subset.getNodes())
173 {
175 node->getID());
176
177 pos.setNodeID(node->getID());
178 pos.setCoordinates(*node);
179 auto const& ic_value = ic(t, pos);
180
181 auto global_index =
182 std::abs(dof_table_of_process.getGlobalIndex(l, variable_id,
183 component_id));
184#ifdef USE_PETSC
185 // The global indices of the ghost entries of the global
186 // matrix or the global vectors need to be set as negative
187 // values for equation assembly, however the global indices
188 // start from zero. Therefore, any ghost entry with zero
189 // index is assigned an negative value of the vector size
190 // or the matrix dimension. To assign the initial value for
191 // the ghost entries, the negative indices of the ghost
192 // entries are restored to zero.
193 if (global_index == x.size())
194 global_index = 0;
195#endif
196 x.set(global_index, ic_value[component_id]);
197 }
198 }
199 }
200
202 MathLib::LinAlg::copy(x, x_prev); // pushState
203
206
207 setInitialConditionsConcreteProcess(process_solutions, t, process_id);
208}
209
211 const int /*process_id*/) const
212{
213 auto const& l = *_local_to_global_index_map;
214 return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
215 &l.getGhostIndices(), &_sparsity_pattern};
216}
217
219 const int process_id)
220{
221 auto const& variables_per_process = getProcessVariables(process_id);
222 for (auto const& variable : variables_per_process)
223 {
224 variable.get().updateDeactivatedSubdomains(time);
225 }
226
228
229 auto active_elements_ids = ranges::views::transform(
230 [](auto const& variable)
231 { return variable.get().getActiveElementIDs(); });
232
233 // Early return if there's any process variable with all elements active.
234 if (ranges::any_of(variables_per_process | active_elements_ids,
235 [](auto const& vector) { return vector.empty(); }))
236 {
237 return;
238 }
239
240 // Some process variable has deactivated elements. Create union of active
241 // ids.
242
243 _ids_of_active_elements = // there is at least one process variable.
244 variables_per_process[0].get().getActiveElementIDs();
245
246 for (auto const& pv_active_element_ids :
247 variables_per_process | ranges::views::drop(1) | active_elements_ids)
248 {
249 std::vector<std::size_t> new_active_elements;
250 new_active_elements.reserve(_ids_of_active_elements.size() +
251 pv_active_element_ids.size());
252 ranges::set_union(_ids_of_active_elements, pv_active_element_ids,
253 std::back_inserter(new_active_elements));
254
255 _ids_of_active_elements = std::move(new_active_elements);
256 }
257}
258
259void Process::preAssemble(const double t, double const dt,
260 GlobalVector const& x)
261{
263}
264
265void Process::assemble(const double t, double const dt,
266 std::vector<GlobalVector*> const& x,
267 std::vector<GlobalVector*> const& x_prev,
268 int const process_id, GlobalMatrix& M, GlobalMatrix& K,
269 GlobalVector& b)
270{
271 assert(x.size() == x_prev.size());
272
273 setLocalAccessibleVectors(x, x_prev);
274
275 assembleConcreteProcess(t, dt, x, x_prev, process_id, M, K, b);
276
277 // the last argument is for the jacobian, nullptr is for a unused jacobian
278 _boundary_conditions[process_id].applyNaturalBC(t, x, process_id, K, b,
279 nullptr);
280
281 // the last argument is for the jacobian, nullptr is for a unused jacobian
282 _source_term_collections[process_id].integrate(t, *x[process_id], b,
283 nullptr);
284}
285
286void Process::assembleWithJacobian(const double t, double const dt,
287 std::vector<GlobalVector*> const& x,
288 std::vector<GlobalVector*> const& x_prev,
289 int const process_id, GlobalMatrix& M,
291 GlobalMatrix& Jac)
292{
293 assert(x.size() == x_prev.size());
294
295 setLocalAccessibleVectors(x, x_prev);
296
297 assembleWithJacobianConcreteProcess(t, dt, x, x_prev, process_id, M, K, b,
298 Jac);
299
300 // TODO: apply BCs to Jacobian.
301 _boundary_conditions[process_id].applyNaturalBC(t, x, process_id, K, b,
302 &Jac);
303
304 // the last argument is for the jacobian, nullptr is for a unused jacobian
305 _source_term_collections[process_id].integrate(t, *x[process_id], b, &Jac);
306}
307
309{
311 {
313
314 return;
315 }
316
317 // For staggered scheme:
318 const int specified_process_id = 0;
320}
321
323{
324 // Create single component dof in every of the mesh nodes.
326 std::make_unique<MeshLib::MeshSubset>(_mesh, _mesh.getNodes());
327
328 // Vector of mesh subsets.
329 std::vector<MeshLib::MeshSubset> all_mesh_subsets;
330
331 // Collect the mesh subsets in a vector for the components of each
332 // variables.
333 for (ProcessVariable const& pv : _process_variables[0])
334 {
335 std::generate_n(std::back_inserter(all_mesh_subsets),
336 pv.getNumberOfGlobalComponents(),
337 [&]() { return *_mesh_subset_all_nodes; });
338 }
339
340 // Create a vector of the number of variable components
341 std::vector<int> vec_var_n_components;
342 transform(cbegin(_process_variables[0]), cend(_process_variables[0]),
343 back_inserter(vec_var_n_components),
344 [](ProcessVariable const& pv)
345 { return pv.getNumberOfGlobalComponents(); });
346
348 std::make_unique<NumLib::LocalToGlobalIndexMap>(
349 std::move(all_mesh_subsets), vec_var_n_components,
351
353}
354
356 const int specified_process_id)
357{
358 // Create single component dof in every of the mesh nodes.
360 std::make_unique<MeshLib::MeshSubset>(_mesh, _mesh.getNodes());
361
362 // Vector of mesh subsets.
363 std::vector<MeshLib::MeshSubset> all_mesh_subsets;
364
365 // Vector of the number of variable components
366 std::vector<int> vec_var_n_components;
367 // Collect the mesh subsets in a vector for each variables' components.
368 std::generate_n(std::back_inserter(all_mesh_subsets),
369 _process_variables[specified_process_id][0]
370 .get()
371 .getNumberOfGlobalComponents(),
372 [&]() { return *_mesh_subset_all_nodes; });
373
374 // Create a vector of the number of variable components.
375 vec_var_n_components.push_back(_process_variables[specified_process_id][0]
376 .get()
377 .getNumberOfGlobalComponents());
379 std::make_unique<NumLib::LocalToGlobalIndexMap>(
380 std::move(all_mesh_subsets), vec_var_n_components,
382
384}
385
386std::vector<NumLib::LocalToGlobalIndexMap const*> Process::getDOFTables(
387 int const number_of_processes) const
388{
389 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
390 dof_tables.reserve(number_of_processes);
391 std::generate_n(std::back_inserter(dof_tables), number_of_processes,
392 [&]() { return &getDOFTable(dof_tables.size()); });
393 return dof_tables;
394}
395
396std::tuple<NumLib::LocalToGlobalIndexMap*, bool>
398{
399 if (_local_to_global_index_map->getNumberOfGlobalComponents() == 1)
400 {
401 // For single-variable-single-component processes reuse the existing DOF
402 // table.
403 const bool manage_storage = false;
404 return std::make_tuple(_local_to_global_index_map.get(),
405 manage_storage);
406 }
407
408 // Otherwise construct a new DOF table.
409 std::vector<MeshLib::MeshSubset> all_mesh_subsets_single_component;
410 all_mesh_subsets_single_component.emplace_back(*_mesh_subset_all_nodes);
411
412 const bool manage_storage = true;
413
414 return std::make_tuple(new NumLib::LocalToGlobalIndexMap(
415 std::move(all_mesh_subsets_single_component),
416 // by location order is needed for output
418 manage_storage);
419}
420
422{
423 NumLib::LocalToGlobalIndexMap* dof_table_single_component;
424 bool manage_storage;
425
426 std::tie(dof_table_single_component, manage_storage) =
428
429 std::unique_ptr<NumLib::Extrapolator> extrapolator(
431 *dof_table_single_component));
432
433 // TODO: Later on the DOF table can change during the simulation!
435 std::move(extrapolator), dof_table_single_component, manage_storage);
436}
437
443
444void Process::preTimestep(std::vector<GlobalVector*> const& x, const double t,
445 const double delta_t, const int process_id)
446{
447 for (auto* const solution : x)
449 preTimestepConcreteProcess(x, t, delta_t, process_id);
450
451 _boundary_conditions[process_id].preTimestep(t, x, process_id);
452}
453
454void Process::postTimestep(std::vector<GlobalVector*> const& x,
455 std::vector<GlobalVector*> const& x_prev,
456 const double t, const double delta_t,
457 int const process_id)
458{
459 setLocalAccessibleVectors(x, x_prev);
460
461 postTimestepConcreteProcess(x, x_prev, t, delta_t, process_id);
462
463 _boundary_conditions[process_id].postTimestep(t, x, process_id);
464}
465
466void Process::postNonLinearSolver(std::vector<GlobalVector*> const& x,
467 std::vector<GlobalVector*> const& x_prev,
468 const double t, double const dt,
469 int const process_id)
470{
471 setLocalAccessibleVectors(x, x_prev);
472
473 postNonLinearSolverConcreteProcess(x, x_prev, t, dt, process_id);
474}
475
477 double const dt,
478 std::vector<GlobalVector*> const& x,
479 GlobalVector const& x_prev,
480 int const process_id)
481{
482 for (auto const* solution : x)
485
486 computeSecondaryVariableConcrete(t, dt, x, x_prev, process_id);
487}
488
489void Process::preIteration(const unsigned iter, const GlobalVector& x)
490{
493}
494
500
501void Process::preOutput(const double t, double const dt,
502 std::vector<GlobalVector*> const& x,
503 std::vector<GlobalVector*> const& x_prev,
504 int const process_id)
505{
506 for (auto const* solution : x)
508
509 preOutputConcreteProcess(t, dt, x, x_prev, process_id);
510}
511
512std::vector<GlobalIndexType>
514{
515 std::vector<GlobalIndexType> indices;
516
517 for (std::size_t process_id = 0; process_id < _process_variables.size();
518 process_id++)
519 {
520 auto const& dof_table_of_process = getDOFTable(process_id);
521
522 auto const& per_process_variables = _process_variables[process_id];
523 for (std::size_t variable_id = 0;
524 variable_id < per_process_variables.size();
525 variable_id++)
526 {
527 auto const& pv = per_process_variables[variable_id];
528 DBUG("Set the initial condition of variable {:s} of process {:d}.",
529 pv.get().getName().data(), process_id);
530
531 if ((pv.get().compensateNonEquilibriumInitialResiduum()))
532 {
533 continue;
534 }
535
536 auto const num_comp = pv.get().getNumberOfGlobalComponents();
537
538 for (int component_id = 0; component_id < num_comp; ++component_id)
539 {
540 auto const& mesh_subset = dof_table_of_process.getMeshSubset(
541 variable_id, component_id);
542 auto const mesh_id = mesh_subset.getMeshID();
543 for (auto const* node : mesh_subset.getNodes())
544 {
545 MeshLib::Location const l(
546 mesh_id, MeshLib::MeshItemType::Node, node->getID());
547
548 auto global_index =
549 std::abs(dof_table_of_process.getGlobalIndex(
550 l, variable_id, component_id));
551
552 indices.push_back(global_index);
553 }
554 }
555 }
556 }
557
558 return indices;
559}
560
561} // namespace ProcessLib
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
Global vector based on Eigen vector.
Definition EigenVector.h:25
std::vector< Node * > const & getNodes() const
Get the nodes-vector for the mesh.
Definition Mesh.h:106
void setNodeID(std::size_t node_id)
void setCoordinates(MathLib::Point3d const &coordinates)
int getNumberOfGlobalComponents() const
Returns the number of components of the process variable.
virtual NumLib::IterationResult postIterationConcreteProcess(GlobalVector const &)
Definition Process.h:294
virtual 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)=0
virtual void preTimestepConcreteProcess(std::vector< GlobalVector * > const &, const double, const double, const int)
Definition Process.h:256
std::vector< BoundaryConditionCollection > _boundary_conditions
Definition Process.h:395
std::vector< NumLib::LocalToGlobalIndexMap const * > getDOFTables(int const number_of_processes) const
Definition Process.cpp:386
virtual void preIterationConcreteProcess(const unsigned, GlobalVector const &)
Definition Process.h:280
void constructDofTableOfSpecifiedProcessStaggeredScheme(const int specified_process_id)
Definition Process.cpp:355
static PROCESSLIB_EXPORT const std::string constant_one_parameter_name
Definition Process.h:46
void preAssemble(const double t, double const dt, GlobalVector const &x) final
Definition Process.cpp:259
virtual void initializeBoundaryConditions(std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
Definition Process.cpp:104
void initialize(std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
Definition Process.cpp:118
void initializeProcessBoundaryConditionsAndSourceTerms(const NumLib::LocalToGlobalIndexMap &dof_table, const int process_id, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
Definition Process.cpp:89
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:406
void assembleWithJacobian(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) final
Definition Process.cpp:286
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:265
virtual void postTimestepConcreteProcess(std::vector< GlobalVector * > const &, std::vector< GlobalVector * > const &, const double, const double, int const)
Definition Process.h:264
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:285
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_all_nodes
Definition Process.h:358
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:501
MeshLib::Mesh & _mesh
Definition Process.h:357
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:390
std::vector< std::size_t > const & getActiveElementIDs() const
Definition Process.h:160
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:495
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:454
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:444
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:476
void initializeExtrapolator()
Definition Process.cpp:421
virtual void postNonLinearSolverConcreteProcess(std::vector< GlobalVector * > const &, std::vector< GlobalVector * > const &, const double, double const, int const)
Definition Process.h:273
std::vector< std::reference_wrapper< ProcessVariable > > const & getProcessVariables(const int process_id) const
Definition Process.h:155
virtual void constructDofTable()
Definition Process.cpp:308
void constructMonolithicProcessDofTable()
Definition Process.cpp:322
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:466
std::vector< GlobalIndexType > getIndicesOfResiduumWithoutInitialCompensation() const override
Definition Process.cpp:513
unsigned const _integration_order
Definition Process.h:374
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
Definition Process.h:360
std::vector< SourceTermCollection > _source_term_collections
Definition Process.h:401
void updateDeactivatedSubdomains(double const time, const int process_id)
Definition Process.cpp:218
void computeSparsityPattern()
Definition Process.cpp:438
MathLib::MatrixSpecifications getMatrixSpecifications(const int process_id) const override
Definition Process.cpp:210
void setInitialConditions(std::vector< GlobalVector * > &process_solutions, std::vector< GlobalVector * > const &process_solutions_prev, double const t, int const process_id)
Definition Process.cpp:139
virtual void preOutputConcreteProcess(const double, double const, std::vector< GlobalVector * > const &, std::vector< GlobalVector * > const &, int const)
Definition Process.h:300
ExtrapolatorData _extrapolator_data
Definition Process.h:403
GlobalSparsityPattern _sparsity_pattern
Definition Process.h:382
void preIteration(const unsigned iter, GlobalVector const &x) final
Definition Process.cpp:489
const bool _use_monolithic_scheme
Definition Process.h:369
virtual std::tuple< NumLib::LocalToGlobalIndexMap *, bool > getDOFTableForExtrapolatorData() const
Definition Process.cpp:397
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:170
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.
void setLocalAccessibleVectors(std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev)
Definition Process.cpp:26