WellsManager.hpp
Go to the documentation of this file.
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 
24 #include <opm/parser/eclipse/EclipseState/EclipseState.hpp>
25 
28 #include <opm/parser/eclipse/EclipseState/Schedule/GroupTree.hpp>
29 
31 
32 struct Wells;
33 struct UnstructuredGrid;
34 
35 
36 namespace Opm
37 {
38  struct WellData
39  {
42  // WellControlType control;
43  // double target;
45  // Opm::InjectionSpecification::InjectorType injected_phase;
47  };
48 
49 
50  struct PerfData
51  {
52  int cell;
53  double well_index;
54  };
60  {
61  public:
63  WellsManager();
64 
70  explicit WellsManager(struct Wells* W);
71 
76  template<class F2C, class FC>
77  WellsManager(const Opm::EclipseStateConstPtr eclipseState,
78  const size_t timeStep,
79  int num_cells,
80  const int* global_cell,
81  const int* cart_dims,
82  int dimensions,
83  const F2C& f2c,
84  FC begin_face_centroids,
85  const double* permeability,
86  bool is_parallel_run=false);
87 
88  WellsManager(const Opm::EclipseStateConstPtr eclipseState,
89  const size_t timeStep,
90  const UnstructuredGrid& grid,
91  const double* permeability);
93  ~WellsManager();
94 
96  bool empty() const;
97 
101  const Wells* c_wells() const;
102 
104  const WellCollection& wellCollection() const;
105 
126  bool conditionsMet(const std::vector<double>& well_bhp,
127  const std::vector<double>& well_reservoirrates_phase,
128  const std::vector<double>& well_surfacerates_phase);
129 
139  void applyExplicitReinjectionControls(const std::vector<double>& well_reservoirrates_phase,
140  const std::vector<double>& well_surfacerates_phase);
141 
142 
143  private:
144  template<class C2F, class FC>
145  void init(const Opm::EclipseStateConstPtr eclipseState,
146  const size_t timeStep,
147  int num_cells,
148  const int* global_cell,
149  const int* cart_dims,
150  int dimensions,
151  const C2F& cell_to_faces,
152  FC begin_face_centroids,
153  const double* permeability);
154  // Disable copying and assignment.
155  WellsManager(const WellsManager& other);
156  WellsManager& operator=(const WellsManager& other);
157  static void setupCompressedToCartesian(const int* global_cell, int number_of_cells, std::map<int,int>& cartesian_to_compressed );
158  void setupWellControls(std::vector<WellConstPtr>& wells, size_t timeStep,
159  std::vector<std::string>& well_names, const PhaseUsage& phaseUsage,
160  const std::vector<int>& wells_on_proc);
161 
162  template<class C2F, class FC, class NTG>
163  void createWellsFromSpecs( std::vector<WellConstPtr>& wells, size_t timeStep,
164  const C2F& cell_to_faces,
165  const int* cart_dims,
166  FC begin_face_centroids,
167  int dimensions,
168  std::vector<double>& dz,
169  std::vector<std::string>& well_names,
170  std::vector<WellData>& well_data,
171  std::map<std::string, int> & well_names_to_index,
172  const PhaseUsage& phaseUsage,
173  const std::map<int,int>& cartesian_to_compressed,
174  const double* permeability,
175  const NTG& ntg,
176  std::vector<int>& wells_on_proc);
177 
178  void addChildGroups(GroupTreeNodeConstPtr parentNode, ScheduleConstPtr schedule, size_t timeStep, const PhaseUsage& phaseUsage);
179  void setupGuideRates(std::vector<WellConstPtr>& wells, const size_t timeStep, std::vector<WellData>& well_data, std::map<std::string, int>& well_names_to_index);
180 
181 
182  // Data
183  Wells* w_;
184  WellCollection well_collection_;
185  // Whether this is a parallel simulation
186  bool is_parallel_run_;
187  };
188 
189 } // namespace Opm
190 
191 #include "WellsManager_impl.hpp"
192 #endif // OPM_WELLSMANAGER_HEADER_INCLUDED
Definition: wells.h:50
Definition: grid.h:98
Definition: AnisotropicEikonal.hpp:43
const WellCollection & wellCollection() const
Access the well group hierarchy.
double reference_bhp_depth
Definition: WellsManager.hpp:44
bool empty() const
Does the "deck" define any wells?
~WellsManager()
Destructor.
WellsManager()
Default constructor – no wells.
int welspecsline
Definition: WellsManager.hpp:46
Definition: WellsManager.hpp:59
bool allowCrossFlow
Definition: WellsManager.hpp:41
int dimensions(const UnstructuredGrid &grid)
Get the dimensions of a grid.
bool conditionsMet(const std::vector< double > &well_bhp, const std::vector< double > &well_reservoirrates_phase, const std::vector< double > &well_surfacerates_phase)
void applyExplicitReinjectionControls(const std::vector< double > &well_reservoirrates_phase, const std::vector< double > &well_surfacerates_phase)
Definition: WellCollection.hpp:35
WellType type
Definition: WellsManager.hpp:40
Definition: WellsManager.hpp:50
int cell
Definition: WellsManager.hpp:52
WellType
Definition: wells.h:41
const UnstructuredGrid & grid
Definition: ColumnExtract.hpp:31
const Wells * c_wells() const
Definition: WellsManager.hpp:38
Definition: BlackoilPhases.hpp:36
double well_index
Definition: WellsManager.hpp:53