My Project
AggregateAquiferData.hpp
1 /*
2  Copyright (c) 2021 Equinor 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_AGGREGATE_AQUIFER_DATA_HPP
21 #define OPM_AGGREGATE_AQUIFER_DATA_HPP
22 
23 #include <opm/output/eclipse/InteHEAD.hpp>
25 
26 #include <opm/output/data/Aquifer.hpp>
27 
28 #include <vector>
29 
30 namespace Opm {
31  class AquiferConfig;
32  class EclipseGrid;
33  class SummaryState;
34  class UnitSystem;
35 } // Opm
36 
37 namespace Opm { namespace RestartIO { namespace Helpers {
38 
40  {
41  public:
56  const AquiferConfig& aqConfig,
57  const EclipseGrid& grid);
58 
78  const data::Aquifers& aquData,
79  const SummaryState& summaryState,
80  const UnitSystem& usys);
81 
88  {
89  return this->maxActiveAnalyticAquiferID_;
90  }
91 
93  const std::vector<int>& getIntegerAquiferData() const
94  {
95  return this->integerAnalyticAq_.data();
96  }
97 
99  const std::vector<float>& getSinglePrecAquiferData() const
100  {
101  return this->singleprecAnalyticAq_.data();
102  }
103 
105  const std::vector<double>& getDoublePrecAquiferData() const
106  {
107  return this->doubleprecAnalyticAq_.data();
108  }
109 
111  const std::vector<int>& getNumericAquiferIntegerData() const
112  {
113  return this->integerNumericAq_.data();
114  }
115 
117  const std::vector<double>& getNumericAquiferDoublePrecData() const
118  {
119  return this->doubleprecNumericAq_.data();
120  }
121 
127  const std::vector<int>& getIntegerAquiferConnectionData(const int aquiferID) const
128  {
129  return this->integerAnalyticAquiferConn_[aquiferID - 1].data();
130  }
131 
138  const std::vector<float>& getSinglePrecAquiferConnectionData(const int aquiferID) const
139  {
140  return this->singleprecAnalyticAquiferConn_[aquiferID - 1].data();
141  }
142 
150  const std::vector<double>& getDoublePrecAquiferConnectionData(const int aquiferID) const
151  {
152  return this->doubleprecAnalyticAquiferConn_[aquiferID - 1].data();
153  }
154 
155  private:
156  int maxActiveAnalyticAquiferID_{0};
157 
158  std::vector<int> numActiveConn_{};
159  std::vector<double> totalInflux_{};
160 
162  WindowedArray<int> integerAnalyticAq_;
163 
165  WindowedArray<float> singleprecAnalyticAq_;
166 
168  WindowedArray<double> doubleprecAnalyticAq_;
169 
171  WindowedArray<int> integerNumericAq_;
172 
174  WindowedArray<double> doubleprecNumericAq_;
175 
178  std::vector<WindowedArray<int>> integerAnalyticAquiferConn_;
179 
182  std::vector<WindowedArray<float>> singleprecAnalyticAquiferConn_;
183 
186  std::vector<WindowedArray<double>> doubleprecAnalyticAquiferConn_;
187  };
188 
189 }}} // Opm::RestartIO::Helpers
190 
191 #endif // OPM_AGGREGATE_WELL_DATA_HPP
Provide facilities to simplify constructing restart vectors such as IWEL or RSEG.
Definition: AquiferConfig.hpp:44
About cell information and dimension: The actual grid information is held in a pointer to an ERT ecl_...
Definition: EclipseGrid.hpp:54
Definition: AggregateAquiferData.hpp:40
AggregateAquiferData(const InteHEAD::AquiferDims &aqDims, const AquiferConfig &aqConfig, const EclipseGrid &grid)
Constructor.
const std::vector< double > & getNumericAquiferDoublePrecData() const
Retrieve Double Precision Aquifer Data Array for Numeric Aquifers.
Definition: AggregateAquiferData.hpp:117
const std::vector< double > & getDoublePrecAquiferData() const
Retrieve Floating-Point (Double Precision) Aquifer Data Array.
Definition: AggregateAquiferData.hpp:105
const std::vector< int > & getIntegerAquiferConnectionData(const int aquiferID) const
Retrieve Integer Aquifer Connection Data Array (analytic aquifers)
Definition: AggregateAquiferData.hpp:127
int maximumActiveAnalyticAquiferID() const
Retrieve the maximum active aquifer ID over all analytic aquifers.
Definition: AggregateAquiferData.hpp:87
void captureDynamicdAquiferData(const AquiferConfig &aqConfig, const data::Aquifers &aquData, const SummaryState &summaryState, const UnitSystem &usys)
Linearise dynamic information pertinent to analytic aquifers into internal arrays.
const std::vector< double > & getDoublePrecAquiferConnectionData(const int aquiferID) const
Retrieve Floating-Point (Double Precision) Aquifer Connection Data Array (analytic aquifers)
Definition: AggregateAquiferData.hpp:150
const std::vector< int > & getIntegerAquiferData() const
Retrieve Integer Aquifer Data Array.
Definition: AggregateAquiferData.hpp:93
const std::vector< float > & getSinglePrecAquiferConnectionData(const int aquiferID) const
Retrieve Floating-Point (Real) Aquifer Connection Data Array (analytic aquifers)
Definition: AggregateAquiferData.hpp:138
const std::vector< int > & getNumericAquiferIntegerData() const
Retrieve Integer Aquifer Data Array for Numeric Aquifers.
Definition: AggregateAquiferData.hpp:111
const std::vector< float > & getSinglePrecAquiferData() const
Retrieve Floating-Point (Real) Aquifer Data Array.
Definition: AggregateAquiferData.hpp:99
const std::vector< T > & data() const
Get read-only access to full, linearised data items for all windows.
Definition: WindowedArray.hpp:131
Definition: SummaryState.hpp:69
Definition: UnitSystem.hpp:33
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition: Exceptions.hpp:29
Definition: InteHEAD.hpp:149