opm-common
Loading...
Searching...
No Matches
Summary.hpp
1/*
2 Copyright 2016 Statoil ASA.
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 3 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
20#ifndef OPM_OUTPUT_SUMMARY_HPP
21#define OPM_OUTPUT_SUMMARY_HPP
22
23#include <opm/output/data/Aquifer.hpp>
25
26#include <map>
27#include <memory>
28#include <optional>
29#include <string>
30#include <unordered_map>
31#include <utility>
32#include <vector>
33
34namespace Opm {
35 class EclipseGrid;
36 class EclipseState;
37 class Inplace;
38 class Schedule;
39 class SummaryConfig;
40 class SummaryState;
41} // namespace Opm
42
43namespace Opm::data {
45 class InterRegFlowMap;
47 class Wells;
48} // namespace Opm::data
49
50namespace Opm::out {
51
57{
58public:
65 {
69 using GlobalProcessParameters = std::map<std::string, double>;
70
73 using RegionParameters = std::map<std::string, std::vector<double>>;
74
79 using BlockValues = std::map<std::pair<std::string, int>, double>;
80
84 using InterRegFlowValues = std::unordered_map<std::string, data::InterRegFlowMap>;
85
88 {
92 const Inplace* current {nullptr};
93
97 const Inplace* initial {nullptr};
98 };
99
103 const data::Wells* well_solution {nullptr};
104
111
117
123
128
133 const BlockValues* block_values {nullptr};
134
138 const data::Aquifers* aquifer_values {nullptr};
139
144
149 };
150
180 Summary(SummaryConfig& sumcfg,
181 const EclipseState& es,
182 const EclipseGrid& grid,
183 const Schedule& sched,
184 const std::string& basename = "",
185 const bool writeEsmry = false);
186
190 ~Summary();
191
207 void add_timestep(const SummaryState& st,
208 const int report_step,
209 const int ministep_id,
210 const bool isSubstep);
211
228 void eval(const int report_step,
229 const double secs_elapsed,
230 const DynamicSimulatorState& values,
231 SummaryState& summary_state) const;
232
239 void write(const bool is_final_summary = false) const;
240
241private:
243 class SummaryImplementation;
244
246 std::unique_ptr<SummaryImplementation> pImpl_;
247};
248
249} // namespace Opm::out
250
251#endif // OPM_OUTPUT_SUMMARY_HPP
Facility for converting collection of region ID pairs into a sparse (CSR) adjacency matrix representa...
About cell information and dimension: The actual grid information is held in a pointer to an ERT ecl_...
Definition EclipseGrid.hpp:62
Definition EclipseState.hpp:66
Definition Inplace.hpp:36
Definition Schedule.hpp:101
Collection of run's summary vectors.
Definition SummaryConfig.hpp:299
Definition SummaryState.hpp:73
Definition Groups.hpp:183
Form CSR adjacency matrix representation of inter-region flow rate graph provided as a list of connec...
Definition InterRegFlowMap.hpp:47
Definition Wells.hpp:1197
void add_timestep(const SummaryState &st, const int report_step, const int ministep_id, const bool isSubstep)
Linearise summary values into internal buffer for output purposes.
Definition Summary.cpp:5902
void write(const bool is_final_summary=false) const
Write all current summary vector buffers to output files.
Definition Summary.cpp:5910
~Summary()
Destructor.
Definition Summary.cpp:5915
Summary(SummaryConfig &sumcfg, const EclipseState &es, const EclipseGrid &grid, const Schedule &sched, const std::string &basename="", const bool writeEsmry=false)
Constructor.
Definition Summary.cpp:5874
void eval(const int report_step, const double secs_elapsed, const DynamicSimulatorState &values, SummaryState &summary_state) const
Calculate summary vector values.
Definition Summary.cpp:5883
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30
Definition Wells.hpp:1279
Volumes of fluids-in-place.
Definition Summary.hpp:88
const Inplace * current
Current fluids-in-place.
Definition Summary.hpp:92
const Inplace * initial
Initial fluids-in-place.
Definition Summary.hpp:97
Layer of indirection for transferring dynamic state objects into vector calculation engine.
Definition Summary.hpp:65
std::unordered_map< std::string, data::InterRegFlowMap > InterRegFlowValues
Collection of named inter-region flows (rates and cumulatives).
Definition Summary.hpp:84
const data::WellBlockAveragePressures * wbp
Well-block averaged pressures.
Definition Summary.hpp:110
const GlobalProcessParameters * single_values
Aggregate information about the simulation process such as number of linear and non-linear iterations...
Definition Summary.hpp:122
const data::Wells * well_solution
Dynamic state variables at the well, connection, and segment levels.
Definition Summary.hpp:103
const data::GroupAndNetworkValues * group_and_nwrk_solution
Dynamic state at the group and network levels (e.g., mode of control and node pressures).
Definition Summary.hpp:116
const InterRegFlowValues * interreg_flows
Inter-region flows (rates and cumulatives).
Definition Summary.hpp:143
std::map< std::string, std::vector< double > > RegionParameters
Collection of named per-region quantities.
Definition Summary.hpp:73
const RegionParameters * region_values
Per region dynamic state such as pressures.
Definition Summary.hpp:127
const BlockValues * block_values
Block (cell) level dynamic state values.
Definition Summary.hpp:133
std::map< std::pair< std::string, int >, double > BlockValues
Collection of per-block (cell) quantities.
Definition Summary.hpp:79
std::map< std::string, double > GlobalProcessParameters
Collection of named scalar quantities such as field-wide pressures, rates, and volumes,...
Definition Summary.hpp:69
VolumeInPlace inplace
Fluid phase volumes in place at the field and region levels.
Definition Summary.hpp:148
const data::Aquifers * aquifer_values
Aquifer level dynamic state values.
Definition Summary.hpp:138