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 cell_average_data_(mesh),
58 _jacobian_assembler(std::move(jacobian_assembler)),
59 _global_assembler(*_jacobian_assembler),
60 _use_monolithic_scheme(use_monolithic_scheme),
61 _integration_order(integration_order),
62 _process_variables(std::move(process_variables)),
63 _boundary_conditions(
64 [&](const std::size_t number_of_process_variables)
65 -> std::vector<BoundaryConditionCollection>
66 {
67 std::vector<BoundaryConditionCollection> pcs_BCs;
68 pcs_BCs.reserve(number_of_process_variables);
69 for (std::size_t i = 0; i < number_of_process_variables; i++)
70 {
71 pcs_BCs.emplace_back(BoundaryConditionCollection(parameters));
72 }
73 return pcs_BCs;
74 }(_process_variables.size())),
75 _source_term_collections(
76 [&](const std::size_t number_of_processes)
77 -> std::vector<SourceTermCollection>
78 {
79 std::vector<SourceTermCollection> pcs_sts;
80 pcs_sts.reserve(number_of_processes);
81 for (std::size_t i = 0; i < number_of_processes; i++)
82 {
83 pcs_sts.emplace_back(SourceTermCollection(parameters));
84 }
85 return pcs_sts;
86 }(_process_variables.size()))
87{
88}
89
91 const NumLib::LocalToGlobalIndexMap& dof_table, const int process_id,
92 std::map<int, std::shared_ptr<MaterialPropertyLib::Medium>> const& media)
93{
94 auto const& per_process_variables = _process_variables[process_id];
95 auto& per_process_BCs = _boundary_conditions[process_id];
96
97 per_process_BCs.addBCsForProcessVariables(per_process_variables, dof_table,
98 _integration_order, *this, media);
99
100 auto& per_process_sts = _source_term_collections[process_id];
101 per_process_sts.addSourceTermsForProcessVariables(
102 per_process_variables, dof_table, _integration_order);
103}
104
106 std::map<int, std::shared_ptr<MaterialPropertyLib::Medium>> const& media)
107{
108 // The number of processes is identical to the size of _process_variables,
109 // the vector contains variables for different processes. See the
110 // documentation of _process_variables.
111 const std::size_t number_of_processes = _process_variables.size();
112 for (std::size_t pcs_id = 0; pcs_id < number_of_processes; pcs_id++)
113 {
115 *_local_to_global_index_map, pcs_id, media);
116 }
117}
118
120 std::map<int, std::shared_ptr<MaterialPropertyLib::Medium>> const& media)
121{
122 DBUG("Initialize process.");
123
124 DBUG("Construct dof mappings.");
126
127 DBUG("Compute sparsity pattern");
129
130 DBUG("Initialize the extrapolator");
132
135
136 DBUG("Initialize boundary conditions.");
138}
139
141 std::vector<GlobalVector*>& process_solutions,
142 std::vector<GlobalVector*> const& process_solutions_prev,
143 double const t,
144 int const process_id)
145{
146 auto& x = *process_solutions[process_id];
147 auto& x_prev = *process_solutions_prev[process_id];
148
149 // getDOFTableOfProcess can be overloaded by the specific process.
150 auto const& dof_table_of_process = getDOFTable(process_id);
151
152 auto const& per_process_variables = _process_variables[process_id];
153 for (std::size_t variable_id = 0;
154 variable_id < per_process_variables.size();
155 variable_id++)
156 {
159
160 auto const& pv = per_process_variables[variable_id];
161 DBUG("Set the initial condition of variable {:s} of process {:d}.",
162 pv.get().getName().data(), process_id);
163
164 auto const& ic = pv.get().getInitialCondition();
165
166 auto const num_comp = pv.get().getNumberOfGlobalComponents();
167
168 for (int component_id = 0; component_id < num_comp; ++component_id)
169 {
170 auto const& mesh_subset =
171 dof_table_of_process.getMeshSubset(variable_id, component_id);
172 auto const mesh_id = mesh_subset.getMeshID();
173 for (auto const* node : mesh_subset.getNodes())
174 {
176 node->getID());
177
178 pos.setNodeID(node->getID());
179 pos.setCoordinates(*node);
180 auto const& ic_value = ic(t, pos);
181
182 auto global_index =
183 std::abs(dof_table_of_process.getGlobalIndex(l, variable_id,
184 component_id));
185#ifdef USE_PETSC
186 // The global indices of the ghost entries of the global
187 // matrix or the global vectors need to be set as negative
188 // values for equation assembly, however the global indices
189 // start from zero. Therefore, any ghost entry with zero
190 // index is assigned an negative value of the vector size
191 // or the matrix dimension. To assign the initial value for
192 // the ghost entries, the negative indices of the ghost
193 // entries are restored to zero.
194 if (global_index == x.size())
195 global_index = 0;
196#endif
197 x.set(global_index, ic_value[component_id]);
198 }
199 }
200 }
201
203 MathLib::LinAlg::copy(x, x_prev); // pushState
204
207
208 setInitialConditionsConcreteProcess(process_solutions, t, process_id);
209}
210
212 const int /*process_id*/) const
213{
214 auto const& l = *_local_to_global_index_map;
215 return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
216 &l.getGhostIndices(), &_sparsity_pattern};
217}
218
220 const int process_id)
221{
222 auto const& variables_per_process = getProcessVariables(process_id);
223 for (auto const& variable : variables_per_process)
224 {
225 variable.get().updateDeactivatedSubdomains(time);
226 }
227
229
230 auto active_elements_ids = ranges::views::transform(
231 [](auto const& variable)
232 { return variable.get().getActiveElementIDs(); });
233
234 // Early return if there's any process variable with all elements active.
235 if (ranges::any_of(variables_per_process | active_elements_ids,
236 [](auto const& vector) { return vector.empty(); }))
237 {
238 return;
239 }
240
241 // Some process variable has deactivated elements. Create union of active
242 // ids.
243
244 _ids_of_active_elements = // there is at least one process variable.
245 variables_per_process[0].get().getActiveElementIDs();
246
247 for (auto const& pv_active_element_ids :
248 variables_per_process | ranges::views::drop(1) | active_elements_ids)
249 {
250 std::vector<std::size_t> new_active_elements;
251 new_active_elements.reserve(_ids_of_active_elements.size() +
252 pv_active_element_ids.size());
253 ranges::set_union(_ids_of_active_elements, pv_active_element_ids,
254 std::back_inserter(new_active_elements));
255
256 _ids_of_active_elements = std::move(new_active_elements);
257 }
258}
259
260void Process::preAssemble(const double t, double const dt,
261 GlobalVector const& x)
262{
264}
265
266void Process::assemble(const double t, double const dt,
267 std::vector<GlobalVector*> const& x,
268 std::vector<GlobalVector*> const& x_prev,
269 int const process_id, GlobalMatrix& M, GlobalMatrix& K,
270 GlobalVector& b)
271{
272 assert(x.size() == x_prev.size());
273
274 setLocalAccessibleVectors(x, x_prev);
275
276 assembleConcreteProcess(t, dt, x, x_prev, process_id, M, K, b);
277
278 // the last argument is for the jacobian, nullptr is for a unused jacobian
279 _boundary_conditions[process_id].applyNaturalBC(t, x, process_id, &K, b,
280 nullptr);
281
282 // the last argument is for the jacobian, nullptr is for a unused jacobian
283 _source_term_collections[process_id].integrate(t, *x[process_id], b,
284 nullptr);
285}
286
287void Process::assembleWithJacobian(const double t, double const dt,
288 std::vector<GlobalVector*> const& x,
289 std::vector<GlobalVector*> const& x_prev,
290 int const process_id, GlobalVector& b,
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, b, Jac);
298
299 // TODO: apply BCs to Jacobian.
300 _boundary_conditions[process_id].applyNaturalBC(t, x, process_id, nullptr,
301 b, &Jac);
302
303 // the last argument is for the jacobian, nullptr is for a unused jacobian
304 _source_term_collections[process_id].integrate(t, *x[process_id], b, &Jac);
305}
306
308{
310 {
312
313 return;
314 }
315
316 // For staggered scheme:
317 const int specified_process_id = 0;
319}
320
322{
323 // Create single component dof in every of the mesh nodes.
325 std::make_unique<MeshLib::MeshSubset>(_mesh, _mesh.getNodes());
326
327 // Vector of mesh subsets.
328 std::vector<MeshLib::MeshSubset> all_mesh_subsets;
329
330 // Collect the mesh subsets in a vector for the components of each
331 // variables.
332 for (ProcessVariable const& pv : _process_variables[0])
333 {
334 std::generate_n(std::back_inserter(all_mesh_subsets),
335 pv.getNumberOfGlobalComponents(),
336 [&]() { return *_mesh_subset_all_nodes; });
337 }
338
339 // Create a vector of the number of variable components
340 std::vector<int> vec_var_n_components;
341 transform(cbegin(_process_variables[0]), cend(_process_variables[0]),
342 back_inserter(vec_var_n_components),
343 [](ProcessVariable const& pv)
344 { return pv.getNumberOfGlobalComponents(); });
345
347 std::make_unique<NumLib::LocalToGlobalIndexMap>(
348 std::move(all_mesh_subsets), vec_var_n_components,
350
352}
353
355 const int specified_process_id)
356{
357 // Create single component dof in every of the mesh nodes.
359 std::make_unique<MeshLib::MeshSubset>(_mesh, _mesh.getNodes());
360
361 // Vector of mesh subsets.
362 std::vector<MeshLib::MeshSubset> all_mesh_subsets;
363
364 // Vector of the number of variable components
365 std::vector<int> vec_var_n_components;
366 // Collect the mesh subsets in a vector for each variables' components.
367 std::generate_n(std::back_inserter(all_mesh_subsets),
368 _process_variables[specified_process_id][0]
369 .get()
370 .getNumberOfGlobalComponents(),
371 [&]() { return *_mesh_subset_all_nodes; });
372
373 // Create a vector of the number of variable components.
374 vec_var_n_components.push_back(_process_variables[specified_process_id][0]
375 .get()
376 .getNumberOfGlobalComponents());
378 std::make_unique<NumLib::LocalToGlobalIndexMap>(
379 std::move(all_mesh_subsets), vec_var_n_components,
381
383}
384
385std::vector<NumLib::LocalToGlobalIndexMap const*> Process::getDOFTables(
386 int const number_of_processes) const
387{
388 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
389 dof_tables.reserve(number_of_processes);
390 std::generate_n(std::back_inserter(dof_tables), number_of_processes,
391 [&]() { return &getDOFTable(dof_tables.size()); });
392 return dof_tables;
393}
394
395std::tuple<NumLib::LocalToGlobalIndexMap*, bool>
397{
398 if (_local_to_global_index_map->getNumberOfGlobalComponents() == 1)
399 {
400 // For single-variable-single-component processes reuse the existing DOF
401 // table.
402 const bool manage_storage = false;
403 return std::make_tuple(_local_to_global_index_map.get(),
404 manage_storage);
405 }
406
407 // Otherwise construct a new DOF table.
408 std::vector<MeshLib::MeshSubset> all_mesh_subsets_single_component;
409 all_mesh_subsets_single_component.emplace_back(*_mesh_subset_all_nodes);
410
411 const bool manage_storage = true;
412
413 return std::make_tuple(new NumLib::LocalToGlobalIndexMap(
414 std::move(all_mesh_subsets_single_component),
415 // by location order is needed for output
417 manage_storage);
418}
419
421{
422 NumLib::LocalToGlobalIndexMap* dof_table_single_component;
423 bool manage_storage;
424
425 std::tie(dof_table_single_component, manage_storage) =
427
428 std::unique_ptr<NumLib::Extrapolator> extrapolator(
430 *dof_table_single_component));
431
432 // TODO: Later on the DOF table can change during the simulation!
434 std::move(extrapolator), dof_table_single_component, manage_storage);
435}
436
442
443void Process::preTimestep(std::vector<GlobalVector*> const& x, const double t,
444 const double delta_t, const int process_id)
445{
446 for (auto* const solution : x)
448 preTimestepConcreteProcess(x, t, delta_t, process_id);
449
450 _boundary_conditions[process_id].preTimestep(t, x, process_id);
451}
452
453void Process::postTimestep(std::vector<GlobalVector*> const& x,
454 std::vector<GlobalVector*> const& x_prev,
455 const double t, const double delta_t,
456 int const process_id)
457{
458 setLocalAccessibleVectors(x, x_prev);
459
460 postTimestepConcreteProcess(x, x_prev, t, delta_t, process_id);
461
462 _boundary_conditions[process_id].postTimestep(t, x, process_id);
463}
464
465void Process::postNonLinearSolver(std::vector<GlobalVector*> const& x,
466 std::vector<GlobalVector*> const& x_prev,
467 const double t, double const dt,
468 int const process_id)
469{
470 setLocalAccessibleVectors(x, x_prev);
471
472 postNonLinearSolverConcreteProcess(x, x_prev, t, dt, process_id);
473}
474
476 double const dt,
477 std::vector<GlobalVector*> const& x,
478 GlobalVector const& x_prev,
479 int const process_id)
480{
481 for (auto const* solution : x)
484
485 computeSecondaryVariableConcrete(t, dt, x, x_prev, process_id);
486}
487
488void Process::preIteration(const unsigned iter, const GlobalVector& x)
489{
492}
493
499
500void Process::preOutput(const double t, double const dt,
501 std::vector<GlobalVector*> const& x,
502 std::vector<GlobalVector*> const& x_prev,
503 int const process_id)
504{
505 for (auto const* solution : x)
507
508 preOutputConcreteProcess(t, dt, x, x_prev, process_id);
509}
510
511std::vector<GlobalIndexType>
513{
514 std::vector<GlobalIndexType> indices;
515
516 for (std::size_t process_id = 0; process_id < _process_variables.size();
517 process_id++)
518 {
519 auto const& dof_table_of_process = getDOFTable(process_id);
520
521 auto const& per_process_variables = _process_variables[process_id];
522 for (std::size_t variable_id = 0;
523 variable_id < per_process_variables.size();
524 variable_id++)
525 {
526 auto const& pv = per_process_variables[variable_id];
527 DBUG("Set the initial condition of variable {:s} of process {:d}.",
528 pv.get().getName().data(), process_id);
529
530 if ((pv.get().compensateNonEquilibriumInitialResiduum()))
531 {
532 continue;
533 }
534
535 auto const num_comp = pv.get().getNumberOfGlobalComponents();
536
537 for (int component_id = 0; component_id < num_comp; ++component_id)
538 {
539 auto const& mesh_subset = dof_table_of_process.getMeshSubset(
540 variable_id, component_id);
541 auto const mesh_id = mesh_subset.getMeshID();
542 for (auto const* node : mesh_subset.getNodes())
543 {
544 MeshLib::Location const l(
545 mesh_id, MeshLib::MeshItemType::Node, node->getID());
546
547 auto global_index =
548 std::abs(dof_table_of_process.getGlobalIndex(
549 l, variable_id, component_id));
550
551 indices.push_back(global_index);
552 }
553 }
554 }
555 }
556
557 return indices;
558}
559
560} // 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:302
virtual void preTimestepConcreteProcess(std::vector< GlobalVector * > const &, const double, const double, const int)
Definition Process.h:264
std::vector< BoundaryConditionCollection > _boundary_conditions
Definition Process.h:405
std::vector< NumLib::LocalToGlobalIndexMap const * > getDOFTables(int const number_of_processes) const
Definition Process.cpp:385
virtual void preIterationConcreteProcess(const unsigned, GlobalVector const &)
Definition Process.h:288
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 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:416
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
virtual void postTimestepConcreteProcess(std::vector< GlobalVector * > const &, std::vector< GlobalVector * > const &, const double, const double, int const)
Definition Process.h:272
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:293
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_all_nodes
Definition Process.h:366
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:365
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:400
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:241
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
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:281
virtual void constructDofTable()
Definition Process.cpp:307
void constructMonolithicProcessDofTable()
Definition Process.cpp:321
virtual void preAssembleConcreteProcess(const double, double const, GlobalVector const &)
Definition Process.h:248
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:384
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
Definition Process.h:368
std::vector< SourceTermCollection > _source_term_collections
Definition Process.h:411
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:308
std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > const & getProcessVariables() const
Definition Process.h:156
ExtrapolatorData _extrapolator_data
Definition Process.h:413
GlobalSparsityPattern _sparsity_pattern
Definition Process.h:392
void preIteration(const unsigned iter, GlobalVector const &x) final
Definition Process.cpp:488
const bool _use_monolithic_scheme
Definition Process.h:379
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.
void finalizeAssembly(PETScMatrix &A)
Definition LinAlg.cpp:198
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