27#ifndef EWOMS_OVERLAPPING_BLOCK_VECTOR_HH
28#define EWOMS_OVERLAPPING_BLOCK_VECTOR_HH
33#include <opm/material/common/Valgrind.hpp>
35#include <dune/istl/bvector.hh>
36#include <dune/common/fvector.hh>
48template <
class FieldVector,
class Overlap>
51 using ParentType = Dune::BlockVector<FieldVector>;
52 using BlockVector = Dune::BlockVector<FieldVector>;
60 : ParentType(overlap.numDomestic()), overlap_(&overlap)
68 , numIndicesSendBuff_(obv.numIndicesSendBuff_)
69 , indicesSendBuff_(obv.indicesSendBuff_)
70 , indicesRecvBuff_(obv.indicesRecvBuff_)
71 , valuesSendBuff_(obv.valuesSendBuff_)
72 , valuesRecvBuff_(obv.valuesRecvBuff_)
73 , overlap_(obv.overlap_)
86 using ParentType::operator=;
94 ParentType::operator=(obv);
95 numIndicesSendBuff_ = obv.numIndicesSendBuff_;
96 indicesSendBuff_ = obv.indicesSendBuff_;
97 indicesRecvBuff_ = obv.indicesRecvBuff_;
98 valuesSendBuff_ = obv.valuesSendBuff_;
99 valuesRecvBuff_ = obv.valuesRecvBuff_;
100 overlap_ = obv.overlap_;
108 template <
class BlockVector>
111 size_t numDomestic = overlap_->numDomestic();
114 for (
unsigned domRowIdx = 0; domRowIdx < numDomestic; ++domRowIdx) {
115 Index nativeRowIdx = overlap_->domesticToNative(
static_cast<Index>(domRowIdx));
116 if (nativeRowIdx < 0)
117 (*this)[domRowIdx] = 0.0;
119 (*
this)[domRowIdx] = nativeBlockVector[nativeRowIdx];
132 template <
class NativeBlockVector>
133 void assign(
const NativeBlockVector& nativeBlockVector)
135 Index numDomestic = overlap_->numDomestic();
138 for (
Index domRowIdx = 0; domRowIdx < numDomestic; ++domRowIdx) {
139 Index nativeRowIdx = overlap_->domesticToNative(domRowIdx);
140 if (nativeRowIdx < 0)
141 (*this)[
static_cast<unsigned>(domRowIdx)] = 0.0;
143 (*
this)[
static_cast<unsigned>(domRowIdx)] = nativeBlockVector[
static_cast<unsigned>(nativeRowIdx)];
155 template <
class NativeBlockVector>
156 void assignTo(NativeBlockVector& nativeBlockVector)
const
159 size_t numNative = overlap_->numNative();
160 nativeBlockVector.resize(numNative);
161 for (
unsigned nativeRowIdx = 0; nativeRowIdx < numNative; ++nativeRowIdx) {
162 Index domRowIdx = overlap_->nativeToDomestic(
static_cast<Index>(nativeRowIdx));
165 nativeBlockVector[nativeRowIdx] = 0.0;
167 nativeBlockVector[nativeRowIdx] = (*this)[
static_cast<unsigned>(domRowIdx)];
178 for (
const auto peerRank: overlap_->peerSet())
179 sendEntries_(peerRank);
182 for (
const auto peerRank: overlap_->peerSet())
183 receiveFromMaster_(peerRank);
196 for (
const auto peerRank: overlap_->peerSet())
197 sendEntries_(peerRank);
200 for (
const auto peerRank: overlap_->peerSet())
201 receiveAdd_(peerRank);
209 for (
unsigned i = 0; i < this->size(); ++i) {
210 std::cout <<
"row " << i << (overlap_->isLocal(i) ?
" " :
"*")
211 <<
": " << (*
this)[i] <<
"\n" << std::flush;
216 void createBuffers_()
220 typename PeerSet::const_iterator peerIt;
221 typename PeerSet::const_iterator peerEndIt = overlap_->peerSet().end();
224 peerIt = overlap_->peerSet().begin();
225 for (; peerIt != peerEndIt; ++peerIt) {
228 size_t numEntries = overlap_->foreignOverlapSize(peerRank);
229 numIndicesSendBuff_[peerRank] = std::make_shared<MpiBuffer<unsigned> >(1);
230 indicesSendBuff_[peerRank] = std::make_shared<MpiBuffer<Index> >(numEntries);
231 valuesSendBuff_[peerRank] = std::make_shared<MpiBuffer<FieldVector> >(numEntries);
235 for (
unsigned i = 0; i < numEntries; ++i) {
236 Index domRowIdx = overlap_->foreignOverlapOffsetToDomesticIdx(peerRank, i);
237 indicesSendBuff[i] = overlap_->domesticToGlobal(domRowIdx);
241 (*numIndicesSendBuff_[peerRank])[0] =
static_cast<unsigned>(numEntries);
242 numIndicesSendBuff_[peerRank]->send(peerRank);
245 indicesSendBuff.send(peerRank);
249 peerIt = overlap_->peerSet().begin();
250 for (; peerIt != peerEndIt; ++peerIt) {
254 MpiBuffer<unsigned> numRowsRecvBuff(1);
255 numRowsRecvBuff.receive(peerRank);
256 unsigned numRows = numRowsRecvBuff[0];
259 indicesRecvBuff_[peerRank] = std::shared_ptr<MpiBuffer<Index> >(
260 new MpiBuffer<Index>(numRows));
261 valuesRecvBuff_[peerRank] = std::shared_ptr<MpiBuffer<FieldVector> >(
262 new MpiBuffer<FieldVector>(numRows));
263 MpiBuffer<Index>& indicesRecvBuff = *indicesRecvBuff_[peerRank];
266 indicesRecvBuff.receive(peerRank);
269 for (
unsigned i = 0; i != numRows; ++i) {
270 Index globalRowIdx = indicesRecvBuff[i];
271 Index domRowIdx = overlap_->globalToDomestic(globalRowIdx);
273 indicesRecvBuff[i] = domRowIdx;
278 peerIt = overlap_->peerSet().begin();
279 for (; peerIt != peerEndIt; ++peerIt) {
281 numIndicesSendBuff_[peerRank]->wait();
282 indicesSendBuff_[peerRank]->wait();
286 MpiBuffer<Index>& indicesSendBuff = *indicesSendBuff_[peerRank];
287 for (
unsigned i = 0; i < indicesSendBuff.size(); ++i) {
288 indicesSendBuff[i] = overlap_->globalToDomestic(indicesSendBuff[i]);
297 const MpiBuffer<Index>& indices = *indicesSendBuff_[peerRank];
298 MpiBuffer<FieldVector>& values = *valuesSendBuff_[peerRank];
299 for (
unsigned i = 0; i < indices.size(); ++i)
300 values[i] = (*
this)[
static_cast<unsigned>(indices[i])];
302 values.send(peerRank);
305 void waitSendFinished_()
307 typename PeerSet::const_iterator peerIt;
308 typename PeerSet::const_iterator peerEndIt = overlap_->peerSet().end();
311 peerIt = overlap_->peerSet().begin();
312 for (; peerIt != peerEndIt; ++peerIt) {
314 valuesSendBuff_[peerRank]->wait();
320 const MpiBuffer<Index>& indices = *indicesRecvBuff_[peerRank];
321 MpiBuffer<FieldVector>& values = *valuesRecvBuff_[peerRank];
324 values.receive(peerRank);
327 for (
unsigned j = 0; j < indices.size(); ++j) {
328 Index domRowIdx = indices[j];
329 if (overlap_->masterRank(domRowIdx) == peerRank) {
330 (*this)[
static_cast<unsigned>(domRowIdx)] = values[j];
337 const MpiBuffer<Index>& indices = *indicesRecvBuff_[peerRank];
338 MpiBuffer<FieldVector>& values = *valuesRecvBuff_[peerRank];
341 values.receive(peerRank);
344 for (
unsigned j = 0; j < indices.size(); ++j) {
345 Index domRowIdx = indices[j];
346 (*this)[
static_cast<unsigned>(domRowIdx)] += values[j];
350 std::map<ProcessRank, std::shared_ptr<MpiBuffer<unsigned> > > numIndicesSendBuff_;
351 std::map<ProcessRank, std::shared_ptr<MpiBuffer<Index> > > indicesSendBuff_;
352 std::map<ProcessRank, std::shared_ptr<MpiBuffer<Index> > > indicesRecvBuff_;
353 std::map<ProcessRank, std::shared_ptr<MpiBuffer<FieldVector> > > valuesSendBuff_;
354 std::map<ProcessRank, std::shared_ptr<MpiBuffer<FieldVector> > > valuesRecvBuff_;
356 const Overlap *overlap_;
An overlap aware block vector.
Definition: overlappingblockvector.hh:50
void syncAdd()
Syncronize all values of the block vector by adding up the values of all peer ranks.
Definition: overlappingblockvector.hh:193
OverlappingBlockVector & operator=(const OverlappingBlockVector &obv)
Assignment operator.
Definition: overlappingblockvector.hh:92
OverlappingBlockVector()
Default constructor.
Definition: overlappingblockvector.hh:79
OverlappingBlockVector(const Overlap &overlap)
Given a domestic overlap object, create an overlapping block vector coherent to it.
Definition: overlappingblockvector.hh:59
OverlappingBlockVector(const OverlappingBlockVector &obv)
Copy constructor.
Definition: overlappingblockvector.hh:66
void assign(const NativeBlockVector &nativeBlockVector)
Assign an overlapping block vector from a non-overlapping one, border entries are assigned using thei...
Definition: overlappingblockvector.hh:133
void print() const
Definition: overlappingblockvector.hh:207
void assignTo(NativeBlockVector &nativeBlockVector) const
Assign the local values to a non-overlapping block vector.
Definition: overlappingblockvector.hh:156
void assignAddBorder(const BlockVector &nativeBlockVector)
Assign an overlapping block vector from a non-overlapping one, border entries are added.
Definition: overlappingblockvector.hh:109
void sync()
Syncronize all values of the block vector from their master process.
Definition: overlappingblockvector.hh:175
unsigned ProcessRank
The type of the rank of a process.
Definition: overlaptypes.hh:49
int Index
The type of an index of a degree of freedom.
Definition: overlaptypes.hh:44
Definition: blackoilboundaryratevector.hh:37
This files provides several data structures for storing tuples of indices of remote and/or local proc...