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