CollectDataOnIORank_impl.hpp
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 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 2 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
23
24#ifndef OPM_COLLECT_DATA_ON_IO_RANK_IMPL_HPP
25#define OPM_COLLECT_DATA_ON_IO_RANK_IMPL_HPP
26
28
29#include <dune/common/version.hh>
30#include <dune/grid/common/gridenums.hh>
31#include <dune/grid/common/mcmgmapper.hh>
32#include <dune/grid/common/partitionset.hh>
33
34#include <opm/grid/common/CartesianIndexMapper.hpp>
35
36#include <algorithm>
37#include <cassert>
38#include <stdexcept>
39#include <string>
40#include <vector>
41
42namespace {
43 std::vector<std::string> toVector(const std::set<std::string>& string_set)
44 {
45 return { string_set.begin(), string_set.end() };
46 }
47}
48
49namespace Opm {
50
51// global id
53{
54 int globalId_;
55 int localIndex_;
56 bool isInterior_;
57
58public:
60 : globalId_(-1)
61 , localIndex_(-1)
62 , isInterior_(true)
63 {}
64 void setGhost()
65 { isInterior_ = false; }
66
67 void setId(int globalId)
68 { globalId_ = globalId; }
70 { localIndex_ = localIndex; }
71
72 int localIndex () const
73 { return localIndex_; }
74 int id () const
75 { return globalId_; }
76 bool isInterior() const
77 { return isInterior_; }
78};
79
80using IndexMapType = std::vector<int>;
81using IndexMapStorageType = std::vector<IndexMapType>;
82using P2PCommunicatorType = Dune::Point2PointCommunicator<Dune::SimpleMessageBuffer>;
84class DistributeIndexMapping : public P2PCommunicatorType::DataHandleInterface
85{
86protected:
87 const std::vector<int>& distributedGlobalIndex_;
90 std::map<int, int> globalPosition_;
91 std::set<int>& recv_;
92 std::vector<int>& ranks_;
93
94public:
95 DistributeIndexMapping(const std::vector<int>& globalIndex,
96 const std::vector<int>& distributedGlobalIndex,
97 IndexMapType& localIndexMap,
98 IndexMapStorageType& indexMaps,
99 std::vector<int>& ranks,
100 std::set<int>& recv,
101 bool isIORank)
102 : distributedGlobalIndex_(distributedGlobalIndex)
103 , localIndexMap_(localIndexMap)
104 , indexMaps_(indexMaps)
106 , recv_(recv)
107 , ranks_(ranks)
108 {
109 std::size_t size = globalIndex.size();
110 // create mapping globalIndex --> localIndex
111 if ( isIORank ) // ioRank
112 for (std::size_t index = 0; index < size; ++index)
113 globalPosition_.insert(std::make_pair(globalIndex[index], index));
114
115 // we need to create a mapping from local to global
116 if (!indexMaps_.empty()) {
117 if (isIORank)
118 ranks_.resize(size, -1);
119 // for the ioRank create a localIndex to index in global state map
120 IndexMapType& indexMap = indexMaps_.back();
121 std::size_t localSize = localIndexMap_.size();
122 indexMap.resize(localSize);
123 for (std::size_t i = 0; i < localSize; ++i)
124 {
126 indexMap[i] = id;
127 }
128 }
129 }
130
132 {
133 if (!indexMaps_.size())
134 return;
135
136 if ( ranks_.size() )
137 {
138 auto rankIt = recv_.begin();
139 // translate index maps from global cartesian to index
140 for (auto& indexMap: indexMaps_)
141 {
142 int rank = 0;
143 if (rankIt != recv_.end())
144 rank = *rankIt;
145
146 for (auto&& entry: indexMap)
148 auto candidate = globalPosition_.find(entry);
149 assert(candidate != globalPosition_.end());
150 entry = candidate->second;
151 // Using max should be backwards compatible
152 ranks_[entry] = std::max(ranks_[entry], rank);
153 }
154 if (rankIt != recv_.end())
155 ++rankIt;
156 }
157 #ifndef NDEBUG
158 for (const auto& rank: ranks_)
159 assert(rank>=0);
160 #endif
161 }
162 }
163
164 void pack(int link, MessageBufferType& buffer)
165 {
166 // we should only get one link
167 if (link != 0)
168 throw std::logic_error("link in method pack is not 0 as execpted");
169
170 // pack all interior global cell id's
171 int size = localIndexMap_.size();
172 buffer.write(size);
173
174 for (int index = 0; index < size; ++index) {
175 int globalIdx = distributedGlobalIndex_[localIndexMap_[index]];
176 buffer.write(globalIdx);
177 }
178 }
179
180 void unpack(int link, MessageBufferType& buffer)
181 {
182 // get index map for current link
183 IndexMapType& indexMap = indexMaps_[link];
184 assert(!globalPosition_.empty());
185
186 // unpack all interior global cell id's
187 int numCells = 0;
188 buffer.read(numCells);
189 indexMap.resize(numCells);
190 for (int index = 0; index < numCells; ++index) {
191 buffer.read(indexMap[index]);
192 }
193 }
194};
195
196
198template<class EquilMapper, class Mapper>
200{
201public:
202 ElementIndexScatterHandle(const EquilMapper& sendMapper, const Mapper& recvMapper, std::vector<int>& elementIndices)
203 : sendMapper_(sendMapper), recvMapper_(recvMapper), elementIndices_(elementIndices)
204 {}
205 using DataType = int;
206 bool fixedSize(int /*dim*/, int /*codim*/)
207 {
208 return true;
209 }
210
211 template<class T>
212 std::size_t size(const T&)
213 {
214 return 1;
215 }
216 template<class B, class T>
217 void gather(B& buffer, const T& t)
218 {
219 buffer.write(sendMapper_.index(t));
220 }
221 template<class B, class T>
222 void scatter(B& buffer, const T& t, std::size_t)
223 {
224 buffer.read(elementIndices_[recvMapper_.index(t)]);
225 }
226
227 bool contains(int dim, int codim)
228 {
229 return dim==3 && codim==0;
230 }
231private:
232 const EquilMapper& sendMapper_;
233 const Mapper& recvMapper_;
234 std::vector<int>& elementIndices_;
235};
236
238template<class Mapper>
240{
241public:
242 ElementIndexHandle(const Mapper& mapper, std::vector<int>& elementIndices)
243 : mapper_(mapper), elementIndices_(elementIndices)
244 {}
245 using DataType = int;
246 bool fixedSize(int /*dim*/, int /*codim*/)
247 {
248 return true;
249 }
250
251 template<class T>
252 std::size_t size(const T&)
253 {
254 return 1;
255 }
256 template<class B, class T>
257 void gather(B& buffer, const T& t)
258 {
259 buffer.write(elementIndices_[mapper_.index(t)]);
260 }
261 template<class B, class T>
262 void scatter(B& buffer, const T& t, std::size_t)
263 {
264 buffer.read(elementIndices_[mapper_.index(t)]);
265 }
266
267 bool contains(int dim, int codim)
268 {
269 return dim==3 && codim==0;
270 }
271private:
272 const Mapper& mapper_;
273 std::vector<int>& elementIndices_;
274};
275
276class PackUnPackCellData : public P2PCommunicatorType::DataHandleInterface
277{
278 const data::Solution& localCellData_;
279 data::Solution& globalCellData_;
280
281 const IndexMapType& localIndexMap_;
282 const IndexMapStorageType& indexMaps_;
283
284public:
285 PackUnPackCellData(const data::Solution& localCellData,
286 data::Solution& globalCellData,
287 const IndexMapType& localIndexMap,
288 const IndexMapStorageType& indexMaps,
289 std::size_t globalSize,
290 bool isIORank)
291 : localCellData_(localCellData)
292 , globalCellData_(globalCellData)
293 , localIndexMap_(localIndexMap)
294 , indexMaps_(indexMaps)
295 {
296 if (isIORank) {
297 // add missing data to global cell data
298 for (const auto& pair : localCellData_) {
299 const std::string& key = pair.first;
300 std::size_t containerSize = globalSize;
301 [[maybe_unused]] auto ret = globalCellData_.insert(key, pair.second.dim,
302 std::vector<double>(containerSize),
303 pair.second.target);
304 assert(ret.second);
305 }
306
307 MessageBufferType buffer;
308 pack(0, buffer);
309
310 // the last index map is the local one
311 doUnpack(indexMaps.back(), buffer);
312 }
313 }
314
315 // pack all data associated with link
316 void pack(int link, MessageBufferType& buffer)
317 {
318 // we should only get one link
319 if (link != 0)
320 throw std::logic_error("link in method pack is not 0 as expected");
321
322 // write all cell data registered in local state
323 for (const auto& pair : localCellData_) {
324 const auto& data = pair.second.data<double>();
325
326 // write all data from local data to buffer
327 write(buffer, localIndexMap_, data);
328 }
329 }
330
331 void doUnpack(const IndexMapType& indexMap, MessageBufferType& buffer)
332 {
333 // we loop over the data as
334 // its order governs the order the data got received.
335 for (auto& pair : localCellData_) {
336 const std::string& key = pair.first;
337 auto& data = globalCellData_.data<double>(key);
338
339 //write all data from local cell data to buffer
340 read(buffer, indexMap, data);
341 }
342 }
343
344 // unpack all data associated with link
345 void unpack(int link, MessageBufferType& buffer)
346 { doUnpack(indexMaps_[link], buffer); }
347
348protected:
349 template <class Vector>
351 const IndexMapType& localIndexMap,
352 const Vector& vector,
353 unsigned int offset = 0,
354 unsigned int stride = 1) const
355 {
356 unsigned int size = localIndexMap.size();
357 buffer.write(size);
358 assert(vector.size() >= stride * size);
359 for (unsigned int i=0; i<size; ++i)
360 {
361 unsigned int index = localIndexMap[i] * stride + offset;
362 assert(index < vector.size());
363 buffer.write(vector[index]);
364 }
365 }
366
367 template <class Vector>
369 const IndexMapType& indexMap,
370 Vector& vector,
371 unsigned int offset = 0,
372 unsigned int stride = 1) const
373 {
374 unsigned int size = 0;
375 buffer.read(size);
376 assert(size == indexMap.size());
377 for (unsigned int i=0; i<size; ++i) {
378 unsigned int index = indexMap[i] * stride + offset;
379 assert(index < vector.size());
380 buffer.read(vector[index]);
381 }
382 }
383};
384
385class PackUnPackWellData : public P2PCommunicatorType::DataHandleInterface
386{
387 const data::Wells& localWellData_;
388 data::Wells& globalWellData_;
389
390public:
391 PackUnPackWellData(const data::Wells& localWellData,
392 data::Wells& globalWellData,
393 bool isIORank)
394 : localWellData_(localWellData)
395 , globalWellData_(globalWellData)
396 {
397 if (isIORank) {
398 MessageBufferType buffer;
399 pack(0, buffer);
400
401 // pass a dummy_link to satisfy virtual class
402 int dummyLink = -1;
403 unpack(dummyLink, buffer);
404 }
405 }
406
407 // pack all data associated with link
408 void pack(int link, MessageBufferType& buffer)
409 {
410 // we should only get one link
411 if (link != 0)
412 throw std::logic_error("link in method pack is not 0 as expected");
413
414 localWellData_.write(buffer);
415 }
416
417 // unpack all data associated with link
418 void unpack(int /*link*/, MessageBufferType& buffer)
419 { globalWellData_.read(buffer); }
420
421};
422
423class PackUnPackGroupAndNetworkValues : public P2PCommunicatorType::DataHandleInterface
424{
425 const data::GroupAndNetworkValues& localGroupAndNetworkData_;
426 data::GroupAndNetworkValues& globalGroupAndNetworkData_;
427
428public:
429 PackUnPackGroupAndNetworkValues(const data::GroupAndNetworkValues& localGroupAndNetworkData,
430 data::GroupAndNetworkValues& globalGroupAndNetworkData,
431 const bool isIORank)
432 : localGroupAndNetworkData_ (localGroupAndNetworkData)
433 , globalGroupAndNetworkData_(globalGroupAndNetworkData)
434 {
435 if (! isIORank) { return; }
436
437 MessageBufferType buffer;
438 this->pack(0, buffer);
439
440 // pass a dummy_link to satisfy virtual class
441 const int dummyLink = -1;
442 this->unpack(dummyLink, buffer);
443 }
444
445 // pack all data associated with link
446 void pack(int link, MessageBufferType& buffer)
447 {
448 // we should only get one link
449 if (link != 0) {
450 throw std::logic_error {
451 "link in method pack is not 0 as expected"
452 };
453 }
454
455 // write all group and network (node/branch) data
456 this->localGroupAndNetworkData_.write(buffer);
457 }
458
459 // unpack all data associated with link
460 void unpack(int /*link*/, MessageBufferType& buffer)
461 { this->globalGroupAndNetworkData_.read(buffer); }
462};
463
464class PackUnPackBlockData : public P2PCommunicatorType::DataHandleInterface
465{
466 const std::map<std::pair<std::string, int>, double>& localBlockData_;
467 std::map<std::pair<std::string, int>, double>& globalBlockValues_;
468
469public:
470 PackUnPackBlockData(const std::map<std::pair<std::string, int>, double>& localBlockData,
471 std::map<std::pair<std::string, int>, double>& globalBlockValues,
472 bool isIORank)
473 : localBlockData_(localBlockData)
474 , globalBlockValues_(globalBlockValues)
475 {
476 if (isIORank) {
477 MessageBufferType buffer;
478 pack(0, buffer);
479
480 // pass a dummyLink to satisfy virtual class
481 int dummyLink = -1;
482 unpack(dummyLink, buffer);
483 }
484 }
485
486 // pack all data associated with link
487 void pack(int link, MessageBufferType& buffer)
488 {
489 // we should only get one link
490 if (link != 0)
491 throw std::logic_error("link in method pack is not 0 as expected");
492
493 // write all block data
494 unsigned int size = localBlockData_.size();
495 buffer.write(size);
496 for (const auto& map : localBlockData_) {
497 buffer.write(map.first.first);
498 buffer.write(map.first.second);
499 buffer.write(map.second);
500 }
501 }
502
503 // unpack all data associated with link
504 void unpack(int /*link*/, MessageBufferType& buffer)
505 {
506 // read all block data
507 unsigned int size = 0;
508 buffer.read(size);
509 for (std::size_t i = 0; i < size; ++i) {
510 std::string name;
511 int idx;
512 double data;
513 buffer.read(name);
514 buffer.read(idx);
515 buffer.read(data);
516 globalBlockValues_[std::make_pair(name, idx)] = data;
517 }
518 }
519};
520
521class PackUnPackWBPData : public P2PCommunicatorType::DataHandleInterface
522{
523 const data::WellBlockAveragePressures& localWBPData_;
524 data::WellBlockAveragePressures& globalWBPValues_;
525
526public:
527 PackUnPackWBPData(const data::WellBlockAveragePressures& localWBPData,
528 data::WellBlockAveragePressures& globalWBPValues,
529 const bool isIORank)
530 : localWBPData_ (localWBPData)
531 , globalWBPValues_(globalWBPValues)
532 {
533 if (! isIORank) {
534 return;
535 }
536
537 MessageBufferType buffer;
538 pack(0, buffer);
539
540 // Pass a dummy link to satisfy base class API requirement
541 const int dummyLink = -1;
542 unpack(dummyLink, buffer);
543 }
544
545 // Pack all data associated with link
546 void pack(int link, MessageBufferType& buffer)
547 {
548 // We should only get one link
549 if (link != 0) {
550 throw std::logic_error {
551 "link (" + std::to_string(link) +
552 ") in method pack() is not 0 as expected"
553 };
554 }
555
556 // Write all local, per-well, WBP data
557 this->localWBPData_.write(buffer);
558 }
559
560 // Unpack all data associated with link
561 void unpack([[maybe_unused]] const int link,
562 MessageBufferType& buffer)
563 {
564 this->globalWBPValues_.read(buffer);
565 }
566};
567
568class PackUnPackWellTestState : public P2PCommunicatorType::DataHandleInterface
569{
570public:
571 PackUnPackWellTestState(const WellTestState& localWTestState,
572 WellTestState& globalWTestState,
573 bool isIORank)
574 : local_(localWTestState)
575 , global_(globalWTestState)
576 {
577 if (isIORank) {
578 MessageBufferType buffer;
579 pack(0, buffer);
580
581 // pass a dummyLink to satisfy virtual class
582 int dummyLink = -1;
583 unpack(dummyLink, buffer);
584 }
585 }
586
587 void pack(int link, MessageBufferType& buffer) {
588 if (link != 0)
589 throw std::logic_error("link in method pack is not 0 as expected");
590 this->local_.pack(buffer);
591 }
592
593 void unpack(int, MessageBufferType& buffer) {
594 this->global_.unpack(buffer);
595 }
596
597private:
598 const WellTestState& local_;
599 WellTestState& global_;
600};
601
602class PackUnPackAquiferData : public P2PCommunicatorType::DataHandleInterface
603{
604 const data::Aquifers& localAquiferData_;
605 data::Aquifers& globalAquiferData_;
606
607public:
608 PackUnPackAquiferData(const data::Aquifers& localAquiferData,
609 data::Aquifers& globalAquiferData,
610 bool isIORank)
611 : localAquiferData_(localAquiferData)
612 , globalAquiferData_(globalAquiferData)
613 {
614 if (isIORank) {
615 MessageBufferType buffer;
616 pack(0, buffer);
617
618 // pass a dummyLink to satisfy virtual class
619 int dummyLink = -1;
620 unpack(dummyLink, buffer);
621 }
622 }
623
624 // pack all data associated with link
625 void pack(int link, MessageBufferType& buffer)
626 {
627 // we should only get one link
628 if (link != 0)
629 throw std::logic_error("link in method pack is not 0 as expected");
630
631 int size = localAquiferData_.size();
632 buffer.write(size);
633 for (const auto& [key, data] : localAquiferData_) {
634 buffer.write(key);
635 data.write(buffer);
636 }
637 }
638
639 // unpack all data associated with link
640 void unpack(int /*link*/, MessageBufferType& buffer)
641 {
642 int size;
643 buffer.read(size);
644 for (int i = 0; i < size; ++i) {
645 int key;
646 buffer.read(key);
647 data::AquiferData data;
648 data.read(buffer);
649
650 auto& aq = this->globalAquiferData_[key];
651
652 this->unpackCommon(data, aq);
653
654 if (auto const* aquFet = data.typeData.get<data::AquiferType::Fetkovich>();
655 aquFet != nullptr)
656 {
657 this->unpackAquFet(*aquFet, aq);
658 }
659 else if (auto const* aquCT = data.typeData.get<data::AquiferType::CarterTracy>();
660 aquCT != nullptr)
661 {
662 this->unpackCarterTracy(*aquCT, aq);
663 }
664 else if (auto const* aquNum = data.typeData.get<data::AquiferType::Numerical>();
665 aquNum != nullptr)
666 {
667 this->unpackNumericAquifer(*aquNum, aq);
668 }
669 }
670 }
671
672private:
673 void unpackCommon(const data::AquiferData& data, data::AquiferData& aq)
674 {
675 aq.aquiferID = std::max(aq.aquiferID, data.aquiferID);
676 aq.pressure = std::max(aq.pressure, data.pressure);
677 aq.initPressure = std::max(aq.initPressure, data.initPressure);
678 aq.datumDepth = std::max(aq.datumDepth, data.datumDepth);
679 aq.fluxRate += data.fluxRate;
680 aq.volume += data.volume;
681 }
682
683 void unpackAquFet(const data::FetkovichData& aquFet, data::AquiferData& aq)
684 {
685 if (! aq.typeData.is<data::AquiferType::Fetkovich>()) {
686 auto* tData = aq.typeData.create<data::AquiferType::Fetkovich>();
687 *tData = aquFet;
688 }
689 else {
690 const auto& src = aquFet;
691 auto& dst = *aq.typeData.getMutable<data::AquiferType::Fetkovich>();
692
693 dst.initVolume = std::max(dst.initVolume , src.initVolume);
694 dst.prodIndex = std::max(dst.prodIndex , src.prodIndex);
695 dst.timeConstant = std::max(dst.timeConstant, src.timeConstant);
696 }
697 }
698
699 void unpackCarterTracy(const data::CarterTracyData& aquCT, data::AquiferData& aq)
700 {
701 if (! aq.typeData.is<data::AquiferType::CarterTracy>()) {
702 auto* tData = aq.typeData.create<data::AquiferType::CarterTracy>();
703 *tData = aquCT;
704 }
705 else {
706 const auto& src = aquCT;
707 auto& dst = *aq.typeData.getMutable<data::AquiferType::CarterTracy>();
708
709 dst.timeConstant = std::max(dst.timeConstant , src.timeConstant);
710 dst.influxConstant = std::max(dst.influxConstant, src.influxConstant);
711 dst.waterDensity = std::max(dst.waterDensity , src.waterDensity);
712 dst.waterViscosity = std::max(dst.waterViscosity, src.waterViscosity);
713
714 dst.dimensionless_time = std::max(dst.dimensionless_time , src.dimensionless_time);
715 dst.dimensionless_pressure = std::max(dst.dimensionless_pressure, src.dimensionless_pressure);
716 }
717 }
718
719 void unpackNumericAquifer(const data::NumericAquiferData& aquNum, data::AquiferData& aq)
720 {
721 if (! aq.typeData.is<data::AquiferType::Numerical>()) {
722 auto* tData = aq.typeData.create<data::AquiferType::Numerical>();
723 *tData = aquNum;
724 }
725 else {
726 const auto& src = aquNum;
727 auto& dst = *aq.typeData.getMutable<data::AquiferType::Numerical>();
728
729 std::transform(src.initPressure.begin(),
730 src.initPressure.end(),
731 dst.initPressure.begin(),
732 dst.initPressure.begin(),
733 [](const double p0_1, const double p0_2)
734 {
735 return std::max(p0_1, p0_2);
736 });
737 }
738 }
739};
740
741class PackUnpackInterRegFlows : public P2PCommunicatorType::DataHandleInterface
742{
743 const InterRegFlowMap& localInterRegFlows_;
744 InterRegFlowMap& globalInterRegFlows_;
745
746public:
747 PackUnpackInterRegFlows(const InterRegFlowMap& localInterRegFlows,
748 InterRegFlowMap& globalInterRegFlows,
749 const bool isIORank)
750 : localInterRegFlows_(localInterRegFlows)
751 , globalInterRegFlows_(globalInterRegFlows)
752 {
753 if (! isIORank) { return; }
754
755 MessageBufferType buffer;
756 this->pack(0, buffer);
757
758 // pass a dummy_link to satisfy virtual class
759 const int dummyLink = -1;
760 this->unpack(dummyLink, buffer);
761 }
762
763 // pack all data associated with link
764 void pack(int link, MessageBufferType& buffer)
765 {
766 // we should only get one link
767 if (link != 0) {
768 throw std::logic_error {
769 "link in method pack is not 0 as expected"
770 };
771 }
772
773 // write all inter-region flow data
774 this->localInterRegFlows_.write(buffer);
775 }
776
777 // unpack all data associated with link
778 void unpack(int /*link*/, MessageBufferType& buffer)
779 { this->globalInterRegFlows_.read(buffer); }
780};
781
782class PackUnpackFlows : public P2PCommunicatorType::DataHandleInterface
783{
784 const std::array<FlowsData<double>, 3>& localFlows_;
785 std::array<FlowsData<double>, 3>& globalFlows_;
786
787public:
788 PackUnpackFlows(const std::array<FlowsData<double>, 3> & localFlows,
789 std::array<FlowsData<double>, 3>& globalFlows,
790 const bool isIORank)
791 : localFlows_(localFlows)
792 , globalFlows_(globalFlows)
793 {
794 if (! isIORank) { return; }
795
796 MessageBufferType buffer;
797 this->pack(0, buffer);
798
799 // pass a dummy_link to satisfy virtual class
800 const int dummyLink = -1;
801 this->unpack(dummyLink, buffer);
802 }
803
804 void pack(int link, MessageBufferType& buffer)
805 {
806 if (link != 0)
807 throw std::logic_error("link in method pack is not 0 as expected");
808 for (int i = 0; i < 3; ++i) {
809 const auto& name = localFlows_[i].name;
810 buffer.write(name);
811 const std::size_t size = localFlows_[i].indices.size();
812 buffer.write(size);
813 for (std::size_t j = 0; j < size; ++j) {
814 const auto& nncIdx = localFlows_[i].indices[j];
815 buffer.write(nncIdx);
816 const auto& flows = localFlows_[i].values[j];
817 buffer.write(flows);
818 }
819 }
820 }
821
822 void unpack(int /*link*/, MessageBufferType& buffer)
823 {
824 for (int i = 0; i < 3; ++i) {
825 std::string name;
826 buffer.read(name);
827 globalFlows_[i].name = name;
828 std::size_t size = 0;
829 buffer.read(size);
830 for (unsigned int j = 0; j < size; ++j) {
831 int nncIdx;
832 double data;
833 buffer.read(nncIdx);
834 buffer.read(data);
835 if (nncIdx < 0)
836 continue;
837 // This array is initialized in the collect(...) method below
838 globalFlows_[i].values[nncIdx] = data;
839 }
840 }
841 }
842};
843
844template <class Grid, class EquilGrid, class GridView>
846CollectDataOnIORank(const Grid& grid, const EquilGrid* equilGrid,
847 const GridView& localGridView,
848 const Dune::CartesianIndexMapper<Grid>& cartMapper,
849 const Dune::CartesianIndexMapper<EquilGrid>* equilCartMapper,
850 const std::set<std::string>& fipRegionsInterregFlow)
851 : toIORankComm_(grid.comm())
852 , globalInterRegFlows_(InterRegFlowMap::createMapFromNames(toVector(fipRegionsInterregFlow)))
853{
854 // index maps only have to be build when reordering is needed
855 if (!needsReordering && !isParallel())
856 return;
857
858 const CollectiveCommunication& comm = grid.comm();
859
860 {
861 std::set<int> send, recv;
862 using EquilGridView = typename EquilGrid::LeafGridView;
863 typename std::is_same<Grid, EquilGrid>::type isSameGrid;
864
865 typedef Dune::MultipleCodimMultipleGeomTypeMapper<GridView> ElementMapper;
866 ElementMapper elemMapper(localGridView, Dune::mcmgElementLayout());
867 sortedCartesianIdx_.reserve(localGridView.size(0));
868
869 for (const auto& elem : elements(localGridView, Dune::Partitions::interior))
870 {
871 auto idx = elemMapper.index(elem);
872 sortedCartesianIdx_.push_back(cartMapper.cartesianIndex(idx));
873 }
874
875 std::sort(sortedCartesianIdx_.begin(), sortedCartesianIdx_.end());
876 localIdxToGlobalIdx_.resize(localGridView.size(0), -1);
877
878 // the I/O rank receives from all other ranks
879 if (isIORank()) {
880 // We need a mapping from local to global grid, here we
881 // use equilGrid which represents a view on the global grid
882 // reserve memory
883 const std::size_t globalSize = equilGrid->leafGridView().size(0);
884 globalCartesianIndex_.resize(globalSize, -1);
885 const EquilGridView equilGridView = equilGrid->leafGridView();
886
887 using EquilElementMapper = Dune::MultipleCodimMultipleGeomTypeMapper<EquilGridView>;
888 EquilElementMapper equilElemMapper(equilGridView, Dune::mcmgElementLayout());
889
890 // Scatter the global index to local index for lookup during restart
891 if constexpr (isSameGrid) {
893 grid.scatterData(handle);
894 }
895
896 // loop over all elements (global grid) and store Cartesian index
897 for (const auto& elem : elements(equilGrid->leafGridView())) {
898 int elemIdx = equilElemMapper.index(elem);
899 int cartElemIdx = equilCartMapper->cartesianIndex(elemIdx);
900 globalCartesianIndex_[elemIdx] = cartElemIdx;
901 }
902
903 for (int i = 0; i < comm.size(); ++i) {
904 if (i != ioRank)
905 recv.insert(i);
906 }
907 }
908 else
909 {
910 // all other simply send to the I/O rank
911 send.insert(ioRank);
912
913 // Scatter the global index to local index for lookup during restart
914 // This is a bit hacky since the type differs from the iorank.
915 // But should work since we only receive, i.e. use the second parameter.
916 if constexpr (isSameGrid) {
918 grid.scatterData(handle);
919 }
920 }
921
922 // Sync the global element indices
923 if constexpr (isSameGrid) {
925 grid.communicate(handle, Dune::InteriorBorder_All_Interface,
926 Dune::ForwardCommunication);
927 }
928
929 localIndexMap_.clear();
930 const std::size_t gridSize = grid.size(0);
931 localIndexMap_.reserve(gridSize);
932
933 // store the local Cartesian index
934 IndexMapType distributedCartesianIndex;
935 distributedCartesianIndex.resize(gridSize, -1);
936
937 // A mapping for the whole grid (including the ghosts) is needed for restarts
938 for (const auto& elem : elements(localGridView, Dune::Partitions::interior)) {
939 int elemIdx = elemMapper.index(elem);
940 distributedCartesianIndex[elemIdx] = cartMapper.cartesianIndex(elemIdx);
941
942 // only store interior element for collection
943 assert(elem.partitionType() == Dune::InteriorEntity);
944
945 localIndexMap_.push_back(elemIdx);
946 }
947
948 // insert send and recv linkage to communicator
949 toIORankComm_.insertRequest(send, recv);
950
951 // need an index map for each rank
952 indexMaps_.clear();
953 indexMaps_.resize(comm.size());
954
955 // distribute global id's to io rank for later association of dof's
957 distributedCartesianIndex,
961 recv,
962 isIORank());
963 toIORankComm_.exchange(distIndexMapping);
964 }
965}
966
967template <class Grid, class EquilGrid, class GridView>
969collect(const data::Solution& localCellData,
970 const std::map<std::pair<std::string, int>, double>& localBlockData,
971 const data::Wells& localWellData,
972 const data::WellBlockAveragePressures& localWBPData,
973 const data::GroupAndNetworkValues& localGroupAndNetworkData,
974 const data::Aquifers& localAquiferData,
975 const WellTestState& localWellTestState,
976 const InterRegFlowMap& localInterRegFlows,
977 const std::array<FlowsData<double>, 3>& localFlowsn,
978 const std::array<FlowsData<double>, 3>& localFloresn)
979{
980 globalCellData_ = {};
981 globalBlockData_.clear();
982 globalWellData_.clear();
983 globalWBPData_.values.clear();
984 globalGroupAndNetworkData_.clear();
985 globalAquiferData_.clear();
986 globalWellTestState_.clear();
987 this->globalInterRegFlows_.clear();
988 globalFlowsn_ = {};
989 globalFloresn_ = {};
990
991 // index maps only have to be build when reordering is needed
992 if(!needsReordering && !isParallel())
993 return;
994
995 // this also linearises the local buffers on ioRank
996 PackUnPackCellData packUnpackCellData {
997 localCellData,
998 this->globalCellData_,
999 this->localIndexMap_,
1000 this->indexMaps_,
1001 this->numCells(),
1002 this->isIORank()
1003 };
1004
1005 if (! isParallel()) {
1006 // no need to collect anything.
1007 return;
1008 }
1009
1010 // Set the right sizes for Flowsn and Floresn
1011 for (int i = 0; i < 3; ++i) {
1012 const std::size_t sizeFlr = localFloresn[i].indices.size();
1013 globalFloresn_[i].resize(sizeFlr);
1014 const std::size_t sizeFlo = localFlowsn[i].indices.size();
1015 globalFlowsn_[i].resize(sizeFlo);
1016 }
1017
1018 PackUnPackWellData packUnpackWellData {
1019 localWellData,
1020 this->globalWellData_,
1021 this->isIORank()
1022 };
1023
1024 PackUnPackGroupAndNetworkValues packUnpackGroupAndNetworkData {
1025 localGroupAndNetworkData,
1026 this->globalGroupAndNetworkData_,
1027 this->isIORank()
1028 };
1029
1030 PackUnPackBlockData packUnpackBlockData {
1031 localBlockData,
1032 this->globalBlockData_,
1033 this->isIORank()
1034 };
1035
1036 PackUnPackWBPData packUnpackWBPData {
1037 localWBPData,
1038 this->globalWBPData_,
1039 this->isIORank()
1040 };
1041
1042 PackUnPackAquiferData packUnpackAquiferData {
1043 localAquiferData,
1044 this->globalAquiferData_,
1045 this->isIORank()
1046 };
1047
1048 PackUnPackWellTestState packUnpackWellTestState {
1049 localWellTestState,
1050 this->globalWellTestState_,
1051 this->isIORank()
1052 };
1053
1054 PackUnpackInterRegFlows packUnpackInterRegFlows {
1055 localInterRegFlows,
1056 this->globalInterRegFlows_,
1057 this->isIORank()
1058 };
1059
1060 PackUnpackFlows packUnpackFlowsn {
1061 localFlowsn,
1062 this->globalFlowsn_,
1063 this->isIORank()
1064 };
1065
1066 PackUnpackFlows packUnpackFloresn {
1067 localFloresn,
1068 this->globalFloresn_,
1069 this->isIORank()
1070 };
1071
1072 toIORankComm_.exchange(packUnpackCellData);
1073 toIORankComm_.exchange(packUnpackWellData);
1074 toIORankComm_.exchange(packUnpackGroupAndNetworkData);
1075 toIORankComm_.exchange(packUnpackBlockData);
1076 toIORankComm_.exchange(packUnpackWBPData);
1077 toIORankComm_.exchange(packUnpackAquiferData);
1078 toIORankComm_.exchange(packUnpackWellTestState);
1079 toIORankComm_.exchange(packUnpackInterRegFlows);
1080 toIORankComm_.exchange(packUnpackFlowsn);
1081 toIORankComm_.exchange(packUnpackFloresn);
1082
1083#ifndef NDEBUG
1084 // make sure every process is on the same page
1085 toIORankComm_.barrier();
1086#endif
1087}
1088
1089template <class Grid, class EquilGrid, class GridView>
1091localIdxToGlobalIdx(unsigned localIdx) const
1092{
1093 if (!isParallel()) {
1094 return localIdx;
1095 }
1096
1097 if (this->localIdxToGlobalIdx_.empty()) {
1098 throw std::logic_error("index map is not created on this rank");
1099 }
1100
1101 if (localIdx > this->localIdxToGlobalIdx_.size()) {
1102 throw std::logic_error("local index is outside map range");
1103 }
1104
1105 return this->localIdxToGlobalIdx_[localIdx];
1106}
1107
1108template <class Grid, class EquilGrid, class GridView>
1110isCartIdxOnThisRank(int cartIdx) const
1111{
1112 if (! this->isParallel()) {
1113 return true;
1114 }
1115
1116 assert (! needsReordering);
1117
1118 auto candidate = std::lower_bound(this->sortedCartesianIdx_.begin(),
1119 this->sortedCartesianIdx_.end(),
1120 cartIdx);
1121
1122 return (candidate != sortedCartesianIdx_.end())
1123 && (*candidate == cartIdx);
1124}
1125
1126} // end namespace Opm
1127
1128#endif // OPM_COLLECT_DATA_ON_IO_RANK_IMPL_HPP
Definition: CollectDataOnIORank.hpp:49
void collect(const data::Solution &localCellData, const std::map< std::pair< std::string, int >, double > &localBlockData, const data::Wells &localWellData, const data::WellBlockAveragePressures &localWBPData, const data::GroupAndNetworkValues &localGroupAndNetworkData, const data::Aquifers &localAquiferData, const WellTestState &localWellTestState, const InterRegFlowMap &interRegFlows, const std::array< FlowsData< double >, 3 > &localFlowsn, const std::array< FlowsData< double >, 3 > &localFloresn)
Definition: CollectDataOnIORank_impl.hpp:969
@ ioRank
Definition: CollectDataOnIORank.hpp:65
std::vector< int > localIdxToGlobalIdx_
Definition: CollectDataOnIORank.hpp:163
IndexMapStorageType indexMaps_
Definition: CollectDataOnIORank.hpp:154
std::vector< int > sortedCartesianIdx_
sorted list of cartesian indices present-
Definition: CollectDataOnIORank.hpp:169
static const bool needsReordering
Definition: CollectDataOnIORank.hpp:67
int localIdxToGlobalIdx(unsigned localIdx) const
Definition: CollectDataOnIORank_impl.hpp:1091
P2PCommunicatorType toIORankComm_
Definition: CollectDataOnIORank.hpp:150
bool isParallel() const
Definition: CollectDataOnIORank.hpp:128
bool isIORank() const
Definition: CollectDataOnIORank.hpp:125
typename Grid::CollectiveCommunication CollectiveCommunication
Definition: CollectDataOnIORank.hpp:58
bool isCartIdxOnThisRank(int cartIdx) const
Definition: CollectDataOnIORank_impl.hpp:1110
std::vector< int > globalRanks_
Definition: CollectDataOnIORank.hpp:155
CollectDataOnIORank(const Grid &grid, const EquilGrid *equilGrid, const GridView &gridView, const Dune::CartesianIndexMapper< Grid > &cartMapper, const Dune::CartesianIndexMapper< EquilGrid > *equilCartMapper, const std::set< std::string > &fipRegionsInterregFlow={})
Definition: CollectDataOnIORank_impl.hpp:846
IndexMapType localIndexMap_
Definition: CollectDataOnIORank.hpp:153
std::vector< int > IndexMapType
Definition: CollectDataOnIORank.hpp:60
IndexMapType globalCartesianIndex_
Definition: CollectDataOnIORank.hpp:152
Definition: CollectDataOnIORank_impl.hpp:85
void pack(int link, MessageBufferType &buffer)
Definition: CollectDataOnIORank_impl.hpp:164
std::map< int, int > globalPosition_
Definition: CollectDataOnIORank_impl.hpp:90
std::vector< int > & ranks_
Definition: CollectDataOnIORank_impl.hpp:92
IndexMapStorageType & indexMaps_
Definition: CollectDataOnIORank_impl.hpp:89
DistributeIndexMapping(const std::vector< int > &globalIndex, const std::vector< int > &distributedGlobalIndex, IndexMapType &localIndexMap, IndexMapStorageType &indexMaps, std::vector< int > &ranks, std::set< int > &recv, bool isIORank)
Definition: CollectDataOnIORank_impl.hpp:95
IndexMapType & localIndexMap_
Definition: CollectDataOnIORank_impl.hpp:88
const std::vector< int > & distributedGlobalIndex_
Definition: CollectDataOnIORank_impl.hpp:87
void unpack(int link, MessageBufferType &buffer)
Definition: CollectDataOnIORank_impl.hpp:180
~DistributeIndexMapping()
Definition: CollectDataOnIORank_impl.hpp:131
std::set< int > & recv_
Definition: CollectDataOnIORank_impl.hpp:91
Communication handle to scatter the global index.
Definition: CollectDataOnIORank_impl.hpp:240
std::size_t size(const T &)
Definition: CollectDataOnIORank_impl.hpp:252
bool fixedSize(int, int)
Definition: CollectDataOnIORank_impl.hpp:246
void scatter(B &buffer, const T &t, std::size_t)
Definition: CollectDataOnIORank_impl.hpp:262
void gather(B &buffer, const T &t)
Definition: CollectDataOnIORank_impl.hpp:257
int DataType
Definition: CollectDataOnIORank_impl.hpp:245
bool contains(int dim, int codim)
Definition: CollectDataOnIORank_impl.hpp:267
ElementIndexHandle(const Mapper &mapper, std::vector< int > &elementIndices)
Definition: CollectDataOnIORank_impl.hpp:242
Communication handle to scatter the global index.
Definition: CollectDataOnIORank_impl.hpp:200
std::size_t size(const T &)
Definition: CollectDataOnIORank_impl.hpp:212
bool contains(int dim, int codim)
Definition: CollectDataOnIORank_impl.hpp:227
void gather(B &buffer, const T &t)
Definition: CollectDataOnIORank_impl.hpp:217
ElementIndexScatterHandle(const EquilMapper &sendMapper, const Mapper &recvMapper, std::vector< int > &elementIndices)
Definition: CollectDataOnIORank_impl.hpp:202
void scatter(B &buffer, const T &t, std::size_t)
Definition: CollectDataOnIORank_impl.hpp:222
bool fixedSize(int, int)
Definition: CollectDataOnIORank_impl.hpp:206
int DataType
Definition: CollectDataOnIORank_impl.hpp:205
Definition: CollectDataOnIORank_impl.hpp:53
int localIndex() const
Definition: CollectDataOnIORank_impl.hpp:72
void setIndex(int localIndex)
Definition: CollectDataOnIORank_impl.hpp:69
int id() const
Definition: CollectDataOnIORank_impl.hpp:74
GlobalCellIndex()
Definition: CollectDataOnIORank_impl.hpp:59
void setGhost()
Definition: CollectDataOnIORank_impl.hpp:64
bool isInterior() const
Definition: CollectDataOnIORank_impl.hpp:76
void setId(int globalId)
Definition: CollectDataOnIORank_impl.hpp:67
Inter-region flow accumulation maps for all region definition arrays.
Definition: InterRegFlows.hpp:179
void read(MessageBufferType &buffer)
Definition: InterRegFlows.hpp:323
void write(MessageBufferType &buffer) const
Definition: InterRegFlows.hpp:296
Definition: CollectDataOnIORank_impl.hpp:603
void unpack(int, MessageBufferType &buffer)
Definition: CollectDataOnIORank_impl.hpp:640
void pack(int link, MessageBufferType &buffer)
Definition: CollectDataOnIORank_impl.hpp:625
PackUnPackAquiferData(const data::Aquifers &localAquiferData, data::Aquifers &globalAquiferData, bool isIORank)
Definition: CollectDataOnIORank_impl.hpp:608
Definition: CollectDataOnIORank_impl.hpp:465
void pack(int link, MessageBufferType &buffer)
Definition: CollectDataOnIORank_impl.hpp:487
PackUnPackBlockData(const std::map< std::pair< std::string, int >, double > &localBlockData, std::map< std::pair< std::string, int >, double > &globalBlockValues, bool isIORank)
Definition: CollectDataOnIORank_impl.hpp:470
void unpack(int, MessageBufferType &buffer)
Definition: CollectDataOnIORank_impl.hpp:504
Definition: CollectDataOnIORank_impl.hpp:277
void unpack(int link, MessageBufferType &buffer)
Definition: CollectDataOnIORank_impl.hpp:345
PackUnPackCellData(const data::Solution &localCellData, data::Solution &globalCellData, const IndexMapType &localIndexMap, const IndexMapStorageType &indexMaps, std::size_t globalSize, bool isIORank)
Definition: CollectDataOnIORank_impl.hpp:285
void pack(int link, MessageBufferType &buffer)
Definition: CollectDataOnIORank_impl.hpp:316
void read(MessageBufferType &buffer, const IndexMapType &indexMap, Vector &vector, unsigned int offset=0, unsigned int stride=1) const
Definition: CollectDataOnIORank_impl.hpp:368
void write(MessageBufferType &buffer, const IndexMapType &localIndexMap, const Vector &vector, unsigned int offset=0, unsigned int stride=1) const
Definition: CollectDataOnIORank_impl.hpp:350
void doUnpack(const IndexMapType &indexMap, MessageBufferType &buffer)
Definition: CollectDataOnIORank_impl.hpp:331
Definition: CollectDataOnIORank_impl.hpp:424
void unpack(int, MessageBufferType &buffer)
Definition: CollectDataOnIORank_impl.hpp:460
PackUnPackGroupAndNetworkValues(const data::GroupAndNetworkValues &localGroupAndNetworkData, data::GroupAndNetworkValues &globalGroupAndNetworkData, const bool isIORank)
Definition: CollectDataOnIORank_impl.hpp:429
void pack(int link, MessageBufferType &buffer)
Definition: CollectDataOnIORank_impl.hpp:446
Definition: CollectDataOnIORank_impl.hpp:522
void pack(int link, MessageBufferType &buffer)
Definition: CollectDataOnIORank_impl.hpp:546
PackUnPackWBPData(const data::WellBlockAveragePressures &localWBPData, data::WellBlockAveragePressures &globalWBPValues, const bool isIORank)
Definition: CollectDataOnIORank_impl.hpp:527
void unpack(const int link, MessageBufferType &buffer)
Definition: CollectDataOnIORank_impl.hpp:561
Definition: CollectDataOnIORank_impl.hpp:386
PackUnPackWellData(const data::Wells &localWellData, data::Wells &globalWellData, bool isIORank)
Definition: CollectDataOnIORank_impl.hpp:391
void pack(int link, MessageBufferType &buffer)
Definition: CollectDataOnIORank_impl.hpp:408
void unpack(int, MessageBufferType &buffer)
Definition: CollectDataOnIORank_impl.hpp:418
Definition: CollectDataOnIORank_impl.hpp:569
void pack(int link, MessageBufferType &buffer)
Definition: CollectDataOnIORank_impl.hpp:587
PackUnPackWellTestState(const WellTestState &localWTestState, WellTestState &globalWTestState, bool isIORank)
Definition: CollectDataOnIORank_impl.hpp:571
void unpack(int, MessageBufferType &buffer)
Definition: CollectDataOnIORank_impl.hpp:593
Definition: CollectDataOnIORank_impl.hpp:783
void pack(int link, MessageBufferType &buffer)
Definition: CollectDataOnIORank_impl.hpp:804
void unpack(int, MessageBufferType &buffer)
Definition: CollectDataOnIORank_impl.hpp:822
PackUnpackFlows(const std::array< FlowsData< double >, 3 > &localFlows, std::array< FlowsData< double >, 3 > &globalFlows, const bool isIORank)
Definition: CollectDataOnIORank_impl.hpp:788
Definition: CollectDataOnIORank_impl.hpp:742
void pack(int link, MessageBufferType &buffer)
Definition: CollectDataOnIORank_impl.hpp:764
void unpack(int, MessageBufferType &buffer)
Definition: CollectDataOnIORank_impl.hpp:778
PackUnpackInterRegFlows(const InterRegFlowMap &localInterRegFlows, InterRegFlowMap &globalInterRegFlows, const bool isIORank)
Definition: CollectDataOnIORank_impl.hpp:747
Definition: BlackoilPhases.hpp:27
std::vector< int > IndexMapType
Definition: CollectDataOnIORank_impl.hpp:80
Dune::Point2PointCommunicator< Dune::SimpleMessageBuffer > P2PCommunicatorType
Definition: CollectDataOnIORank_impl.hpp:82
std::vector< IndexMapType > IndexMapStorageType
Definition: CollectDataOnIORank_impl.hpp:81
std::string to_string(const ConvergenceReport::ReservoirFailure::Type t)
typename P2PCommunicatorType::MessageBufferType MessageBufferType
Definition: CollectDataOnIORank_impl.hpp:83