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