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)
239 this->eclIO_ = std::make_unique<EclipseIO>
241 UgGridHelpers::createEclipseGrid(*equilGrid,
eclState_.getInputGrid()),
242 this->schedule_, summaryConfig,
"", enableEsmry);
247 int numWorkerThreads = 0;
249 numWorkerThreads = 1;
255template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
263template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
267 if (collectOnIORank_.isIORank()) {
268 std::map<std::string, std::vector<int>> integerVectors;
269 if (collectOnIORank_.isParallel()) {
270 integerVectors.emplace(
"MPI_RANK", collectOnIORank_.globalRanks());
273 eclIO_->writeInitial(*this->outputTrans_,
276 this->outputTrans_.reset();
279template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
284 if (collectOnIORank_.isIORank()) {
285 auto cartMap = cartesianToCompressed(equilGrid_->size(0), UgGridHelpers::globalCell(*equilGrid_));
286 computeTrans_(cartMap, map);
287 exportNncStructure_(cartMap, map);
291 if (collectOnIORank_.isParallel()) {
292 const auto& comm = grid_.comm();
299template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
302computeTrans_(
const std::unordered_map<int,int>& cartesianToActive,
303 const std::function<
unsigned int(
unsigned int)>& map)
const
306 outputTrans_ = std::make_unique<data::Solution>();
309 const auto& cartMapper = *equilCartMapper_;
310 const auto& cartDims = cartMapper.cartesianDimensions();
312 auto createCellData = [&cartDims]() {
313 return data::CellData{
314 UnitSystem::measure::transmissibility,
315 std::vector<double>(cartDims[0] * cartDims[1] * cartDims[2], 0.0),
316 data::TargetType::INIT
320 outputTrans_->clear();
321 outputTrans_->emplace(
"TRANX", createCellData());
322 outputTrans_->emplace(
"TRANY", createCellData());
323 outputTrans_->emplace(
"TRANZ", createCellData());
325 auto& tranx = this->outputTrans_->at(
"TRANX");
326 auto& trany = this->outputTrans_->at(
"TRANY");
327 auto& tranz = this->outputTrans_->at(
"TRANZ");
329 using GlobalGridView =
typename EquilGrid::LeafGridView;
330 using GlobElementMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GlobalGridView>;
331 const GlobalGridView& globalGridView = this->equilGrid_->leafGridView();
332 const GlobElementMapper globalElemMapper { globalGridView, Dune::mcmgElementLayout() };
334 auto isNumAquCell = [numAquCell = this->eclState_.aquifer().hasNumericalAquifer()
335 ? this->eclState_.aquifer().numericalAquifers().allAquiferCellIds()
336 : std::vector<std::size_t>{}]
337 (
const std::size_t cellIdx)
339 return std::binary_search(numAquCell.begin(), numAquCell.end(), cellIdx);
342 for (
const auto& elem : elements(globalGridView)) {
343 for (
const auto& is : intersections(globalGridView, elem)) {
348 unsigned c1 = globalElemMapper.index(is.inside());
349 unsigned c2 = globalElemMapper.index(is.outside());
355 const int cartIdx1 = cartMapper.cartesianIndex( c1 );
356 const int cartIdx2 = cartMapper.cartesianIndex( c2 );
358 if (isNumAquCell(cartIdx1) || isNumAquCell(cartIdx2)) {
367 assert(cartIdx1 <= cartIdx2);
368 int gc1 = std::min(cartIdx1, cartIdx2);
369 int gc2 = std::max(cartIdx1, cartIdx2);
377 if (gc2 - gc1 == 1 && cartDims[0] > 1 ) {
378 tranx.template data<double>()[gc1] = globalTrans().transmissibility(c1, c2);
382 if (gc2 - gc1 == cartDims[0] && cartDims[1] > 1) {
383 trany.template data<double>()[gc1] = globalTrans().transmissibility(c1, c2);
387 if ( gc2 - gc1 == cartDims[0]*cartDims[1] ||
388 directVerticalNeighbors(cartDims, cartesianToActive, gc1, gc2))
389 tranz.template data<double>()[gc1] = globalTrans().transmissibility(c1, c2);
394template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
396EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
397exportNncStructure_(
const std::unordered_map<int,int>& cartesianToActive,
398 const std::function<
unsigned int(
unsigned int)>& map)
const
400 auto isNumAquCell = [numAquCell = this->eclState_.aquifer().hasNumericalAquifer()
401 ? this->eclState_.aquifer().numericalAquifers().allAquiferCellIds()
402 : std::vector<std::size_t>{}]
403 (
const std::size_t cellIdx)
405 return std::binary_search(numAquCell.begin(), numAquCell.end(), cellIdx);
408 auto isNumAquConn = [&isNumAquCell](
const std::size_t cellIdx1,
409 const std::size_t cellIdx2)
411 return isNumAquCell(cellIdx1) || isNumAquCell(cellIdx2);
414 auto isCartesianNeighbour = [nx = this->eclState_.getInputGrid().getNX(),
415 ny = this->eclState_.getInputGrid().getNY()]
416 (
const std::size_t cellIdx1,
const std::size_t cellIdx2)
418 const auto cellDiff = cellIdx2 - cellIdx1;
420 return (cellDiff == 1)
422 || (cellDiff == nx * ny);
425 auto activeCell = [&cartesianToActive](
const std::size_t cellIdx)
427 auto pos = cartesianToActive.find(cellIdx);
428 return (pos == cartesianToActive.end()) ? -1 : pos->second;
431 const auto& nncData = this->eclState_.getInputNNC().input();
432 const auto& unitSystem = this->eclState_.getDeckUnitSystem();
434 for (
const auto& entry : nncData) {
443 assert (entry.cell2 >= entry.cell1);
445 if (! isCartesianNeighbour(entry.cell1, entry.cell2) ||
446 isNumAquConn(entry.cell1, entry.cell2))
451 const auto c1 = activeCell(entry.cell1);
452 const auto c2 = activeCell(entry.cell2);
454 if ((c1 < 0) || (c2 < 0)) {
460 const auto trans = this->globalTrans().transmissibility(c1, c2);
461 const auto tt = unitSystem
462 .from_si(UnitSystem::measure::transmissibility, trans);
467 if (std::isnormal(tt) && ! (tt < 1.0e-6)) {
468 this->outputNnc_.emplace_back(entry.cell1, entry.cell2, trans);
473 auto isDirectNeighbours = [&isCartesianNeighbour, &cartesianToActive,
474 cartDims = &this->cartMapper_.cartesianDimensions()]
475 (
const std::size_t cellIdx1,
const std::size_t cellIdx2)
477 return isCartesianNeighbour(cellIdx1, cellIdx2)
478 || directVerticalNeighbors(*cartDims, cartesianToActive, cellIdx1, cellIdx2);
481 using GlobalGridView =
typename EquilGrid::LeafGridView;
482 using GlobElementMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GlobalGridView>;
483 const GlobalGridView& globalGridView = this->equilGrid_->leafGridView();
484 const GlobElementMapper globalElemMapper { globalGridView, Dune::mcmgElementLayout() };
487 const auto& equilCartMapper = *equilCartMapper_;
488 for (
const auto& elem : elements(globalGridView)) {
489 for (
const auto& is : intersections(globalGridView, elem)) {
494 unsigned c1 = globalElemMapper.index(is.inside());
495 unsigned c2 = globalElemMapper.index(is.outside());
500 std::size_t cc1 = equilCartMapper.cartesianIndex( c1 );
501 std::size_t cc2 = equilCartMapper.cartesianIndex( c2 );
512 if (isNumAquConn(cc1, cc2) || ! isDirectNeighbours(cc1, cc2)) {
515 auto t = this->globalTrans().transmissibility(c1, c2);
516 auto candidate = std::lower_bound(nncData.begin(), nncData.end(),
517 NNCdata { cc1, cc2, 0.0 });
519 while ((candidate != nncData.end()) &&
520 (candidate->cell1 == cc1) &&
521 (candidate->cell2 == cc2))
523 t -= candidate->trans;
532 const auto tt = unitSystem
533 .from_si(UnitSystem::measure::transmissibility, t);
535 if (std::isnormal(tt) && (tt > 1.0e-12)) {
536 this->outputNnc_.emplace_back(cc1, cc2, t);
542 return this->outputNnc_;
545template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
548 const std::optional<int> timeStepNum,
549 const bool isSubStep,
550 data::Solution&& localCellData,
551 data::Wells&& localWellData,
552 data::GroupAndNetworkValues&& localGroupAndNetworkData,
553 data::Aquifers&& localAquiferData,
554 WellTestState&& localWTestState,
555 const Action::State& actionState,
556 const UDQState& udqState,
557 const SummaryState& summaryState,
558 const std::vector<Scalar>& thresholdPressure,
561 bool doublePrecision,
567 const auto isParallel = this->collectOnIORank_.isParallel();
568 const bool needsReordering = this->collectOnIORank_.doesNeedReordering();
570 RestartValue restartValue {
571 (isParallel || needsReordering)
572 ? this->collectOnIORank_.globalCellData()
573 : std::move(localCellData),
575 isParallel ? this->collectOnIORank_.globalWellData()
576 : std::move(localWellData),
578 isParallel ? this->collectOnIORank_.globalGroupAndNetworkData()
579 : std::move(localGroupAndNetworkData),
581 isParallel ? this->collectOnIORank_.globalAquiferData()
582 : std::move(localAquiferData)
585 if (eclState_.getSimulationConfig().useThresholdPressure()) {
586 restartValue.addExtra(
"THRESHPR", UnitSystem::measure::pressure,
592 restartValue.addExtra(
"OPMEXTRA", std::vector<double>(1, nextStepSize));
597 const auto flowsn_global = isParallel ? this->collectOnIORank_.globalFlowsn() : std::move(flowsn);
598 for (
const auto& flows : flowsn_global) {
599 if (flows.name.empty())
601 if (flows.name ==
"FLOGASN+") {
602 restartValue.addExtra(flows.name, UnitSystem::measure::gas_surface_rate, flows.values);
604 restartValue.addExtra(flows.name, UnitSystem::measure::liquid_surface_rate, flows.values);
609 const auto floresn_global = isParallel ? this->collectOnIORank_.globalFloresn() : std::move(floresn);
610 for (
const auto& flores : floresn_global) {
611 if (flores.name.empty()) {
614 restartValue.addExtra(flores.name, UnitSystem::measure::rate, flores.values);
618 std::vector<Opm::RestartValue> restartValues{};
620 if ( !isParallel && !needsReordering && (this->eclState_.getLgrs().size()>0) && (this->grid_.maxLevel()>0) ) {
625 Opm::Lgr::extractRestartValueLevelGrids<Grid>(this->grid_, restartValue, restartValues);
628 restartValues.reserve(1);
629 restartValues.push_back(std::move(restartValue));
635 this->taskletRunner_->barrier();
638 if (this->taskletRunner_->failure()) {
639 throw std::runtime_error(
"Failure in the TaskletRunner while writing output.");
643 auto eclWriteTasklet = std::make_shared<EclWriteTasklet>(
645 isParallel ? this->collectOnIORank_.globalWellTestState() : std::move(localWTestState),
646 summaryState, udqState, *this->eclIO_,
647 reportStepNum, timeStepNum, isSubStep, curTime, std::move(restartValues), doublePrecision);
650 this->taskletRunner_->dispatch(std::move(eclWriteTasklet));
653template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
656 const Scalar curTime,
657 const data::Wells& localWellData,
658 const data::WellBlockAveragePressures& localWBPData,
659 const data::GroupAndNetworkValues& localGroupAndNetworkData,
660 const std::map<int,data::AquiferData>& localAquiferData,
661 const std::map<std::pair<std::string, int>,
double>& blockData,
662 const std::map<std::string, double>& miscSummaryData,
663 const std::map<std::string, std::vector<double>>& regionData,
664 const Inplace& inplace,
665 const std::optional<Inplace>& initialInPlace,
667 SummaryState& summaryState,
670 if (collectOnIORank_.isIORank()) {
671 const auto& summary = eclIO_->summary();
673 const auto& wellData = this->collectOnIORank_.isParallel()
674 ? this->collectOnIORank_.globalWellData()
677 const auto& wbpData = this->collectOnIORank_.isParallel()
678 ? this->collectOnIORank_.globalWBPData()
681 const auto& groupAndNetworkData = this->collectOnIORank_.isParallel()
682 ? this->collectOnIORank_.globalGroupAndNetworkData()
683 : localGroupAndNetworkData;
685 const auto& aquiferData = this->collectOnIORank_.isParallel()
686 ? this->collectOnIORank_.globalAquiferData()
689 summary.eval(summaryState,
701 getInterRegFlowsAsMap(interRegFlows));
707 const auto udq_step = reportStepNum - 1;
709 this->schedule_[udq_step].udq()
711 this->schedule_.wellMatcher(udq_step),
712 this->schedule_[udq_step].group_order(),
713 this->schedule_.segmentMatcherFactory(udq_step),
714 [es = std::cref(this->eclState_)]() {
715 return std::make_unique<RegionSetMatcher>
716 (es.get().fipRegionStatistics());
718 summaryState, udqState);
722 if (collectOnIORank_.isParallel()) {
729template<
class Gr
id,
class EquilGr
id,
class Gr
idView,
class ElementMapper,
class Scalar>
734 assert (globalTrans_);
735 return *globalTrans_;
Definition: CollectDataOnIORank.hpp:49
bool isIORank() const
Definition: CollectDataOnIORank.hpp:126
Definition: EclGenericWriter.hpp:65
const EclipseState & eclState_
Definition: EclGenericWriter.hpp:156
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:547
CollectDataOnIORankType collectOnIORank_
Definition: EclGenericWriter.hpp:152
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 std::optional< Inplace > &initialInPlace, const InterRegFlowMap &interRegFlows, SummaryState &summaryState, UDQState &udqState)
Definition: EclGenericWriter_impl.hpp:655
void extractOutputTransAndNNC(const std::function< unsigned int(unsigned int)> &map)
Definition: EclGenericWriter_impl.hpp:282
std::unique_ptr< TaskletRunner > taskletRunner_
Definition: EclGenericWriter.hpp:158
void writeInit()
Definition: EclGenericWriter_impl.hpp:265
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:257
const TransmissibilityType & globalTrans() const
Definition: EclGenericWriter_impl.hpp:732
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
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:43
Avoid mistakes in calls to broadcast() by wrapping the root argument in an explicit type.
Definition: MPISerializer.hpp:33