EclGenericWriter.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*/
28#ifndef OPM_ECL_GENERIC_WRITER_HPP
29#define OPM_ECL_GENERIC_WRITER_HPP
30
32
33#include <opm/output/data/Groups.hpp>
34
38
39#include <map>
40#include <memory>
41#include <optional>
42#include <string>
43#include <utility>
44#include <vector>
45
46namespace Opm {
47
48class EclipseIO;
49class EclipseState;
50class InterRegFlowMap;
51class Inplace;
52template <class Grid> class LevelCartesianIndexMapper;
53struct NNCdata;
54class Schedule;
55class SummaryConfig;
56class SummaryState;
57class UDQState;
58
59} // namespace Opm
60
61namespace Opm { namespace Action {
62class State;
63}} // namespace Opm::Action
64
65namespace Opm {
66
67template <class Grid, class EquilGrid, class GridView, class ElementMapper, class Scalar>
69{
74
75public:
76 EclGenericWriter(const Schedule& schedule,
77 const EclipseState& eclState,
78 const SummaryConfig& summaryConfig,
79 const Grid& grid,
80 const EquilGrid* equilGrid,
81 const GridView& gridView,
82 const Dune::CartesianIndexMapper<Grid>& cartMapper,
83 const Dune::CartesianIndexMapper<EquilGrid>* equilCartMapper,
84 bool enableAsyncOutput,
85 bool enableEsmry);
86
87 const EclipseIO& eclIO() const;
88
89 void writeInit();
90
92 {
94 }
95
97 {
98 sub_step_report_ = report;
99 }
101 {
102 simulation_report_ = report;
103 }
104
105 const std::vector<std::vector<NNCdata>>& getOutputNnc() const
106 {
107 return outputNnc_;
108 }
109
111 {
112 return collectOnIORank_;
113 }
114
115 void extractOutputTransAndNNC(const std::function<unsigned int(unsigned int)>& map);
116
117protected:
118 const TransmissibilityType& globalTrans() const;
119 unsigned int gridEquilIdxToGridIdx(unsigned int elemIndex) const;
120
121 void doWriteOutput(const int reportStepNum,
122 const std::optional<int> timeStepNum,
123 const bool isSubStep,
124 const bool forcedSimulationFinished,
125 data::Solution&& localCellData,
126 data::Wells&& localWellData,
127 data::GroupAndNetworkValues&& localGroupAndNetworkData,
128 data::Aquifers&& localAquiferData,
129 WellTestState&& localWTestState,
130 const Action::State& actionState,
131 const UDQState& udqState,
132 const SummaryState& summaryState,
133 const std::vector<Scalar>& thresholdPressure,
134 Scalar curTime,
135 Scalar nextStepSize,
136 bool doublePrecision,
137 bool isFlowsn,
138 std::array<FlowsData<double>,3>&& flowsn,
139 bool isFloresn,
140 std::array<FlowsData<double>, 3>&& floresn);
141
142 void evalSummary(int reportStepNum,
143 Scalar curTime,
144 const data::Wells& localWellData,
145 const data::WellBlockAveragePressures& localWBPData,
146 const data::GroupAndNetworkValues& localGroupAndNetworkData,
147 const std::map<int,data::AquiferData>& localAquiferData,
148 const std::map<std::pair<std::string, int>, double>& blockData,
149 const std::map<std::string, double>& miscSummaryData,
150 const std::map<std::string, std::vector<double>>& regionData,
151 const Inplace& inplace,
152 const Inplace* initialInPlace,
153 const InterRegFlowMap& interRegFlows,
154 SummaryState& summaryState,
155 UDQState& udqState,
156 const data::ReservoirCouplingGroupRates* rcGroupRates = nullptr);
157
159 const Grid& grid_;
160 const GridView& gridView_;
161 const Schedule& schedule_;
162 const EclipseState& eclState_;
163 std::unique_ptr<EclipseIO> eclIO_;
164 std::unique_ptr<TaskletRunner> taskletRunner_;
169 const EquilGrid* equilGrid_;
172
173 // Regular NNCs per grid: internal to a grid.
174 // Both cells belong to the same level grid, either the main grid or a level/local grid.
175 // nnc.cell1 (NNC1 in *.EGRID) Level/Local Cartesian index of cell1
176 // nnc.cell2 (NNC2 in *.EGRID) Level/Local Cartesian index of cell2
177 // Equivalent to TRANNNC in *.INIT
178 mutable std::vector<std::vector<NNCdata>> outputNnc_;
179
180 // NNCs between main (level zero) grid and local grids:
181 // nnc.cell1 (NNCG in *.EGRID) Cartesian index of cell1 in the main grid where the cell belongs to.
182 // nnc.cell2 (NNCL in *.EGRID) Level/Local Cartesian index of cell2 in the refined level grid where the cell belongs to.
183 // Equivalent to TRANGL in *.INIT
184 mutable std::vector<std::vector<NNCdata>> outputNncGlobalLocal_; // here GLOBAL refers to level 0 grid, local to any LGR (refined grid)
185
186 // Amalgamated NNCs: nncs between different LGRs. For example, nested refinement or neighboring LGRs.
187 // The cells belong to different refined level grids.
188 // nnc.cell1 (NNA1 in *.EGRID) Level/Local Cartesian index of cell1 (in its level grid: level1)
189 // nnc.cell2 (NNA2 in *.EGRID) Level/Local Cartesian index of cell2 (in its level grid: level2, with level2 > level1).
190 // Equivalent to TRANLL in *.INIT
191 mutable std::vector<std::vector<std::vector<NNCdata>>> outputAmalgamatedNnc_;
192
193 mutable std::unique_ptr<std::vector<data::Solution>> outputTrans_;
194
195private:
196 template<typename LevelIndicesFunction, typename OriginIndicesFunction>
197 void computeTrans_(const std::vector<std::unordered_map<int,int>>& levelCartToLevelCompressed,
198 const std::function<unsigned int(unsigned int)>& map,
199 const LevelIndicesFunction& computeLevelIndices,
200 const std::function<int(int, int)>& computeLevelCartIdx,
201 const std::function<std::array<int,3>(int)>& computeLevelCartDimensions,
202 const OriginIndicesFunction& computeOriginIndices) const;
203
204 template<typename LevelIndicesFunction, typename OriginIndicesFunction>
205 std::vector<std::vector<NNCdata>> exportNncStructure_(const std::vector<std::unordered_map<int,int>>& levelCartToLevelCompressed,
206 const std::function<unsigned int(unsigned int)>& map,
207 const LevelIndicesFunction& computeLevelIndices,
208 const std::function<int(int, int)>& computeLevelCartIdx,
209 const std::function<std::array<int,3>(int)>& computeLevelCartDimensions,
210 const OriginIndicesFunction& computeOriginIndices) const;
211
213 bool isCartesianNeighbour_(const std::array<int,3>& levelCartDims,
214 const std::size_t levelCartIdx1,
215 const std::size_t levelCartIdx2) const;
216
217 bool isDirectNeighbours_(const std::unordered_map<int,int>& levelCartesianToActive,
218 const std::array<int,3>& levelCartDims,
219 const std::size_t levelCartIdx1,
220 const std::size_t levelCartIdx2) const;
221
222 auto activeCell_(const std::unordered_map<int,int>& levelCartToLevelCompressed,
223 const std::size_t levelCartIdx) const;
224
225
226
229 bool isNumAquCell_(const std::size_t cartIdx) const;
231 bool isNumAquConn_(const std::size_t cartIdx1, const std::size_t cartIdx2) const;
232
241 template<bool equilGridIsCpGrid>
242 Opm::LevelCartesianIndexMapper<EquilGrid> createLevelCartMapp_() const;
243
251 template<bool equilGridIsCpGrid>
252 std::vector<std::unordered_map<int,int>> createCartesianToActiveMaps_(const Opm::LevelCartesianIndexMapper<EquilGrid>& levelCartMapp) const;
253
259 template<bool equilGridIsCpGrid>
260 std::function<std::array<int,3>(int)> computeLevelCartDimensions_(const Opm::LevelCartesianIndexMapper<EquilGrid>& levelCartMapp,
261 const Dune::CartesianIndexMapper<EquilGrid>& equilCartMapp) const;
262
267 template<bool equilGridIsCpGrid>
268 std::function<int(int, int)> computeLevelCartIdx_(const Opm::LevelCartesianIndexMapper<EquilGrid>& levelCartMapp,
269 const Dune::CartesianIndexMapper<EquilGrid>& equilCartMapp) const;
270
275 template <bool equilGridIsCpGrid>
276 auto computeLevelIndices_() const;
277
282 template <bool equilGridIsCpGrid>
283 auto computeOriginIndices_() const;
284
287 void allocateLevelTrans_(const std::array<int,3>& levelCartDims,
288 data::Solution& levelTrans) const;
294 void allocateAllNncs_(int maxLevel) const;
295};
296
297} // namespace Opm
298
299#endif // OPM_ECL_GENERIC_WRITER_HPP
Definition: CollectDataOnIORank.hpp:49
Definition: CollectDataOnIORank.hpp:56
Definition: EclGenericWriter.hpp:69
std::vector< std::vector< NNCdata > > outputNnc_
Definition: EclGenericWriter.hpp:178
const Grid & grid_
Definition: EclGenericWriter.hpp:159
const EclipseState & eclState_
Definition: EclGenericWriter.hpp:162
void setSubStepReport(const SimulatorReportSingle &report)
Definition: EclGenericWriter.hpp:96
void setSimulationReport(const SimulatorReport &report)
Definition: EclGenericWriter.hpp:100
std::vector< std::vector< NNCdata > > outputNncGlobalLocal_
Definition: EclGenericWriter.hpp:184
std::unique_ptr< std::vector< data::Solution > > outputTrans_
Definition: EclGenericWriter.hpp:193
unsigned int gridEquilIdxToGridIdx(unsigned int elemIndex) const
CollectDataOnIORankType collectOnIORank_
Definition: EclGenericWriter.hpp:158
const Schedule & schedule_
Definition: EclGenericWriter.hpp:161
SimulatorReport simulation_report_
Definition: EclGenericWriter.hpp:171
void setTransmissibilities(const TransmissibilityType *globalTrans)
Definition: EclGenericWriter.hpp:91
SimulatorReportSingle sub_step_report_
Definition: EclGenericWriter.hpp:170
const EquilGrid * equilGrid_
Definition: EclGenericWriter.hpp:169
void extractOutputTransAndNNC(const std::function< unsigned int(unsigned int)> &map)
Definition: EclGenericWriter_impl.hpp:303
Scalar restartTimeStepSize_
Definition: EclGenericWriter.hpp:165
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::vector< std::vector< std::vector< NNCdata > > > outputAmalgamatedNnc_
Definition: EclGenericWriter.hpp:191
std::unique_ptr< TaskletRunner > taskletRunner_
Definition: EclGenericWriter.hpp:164
void writeInit()
Definition: EclGenericWriter_impl.hpp:275
const Dune::CartesianIndexMapper< Grid > & cartMapper_
Definition: EclGenericWriter.hpp:167
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 GridView & gridView_
Definition: EclGenericWriter.hpp:160
const std::vector< std::vector< NNCdata > > & getOutputNnc() const
Definition: EclGenericWriter.hpp:105
const CollectDataOnIORankType & collectOnIORank() const
Definition: EclGenericWriter.hpp:110
const TransmissibilityType * globalTrans_
Definition: EclGenericWriter.hpp:166
const TransmissibilityType & globalTrans() const
Definition: EclGenericWriter_impl.hpp:1105
std::unique_ptr< EclipseIO > eclIO_
Definition: EclGenericWriter.hpp:163
const Dune::CartesianIndexMapper< EquilGrid > * equilCartMapper_
Definition: EclGenericWriter.hpp:168
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
Definition: EclGenericWriter.hpp:52
Definition: Transmissibility.hpp:54
Definition: blackoilbioeffectsmodules.hh:45
Definition: SimulatorReport.hpp:122
A struct for returning timing data from a simulator to its caller.
Definition: SimulatorReport.hpp:34
Provides a mechanism to dispatch work to separate threads.