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>
54#include <opm/simulators/utils/MPISerializer.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 bool forcedSimulationFinished_;
156 explicit EclWriteTasklet(
const Opm::Action::State& actionState,
157 const Opm::WellTestState& wtestState,
158 const Opm::SummaryState& summaryState,
159 const Opm::UDQState& udqState,
160 Opm::EclipseIO& eclIO,
162 std::optional<int> timeStepNum,
164 double secondsElapsed,
165 std::vector<Opm::RestartValue> restartValue,
166 bool writeDoublePrecision,
167 bool forcedSimulationFinished)
168 : actionState_(actionState)
169 , wtestState_(wtestState)
170 , summaryState_(summaryState)
171 , udqState_(udqState)
173 , reportStepNum_(reportStepNum)
174 , timeStepNum_(timeStepNum)
175 , isSubStep_(isSubStep)
176 , secondsElapsed_(secondsElapsed)
177 , restartValue_(std::move(restartValue))
178 , writeDoublePrecision_(writeDoublePrecision)
179 , forcedSimulationFinished_(forcedSimulationFinished)
185 if (this->restartValue_.size() == 1) {
186 this->eclIO_.writeTimeStep(this->actionState_,
190 this->reportStepNum_,
192 this->secondsElapsed_,
193 std::move(this->restartValue_.back()),
194 this->writeDoublePrecision_,
196 forcedSimulationFinished_);
199 this->eclIO_.writeTimeStep(this->actionState_,
203 this->reportStepNum_,
205 this->secondsElapsed_,
206 std::move(this->restartValue_),
207 this->writeDoublePrecision_,
209 forcedSimulationFinished_);
218template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
221 const EclipseState& eclState,
222 const SummaryConfig& summaryConfig,
224 const EquilGrid* equilGrid,
225 const GridView& gridView,
226 const Dune::CartesianIndexMapper<Grid>& cartMapper,
227 const Dune::CartesianIndexMapper<EquilGrid>* equilCartMapper,
228 bool enableAsyncOutput,
230 : collectOnIORank_(grid,
235 summaryConfig.fip_regions_interreg_flow())
237 , gridView_ (gridView)
238 , schedule_ (schedule)
239 , eclState_ (eclState)
240 , cartMapper_ (cartMapper)
241 , equilCartMapper_(equilCartMapper)
242 , equilGrid_ (equilGrid)
245 outputNnc_.resize(1);
247 if (this->collectOnIORank_.isIORank()) {
248 this->eclIO_ = std::make_unique<EclipseIO>
250 UgGridHelpers::createEclipseGrid(*equilGrid, eclState_.getInputGrid()),
251 this->schedule_, summaryConfig,
"", enableEsmry);
256 int numWorkerThreads = 0;
257 if (enableAsyncOutput && collectOnIORank_.isIORank()) {
258 numWorkerThreads = 1;
261 this->taskletRunner_.reset(
new TaskletRunner(numWorkerThreads));
264template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
265const EclipseIO& EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
272template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
273void EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
276 if (collectOnIORank_.isIORank()) {
277 std::map<std::string, std::vector<int>> integerVectors;
278 if (collectOnIORank_.isParallel()) {
279 integerVectors.emplace(
"MPI_RANK", collectOnIORank_.globalRanks());
282 eclIO_->writeInitial(*this->outputTrans_,
284 this->outputNnc_.front());
285 this->outputTrans_.reset();
289template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
291EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
292extractOutputTransAndNNC(
const std::function<
unsigned int(
unsigned int)>& map)
294 if (collectOnIORank_.isIORank()) {
295 constexpr bool equilGridIsCpGrid = std::is_same_v<EquilGrid, Dune::CpGrid>;
297 const auto levelCartMapp = this->createLevelCartMapp_<equilGridIsCpGrid>();
298 const auto levelCartToLevelCompressed = this->createCartesianToActiveMaps_<equilGridIsCpGrid>(levelCartMapp);
299 auto computeLevelIndices = this->computeLevelIndices_<equilGridIsCpGrid>();
300 auto computeLevelCartIdx = this->computeLevelCartIdx_<equilGridIsCpGrid>(levelCartMapp, *(this->equilCartMapper_));
301 auto computeLevelCartDimensions = this->computeLevelCartDimensions_<equilGridIsCpGrid>(levelCartMapp, *(this->equilCartMapper_));
302 auto computeOriginIndices = this->computeOriginIndices_<equilGridIsCpGrid>();
304 computeTrans_(levelCartToLevelCompressed, map, computeLevelIndices,
305 computeLevelCartIdx, computeLevelCartDimensions, computeOriginIndices);
306 exportNncStructure_(levelCartToLevelCompressed, map, computeLevelIndices, computeLevelCartIdx,
307 computeLevelCartDimensions, computeOriginIndices);
311 if (collectOnIORank_.isParallel()) {
312 const auto& comm = grid_.comm();
313 Parallel::MpiSerializer ser(comm);
314 ser.broadcast(Parallel::RootRank{0}, outputNnc_);
319template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
321EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
322isNumAquCell_(
const std::size_t cartIdx)
const
324 const auto& numAquCell = this->eclState_.aquifer().hasNumericalAquifer()
325 ? this->eclState_.aquifer().numericalAquifers().allAquiferCellIds()
326 : std::vector<std::size_t>{};
328 return std::ranges::binary_search(numAquCell.begin(), numAquCell.end(), cartIdx);
331template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
333EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
334isNumAquConn_(
const std::size_t cartIdx1,
335 const std::size_t cartIdx2)
const
337 return isNumAquCell_(cartIdx1) || isNumAquCell_(cartIdx2);
340template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
341template<
bool equilGr
idIsCpGr
id>
343EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
344createLevelCartMapp_()
const
346 if constexpr (equilGridIsCpGrid) {
352template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
353template<
bool equilGr
idIsCpGr
id>
354std::vector<std::unordered_map<int,int>>
355EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
358 if constexpr (equilGridIsCpGrid) {
359 if (this->equilGrid_->maxLevel()) {
360 return Opm::Lgr::levelCartesianToLevelCompressedMaps(*this->equilGrid_, levelCartMapp); }
362 return std::vector<std::unordered_map<int,int>>{ cartesianToCompressed(equilGrid_->size(0), UgGridHelpers::globalCell(*equilGrid_)) };
365 return std::vector<std::unordered_map<int,int>>{ cartesianToCompressed(equilGrid_->size(0), UgGridHelpers::globalCell(*equilGrid_)) };
368template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
369template<
bool equilGr
idIsCpGr
id>
370std::function<std::array<int,3>(
int)>
371EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
375 if constexpr (equilGridIsCpGrid) {
376 return [&](
int level)
378 return levelCartMapp.cartesianDimensions(level);
382 return [&]([[maybe_unused]]
int level)
385 return equilCartMapp.cartesianDimensions();
390template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
391template<
bool equilGr
idIsCpGr
id>
392std::function<int(
int,
int)>
393EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
397 if constexpr (equilGridIsCpGrid) {
398 return [&](
int levelCompressedIdx,
401 return levelCartMapp.cartesianIndex(levelCompressedIdx, level);
405 return [&](
int levelCompressedIdx,
406 [[maybe_unused]]
int level)
409 return equilCartMapp.cartesianIndex(levelCompressedIdx);
414template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
415template <
bool equilGr
idIsCpGr
id>
417EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
418computeLevelIndices_()
const
420 if constexpr (equilGridIsCpGrid) {
421 return [](
const auto& intersection,
425 return std::pair{intersection.inside().getLevelElem().index(), intersection.outside().getLevelElem().index()};
429 return [](
const auto&,
430 const auto& intersectionInsideLeafIdx,
431 const auto& intersectionOutsideLeafIdx)
433 return std::pair{intersectionInsideLeafIdx, intersectionOutsideLeafIdx};
438template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
439template <
bool equilGr
idIsCpGr
id>
441EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
442computeOriginIndices_()
const
444 if constexpr (equilGridIsCpGrid) {
445 return [](
const auto& intersection,
449 return std::pair{intersection.inside().getOrigin().index(), intersection.outside().getOrigin().index()};
453 return [](
const auto&,
454 const auto& intersectionInsideLeafIdx,
455 const auto& intersectionOutsideLeafIdx)
457 return std::pair{intersectionInsideLeafIdx, intersectionOutsideLeafIdx};
462template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
464EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
465allocateLevelTrans_(
const std::array<int,3>& levelCartDims,
466 data::Solution& levelTrans)
const
468 auto createLevelCellData = [&levelCartDims]() {
469 return Opm::data::CellData{
470 Opm::UnitSystem::measure::transmissibility,
471 std::vector<double>(levelCartDims[0] * levelCartDims[1] * levelCartDims[2], 0.0),
472 Opm::data::TargetType::INIT
477 levelTrans.emplace(
"TRANX", createLevelCellData());
478 levelTrans.emplace(
"TRANY", createLevelCellData());
479 levelTrans.emplace(
"TRANZ", createLevelCellData());
482template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
483template<
typename LevelIndicesFunction,
typename OriginIndicesFunction>
485EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
486computeTrans_(
const std::vector<std::unordered_map<int,int>>& levelCartToLevelCompressed,
487 const std::function<
unsigned int(
unsigned int)>& map,
488 const LevelIndicesFunction& computeLevelIndices,
489 const std::function<
int(
int,
int)>& computeLevelCartIdx,
490 const std::function<std::array<int,3>(
int)>& computeLevelCartDims,
491 const OriginIndicesFunction& computeOriginIndices)
const
494 outputTrans_ = std::make_unique<std::vector<data::Solution>>(std::vector<data::Solution>{});
497 using GlobalGridView =
typename EquilGrid::LeafGridView;
498 using GlobElementMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GlobalGridView>;
499 const GlobalGridView& globalGridView = this->equilGrid_->leafGridView();
500 const GlobElementMapper globalElemMapper { globalGridView, Dune::mcmgElementLayout() };
503 int maxLevel = this->equilGrid_->maxLevel();
505 outputTrans_->resize(maxLevel+1);
507 for (
int level = 0; level <= maxLevel; ++level) {
508 allocateLevelTrans_(computeLevelCartDims(level), this->outputTrans_->at(level));
511 for (
const auto& elem : elements(globalGridView)) {
512 for (
const auto& is : intersections(globalGridView, elem)) {
516 if ( is.inside().level() != is.outside().level() )
520 unsigned c1 = globalElemMapper.index(is.inside());
521 unsigned c2 = globalElemMapper.index(is.outside());
526 int level = is.inside().level();
529 const auto& [levelInIdx, levelOutIdx] = computeLevelIndices(is, c1, c2);
531 const int levelCartIdxIn = computeLevelCartIdx(levelInIdx, level);
532 const int levelCartIdxOut = computeLevelCartIdx(levelOutIdx, level);
538 const auto [originInIdx, originOutIdx] = computeOriginIndices(is, c1, c2);
540 const auto originCartIdxIn = computeLevelCartIdx(originInIdx, 0);
541 const auto originCartIdxOut = computeLevelCartIdx(originOutIdx, 0);
544 if (isNumAquCell_(originCartIdxIn) || isNumAquCell_(originCartIdxOut)) {
554 const auto minLevelCartIdx = std::min(levelCartIdxIn, levelCartIdxOut);
555 const auto maxLevelCartIdx = std::max(levelCartIdxIn, levelCartIdxOut);
557 const auto& levelCartDims = computeLevelCartDims(level);
565 if (maxLevelCartIdx - minLevelCartIdx == 1 && levelCartDims[0] > 1 ) {
566 outputTrans_->at(level).at(
"TRANX").template data<double>()[minLevelCartIdx] = globalTrans().transmissibility(c1, c2);
570 if (maxLevelCartIdx - minLevelCartIdx == levelCartDims[0] && levelCartDims[1] > 1) {
571 outputTrans_->at(level).at(
"TRANY").template data<double>()[minLevelCartIdx] = globalTrans().transmissibility(c1, c2);
575 if ( maxLevelCartIdx - minLevelCartIdx == levelCartDims[0]*levelCartDims[1] ||
576 directVerticalNeighbors(levelCartDims,
577 levelCartToLevelCompressed[level],
580 outputTrans_->at(level).at(
"TRANZ").template data<double>()[minLevelCartIdx] = globalTrans().transmissibility(c1, c2);
586template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
588EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
589isCartesianNeighbour_(
const std::array<int,3>& levelCartDims,
590 const std::size_t levelCartIdx1,
591 const std::size_t levelCartIdx2)
const
593 const int diff = levelCartIdx2 - levelCartIdx1;
596 || (diff == levelCartDims[0])
597 || (diff == (levelCartDims[0] * levelCartDims[1]));
600template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
602EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
603isDirectNeighbours_(
const std::unordered_map<int,int>& levelCartesianToActive,
604 const std::array<int,3>& levelCartDims,
605 const std::size_t levelCartIdx1,
606 const std::size_t levelCartIdx2)
const
608 return isCartesianNeighbour_(levelCartDims, levelCartIdx1, levelCartIdx2)
609 || directVerticalNeighbors(levelCartDims, levelCartesianToActive, levelCartIdx1, levelCartIdx2);
612template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
614EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
615activeCell_(
const std::unordered_map<int,int>& levelCartToLevelCompressed,
616 const std::size_t levelCartIdx)
const
618 auto pos = levelCartToLevelCompressed.find(levelCartIdx);
619 return (pos == levelCartToLevelCompressed.end()) ? -1 : pos->second;
622template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
624EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
625allocateAllNncs_(
int maxLevel)
const
627 this->outputNnc_.resize(maxLevel+1);
635 this->outputNncGlobalLocal_.resize(maxLevel);
643 this->outputAmalgamatedNnc_.resize(maxLevel-1);
644 for (
int i = 0; i < maxLevel-1; ++i) {
645 this->outputAmalgamatedNnc_[i].resize(maxLevel-1-i);
650template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
651template<
typename LevelIndicesFunction,
typename OriginIndicesFunction>
652std::vector<std::vector<NNCdata>>
653EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
654exportNncStructure_(
const std::vector<std::unordered_map<int,int>>& levelCartToLevelCompressed,
655 const std::function<
unsigned int(
unsigned int)>& map,
656 const LevelIndicesFunction& computeLevelIndices,
657 const std::function<
int(
int,
int)>& computeLevelCartIdx,
658 const std::function<std::array<int,3>(
int)>& computeLevelCartDims,
659 const OriginIndicesFunction& computeOriginIndices)
const
661 const auto& nncData = this->eclState_.getInputNNC().input();
662 const auto& nncEdit = this->eclState_.getInputNNC().edit();
663 const auto& nncEditr = this->eclState_.getInputNNC().editr();
664 const auto& unitSystem = this->eclState_.getDeckUnitSystem();
665 const auto& transMult = this->eclState_.getTransMult();
668 const auto& equilCartMapper = *equilCartMapper_;
670 const auto& level0CartDims = equilCartMapper.cartesianDimensions();
672 int maxLevel = this->equilGrid_->maxLevel();
673 allocateAllNncs_(maxLevel);
675 using GlobalGridView =
typename EquilGrid::LeafGridView;
676 using GlobElementMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GlobalGridView>;
677 const GlobalGridView& globalGridView = this->equilGrid_->leafGridView();
678 const GlobElementMapper globalElemMapper { globalGridView, Dune::mcmgElementLayout() };
680 for (
const auto& elem : elements(globalGridView)) {
681 for (
const auto& is : intersections(globalGridView, elem)) {
686 unsigned c1 = globalElemMapper.index(is.inside());
687 unsigned c2 = globalElemMapper.index(is.outside());
692 if ( is.inside().level() != is.outside().level() ) {
694 const auto& [levelInIdx, levelOutIdx] = computeLevelIndices(is, c1, c2);
696 const int levelIn = is.inside().level();
697 const int levelOut = is.outside().level();
699 auto levelCartIdxIn = computeLevelCartIdx(levelInIdx, levelIn);
700 auto levelCartIdxOut = computeLevelCartIdx(levelOutIdx, levelOut);
703 std::pair<int,int> smallerPair = {levelIn, levelCartIdxIn},
704 largerPair = {levelOut, levelCartIdxOut};
705 if (smallerPair.first > largerPair.first) {
706 std::swap(smallerPair, largerPair);
709 const auto& [smallerLevel, smallerLevelCartIdx] = smallerPair;
710 const auto& [largerLevel, largerLevelCartIdx] = largerPair;
712 auto t = this->globalTrans().transmissibility(c1, c2);
719 const auto tt = unitSystem
720 .from_si(UnitSystem::measure::transmissibility, t);
722 if (std::isnormal(tt) && (tt > 1.0e-12)) {
724 if (smallerLevel == 0) {
725 this->outputNncGlobalLocal_[largerLevel-1].emplace_back(smallerLevelCartIdx, largerLevelCartIdx, t);
728 assert(smallerLevel >= 1);
729 this->outputAmalgamatedNnc_[smallerLevel-1][largerLevel-smallerLevel-1].emplace_back(smallerLevelCartIdx, largerLevelCartIdx, t);
735 assert(is.inside().level() == is.outside().level());
736 const int level = is.inside().level();
742 const auto [originInIdx, originOutIdx] = computeOriginIndices(is, c1, c2);
744 const std::size_t originCartIdxIn = computeLevelCartIdx(originInIdx, 0);
745 const std::size_t originCartIdxOut = computeLevelCartIdx(originOutIdx, 0);
748 const auto& [levelInIdx, levelOutIdx] = computeLevelIndices(is, c1, c2);
750 auto levelCartIdxIn = computeLevelCartIdx(levelInIdx, level);
751 auto levelCartIdxOut = computeLevelCartIdx(levelOutIdx, level);
753 if ( levelCartIdxOut < levelCartIdxIn )
754 std::swap(levelCartIdxIn, levelCartIdxOut);
762 const auto& levelCartDims = computeLevelCartDims(level);
765 assert(!isNumAquConn_(originCartIdxIn, originCartIdxOut) || level == 0);
767 if (isNumAquConn_(originCartIdxIn, originCartIdxOut) ||
768 ! isDirectNeighbours_(levelCartToLevelCompressed[level],
770 levelCartIdxIn, levelCartIdxOut)) {
773 auto t = this->globalTrans().transmissibility(c1, c2);
776 auto candidate = std::lower_bound(nncData.begin(), nncData.end(),
777 NNCdata { originCartIdxIn, originCartIdxOut, 0.0 });
778 const auto transMlt = transMult.getRegionMultiplierNNC(originCartIdxIn, originCartIdxOut);
779 bool foundNncEditr =
false;
781 while ((candidate != nncData.end()) &&
782 (candidate->cell1 == originCartIdxIn) &&
783 (candidate->cell2 == originCartIdxOut))
785 auto trans = candidate->trans;
787 if (! nncEditr.empty()) {
788 auto it = std::lower_bound(nncEditr.begin(), nncEditr.end(),
789 NNCdata { originCartIdxIn, originCartIdxOut, 0.0 });
790 foundNncEditr = it != nncEditr.end() && it->cell1 == originCartIdxIn && it->cell2 == originCartIdxOut;
796 if (! nncEdit.empty()) {
797 auto it = std::lower_bound(nncEdit.begin(), nncEdit.end(),
798 NNCdata { originCartIdxIn, originCartIdxOut, 0.0 });
799 if (it != nncEdit.end() && it->cell1 == originCartIdxIn && it->cell2 == originCartIdxOut) {
817 const auto tt = unitSystem
818 .from_si(UnitSystem::measure::transmissibility, t);
820 if (std::isnormal(tt) && (tt > 1.0e-12)) {
821 this->outputNnc_[level].emplace_back(levelCartIdxIn, levelCartIdxOut, t);
829 std::vector<NNCdata> inputedNnc{};
830 const auto generatedNnc = outputNnc_[0];
834 for (
const auto& entry : nncData) {
843 assert (entry.cell2 >= entry.cell1);
845 if (! isCartesianNeighbour_(level0CartDims, entry.cell1, entry.cell2) ||
846 isNumAquConn_(entry.cell1, entry.cell2))
848 bool foundNncEdit =
false;
849 auto trans = entry.trans;
850 if (! nncEdit.empty()) {
851 auto it = std::lower_bound(nncEdit.begin(), nncEdit.end(),
852 NNCdata {entry.cell1, entry.cell2, 0.0 });
853 if (it != nncEdit.end() && it->cell1 == entry.cell1 && it->cell2 == entry.cell2) {
858 if (! foundNncEdit) {
862 const auto c1 = activeCell_(levelCartToLevelCompressed[0], entry.cell1);
863 const auto c2 = activeCell_(levelCartToLevelCompressed[0], entry.cell2);
865 if ((c1 < 0) || (c2 < 0)) {
871 trans = this->globalTrans().transmissibility(c1, c2);
873 if (! generatedNnc.empty()) {
874 for (
const auto& generated : generatedNnc) {
875 if (entry.cell1 == generated.cell1 && entry.cell2 == generated.cell2) {
876 trans -= generated.trans;
882 const auto tt = unitSystem
883 .from_si(UnitSystem::measure::transmissibility, trans);
888 if (std::isnormal(tt) && ! (tt < 1.0e-6)) {
889 inputedNnc.emplace_back(entry.cell1, entry.cell2, trans);
894 this->outputNnc_[0].insert(this->outputNnc_[0].begin(), inputedNnc.begin(), inputedNnc.end());
895 return this->outputNnc_;
898template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
899void EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
900doWriteOutput(
const int reportStepNum,
901 const std::optional<int> timeStepNum,
902 const bool isSubStep,
903 const bool isForcedFinalOutput,
904 data::Solution&& localCellData,
905 data::Wells&& localWellData,
906 data::GroupAndNetworkValues&& localGroupAndNetworkData,
907 data::Aquifers&& localAquiferData,
908 WellTestState&& localWTestState,
909 const Action::State& actionState,
910 const UDQState& udqState,
911 const SummaryState& summaryState,
912 const std::vector<Scalar>& thresholdPressure,
915 bool doublePrecision,
917 std::array<FlowsData<double>, 3>&& flowsn,
919 std::array<FlowsData<double>, 3>&& floresn)
921 const auto isParallel = this->collectOnIORank_.isParallel();
922 const bool needsReordering = this->collectOnIORank_.doesNeedReordering();
924 RestartValue restartValue {
925 (isParallel || needsReordering)
926 ? this->collectOnIORank_.globalCellData()
927 : std::move(localCellData),
929 isParallel ? this->collectOnIORank_.globalWellData()
930 : std::move(localWellData),
932 isParallel ? this->collectOnIORank_.globalGroupAndNetworkData()
933 : std::move(localGroupAndNetworkData),
935 isParallel ? this->collectOnIORank_.globalAquiferData()
936 : std::move(localAquiferData)
939 if (eclState_.getSimulationConfig().useThresholdPressure()) {
940 restartValue.addExtra(
"THRESHPR", UnitSystem::measure::pressure,
946 restartValue.addExtra(
"OPMEXTRA", std::vector<double>(1, nextStepSize));
951 const auto flowsn_global = isParallel ? this->collectOnIORank_.globalFlowsn() : std::move(flowsn);
952 for (
const auto& flows : flowsn_global) {
953 if (flows.name.empty())
955 if (flows.name ==
"FLOGASN+") {
956 restartValue.addExtra(flows.name, UnitSystem::measure::gas_surface_rate, flows.values);
958 restartValue.addExtra(flows.name, UnitSystem::measure::liquid_surface_rate, flows.values);
963 const auto floresn_global = isParallel ? this->collectOnIORank_.globalFloresn() : std::move(floresn);
964 for (
const auto& flores : floresn_global) {
965 if (flores.name.empty()) {
968 restartValue.addExtra(flores.name, UnitSystem::measure::rate, flores.values);
972 std::vector<Opm::RestartValue> restartValues{};
974 if ( !isParallel && !needsReordering && (this->eclState_.getLgrs().size()>0) && (this->grid_.maxLevel()>0) ) {
979 Opm::Lgr::extractRestartValueLevelGrids<Grid>(this->grid_, restartValue, restartValues);
982 restartValues.reserve(1);
983 restartValues.push_back(std::move(restartValue));
989 this->taskletRunner_->barrier();
992 if (this->taskletRunner_->failure()) {
993 throw std::runtime_error(
"Failure in the TaskletRunner while writing output.");
997 auto eclWriteTasklet = std::make_shared<EclWriteTasklet>(
999 isParallel ? this->collectOnIORank_.globalWellTestState() : std::move(localWTestState),
1000 summaryState, udqState, *this->eclIO_,
1001 reportStepNum, timeStepNum, isSubStep, curTime, std::move(restartValues), doublePrecision,
1002 isForcedFinalOutput);
1005 this->taskletRunner_->dispatch(std::move(eclWriteTasklet));
1008template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
1009void EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
1010evalSummary(
const int reportStepNum,
1011 const Scalar curTime,
1012 const data::Wells& localWellData,
1013 const data::WellBlockAveragePressures& localWBPData,
1014 const data::GroupAndNetworkValues& localGroupAndNetworkData,
1015 const std::map<int,data::AquiferData>& localAquiferData,
1016 const std::map<std::pair<std::string, int>,
double>& blockData,
1017 const std::map<std::string, double>& miscSummaryData,
1018 const std::map<std::string, std::vector<double>>& regionData,
1019 const Inplace& inplace,
1020 const Inplace* initialInPlace,
1021 const InterRegFlowMap& interRegFlows,
1022 SummaryState& summaryState,
1025 if (collectOnIORank_.isIORank()) {
1026 const auto& summary = eclIO_->summary();
1028 const auto& wellData = this->collectOnIORank_.isParallel()
1029 ? this->collectOnIORank_.globalWellData()
1032 const auto& wbpData = this->collectOnIORank_.isParallel()
1033 ? this->collectOnIORank_.globalWBPData()
1036 const auto& groupAndNetworkData = this->collectOnIORank_.isParallel()
1037 ? this->collectOnIORank_.globalGroupAndNetworkData()
1038 : localGroupAndNetworkData;
1040 const auto& aquiferData = this->collectOnIORank_.isParallel()
1041 ? this->collectOnIORank_.globalAquiferData()
1044 const auto interreg_flows = getInterRegFlowsAsMap(interRegFlows);
1046 const auto values = out::Summary::DynamicSimulatorState {
1049 &groupAndNetworkData,
1061 summary.eval(reportStepNum, curTime, values, summaryState);
1067 const auto udq_step = reportStepNum - 1;
1069 this->schedule_[udq_step].udq()
1071 this->schedule_.wellMatcher(udq_step),
1072 this->schedule_[udq_step].group_order(),
1073 this->schedule_.segmentMatcherFactory(udq_step),
1074 [es = std::cref(this->eclState_)]() {
1075 return std::make_unique<RegionSetMatcher>
1076 (es.get().fipRegionStatistics());
1078 summaryState, udqState);
1082 if (collectOnIORank_.isParallel()) {
1083 Parallel::MpiSerializer ser(grid_.comm());
1084 ser.append(summaryState);
1089template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
1090const typename EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::TransmissibilityType&
1091EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
1094 assert (globalTrans_);
1095 return *globalTrans_;
Collects necessary output values and pass it to opm-common's ECL output.
Definition CollectDataOnIORank.hpp:49
Definition EclGenericWriter.hpp:67
Inter-region flow accumulation maps for all region definition arrays.
Definition InterRegFlows.hpp:179
const std::vector< std::string > & names() const
Names of all applicable region definition arrays.
Definition InterRegFlows.cpp:181
std::vector< data::InterRegFlowMap > getInterRegFlows() const
Get read-only access to the underlying CSR representation.
Definition InterRegFlows.cpp:187
Definition EclGenericWriter.hpp:50
The base class for tasklets.
Definition tasklets.hpp:45
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:45