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;
50struct NNCdata;
51class Schedule;
52class SummaryConfig;
53class SummaryState;
54class UDQState;
55
56} // namespace Opm
57
58namespace Opm { namespace Action {
59class State;
60}} // namespace Opm::Action
61
62namespace Opm {
63
64template <class Grid, class EquilGrid, class GridView, class ElementMapper, class Scalar>
66{
71
72public:
73 EclGenericWriter(const Schedule& schedule,
74 const EclipseState& eclState,
75 const SummaryConfig& summaryConfig,
76 const Grid& grid,
77 const EquilGrid* equilGrid,
78 const GridView& gridView,
79 const Dune::CartesianIndexMapper<Grid>& cartMapper,
80 const Dune::CartesianIndexMapper<EquilGrid>* equilCartMapper,
81 bool enableAsyncOutput,
82 bool enableEsmry);
83
84 const EclipseIO& eclIO() const;
85
86 void writeInit();
87
89 {
91 }
92
94 {
95 sub_step_report_ = report;
96 }
98 {
99 simulation_report_ = report;
100 }
101
102 const std::vector<NNCdata>& getOutputNnc() const
103 {
104 return outputNnc_;
105 }
106
108 {
109 return collectOnIORank_;
110 }
111
112 void extractOutputTransAndNNC(const std::function<unsigned int(unsigned int)>& map);
113
114protected:
115 const TransmissibilityType& globalTrans() const;
116 unsigned int gridEquilIdxToGridIdx(unsigned int elemIndex) const;
117
118 void doWriteOutput(const int reportStepNum,
119 const std::optional<int> timeStepNum,
120 const bool isSubStep,
121 data::Solution&& localCellData,
122 data::Wells&& localWellData,
123 data::GroupAndNetworkValues&& localGroupAndNetworkData,
124 data::Aquifers&& localAquiferData,
125 WellTestState&& localWTestState,
126 const Action::State& actionState,
127 const UDQState& udqState,
128 const SummaryState& summaryState,
129 const std::vector<Scalar>& thresholdPressure,
130 Scalar curTime,
131 Scalar nextStepSize,
132 bool doublePrecision,
133 bool isFlowsn,
134 std::array<FlowsData<double>,3>&& flowsn,
135 bool isFloresn,
136 std::array<FlowsData<double>, 3>&& floresn);
137
138 void evalSummary(int reportStepNum,
139 Scalar curTime,
140 const data::Wells& localWellData,
141 const data::WellBlockAveragePressures& localWBPData,
142 const data::GroupAndNetworkValues& localGroupAndNetworkData,
143 const std::map<int,data::AquiferData>& localAquiferData,
144 const std::map<std::pair<std::string, int>, double>& blockData,
145 const std::map<std::string, double>& miscSummaryData,
146 const std::map<std::string, std::vector<double>>& regionData,
147 const Inplace& inplace,
148 const Inplace* initialInPlace,
149 const InterRegFlowMap& interRegFlows,
150 SummaryState& summaryState,
151 UDQState& udqState);
152
154 const Grid& grid_;
155 const GridView& gridView_;
156 const Schedule& schedule_;
157 const EclipseState& eclState_;
158 std::unique_ptr<EclipseIO> eclIO_;
159 std::unique_ptr<TaskletRunner> taskletRunner_;
164 const EquilGrid* equilGrid_;
167 mutable std::vector<NNCdata> outputNnc_;
168 mutable std::unique_ptr<data::Solution> outputTrans_;
169
170private:
171 void computeTrans_(const std::unordered_map<int,int>& cartesianToActive, const std::function<unsigned int(unsigned int)>& map) const;
172 std::vector<NNCdata> exportNncStructure_(const std::unordered_map<int,int>& cartesianToActive, const std::function<unsigned int(unsigned int)>& map) const;
173};
174
175} // namespace Opm
176
177#endif // OPM_ECL_GENERIC_WRITER_HPP
Definition: CollectDataOnIORank.hpp:49
Definition: CollectDataOnIORank.hpp:56
Definition: EclGenericWriter.hpp:66
const Grid & grid_
Definition: EclGenericWriter.hpp:154
const EclipseState & eclState_
Definition: EclGenericWriter.hpp:157
void setSubStepReport(const SimulatorReportSingle &report)
Definition: EclGenericWriter.hpp:93
void setSimulationReport(const SimulatorReport &report)
Definition: EclGenericWriter.hpp:97
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:553
unsigned int gridEquilIdxToGridIdx(unsigned int elemIndex) const
CollectDataOnIORankType collectOnIORank_
Definition: EclGenericWriter.hpp:153
const Schedule & schedule_
Definition: EclGenericWriter.hpp:156
std::unique_ptr< data::Solution > outputTrans_
Definition: EclGenericWriter.hpp:168
SimulatorReport simulation_report_
Definition: EclGenericWriter.hpp:166
void setTransmissibilities(const TransmissibilityType *globalTrans)
Definition: EclGenericWriter.hpp:88
SimulatorReportSingle sub_step_report_
Definition: EclGenericWriter.hpp:165
const EquilGrid * equilGrid_
Definition: EclGenericWriter.hpp:164
void extractOutputTransAndNNC(const std::function< unsigned int(unsigned int)> &map)
Definition: EclGenericWriter_impl.hpp:282
Scalar restartTimeStepSize_
Definition: EclGenericWriter.hpp:160
std::unique_ptr< TaskletRunner > taskletRunner_
Definition: EclGenericWriter.hpp:159
void writeInit()
Definition: EclGenericWriter_impl.hpp:265
std::vector< NNCdata > outputNnc_
Definition: EclGenericWriter.hpp:167
const Dune::CartesianIndexMapper< Grid > & cartMapper_
Definition: EclGenericWriter.hpp:162
const std::vector< NNCdata > & getOutputNnc() const
Definition: EclGenericWriter.hpp:102
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:257
const GridView & gridView_
Definition: EclGenericWriter.hpp:155
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:661
const CollectDataOnIORankType & collectOnIORank() const
Definition: EclGenericWriter.hpp:107
const TransmissibilityType * globalTrans_
Definition: EclGenericWriter.hpp:161
const TransmissibilityType & globalTrans() const
Definition: EclGenericWriter_impl.hpp:743
std::unique_ptr< EclipseIO > eclIO_
Definition: EclGenericWriter.hpp:158
const Dune::CartesianIndexMapper< EquilGrid > * equilCartMapper_
Definition: EclGenericWriter.hpp:163
Inter-region flow accumulation maps for all region definition arrays.
Definition: InterRegFlows.hpp:179
Definition: Transmissibility.hpp:54
Definition: blackoilbioeffectsmodules.hh:43
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.