23#ifndef OPM_ECL_GENERIC_WRITER_IMPL_HPP
24#define OPM_ECL_GENERIC_WRITER_IMPL_HPP
26#include <dune/grid/common/mcmgmapper.hh>
28#include <opm/grid/cpgrid/LgrOutputHelpers.hpp>
29#include <opm/grid/GridHelpers.hpp>
30#include <opm/grid/utility/cartesianToCompressed.hpp>
32#include <opm/input/eclipse/EclipseState/EclipseState.hpp>
33#include <opm/input/eclipse/EclipseState/Grid/RegionSetMatcher.hpp>
34#include <opm/input/eclipse/EclipseState/Grid/NNC.hpp>
35#include <opm/input/eclipse/EclipseState/SummaryConfig/SummaryConfig.hpp>
37#include <opm/input/eclipse/Schedule/Action/State.hpp>
38#include <opm/input/eclipse/Schedule/RPTConfig.hpp>
39#include <opm/input/eclipse/Schedule/Schedule.hpp>
40#include <opm/input/eclipse/Schedule/SummaryState.hpp>
41#include <opm/input/eclipse/Schedule/UDQ/UDQConfig.hpp>
42#include <opm/input/eclipse/Schedule/UDQ/UDQState.hpp>
43#include <opm/input/eclipse/Schedule/Well/WellConnections.hpp>
44#include <opm/input/eclipse/Schedule/Well/WellMatcher.hpp>
46#include <opm/input/eclipse/Units/UnitSystem.hpp>
48#include <opm/output/eclipse/EclipseIO.hpp>
49#include <opm/output/eclipse/RestartValue.hpp>
50#include <opm/output/eclipse/Summary.hpp>
70#include <unordered_map>
90bool directVerticalNeighbors(
const std::array<int, 3>& cartDims,
91 const std::unordered_map<int,int>& cartesianToActive,
92 int smallGlobalIndex,
int largeGlobalIndex)
94 assert(smallGlobalIndex <= largeGlobalIndex);
95 std::array<int, 3> ijk1, ijk2;
96 auto globalToIjk = [cartDims](
int gc) {
97 std::array<int, 3> ijk;
98 ijk[0] = gc % cartDims[0];
100 ijk[1] = gc % cartDims[1];
101 ijk[2] = gc / cartDims[1];
104 ijk1 = globalToIjk(smallGlobalIndex);
105 ijk2 = globalToIjk(largeGlobalIndex);
106 assert(ijk2[2]>=ijk1[2]);
108 if ( ijk1[0] == ijk2[0] && ijk1[1] == ijk2[1] && (ijk2[2] - ijk1[2]) > 1)
110 assert((largeGlobalIndex-smallGlobalIndex)%(cartDims[0]*cartDims[1])==0);
111 for (
int gi = smallGlobalIndex + cartDims[0] * cartDims[1]; gi < largeGlobalIndex;
112 gi += cartDims[0] * cartDims[1] )
114 if ( cartesianToActive.find( gi ) != cartesianToActive.end() )
124std::unordered_map<std::string, Opm::data::InterRegFlowMap>
127 auto maps = std::unordered_map<std::string, Opm::data::InterRegFlowMap>{};
129 const auto& regionNames = map.
names();
131 const auto nmap = regionNames.size();
134 for (
auto mapID = 0*nmap; mapID < nmap; ++mapID) {
135 maps.emplace(regionNames[mapID], std::move(flows[mapID]));
143 Opm::Action::State actionState_;
144 Opm::WellTestState wtestState_;
145 Opm::SummaryState summaryState_;
146 Opm::UDQState udqState_;
147 Opm::EclipseIO& eclIO_;
149 std::optional<int> timeStepNum_;
151 double secondsElapsed_;
152 std::vector<Opm::RestartValue> restartValue_;
153 bool writeDoublePrecision_;
155 bool forcedSimulationFinished_;
157 explicit EclWriteTasklet(
const Opm::Action::State& actionState,
158 const Opm::WellTestState& wtestState,
159 const Opm::SummaryState& summaryState,
160 const Opm::UDQState& udqState,
161 Opm::EclipseIO& eclIO,
163 std::optional<int> timeStepNum,
165 double secondsElapsed,
166 std::vector<Opm::RestartValue> restartValue,
167 bool writeDoublePrecision,
168 bool forcedSimulationFinished)
169 : actionState_(actionState)
170 , wtestState_(wtestState)
171 , summaryState_(summaryState)
172 , udqState_(udqState)
174 , reportStepNum_(reportStepNum)
175 , timeStepNum_(timeStepNum)
176 , isSubStep_(isSubStep)
177 , secondsElapsed_(secondsElapsed)
178 , restartValue_(std::move(restartValue))
179 , writeDoublePrecision_(writeDoublePrecision)
180 , forcedSimulationFinished_(forcedSimulationFinished)
186 if (this->restartValue_.size() == 1) {
187 this->eclIO_.writeTimeStep(this->actionState_,
191 this->reportStepNum_,
193 this->secondsElapsed_,
194 std::move(this->restartValue_.back()),
195 this->writeDoublePrecision_,
197 forcedSimulationFinished_);
200 this->eclIO_.writeTimeStep(this->actionState_,
204 this->reportStepNum_,
206 this->secondsElapsed_,
207 std::move(this->restartValue_),
208 this->writeDoublePrecision_,
210 forcedSimulationFinished_);
219template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
222 const EclipseState& eclState,
223 const SummaryConfig& summaryConfig,
225 const EquilGrid* equilGrid,
226 const GridView& gridView,
229 bool enableAsyncOutput,
231 : collectOnIORank_(grid,
236 summaryConfig.fip_regions_interreg_flow())
238 , gridView_ (gridView)
239 , schedule_ (schedule)
240 , eclState_ (eclState)
241 , cartMapper_ (cartMapper)
242 , equilCartMapper_(equilCartMapper)
243 , equilGrid_ (equilGrid)
249 this->
eclIO_ = std::make_unique<EclipseIO>
251 UgGridHelpers::createEclipseGrid(*equilGrid,
eclState_.getInputGrid()),
252 this->schedule_, summaryConfig,
"", enableEsmry);
257 int numWorkerThreads = 0;
259 numWorkerThreads = 1;
265template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
273template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
277 if (collectOnIORank_.isIORank()) {
278 std::map<std::string, std::vector<int>> integerVectors;
279 if (collectOnIORank_.isParallel()) {
280 integerVectors.emplace(
"MPI_RANK", collectOnIORank_.globalRanks());
283 if (
const auto& lgrs = this->eclState_.getLgrs(); lgrs.size() > 0) {
285 const auto nncCollection = Opm::NNCCollection::fromLGROutputContainers(this->outputNnc_,
286 this->outputNncGlobalLocal_,
287 this->outputAmalgamatedNnc_);
288 eclIO_->writeInitial(*this->outputTrans_,
292 eclIO_->writeInitial(*this->outputTrans_,
294 this->outputNnc_.front());
296 this->outputTrans_.reset();
300template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
305 if (collectOnIORank_.isIORank()) {
306 constexpr bool equilGridIsCpGrid = std::is_same_v<EquilGrid, Dune::CpGrid>;
308 const auto levelCartMapp = this->createLevelCartMapp_<equilGridIsCpGrid>();
309 const auto levelCartToLevelCompressed = this->createCartesianToActiveMaps_<equilGridIsCpGrid>(levelCartMapp);
310 auto computeLevelIndices = this->computeLevelIndices_<equilGridIsCpGrid>();
311 auto computeLevelCartIdx = this->computeLevelCartIdx_<equilGridIsCpGrid>(levelCartMapp, *(this->equilCartMapper_));
312 auto computeLevelCartDimensions = this->computeLevelCartDimensions_<equilGridIsCpGrid>(levelCartMapp, *(this->equilCartMapper_));
313 auto computeOriginIndices = this->computeOriginIndices_<equilGridIsCpGrid>();
315 computeTrans_(levelCartToLevelCompressed, map, computeLevelIndices,
316 computeLevelCartIdx, computeLevelCartDimensions, computeOriginIndices);
317 exportNncStructure_(levelCartToLevelCompressed, map, computeLevelIndices, computeLevelCartIdx,
318 computeLevelCartDimensions, computeOriginIndices);
322 if (collectOnIORank_.isParallel()) {
323 const auto& comm = grid_.comm();
330template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
335 const auto& numAquCell = this->eclState_.aquifer().hasNumericalAquifer()
336 ? this->eclState_.aquifer().numericalAquifers().allAquiferCellIds()
337 : std::vector<std::size_t>{};
339 return std::ranges::binary_search(numAquCell.begin(), numAquCell.end(), cartIdx);
342template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
344EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
345isNumAquConn_(
const std::size_t cartIdx1,
346 const std::size_t cartIdx2)
const
348 return isNumAquCell_(cartIdx1) || isNumAquCell_(cartIdx2);
351template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
352template<
bool equilGr
idIsCpGr
id>
354EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
355createLevelCartMapp_()
const
357 if constexpr (equilGridIsCpGrid) {
363template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
364template<
bool equilGr
idIsCpGr
id>
365std::vector<std::unordered_map<int,int>>
366EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
369 if constexpr (equilGridIsCpGrid) {
370 if (this->equilGrid_->maxLevel()) {
371 return Opm::Lgr::levelCartesianToLevelCompressedMaps(*this->equilGrid_, levelCartMapp); }
373 return std::vector<std::unordered_map<int,int>>{ cartesianToCompressed(equilGrid_->size(0), UgGridHelpers::globalCell(*equilGrid_)) };
376 return std::vector<std::unordered_map<int,int>>{ cartesianToCompressed(equilGrid_->size(0), UgGridHelpers::globalCell(*equilGrid_)) };
379template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
380template<
bool equilGr
idIsCpGr
id>
381std::function<std::array<int,3>(
int)>
382EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
386 if constexpr (equilGridIsCpGrid) {
387 return [&](
int level)
389 return levelCartMapp.cartesianDimensions(level);
393 return [&]([[maybe_unused]]
int level)
396 return equilCartMapp.cartesianDimensions();
401template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
402template<
bool equilGr
idIsCpGr
id>
403std::function<int(
int,
int)>
404EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
408 if constexpr (equilGridIsCpGrid) {
409 return [&](
int levelCompressedIdx,
412 return levelCartMapp.cartesianIndex(levelCompressedIdx, level);
416 return [&](
int levelCompressedIdx,
417 [[maybe_unused]]
int level)
420 return equilCartMapp.cartesianIndex(levelCompressedIdx);
425template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
426template <
bool equilGr
idIsCpGr
id>
428EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
429computeLevelIndices_()
const
431 if constexpr (equilGridIsCpGrid) {
432 return [](
const auto& intersection,
436 return std::pair{intersection.inside().getLevelElem().index(), intersection.outside().getLevelElem().index()};
440 return [](
const auto&,
441 const auto& intersectionInsideLeafIdx,
442 const auto& intersectionOutsideLeafIdx)
444 return std::pair{intersectionInsideLeafIdx, intersectionOutsideLeafIdx};
449template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
450template <
bool equilGr
idIsCpGr
id>
452EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
453computeOriginIndices_()
const
455 if constexpr (equilGridIsCpGrid) {
456 return [](
const auto& intersection,
460 return std::pair{intersection.inside().getOrigin().index(), intersection.outside().getOrigin().index()};
464 return [](
const auto&,
465 const auto& intersectionInsideLeafIdx,
466 const auto& intersectionOutsideLeafIdx)
468 return std::pair{intersectionInsideLeafIdx, intersectionOutsideLeafIdx};
473template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
475EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
476allocateLevelTrans_(
const std::array<int,3>& levelCartDims,
477 data::Solution& levelTrans)
const
479 auto createLevelCellData = [&levelCartDims]() {
480 return Opm::data::CellData{
481 Opm::UnitSystem::measure::transmissibility,
482 std::vector<double>(levelCartDims[0] * levelCartDims[1] * levelCartDims[2], 0.0),
483 Opm::data::TargetType::INIT
488 levelTrans.emplace(
"TRANX", createLevelCellData());
489 levelTrans.emplace(
"TRANY", createLevelCellData());
490 levelTrans.emplace(
"TRANZ", createLevelCellData());
493template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
494template<
typename LevelIndicesFunction,
typename OriginIndicesFunction>
496EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
497computeTrans_(
const std::vector<std::unordered_map<int,int>>& levelCartToLevelCompressed,
498 const std::function<
unsigned int(
unsigned int)>& map,
499 const LevelIndicesFunction& computeLevelIndices,
500 const std::function<
int(
int,
int)>& computeLevelCartIdx,
501 const std::function<std::array<int,3>(
int)>& computeLevelCartDims,
502 const OriginIndicesFunction& computeOriginIndices)
const
505 outputTrans_ = std::make_unique<std::vector<data::Solution>>(std::vector<data::Solution>{});
508 using GlobalGridView =
typename EquilGrid::LeafGridView;
509 using GlobElementMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GlobalGridView>;
510 const GlobalGridView& globalGridView = this->equilGrid_->leafGridView();
511 const GlobElementMapper globalElemMapper { globalGridView, Dune::mcmgElementLayout() };
514 int maxLevel = this->equilGrid_->maxLevel();
516 outputTrans_->resize(maxLevel+1);
518 for (
int level = 0; level <= maxLevel; ++level) {
519 allocateLevelTrans_(computeLevelCartDims(level), this->outputTrans_->at(level));
522 for (
const auto& elem : elements(globalGridView)) {
523 for (
const auto& is : intersections(globalGridView, elem)) {
527 if ( is.inside().level() != is.outside().level() )
531 unsigned c1 = globalElemMapper.index(is.inside());
532 unsigned c2 = globalElemMapper.index(is.outside());
537 int level = is.inside().level();
540 const auto& [levelInIdx, levelOutIdx] = computeLevelIndices(is, c1, c2);
542 const int levelCartIdxIn = computeLevelCartIdx(levelInIdx, level);
543 const int levelCartIdxOut = computeLevelCartIdx(levelOutIdx, level);
549 const auto [originInIdx, originOutIdx] = computeOriginIndices(is, c1, c2);
551 const auto originCartIdxIn = computeLevelCartIdx(originInIdx, 0);
552 const auto originCartIdxOut = computeLevelCartIdx(originOutIdx, 0);
555 if (isNumAquCell_(originCartIdxIn) || isNumAquCell_(originCartIdxOut)) {
565 const auto minLevelCartIdx = std::min(levelCartIdxIn, levelCartIdxOut);
566 const auto maxLevelCartIdx = std::max(levelCartIdxIn, levelCartIdxOut);
568 const auto& levelCartDims = computeLevelCartDims(level);
576 if (maxLevelCartIdx - minLevelCartIdx == 1 && levelCartDims[0] > 1 ) {
577 outputTrans_->at(level).at(
"TRANX").template data<double>()[minLevelCartIdx] = globalTrans().transmissibility(c1, c2);
581 if (maxLevelCartIdx - minLevelCartIdx == levelCartDims[0] && levelCartDims[1] > 1) {
582 outputTrans_->at(level).at(
"TRANY").template data<double>()[minLevelCartIdx] = globalTrans().transmissibility(c1, c2);
586 if ( maxLevelCartIdx - minLevelCartIdx == levelCartDims[0]*levelCartDims[1] ||
587 directVerticalNeighbors(levelCartDims,
588 levelCartToLevelCompressed[level],
591 outputTrans_->at(level).at(
"TRANZ").template data<double>()[minLevelCartIdx] = globalTrans().transmissibility(c1, c2);
597template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
599EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
600isCartesianNeighbour_(
const std::array<int,3>& levelCartDims,
601 const std::size_t levelCartIdx1,
602 const std::size_t levelCartIdx2)
const
604 const int diff = levelCartIdx2 - levelCartIdx1;
607 || (diff == levelCartDims[0])
608 || (diff == (levelCartDims[0] * levelCartDims[1]));
611template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
613EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
614isDirectNeighbours_(
const std::unordered_map<int,int>& levelCartesianToActive,
615 const std::array<int,3>& levelCartDims,
616 const std::size_t levelCartIdx1,
617 const std::size_t levelCartIdx2)
const
619 return isCartesianNeighbour_(levelCartDims, levelCartIdx1, levelCartIdx2)
620 || directVerticalNeighbors(levelCartDims, levelCartesianToActive, levelCartIdx1, levelCartIdx2);
623template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
625EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
626activeCell_(
const std::unordered_map<int,int>& levelCartToLevelCompressed,
627 const std::size_t levelCartIdx)
const
629 auto pos = levelCartToLevelCompressed.find(levelCartIdx);
630 return (pos == levelCartToLevelCompressed.end()) ? -1 : pos->second;
633template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
635EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
636allocateAllNncs_(
int maxLevel)
const
638 this->outputNnc_.resize(maxLevel+1);
646 this->outputNncGlobalLocal_.resize(maxLevel);
654 this->outputAmalgamatedNnc_.resize(maxLevel-1);
655 for (
int i = 0; i < maxLevel-1; ++i) {
656 this->outputAmalgamatedNnc_[i].resize(maxLevel-1-i);
661template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
662template<
typename LevelIndicesFunction,
typename OriginIndicesFunction>
663std::vector<std::vector<NNCdata>>
664EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
665exportNncStructure_(
const std::vector<std::unordered_map<int,int>>& levelCartToLevelCompressed,
666 const std::function<
unsigned int(
unsigned int)>& map,
667 const LevelIndicesFunction& computeLevelIndices,
668 const std::function<
int(
int,
int)>& computeLevelCartIdx,
669 const std::function<std::array<int,3>(
int)>& computeLevelCartDims,
670 const OriginIndicesFunction& computeOriginIndices)
const
672 const auto& nncData = this->eclState_.getInputNNC().input();
673 const auto& nncEdit = this->eclState_.getInputNNC().edit();
674 const auto& nncEditr = this->eclState_.getInputNNC().editr();
675 const auto& unitSystem = this->eclState_.getDeckUnitSystem();
676 const auto& transMult = this->eclState_.getTransMult();
679 const auto& equilCartMapper = *equilCartMapper_;
681 const auto& level0CartDims = equilCartMapper.cartesianDimensions();
683 int maxLevel = this->equilGrid_->maxLevel();
684 allocateAllNncs_(maxLevel);
686 using GlobalGridView =
typename EquilGrid::LeafGridView;
687 using GlobElementMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GlobalGridView>;
688 const GlobalGridView& globalGridView = this->equilGrid_->leafGridView();
689 const GlobElementMapper globalElemMapper { globalGridView, Dune::mcmgElementLayout() };
691 for (
const auto& elem : elements(globalGridView)) {
692 for (
const auto& is : intersections(globalGridView, elem)) {
697 unsigned c1 = globalElemMapper.index(is.inside());
698 unsigned c2 = globalElemMapper.index(is.outside());
703 if ( is.inside().level() != is.outside().level() ) {
705 const auto& [levelInIdx, levelOutIdx] = computeLevelIndices(is, c1, c2);
707 const int levelIn = is.inside().level();
708 const int levelOut = is.outside().level();
710 auto levelCartIdxIn = computeLevelCartIdx(levelInIdx, levelIn);
711 auto levelCartIdxOut = computeLevelCartIdx(levelOutIdx, levelOut);
714 std::pair<int,int> smallerPair = {levelIn, levelCartIdxIn},
715 largerPair = {levelOut, levelCartIdxOut};
716 if (smallerPair.first > largerPair.first) {
717 std::swap(smallerPair, largerPair);
720 const auto& [smallerLevel, smallerLevelCartIdx] = smallerPair;
721 const auto& [largerLevel, largerLevelCartIdx] = largerPair;
723 auto t = this->globalTrans().transmissibility(c1, c2);
730 const auto tt = unitSystem
731 .from_si(UnitSystem::measure::transmissibility, t);
733 if (std::isnormal(tt) && (tt > 1.0e-12)) {
735 if (smallerLevel == 0) {
736 this->outputNncGlobalLocal_[largerLevel-1].emplace_back(smallerLevelCartIdx, largerLevelCartIdx, t);
739 assert(smallerLevel >= 1);
740 this->outputAmalgamatedNnc_[smallerLevel-1][largerLevel-smallerLevel-1].emplace_back(smallerLevelCartIdx, largerLevelCartIdx, t);
746 assert(is.inside().level() == is.outside().level());
747 const int level = is.inside().level();
753 const auto [originInIdx, originOutIdx] = computeOriginIndices(is, c1, c2);
755 const std::size_t originCartIdxIn = computeLevelCartIdx(originInIdx, 0);
756 const std::size_t originCartIdxOut = computeLevelCartIdx(originOutIdx, 0);
759 const auto& [levelInIdx, levelOutIdx] = computeLevelIndices(is, c1, c2);
761 auto levelCartIdxIn = computeLevelCartIdx(levelInIdx, level);
762 auto levelCartIdxOut = computeLevelCartIdx(levelOutIdx, level);
764 if ( levelCartIdxOut < levelCartIdxIn )
765 std::swap(levelCartIdxIn, levelCartIdxOut);
773 const auto& levelCartDims = computeLevelCartDims(level);
776 assert(!isNumAquConn_(originCartIdxIn, originCartIdxOut) || level == 0);
778 if (isNumAquConn_(originCartIdxIn, originCartIdxOut) ||
779 ! isDirectNeighbours_(levelCartToLevelCompressed[level],
781 levelCartIdxIn, levelCartIdxOut)) {
784 auto t = this->globalTrans().transmissibility(c1, c2);
787 auto candidate = std::lower_bound(nncData.begin(), nncData.end(),
788 NNCdata { originCartIdxIn, originCartIdxOut, 0.0 });
789 const auto transMlt = transMult.getRegionMultiplierNNC(originCartIdxIn, originCartIdxOut);
790 bool foundNncEditr =
false;
792 while ((candidate != nncData.end()) &&
793 (candidate->cell1 == originCartIdxIn) &&
794 (candidate->cell2 == originCartIdxOut))
796 auto trans = candidate->trans;
798 if (! nncEditr.empty()) {
799 auto it = std::lower_bound(nncEditr.begin(), nncEditr.end(),
800 NNCdata { originCartIdxIn, originCartIdxOut, 0.0 });
801 foundNncEditr = it != nncEditr.end() && it->cell1 == originCartIdxIn && it->cell2 == originCartIdxOut;
807 if (! nncEdit.empty()) {
808 auto it = std::lower_bound(nncEdit.begin(), nncEdit.end(),
809 NNCdata { originCartIdxIn, originCartIdxOut, 0.0 });
810 if (it != nncEdit.end() && it->cell1 == originCartIdxIn && it->cell2 == originCartIdxOut) {
828 const auto tt = unitSystem
829 .from_si(UnitSystem::measure::transmissibility, t);
831 if (std::isnormal(tt) && (tt > 1.0e-12)) {
832 this->outputNnc_[level].emplace_back(levelCartIdxIn, levelCartIdxOut, t);
840 std::vector<NNCdata> inputedNnc{};
841 const auto generatedNnc = outputNnc_[0];
845 for (
const auto& entry : nncData) {
854 assert (entry.cell2 >= entry.cell1);
856 if (! isCartesianNeighbour_(level0CartDims, entry.cell1, entry.cell2) ||
857 isNumAquConn_(entry.cell1, entry.cell2))
859 bool foundNncEdit =
false;
860 auto trans = entry.trans;
861 if (! nncEdit.empty()) {
862 auto it = std::lower_bound(nncEdit.begin(), nncEdit.end(),
863 NNCdata {entry.cell1, entry.cell2, 0.0 });
864 if (it != nncEdit.end() && it->cell1 == entry.cell1 && it->cell2 == entry.cell2) {
869 if (! foundNncEdit) {
873 const auto c1 = activeCell_(levelCartToLevelCompressed[0], entry.cell1);
874 const auto c2 = activeCell_(levelCartToLevelCompressed[0], entry.cell2);
876 if ((c1 < 0) || (c2 < 0)) {
882 trans = this->globalTrans().transmissibility(c1, c2);
884 if (! generatedNnc.empty()) {
885 for (
const auto& generated : generatedNnc) {
886 if (entry.cell1 == generated.cell1 && entry.cell2 == generated.cell2) {
887 trans -= generated.trans;
893 const auto tt = unitSystem
894 .from_si(UnitSystem::measure::transmissibility, trans);
899 if (std::isnormal(tt) && ! (tt < 1.0e-6)) {
900 inputedNnc.emplace_back(entry.cell1, entry.cell2, trans);
905 this->outputNnc_[0].insert(this->outputNnc_[0].begin(), inputedNnc.begin(), inputedNnc.end());
906 return this->outputNnc_;
909template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
912 const std::optional<int> timeStepNum,
913 const bool isSubStep,
914 const bool isForcedFinalOutput,
915 data::Solution&& localCellData,
916 data::Wells&& localWellData,
917 data::GroupAndNetworkValues&& localGroupAndNetworkData,
918 data::Aquifers&& localAquiferData,
919 WellTestState&& localWTestState,
920 const Action::State& actionState,
921 const UDQState& udqState,
922 const SummaryState& summaryState,
923 const std::vector<Scalar>& thresholdPressure,
926 bool doublePrecision,
932 const auto isParallel = this->collectOnIORank_.isParallel();
933 const bool needsReordering = this->collectOnIORank_.doesNeedReordering();
935 RestartValue restartValue {
936 (isParallel || needsReordering)
937 ? this->collectOnIORank_.globalCellData()
938 : std::move(localCellData),
940 isParallel ? this->collectOnIORank_.globalWellData()
941 : std::move(localWellData),
943 isParallel ? this->collectOnIORank_.globalGroupAndNetworkData()
944 : std::move(localGroupAndNetworkData),
946 isParallel ? this->collectOnIORank_.globalAquiferData()
947 : std::move(localAquiferData)
950 if (eclState_.getSimulationConfig().useThresholdPressure()) {
951 restartValue.addExtra(
"THRESHPR", UnitSystem::measure::pressure,
957 restartValue.addExtra(
"OPMEXTRA", std::vector<double>(1, nextStepSize));
962 const auto flowsn_global = isParallel ? this->collectOnIORank_.globalFlowsn() : std::move(flowsn);
963 for (
const auto& flows : flowsn_global) {
964 if (flows.name.empty())
966 if (flows.name ==
"FLOGASN+") {
967 restartValue.addExtra(flows.name, UnitSystem::measure::gas_surface_rate, flows.values);
969 restartValue.addExtra(flows.name, UnitSystem::measure::liquid_surface_rate, flows.values);
974 const auto floresn_global = isParallel ? this->collectOnIORank_.globalFloresn() : std::move(floresn);
975 for (
const auto& flores : floresn_global) {
976 if (flores.name.empty()) {
979 restartValue.addExtra(flores.name, UnitSystem::measure::rate, flores.values);
983 std::vector<Opm::RestartValue> restartValues{};
985 if ( !isParallel && !needsReordering && (this->eclState_.getLgrs().size()>0) && (this->grid_.maxLevel()>0) ) {
990 Opm::Lgr::extractRestartValueLevelGrids<Grid>(this->grid_, restartValue, restartValues);
993 restartValues.reserve(1);
994 restartValues.push_back(std::move(restartValue));
1000 this->taskletRunner_->barrier();
1003 if (this->taskletRunner_->failure()) {
1004 throw std::runtime_error(
"Failure in the TaskletRunner while writing output.");
1008 auto eclWriteTasklet = std::make_shared<EclWriteTasklet>(
1010 isParallel ? this->collectOnIORank_.globalWellTestState() : std::move(localWTestState),
1011 summaryState, udqState, *this->eclIO_,
1012 reportStepNum, timeStepNum, isSubStep, curTime, std::move(restartValues), doublePrecision,
1013 isForcedFinalOutput);
1016 this->taskletRunner_->dispatch(std::move(eclWriteTasklet));
1019template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
1022 const Scalar curTime,
1023 const data::Wells& localWellData,
1024 const data::WellBlockAveragePressures& localWBPData,
1025 const data::GroupAndNetworkValues& localGroupAndNetworkData,
1026 const std::map<int,data::AquiferData>& localAquiferData,
1027 const std::map<std::pair<std::string, int>,
double>& blockData,
1028 const std::map<std::string, double>& miscSummaryData,
1029 const std::map<std::string, std::vector<double>>& regionData,
1030 const Inplace& inplace,
1031 const Inplace* initialInPlace,
1033 SummaryState& summaryState,
1035 const data::ReservoirCouplingGroupRates* rcGroupRates)
1037 if (collectOnIORank_.isIORank()) {
1038 const auto& summary = eclIO_->summary();
1040 const auto& wellData = this->collectOnIORank_.isParallel()
1041 ? this->collectOnIORank_.globalWellData()
1044 const auto& wbpData = this->collectOnIORank_.isParallel()
1045 ? this->collectOnIORank_.globalWBPData()
1048 const auto& groupAndNetworkData = this->collectOnIORank_.isParallel()
1049 ? this->collectOnIORank_.globalGroupAndNetworkData()
1050 : localGroupAndNetworkData;
1052 const auto& aquiferData = this->collectOnIORank_.isParallel()
1053 ? this->collectOnIORank_.globalAquiferData()
1056 const auto interreg_flows = getInterRegFlowsAsMap(interRegFlows);
1058 const auto values = out::Summary::DynamicSimulatorState {
1061 &groupAndNetworkData,
1074 summary.eval(reportStepNum, curTime, values, summaryState);
1080 const auto udq_step = reportStepNum - 1;
1082 this->schedule_[udq_step].udq()
1084 this->schedule_.wellMatcher(udq_step),
1085 this->schedule_[udq_step].group_order(),
1086 this->schedule_.segmentMatcherFactory(udq_step),
1087 [es = std::cref(this->eclState_)]() {
1088 return std::make_unique<RegionSetMatcher>
1089 (es.get().fipRegionStatistics());
1091 summaryState, udqState);
1095 if (collectOnIORank_.isParallel()) {
1097 ser.
append(summaryState);
1102template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
1107 assert (globalTrans_);
1108 return *globalTrans_;
Definition: CollectDataOnIORank.hpp:49
bool isIORank() const
Definition: CollectDataOnIORank.hpp:126
Definition: EclGenericWriter.hpp:69
std::vector< std::vector< NNCdata > > outputNnc_
Definition: EclGenericWriter.hpp:178
const EclipseState & eclState_
Definition: EclGenericWriter.hpp:162
CollectDataOnIORankType collectOnIORank_
Definition: EclGenericWriter.hpp:158
void extractOutputTransAndNNC(const std::function< unsigned int(unsigned int)> &map)
Definition: EclGenericWriter_impl.hpp:303
void evalSummary(int reportStepNum, Scalar curTime, const data::Wells &localWellData, const data::WellBlockAveragePressures &localWBPData, const data::GroupAndNetworkValues &localGroupAndNetworkData, const std::map< int, data::AquiferData > &localAquiferData, const std::map< std::pair< std::string, int >, double > &blockData, const std::map< std::string, double > &miscSummaryData, const std::map< std::string, std::vector< double > > ®ionData, const Inplace &inplace, const Inplace *initialInPlace, const InterRegFlowMap &interRegFlows, SummaryState &summaryState, UDQState &udqState, const data::ReservoirCouplingGroupRates *rcGroupRates=nullptr)
Definition: EclGenericWriter_impl.hpp:1021
std::unique_ptr< TaskletRunner > taskletRunner_
Definition: EclGenericWriter.hpp:164
void writeInit()
Definition: EclGenericWriter_impl.hpp:275
EclGenericWriter(const Schedule &schedule, const EclipseState &eclState, const SummaryConfig &summaryConfig, const Grid &grid, const EquilGrid *equilGrid, const GridView &gridView, const Dune::CartesianIndexMapper< Grid > &cartMapper, const Dune::CartesianIndexMapper< EquilGrid > *equilCartMapper, bool enableAsyncOutput, bool enableEsmry)
Definition: EclGenericWriter_impl.hpp:221
const EclipseIO & eclIO() const
Definition: EclGenericWriter_impl.hpp:267
const TransmissibilityType & globalTrans() const
Definition: EclGenericWriter_impl.hpp:1105
std::unique_ptr< EclipseIO > eclIO_
Definition: EclGenericWriter.hpp:163
void doWriteOutput(const int reportStepNum, const std::optional< int > timeStepNum, const bool isSubStep, const bool forcedSimulationFinished, data::Solution &&localCellData, data::Wells &&localWellData, data::GroupAndNetworkValues &&localGroupAndNetworkData, data::Aquifers &&localAquiferData, WellTestState &&localWTestState, const Action::State &actionState, const UDQState &udqState, const SummaryState &summaryState, const std::vector< Scalar > &thresholdPressure, Scalar curTime, Scalar nextStepSize, bool doublePrecision, bool isFlowsn, std::array< FlowsData< double >, 3 > &&flowsn, bool isFloresn, std::array< FlowsData< double >, 3 > &&floresn)
Definition: EclGenericWriter_impl.hpp:911
Inter-region flow accumulation maps for all region definition arrays.
Definition: InterRegFlows.hpp:179
std::vector< data::InterRegFlowMap > getInterRegFlows() const
const std::vector< std::string > & names() const
Definition: EclGenericWriter.hpp:52
Class for serializing and broadcasting data using MPI.
Definition: MPISerializer.hpp:38
void append(T &data, int root=0)
Serialize and broadcast on root process, de-serialize and append on others.
Definition: MPISerializer.hpp:82
void broadcast(RootRank rootrank, Args &&... args)
Definition: MPISerializer.hpp:47
The base class for tasklets.
Definition: tasklets.hpp:45
Handles where a given tasklet is run.
Definition: tasklets.hpp:93
Definition: Transmissibility.hpp:54
Definition: blackoilbioeffectsmodules.hh:45
Avoid mistakes in calls to broadcast() by wrapping the root argument in an explicit type.
Definition: MPISerializer.hpp:33