OGS
PETScVector.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 "PETScVector.h"
5
6#include <algorithm>
7#include <cassert>
8
9#include "BaseLib/Algorithm.h"
10#include "BaseLib/Error.h"
11#include "BaseLib/MPI.h"
12
13namespace MathLib
14{
15PETScVector::PETScVector(const PetscInt vec_size, const bool is_global_size)
16{
17 if (is_global_size)
18 {
19 VecCreate(PETSC_COMM_WORLD, &v_);
20 VecSetSizes(v_, PETSC_DECIDE, vec_size);
21 }
22 else
23 {
24 // Fix size partitioning
25 // the size can be associated to specific memory allocation of a matrix
26 VecCreateMPI(PETSC_COMM_WORLD, vec_size, PETSC_DECIDE, &v_);
27 }
28
29 config();
30}
31
32PETScVector::PETScVector(const PetscInt vec_size,
33 const std::vector<PetscInt>& ghost_ids,
34 const bool is_global_size)
35 : size_ghosts_{static_cast<PetscInt>(ghost_ids.size())},
37{
38 PetscInt nghosts = static_cast<PetscInt>(ghost_ids.size());
39 if (is_global_size)
40 {
41 VecCreateGhost(PETSC_COMM_WORLD, PETSC_DECIDE, vec_size, nghosts,
42 ghost_ids.data(), &v_);
43 }
44 else
45 {
46 VecCreate(PETSC_COMM_WORLD, &v_);
47 VecSetType(v_, VECMPI);
48 VecSetSizes(v_, vec_size, PETSC_DECIDE);
49 VecMPISetGhost(v_, nghosts, ghost_ids.data());
50 }
51
52 config();
53
54 for (PetscInt i = 0; i < nghosts; i++)
55 {
56 global_ids2local_ids_ghost_.emplace(ghost_ids[i], size_loc_ + i);
57 }
58}
59
60PETScVector::PETScVector(const PETScVector& existing_vec, const bool deep_copy)
61{
62 shallowCopy(existing_vec);
63
64 // Copy values
65 if (deep_copy)
66 {
67 VecCopy(existing_vec.v_, v_);
68 }
69}
70
83
85{
86 VecSetFromOptions(v_);
87 VecGetOwnershipRange(v_, &start_rank_, &end_rank_);
88
89 VecGetLocalSize(v_, &size_loc_);
90 VecGetSize(v_, &size_);
91
92 VecSetOption(v_, VEC_IGNORE_NEGATIVE_INDICES, PETSC_TRUE);
93}
94
96{
97 VecAssemblyBegin(v_);
98 VecAssemblyEnd(v_);
99}
100
101void PETScVector::getGlobalVector(std::vector<PetscScalar>& u) const
102{
103#ifdef TEST_MEM_PETSC
104 PetscLogDouble mem1, mem2;
105 PetscMemoryGetCurrentUsage(&mem1);
106#endif
107
108 assert(static_cast<PetscInt>(u.size()) == size_);
109
110 PetscInt state;
111 VecLockGet(v_, &state);
112 if (state != 0)
113 {
114 OGS_FATAL("PETSc vector is already locked for {:s} access.",
115 state > 0 ? "read" : "write");
116 }
117 PetscScalar* xp = nullptr;
118 VecGetArray(v_, &xp);
119
120 BaseLib::MPI::Mpi mpi{PETSC_COMM_WORLD};
121 BaseLib::MPI::allgatherv(std::span(xp, size_loc_), u, mpi);
122
123 // This following line may be needed late on
124 // for a communication load balance:
125 // MPI_Barrier(PETSC_COMM_WORLD);
126
127 VecRestoreArray(v_, &xp);
128
129// TEST
130#ifdef TEST_MEM_PETSC
131 PetscMemoryGetCurrentUsage(&mem2);
132 PetscPrintf(
133 PETSC_COMM_WORLD,
134 "### Memory usage by Updating. Before :%f After:%f Increase:%d\n", mem1,
135 mem2, (int)(mem2 - mem1));
136#endif
137}
138
140{
142 {
144 return;
145 }
146
147 entry_array_.resize(size_);
149}
150
151void PETScVector::copyValues(std::vector<PetscScalar>& u) const
152{
153 u.resize(getLocalSize() + getGhostSize());
154
155 PetscScalar* loc_x = getLocalVector();
156 std::copy_n(loc_x, getLocalSize() + getGhostSize(), u.begin());
157 restoreArray(loc_x);
158}
159
160void PETScVector::copyValues(std::span<PetscScalar> u) const
161{
162 if (u.size() != static_cast<std::size_t>(getLocalSize() + getGhostSize()))
163 {
164 OGS_FATAL(
165 "PETScVector::copyValues() size mismatch. Trying to copy a vector "
166 "of size {:d} to a span of size {:d}.",
167 getLocalSize() + getGhostSize(), u.size());
168 }
169
170 PetscScalar* loc_x = getLocalVector();
171 std::copy_n(loc_x, getLocalSize() + getGhostSize(), u.begin());
172 restoreArray(loc_x);
173}
174
175PetscScalar PETScVector::get(const PetscInt idx) const
176{
178 {
179 return entry_array_[getLocalIndex(idx)];
180 }
181
182 return entry_array_[idx];
183}
184
185std::vector<PetscScalar> PETScVector::get(
186 std::vector<IndexType> const& indices) const
187{
188 std::vector<PetscScalar> local_x(indices.size());
189 std::transform(indices.begin(), indices.end(), local_x.begin(),
190 [this](IndexType index) { return get(index); });
191 return local_x;
192}
193
194PetscScalar* PETScVector::getLocalVector() const
195{
196 PetscScalar* loc_array;
198 {
199 VecGhostUpdateBegin(v_, INSERT_VALUES, SCATTER_FORWARD);
200 VecGhostUpdateEnd(v_, INSERT_VALUES, SCATTER_FORWARD);
201 VecGhostGetLocalForm(v_, &v_loc_);
202 VecGetArray(v_loc_, &loc_array);
203 }
204 else
205 {
206 VecGetArray(v_, &loc_array);
207 }
208
209 return loc_array;
210}
211
212PetscInt PETScVector::getLocalIndex(const PetscInt global_index) const
213{
214 if (global_index >= 0) // non-ghost entry.
215 {
216#ifndef NDEBUG
217 if (global_index < start_rank_ || global_index >= end_rank_)
218 {
219 OGS_FATAL(
220 "The global index {:d} is out of the range `[`{:d}, {:d}`)` of "
221 "the current rank.",
222 global_index, start_rank_, end_rank_);
223 }
224#endif
225 return global_index - start_rank_;
226 }
227
228 // A special case for a ghost location with global index equal to
229 // the size of the local vector:
230 PetscInt real_global_index = (-global_index == size_) ? 0 : -global_index;
231
232#ifndef NDEBUG
233 if (global_ids2local_ids_ghost_.find(real_global_index) ==
236 {
237 OGS_FATAL("The global index {:d} is not found as a ghost ID",
238 global_index);
239 }
240#endif
241
242 return global_ids2local_ids_ghost_.at(real_global_index);
243}
244
245void PETScVector::restoreArray(PetscScalar* array) const
246{
248 {
249 VecRestoreArray(v_loc_, &array);
250 VecGhostRestoreLocalForm(v_, &v_loc_);
251 }
252 else
253 {
254 VecRestoreArray(v_, &array);
255 }
256}
257
258void PETScVector::viewer(const std::string& file_name,
259 const PetscViewerFormat vw_format) const
260{
261 PetscViewer viewer;
262 PetscViewerASCIIOpen(PETSC_COMM_WORLD, file_name.c_str(), &viewer);
263 PetscViewerPushFormat(viewer, vw_format);
264
265 PetscObjectSetName((PetscObject)v_, file_name.c_str());
266 VecView(v_, viewer);
267
268#define nEXIT_TEST
269#ifdef EXIT_TEST
270 VecDestroy(v_);
271 PetscFinalize();
272 exit(0);
273#endif
274}
275
277{
278 destroy();
279
280 VecDuplicate(v.getRawVector(), &v_);
281
282 start_rank_ = v.start_rank_;
283 end_rank_ = v.end_rank_;
284 size_ = v.size_;
285 size_loc_ = v.size_loc_;
286 size_ghosts_ = v.size_ghosts_;
287 created_with_ghost_id_ = v.created_with_ghost_id_;
288 global_ids2local_ids_ghost_ = v.global_ids2local_ids_ghost_;
289
290 config();
291}
292
297
298} // namespace MathLib
#define OGS_FATAL(...)
Definition Error.h:19
PetscInt start_rank_
Starting index in a rank.
void finalizeAssembly()
Perform MPI collection of assembled entries in buffer.
std::vector< PetscScalar > entry_array_
Array containing the entries of the vector. If the vector is created without given ghost IDs,...
void config()
A function called by constructors to configure members.
void setLocalAccessibleVector() const
void restoreArray(PetscScalar *array) const
PetscInt getLocalSize() const
Get the number of entries in the same rank.
Definition PETScVector.h:79
PetscInt size_loc_
Size of local entries.
bool created_with_ghost_id_
Flag to indicate whether the vector is created with ghost entry indices.
PetscScalar * getLocalVector() const
void getGlobalVector(std::vector< PetscScalar > &u) const
std::vector< PetscScalar > get(std::vector< IndexType > const &indices) const
PetscInt size_ghosts_
Size of local ghost entries.
void copyValues(std::vector< PetscScalar > &u) const
PetscInt size_
Size of the vector.
PetscInt size() const
Get the global size of the vector.
Definition PETScVector.h:77
void shallowCopy(const PETScVector &v)
PetscInt getGhostSize() const
Get the number of ghost entries in the same rank.
Definition PETScVector.h:81
std::map< PetscInt, PetscInt > global_ids2local_ids_ghost_
Map global indices of ghost entries to local indices.
PetscInt getLocalIndex(const PetscInt global_index) const
Get local index by a global index.
PetscInt end_rank_
Ending index in a rank.
void viewer(const std::string &file_name, const PetscViewerFormat vw_format=PETSC_VIEWER_ASCII_MATLAB) const
static std::vector< int > allgatherv(std::span< T > const send_buffer, std::vector< std::remove_const_t< T > > &receive_buffer, Mpi const &mpi)
Definition MPI.h:158
void finalizeVectorAssembly(VEC_T &)
General function to finalize the vector assembly.
static const double u
static const double v