overlappingblockvector.hh
Go to the documentation of this file.
1 // -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2 // vi: set et ts=4 sw=4 sts=4:
3 /*
4  Copyright (C) 2011-2013 by Andreas Lauser
5 
6  This file is part of the Open Porous Media project (OPM).
7 
8  OPM is free software: you can redistribute it and/or modify
9  it under the terms of the GNU General Public License as published by
10  the Free Software Foundation, either version 2 of the License, or
11  (at your option) any later version.
12 
13  OPM is distributed in the hope that it will be useful,
14  but WITHOUT ANY WARRANTY; without even the implied warranty of
15  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16  GNU General Public License for more details.
17 
18  You should have received a copy of the GNU General Public License
19  along with OPM. If not, see <http://www.gnu.org/licenses/>.
20 */
25 #ifndef EWOMS_OVERLAPPING_BLOCK_VECTOR_HH
26 #define EWOMS_OVERLAPPING_BLOCK_VECTOR_HH
27 
28 #include "overlaptypes.hh"
29 
31 #include <opm/material/common/Valgrind.hpp>
32 
33 #include <dune/istl/bvector.hh>
34 #include <dune/common/fvector.hh>
35 
36 #include <memory>
37 #include <map>
38 #include <iostream>
39 
40 namespace Ewoms {
41 namespace Linear {
42 
46 template <class FieldVector, class Overlap>
47 class OverlappingBlockVector : public Dune::BlockVector<FieldVector>
48 {
49  typedef Dune::BlockVector<FieldVector> ParentType;
50  typedef Dune::BlockVector<FieldVector> BlockVector;
51 
52 public:
57  OverlappingBlockVector(const Overlap &overlap)
58  : ParentType(overlap.numDomestic()), overlap_(&overlap)
59  { createBuffers_(); }
60 
65  : ParentType(obv)
66  , numIndicesSendBuff_(obv.numIndicesSendBuff_)
67  , indicesSendBuff_(obv.indicesSendBuff_)
68  , indicesRecvBuff_(obv.indicesRecvBuff_)
69  , valuesSendBuff_(obv.valuesSendBuff_)
70  , valuesRecvBuff_(obv.valuesRecvBuff_)
71  , overlap_(obv.overlap_)
72  {}
73 
78  {}
79 
81 
84  using ParentType::operator=;
86 
91  {
92  ParentType::operator=(obv);
93  numIndicesSendBuff_ = obv.numIndicesSendBuff_;
94  indicesSendBuff_ = obv.indicesSendBuff_;
95  indicesRecvBuff_ = obv.indicesRecvBuff_;
96  valuesSendBuff_ = obv.valuesSendBuff_;
97  valuesRecvBuff_ = obv.valuesRecvBuff_;
98  overlap_ = obv.overlap_;
99  return *this;
100  }
101 
106  void assignAddBorder(const BlockVector &nativeBlockVector)
107  {
108  int numDomestic = overlap_->numDomestic();
109 
110  // assign the local rows from the non-overlapping block vector
111  for (int domRowIdx = 0; domRowIdx < numDomestic; ++domRowIdx) {
112  int nativeRowIdx = overlap_->domesticToNative(domRowIdx);
113  if (nativeRowIdx < 0)
114  (*this)[domRowIdx] = 0.0;
115  else
116  (*this)[domRowIdx] = nativeBlockVector[nativeRowIdx];
117  }
118 
119  // add up the contents of border rows, for the remaining rows,
120  // get the values from their respective master process.
121  syncAddBorder();
122  }
123 
128  void assign(const BlockVector &nativeBlockVector)
129  {
130  int numDomestic = overlap_->numDomestic();
131 
132  // assign the local rows from the non-overlapping block vector
133  for (int domRowIdx = 0; domRowIdx < numDomestic; ++domRowIdx) {
134  int nativeRowIdx = overlap_->domesticToNative(domRowIdx);
135  if (nativeRowIdx < 0)
136  (*this)[domRowIdx] = 0.0;
137  else
138  (*this)[domRowIdx] = nativeBlockVector[nativeRowIdx];
139  }
140 
141  // add up the contents of border rows, for the remaining rows,
142  // get the values from their respective master process.
143  sync();
144  }
145 
150  void assignTo(BlockVector &nativeBlockVector) const
151  {
152  // assign the local rows
153  int numNative = overlap_->numNative();
154  nativeBlockVector.resize(numNative);
155  for (int nativeRowIdx = 0; nativeRowIdx < numNative; ++nativeRowIdx) {
156  int domRowIdx = overlap_->nativeToDomestic(nativeRowIdx);
157 
158  if (domRowIdx < 0)
159  nativeBlockVector[nativeRowIdx] = 0.0;
160  else
161  nativeBlockVector[nativeRowIdx] = (*this)[domRowIdx];
162  }
163  }
164 
169  void sync()
170  {
171  typename PeerSet::const_iterator peerIt;
172  typename PeerSet::const_iterator peerEndIt = overlap_->peerSet().end();
173 
174  // send all entries to all peers
175  peerIt = overlap_->peerSet().begin();
176  for (; peerIt != peerEndIt; ++peerIt) {
177  int peerRank = *peerIt;
178  sendEntries_(peerRank);
179  }
180 
181  // recieve all entries to the peers
182  peerIt = overlap_->peerSet().begin();
183  for (; peerIt != peerEndIt; ++peerIt) {
184  int peerRank = *peerIt;
185  receiveFromMaster_(peerRank);
186  }
187 
188  // wait until we have send everything
189  waitSendFinished_();
190  }
191 
196  void syncAdd()
197  {
198  typename PeerSet::const_iterator peerIt;
199  typename PeerSet::const_iterator peerEndIt = overlap_->peerSet().end();
200 
201  // send all entries to all peers
202  peerIt = overlap_->peerSet().begin();
203  for (; peerIt != peerEndIt; ++peerIt) {
204  int peerRank = *peerIt;
205  sendEntries_(peerRank);
206  }
207 
208  // recieve all entries to the peers
209  peerIt = overlap_->peerSet().begin();
210  for (; peerIt != peerEndIt; ++peerIt) {
211  int peerRank = *peerIt;
212  receiveAdd_(peerRank);
213  }
214 
215  // wait until we have send everything
216  waitSendFinished_();
217  }
218 
224  {
225  typename PeerSet::const_iterator peerIt;
226  typename PeerSet::const_iterator peerEndIt = overlap_->peerSet().end();
227 
228  // send all entries to all peers
229  peerIt = overlap_->peerSet().begin();
230  for (; peerIt != peerEndIt; ++peerIt) {
231  int peerRank = *peerIt;
232  sendEntries_(peerRank);
233  }
234 
235  // recieve all entries to the peers
236  peerIt = overlap_->peerSet().begin();
237  for (; peerIt != peerEndIt; ++peerIt) {
238  int peerRank = *peerIt;
239  receiveAddBorder_(peerRank);
240  }
241 
242  // wait until we have send everything
243  waitSendFinished_();
244 
245  }
246 
247  void print() const
248  {
249  for (int i = 0; i < this->size(); ++i) {
250  std::cout << "row " << i << (overlap_->isLocal(i) ? " " : "*")
251  << ": " << (*this)[i] << "\n" << std::flush;
252  }
253  }
254 
255 private:
256  void createBuffers_()
257  {
258 #if HAVE_MPI
259  // create array for the front indices
260  typename PeerSet::const_iterator peerIt;
261  typename PeerSet::const_iterator peerEndIt = overlap_->peerSet().end();
262 
263  // send all indices to the peers
264  peerIt = overlap_->peerSet().begin();
265  for (; peerIt != peerEndIt; ++peerIt) {
266  int peerRank = *peerIt;
267 
268  int numEntries = overlap_->foreignOverlapSize(peerRank);
269  numIndicesSendBuff_[peerRank] = std::make_shared<MpiBuffer<int> >(1);
270  indicesSendBuff_[peerRank] = std::make_shared<MpiBuffer<Index> >(numEntries);
271  valuesSendBuff_[peerRank] = std::make_shared<MpiBuffer<FieldVector> >(numEntries);
272 
273  // fill the indices buffer with global indices
274  MpiBuffer<Index> &indicesSendBuff = *indicesSendBuff_[peerRank];
275  for (int i = 0; i < numEntries; ++i) {
276  int domRowIdx = overlap_->foreignOverlapOffsetToDomesticIdx(peerRank, i);
277  indicesSendBuff[i] = overlap_->domesticToGlobal(domRowIdx);
278  }
279 
280  // first, send the number of indices
281  (*numIndicesSendBuff_[peerRank])[0] = numEntries;
282  numIndicesSendBuff_[peerRank]->send(peerRank);
283 
284  // then, send the indices themselfs
285  indicesSendBuff.send(peerRank);
286  }
287 
288  // receive the indices from the peers
289  peerIt = overlap_->peerSet().begin();
290  for (; peerIt != peerEndIt; ++peerIt) {
291  int peerRank = *peerIt;
292 
293  // receive size of overlap to peer
294  MpiBuffer<int> numRowsRecvBuff(1);
295  numRowsRecvBuff.receive(peerRank);
296  int numRows = numRowsRecvBuff[0];
297 
298  // then, create the MPI buffers
299  indicesRecvBuff_[peerRank] = std::shared_ptr<MpiBuffer<Index> >(
300  new MpiBuffer<Index>(numRows));
301  valuesRecvBuff_[peerRank] = std::shared_ptr<MpiBuffer<FieldVector> >(
302  new MpiBuffer<FieldVector>(numRows));
303  MpiBuffer<Index> &indicesRecvBuff = *indicesRecvBuff_[peerRank];
304 
305  // next, receive the actual indices
306  indicesRecvBuff.receive(peerRank);
307 
308  // finally, translate the global indices to domestic ones
309  for (int i = 0; i != numRows; ++i) {
310  int globalRowIdx = indicesRecvBuff[i];
311  int domRowIdx = overlap_->globalToDomestic(globalRowIdx);
312 
313  indicesRecvBuff[i] = domRowIdx;
314  }
315  }
316 
317  // wait for all send operations to complete
318  peerIt = overlap_->peerSet().begin();
319  for (; peerIt != peerEndIt; ++peerIt) {
320  int peerRank = *peerIt;
321  numIndicesSendBuff_[peerRank]->wait();
322  indicesSendBuff_[peerRank]->wait();
323 
324  // convert the global indices of the send buffer to
325  // domestic ones
326  MpiBuffer<Index> &indicesSendBuff = *indicesSendBuff_[peerRank];
327  for (unsigned i = 0; i < indicesSendBuff.size(); ++i) {
328  indicesSendBuff[i] = overlap_->globalToDomestic(indicesSendBuff[i]);
329  }
330  }
331 #endif // HAVE_MPI
332  }
333 
334  void sendEntries_(int peerRank)
335  {
336  // copy the values into the send buffer
337  const MpiBuffer<Index> &indices = *indicesSendBuff_[peerRank];
338  MpiBuffer<FieldVector> &values = *valuesSendBuff_[peerRank];
339  for (unsigned i = 0; i < indices.size(); ++i)
340  values[i] = (*this)[indices[i]];
341 
342  values.send(peerRank);
343  }
344 
345  void waitSendFinished_()
346  {
347  typename PeerSet::const_iterator peerIt;
348  typename PeerSet::const_iterator peerEndIt = overlap_->peerSet().end();
349 
350  // send all entries to all peers
351  peerIt = overlap_->peerSet().begin();
352  for (; peerIt != peerEndIt; ++peerIt) {
353  int peerRank = *peerIt;
354  valuesSendBuff_[peerRank]->wait();
355  }
356  }
357 
358  void receiveFromMaster_(int peerRank)
359  {
360  const MpiBuffer<Index> &indices = *indicesRecvBuff_[peerRank];
361  MpiBuffer<FieldVector> &values = *valuesRecvBuff_[peerRank];
362 
363  // receive the values from the peer
364  values.receive(peerRank);
365 
366  // copy them into the block vector
367  for (unsigned j = 0; j < indices.size(); ++j) {
368  Index domRowIdx = indices[j];
369  if (overlap_->masterRank(domRowIdx) == peerRank) {
370  (*this)[domRowIdx] = values[j];
371  }
372  }
373  }
374 
375  void receiveAddBorder_(int peerRank)
376  {
377  const MpiBuffer<Index> &indices = *indicesRecvBuff_[peerRank];
378  MpiBuffer<FieldVector> &values = *valuesRecvBuff_[peerRank];
379 
380  // receive the values from the peer
381  values.receive(peerRank);
382 
383  // add up the values of rows on the shared boundary
384  for (unsigned j = 0; j < indices.size(); ++j) {
385  int domRowIdx = indices[j];
386  if (overlap_->isBorderWith(domRowIdx, peerRank))
387  (*this)[domRowIdx] += values[j];
388  else
389  (*this)[domRowIdx] = values[j];
390  }
391  }
392 
393  void receiveAdd_(int peerRank)
394  {
395  const MpiBuffer<Index> &indices = *indicesRecvBuff_[peerRank];
396  MpiBuffer<FieldVector> &values = *valuesRecvBuff_[peerRank];
397 
398  // receive the values from the peer
399  values.receive(peerRank);
400 
401  // add up the values of rows on the shared boundary
402  for (int j = 0; j < indices.size(); ++j) {
403  int domRowIdx = indices[j];
404  (*this)[domRowIdx] += values[j];
405  }
406  }
407 
408  std::map<ProcessRank, std::shared_ptr<MpiBuffer<int> > > numIndicesSendBuff_;
409  std::map<ProcessRank, std::shared_ptr<MpiBuffer<Index> > > indicesSendBuff_;
410  std::map<ProcessRank, std::shared_ptr<MpiBuffer<Index> > > indicesRecvBuff_;
411  std::map<ProcessRank, std::shared_ptr<MpiBuffer<FieldVector> > > valuesSendBuff_;
412  std::map<ProcessRank, std::shared_ptr<MpiBuffer<FieldVector> > > valuesRecvBuff_;
413 
414  const Overlap *overlap_;
415 };
416 
417 } // namespace Linear
418 } // namespace Ewoms
419 
420 #endif
void assignTo(BlockVector &nativeBlockVector) const
Assign the local values to a non-overlapping block vector.
Definition: overlappingblockvector.hh:150
This files provides several data structures for storing tuples of indices of remote and/or local proc...
Simplifies handling of buffers to be used in conjunction with MPI.
OverlappingBlockVector(const OverlappingBlockVector &obv)
Copy constructor.
Definition: overlappingblockvector.hh:64
An overlap aware block vector.
Definition: overlappingblockvector.hh:47
void print() const
Definition: overlappingblockvector.hh:247
OverlappingBlockVector()
Default constructor.
Definition: overlappingblockvector.hh:77
void syncAdd()
Syncronize all values of the block vector by adding up the values of all peer ranks.
Definition: overlappingblockvector.hh:196
void assignAddBorder(const BlockVector &nativeBlockVector)
Assign an overlapping block vector from a non-overlapping one, border entries are added...
Definition: overlappingblockvector.hh:106
Definition: baseauxiliarymodule.hh:35
OverlappingBlockVector(const Overlap &overlap)
Given a domestic overlap object, create an overlapping block vector coherent to it.
Definition: overlappingblockvector.hh:57
OverlappingBlockVector & operator=(const OverlappingBlockVector &obv)
Assignment operator.
Definition: overlappingblockvector.hh:90
void sync()
Syncronize all values of the block vector from their master process.
Definition: overlappingblockvector.hh:169
int Index
The type of an index of a degree of freedom.
Definition: overlaptypes.hh:43
void syncAddBorder()
Syncronize all values of the block vector from the master rank, but add up the entries on the border...
Definition: overlappingblockvector.hh:223
void assign(const BlockVector &nativeBlockVector)
Assign an overlapping block vector from a non-overlapping one, border entries are assigned using thei...
Definition: overlappingblockvector.hh:128