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