AMDiS 2.11-git
The Adaptive Multi-Dimensional Simulation Toolbox
 
Loading...
Searching...
No Matches
VectorBackend.hpp
1#pragma once
2
3#include <petscvec.h>
4
5#include <dune/common/ftraits.hh>
6#include <dune/common/rangeutilities.hh>
7#include <dune/common/shared_ptr.hh>
8#include <dune/functions/common/multiindex.hh>
9
10#include <amdis/Output.hpp>
11#include <amdis/algorithm/ForEach.hpp>
12#include <amdis/algorithm/InnerProduct.hpp>
13#include <amdis/algorithm/Transform.hpp>
14#include <amdis/common/FakeContainer.hpp>
15#include <amdis/operations/Assigner.hpp>
16#include <amdis/typetree/MultiIndex.hpp>
17
18namespace AMDiS
19{
21 template <class DofMap>
23 {
24 using Self = PetscVector;
25
26 template <class A>
27 using AssignMode = std::conditional_t<std::is_same_v<A,Assigner::plus_assign>,
28 std::integral_constant<::InsertMode, ADD_VALUES>,
29 std::integral_constant<::InsertMode, INSERT_VALUES>>;
30
31 public:
33 using BaseVector = ::Vec;
34
36 using value_type = PetscScalar;
37
39 using size_type = PetscInt;
40
41 public:
43 explicit PetscVector(DofMap const& dofMap)
44 : dofMap_(dofMap)
45 {}
46
49 : dofMap_(other.dofMap_)
50 {
51 swap(*this, other);
52 }
53
56 : dofMap_(other.dofMap_)
57 , initialized_(other.initialized_)
58 {
59 if (other.initialized_) {
60 VecDuplicate(other.vector_, &vector_);
61
62 // copy the data including the ghost values
63 Vec local, otherLocal;
64 VecGhostGetLocalForm(vector_, &local);
65 VecGhostGetLocalForm(other.vector_, &otherLocal);
66 VecCopy(otherLocal, local);
67 VecGhostRestoreLocalForm(vector_, &local);
68 VecGhostRestoreLocalForm(other.vector_, &otherLocal);
69 }
70 }
71
74 {
75 destroy();
76 }
77
80 {
81 swap(*this, other);
82 return *this;
83 }
84
87 {
88 return *this = PetscVector(other);
89 }
90
92 BaseVector const& vector() const
93 {
94 return vector_;
95 }
96
99 {
100 return vector_;
101 }
102
104 std::size_t localSize() const
105 {
106 return dofMap_.localSize();
107 }
108
110 std::size_t globalSize() const
111 {
112 return dofMap_.globalSize();
113 }
114
116 // [[possibly collective]]
117 template <class Basis>
118 void init(Basis const&, bool clear)
119 {
120 Dune::Timer t;
121 if (initialized_)
122 destroy();
123
124 // create the vector
125 create();
126 info(3, " create new vector needed {} seconds", t.elapsed());
127 t.reset();
128
129 if (clear)
130 setZero();
131
132 initialized_ = true;
133 }
134
136 // [[collective]]
137 void finish()
138 {
139 Dune::Timer t;
140 VecAssemblyBegin(vector_);
141 VecAssemblyEnd(vector_);
142 info(3, " finish vector assembly needed {} seconds", t.elapsed());
143 }
144
146 // [[collective]]
148 {
149 Dune::Timer t;
150 VecGhostUpdateBegin(vector_, INSERT_VALUES, SCATTER_FORWARD);
151 VecGhostUpdateEnd(vector_, INSERT_VALUES, SCATTER_FORWARD);
152 info(3, " synchronize vector needed {} seconds", t.elapsed());
153 }
154
156 // [[not collective]]
157 template <class MultiIndex>
158 PetscScalar at(MultiIndex const& idx) const
159 {
160 PetscScalar y = 0;
161 gather(std::array<MultiIndex,1>{idx}, &y);
162 return y;
163 }
164
166 // [[not collective]]
167 template <class MultiIndex, class Assign>
168 void insert(MultiIndex const& dof, PetscScalar value, Assign)
169 {
170 assert(initialized_);
171 VecSetValue(vector_, dofMap_.global(flatMultiIndex(dof)), value, AssignMode<Assign>::value);
172 }
173
175 // [[not collective]]
176 template <class IndexRange, class OutputIterator>
177 void gather(IndexRange const& localInd, OutputIterator buffer) const
178 {
179 forEach(localInd, [&buffer](auto&&, auto const& value) { *buffer = value; ++buffer; });
180 }
181
183
187 // [[not collective]]
188 template <class IndexRange, class LocalValues, class MaskRange, class Assign>
189 void scatter(IndexRange const& localInd, LocalValues const& localVal, MaskRange const& mask, Assign)
190 {
191 if constexpr(not std::is_same_v<typename LocalValues::value_type, PetscScalar>)
192 {
193 std::vector<PetscScalar> newLocalVal(localVal.begin(), localVal.end());
194 scatter(localInd, newLocalVal, mask, Assign{});
195 return;
196 }
197
198 assert(initialized_);
199 assert(localInd.size() == localVal.size());
200
201 // map local to global indices
202 std::vector<PetscInt> globalInd;
203 globalInd.reserve(localInd.size());
204
205 auto ind_it = std::begin(localInd);
206 auto mask_it = std::begin(mask);
207 for (; ind_it != std::end(localInd); ++mask_it, ++ind_it)
208 globalInd.push_back(*mask_it ? dofMap_.global(flatMultiIndex(*ind_it)) : -1);
209
210 if (!std::is_same_v<MaskRange, FakeContainer<bool,true>>)
211 VecSetOption(vector_, VEC_IGNORE_NEGATIVE_INDICES, PETSC_TRUE);
212
213 VecSetValues(vector_, globalInd.size(), globalInd.data(), localVal.data(), AssignMode<Assign>::value);
214 VecSetOption(vector_, VEC_IGNORE_NEGATIVE_INDICES, PETSC_FALSE);
215 }
216
219
221 void axpy (const PetscScalar& alpha, const Self& x)
222 {
223 Vec localGhostedVec = nullptr, x_localGhostedVec = nullptr;
224 VecGhostGetLocalForm(vector_, &localGhostedVec);
225 VecGhostGetLocalForm(x.vector_, &x_localGhostedVec);
226 VecAXPY(localGhostedVec, alpha, x_localGhostedVec);
227 VecGhostRestoreLocalForm(vector_, &localGhostedVec);
228 VecGhostRestoreLocalForm(x.vector_, &x_localGhostedVec);
229 }
230
232 void aypx (const PetscScalar& alpha, const Self& x)
233 {
234 Vec localGhostedVec = nullptr, x_localGhostedVec = nullptr;
235 VecGhostGetLocalForm(vector_, &localGhostedVec);
236 VecGhostGetLocalForm(x.vector_, &x_localGhostedVec);
237 VecAYPX(localGhostedVec, alpha, x_localGhostedVec);
238 VecGhostRestoreLocalForm(vector_, &localGhostedVec);
239 VecGhostRestoreLocalForm(x.vector_, &x_localGhostedVec);
240 }
241
243 void scale (const PetscScalar& alpha)
244 {
245 Vec localGhostedVec = nullptr;
246 VecGhostGetLocalForm(vector_, &localGhostedVec);
247 VecScale(localGhostedVec, alpha);
248 VecGhostRestoreLocalForm(vector_, &localGhostedVec);
249 }
250
252
254 // [[not collective]]
255 template <class IndexRange, class Func>
256 void forEach(IndexRange const& localInd, Func&& f) const
257 {
258 // Obtains the local ghosted representation of a parallel vector
259 Vec localGhostedVec = nullptr;
260 VecGhostGetLocalForm(vector_, &localGhostedVec);
261 assert(localGhostedVec != nullptr);
262
263 // A pointer to a contiguous array that contains this processor's portion of the vector data.
264 PetscScalar* ptr;
265 VecGetArray(localGhostedVec, &ptr);
266
267 for (auto const& dof : localInd)
268 {
269 size_type i = flatMultiIndex(dof);
270 if (dofMap_.owner(i)) {
271 // the index is a purely local entry
272 f(dof, ptr[dofMap_.dofToLocal(i)]);
273 }
274 else {
275 // ghost entries are stored at the end of the array
276 f(dof, ptr[dofMap_.localSize() + dofMap_.dofToGhost(i)]);
277 }
278 }
279
280 VecRestoreArray(localGhostedVec, &ptr);
281 VecGhostRestoreLocalForm(vector_, &localGhostedVec);
282 }
283
285 // [[not collective]]
286 template <class Func>
287 void forEachLocal(Func&& f) const
288 {
289 // A pointer to a contiguous array that contains this processor's portion of the vector data.
290 PetscScalar const* ptr;
291 VecGetArrayRead(vector_, &ptr);
292
293 for (size_type i = 0; i < size_type(localSize()); ++i)
294 f(ptr[dofMap_.dofToLocal(i)]);
295
296 VecRestoreArrayRead(vector_, &ptr);
297 }
298
300 // [[mutable]]
301 // [[not collective]]
302 template <class Func>
303 void forEachLocal(Func&& f)
304 {
305 // A pointer to a contiguous array that contains this processor's portion of the vector data.
306 PetscScalar *ptr;
307 VecGetArray(vector_, &ptr);
308
309 for (size_type i = 0; i < size_type(localSize()); ++i)
310 f(ptr[dofMap_.dofToLocal(i)]);
311
312 VecRestoreArray(vector_, &ptr);
313 }
314
316 // [[not collective]]
317 template <class Operation, class... VecIns>
318 static void transformLocal(PetscVector& vecOut, Operation op, VecIns const&... vecIn)
319 {
320 PetscScalar *ptrOut;
321 std::array<PetscScalar const*,sizeof...(VecIns)> ptrsIn{};
322
323 // Obtain ptrs to internal local arrays
324 VecGetArray(vecOut.vector_, &ptrOut);
325 Ranges::forIndices<sizeof...(VecIns)>([&](auto ii) {
326 VecGetArrayRead(std::get<ii>(std::tie(vecIn...)).vector_, &ptrsIn[ii]);
327 });
328
329 for (size_type i = 0; i < size_type(vecOut.localSize()); ++i) {
330 auto idx = vecOut.dofMap_.dofToLocal(i);
331 Ranges::applyIndices<sizeof...(VecIns)>([&](auto... ii) {
332 ptrOut[idx] = op(ptrsIn[ii][idx]...);
333 });
334 }
335
336 // Restore array obtained with VecGetArray[Read]()
337 VecRestoreArray(vecOut.vector_, &ptrOut);
338 Ranges::forIndices<sizeof...(VecIns)>([&](auto ii) {
339 VecRestoreArrayRead(std::get<ii>(std::tie(vecIn...)).vector_, &ptrsIn[ii]);
340 });
341 }
342
344 // [[not collective]]
345 template <class T, class BinOp1, class BinOp2>
346 static T innerProductLocal(PetscVector const& in1, PetscVector const& in2, T init, BinOp1 op1, BinOp2 op2)
347 {
348 PetscScalar const* ptr1;
349 PetscScalar const* ptr2;
350
351 // Obtain ptrs to internal local arrays
352 VecGetArrayRead(in1.vector_, &ptr1);
353 VecGetArrayRead(in2.vector_, &ptr2);
354
355 for (size_type i = 0; i < size_type(in1.localSize()); ++i) {
356 auto idx = in1.dofMap_.dofToLocal(i);
357 init = op1(std::move(init), op2(ptr1[idx], ptr2[idx]));
358 }
359
360 // Restore array obtained with VecGetArrayRead()
361 VecRestoreArrayRead(in1.vector_, &ptr1);
362 VecRestoreArrayRead(in2.vector_, &ptr2);
363
364 return init;
365 }
366
368 // [[not collective]]
369 void set(PetscScalar value)
370 {
371 Vec localGhostedVec = nullptr;
372 VecGhostGetLocalForm(vector_, &localGhostedVec);
373 assert(localGhostedVec != nullptr);
374
375 VecSet(localGhostedVec, value);
376 VecGhostRestoreLocalForm(vector_, &localGhostedVec);
377 }
378
380 // [[not collective]]
381 void setZero()
382 {
383 Vec localGhostedVec = nullptr;
384 VecGhostGetLocalForm(vector_, &localGhostedVec);
385 assert(localGhostedVec != nullptr);
386
387 VecZeroEntries(localGhostedVec);
388 VecGhostRestoreLocalForm(vector_, &localGhostedVec);
389 }
390
391 friend void swap(PetscVector& lhs, PetscVector& rhs)
392 {
393 assert(&lhs.dofMap_ == &rhs.dofMap_);
394 std::swap(lhs.vector_, rhs.vector_);
395 std::swap(lhs.initialized_, rhs.initialized_);
396 }
397
398 // An MPI Communicator or PETSC_COMM_SELF
399 MPI_Comm comm() const
400 {
401 return dofMap_.comm();
402 }
403
404 private:
405 // Create a new Vector
406 void create()
407 {
408 VecCreateGhost(comm(),
409 dofMap_.localSize(), dofMap_.globalSize(), dofMap_.ghostSize(), dofMap_.ghostIndices().data(),
410 &vector_);
411 }
412
413 // Destroy an old vector if created before
414 void destroy()
415 {
416 if (initialized_)
417 VecDestroy(&vector_);
418
419 initialized_ = false;
420 }
421
422 private:
423 // The local-to-global map
424 DofMap const& dofMap_;
425
427 mutable Vec vector_ = nullptr;
428 bool initialized_ = false;
429 };
430
431
432 namespace Recursive
433 {
434 template <class DofMap>
435 struct ForEach<PetscVector<DofMap>>
436 {
437 template <class Vec, class UnaryFunction>
438 static void impl(Vec&& vec, UnaryFunction f)
439 {
440 vec.forEachLocal(f);
441 }
442 };
443
444 template <class DofMap>
445 struct Transform<PetscVector<DofMap>>
446 {
447 template <class Operation, class... VecIns>
448 static void impl(PetscVector<DofMap>& vecOut, Operation op, VecIns const&... vecIn)
449 {
450 PetscVector<DofMap>::transformLocal(vecOut, op, vecIn...);
451 }
452 };
453
454 template <class DofMap>
455 struct InnerProduct<PetscVector<DofMap>>
456 {
457 template <class T, class BinOp1, class BinOp2>
458 static T impl(PetscVector<DofMap> const& in1, PetscVector<DofMap> const& in2, T init, BinOp1 op1, BinOp2 op2)
459 {
460 return PetscVector<DofMap>::innerProductLocal(in1, in2, std::move(init), op1, op2);
461 }
462 };
463
464 } // end namespace Recursive
465} // end namespace AMDiS
A container-like data-structure not storing anything and with empty implementations in many container...
Definition FakeContainer.hpp:36
The basic container that stores a base vector data.
Definition VectorBackend.hpp:23
PetscVector & operator=(PetscVector const &other)
Copy assignment operator.
Definition VectorBackend.hpp:86
void insert(MultiIndex const &dof, PetscScalar value, Assign)
Access the entry i of the vector with write-access.
Definition VectorBackend.hpp:168
void forEachLocal(Func &&f)
Perform the f(value) on the local elements of this vector.
Definition VectorBackend.hpp:303
PetscVector & operator=(PetscVector &&other)
Move assignment operator.
Definition VectorBackend.hpp:79
::Vec BaseVector
The type of the base vector.
Definition VectorBackend.hpp:33
void forEachLocal(Func &&f) const
Perform the f(value) on the local elements of this vector.
Definition VectorBackend.hpp:287
void axpy(const PetscScalar &alpha, const Self &x)
Implements *this = *this + alpha * x.
Definition VectorBackend.hpp:221
void setZero()
Zero all entries, including the ghost entries.
Definition VectorBackend.hpp:381
void synchronize()
Update the ghost regions with values from the owning process.
Definition VectorBackend.hpp:147
void finish()
Finish assembly. Must be called before vector can be used in a KSP.
Definition VectorBackend.hpp:137
PetscVector(PetscVector &&other)
Move constructor.
Definition VectorBackend.hpp:48
PetscScalar at(MultiIndex const &idx) const
Access the entry i of the vector with read-only-access.
Definition VectorBackend.hpp:158
static T innerProductLocal(PetscVector const &in1, PetscVector const &in2, T init, BinOp1 op1, BinOp2 op2)
Perform the op1(init, op2(value1, value2)) on the local elements of [in1] and [in2].
Definition VectorBackend.hpp:346
PetscVector(PetscVector const &other)
Copy constructor.
Definition VectorBackend.hpp:55
~PetscVector()
Destructor.
Definition VectorBackend.hpp:73
BaseVector & vector()
Return the data-vector vector_.
Definition VectorBackend.hpp:98
std::size_t globalSize() const
Return the number of entries in the global vector.
Definition VectorBackend.hpp:110
static void transformLocal(PetscVector &vecOut, Operation op, VecIns const &... vecIn)
Perform the valueOut = op(valueIn...) on the local elements of [vecOut] and [vecIn....
Definition VectorBackend.hpp:318
void init(Basis const &, bool clear)
Resize the vector_ to the size s.
Definition VectorBackend.hpp:118
void scatter(IndexRange const &localInd, LocalValues const &localVal, MaskRange const &mask, Assign)
Store all values given by localVal with indices localInd in the vector.
Definition VectorBackend.hpp:189
PetscVector(DofMap const &dofMap)
Constructor. Constructs new BaseVector.
Definition VectorBackend.hpp:43
std::size_t localSize() const
Return the number of entries in the local part of the vector.
Definition VectorBackend.hpp:104
PetscInt size_type
The index/size - type.
Definition VectorBackend.hpp:39
void aypx(const PetscScalar &alpha, const Self &x)
Implements *this = alpha * (*this) + x.
Definition VectorBackend.hpp:232
BaseVector const & vector() const
Return the data-vector vector_.
Definition VectorBackend.hpp:92
void set(PetscScalar value)
Set all entries to value, including the ghost entries.
Definition VectorBackend.hpp:369
void gather(IndexRange const &localInd, OutputIterator buffer) const
Collect all values from this vector to indices given by localInd and store it into the output buffer.
Definition VectorBackend.hpp:177
void scale(const PetscScalar &alpha)
scale all entries by the factor alpha
Definition VectorBackend.hpp:243
PetscScalar value_type
The type of the elements of the DOFVector.
Definition VectorBackend.hpp:36
void forEach(IndexRange const &localInd, Func &&f) const
Apply the functor f to each vector entry to local index in localInd range.
Definition VectorBackend.hpp:256