All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Pages
WellsManager.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_WELLSMANAGER_HEADER_INCLUDED
21 #define OPM_WELLSMANAGER_HEADER_INCLUDED
22 
23 #include <unordered_set>
24 
25 #include <opm/parser/eclipse/EclipseState/EclipseState.hpp>
26 
27 #include <opm/core/wells/WellCollection.hpp>
28 #include <opm/core/wells/WellsGroup.hpp>
29 #include <opm/core/wells/DynamicListEconLimited.hpp>
30 #include <opm/parser/eclipse/EclipseState/Schedule/GroupTree.hpp>
31 
33 
34 struct Wells;
35 struct UnstructuredGrid;
36 
37 
38 namespace Opm
39 {
40 
41  class Schedule;
42 
43  struct WellData
44  {
45  WellType type;
46  bool allowCrossFlow;
47  // WellControlType control;
48  // double target;
49  double reference_bhp_depth;
50  // Opm::InjectionSpecification::InjectorType injected_phase;
51  int welspecsline;
52  };
53 
54 
55  struct PerfData
56  {
57  int cell;
58  double well_index;
59  int satnumid;
60  };
66  {
67  public:
69  WellsManager();
70 
76  explicit WellsManager(struct Wells* W);
77 
86  template<class F2C, class FC>
87  WellsManager(const Opm::EclipseState& eclipseState,
88  const size_t timeStep,
89  int num_cells,
90  const int* global_cell,
91  const int* cart_dims,
92  int dimensions,
93  const F2C& f2c,
94  FC begin_face_centroids,
95  const DynamicListEconLimited& list_econ_limited,
96  bool is_parallel_run=false,
97  const std::unordered_set<std::string>& deactivated_wells = std::unordered_set<std::string> ());
98 
99  WellsManager(const Opm::EclipseState& eclipseState,
100  const size_t timeStep,
101  const UnstructuredGrid& grid);
103  ~WellsManager();
104 
106  bool empty() const;
107 
111  const Wells* c_wells() const;
112 
114  const WellCollection& wellCollection() const;
116 
137  bool conditionsMet(const std::vector<double>& well_bhp,
138  const std::vector<double>& well_reservoirrates_phase,
139  const std::vector<double>& well_surfacerates_phase);
140 
150  void applyExplicitReinjectionControls(const std::vector<double>& well_reservoirrates_phase,
151  const std::vector<double>& well_surfacerates_phase);
152 
153 
154  private:
155  template<class C2F, class FC>
156  void init(const Opm::EclipseState& eclipseState,
157  const size_t timeStep,
158  int num_cells,
159  const int* global_cell,
160  const int* cart_dims,
161  int dimensions,
162  const C2F& cell_to_faces,
163  FC begin_face_centroids,
164  const DynamicListEconLimited& list_econ_limited,
165  const std::unordered_set<std::string>& deactivated_wells);
166  // Disable copying and assignment.
167  WellsManager(const WellsManager& other);
168  WellsManager& operator=(const WellsManager& other);
169  static void setupCompressedToCartesian(const int* global_cell, int number_of_cells, std::map<int,int>& cartesian_to_compressed );
170  void setupWellControls(std::vector<const Well*>& wells, size_t timeStep,
171  std::vector<std::string>& well_names, const PhaseUsage& phaseUsage,
172  const std::vector<int>& wells_on_proc,
173  const DynamicListEconLimited& list_econ_limited);
174 
175  template<class C2F, class FC, class NTG>
176  void createWellsFromSpecs( std::vector<const Well*>& wells, size_t timeStep,
177  const C2F& cell_to_faces,
178  const int* cart_dims,
179  FC begin_face_centroids,
180  int dimensions,
181  std::vector<double>& dz,
182  std::vector<std::string>& well_names,
183  std::vector<WellData>& well_data,
184  std::map<std::string, int> & well_names_to_index,
185  const PhaseUsage& phaseUsage,
186  const std::map<int,int>& cartesian_to_compressed,
187  const double* permeability,
188  const NTG& ntg,
189  std::vector<int>& wells_on_proc,
190  const std::unordered_set<std::string>& deactivated_wells,
191  const DynamicListEconLimited& list_econ_limited);
192 
193  void setupGuideRates(std::vector<const Well*>& wells, const size_t timeStep, std::vector<WellData>& well_data, std::map<std::string, int>& well_names_to_index);
194 
195  // Data
196  Wells* w_;
197  WellCollection well_collection_;
198  // Whether this is a parallel simulation
199  bool is_parallel_run_;
200  };
201 
202 } // namespace Opm
203 
204 #include "WellsManager_impl.hpp"
205 #endif // OPM_WELLSMANAGER_HEADER_INCLUDED
Facility for accessing active subset of data arrays defined for all global cells. ...
void applyExplicitReinjectionControls(const std::vector< double > &well_reservoirrates_phase, const std::vector< double > &well_surfacerates_phase)
Applies explicit reinjection controls.
Definition: WellsManager.cpp:399
Definition: WellsManager.hpp:43
const WellCollection & wellCollection() const
Access the well group hierarchy.
Definition: WellsManager.cpp:369
~WellsManager()
Destructor.
Definition: WellsManager.cpp:347
Definition: WellsManager.hpp:55
Data structure aggregating static information about all wells in a scenario.
Definition: wells.h:50
bool conditionsMet(const std::vector< double > &well_bhp, const std::vector< double > &well_reservoirrates_phase, const std::vector< double > &well_surfacerates_phase)
Checks if each condition is met, applies well controls where needed (that is, it either changes the a...
Definition: WellsManager.cpp:380
to handle the wells and connections violating economic limits.
Definition: DynamicListEconLimited.hpp:32
WellsManager()
Default constructor – no wells.
Definition: WellsManager.cpp:317
WellType
Well type indicates desired/expected well behaviour.
Definition: wells.h:41
bool empty() const
Does the &quot;deck&quot; define any wells?
Definition: WellsManager.cpp:354
const Wells * c_wells() const
Access the managed Wells.
Definition: WellsManager.cpp:364
Definition: BlackoilPhases.hpp:36
This class manages a Wells struct in the sense that it encapsulates creation and destruction of the w...
Definition: WellsManager.hpp:65
Definition: WellCollection.hpp:35