24#ifndef OPM_COLLECT_DATA_ON_IO_RANK_IMPL_HPP
25#define OPM_COLLECT_DATA_ON_IO_RANK_IMPL_HPP
27#include <opm/simulators/flow/CollectDataOnIORank.hpp>
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>
34#include <opm/grid/common/CartesianIndexMapper.hpp>
43 std::vector<std::string> toVector(
const std::set<std::string>& string_set)
45 return { string_set.begin(), string_set.end() };
65 { isInterior_ =
false; }
67 void setId(
int globalId)
68 { globalId_ = globalId; }
69 void setIndex(
int localIndex)
70 { localIndex_ = localIndex; }
72 int localIndex ()
const
73 {
return localIndex_; }
76 bool isInterior()
const
77 {
return isInterior_; }
80using IndexMapType = std::vector<int>;
81using IndexMapStorageType = std::vector<IndexMapType>;
82using P2PCommunicatorType = Dune::Point2PointCommunicator<Dune::SimpleMessageBuffer>;
83using MessageBufferType =
typename P2PCommunicatorType::MessageBufferType;
84class DistributeIndexMapping :
public P2PCommunicatorType::DataHandleInterface
87 const std::vector<int>& distributedGlobalIndex_;
88 IndexMapType& localIndexMap_;
89 IndexMapStorageType& indexMaps_;
90 std::map<int, int> globalPosition_;
92 std::vector<int>& ranks_;
95 DistributeIndexMapping(
const std::vector<int>& globalIndex,
96 const std::vector<int>& distributedGlobalIndex,
97 IndexMapType& localIndexMap,
98 IndexMapStorageType& indexMaps,
99 std::vector<int>& ranks,
102 : distributedGlobalIndex_(distributedGlobalIndex)
103 , localIndexMap_(localIndexMap)
104 , indexMaps_(indexMaps)
109 std::size_t size = globalIndex.size();
112 for (std::size_t index = 0; index < size; ++index)
113 globalPosition_.insert(std::make_pair(globalIndex[index], index));
116 if (!indexMaps_.empty()) {
118 ranks_.resize(size, -1);
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)
125 int id = distributedGlobalIndex_[localIndexMap_[i]];
131 ~DistributeIndexMapping()
133 if (!indexMaps_.size())
138 auto rankIt = recv_.begin();
140 for (
auto& indexMap: indexMaps_)
143 if (rankIt != recv_.end())
146 for (
auto&& entry: indexMap)
148 auto candidate = globalPosition_.find(entry);
149 assert(candidate != globalPosition_.end());
150 entry = candidate->second;
152 ranks_[entry] = std::max(ranks_[entry], rank);
154 if (rankIt != recv_.end())
158 for (
const auto& rank: ranks_)
164 void pack(
int link, MessageBufferType& buffer)
168 throw std::logic_error(
"link in method pack is not 0 as execpted");
171 int size = localIndexMap_.size();
174 for (
int index = 0; index < size; ++index) {
175 int globalIdx = distributedGlobalIndex_[localIndexMap_[index]];
176 buffer.write(globalIdx);
180 void unpack(
int link, MessageBufferType& buffer)
183 IndexMapType& indexMap = indexMaps_[link];
184 assert(!globalPosition_.empty());
188 buffer.read(numCells);
189 indexMap.resize(numCells);
190 for (
int index = 0; index < numCells; ++index) {
191 buffer.read(indexMap[index]);
198template<
class EquilMapper,
class Mapper>
199class ElementIndexScatterHandle
202 ElementIndexScatterHandle(
const EquilMapper& sendMapper,
const Mapper& recvMapper, std::vector<int>& elementIndices)
203 : sendMapper_(sendMapper), recvMapper_(recvMapper), elementIndices_(elementIndices)
205 using DataType = int;
206 bool fixedSize(
int ,
int )
212 std::size_t size(
const T&)
216 template<
class B,
class T>
217 void gather(B& buffer,
const T& t)
219 buffer.write(sendMapper_.index(t));
221 template<
class B,
class T>
222 void scatter(B& buffer,
const T& t, std::size_t)
224 buffer.read(elementIndices_[recvMapper_.index(t)]);
227 bool contains(
int dim,
int codim)
229 return dim==3 && codim==0;
232 const EquilMapper& sendMapper_;
233 const Mapper& recvMapper_;
234 std::vector<int>& elementIndices_;
238template<
class Mapper>
239class ElementIndexHandle
242 ElementIndexHandle(
const Mapper& mapper, std::vector<int>& elementIndices)
243 : mapper_(mapper), elementIndices_(elementIndices)
245 using DataType = int;
246 bool fixedSize(
int ,
int )
252 std::size_t size(
const T&)
256 template<
class B,
class T>
257 void gather(B& buffer,
const T& t)
259 buffer.write(elementIndices_[mapper_.index(t)]);
261 template<
class B,
class T>
262 void scatter(B& buffer,
const T& t, std::size_t)
264 buffer.read(elementIndices_[mapper_.index(t)]);
267 bool contains(
int dim,
int codim)
269 return dim==3 && codim==0;
272 const Mapper& mapper_;
273 std::vector<int>& elementIndices_;
276class PackUnPackCellData :
public P2PCommunicatorType::DataHandleInterface
278 const data::Solution& localCellData_;
279 data::Solution& globalCellData_;
281 const IndexMapType& localIndexMap_;
282 const IndexMapStorageType& indexMaps_;
285 PackUnPackCellData(
const data::Solution& localCellData,
286 data::Solution& globalCellData,
287 const IndexMapType& localIndexMap,
288 const IndexMapStorageType& indexMaps,
289 std::size_t globalSize,
291 : localCellData_(localCellData)
292 , globalCellData_(globalCellData)
293 , localIndexMap_(localIndexMap)
294 , indexMaps_(indexMaps)
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),
307 MessageBufferType buffer;
311 doUnpack(indexMaps.back(), buffer);
316 void pack(
int link, MessageBufferType& buffer)
320 throw std::logic_error(
"link in method pack is not 0 as expected");
323 for (
const auto& pair : localCellData_) {
324 const auto& data = pair.second.data<
double>();
327 write(buffer, localIndexMap_, data);
331 void doUnpack(
const IndexMapType& indexMap, MessageBufferType& buffer)
335 for (
auto& pair : localCellData_) {
336 const std::string& key = pair.first;
337 auto& data = globalCellData_.data<
double>(key);
340 read(buffer, indexMap, data);
345 void unpack(
int link, MessageBufferType& buffer)
346 { doUnpack(indexMaps_[link], buffer); }
349 template <
class Vector>
350 void write(MessageBufferType& buffer,
351 const IndexMapType& localIndexMap,
352 const Vector& vector,
353 unsigned int offset = 0,
354 unsigned int stride = 1)
const
356 unsigned int size = localIndexMap.size();
358 assert(vector.size() >= stride * size);
359 for (
unsigned int i=0; i<size; ++i)
361 unsigned int index = localIndexMap[i] * stride + offset;
362 assert(index < vector.size());
363 buffer.write(vector[index]);
367 template <
class Vector>
368 void read(MessageBufferType& buffer,
369 const IndexMapType& indexMap,
371 unsigned int offset = 0,
372 unsigned int stride = 1)
const
374 unsigned int size = 0;
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]);
385class PackUnPackWellData :
public P2PCommunicatorType::DataHandleInterface
387 const data::Wells& localWellData_;
388 data::Wells& globalWellData_;
391 PackUnPackWellData(
const data::Wells& localWellData,
392 data::Wells& globalWellData,
394 : localWellData_(localWellData)
395 , globalWellData_(globalWellData)
398 MessageBufferType buffer;
403 unpack(dummyLink, buffer);
408 void pack(
int link, MessageBufferType& buffer)
412 throw std::logic_error(
"link in method pack is not 0 as expected");
414 localWellData_.write(buffer);
418 void unpack(
int , MessageBufferType& buffer)
419 { globalWellData_.read(buffer); }
423class PackUnPackGroupAndNetworkValues :
public P2PCommunicatorType::DataHandleInterface
425 const data::GroupAndNetworkValues& localGroupAndNetworkData_;
426 data::GroupAndNetworkValues& globalGroupAndNetworkData_;
429 PackUnPackGroupAndNetworkValues(
const data::GroupAndNetworkValues& localGroupAndNetworkData,
430 data::GroupAndNetworkValues& globalGroupAndNetworkData,
432 : localGroupAndNetworkData_ (localGroupAndNetworkData)
433 , globalGroupAndNetworkData_(globalGroupAndNetworkData)
435 if (! isIORank) {
return; }
437 MessageBufferType buffer;
438 this->pack(0, buffer);
441 const int dummyLink = -1;
442 this->unpack(dummyLink, buffer);
446 void pack(
int link, MessageBufferType& buffer)
450 throw std::logic_error {
451 "link in method pack is not 0 as expected"
456 this->localGroupAndNetworkData_.write(buffer);
460 void unpack(
int , MessageBufferType& buffer)
461 { this->globalGroupAndNetworkData_.read(buffer); }
464class PackUnPackBlockData :
public P2PCommunicatorType::DataHandleInterface
466 const std::map<std::pair<std::string, int>,
double>& localBlockData_;
467 std::map<std::pair<std::string, int>,
double>& globalBlockValues_;
470 PackUnPackBlockData(
const std::map<std::pair<std::string, int>,
double>& localBlockData,
471 std::map<std::pair<std::string, int>,
double>& globalBlockValues,
473 : localBlockData_(localBlockData)
474 , globalBlockValues_(globalBlockValues)
477 MessageBufferType buffer;
482 unpack(dummyLink, buffer);
487 void pack(
int link, MessageBufferType& buffer)
491 throw std::logic_error(
"link in method pack is not 0 as expected");
494 unsigned int size = localBlockData_.size();
496 for (
const auto& map : localBlockData_) {
497 buffer.write(map.first.first);
498 buffer.write(map.first.second);
499 buffer.write(map.second);
504 void unpack(
int , MessageBufferType& buffer)
507 unsigned int size = 0;
509 for (std::size_t i = 0; i < size; ++i) {
516 globalBlockValues_[std::make_pair(name, idx)] = data;
521class PackUnPackWBPData :
public P2PCommunicatorType::DataHandleInterface
523 const data::WellBlockAveragePressures& localWBPData_;
524 data::WellBlockAveragePressures& globalWBPValues_;
527 PackUnPackWBPData(
const data::WellBlockAveragePressures& localWBPData,
528 data::WellBlockAveragePressures& globalWBPValues,
530 : localWBPData_ (localWBPData)
531 , globalWBPValues_(globalWBPValues)
537 MessageBufferType buffer;
541 const int dummyLink = -1;
542 unpack(dummyLink, buffer);
546 void pack(
int link, MessageBufferType& buffer)
550 throw std::logic_error {
551 "link (" + std::to_string(link) +
552 ") in method pack() is not 0 as expected"
557 this->localWBPData_.write(buffer);
561 void unpack([[maybe_unused]]
const int link,
562 MessageBufferType& buffer)
564 this->globalWBPValues_.read(buffer);
568class PackUnPackWellTestState :
public P2PCommunicatorType::DataHandleInterface
571 PackUnPackWellTestState(
const WellTestState& localWTestState,
572 WellTestState& globalWTestState,
574 : local_(localWTestState)
575 , global_(globalWTestState)
578 MessageBufferType buffer;
583 unpack(dummyLink, buffer);
587 void pack(
int link, MessageBufferType& buffer) {
589 throw std::logic_error(
"link in method pack is not 0 as expected");
590 this->local_.pack(buffer);
593 void unpack(
int, MessageBufferType& buffer) {
594 this->global_.unpack(buffer);
598 const WellTestState& local_;
599 WellTestState& global_;
602class PackUnPackAquiferData :
public P2PCommunicatorType::DataHandleInterface
604 const data::Aquifers& localAquiferData_;
605 data::Aquifers& globalAquiferData_;
608 PackUnPackAquiferData(
const data::Aquifers& localAquiferData,
609 data::Aquifers& globalAquiferData,
611 : localAquiferData_(localAquiferData)
612 , globalAquiferData_(globalAquiferData)
615 MessageBufferType buffer;
620 unpack(dummyLink, buffer);
625 void pack(
int link, MessageBufferType& buffer)
629 throw std::logic_error(
"link in method pack is not 0 as expected");
631 int size = localAquiferData_.size();
633 for (
const auto& [key, data] : localAquiferData_) {
640 void unpack(
int , MessageBufferType& buffer)
644 for (
int i = 0; i < size; ++i) {
647 data::AquiferData data;
650 auto& aq = this->globalAquiferData_[key];
652 this->unpackCommon(data, aq);
654 if (
auto const* aquFet = data.typeData.get<data::AquiferType::Fetkovich>();
657 this->unpackAquFet(*aquFet, aq);
659 else if (
auto const* aquCT = data.typeData.get<data::AquiferType::CarterTracy>();
662 this->unpackCarterTracy(*aquCT, aq);
664 else if (
auto const* aquNum = data.typeData.get<data::AquiferType::Numerical>();
667 this->unpackNumericAquifer(*aquNum, aq);
673 void unpackCommon(
const data::AquiferData& data, data::AquiferData& aq)
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;
683 void unpackAquFet(
const data::FetkovichData& aquFet, data::AquiferData& aq)
685 if (! aq.typeData.is<data::AquiferType::Fetkovich>()) {
686 auto* tData = aq.typeData.create<data::AquiferType::Fetkovich>();
690 const auto& src = aquFet;
691 auto& dst = *aq.typeData.getMutable<data::AquiferType::Fetkovich>();
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);
699 void unpackCarterTracy(
const data::CarterTracyData& aquCT, data::AquiferData& aq)
701 if (! aq.typeData.is<data::AquiferType::CarterTracy>()) {
702 auto* tData = aq.typeData.create<data::AquiferType::CarterTracy>();
706 const auto& src = aquCT;
707 auto& dst = *aq.typeData.getMutable<data::AquiferType::CarterTracy>();
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);
714 dst.dimensionless_time = std::max(dst.dimensionless_time , src.dimensionless_time);
715 dst.dimensionless_pressure = std::max(dst.dimensionless_pressure, src.dimensionless_pressure);
719 void unpackNumericAquifer(
const data::NumericAquiferData& aquNum, data::AquiferData& aq)
721 if (! aq.typeData.is<data::AquiferType::Numerical>()) {
722 auto* tData = aq.typeData.create<data::AquiferType::Numerical>();
726 const auto& src = aquNum;
727 auto& dst = *aq.typeData.getMutable<data::AquiferType::Numerical>();
729 std::ranges::transform(src.initPressure, dst.initPressure, dst.initPressure.begin(),
730 [](
const double p0_1,
const double p0_2)
731 { return std::max(p0_1, p0_2); });
736class PackUnpackInterRegFlows :
public P2PCommunicatorType::DataHandleInterface
745 : localInterRegFlows_(localInterRegFlows)
746 , globalInterRegFlows_(globalInterRegFlows)
748 if (! isIORank) {
return; }
750 MessageBufferType buffer;
751 this->pack(0, buffer);
754 const int dummyLink = -1;
755 this->unpack(dummyLink, buffer);
759 void pack(
int link, MessageBufferType& buffer)
763 throw std::logic_error {
764 "link in method pack is not 0 as expected"
769 this->localInterRegFlows_.
write(buffer);
773 void unpack(
int , MessageBufferType& buffer)
774 { this->globalInterRegFlows_.
read(buffer); }
777class PackUnpackFlows :
public P2PCommunicatorType::DataHandleInterface
779 const std::array<FlowsData<double>, 3>& localFlows_;
780 std::array<FlowsData<double>, 3>& globalFlows_;
786 : localFlows_(localFlows)
787 , globalFlows_(globalFlows)
789 if (! isIORank) {
return; }
791 MessageBufferType buffer;
792 this->pack(0, buffer);
795 const int dummyLink = -1;
796 this->unpack(dummyLink, buffer);
799 void pack(
int link, MessageBufferType& buffer)
802 throw std::logic_error(
"link in method pack is not 0 as expected");
803 for (
int i = 0; i < 3; ++i) {
804 const auto& name = localFlows_[i].name;
806 const std::size_t size = localFlows_[i].indices.size();
808 for (std::size_t j = 0; j < size; ++j) {
809 const auto& nncIdx = localFlows_[i].indices[j];
810 buffer.write(nncIdx);
811 const auto& flows = localFlows_[i].values[j];
817 void unpack(
int , MessageBufferType& buffer)
819 for (
int i = 0; i < 3; ++i) {
822 globalFlows_[i].name = name;
823 std::size_t size = 0;
825 for (
unsigned int j = 0; j < size; ++j) {
833 globalFlows_[i].values[nncIdx] = data;
839template <
class Gr
id,
class EquilGr
id,
class Gr
idView>
842 const GridView& localGridView,
845 const std::set<std::string>& fipRegionsInterregFlow)
846 : toIORankComm_(grid.comm())
847 , globalInterRegFlows_(
InterRegFlowMap::createMapFromNames(toVector(fipRegionsInterregFlow)))
850 if ((!needsReordering && !isParallel()) || (isParallel() && (grid.maxLevel()>0)))
853 const CollectiveCommunication& comm = grid.comm();
856 std::set<int> send, recv;
857 using EquilGridView = typename EquilGrid::LeafGridView;
858 typename std::is_same<Grid, EquilGrid>::type isSameGrid;
860 typedef Dune::MultipleCodimMultipleGeomTypeMapper<GridView> ElementMapper;
861 ElementMapper elemMapper(localGridView, Dune::mcmgElementLayout());
862 sortedCartesianIdx_.reserve(localGridView.size(0));
864 for (const auto& elem : elements(localGridView, Dune::Partitions::interior))
866 auto idx = elemMapper.index(elem);
867 sortedCartesianIdx_.push_back(cartMapper.cartesianIndex(idx));
870 std::ranges::sort(sortedCartesianIdx_);
871 localIdxToGlobalIdx_.resize(localGridView.size(0), -1);
878 const std::size_t globalSize = equilGrid->leafGridView().size(0);
879 globalCartesianIndex_.resize(globalSize, -1);
880 const EquilGridView equilGridView = equilGrid->leafGridView();
882 using EquilElementMapper = Dune::MultipleCodimMultipleGeomTypeMapper<EquilGridView>;
883 EquilElementMapper equilElemMapper(equilGridView, Dune::mcmgElementLayout());
886 if constexpr (isSameGrid) {
888 grid.scatterData(handle);
892 for (
const auto& elem : elements(equilGrid->leafGridView())) {
893 int elemIdx = equilElemMapper.index(elem);
894 int cartElemIdx = equilCartMapper->cartesianIndex(elemIdx);
895 globalCartesianIndex_[elemIdx] = cartElemIdx;
898 for (
int i = 0; i < comm.size(); ++i) {
911 if constexpr (isSameGrid) {
912 ElementIndexScatterHandle<ElementMapper, ElementMapper> handle(elemMapper, elemMapper, localIdxToGlobalIdx_);
913 grid.scatterData(handle);
918 if constexpr (isSameGrid) {
919 ElementIndexHandle<ElementMapper> handle(elemMapper, localIdxToGlobalIdx_);
920 grid.communicate(handle, Dune::InteriorBorder_All_Interface,
921 Dune::ForwardCommunication);
924 localIndexMap_.clear();
925 const std::size_t gridSize = grid.size(0);
926 localIndexMap_.reserve(gridSize);
929 IndexMapType distributedCartesianIndex;
930 distributedCartesianIndex.resize(gridSize, -1);
933 for (
const auto& elem : elements(localGridView, Dune::Partitions::interior)) {
934 int elemIdx = elemMapper.index(elem);
935 distributedCartesianIndex[elemIdx] = cartMapper.cartesianIndex(elemIdx);
938 assert(elem.partitionType() == Dune::InteriorEntity);
940 localIndexMap_.push_back(elemIdx);
944 toIORankComm_.insertRequest(send, recv);
948 indexMaps_.resize(comm.size());
951 DistributeIndexMapping distIndexMapping(globalCartesianIndex_,
952 distributedCartesianIndex,
958 toIORankComm_.exchange(distIndexMapping);
962template <
class Gr
id,
class EquilGr
id,
class Gr
idView>
963void CollectDataOnIORank<Grid,EquilGrid,GridView>::
964collect(
const data::Solution& localCellData,
965 const std::map<std::pair<std::string, int>,
double>& localBlockData,
966 std::map<std::pair<std::string, int>,
double>& extraBlockData,
967 const data::Wells& localWellData,
968 const data::WellBlockAveragePressures& localWBPData,
969 const data::GroupAndNetworkValues& localGroupAndNetworkData,
970 const data::Aquifers& localAquiferData,
971 const WellTestState& localWellTestState,
972 const InterRegFlowMap& localInterRegFlows,
973 const std::array<FlowsData<double>, 3>& localFlowsn,
974 const std::array<FlowsData<double>, 3>& localFloresn)
976 globalCellData_ = {};
977 globalBlockData_.clear();
978 std::map<std::pair<std::string,int>,
double> globalExtraBlockData;
979 globalWellData_.clear();
980 globalWBPData_.values.clear();
981 globalGroupAndNetworkData_.clear();
982 globalAquiferData_.clear();
983 globalWellTestState_.clear();
984 this->globalInterRegFlows_.clear();
989 if(!needsReordering && !isParallel())
993 PackUnPackCellData packUnpackCellData {
995 this->globalCellData_,
996 this->localIndexMap_,
1002 if (! isParallel()) {
1008 for (
int i = 0; i < 3; ++i) {
1009 const std::size_t sizeFlr = localFloresn[i].indices.size();
1010 globalFloresn_[i].resize(sizeFlr);
1011 const std::size_t sizeFlo = localFlowsn[i].indices.size();
1012 globalFlowsn_[i].resize(sizeFlo);
1015 PackUnPackWellData packUnpackWellData {
1017 this->globalWellData_,
1021 PackUnPackGroupAndNetworkValues packUnpackGroupAndNetworkData {
1022 localGroupAndNetworkData,
1023 this->globalGroupAndNetworkData_,
1027 PackUnPackBlockData packUnpackBlockData {
1029 this->globalBlockData_,
1033 PackUnPackBlockData packUnpackExtraBlockData {
1035 globalExtraBlockData,
1039 PackUnPackWBPData packUnpackWBPData {
1041 this->globalWBPData_,
1045 PackUnPackAquiferData packUnpackAquiferData {
1047 this->globalAquiferData_,
1051 PackUnPackWellTestState packUnpackWellTestState {
1053 this->globalWellTestState_,
1057 PackUnpackInterRegFlows packUnpackInterRegFlows {
1059 this->globalInterRegFlows_,
1063 PackUnpackFlows packUnpackFlowsn {
1065 this->globalFlowsn_,
1069 PackUnpackFlows packUnpackFloresn {
1071 this->globalFloresn_,
1075 toIORankComm_.exchange(packUnpackCellData);
1076 toIORankComm_.exchange(packUnpackWellData);
1077 toIORankComm_.exchange(packUnpackGroupAndNetworkData);
1078 toIORankComm_.exchange(packUnpackBlockData);
1079 toIORankComm_.exchange(packUnpackExtraBlockData);
1080 toIORankComm_.exchange(packUnpackWBPData);
1081 toIORankComm_.exchange(packUnpackAquiferData);
1082 toIORankComm_.exchange(packUnpackWellTestState);
1083 toIORankComm_.exchange(packUnpackInterRegFlows);
1084 toIORankComm_.exchange(packUnpackFlowsn);
1085 toIORankComm_.exchange(packUnpackFloresn);
1089 toIORankComm_.barrier();
1092 extraBlockData = globalExtraBlockData;
1095template <
class Gr
id,
class EquilGr
id,
class Gr
idView>
1096int CollectDataOnIORank<Grid,EquilGrid,GridView>::
1097localIdxToGlobalIdx(
unsigned localIdx)
const
1099 if (!isParallel()) {
1103 if (this->localIdxToGlobalIdx_.empty()) {
1104 throw std::logic_error(
"index map is not created on this rank");
1107 if (localIdx >= this->localIdxToGlobalIdx_.size()) {
1108 throw std::logic_error(
"local index is outside map range");
1111 return this->localIdxToGlobalIdx_[localIdx];
1114template <
class Gr
id,
class EquilGr
id,
class Gr
idView>
1115bool CollectDataOnIORank<Grid,EquilGrid,GridView>::
1116isCartIdxOnThisRank(
int cartIdx)
const
1118 if (! this->isParallel()) {
1122 assert (! needsReordering);
1124 auto candidate = std::lower_bound(this->sortedCartesianIdx_.begin(),
1125 this->sortedCartesianIdx_.end(),
1128 return (candidate != sortedCartesianIdx_.end())
1129 && (*candidate == cartIdx);
Definition CollectDataOnIORank.hpp:49
Definition CollectDataOnIORank.hpp:56
Communication handle to scatter the global index.
Definition CollectDataOnIORank_impl.hpp:200
Inter-region flow accumulation maps for all region definition arrays.
Definition InterRegFlows.hpp:179
void read(MessageBufferType &buffer)
Reconstitute internal object representation from MPI message buffer.
Definition InterRegFlows.hpp:323
void write(MessageBufferType &buffer) const
Serialise internal representation to MPI message buffer.
Definition InterRegFlows.hpp:296
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:45
Simple container for FLOWS data.
Definition FlowsData.hpp:32