12#ifndef DUMUX_LINEAR_PARALLELHELPERS_HH
13#define DUMUX_LINEAR_PARALLELHELPERS_HH
15#include <dune/common/exceptions.hh>
16#include <dune/geometry/dimension.hh>
17#include <dune/grid/common/datahandleif.hh>
18#include <dune/grid/common/partitionset.hh>
19#include <dune/istl/owneroverlapcopy.hh>
20#include <dune/istl/paamg/pinfo.hh>
21#include <dune/istl/bvector.hh>
22#include <dune/istl/multitypeblockvector.hh>
28template<
class LinearSolverTraits,
bool canCommunicate = false>
31template<
class LinearSolverTraits>
34 using GridView =
typename LinearSolverTraits::GridView;
35 using DofMapper =
typename LinearSolverTraits::DofMapper;
36 static constexpr int dofCodim = LinearSolverTraits::dofCodim;
39 static constexpr std::size_t ghostMarker_ = 1<<24;
41 class BaseGatherScatter
44 BaseGatherScatter(
const DofMapper& mapper)
47 template<
class EntityType>
48 int index(
const EntityType& e)
const
49 {
return mapper_.index(e); }
51 bool contains(
int dim,
int codim)
const
52 {
return dofCodim == codim; }
55 bool fixedSize(
int dim,
int codim)
const
58 template<
class EntityType>
59 std::size_t size(EntityType& e)
const
62 template<
class EntityType>
63 bool isNeitherInteriorNorBorderEntity(EntityType& e)
const
64 {
return e.partitionType() != Dune::InteriorEntity && e.partitionType() != Dune::BorderEntity; }
67 const DofMapper& mapper_;
76 class GhostGatherScatter
77 :
public BaseGatherScatter
78 ,
public Dune::CommDataHandleIF<GhostGatherScatter, std::size_t>
81 using DataType = std::size_t;
82 using BaseGatherScatter::contains;
83 using BaseGatherScatter::fixedSize;
84 using BaseGatherScatter::size;
86 GhostGatherScatter(std::vector<std::size_t>& ranks,
const DofMapper& mapper)
87 : BaseGatherScatter(mapper)
91 template<
class MessageBuffer,
class EntityType>
92 void gather(MessageBuffer& buff,
const EntityType& e)
const
94 auto& data = ranks_[this->index(e)];
95 if (this->isNeitherInteriorNorBorderEntity(e))
100 template<
class MessageBuffer,
class EntityType>
101 void scatter(MessageBuffer& buff,
const EntityType& e, std::size_t n)
105 auto& data = ranks_[this->index(e)];
106 if (this->isNeitherInteriorNorBorderEntity(e))
110 std::vector<std::size_t>& ranks_;
119 class InteriorBorderGatherScatter
120 :
public BaseGatherScatter
121 ,
public Dune::CommDataHandleIF<InteriorBorderGatherScatter, std::size_t>
124 using DataType = std::size_t;
125 using BaseGatherScatter::contains;
126 using BaseGatherScatter::fixedSize;
127 using BaseGatherScatter::size;
129 InteriorBorderGatherScatter(std::vector<std::size_t>& ranks,
const DofMapper& mapper)
130 : BaseGatherScatter(mapper)
137 template<
class MessageBuffer,
class EntityType>
138 void gather(MessageBuffer& buff,
const EntityType& e)
const
140 auto& data = ranks_[this->index(e)];
141 if (this->isNeitherInteriorNorBorderEntity(e))
149 template<
class MessageBuffer,
class EntityType>
150 void scatter(MessageBuffer& buff,
const EntityType& e, std::size_t n)
154 auto& data = ranks_[this->index(e)];
160 data = this->isNeitherInteriorNorBorderEntity(e) ? x : min(data, x);
163 std::vector<std::size_t>& ranks_;
170 struct NeighbourGatherScatter
171 :
public BaseGatherScatter
172 ,
public Dune::CommDataHandleIF<NeighbourGatherScatter, int>
174 using DataType = int;
175 using BaseGatherScatter::contains;
176 using BaseGatherScatter::fixedSize;
177 using BaseGatherScatter::size;
179 NeighbourGatherScatter(
const DofMapper& mapper,
int rank, std::set<int>& neighbours)
180 : BaseGatherScatter(mapper)
182 , neighbours_(neighbours)
185 template<
class MessageBuffer,
class EntityType>
186 void gather(MessageBuffer& buff,
const EntityType& e)
const
191 template<
class MessageBuffer,
class EntityType>
192 void scatter(MessageBuffer& buff,
const EntityType& e, std::size_t n)
196 neighbours_.insert(x);
200 std::set<int>& neighbours_;
208 struct SharedGatherScatter
209 :
public BaseGatherScatter
210 ,
public Dune::CommDataHandleIF<SharedGatherScatter, int>
212 using DataType = int;
213 using BaseGatherScatter::contains;
214 using BaseGatherScatter::fixedSize;
215 using BaseGatherScatter::size;
217 SharedGatherScatter(std::vector<int>& shared,
const DofMapper& mapper)
218 : BaseGatherScatter(mapper)
222 template<
class MessageBuffer,
class EntityType>
223 void gather(MessageBuffer& buff, EntityType& e)
const
229 template<
class MessageBuffer,
class EntityType>
230 void scatter(MessageBuffer& buff,
const EntityType &e, std::size_t n)
234 auto& data = shared_[this->index(e)];
238 std::vector<int>& shared_;
245 template<
typename GlobalIndex>
246 struct GlobalIndexGatherScatter
247 :
public BaseGatherScatter
248 ,
public Dune::CommDataHandleIF<GlobalIndexGatherScatter<GlobalIndex>, GlobalIndex>
250 using DataType = GlobalIndex;
251 using BaseGatherScatter::contains;
252 using BaseGatherScatter::fixedSize;
253 using BaseGatherScatter::size;
255 GlobalIndexGatherScatter(std::vector<GlobalIndex>& globalIndices,
const DofMapper& mapper)
256 : BaseGatherScatter(mapper)
257 , globalIndices_(globalIndices)
260 template<
class MessageBuffer,
class EntityType>
261 void gather(MessageBuffer& buff,
const EntityType& e)
const
263 buff.write(globalIndices_[this->index(e)]);
266 template<
class MessageBuffer,
class EntityType>
267 void scatter(MessageBuffer& buff,
const EntityType& e, std::size_t n)
272 DataType& data = globalIndices_[this->index(e)];
276 std::vector<GlobalIndex>& globalIndices_;
282 : gridView_(gridView), mapper_(mapper)
284 if constexpr (Detail::canCommunicate<typename GridView::Traits::Grid, dofCodim>)
285 initGhostsAndOwners_();
287 DUNE_THROW(Dune::InvalidStateException,
288 "Cannot initialize parallel helper for a grid that cannot communicate codim-" << dofCodim <<
"-entities."
293 {
return isGhost_[i] == ghostMarker_; }
296 {
return isOwned_[i] == 1; }
312 if constexpr (Detail::canCommunicate<typename GridView::Traits::Grid, dofCodim>)
314 if (gridView_.comm().size() <= 1)
316 comm.remoteIndices().template rebuild<false>();
321 std::vector<int> isShared(mapper_.size(),
false);
322 SharedGatherScatter sgs(isShared, mapper_);
323 gridView_.communicate(sgs, Dune::All_All_Interface, Dune::ForwardCommunication);
326 using GlobalIndex =
typename Comm::ParallelIndexSet::GlobalIndex;
327 GlobalIndex count = 0;
328 for (std::size_t i = 0; i < isShared.size(); ++i)
329 if (isShared[i] && isOwned_[i] == 1)
332 std::vector<GlobalIndex> counts(gridView_.comm().size());
333 gridView_.comm().allgather(&count, 1, counts.data());
336 const int rank = gridView_.comm().rank();
337 auto start = std::accumulate(counts.begin(), counts.begin() + rank, GlobalIndex(0));
340 std::vector<GlobalIndex> globalIndices(mapper_.size(), std::numeric_limits<GlobalIndex>::max());
341 for (std::size_t i = 0; i < globalIndices.size(); ++i)
343 if (isOwned_[i] == 1 && isShared[i])
345 globalIndices[i] = start;
351 GlobalIndexGatherScatter<GlobalIndex> gigs(globalIndices, mapper_);
352 gridView_.communicate(gigs, Dune::All_All_Interface, Dune::ForwardCommunication);
355 resizeAndFillIndexSet_(comm, globalIndices);
358 std::set<int> neighbours;
359 NeighbourGatherScatter ngs(mapper_, gridView_.comm().rank(), neighbours);
360 gridView_.communicate(ngs, Dune::All_All_Interface, Dune::ForwardCommunication);
361 comm.remoteIndices().setNeighbours(neighbours);
363 comm.remoteIndices().template rebuild<false>();
366 DUNE_THROW(Dune::InvalidStateException,
367 "Cannot build parallel index set for a grid that cannot communicate codim-" << dofCodim <<
"-entities."
377 {
return gridView_; }
388 void initGhostsAndOwners_()
390 const auto rank = gridView_.comm().rank();
391 isOwned_.assign(mapper_.size(), rank);
394 GhostGatherScatter ggs(isOwned_, mapper_);
396 if (gridView_.comm().size() > 1)
397 gridView_.communicate(ggs, Dune::InteriorBorder_All_Interface, Dune::ForwardCommunication);
405 InteriorBorderGatherScatter dh(isOwned_, mapper_);
407 if (gridView_.comm().size() > 1)
408 gridView_.communicate(dh, Dune::InteriorBorder_InteriorBorder_Interface, Dune::ForwardCommunication);
411 for (
auto& v : isOwned_)
412 v = (v == rank) ? 1 : 0;
415 template<
class Comm,
class GlobalIndices>
416 void resizeAndFillIndexSet_(Comm& comm,
const GlobalIndices& globalIndices)
const
418 comm.indexSet().beginResize();
422 for (std::size_t localIndex = 0; localIndex < globalIndices.size(); ++localIndex)
424 const auto globalIndex = globalIndices[localIndex];
425 if (globalIndex != std::numeric_limits<typename GlobalIndices::value_type>::max())
427 const bool isOwned = isOwned_[localIndex] > 0;
428 const auto attr = getAttribute_(comm, isOwned, isGhost(localIndex));
429 using LocalIndex =
typename Comm::ParallelIndexSet::LocalIndex;
430 comm.indexSet().add(globalIndex, LocalIndex{localIndex, attr});
434 comm.indexSet().endResize();
438 Dune::OwnerOverlapCopyAttributeSet::AttributeSet
439 getAttribute_(
const Comm& comm,
const bool isOwned,
const bool isGhost)
const
442 return Dune::OwnerOverlapCopyAttributeSet::owner;
443 else if (isGhost && (comm.category() == Dune::SolverCategory::nonoverlapping) )
444 return Dune::OwnerOverlapCopyAttributeSet::overlap;
446 return Dune::OwnerOverlapCopyAttributeSet::copy;
449 const GridView gridView_;
450 const DofMapper& mapper_;
452 std::vector<std::size_t> isOwned_;
454 std::vector<std::size_t> isGhost_;
467template<
class LinearSolverTraits>
474template<
class Gr
idView,
class DofMapper,
int dofCodim>
479 : gridView_(gridView), mapper_(mapper)
483 template<
class Block,
class Alloc>
486 if constexpr (Detail::canCommunicate<typename GridView::Traits::Grid, dofCodim>)
489 if (gridView_.comm().size() > 1)
490 gridView_.communicate(gs, Dune::InteriorBorder_InteriorBorder_Interface,
491 Dune::ForwardCommunication);
494 DUNE_THROW(Dune::InvalidStateException,
"Cannot call makeNonOverlappingConsistent for a grid that cannot communicate codim-" << dofCodim <<
"-entities.");
498 template<
class... Blocks>
501 DUNE_THROW(Dune::NotImplemented,
"makeNonOverlappingConsistent for Dune::MultiTypeBlockVector");
505 const GridView gridView_;
506 const DofMapper& mapper_;
517template<
class Matrix,
class GridView,
518 class RowDofMapper,
int rowDofCodim>
521 using IdType =
typename GridView::Traits::Grid::Traits::GlobalIdSet::IdType;
544 template<
class ColIsGhostFunc>
545 struct MatrixPatternExchange
546 :
public Dune::CommDataHandleIF<MatrixPatternExchange<ColIsGhostFunc>, IdType>
549 using DataType = IdType;
561 MatrixPatternExchange(
const RowDofMapper& rowEntityMapper,
562 const std::map<IdType,int>& globalToLocal,
563 const std::map<int,IdType>& localToGlobal,
565 std::vector<std::set<int>>& sparsityPattern,
566 const ColIsGhostFunc& isGhostColumDof)
567 : rowEntityMapper_(rowEntityMapper), idToIndex_(globalToLocal), indexToID_(localToGlobal)
568 , sparsityPattern_(sparsityPattern), A_(A), isGhostColumDof_(isGhostColumDof)
570 sparsityPattern_.resize(A.N());
576 bool contains (
int dim,
int codim)
const
577 {
return (codim == rowDofCodim); }
580 bool fixedSize(
int dim,
int codim)
const
586 template<
class EntityType>
587 std::size_t size(EntityType& e)
const
589 const auto rowIdx = rowEntityMapper_.index(e);
593 for (
auto colIt = A_[rowIdx].begin(); colIt != A_[rowIdx].end(); ++colIt)
594 if (indexToID_.count(colIt.index()))
603 template<
class MessageBuffer,
class EntityType>
604 void gather(MessageBuffer& buff,
const EntityType& e)
const
606 const auto rowIdx = rowEntityMapper_.index(e);
609 for (
auto colIt = A_[rowIdx].begin(); colIt != A_[rowIdx].end(); ++colIt)
610 if (
auto it = indexToID_.find(colIt.index()); it != indexToID_.end())
611 buff.write(it->second);
617 template<
class MessageBuffer,
class EntityType>
618 void scatter(MessageBuffer& buff,
const EntityType& e, std::size_t n)
620 const auto rowIdx = rowEntityMapper_.index(e);
621 for (std::size_t k = 0; k < n; k++)
628 if (
const auto it = idToIndex_.find(
id); it != idToIndex_.end())
630 const auto colIdx = it->second;
633 if (!isGhostColumDof_(colIdx))
635 sparsityPattern_[rowIdx].insert(colIdx);
641 const RowDofMapper& rowEntityMapper_;
642 const std::map<IdType, int>& idToIndex_;
643 const std::map<int, IdType>& indexToID_;
644 std::vector<std::set<int>>& sparsityPattern_;
646 const ColIsGhostFunc& isGhostColumDof_;
654 typename Matrix::block_type second;
658 struct MatrixEntryExchange
659 :
public Dune::CommDataHandleIF<MatrixEntryExchange, MatrixEntry>
662 using DataType = MatrixEntry;
664 MatrixEntryExchange(
const RowDofMapper& rowEntityMapper,
665 const std::map<IdType, int>& globalToLocal,
666 const std::map<int, IdType>& localToGlobal,
668 : rowEntityMapper_(rowEntityMapper), idToIndex_(globalToLocal), indexToID_(localToGlobal), A_(A)
674 bool contains(
int dim,
int codim)
const
675 {
return (codim == rowDofCodim); }
678 bool fixedSize(
int dim,
int codim)
const
684 template<
class EntityType>
685 std::size_t size(EntityType& e)
const
687 const auto rowIdx = rowEntityMapper_.index(e);
689 for (
auto colIt = A_[rowIdx].begin(); colIt != A_[rowIdx].end(); ++colIt)
690 if (indexToID_.count(colIt.index()))
699 template<
class MessageBuffer,
class EntityType>
700 void gather(MessageBuffer& buff,
const EntityType& e)
const
702 const auto rowIdx = rowEntityMapper_.index(e);
704 for (
auto colIt = A_[rowIdx].begin(); colIt != A_[rowIdx].end(); ++colIt)
705 if (
auto it = indexToID_.find(colIt.index()); it != indexToID_.end())
706 buff.write(MatrixEntry{ it->second, *colIt });
712 template<
class MessageBuffer,
class EntityType>
713 void scatter(MessageBuffer& buff,
const EntityType& e, std::size_t n)
715 const auto rowIdx = rowEntityMapper_.index(e);
716 for (std::size_t k = 0; k < n; k++)
720 const auto& [colDofID, matrixBlock] = m;
723 if (
auto it = idToIndex_.find(colDofID); it != idToIndex_.end())
724 if (A_[rowIdx].find(it->second) != A_[rowIdx].end())
725 A_[rowIdx][it->second] += matrixBlock;
730 const RowDofMapper& rowEntityMapper_;
731 const std::map<IdType, int>& idToIndex_;
732 const std::map<int, IdType>& indexToID_;
740 : gridView_(gridView), mapper_(mapper)
745 std::vector<bool> handledDof(gridView_.size(rowDofCodim),
false);
746 for (
const auto& element : elements(gridView_))
748 for (
int i = 0; i < element.subEntities(rowDofCodim); ++i)
750 const auto entity = element.template subEntity<rowDofCodim>(i);
751 if (entity.partitionType() == Dune::BorderEntity)
753 const auto localRowIdx = mapper_.index(entity);
754 if (!handledDof[localRowIdx])
756 IdType dofIdxGlobal = gridView_.grid().globalIdSet().id(entity);
757 idToIndex_.emplace(dofIdxGlobal, localRowIdx);
758 indexToID_.emplace(localRowIdx, dofIdxGlobal);
770 template<
class IsGhostFunc>
773 if constexpr (Detail::canCommunicate<typename GridView::Grid, rowDofCodim>)
775 if (gridView_.comm().size() <= 1)
779 Matrix matrixAsAssembled(A);
782 std::size_t numNonZeroEntries = 0;
783 std::vector<std::set<int>> sparsityPattern;
784 MatrixPatternExchange<IsGhostFunc> dataHandle(mapper_, idToIndex_, indexToID_, A, sparsityPattern, isGhost);
785 gridView_.communicate(
786 dataHandle, Dune::InteriorBorder_InteriorBorder_Interface, Dune::ForwardCommunication
790 for (
auto rowIt = A.begin(); rowIt != A.end(); ++rowIt)
792 const auto colEndIt = A[rowIt.index()].end();
793 for (
auto colIt = rowIt->begin(); colIt != colEndIt; ++colIt)
794 if (!sparsityPattern[rowIt.index()].count(colIt.index()))
795 sparsityPattern[rowIt.index()].insert(colIt.index());
797 numNonZeroEntries += sparsityPattern[rowIt.index()].size();
801 A.setSize(matrixAsAssembled.N(), matrixAsAssembled.M(), numNonZeroEntries);
802 A.setBuildMode(Matrix::row_wise);
803 auto citer = A.createbegin();
804 for (
auto i = sparsityPattern.begin(), end = sparsityPattern.end(); i!=end; ++i, ++citer)
805 for (
auto si = i->begin(), send = i->end(); si!=send; ++si)
810 const auto rowEndIt = matrixAsAssembled.end();
811 for (
auto rowIt = matrixAsAssembled.begin(); rowIt != rowEndIt; ++rowIt)
812 for (
auto colIt = matrixAsAssembled[rowIt.index()].begin(); colIt != matrixAsAssembled[rowIt.index()].end(); ++colIt)
813 A[rowIt.index()][colIt.index()] = *colIt;
819 DUNE_THROW(Dune::InvalidStateException,
"Cannot call extendMatrix for a grid that cannot communicate codim-" << rowDofCodim <<
"-entities.");
837 if constexpr (Detail::canCommunicate<typename GridView::Grid, rowDofCodim>)
839 if (gridView_.comm().size() <= 1)
842 MatrixEntryExchange dataHandle(mapper_, idToIndex_, indexToID_, A);
843 gridView_.communicate(
844 dataHandle, Dune::InteriorBorder_InteriorBorder_Interface, Dune::ForwardCommunication
848 DUNE_THROW(Dune::InvalidStateException,
849 "Cannot call sumEntries for a grid that cannot communicate codim-" << rowDofCodim <<
"-entities."
854 const GridView gridView_;
855 const RowDofMapper& mapper_;
856 std::map<IdType, int> idToIndex_;
857 std::map<int, IdType> indexToID_;
865 class Matrix,
class ParallelHelper>
868 if constexpr (ParallelTraits::isNonOverlapping)
872 using GridView =
typename LinearSolverTraits::GridView;
873 using DofMapper =
typename LinearSolverTraits::DofMapper;
874 static constexpr int dofCodim = LinearSolverTraits::dofCodim;
876 matrixHelper.
extendMatrix(A, [&pHelper](
auto idx){
return pHelper.isGhost(idx); });
885 class Vector,
class ParallelHelper>
888 if constexpr (ParallelTraits::isNonOverlapping)
892 using GridView =
typename LinearSolverTraits::GridView;
893 using DofMapper =
typename LinearSolverTraits::DofMapper;
894 static constexpr int dofCodim = LinearSolverTraits::dofCodim;
904 class Matrix,
class Vector,
class ParallelHelper>
907 prepareMatrixParallel<LinearSolverTraits, ParallelTraits>(A, pHelper);
908 prepareVectorParallel<LinearSolverTraits, ParallelTraits>(b, pHelper);
915 class Matrix,
class Vector,
class ParallelHelper>
917 std::shared_ptr<typename ParallelTraits::Comm>& comm,
918 std::shared_ptr<typename ParallelTraits::LinearOperator>& fop,
919 std::shared_ptr<typename ParallelTraits::ScalarProduct>& sp,
920 ParallelHelper& pHelper)
922 prepareLinearAlgebraParallel<LinearSolverTraits, ParallelTraits>(A, b, pHelper);
923 const auto category = ParallelTraits::isNonOverlapping ?
924 Dune::SolverCategory::nonoverlapping : Dune::SolverCategory::overlapping;
926 comm = std::make_shared<typename ParallelTraits::Comm>(pHelper.gridView().comm(), category);
927 pHelper.createParallelIndexSet(*comm);
928 fop = std::make_shared<typename ParallelTraits::LinearOperator>(A, *comm);
929 sp = std::make_shared<typename ParallelTraits::ScalarProduct>(*comm);
const DofMapper & dofMapper() const
Return the dofMapper.
Definition: parallelhelpers.hh:372
ParallelISTLHelperImpl(const GridView &gridView, const DofMapper &mapper)
Definition: parallelhelpers.hh:281
bool isOwned(std::size_t i) const
Definition: parallelhelpers.hh:295
const GridView & gridView() const
Return the gridView.
Definition: parallelhelpers.hh:376
bool isGhost(std::size_t i) const
Definition: parallelhelpers.hh:292
void createParallelIndexSet(Comm &comm) const
Creates a parallel index set for a communicator.
Definition: parallelhelpers.hh:310
Definition: parallelhelpers.hh:29
Helper class for adding up matrix entries for border entities.
Definition: parallelhelpers.hh:520
ParallelMatrixHelper(const GridView &gridView, const RowDofMapper &mapper)
Definition: parallelhelpers.hh:739
void sumEntries(Matrix &A)
Sums up the entries corresponding to border entities (usually vertices or faces)
Definition: parallelhelpers.hh:835
void extendMatrix(Matrix &A, const IsGhostFunc &isGhost)
communicates values for the sparsity pattern of the new matrix.
Definition: parallelhelpers.hh:771
Definition: parallelhelpers.hh:476
void makeNonOverlappingConsistent(Dune::MultiTypeBlockVector< Blocks... > &v) const
Make a vector consistent for non-overlapping domain decomposition methods.
Definition: parallelhelpers.hh:499
ParallelVectorHelper(const GridView &gridView, const DofMapper &mapper)
Definition: parallelhelpers.hh:478
void makeNonOverlappingConsistent(Dune::BlockVector< Block, Alloc > &v) const
Make a vector consistent for non-overlapping domain decomposition methods.
Definition: parallelhelpers.hh:484
A data handle class to exchange entries of a vector.
Definition: vectorcommdatahandle.hh:66
Definition: variablesbackend.hh:34
dune-grid capabilities compatibility layer
Distance implementation details.
Definition: cvfelocalresidual.hh:25
static constexpr bool canCommunicate
Definition: gridcapabilities.hh:51
void prepareMatrixParallel(Matrix &A, ParallelHelper &pHelper)
Prepare a matrix for parallel solvers.
Definition: parallelhelpers.hh:866
void prepareLinearAlgebraParallel(Matrix &A, Vector &b, ParallelHelper &pHelper)
Prepare linear algebra variables for parallel solvers.
Definition: parallelhelpers.hh:905
LinearSolverTraitsImpl< GridGeometry, typename GridGeometry::DiscretizationMethod > LinearSolverTraits
The type traits required for using the IstlFactoryBackend.
Definition: linearsolvertraits.hh:37
void prepareVectorParallel(Vector &b, ParallelHelper &pHelper)
Prepare a vector for parallel solvers.
Definition: parallelhelpers.hh:886
The implementation is specialized for the different discretizations.
Definition: linearsolvertraits.hh:31
Contains a class to exchange entries of a vector.