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