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>
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>
21 template <
class DofMap>
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>>;
49 : dofMap_(other.dofMap_)
56 : dofMap_(other.dofMap_)
57 , initialized_(other.initialized_)
59 if (other.initialized_) {
60 VecDuplicate(other.vector_, &vector_);
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);
106 return dofMap_.localSize();
112 return dofMap_.globalSize();
117 template <
class Basis>
118 void init(Basis
const&,
bool clear)
126 info(3,
" create new vector needed {} seconds", t.elapsed());
140 VecAssemblyBegin(vector_);
141 VecAssemblyEnd(vector_);
142 info(3,
" finish vector assembly needed {} seconds", t.elapsed());
150 VecGhostUpdateBegin(vector_, INSERT_VALUES, SCATTER_FORWARD);
151 VecGhostUpdateEnd(vector_, INSERT_VALUES, SCATTER_FORWARD);
152 info(3,
" synchronize vector needed {} seconds", t.elapsed());
157 template <
class MultiIndex>
158 PetscScalar
at(MultiIndex
const& idx)
const
161 gather(std::array<MultiIndex,1>{idx}, &y);
167 template <
class MultiIndex,
class Assign>
168 void insert(MultiIndex
const& dof, PetscScalar value, Assign)
170 assert(initialized_);
171 VecSetValue(vector_, dofMap_.global(flatMultiIndex(dof)), value, AssignMode<Assign>::value);
176 template <
class IndexRange,
class OutputIterator>
177 void gather(IndexRange
const& localInd, OutputIterator buffer)
const
179 forEach(localInd, [&buffer](
auto&&,
auto const& value) { *buffer = value; ++buffer; });
188 template <
class IndexRange,
class LocalValues,
class MaskRange,
class Assign>
189 void scatter(IndexRange
const& localInd, LocalValues
const& localVal, MaskRange
const& mask, Assign)
191 if constexpr(not std::is_same_v<typename LocalValues::value_type, PetscScalar>)
193 std::vector<PetscScalar> newLocalVal(localVal.begin(), localVal.end());
194 scatter(localInd, newLocalVal, mask, Assign{});
198 assert(initialized_);
199 assert(localInd.size() == localVal.size());
202 std::vector<PetscInt> globalInd;
203 globalInd.reserve(localInd.size());
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);
211 VecSetOption(vector_, VEC_IGNORE_NEGATIVE_INDICES, PETSC_TRUE);
213 VecSetValues(vector_, globalInd.size(), globalInd.data(), localVal.data(), AssignMode<Assign>::value);
214 VecSetOption(vector_, VEC_IGNORE_NEGATIVE_INDICES, PETSC_FALSE);
221 void axpy (
const PetscScalar& alpha,
const Self& x)
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);
232 void aypx (
const PetscScalar& alpha,
const Self& x)
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);
243 void scale (
const PetscScalar& alpha)
245 Vec localGhostedVec =
nullptr;
246 VecGhostGetLocalForm(vector_, &localGhostedVec);
247 VecScale(localGhostedVec, alpha);
248 VecGhostRestoreLocalForm(vector_, &localGhostedVec);
255 template <
class IndexRange,
class Func>
256 void forEach(IndexRange
const& localInd, Func&& f)
const
259 Vec localGhostedVec =
nullptr;
260 VecGhostGetLocalForm(vector_, &localGhostedVec);
261 assert(localGhostedVec !=
nullptr);
265 VecGetArray(localGhostedVec, &ptr);
267 for (
auto const& dof : localInd)
270 if (dofMap_.owner(i)) {
272 f(dof, ptr[dofMap_.dofToLocal(i)]);
276 f(dof, ptr[dofMap_.localSize() + dofMap_.dofToGhost(i)]);
280 VecRestoreArray(localGhostedVec, &ptr);
281 VecGhostRestoreLocalForm(vector_, &localGhostedVec);
286 template <
class Func>
290 PetscScalar
const* ptr;
291 VecGetArrayRead(vector_, &ptr);
294 f(ptr[dofMap_.dofToLocal(i)]);
296 VecRestoreArrayRead(vector_, &ptr);
302 template <
class Func>
307 VecGetArray(vector_, &ptr);
310 f(ptr[dofMap_.dofToLocal(i)]);
312 VecRestoreArray(vector_, &ptr);
317 template <
class Operation,
class... VecIns>
321 std::array<PetscScalar
const*,
sizeof...(VecIns)> ptrsIn{};
324 VecGetArray(vecOut.vector_, &ptrOut);
325 Ranges::forIndices<
sizeof...(VecIns)>([&](
auto ii) {
326 VecGetArrayRead(std::get<ii>(std::tie(vecIn...)).vector_, &ptrsIn[ii]);
330 auto idx = vecOut.dofMap_.dofToLocal(i);
331 Ranges::applyIndices<
sizeof...(VecIns)>([&](
auto... ii) {
332 ptrOut[idx] = op(ptrsIn[ii][idx]...);
337 VecRestoreArray(vecOut.vector_, &ptrOut);
338 Ranges::forIndices<
sizeof...(VecIns)>([&](
auto ii) {
339 VecRestoreArrayRead(std::get<ii>(std::tie(vecIn...)).vector_, &ptrsIn[ii]);
345 template <
class T,
class BinOp1,
class BinOp2>
348 PetscScalar
const* ptr1;
349 PetscScalar
const* ptr2;
352 VecGetArrayRead(in1.vector_, &ptr1);
353 VecGetArrayRead(in2.vector_, &ptr2);
356 auto idx = in1.dofMap_.dofToLocal(i);
357 init = op1(std::move(
init), op2(ptr1[idx], ptr2[idx]));
361 VecRestoreArrayRead(in1.vector_, &ptr1);
362 VecRestoreArrayRead(in2.vector_, &ptr2);
369 void set(PetscScalar value)
371 Vec localGhostedVec =
nullptr;
372 VecGhostGetLocalForm(vector_, &localGhostedVec);
373 assert(localGhostedVec !=
nullptr);
375 VecSet(localGhostedVec, value);
376 VecGhostRestoreLocalForm(vector_, &localGhostedVec);
383 Vec localGhostedVec =
nullptr;
384 VecGhostGetLocalForm(vector_, &localGhostedVec);
385 assert(localGhostedVec !=
nullptr);
387 VecZeroEntries(localGhostedVec);
388 VecGhostRestoreLocalForm(vector_, &localGhostedVec);
393 assert(&lhs.dofMap_ == &rhs.dofMap_);
394 std::swap(lhs.vector_, rhs.vector_);
395 std::swap(lhs.initialized_, rhs.initialized_);
399 MPI_Comm comm()
const
401 return dofMap_.comm();
408 VecCreateGhost(comm(),
409 dofMap_.localSize(), dofMap_.globalSize(), dofMap_.ghostSize(), dofMap_.ghostIndices().data(),
417 VecDestroy(&vector_);
419 initialized_ =
false;
424 DofMap
const& dofMap_;
427 mutable Vec vector_ =
nullptr;
428 bool initialized_ =
false;
434 template <
class DofMap>
435 struct ForEach<PetscVector<DofMap>>
437 template <
class Vec,
class UnaryFunction>
438 static void impl(Vec&& vec, UnaryFunction f)
444 template <
class DofMap>
445 struct Transform<PetscVector<DofMap>>
447 template <
class Operation,
class... VecIns>
448 static void impl(PetscVector<DofMap>& vecOut, Operation op, VecIns
const&... vecIn)
454 template <
class DofMap>
455 struct InnerProduct<PetscVector<DofMap>>
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)
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