WellState.hpp
1 /*
2  Copyright 2012 SINTEF ICT, Applied Mathematics.
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_WELLSTATE_HEADER_INCLUDED
21 #define OPM_WELLSTATE_HEADER_INCLUDED
22 
23 #include <opm/core/props/BlackoilPhases.hpp>
24 #include <opm/core/wells.h>
25 #include <opm/core/well_controls.h>
26 #include <opm/output/data/Wells.hpp>
27 
28 #include <array>
29 #include <map>
30 #include <memory>
31 #include <string>
32 #include <vector>
33 #include <cassert>
34 #include <cstddef>
35 
36 namespace Opm
37 {
39  class WellState
40  {
41  public:
42  typedef std::array< int, 3 > mapentry_t;
43  typedef std::map< std::string, mapentry_t > WellMapType;
44 
45  template <class State>
46  void init(const Wells* wells, const State& state)
47  {
48  init(wells, state.pressure());
49  }
50 
56  void init(const Wells* wells, const std::vector<double>& cellPressures)
57  {
58  // clear old name mapping
59  wellMap_.clear();
60  wells_.reset( clone_wells( wells ) );
61 
62  if (wells) {
63  const int nw = wells->number_of_wells;
64  const int np = wells->number_of_phases;
65  bhp_.resize(nw);
66  thp_.resize(nw);
67  temperature_.resize(nw, 273.15 + 20); // standard temperature for now
68  wellrates_.resize(nw * np, 0.0);
69  for (int w = 0; w < nw; ++w) {
70  assert((wells->type[w] == INJECTOR) || (wells->type[w] == PRODUCER));
71  const WellControls* ctrl = wells->ctrls[w];
72  const int num_perf_this_well = wells->well_connpos[w + 1] - wells->well_connpos[w];
73 
74  // setup wellname -> well index mapping
75  {
76  assert( wells->name[ w ] );
77  std::string name( wells->name[ w ] );
78  assert( name.size() > 0 );
79  mapentry_t& wellMapEntry = wellMap_[name];
80  wellMapEntry[ 0 ] = w;
81  wellMapEntry[ 1 ] = wells->well_connpos[w];
82  // also store the number of perforations in this well
83  wellMapEntry[ 2 ] = num_perf_this_well;
84  }
85 
86  if ( num_perf_this_well == 0 )
87  {
88  // No perforations of the well. Initialize to zero.
89  for (int p = 0; p < np; ++p) {
90  wellrates_[np*w + p] = 0.0;
91  }
92  bhp_[w] = 0.;
93  thp_[w] = 0.;
94  continue;
95  }
96 
97  if (well_controls_well_is_stopped(ctrl)) {
98  // Stopped well:
99  // 1. Rates: assign zero well rates.
100  for (int p = 0; p < np; ++p) {
101  wellrates_[np*w + p] = 0.0;
102  }
103  // 2. Bhp: assign bhp equal to bhp control, if
104  // applicable, otherwise assign equal to
105  // first perforation cell pressure.
106  if (well_controls_get_current_type(ctrl) == BHP) {
107  bhp_[w] = well_controls_get_current_target( ctrl );
108  } else {
109  const int first_cell = wells->well_cells[wells->well_connpos[w]];
110  bhp_[w] = cellPressures[first_cell];
111  }
112  } else {
113  // Open well:
114  // 1. Rates: initialize well rates to match controls
115  // if type is SURFACE_RATE. Otherwise, we
116  // cannot set the correct value here, so we
117  // assign a small rate with the correct
118  // sign so that any logic depending on that
119  // sign will work as expected.
120  if (well_controls_get_current_type(ctrl) == SURFACE_RATE) {
121  const double rate_target = well_controls_get_current_target(ctrl);
122  const double * distr = well_controls_get_current_distr( ctrl );
123  for (int p = 0; p < np; ++p) {
124  wellrates_[np*w + p] = rate_target * distr[p];
125  }
126  } else {
127  const double small_rate = 1e-14;
128  const double sign = (wells->type[w] == INJECTOR) ? 1.0 : -1.0;
129  for (int p = 0; p < np; ++p) {
130  wellrates_[np*w + p] = small_rate * sign;
131  }
132  }
133 
134  // 2. Bhp: initialize bhp to be target pressure if
135  // bhp-controlled well, otherwise set to a
136  // little above or below (depending on if
137  // the well is an injector or producer)
138  // pressure in first perforation cell.
139  if (well_controls_get_current_type(ctrl) == BHP) {
140  bhp_[w] = well_controls_get_current_target( ctrl );
141  } else {
142  const int first_cell = wells->well_cells[wells->well_connpos[w]];
143  const double safety_factor = (wells->type[w] == INJECTOR) ? 1.01 : 0.99;
144  bhp_[w] = safety_factor*cellPressures[first_cell];
145  }
146  }
147 
148  // 3. Thp: assign thp equal to thp target/limit, if applicable,
149  // otherwise keep it zero. Basically, the value should not be used
150  // in the simulation at all.
151  const int nwc = well_controls_get_num(ctrl);
152  for (int ctrl_index = 0; ctrl_index < nwc; ++ctrl_index) {
153  if (well_controls_iget_type(ctrl, ctrl_index) == THP) {
154  thp_[w] = well_controls_iget_target(ctrl, ctrl_index);
155  break;
156  }
157  }
158  }
159 
160  // The perforation rates and perforation pressures are
161  // not expected to be consistent with bhp_ and wellrates_
162  // after init().
163  perfrates_.resize(wells->well_connpos[nw], 0.0);
164  perfpress_.resize(wells->well_connpos[nw], -1e100);
165  }
166  }
167 
169  std::vector<double>& bhp() { return bhp_; }
170  const std::vector<double>& bhp() const { return bhp_; }
171 
173  std::vector<double>& thp() { return thp_; }
174  const std::vector<double>& thp() const { return thp_; }
175 
177  std::vector<double>& temperature() { return temperature_; }
178  const std::vector<double>& temperature() const { return temperature_; }
179 
181  std::vector<double>& wellRates() { return wellrates_; }
182  const std::vector<double>& wellRates() const { return wellrates_; }
183 
185  std::vector<double>& perfRates() { return perfrates_; }
186  const std::vector<double>& perfRates() const { return perfrates_; }
187 
189  std::vector<double>& perfPress() { return perfpress_; }
190  const std::vector<double>& perfPress() const { return perfpress_; }
191 
192  size_t getRestartBhpOffset() const {
193  return 0;
194  }
195 
196  size_t getRestartPerfPressOffset() const {
197  return bhp_.size();
198  }
199 
200  size_t getRestartPerfRatesOffset() const {
201  return getRestartPerfPressOffset() + perfpress_.size();
202  }
203 
204  size_t getRestartTemperatureOffset() const {
205  return getRestartPerfRatesOffset() + perfrates_.size();
206  }
207 
208  size_t getRestartWellRatesOffset() const {
209  return getRestartTemperatureOffset() + temperature_.size();
210  }
211 
212  const WellMapType& wellMap() const { return wellMap_; }
213  WellMapType& wellMap() { return wellMap_; }
214 
216  int numWells() const
217  {
218  return bhp().size();
219  }
220 
222  int numPhases() const
223  {
224  return wellRates().size() / numWells();
225  }
226 
227  virtual data::Wells report(const PhaseUsage& pu) const
228  {
229  using rt = data::Rates::opt;
230 
231  data::Wells dw;
232  for( const auto& itr : this->wellMap_ ) {
233  const auto well_index = itr.second[ 0 ];
234 
235  auto& well = dw[ itr.first ];
236  well.bhp = this->bhp().at( well_index );
237  well.thp = this->thp().at( well_index );
238  well.temperature = this->temperature().at( well_index );
239 
240  const auto wellrate_index = well_index * pu.num_phases;
241  const auto& wv = this->wellRates();
242  if( pu.phase_used[BlackoilPhases::Aqua] ) {
243  well.rates.set( rt::wat, wv[ wellrate_index + pu.phase_pos[BlackoilPhases::Aqua] ] );
244  }
245 
246  if( pu.phase_used[BlackoilPhases::Liquid] ) {
247  well.rates.set( rt::oil, wv[ wellrate_index + pu.phase_pos[BlackoilPhases::Liquid] ] );
248  }
249 
250  if( pu.phase_used[BlackoilPhases::Vapour] ) {
251  well.rates.set( rt::gas, wv[ wellrate_index + pu.phase_pos[BlackoilPhases::Vapour] ] );
252  }
253 
254  const int num_perf_well = this->wells_->well_connpos[ well_index + 1 ]
255  - this->wells_->well_connpos[ well_index ];
256  well.completions.resize(num_perf_well);
257 
258  for( int i = 0; i < num_perf_well; ++i ) {
259  const auto wi = this->wells_->well_connpos[ well_index ] + i;
260  const auto active_index = this->wells_->well_cells[ wi ];
261 
262  auto& completion = well.completions[ i ];
263  completion.index = active_index;
264  completion.pressure = this->perfPress()[ itr.second[1] + i ];
265  completion.reservoir_rate = this->perfRates()[ itr.second[1] + i ];
266  }
267  assert(num_perf_well == int(well.completions.size()));
268  }
269 
270  return dw;
271 
272  }
273 
274  virtual ~WellState() {}
275 
276  WellState() = default;
277  WellState( const WellState& rhs ) :
278  bhp_( rhs.bhp_ ),
279  thp_( rhs.thp_ ),
280  temperature_( rhs.temperature_ ),
281  wellrates_( rhs.wellrates_ ),
282  perfrates_( rhs.perfrates_ ),
283  perfpress_( rhs.perfpress_ ),
284  wellMap_( rhs.wellMap_ ),
285  wells_( clone_wells( rhs.wells_.get() ) )
286  {}
287 
288  WellState& operator=( const WellState& rhs ) {
289  this->bhp_ = rhs.bhp_;
290  this->thp_ = rhs.thp_;
291  this->temperature_ = rhs.temperature_;
292  this->wellrates_ = rhs.wellrates_;
293  this->perfrates_ = rhs.perfrates_;
294  this->perfpress_ = rhs.perfpress_;
295  this->wellMap_ = rhs.wellMap_;
296  this->wells_.reset( clone_wells( rhs.wells_.get() ) );
297 
298  return *this;
299  }
300 
301  private:
302  std::vector<double> bhp_;
303  std::vector<double> thp_;
304  std::vector<double> temperature_;
305  std::vector<double> wellrates_;
306  std::vector<double> perfrates_;
307  std::vector<double> perfpress_;
308 
309  WellMapType wellMap_;
310 
311  protected:
312  struct wdel {
313  void operator()( Wells* w ) { destroy_wells( w ); }
314  };
315  std::unique_ptr< Wells, wdel > wells_;
316  };
317 
318 } // namespace Opm
319 
320 #endif // OPM_WELLSTATE_HEADER_INCLUDED
Well constrained by BHP target.
Definition: well_controls.h:35
API for managing sets of well controls.
char ** name
Well names.
Definition: wells.h:107
Controls for a single well.
Definition: well_controls.c:72
Well constrained by THP target.
Definition: well_controls.h:36
int * well_cells
Array of perforation cell indices.
Definition: wells.h:87
int number_of_wells
Number of wells.
Definition: wells.h:52
Definition: AnisotropicEikonal.cpp:446
std::vector< double > & wellRates()
One rate per well and phase.
Definition: WellState.hpp:181
Data structure aggregating static information about all wells in a scenario.
Definition: wells.h:50
Well is an injector.
Definition: wells.h:42
std::vector< double > & temperature()
One temperature per well.
Definition: WellState.hpp:177
int * well_connpos
Array of indices into well_cells (and WI).
Definition: wells.h:81
Main OPM-Core well data structure along with functions to create, populate and destroy it...
Definition: WellState.hpp:312
enum WellType * type
Array of well types.
Definition: wells.h:58
int number_of_phases
Number of phases.
Definition: wells.h:53
std::vector< double > & bhp()
One bhp pressure per well.
Definition: WellState.hpp:169
Definition: BlackoilPhases.hpp:36
std::vector< double > & thp()
One thp pressure per well.
Definition: WellState.hpp:173
void init(const Wells *wells, const std::vector< double > &cellPressures)
Allocate and initialize if wells is non-null.
Definition: WellState.hpp:56
Well is a producer.
Definition: wells.h:43
int numWells() const
The number of wells present.
Definition: WellState.hpp:216
int numPhases() const
The number of phases present.
Definition: WellState.hpp:222
struct Wells * clone_wells(const struct Wells *W)
Create a deep-copy (i.e., clone) of an existing Wells object, including its controls.
Definition: wells.c:493
void destroy_wells(struct Wells *W)
Wells object destructor.
Definition: wells.c:310
struct WellControls ** ctrls
Well controls, one set of controls for each well.
Definition: wells.h:102
std::vector< double > & perfRates()
One rate per well connection.
Definition: WellState.hpp:185
std::vector< double > & perfPress()
One pressure per well connection.
Definition: WellState.hpp:189
The state of a set of wells.
Definition: WellState.hpp:39
Well constrained by surface volume flow rate.
Definition: well_controls.h:38