opm-simulators
Loading...
Searching...
No Matches
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/simulators/flow/CollectDataOnIORank.hpp>
35#include <opm/simulators/timestepping/SimulatorReport.hpp>
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>
66class EclGenericWriter
67{
68 using CartesianIndexMapper = Dune::CartesianIndexMapper<Grid>;
69 using EquilCartesianIndexMapper = Dune::CartesianIndexMapper<EquilGrid>;
70 using CollectDataOnIORankType = CollectDataOnIORank<Grid,EquilGrid,GridView>;
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
89 void setTransmissibilities(const TransmissibilityType* globalTrans)
90 {
91 globalTrans_ = globalTrans;
92 }
93
94 void setSubStepReport(const SimulatorReportSingle& report)
95 {
96 sub_step_report_ = report;
97 }
98 void setSimulationReport(const SimulatorReport& report)
99 {
100 simulation_report_ = report;
101 }
102
103 const std::vector<std::vector<NNCdata>>& getOutputNnc() const
104 {
105 return outputNnc_;
106 }
107
108 const CollectDataOnIORankType& collectOnIORank() const
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 const bool forcedSimulationFinished,
123 data::Solution&& localCellData,
124 data::Wells&& localWellData,
125 data::GroupAndNetworkValues&& localGroupAndNetworkData,
126 data::Aquifers&& localAquiferData,
127 WellTestState&& localWTestState,
128 const Action::State& actionState,
129 const UDQState& udqState,
130 const SummaryState& summaryState,
131 const std::vector<Scalar>& thresholdPressure,
132 Scalar curTime,
133 Scalar nextStepSize,
134 bool doublePrecision,
135 bool isFlowsn,
136 std::array<FlowsData<double>,3>&& flowsn,
137 bool isFloresn,
138 std::array<FlowsData<double>, 3>&& floresn);
139
140 void evalSummary(int reportStepNum,
141 Scalar curTime,
142 const data::Wells& localWellData,
143 const data::WellBlockAveragePressures& localWBPData,
144 const data::GroupAndNetworkValues& localGroupAndNetworkData,
145 const std::map<int,data::AquiferData>& localAquiferData,
146 const std::map<std::pair<std::string, int>, double>& blockData,
147 const std::map<std::string, double>& miscSummaryData,
148 const std::map<std::string, std::vector<double>>& regionData,
149 const Inplace& inplace,
150 const Inplace* initialInPlace,
151 const InterRegFlowMap& interRegFlows,
152 SummaryState& summaryState,
153 UDQState& udqState);
154
155 CollectDataOnIORankType collectOnIORank_;
156 const Grid& grid_;
157 const GridView& gridView_;
158 const Schedule& schedule_;
159 const EclipseState& eclState_;
160 std::unique_ptr<EclipseIO> eclIO_;
161 std::unique_ptr<TaskletRunner> taskletRunner_;
162 Scalar restartTimeStepSize_;
163 const TransmissibilityType* globalTrans_ = nullptr;
164 const Dune::CartesianIndexMapper<Grid>& cartMapper_;
165 const Dune::CartesianIndexMapper<EquilGrid>* equilCartMapper_;
166 const EquilGrid* equilGrid_;
167 SimulatorReportSingle sub_step_report_;
168 SimulatorReport simulation_report_;
169
170 // Regular NNCs per grid: internal to a grid.
171 // Both cells belong to the same level grid, either the main grid or a level/local grid.
172 // nnc.cell1 (NNC1 in *.EGRID) Level/Local Cartesian index of cell1
173 // nnc.cell2 (NNC2 in *.EGRID) Level/Local Cartesian index of cell2
174 // Equivalent to TRANNNC in *.INIT
175 mutable std::vector<std::vector<NNCdata>> outputNnc_;
176
177 // NNCs between main (level zero) grid and local grids:
178 // nnc.cell1 (NNCG in *.EGRID) Cartesian index of cell1 in the main grid where the cell belongs to.
179 // nnc.cell2 (NNCL in *.EGRID) Level/Local Cartesian index of cell2 in the refined level grid where the cell belongs to.
180 // Equivalent to TRANGL in *.INIT
181 mutable std::vector<std::vector<NNCdata>> outputNncGlobalLocal_; // here GLOBAL refers to level 0 grid, local to any LGR (refined grid)
182
183 // Amalgamated NNCs: nncs between different LGRs. For example, nested refinement or neighboring LGRs.
184 // The cells belong to different refined level grids.
185 // nnc.cell1 (NNA1 in *.EGRID) Level/Local Cartesian index of cell1 (in its level grid: level1)
186 // nnc.cell2 (NNA2 in *.EGRID) Level/Local Cartesian index of cell2 (in its level grid: level2, with level2 > level1).
187 // Equivalent to TRANLL in *.INIT
188 mutable std::vector<std::vector<std::vector<NNCdata>>> outputAmalgamatedNnc_;
189
190 mutable std::unique_ptr<std::vector<data::Solution>> outputTrans_;
191
192private:
193 template<typename LevelIndicesFunction, typename OriginIndicesFunction>
194 void computeTrans_(const std::vector<std::unordered_map<int,int>>& levelCartToLevelCompressed,
195 const std::function<unsigned int(unsigned int)>& map,
196 const LevelIndicesFunction& computeLevelIndices,
197 const std::function<int(int, int)>& computeLevelCartIdx,
198 const std::function<std::array<int,3>(int)>& computeLevelCartDimensions,
199 const OriginIndicesFunction& computeOriginIndices) const;
200
201 template<typename LevelIndicesFunction, typename OriginIndicesFunction>
202 std::vector<std::vector<NNCdata>> exportNncStructure_(const std::vector<std::unordered_map<int,int>>& levelCartToLevelCompressed,
203 const std::function<unsigned int(unsigned int)>& map,
204 const LevelIndicesFunction& computeLevelIndices,
205 const std::function<int(int, int)>& computeLevelCartIdx,
206 const std::function<std::array<int,3>(int)>& computeLevelCartDimensions,
207 const OriginIndicesFunction& computeOriginIndices) const;
208
210 bool isCartesianNeighbour_(const std::array<int,3>& levelCartDims,
211 const std::size_t levelCartIdx1,
212 const std::size_t levelCartIdx2) const;
213
214 bool isDirectNeighbours_(const std::unordered_map<int,int>& levelCartesianToActive,
215 const std::array<int,3>& levelCartDims,
216 const std::size_t levelCartIdx1,
217 const std::size_t levelCartIdx2) const;
218
219 auto activeCell_(const std::unordered_map<int,int>& levelCartToLevelCompressed,
220 const std::size_t levelCartIdx) const;
221
222
223
226 bool isNumAquCell_(const std::size_t cartIdx) const;
228 bool isNumAquConn_(const std::size_t cartIdx1, const std::size_t cartIdx2) const;
229
238 template<bool equilGridIsCpGrid>
239 Opm::LevelCartesianIndexMapper<EquilGrid> createLevelCartMapp_() const;
240
248 template<bool equilGridIsCpGrid>
249 std::vector<std::unordered_map<int,int>> createCartesianToActiveMaps_(const Opm::LevelCartesianIndexMapper<EquilGrid>& levelCartMapp) const;
250
256 template<bool equilGridIsCpGrid>
257 std::function<std::array<int,3>(int)> computeLevelCartDimensions_(const Opm::LevelCartesianIndexMapper<EquilGrid>& levelCartMapp,
258 const Dune::CartesianIndexMapper<EquilGrid>& equilCartMapp) const;
259
264 template<bool equilGridIsCpGrid>
265 std::function<int(int, int)> computeLevelCartIdx_(const Opm::LevelCartesianIndexMapper<EquilGrid>& levelCartMapp,
266 const Dune::CartesianIndexMapper<EquilGrid>& equilCartMapp) const;
267
272 template <bool equilGridIsCpGrid>
273 auto computeLevelIndices_() const;
274
279 template <bool equilGridIsCpGrid>
280 auto computeOriginIndices_() const;
281
284 void allocateLevelTrans_(const std::array<int,3>& levelCartDims,
285 data::Solution& levelTrans) const;
291 void allocateAllNncs_(int maxLevel) const;
292};
293
294} // namespace Opm
295
296#endif // OPM_ECL_GENERIC_WRITER_HPP
Definition CollectDataOnIORank.hpp:49
Definition CollectDataOnIORank.hpp:56
Inter-region flow accumulation maps for all region definition arrays.
Definition InterRegFlows.hpp:179
Definition EclGenericWriter.hpp:50
Definition Transmissibility.hpp:54
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:45
Simple container for FLOWS data.
Definition FlowsData.hpp:32
A struct for returning timing data from a simulator to its caller.
Definition SimulatorReport.hpp:34
Definition SimulatorReport.hpp:122
Provides a mechanism to dispatch work to separate threads.