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/SummaryConfig/SummaryConfig.hpp>
36#include <opm/input/eclipse/Schedule/Action/State.hpp>
37#include <opm/input/eclipse/Schedule/RPTConfig.hpp>
38#include <opm/input/eclipse/Schedule/Schedule.hpp>
39#include <opm/input/eclipse/Schedule/SummaryState.hpp>
40#include <opm/input/eclipse/Schedule/UDQ/UDQConfig.hpp>
41#include <opm/input/eclipse/Schedule/UDQ/UDQState.hpp>
42#include <opm/input/eclipse/Schedule/Well/WellConnections.hpp>
43#include <opm/input/eclipse/Schedule/Well/WellMatcher.hpp>
45#include <opm/input/eclipse/Units/UnitSystem.hpp>
47#include <opm/output/eclipse/EclipseIO.hpp>
48#include <opm/output/eclipse/RestartValue.hpp>
49#include <opm/output/eclipse/Summary.hpp>
69#include <unordered_map>
89bool directVerticalNeighbors(
const std::array<int, 3>& cartDims,
90 const std::unordered_map<int,int>& cartesianToActive,
91 int smallGlobalIndex,
int largeGlobalIndex)
93 assert(smallGlobalIndex <= largeGlobalIndex);
94 std::array<int, 3> ijk1, ijk2;
95 auto globalToIjk = [cartDims](
int gc) {
96 std::array<int, 3> ijk;
97 ijk[0] = gc % cartDims[0];
99 ijk[1] = gc % cartDims[1];
100 ijk[2] = gc / cartDims[1];
103 ijk1 = globalToIjk(smallGlobalIndex);
104 ijk2 = globalToIjk(largeGlobalIndex);
105 assert(ijk2[2]>=ijk1[2]);
107 if ( ijk1[0] == ijk2[0] && ijk1[1] == ijk2[1] && (ijk2[2] - ijk1[2]) > 1)
109 assert((largeGlobalIndex-smallGlobalIndex)%(cartDims[0]*cartDims[1])==0);
110 for (
int gi = smallGlobalIndex + cartDims[0] * cartDims[1]; gi < largeGlobalIndex;
111 gi += cartDims[0] * cartDims[1] )
113 if ( cartesianToActive.find( gi ) != cartesianToActive.end() )
123std::unordered_map<std::string, Opm::data::InterRegFlowMap>
126 auto maps = std::unordered_map<std::string, Opm::data::InterRegFlowMap>{};
128 const auto& regionNames = map.
names();
130 const auto nmap = regionNames.size();
133 for (
auto mapID = 0*nmap; mapID < nmap; ++mapID) {
134 maps.emplace(regionNames[mapID], std::move(flows[mapID]));
142 Opm::Action::State actionState_;
143 Opm::WellTestState wtestState_;
144 Opm::SummaryState summaryState_;
145 Opm::UDQState udqState_;
146 Opm::EclipseIO& eclIO_;
148 std::optional<int> timeStepNum_;
150 double secondsElapsed_;
151 std::vector<Opm::RestartValue> restartValue_;
152 bool writeDoublePrecision_;
154 explicit EclWriteTasklet(
const Opm::Action::State& actionState,
155 const Opm::WellTestState& wtestState,
156 const Opm::SummaryState& summaryState,
157 const Opm::UDQState& udqState,
158 Opm::EclipseIO& eclIO,
160 std::optional<int> timeStepNum,
162 double secondsElapsed,
163 std::vector<Opm::RestartValue> restartValue,
164 bool writeDoublePrecision)
165 : actionState_(actionState)
166 , wtestState_(wtestState)
167 , summaryState_(summaryState)
168 , udqState_(udqState)
170 , reportStepNum_(reportStepNum)
171 , timeStepNum_(timeStepNum)
172 , isSubStep_(isSubStep)
173 , secondsElapsed_(secondsElapsed)
174 , restartValue_(std::move(restartValue))
175 , writeDoublePrecision_(writeDoublePrecision)
181 if (this->restartValue_.size() == 1) {
182 this->eclIO_.writeTimeStep(this->actionState_,
186 this->reportStepNum_,
188 this->secondsElapsed_,
189 std::move(this->restartValue_.back()),
190 this->writeDoublePrecision_,
194 this->eclIO_.writeTimeStep(this->actionState_,
198 this->reportStepNum_,
200 this->secondsElapsed_,
201 std::move(this->restartValue_),
202 this->writeDoublePrecision_,
212template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
215 const EclipseState& eclState,
216 const SummaryConfig& summaryConfig,
218 const EquilGrid* equilGrid,
219 const GridView& gridView,
222 bool enableAsyncOutput,
224 : collectOnIORank_(grid,
229 summaryConfig.fip_regions_interreg_flow())
231 , gridView_ (gridView)
232 , schedule_ (schedule)
233 , eclState_ (eclState)
234 , cartMapper_ (cartMapper)
235 , equilCartMapper_(equilCartMapper)
236 , equilGrid_ (equilGrid)
242 this->
eclIO_ = std::make_unique<EclipseIO>
244 UgGridHelpers::createEclipseGrid(*equilGrid,
eclState_.getInputGrid()),
245 this->schedule_, summaryConfig,
"", enableEsmry);
250 int numWorkerThreads = 0;
252 numWorkerThreads = 1;
258template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
266template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
270 if (collectOnIORank_.isIORank()) {
271 std::map<std::string, std::vector<int>> integerVectors;
272 if (collectOnIORank_.isParallel()) {
273 integerVectors.emplace(
"MPI_RANK", collectOnIORank_.globalRanks());
276 eclIO_->writeInitial(*this->outputTrans_,
278 this->outputNnc_.front());
279 this->outputTrans_.reset();
283template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
288 if (collectOnIORank_.isIORank()) {
289 constexpr bool equilGridIsCpGrid = std::is_same_v<EquilGrid, Dune::CpGrid>;
291 const auto levelCartMapp = this->createLevelCartMapp_<equilGridIsCpGrid>();
292 const auto levelCartToLevelCompressed = this->createCartesianToActiveMaps_<equilGridIsCpGrid>(levelCartMapp);
293 auto computeLevelIndices = this->computeLevelIndices_<equilGridIsCpGrid>();
294 auto computeLevelCartIdx = this->computeLevelCartIdx_<equilGridIsCpGrid>(levelCartMapp, *(this->equilCartMapper_));
295 auto computeLevelCartDimensions = this->computeLevelCartDimensions_<equilGridIsCpGrid>(levelCartMapp, *(this->equilCartMapper_));
296 auto computeOriginIndices = this->computeOriginIndices_<equilGridIsCpGrid>();
298 computeTrans_(levelCartToLevelCompressed, map, computeLevelIndices,
299 computeLevelCartIdx, computeLevelCartDimensions, computeOriginIndices);
300 exportNncStructure_(levelCartToLevelCompressed, map, computeLevelIndices, computeLevelCartIdx,
301 computeLevelCartDimensions, computeOriginIndices);
305 if (collectOnIORank_.isParallel()) {
306 const auto& comm = grid_.comm();
313template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
318 const auto& numAquCell = this->eclState_.aquifer().hasNumericalAquifer()
319 ? this->eclState_.aquifer().numericalAquifers().allAquiferCellIds()
320 : std::vector<std::size_t>{};
322 return std::ranges::binary_search(numAquCell.begin(), numAquCell.end(), cartIdx);
325template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
327EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
328isNumAquConn_(
const std::size_t cartIdx1,
329 const std::size_t cartIdx2)
const
331 return isNumAquCell_(cartIdx1) || isNumAquCell_(cartIdx2);
334template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
335template<
bool equilGr
idIsCpGr
id>
337EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
338createLevelCartMapp_()
const
340 if constexpr (equilGridIsCpGrid) {
346template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
347template<
bool equilGr
idIsCpGr
id>
348std::vector<std::unordered_map<int,int>>
349EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
352 if constexpr (equilGridIsCpGrid) {
353 if (this->equilGrid_->maxLevel()) {
354 return Opm::Lgr::levelCartesianToLevelCompressedMaps(*this->equilGrid_, levelCartMapp); }
356 return std::vector<std::unordered_map<int,int>>{ cartesianToCompressed(equilGrid_->size(0), UgGridHelpers::globalCell(*equilGrid_)) };
359 return std::vector<std::unordered_map<int,int>>{ cartesianToCompressed(equilGrid_->size(0), UgGridHelpers::globalCell(*equilGrid_)) };
362template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
363template<
bool equilGr
idIsCpGr
id>
364std::function<std::array<int,3>(
int)>
365EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
369 if constexpr (equilGridIsCpGrid) {
370 return [&](
int level)
372 return levelCartMapp.cartesianDimensions(level);
376 return [&]([[maybe_unused]]
int level)
379 return equilCartMapp.cartesianDimensions();
384template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
385template<
bool equilGr
idIsCpGr
id>
386std::function<int(
int,
int)>
387EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
391 if constexpr (equilGridIsCpGrid) {
392 return [&](
int levelCompressedIdx,
395 return levelCartMapp.cartesianIndex(levelCompressedIdx, level);
399 return [&](
int levelCompressedIdx,
400 [[maybe_unused]]
int level)
403 return equilCartMapp.cartesianIndex(levelCompressedIdx);
408template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
409template <
bool equilGr
idIsCpGr
id>
411EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
412computeLevelIndices_()
const
414 if constexpr (equilGridIsCpGrid) {
415 return [](
const auto& intersection,
419 return std::pair{intersection.inside().getLevelElem().index(), intersection.outside().getLevelElem().index()};
423 return [](
const auto&,
424 const auto& intersectionInsideLeafIdx,
425 const auto& intersectionOutsideLeafIdx)
427 return std::pair{intersectionInsideLeafIdx, intersectionOutsideLeafIdx};
432template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
433template <
bool equilGr
idIsCpGr
id>
435EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
436computeOriginIndices_()
const
438 if constexpr (equilGridIsCpGrid) {
439 return [](
const auto& intersection,
443 return std::pair{intersection.inside().getOrigin().index(), intersection.outside().getOrigin().index()};
447 return [](
const auto&,
448 const auto& intersectionInsideLeafIdx,
449 const auto& intersectionOutsideLeafIdx)
451 return std::pair{intersectionInsideLeafIdx, intersectionOutsideLeafIdx};
456template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
458EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
459allocateLevelTrans_(
const std::array<int,3>& levelCartDims,
460 data::Solution& levelTrans)
const
462 auto createLevelCellData = [&levelCartDims]() {
463 return Opm::data::CellData{
464 Opm::UnitSystem::measure::transmissibility,
465 std::vector<double>(levelCartDims[0] * levelCartDims[1] * levelCartDims[2], 0.0),
466 Opm::data::TargetType::INIT
471 levelTrans.emplace(
"TRANX", createLevelCellData());
472 levelTrans.emplace(
"TRANY", createLevelCellData());
473 levelTrans.emplace(
"TRANZ", createLevelCellData());
476template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
477template<
typename LevelIndicesFunction,
typename OriginIndicesFunction>
479EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
480computeTrans_(
const std::vector<std::unordered_map<int,int>>& levelCartToLevelCompressed,
481 const std::function<
unsigned int(
unsigned int)>& map,
482 const LevelIndicesFunction& computeLevelIndices,
483 const std::function<
int(
int,
int)>& computeLevelCartIdx,
484 const std::function<std::array<int,3>(
int)>& computeLevelCartDims,
485 const OriginIndicesFunction& computeOriginIndices)
const
488 outputTrans_ = std::make_unique<std::vector<data::Solution>>(std::vector<data::Solution>{});
491 using GlobalGridView =
typename EquilGrid::LeafGridView;
492 using GlobElementMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GlobalGridView>;
493 const GlobalGridView& globalGridView = this->equilGrid_->leafGridView();
494 const GlobElementMapper globalElemMapper { globalGridView, Dune::mcmgElementLayout() };
497 int maxLevel = this->equilGrid_->maxLevel();
499 outputTrans_->resize(maxLevel+1);
501 for (
int level = 0; level <= maxLevel; ++level) {
502 allocateLevelTrans_(computeLevelCartDims(level), this->outputTrans_->at(level));
505 for (
const auto& elem : elements(globalGridView)) {
506 for (
const auto& is : intersections(globalGridView, elem)) {
510 if ( is.inside().level() != is.outside().level() )
514 unsigned c1 = globalElemMapper.index(is.inside());
515 unsigned c2 = globalElemMapper.index(is.outside());
520 int level = is.inside().level();
523 const auto& [levelInIdx, levelOutIdx] = computeLevelIndices(is, c1, c2);
525 const int levelCartIdxIn = computeLevelCartIdx(levelInIdx, level);
526 const int levelCartIdxOut = computeLevelCartIdx(levelOutIdx, level);
532 const auto [originInIdx, originOutIdx] = computeOriginIndices(is, c1, c2);
534 const auto originCartIdxIn = computeLevelCartIdx(originInIdx, 0);
535 const auto originCartIdxOut = computeLevelCartIdx(originOutIdx, 0);
538 if (isNumAquCell_(originCartIdxIn) || isNumAquCell_(originCartIdxOut)) {
548 const auto minLevelCartIdx = std::min(levelCartIdxIn, levelCartIdxOut);
549 const auto maxLevelCartIdx = std::max(levelCartIdxIn, levelCartIdxOut);
551 const auto& levelCartDims = computeLevelCartDims(level);
559 if (maxLevelCartIdx - minLevelCartIdx == 1 && levelCartDims[0] > 1 ) {
560 outputTrans_->at(level).at(
"TRANX").template data<double>()[minLevelCartIdx] = globalTrans().transmissibility(c1, c2);
564 if (maxLevelCartIdx - minLevelCartIdx == levelCartDims[0] && levelCartDims[1] > 1) {
565 outputTrans_->at(level).at(
"TRANY").template data<double>()[minLevelCartIdx] = globalTrans().transmissibility(c1, c2);
569 if ( maxLevelCartIdx - minLevelCartIdx == levelCartDims[0]*levelCartDims[1] ||
570 directVerticalNeighbors(levelCartDims,
571 levelCartToLevelCompressed[level],
574 outputTrans_->at(level).at(
"TRANZ").template data<double>()[minLevelCartIdx] = globalTrans().transmissibility(c1, c2);
580template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
582EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
583isCartesianNeighbour_(
const std::array<int,3>& levelCartDims,
584 const std::size_t levelCartIdx1,
585 const std::size_t levelCartIdx2)
const
587 const int diff = levelCartIdx2 - levelCartIdx1;
590 || (diff == levelCartDims[0])
591 || (diff == (levelCartDims[0] * levelCartDims[1]));
594template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
596EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
597isDirectNeighbours_(
const std::unordered_map<int,int>& levelCartesianToActive,
598 const std::array<int,3>& levelCartDims,
599 const std::size_t levelCartIdx1,
600 const std::size_t levelCartIdx2)
const
602 return isCartesianNeighbour_(levelCartDims, levelCartIdx1, levelCartIdx2)
603 || directVerticalNeighbors(levelCartDims, levelCartesianToActive, levelCartIdx1, levelCartIdx2);
606template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
608EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
609activeCell_(
const std::unordered_map<int,int>& levelCartToLevelCompressed,
610 const std::size_t levelCartIdx)
const
612 auto pos = levelCartToLevelCompressed.find(levelCartIdx);
613 return (pos == levelCartToLevelCompressed.end()) ? -1 : pos->second;
616template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
618EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
619allocateAllNncs_(
int maxLevel)
const
621 this->outputNnc_.resize(maxLevel+1);
629 this->outputNncGlobalLocal_.resize(maxLevel);
637 this->outputAmalgamatedNnc_.resize(maxLevel-1);
638 for (
int i = 0; i < maxLevel-1; ++i) {
639 this->outputAmalgamatedNnc_[i].resize(maxLevel-1-i);
644template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
645template<
typename LevelIndicesFunction,
typename OriginIndicesFunction>
646std::vector<std::vector<NNCdata>>
647EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
648exportNncStructure_(
const std::vector<std::unordered_map<int,int>>& levelCartToLevelCompressed,
649 const std::function<
unsigned int(
unsigned int)>& map,
650 const LevelIndicesFunction& computeLevelIndices,
651 const std::function<
int(
int,
int)>& computeLevelCartIdx,
652 const std::function<std::array<int,3>(
int)>& computeLevelCartDims,
653 const OriginIndicesFunction& computeOriginIndices)
const
655 const auto& nncData = this->eclState_.getInputNNC().input();
656 const auto& nncEdit = this->eclState_.getInputNNC().edit();
657 const auto& nncEditr = this->eclState_.getInputNNC().editr();
658 const auto& unitSystem = this->eclState_.getDeckUnitSystem();
659 const auto& transMult = this->eclState_.getTransMult();
662 const auto& equilCartMapper = *equilCartMapper_;
664 const auto& level0CartDims = equilCartMapper.cartesianDimensions();
666 int maxLevel = this->equilGrid_->maxLevel();
667 allocateAllNncs_(maxLevel);
669 using GlobalGridView =
typename EquilGrid::LeafGridView;
670 using GlobElementMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GlobalGridView>;
671 const GlobalGridView& globalGridView = this->equilGrid_->leafGridView();
672 const GlobElementMapper globalElemMapper { globalGridView, Dune::mcmgElementLayout() };
674 for (
const auto& elem : elements(globalGridView)) {
675 for (
const auto& is : intersections(globalGridView, elem)) {
680 unsigned c1 = globalElemMapper.index(is.inside());
681 unsigned c2 = globalElemMapper.index(is.outside());
686 if ( is.inside().level() != is.outside().level() ) {
688 const auto& [levelInIdx, levelOutIdx] = computeLevelIndices(is, c1, c2);
690 const int levelIn = is.inside().level();
691 const int levelOut = is.outside().level();
693 auto levelCartIdxIn = computeLevelCartIdx(levelInIdx, levelIn);
694 auto levelCartIdxOut = computeLevelCartIdx(levelOutIdx, levelOut);
697 std::pair<int,int> smallerPair = {levelIn, levelCartIdxIn},
698 largerPair = {levelOut, levelCartIdxOut};
699 if (smallerPair.first > largerPair.first) {
700 std::swap(smallerPair, largerPair);
703 const auto& [smallerLevel, smallerLevelCartIdx] = smallerPair;
704 const auto& [largerLevel, largerLevelCartIdx] = largerPair;
706 auto t = this->globalTrans().transmissibility(c1, c2);
713 const auto tt = unitSystem
714 .from_si(UnitSystem::measure::transmissibility, t);
716 if (std::isnormal(tt) && (tt > 1.0e-12)) {
718 if (smallerLevel == 0) {
719 this->outputNncGlobalLocal_[largerLevel-1].emplace_back(smallerLevelCartIdx, largerLevelCartIdx, t);
722 assert(smallerLevel >= 1);
723 this->outputAmalgamatedNnc_[smallerLevel-1][largerLevel-smallerLevel-1].emplace_back(smallerLevelCartIdx, largerLevelCartIdx, t);
729 assert(is.inside().level() == is.outside().level());
730 const int level = is.inside().level();
736 const auto [originInIdx, originOutIdx] = computeOriginIndices(is, c1, c2);
738 const std::size_t originCartIdxIn = computeLevelCartIdx(originInIdx, 0);
739 const std::size_t originCartIdxOut = computeLevelCartIdx(originOutIdx, 0);
742 const auto& [levelInIdx, levelOutIdx] = computeLevelIndices(is, c1, c2);
744 auto levelCartIdxIn = computeLevelCartIdx(levelInIdx, level);
745 auto levelCartIdxOut = computeLevelCartIdx(levelOutIdx, level);
747 if ( levelCartIdxOut < levelCartIdxIn )
748 std::swap(levelCartIdxIn, levelCartIdxOut);
756 const auto& levelCartDims = computeLevelCartDims(level);
759 assert(!isNumAquConn_(originCartIdxIn, originCartIdxOut) || level == 0);
761 if (isNumAquConn_(originCartIdxIn, originCartIdxOut) ||
762 ! isDirectNeighbours_(levelCartToLevelCompressed[level],
764 levelCartIdxIn, levelCartIdxOut)) {
767 auto t = this->globalTrans().transmissibility(c1, c2);
770 auto candidate = std::lower_bound(nncData.begin(), nncData.end(),
771 NNCdata { originCartIdxIn, originCartIdxOut, 0.0 });
772 const auto transMlt = transMult.getRegionMultiplierNNC(originCartIdxIn, originCartIdxOut);
773 bool foundNncEditr =
false;
775 while ((candidate != nncData.end()) &&
776 (candidate->cell1 == originCartIdxIn) &&
777 (candidate->cell2 == originCartIdxOut))
779 auto trans = candidate->trans;
781 if (! nncEditr.empty()) {
782 auto it = std::lower_bound(nncEditr.begin(), nncEditr.end(),
783 NNCdata { originCartIdxIn, originCartIdxOut, 0.0 });
784 foundNncEditr = it != nncEditr.end() && it->cell1 == originCartIdxIn && it->cell2 == originCartIdxOut;
790 if (! nncEdit.empty()) {
791 auto it = std::lower_bound(nncEdit.begin(), nncEdit.end(),
792 NNCdata { originCartIdxIn, originCartIdxOut, 0.0 });
793 if (it != nncEdit.end() && it->cell1 == originCartIdxIn && it->cell2 == originCartIdxOut) {
811 const auto tt = unitSystem
812 .from_si(UnitSystem::measure::transmissibility, t);
814 if (std::isnormal(tt) && (tt > 1.0e-12)) {
815 this->outputNnc_[level].emplace_back(levelCartIdxIn, levelCartIdxOut, t);
823 std::vector<NNCdata> inputedNnc{};
824 const auto generatedNnc = outputNnc_[0];
828 for (
const auto& entry : nncData) {
837 assert (entry.cell2 >= entry.cell1);
839 if (! isCartesianNeighbour_(level0CartDims, entry.cell1, entry.cell2) ||
840 isNumAquConn_(entry.cell1, entry.cell2))
842 bool foundNncEdit =
false;
843 auto trans = entry.trans;
844 if (! nncEdit.empty()) {
845 auto it = std::lower_bound(nncEdit.begin(), nncEdit.end(),
846 NNCdata {entry.cell1, entry.cell2, 0.0 });
847 if (it != nncEdit.end() && it->cell1 == entry.cell1 && it->cell2 == entry.cell2) {
852 if (! foundNncEdit) {
856 const auto c1 = activeCell_(levelCartToLevelCompressed[0], entry.cell1);
857 const auto c2 = activeCell_(levelCartToLevelCompressed[0], entry.cell2);
859 if ((c1 < 0) || (c2 < 0)) {
865 trans = this->globalTrans().transmissibility(c1, c2);
867 if (! generatedNnc.empty()) {
868 for (
const auto& generated : generatedNnc) {
869 if (entry.cell1 == generated.cell1 && entry.cell2 == generated.cell2) {
870 trans -= generated.trans;
876 const auto tt = unitSystem
877 .from_si(UnitSystem::measure::transmissibility, trans);
882 if (std::isnormal(tt) && ! (tt < 1.0e-6)) {
883 inputedNnc.emplace_back(entry.cell1, entry.cell2, trans);
888 this->outputNnc_[0].insert(this->outputNnc_[0].begin(), inputedNnc.begin(), inputedNnc.end());
889 return this->outputNnc_;
892template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
895 const std::optional<int> timeStepNum,
896 const bool isSubStep,
897 data::Solution&& localCellData,
898 data::Wells&& localWellData,
899 data::GroupAndNetworkValues&& localGroupAndNetworkData,
900 data::Aquifers&& localAquiferData,
901 WellTestState&& localWTestState,
902 const Action::State& actionState,
903 const UDQState& udqState,
904 const SummaryState& summaryState,
905 const std::vector<Scalar>& thresholdPressure,
908 bool doublePrecision,
914 const auto isParallel = this->collectOnIORank_.isParallel();
915 const bool needsReordering = this->collectOnIORank_.doesNeedReordering();
917 RestartValue restartValue {
918 (isParallel || needsReordering)
919 ? this->collectOnIORank_.globalCellData()
920 : std::move(localCellData),
922 isParallel ? this->collectOnIORank_.globalWellData()
923 : std::move(localWellData),
925 isParallel ? this->collectOnIORank_.globalGroupAndNetworkData()
926 : std::move(localGroupAndNetworkData),
928 isParallel ? this->collectOnIORank_.globalAquiferData()
929 : std::move(localAquiferData)
932 if (eclState_.getSimulationConfig().useThresholdPressure()) {
933 restartValue.addExtra(
"THRESHPR", UnitSystem::measure::pressure,
939 restartValue.addExtra(
"OPMEXTRA", std::vector<double>(1, nextStepSize));
944 const auto flowsn_global = isParallel ? this->collectOnIORank_.globalFlowsn() : std::move(flowsn);
945 for (
const auto& flows : flowsn_global) {
946 if (flows.name.empty())
948 if (flows.name ==
"FLOGASN+") {
949 restartValue.addExtra(flows.name, UnitSystem::measure::gas_surface_rate, flows.values);
951 restartValue.addExtra(flows.name, UnitSystem::measure::liquid_surface_rate, flows.values);
956 const auto floresn_global = isParallel ? this->collectOnIORank_.globalFloresn() : std::move(floresn);
957 for (
const auto& flores : floresn_global) {
958 if (flores.name.empty()) {
961 restartValue.addExtra(flores.name, UnitSystem::measure::rate, flores.values);
965 std::vector<Opm::RestartValue> restartValues{};
967 if ( !isParallel && !needsReordering && (this->eclState_.getLgrs().size()>0) && (this->grid_.maxLevel()>0) ) {
972 Opm::Lgr::extractRestartValueLevelGrids<Grid>(this->grid_, restartValue, restartValues);
975 restartValues.reserve(1);
976 restartValues.push_back(std::move(restartValue));
982 this->taskletRunner_->barrier();
985 if (this->taskletRunner_->failure()) {
986 throw std::runtime_error(
"Failure in the TaskletRunner while writing output.");
990 auto eclWriteTasklet = std::make_shared<EclWriteTasklet>(
992 isParallel ? this->collectOnIORank_.globalWellTestState() : std::move(localWTestState),
993 summaryState, udqState, *this->eclIO_,
994 reportStepNum, timeStepNum, isSubStep, curTime, std::move(restartValues), doublePrecision);
997 this->taskletRunner_->dispatch(std::move(eclWriteTasklet));
1000template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
1003 const Scalar curTime,
1004 const data::Wells& localWellData,
1005 const data::WellBlockAveragePressures& localWBPData,
1006 const data::GroupAndNetworkValues& localGroupAndNetworkData,
1007 const std::map<int,data::AquiferData>& localAquiferData,
1008 const std::map<std::pair<std::string, int>,
double>& blockData,
1009 const std::map<std::string, double>& miscSummaryData,
1010 const std::map<std::string, std::vector<double>>& regionData,
1011 const Inplace& inplace,
1012 const Inplace* initialInPlace,
1014 SummaryState& summaryState,
1017 if (collectOnIORank_.isIORank()) {
1018 const auto& summary = eclIO_->summary();
1020 const auto& wellData = this->collectOnIORank_.isParallel()
1021 ? this->collectOnIORank_.globalWellData()
1024 const auto& wbpData = this->collectOnIORank_.isParallel()
1025 ? this->collectOnIORank_.globalWBPData()
1028 const auto& groupAndNetworkData = this->collectOnIORank_.isParallel()
1029 ? this->collectOnIORank_.globalGroupAndNetworkData()
1030 : localGroupAndNetworkData;
1032 const auto& aquiferData = this->collectOnIORank_.isParallel()
1033 ? this->collectOnIORank_.globalAquiferData()
1036 const auto interreg_flows = getInterRegFlowsAsMap(interRegFlows);
1038 const auto values = out::Summary::DynamicSimulatorState {
1041 &groupAndNetworkData,
1053 summary.eval(reportStepNum, curTime, values, summaryState);
1059 const auto udq_step = reportStepNum - 1;
1061 this->schedule_[udq_step].udq()
1063 this->schedule_.wellMatcher(udq_step),
1064 this->schedule_[udq_step].group_order(),
1065 this->schedule_.segmentMatcherFactory(udq_step),
1066 [es = std::cref(this->eclState_)]() {
1067 return std::make_unique<RegionSetMatcher>
1068 (es.get().fipRegionStatistics());
1070 summaryState, udqState);
1074 if (collectOnIORank_.isParallel()) {
1076 ser.
append(summaryState);
1081template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
1086 assert (globalTrans_);
1087 return *globalTrans_;
Definition: CollectDataOnIORank.hpp:49
bool isIORank() const
Definition: CollectDataOnIORank.hpp:126
Definition: EclGenericWriter.hpp:67
std::vector< std::vector< NNCdata > > outputNnc_
Definition: EclGenericWriter.hpp:174
const EclipseState & eclState_
Definition: EclGenericWriter.hpp:158
void doWriteOutput(const int reportStepNum, const std::optional< int > timeStepNum, const bool isSubStep, 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:894
CollectDataOnIORankType collectOnIORank_
Definition: EclGenericWriter.hpp:154
void extractOutputTransAndNNC(const std::function< unsigned int(unsigned int)> &map)
Definition: EclGenericWriter_impl.hpp:286
std::unique_ptr< TaskletRunner > taskletRunner_
Definition: EclGenericWriter.hpp:160
void writeInit()
Definition: EclGenericWriter_impl.hpp:268
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:214
const EclipseIO & eclIO() const
Definition: EclGenericWriter_impl.hpp:260
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)
Definition: EclGenericWriter_impl.hpp:1002
const TransmissibilityType & globalTrans() const
Definition: EclGenericWriter_impl.hpp:1084
std::unique_ptr< EclipseIO > eclIO_
Definition: EclGenericWriter.hpp:159
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:50
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