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