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:
62 using GlobalProcessParameters = std::map<std::string, double>;
63
66 using RegionParameters = std::map<std::string, std::vector<double>>;
67
72 using BlockValues = std::map<std::pair<std::string, int>, double>;
73
77 using InterRegFlowValues = std::unordered_map<std::string, data::InterRegFlowMap>;
78
108 Summary(SummaryConfig& sumcfg,
109 const EclipseState& es,
110 const EclipseGrid& grid,
111 const Schedule& sched,
112 const std::string& basename = "",
113 const bool writeEsmry = false);
114
118 ~Summary();
119
135 void add_timestep(const SummaryState& st,
136 const int report_step,
137 const int ministep_id,
138 const bool isSubstep);
139
190 void eval(SummaryState& summary_state,
191 const int report_step,
192 const double secs_elapsed,
193 const data::Wells& well_solution,
195 const data::GroupAndNetworkValues& group_and_nwrk_solution,
196 const GlobalProcessParameters& single_values,
197 const std::optional<Inplace>& initial_inplace,
198 const Inplace& inplace,
199 const RegionParameters& region_values = {},
200 const BlockValues& block_values = {},
201 const data::Aquifers& aquifers_values = {},
202 const InterRegFlowValues& interreg_flows = {}) const;
203
210 void write(const bool is_final_summary = false) const;
211
212private:
215
217 std::unique_ptr<SummaryImplementation> pImpl_;
218};
219
220} // namespace Opm::out
221
222#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:61
Definition EclipseState.hpp:62
Definition Inplace.hpp:36
Definition Schedule.hpp:101
Definition SummaryConfig.hpp:133
Definition SummaryState.hpp:72
Definition Groups.hpp:212
Form CSR adjacency matrix representation of inter-region flow rate graph provided as a list of connec...
Definition InterRegFlowMap.hpp:47
Definition Wells.hpp:1199
Definition Summary.cpp:4951
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:5637
void write(const bool is_final_summary=false) const
Write all current summary vector buffers to output files.
Definition Summary.cpp:5645
~Summary()
Destructor.
Definition Summary.cpp:5650
std::map< std::string, std::vector< double > > RegionParameters
Collection of named per-region quantities.
Definition Summary.hpp:66
std::map< std::pair< std::string, int >, double > BlockValues
Collection of per-block (cell) quantities.
Definition Summary.hpp:72
void eval(SummaryState &summary_state, const int report_step, const double secs_elapsed, const data::Wells &well_solution, const data::WellBlockAveragePressures &wbp, const data::GroupAndNetworkValues &group_and_nwrk_solution, const GlobalProcessParameters &single_values, const std::optional< Inplace > &initial_inplace, const Inplace &inplace, const RegionParameters &region_values={}, const BlockValues &block_values={}, const data::Aquifers &aquifers_values={}, const InterRegFlowValues &interreg_flows={}) const
Calculate summary vector values.
Definition Summary.cpp:5602
std::unordered_map< std::string, data::InterRegFlowMap > InterRegFlowValues
Collection of named inter-region flows (rates and cumulatives).
Definition Summary.hpp:77
Summary(SummaryConfig &sumcfg, const EclipseState &es, const EclipseGrid &grid, const Schedule &sched, const std::string &basename="", const bool writeEsmry=false)
Constructor.
Definition Summary.cpp:5593
std::map< std::string, double > GlobalProcessParameters
Collection of named scalar quantities such as field-wide pressures, rates, and volumes,...
Definition Summary.hpp:62
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30
Definition Wells.hpp:1294