EclGenericWriter_impl.hpp
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 2 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
23#ifndef OPM_ECL_GENERIC_WRITER_IMPL_HPP
24#define OPM_ECL_GENERIC_WRITER_IMPL_HPP
25
26#include <dune/grid/common/mcmgmapper.hh>
27
28#include <opm/grid/cpgrid/LgrOutputHelpers.hpp>
29#include <opm/grid/GridHelpers.hpp>
30#include <opm/grid/utility/cartesianToCompressed.hpp>
31
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>
36
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>
45
46#include <opm/input/eclipse/Units/UnitSystem.hpp>
47
48#include <opm/output/eclipse/EclipseIO.hpp>
49#include <opm/output/eclipse/RestartValue.hpp>
50#include <opm/output/eclipse/Summary.hpp>
51
53
54#if HAVE_MPI
56#endif
57
58#if HAVE_MPI
59#include <mpi.h>
60#endif
61
62#include <algorithm>
63#include <array>
64#include <cassert>
65#include <cmath>
66#include <functional>
67#include <map>
68#include <memory>
69#include <string>
70#include <unordered_map>
71#include <utility>
72#include <vector>
73
74namespace {
75
90bool directVerticalNeighbors(const std::array<int, 3>& cartDims,
91 const std::unordered_map<int,int>& cartesianToActive,
92 int smallGlobalIndex, int largeGlobalIndex)
93{
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];
99 gc /= cartDims[0];
100 ijk[1] = gc % cartDims[1];
101 ijk[2] = gc / cartDims[1];
102 return ijk;
103 };
104 ijk1 = globalToIjk(smallGlobalIndex);
105 ijk2 = globalToIjk(largeGlobalIndex);
106 assert(ijk2[2]>=ijk1[2]);
107
108 if ( ijk1[0] == ijk2[0] && ijk1[1] == ijk2[1] && (ijk2[2] - ijk1[2]) > 1)
109 {
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] )
113 {
114 if ( cartesianToActive.find( gi ) != cartesianToActive.end() )
116 return false;
117 }
119 return true;
120 } else
121 return false;
122}
123
124std::unordered_map<std::string, Opm::data::InterRegFlowMap>
125getInterRegFlowsAsMap(const Opm::InterRegFlowMap& map)
126{
127 auto maps = std::unordered_map<std::string, Opm::data::InterRegFlowMap>{};
128
129 const auto& regionNames = map.names();
130 auto flows = map.getInterRegFlows();
131 const auto nmap = regionNames.size();
132
133 maps.reserve(nmap);
134 for (auto mapID = 0*nmap; mapID < nmap; ++mapID) {
135 maps.emplace(regionNames[mapID], std::move(flows[mapID]));
136 }
137
138 return maps;
139}
140
141struct EclWriteTasklet : public Opm::TaskletInterface
143 Opm::Action::State actionState_;
144 Opm::WellTestState wtestState_;
145 Opm::SummaryState summaryState_;
146 Opm::UDQState udqState_;
147 Opm::EclipseIO& eclIO_;
148 int reportStepNum_;
149 std::optional<int> timeStepNum_;
150 bool isSubStep_;
151 double secondsElapsed_;
152 std::vector<Opm::RestartValue> restartValue_;
153 bool writeDoublePrecision_;
155 bool forcedSimulationFinished_;
156
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,
162 int reportStepNum,
163 std::optional<int> timeStepNum,
164 bool isSubStep,
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)
173 , eclIO_(eclIO)
174 , reportStepNum_(reportStepNum)
175 , timeStepNum_(timeStepNum)
176 , isSubStep_(isSubStep)
177 , secondsElapsed_(secondsElapsed)
178 , restartValue_(std::move(restartValue))
179 , writeDoublePrecision_(writeDoublePrecision)
180 , forcedSimulationFinished_(forcedSimulationFinished)
181 {}
182
183 // callback to eclIO serial writeTimeStep method
184 void run() override
185 {
186 if (this->restartValue_.size() == 1) {
187 this->eclIO_.writeTimeStep(this->actionState_,
188 this->wtestState_,
189 this->summaryState_,
190 this->udqState_,
191 this->reportStepNum_,
192 this->isSubStep_,
193 this->secondsElapsed_,
194 std::move(this->restartValue_.back()),
195 this->writeDoublePrecision_,
196 this->timeStepNum_,
197 forcedSimulationFinished_);
198 }
199 else{
200 this->eclIO_.writeTimeStep(this->actionState_,
201 this->wtestState_,
202 this->summaryState_,
203 this->udqState_,
204 this->reportStepNum_,
205 this->isSubStep_,
206 this->secondsElapsed_,
207 std::move(this->restartValue_),
208 this->writeDoublePrecision_,
209 this->timeStepNum_,
210 forcedSimulationFinished_);
211 }
212 }
213};
214
215}
216
217namespace Opm {
218
219template<class Grid, class EquilGrid, class GridView, class ElementMapper, class Scalar>
221EclGenericWriter(const Schedule& schedule,
222 const EclipseState& eclState,
223 const SummaryConfig& summaryConfig,
224 const Grid& grid,
225 const EquilGrid* equilGrid,
226 const GridView& gridView,
227 const Dune::CartesianIndexMapper<Grid>& cartMapper,
228 const Dune::CartesianIndexMapper<EquilGrid>* equilCartMapper,
229 bool enableAsyncOutput,
230 bool enableEsmry )
231 : collectOnIORank_(grid,
232 equilGrid,
233 gridView,
234 cartMapper,
235 equilCartMapper,
236 summaryConfig.fip_regions_interreg_flow())
237 , grid_ (grid)
238 , gridView_ (gridView)
239 , schedule_ (schedule)
240 , eclState_ (eclState)
241 , cartMapper_ (cartMapper)
242 , equilCartMapper_(equilCartMapper)
243 , equilGrid_ (equilGrid)
244{
245 // Make sure outputNnc_ vector has at least 1 entry in all ranks.
246 outputNnc_.resize(1);
247
248 if (this->collectOnIORank_.isIORank()) {
249 this->eclIO_ = std::make_unique<EclipseIO>
250 (this->eclState_,
251 UgGridHelpers::createEclipseGrid(*equilGrid, eclState_.getInputGrid()),
252 this->schedule_, summaryConfig, "", enableEsmry);
253 }
254
255 // create output thread if enabled and rank is I/O rank
256 // async output is enabled by default if pthread are enabled
257 int numWorkerThreads = 0;
258 if (enableAsyncOutput && collectOnIORank_.isIORank()) {
259 numWorkerThreads = 1;
260 }
261
262 this->taskletRunner_.reset(new TaskletRunner(numWorkerThreads));
263}
264
265template<class Grid, class EquilGrid, class GridView, class ElementMapper, class Scalar>
267eclIO() const
268{
269 assert(eclIO_);
270 return *eclIO_;
271}
272
273template<class Grid, class EquilGrid, class GridView, class ElementMapper, class Scalar>
275writeInit()
276{
277 if (collectOnIORank_.isIORank()) {
278 std::map<std::string, std::vector<int>> integerVectors;
279 if (collectOnIORank_.isParallel()) {
280 integerVectors.emplace("MPI_RANK", collectOnIORank_.globalRanks());
281 }
282
283 if (const auto& lgrs = this->eclState_.getLgrs(); lgrs.size() > 0) {
284
285 const auto nncCollection = Opm::NNCCollection::fromLGROutputContainers(this->outputNnc_,
286 this->outputNncGlobalLocal_,
287 this->outputAmalgamatedNnc_);
288 eclIO_->writeInitial(*this->outputTrans_,
289 integerVectors,
290 nncCollection);
291 } else {
292 eclIO_->writeInitial(*this->outputTrans_,
293 integerVectors,
294 this->outputNnc_.front());
295 }
296 this->outputTrans_.reset();
297 }
298}
299
300template<class Grid, class EquilGrid, class GridView, class ElementMapper, class Scalar>
301void
303extractOutputTransAndNNC(const std::function<unsigned int(unsigned int)>& map)
304{
305 if (collectOnIORank_.isIORank()) {
306 constexpr bool equilGridIsCpGrid = std::is_same_v<EquilGrid, Dune::CpGrid>;
307
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>();
314
315 computeTrans_(levelCartToLevelCompressed, map, computeLevelIndices,
316 computeLevelCartIdx, computeLevelCartDimensions, computeOriginIndices);
317 exportNncStructure_(levelCartToLevelCompressed, map, computeLevelIndices, computeLevelCartIdx,
318 computeLevelCartDimensions, computeOriginIndices);
319 }
320
321#if HAVE_MPI
322 if (collectOnIORank_.isParallel()) {
323 const auto& comm = grid_.comm();
324 Parallel::MpiSerializer ser(comm);
325 ser.broadcast(Parallel::RootRank{0}, outputNnc_);
326 }
327#endif
328}
329
330template<class Grid, class EquilGrid, class GridView, class ElementMapper, class Scalar>
331bool
333isNumAquCell_(const std::size_t cartIdx) const
334{
335 const auto& numAquCell = this->eclState_.aquifer().hasNumericalAquifer()
336 ? this->eclState_.aquifer().numericalAquifers().allAquiferCellIds()
337 : std::vector<std::size_t>{};
338
339 return std::ranges::binary_search(numAquCell.begin(), numAquCell.end(), cartIdx);
340}
341
342template<class Grid, class EquilGrid, class GridView, class ElementMapper, class Scalar>
343bool
344EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
345isNumAquConn_(const std::size_t cartIdx1,
346 const std::size_t cartIdx2) const
347{
348 return isNumAquCell_(cartIdx1) || isNumAquCell_(cartIdx2);
349}
350
351template<class Grid, class EquilGrid, class GridView, class ElementMapper, class Scalar>
352template<bool equilGridIsCpGrid>
354EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
355createLevelCartMapp_() const
356{
357 if constexpr (equilGridIsCpGrid) {
358 return Opm::LevelCartesianIndexMapper<EquilGrid>(*this->equilGrid_);
359 } else {
360 return Opm::LevelCartesianIndexMapper<EquilGrid>(*equilCartMapper_); }
361}
362
363template<class Grid, class EquilGrid, class GridView, class ElementMapper, class Scalar>
364template<bool equilGridIsCpGrid>
365std::vector<std::unordered_map<int,int>>
366EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
367createCartesianToActiveMaps_(const Opm::LevelCartesianIndexMapper<EquilGrid>& levelCartMapp) const
368{
369 if constexpr (equilGridIsCpGrid) {
370 if (this->equilGrid_->maxLevel()) {
371 return Opm::Lgr::levelCartesianToLevelCompressedMaps(*this->equilGrid_, levelCartMapp); }
372 else {
373 return std::vector<std::unordered_map<int,int>>{ cartesianToCompressed(equilGrid_->size(0), UgGridHelpers::globalCell(*equilGrid_)) };
374 }
375 }
376 return std::vector<std::unordered_map<int,int>>{ cartesianToCompressed(equilGrid_->size(0), UgGridHelpers::globalCell(*equilGrid_)) };
377}
378
379template<class Grid, class EquilGrid, class GridView, class ElementMapper, class Scalar>
380template<bool equilGridIsCpGrid>
381std::function<std::array<int,3>(int)>
382EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
383computeLevelCartDimensions_(const Opm::LevelCartesianIndexMapper<EquilGrid>& levelCartMapp,
384 const Dune::CartesianIndexMapper<EquilGrid>& equilCartMapp) const
385{
386 if constexpr (equilGridIsCpGrid) {
387 return [&](int level)
388 {
389 return levelCartMapp.cartesianDimensions(level);
390 };
391 }
392 else {
393 return [&]([[maybe_unused]] int level)
394 {
395 assert(level == 0); // refinement only supported for CpGrid for now
396 return equilCartMapp.cartesianDimensions();
397 };
398 }
399}
400
401template<class Grid, class EquilGrid, class GridView, class ElementMapper, class Scalar>
402template<bool equilGridIsCpGrid>
403std::function<int(int, int)>
404EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
405computeLevelCartIdx_(const Opm::LevelCartesianIndexMapper<EquilGrid>& levelCartMapp,
406 const Dune::CartesianIndexMapper<EquilGrid>& equilCartMapp) const
407{
408 if constexpr (equilGridIsCpGrid) {
409 return [&](int levelCompressedIdx,
410 int level)
411 {
412 return levelCartMapp.cartesianIndex(levelCompressedIdx, level);
413 };
414 }
415 else {
416 return [&](int levelCompressedIdx,
417 [[maybe_unused]] int level)
418 {
419 assert(level == 0); // refinement only supported for CpGrid for now
420 return equilCartMapp.cartesianIndex(levelCompressedIdx);
421 };
422 }
423}
424
425template<class Grid, class EquilGrid, class GridView, class ElementMapper, class Scalar>
426template <bool equilGridIsCpGrid>
427auto
428EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
429computeLevelIndices_() const
430{
431 if constexpr (equilGridIsCpGrid) {
432 return [](const auto& intersection,
433 const auto&,
434 const auto&)
435 {
436 return std::pair{intersection.inside().getLevelElem().index(), intersection.outside().getLevelElem().index()};
437 };
438 }
439 else {
440 return [](const auto&,
441 const auto& intersectionInsideLeafIdx,
442 const auto& intersectionOutsideLeafIdx)
443 {
444 return std::pair{intersectionInsideLeafIdx, intersectionOutsideLeafIdx};
445 };
446 }
447}
448
449template<class Grid, class EquilGrid, class GridView, class ElementMapper, class Scalar>
450template <bool equilGridIsCpGrid>
451auto
452EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
453computeOriginIndices_() const
454{
455 if constexpr (equilGridIsCpGrid) {
456 return [](const auto& intersection,
457 const auto&,
458 const auto&)
459 {
460 return std::pair{intersection.inside().getOrigin().index(), intersection.outside().getOrigin().index()};
461 };
462 }
463 else {
464 return [](const auto&,
465 const auto& intersectionInsideLeafIdx,
466 const auto& intersectionOutsideLeafIdx)
467 {
468 return std::pair{intersectionInsideLeafIdx, intersectionOutsideLeafIdx};
469 };
470 }
471}
472
473template<class Grid, class EquilGrid, class GridView, class ElementMapper, class Scalar>
474void
475EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
476allocateLevelTrans_(const std::array<int,3>& levelCartDims,
477 data::Solution& levelTrans) const
478{
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
484 };
485 };
486
487 levelTrans.clear();
488 levelTrans.emplace("TRANX", createLevelCellData());
489 levelTrans.emplace("TRANY", createLevelCellData());
490 levelTrans.emplace("TRANZ", createLevelCellData());
491}
492
493template<class Grid, class EquilGrid, class GridView, class ElementMapper, class Scalar>
494template<typename LevelIndicesFunction, typename OriginIndicesFunction>
495void
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
503{
504 if (!outputTrans_) {
505 outputTrans_ = std::make_unique<std::vector<data::Solution>>(std::vector<data::Solution>{});
506 }
507
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() };
512
513 // Refinement supported only for CpGrid for now.
514 int maxLevel = this->equilGrid_->maxLevel();
515
516 outputTrans_->resize(maxLevel+1); // including level zero grid
517
518 for (int level = 0; level <= maxLevel; ++level) {
519 allocateLevelTrans_(computeLevelCartDims(level), this->outputTrans_->at(level));
520 }
521
522 for (const auto& elem : elements(globalGridView)) {
523 for (const auto& is : intersections(globalGridView, elem)) {
524 if (!is.neighbor())
525 continue; // intersection is on the domain boundary
526
527 if ( is.inside().level() != is.outside().level() ) // Those are treated as NNCs
528 continue;
529
530 // Not 'const' because remapped if 'map' is non-null.
531 unsigned c1 = globalElemMapper.index(is.inside());
532 unsigned c2 = globalElemMapper.index(is.outside());
533
534 if (c1 > c2)
535 continue; // we only need to handle each connection once, thank you.
536
537 int level = is.inside().level();
538
539 // For CpGrid with LGRs, level*Idx and c* do not coincide.
540 const auto& [levelInIdx, levelOutIdx] = computeLevelIndices(is, c1, c2);
541
542 const int levelCartIdxIn = computeLevelCartIdx(levelInIdx, level);
543 const int levelCartIdxOut = computeLevelCartIdx(levelOutIdx, level);
544
545 // For CpGrid with LGRs, the origin cell index refers to the coarsest
546 // ancestor cell when the cell is refined. For cells not involved in
547 // any refinement, it corresponds to the geometrically equivalent
548 // cell in the level-zero grid.
549 const auto [originInIdx, originOutIdx] = computeOriginIndices(is, c1, c2);
550
551 const auto originCartIdxIn = computeLevelCartIdx(originInIdx, /* level = */ 0);
552 const auto originCartIdxOut = computeLevelCartIdx(originOutIdx, /* level = */ 0);
553
554 // For level-zero grid, level Cartesian indices coincide with the grid Cartesian indices.
555 if (isNumAquCell_(originCartIdxIn) || isNumAquCell_(originCartIdxOut)) {
556 // Check there are no refined aquifer cells.
557 assert(level == 0);
558 // Connections involving numerical aquifers are always NNCs
559 // for the purpose of file output. This holds even for
560 // connections between cells like (I,J,K) and (I+1,J,K)
561 // which are nominally neighbours in the Cartesian grid.
562 continue;
563 }
564
565 const auto minLevelCartIdx = std::min(levelCartIdxIn, levelCartIdxOut);
566 const auto maxLevelCartIdx = std::max(levelCartIdxIn, levelCartIdxOut);
567
568 const auto& levelCartDims = computeLevelCartDims(level);
569
570 // Re-ordering in case of non-empty mapping between equilGrid to grid
571 if (map) {
572 c1 = map(c1); // equilGridToGrid map
573 c2 = map(c2);
574 }
575
576 if (maxLevelCartIdx - minLevelCartIdx == 1 && levelCartDims[0] > 1 ) {
577 outputTrans_->at(level).at("TRANX").template data<double>()[minLevelCartIdx] = globalTrans().transmissibility(c1, c2);
578 continue; // skip other if clauses as they are false, last one needs some computation
579 }
580
581 if (maxLevelCartIdx - minLevelCartIdx == levelCartDims[0] && levelCartDims[1] > 1) {
582 outputTrans_->at(level).at("TRANY").template data<double>()[minLevelCartIdx] = globalTrans().transmissibility(c1, c2);
583 continue; // skipt next if clause as it needs some computation
584 }
585
586 if ( maxLevelCartIdx - minLevelCartIdx == levelCartDims[0]*levelCartDims[1] ||
587 directVerticalNeighbors(levelCartDims,
588 levelCartToLevelCompressed[level],
589 minLevelCartIdx,
590 maxLevelCartIdx)) {
591 outputTrans_->at(level).at("TRANZ").template data<double>()[minLevelCartIdx] = globalTrans().transmissibility(c1, c2);
592 }
593 }
594 }
595}
596
597template<class Grid, class EquilGrid, class GridView, class ElementMapper, class Scalar>
598bool
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
603{
604 const int diff = levelCartIdx2 - levelCartIdx1;
605
606 return (diff == 1)
607 || (diff == levelCartDims[0])
608 || (diff == (levelCartDims[0] * levelCartDims[1]));
609}
610
611template<class Grid, class EquilGrid, class GridView, class ElementMapper, class Scalar>
612bool
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
618{
619 return isCartesianNeighbour_(levelCartDims, levelCartIdx1, levelCartIdx2)
620 || directVerticalNeighbors(levelCartDims, levelCartesianToActive, levelCartIdx1, levelCartIdx2);
621}
622
623template<class Grid, class EquilGrid, class GridView, class ElementMapper, class Scalar>
624auto
625EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
626activeCell_(const std::unordered_map<int,int>& levelCartToLevelCompressed,
627 const std::size_t levelCartIdx) const
628{
629 auto pos = levelCartToLevelCompressed.find(levelCartIdx);
630 return (pos == levelCartToLevelCompressed.end()) ? -1 : pos->second;
631}
632
633template<class Grid, class EquilGrid, class GridView, class ElementMapper, class Scalar>
634void
635EclGenericWriter<Grid,EquilGrid,GridView,ElementMapper,Scalar>::
636allocateAllNncs_(int maxLevel) const
637{
638 this->outputNnc_.resize(maxLevel+1); // level 0,1,..., maxLevel
639
640 if (maxLevel) {
641 // NNCs between main (level zero) grid and LGRs: level 1, ...., maxLevel.
642 // Example: grid with maxLevel == 3, outputNncGlobalLocal_.size() is maxLevel = 3
643 // outputNncGlobalLocal_[0] -> NNCs between level 0 and level 1
644 // outputNncGlobalLocal_[1] -> NNCs between level 0 and level 2
645 // outputAmalgamatedNnc_[2] -> NNCs between level 0 and level 3
646 this->outputNncGlobalLocal_.resize(maxLevel);
647
648 // NNCs between different refined level grids: (level1, level2)
649 // with 0 < level1 < level2 <= maxLevel
650 // Example: grid with maxLevel == 3, outputAmalgamatedNnc_.size() is maxLevel-1 = 2
651 // outputAmalgamatedNnc_[0][0] -> NNCs between level 1 and level 2
652 // outputAmalgamatedNnc_[0][1] -> NNCs between level 1 and level 3
653 // outputAmalgamatedNnc_[1][2] -> NNCs between level 2 and level 3
654 this->outputAmalgamatedNnc_.resize(maxLevel-1);
655 for (int i = 0; i < maxLevel-1; ++i) {
656 this->outputAmalgamatedNnc_[i].resize(maxLevel-1-i);
657 }
658 }
659}
660
661template<class Grid, class EquilGrid, class GridView, 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
671{
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();
677
678 // Cartesian index mapper for the serial I/O grid
679 const auto& equilCartMapper = *equilCartMapper_;
680
681 const auto& level0CartDims = equilCartMapper.cartesianDimensions();
682
683 int maxLevel = this->equilGrid_->maxLevel();
684 allocateAllNncs_(maxLevel);
685
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() };
690
691 for (const auto& elem : elements(globalGridView)) {
692 for (const auto& is : intersections(globalGridView, elem)) {
693 if (!is.neighbor())
694 continue; // intersection is on the domain boundary
695
696 // Not 'const' because remapped if 'map' is non-null.
697 unsigned c1 = globalElemMapper.index(is.inside());
698 unsigned c2 = globalElemMapper.index(is.outside());
699
700 if (c1 > c2)
701 continue; // we only need to handle each connection once, thank you.
702
703 if ( is.inside().level() != is.outside().level() ) { // TRANGL and TRANLL
704 // For CpGrid with LGRs, level*Idx and c* do not coincide.
705 const auto& [levelInIdx, levelOutIdx] = computeLevelIndices(is, c1, c2);
706
707 const int levelIn = is.inside().level();
708 const int levelOut = is.outside().level();
709
710 auto levelCartIdxIn = computeLevelCartIdx(levelInIdx, levelIn);
711 auto levelCartIdxOut = computeLevelCartIdx(levelOutIdx, levelOut);
712
713 // To store correctly and only once the corresponding NNC
714 std::pair<int,int> smallerPair = {levelIn, levelCartIdxIn},
715 largerPair = {levelOut, levelCartIdxOut};
716 if (smallerPair.first > largerPair.first) {
717 std::swap(smallerPair, largerPair);
718 }
719
720 const auto& [smallerLevel, smallerLevelCartIdx] = smallerPair;
721 const auto& [largerLevel, largerLevelCartIdx] = largerPair;
722
723 auto t = this->globalTrans().transmissibility(c1, c2);
724
725 // ECLIPSE ignores NNCs with zero transmissibility
726 // (different threshold than for NNC with corresponding
727 // EDITNNC above). In addition we do set small
728 // transmissibilities to zero when setting up the simulator.
729 // These will be ignored here, too.
730 const auto tt = unitSystem
731 .from_si(UnitSystem::measure::transmissibility, t);
732
733 if (std::isnormal(tt) && (tt > 1.0e-12)) {
734 // Store always FIRST the level Cartesian index of the cell belonging to the smaller level grid involved.
735 if (smallerLevel == 0) { // NNC between main (level zero) grid and a refined level/local grid
736 this->outputNncGlobalLocal_[largerLevel-1].emplace_back(smallerLevelCartIdx, largerLevelCartIdx, t);
737 }
738 else { // NNC between different refined level/local grids -> amlgamated NNC
739 assert(smallerLevel >= 1);
740 this->outputAmalgamatedNnc_[smallerLevel-1][largerLevel-smallerLevel-1].emplace_back(smallerLevelCartIdx, largerLevelCartIdx, t);
741 }
742 }
743 }
744 else {
745 // the cells sharing the intersection belong to the same level
746 assert(is.inside().level() == is.outside().level());
747 const int level = is.inside().level();
748
749 // For CpGrid with LGRs, the origin cell index refers to the coarsest
750 // ancestor cell when the cell is refined. For cells not involved in
751 // any refinement, it corresponds to the geometrically equivalent
752 // cell in the level-zero grid.
753 const auto [originInIdx, originOutIdx] = computeOriginIndices(is, c1, c2);
754
755 const std::size_t originCartIdxIn = computeLevelCartIdx(originInIdx, /* level = */ 0);
756 const std::size_t originCartIdxOut = computeLevelCartIdx(originOutIdx, /* level = */ 0);
757
758 // For CpGrid with LGRs, level*Idx and c* do not coincide.
759 const auto& [levelInIdx, levelOutIdx] = computeLevelIndices(is, c1, c2);
760
761 auto levelCartIdxIn = computeLevelCartIdx(levelInIdx, level);
762 auto levelCartIdxOut = computeLevelCartIdx(levelOutIdx, level);
763
764 if ( levelCartIdxOut < levelCartIdxIn )
765 std::swap(levelCartIdxIn, levelCartIdxOut);
766
767 // Re-ordering in case of non-empty mapping between equilGrid to grid
768 if (map) {
769 c1 = map(c1); // equilGridToGrid map
770 c2 = map(c2);
771 }
772
773 const auto& levelCartDims = computeLevelCartDims(level);
774
775 // Check there are no refined aquifer connections
776 assert(!isNumAquConn_(originCartIdxIn, originCartIdxOut) || level == 0);
777
778 if (isNumAquConn_(originCartIdxIn, originCartIdxOut) ||
779 ! isDirectNeighbours_(levelCartToLevelCompressed[level],
780 levelCartDims,
781 levelCartIdxIn, levelCartIdxOut)) {
782 // We need to check whether an NNC for this face was also
783 // specified via the NNC keyword in the deck.
784 auto t = this->globalTrans().transmissibility(c1, c2);
785
786 if (level == 0) {
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;
791
792 while ((candidate != nncData.end()) &&
793 (candidate->cell1 == originCartIdxIn) &&
794 (candidate->cell2 == originCartIdxOut))
795 {
796 auto trans = candidate->trans;
797 trans *= transMlt;
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;
802 }
803 if (foundNncEditr) {
804 // Only write one value for EDITNNCR, then skip it here and add it on the second loop below
805 break;
806 }
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) {
811 trans *= it->trans;
812 }
813 }
814 t -= trans;
815 ++candidate;
816 }
817 if (foundNncEditr) {
818 // Only write one value for EDITNNCR, then skip it here and add it on the second loop below
819 continue;
820 }
821 }
822
823 // ECLIPSE ignores NNCs with zero transmissibility
824 // (different threshold than for NNC with corresponding
825 // EDITNNC above). In addition we do set small
826 // transmissibilities to zero when setting up the simulator.
827 // These will be ignored here, too.
828 const auto tt = unitSystem
829 .from_si(UnitSystem::measure::transmissibility, t);
830
831 if (std::isnormal(tt) && (tt > 1.0e-12)) {
832 this->outputNnc_[level].emplace_back(levelCartIdxIn, levelCartIdxOut, t);
833 }
834 }
835 }
836 }
837 }
838
839 // Do not include the generated NNCs transsmisibilities in the input NNCs
840 std::vector<NNCdata> inputedNnc{};
841 const auto generatedNnc = outputNnc_[0];
842
843 // The NNC keyword in the deck is defined only for faces in the level-0 grid.
844 // The same limitation applies to aquifer data.
845 for (const auto& entry : nncData) {
846 // Ignore most explicit NNCs between otherwise neighbouring cells.
847 // We keep NNCs that involve cells with numerical aquifers even if
848 // these might be between neighbouring cells in the Cartesian
849 // grid--e.g., between cells (I,J,K) and (I+1,J,K). All such
850 // connections should be written to NNC output arrays provided the
851 // transmissibility value is sufficiently large.
852 //
853 // The condition cell2 >= cell1 holds by construction of nncData.
854 assert (entry.cell2 >= entry.cell1);
855
856 if (! isCartesianNeighbour_(level0CartDims, entry.cell1, entry.cell2) ||
857 isNumAquConn_(entry.cell1, entry.cell2))
858 {
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) {
865 trans *= it->trans;
866 foundNncEdit = true;
867 }
868 }
869 if (! foundNncEdit) {
870 // Pick up transmissibility value from 'globalTrans()' since
871 // multiplier keywords like MULTREGT might have impacted the
872 // values entered in primary sources like NNC/EDITNNC/EDITNNCR.
873 const auto c1 = activeCell_(levelCartToLevelCompressed[/* level */0], entry.cell1);
874 const auto c2 = activeCell_(levelCartToLevelCompressed[/* level */0], entry.cell2);
875
876 if ((c1 < 0) || (c2 < 0)) {
877 // Connection between inactive cells? Unexpected at this
878 // level. Might consider 'throw'ing if this happens...
879 continue;
880 }
881
882 trans = this->globalTrans().transmissibility(c1, c2);
883
884 if (! generatedNnc.empty()) {
885 for (const auto& generated : generatedNnc) {
886 if (entry.cell1 == generated.cell1 && entry.cell2 == generated.cell2) {
887 trans -= generated.trans;
888 break;
889 }
890 }
891 }
892 }
893 const auto tt = unitSystem
894 .from_si(UnitSystem::measure::transmissibility, trans);
895
896 // ECLIPSE ignores NNCs (with EDITNNC/EDITNNCR applied) with
897 // small transmissibility values. Seems like the threshold is
898 // 1.0e-6 in output units.
899 if (std::isnormal(tt) && ! (tt < 1.0e-6)) {
900 inputedNnc.emplace_back(entry.cell1, entry.cell2, trans);
901 }
902 }
903 }
904 // Write first the inputed NNCs and after the internally computed NNCs
905 this->outputNnc_[0].insert(this->outputNnc_[0].begin(), inputedNnc.begin(), inputedNnc.end());
906 return this->outputNnc_;
907}
908
909template<class Grid, class EquilGrid, class GridView, class ElementMapper, class Scalar>
911doWriteOutput(const int reportStepNum,
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,
924 Scalar curTime,
925 Scalar nextStepSize,
926 bool doublePrecision,
927 bool isFlowsn,
928 std::array<FlowsData<double>, 3>&& flowsn,
929 bool isFloresn,
930 std::array<FlowsData<double>, 3>&& floresn)
931{
932 const auto isParallel = this->collectOnIORank_.isParallel();
933 const bool needsReordering = this->collectOnIORank_.doesNeedReordering();
934
935 RestartValue restartValue {
936 (isParallel || needsReordering)
937 ? this->collectOnIORank_.globalCellData()
938 : std::move(localCellData),
939
940 isParallel ? this->collectOnIORank_.globalWellData()
941 : std::move(localWellData),
942
943 isParallel ? this->collectOnIORank_.globalGroupAndNetworkData()
944 : std::move(localGroupAndNetworkData),
945
946 isParallel ? this->collectOnIORank_.globalAquiferData()
947 : std::move(localAquiferData)
948 };
949
950 if (eclState_.getSimulationConfig().useThresholdPressure()) {
951 restartValue.addExtra("THRESHPR", UnitSystem::measure::pressure,
952 thresholdPressure);
953 }
954
955 // Add suggested next timestep to extra data.
956 if (! isSubStep) {
957 restartValue.addExtra("OPMEXTRA", std::vector<double>(1, nextStepSize));
958 }
959
960 // Add nnc flows and flores.
961 if (isFlowsn) {
962 const auto flowsn_global = isParallel ? this->collectOnIORank_.globalFlowsn() : std::move(flowsn);
963 for (const auto& flows : flowsn_global) {
964 if (flows.name.empty())
965 continue;
966 if (flows.name == "FLOGASN+") {
967 restartValue.addExtra(flows.name, UnitSystem::measure::gas_surface_rate, flows.values);
968 } else {
969 restartValue.addExtra(flows.name, UnitSystem::measure::liquid_surface_rate, flows.values);
970 }
971 }
972 }
973 if (isFloresn) {
974 const auto floresn_global = isParallel ? this->collectOnIORank_.globalFloresn() : std::move(floresn);
975 for (const auto& flores : floresn_global) {
976 if (flores.name.empty()) {
977 continue;
978 }
979 restartValue.addExtra(flores.name, UnitSystem::measure::rate, flores.values);
980 }
981 }
982
983 std::vector<Opm::RestartValue> restartValues{};
984 // only serial, only CpGrid (for now)
985 if ( !isParallel && !needsReordering && (this->eclState_.getLgrs().size()>0) && (this->grid_.maxLevel()>0) ) {
986 // Level cells that appear on the leaf grid view get the data::Solution values from there.
987 // Other cells (i.e., parent cells that vanished due to refinement) get rubbish values for now.
988 // Only data::Solution is restricted to the level grids. Well, GroupAndNetwork, Aquifer are
989 // not modified in this method.
990 Opm::Lgr::extractRestartValueLevelGrids<Grid>(this->grid_, restartValue, restartValues);
991 }
992 else {
993 restartValues.reserve(1); // minimum size
994 restartValues.push_back(std::move(restartValue)); // no LGRs-> only one restart value
995 }
996
997 // make sure that the previous I/O request has been completed
998 // and the number of incomplete tasklets does not increase between
999 // time steps
1000 this->taskletRunner_->barrier();
1001
1002 // check if there might have been a failure in the TaskletRunner
1003 if (this->taskletRunner_->failure()) {
1004 throw std::runtime_error("Failure in the TaskletRunner while writing output.");
1005 }
1006
1007 // create a tasklet to write the data for the current time step to disk
1008 auto eclWriteTasklet = std::make_shared<EclWriteTasklet>(
1009 actionState,
1010 isParallel ? this->collectOnIORank_.globalWellTestState() : std::move(localWTestState),
1011 summaryState, udqState, *this->eclIO_,
1012 reportStepNum, timeStepNum, isSubStep, curTime, std::move(restartValues), doublePrecision,
1013 isForcedFinalOutput);
1014
1015 // finally, start a new output writing job
1016 this->taskletRunner_->dispatch(std::move(eclWriteTasklet));
1017}
1018
1019template<class Grid, class EquilGrid, class GridView, class ElementMapper, class Scalar>
1021evalSummary(const int reportStepNum,
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,
1032 const InterRegFlowMap& interRegFlows,
1033 SummaryState& summaryState,
1034 UDQState& udqState,
1035 const data::ReservoirCouplingGroupRates* rcGroupRates)
1036{
1037 if (collectOnIORank_.isIORank()) {
1038 const auto& summary = eclIO_->summary();
1039
1040 const auto& wellData = this->collectOnIORank_.isParallel()
1041 ? this->collectOnIORank_.globalWellData()
1042 : localWellData;
1043
1044 const auto& wbpData = this->collectOnIORank_.isParallel()
1045 ? this->collectOnIORank_.globalWBPData()
1046 : localWBPData;
1047
1048 const auto& groupAndNetworkData = this->collectOnIORank_.isParallel()
1049 ? this->collectOnIORank_.globalGroupAndNetworkData()
1050 : localGroupAndNetworkData;
1051
1052 const auto& aquiferData = this->collectOnIORank_.isParallel()
1053 ? this->collectOnIORank_.globalAquiferData()
1054 : localAquiferData;
1055
1056 const auto interreg_flows = getInterRegFlowsAsMap(interRegFlows);
1057
1058 const auto values = out::Summary::DynamicSimulatorState {
1059 /* well_solution = */ &wellData,
1060 /* wbp = */ &wbpData,
1061 /* group_and_nwrk_solution = */ &groupAndNetworkData,
1062 /* single_values = */ &miscSummaryData,
1063 /* region_values = */ &regionData,
1064 /* block_values = */ &blockData,
1065 /* aquifer_values = */ &aquiferData,
1066 /* interreg_flows = */ &interreg_flows,
1067 /* rc_group_rates = */ rcGroupRates,
1068 /* inplace = */ {
1069 /* current = */ &inplace,
1070 /* initial = */ initialInPlace
1071 }
1072 };
1073
1074 summary.eval(reportStepNum, curTime, values, summaryState);
1075
1076 // Off-by-one-fun: The reportStepNum argument corresponds to the
1077 // report step these results will be written to, whereas the
1078 // argument to UDQ function evaluation corresponds to the report
1079 // step we are currently on.
1080 const auto udq_step = reportStepNum - 1;
1081
1082 this->schedule_[udq_step].udq()
1083 .eval(udq_step,
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());
1090 },
1091 summaryState, udqState);
1092 }
1093
1094#if HAVE_MPI
1095 if (collectOnIORank_.isParallel()) {
1096 Parallel::MpiSerializer ser(grid_.comm());
1097 ser.append(summaryState);
1098 }
1099#endif
1100}
1101
1102template<class Grid, class EquilGrid, class GridView, class ElementMapper, class Scalar>
1105globalTrans() const
1106{
1107 assert (globalTrans_);
1108 return *globalTrans_;
1109}
1110
1111} // namespace Opm
1112
1113#endif // OPM_ECL_GENERIC_WRITER_IMPL_HPP
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 > > &regionData, 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
virtual void run()=0
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