PvtLiveOil.hpp
Go to the documentation of this file.
1 /*
2  Copyright 2010, 2011, 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_PVTLIVEOIL_HEADER_INCLUDED
21 #define OPM_PVTLIVEOIL_HEADER_INCLUDED
22 
24 
25 #include <opm/parser/eclipse/Deck/Deck.hpp>
26 #include <opm/parser/eclipse/EclipseState/EclipseState.hpp>
27 
28 #include <vector>
29 
30 namespace Opm
31 {
39  class PvtLiveOil : public PvtInterface
40  {
41  public:
42  PvtLiveOil(const std::vector<Opm::PvtoTable>& pvtoTables);
43  virtual ~PvtLiveOil();
44 
46  virtual void mu(const int n,
47  const int* pvtTableIdx,
48  const double* p,
49  const double* T,
50  const double* z,
51  double* output_mu) const;
52 
55  virtual void mu(const int n,
56  const int* pvtTableIdx,
57  const double* p,
58  const double* T,
59  const double* r,
60  double* output_mu,
61  double* output_dmudp,
62  double* output_dmudr) const;
63 
66  virtual void mu(const int n,
67  const int* pvtTableIdx,
68  const double* p,
69  const double* T,
70  const double* r,
71  const PhasePresence* cond,
72  double* output_mu,
73  double* output_dmudp,
74  double* output_dmudr) const;
75 
77  virtual void B(const int n,
78  const int* pvtTableIdx,
79  const double* p,
80  const double* T,
81  const double* z,
82  double* output_B) const;
83 
85  virtual void dBdp(const int n,
86  const int* pvtTableIdx,
87  const double* p,
88  const double* T,
89  const double* z,
90  double* output_B,
91  double* output_dBdp) const;
92 
95  virtual void b(const int n,
96  const int* pvtTableIdx,
97  const double* p,
98  const double* T,
99  const double* r,
100  double* output_b,
101  double* output_dbdp,
102  double* output_dbdr) const;
103 
106  virtual void b(const int n,
107  const int* pvtTableIdx,
108  const double* p,
109  const double* T,
110  const double* r,
111  const PhasePresence* cond,
112  double* output_b,
113  double* output_dbdp,
114  double* output_dbdr) const;
115 
117  virtual void rsSat(const int n,
118  const int* pvtTableIdx,
119  const double* p,
120  double* output_rsSat,
121  double* output_drsSatdp) const;
122 
124  virtual void rvSat(const int n,
125  const int* pvtTableIdx,
126  const double* p,
127  double* output_rvSat,
128  double* output_drvSatdp) const;
129 
131  virtual void R(const int n,
132  const int* pvtTableIdx,
133  const double* p,
134  const double* z,
135  double* output_R) const;
136 
138  virtual void dRdp(const int n,
139  const int* pvtTableIdx,
140  const double* p,
141  const double* z,
142  double* output_R,
143  double* output_dRdp) const;
144 
145  private:
146  int getTableIndex_(const int* pvtTableIdx, int cellIdx) const
147  {
148  if (!pvtTableIdx)
149  return 0;
150  return pvtTableIdx[cellIdx];
151  }
152 
153  double evalB(size_t pvtTableIdx, double press, const double* surfvol) const;
154  void evalBDeriv(size_t pvtTableIdx, double press, const double* surfvol, double& B, double& dBdp) const;
155  double evalR(size_t pvtTableIdx, double press, const double* surfvol) const;
156  void evalRDeriv(size_t pvtTableIdx, double press, const double* surfvol, double& R, double& dRdp) const;
157 
158  // item: 1=>1/B 2=>mu;
159  double miscible_oil(const double press,
160  const double* surfvol,
161  const int pvtTableIdx,
162  const int item,
163  const bool deriv = false) const;
164 
165  double miscible_oil(const double press,
166  const double r,
167  const int pvtTableIdx,
168  const int item,
169  const int deriv = 0) const;
170 
171  double miscible_oil(const double press,
172  const double r,
173  const PhasePresence& cond,
174  const int pvtTableIdx,
175  const int item,
176  const int deriv = 0) const;
177 
178  // PVT properties of live oil (with dissolved gas). We need to
179  // store one table per PVT region.
180  std::vector<std::vector<std::vector<double> > > saturated_oil_table_;
181  std::vector<std::vector<std::vector<std::vector<double> > > > undersat_oil_tables_;
182  };
183 
184 }
185 
186 #endif // OPM_PVTLIVEOIL_HEADER_INCLUDED
187 
Definition: BlackoilPhases.hpp:48
virtual void mu(const int n, const int *pvtTableIdx, const double *p, const double *T, const double *z, double *output_mu) const
Viscosity as a function of p, T and z.
Definition: AnisotropicEikonal.hpp:43
virtual ~PvtLiveOil()
virtual void dBdp(const int n, const int *pvtTableIdx, const double *p, const double *T, const double *z, double *output_B, double *output_dBdp) const
Formation volume factor and p-derivative as functions of p, T and z.
PvtLiveOil(const std::vector< Opm::PvtoTable > &pvtoTables)
virtual void rsSat(const int n, const int *pvtTableIdx, const double *p, double *output_rsSat, double *output_drsSatdp) const
Solution gas/oil ratio and its derivatives at saturated conditions as a function of p...
virtual void rvSat(const int n, const int *pvtTableIdx, const double *p, double *output_rvSat, double *output_drvSatdp) const
Vapor oil/gas ratio and its derivatives at saturated conditions as a function of p.
virtual void dRdp(const int n, const int *pvtTableIdx, const double *p, const double *z, double *output_R, double *output_dRdp) const
Solution factor and p-derivative as functions of p and z.
virtual void R(const int n, const int *pvtTableIdx, const double *p, const double *z, double *output_R) const
Solution factor as a function of p and z.
virtual void B(const int n, const int *pvtTableIdx, const double *p, const double *T, const double *z, double *output_B) const
Formation volume factor as a function of p, T and z.
Definition: PvtInterface.hpp:30
virtual void b(const int n, const int *pvtTableIdx, const double *p, const double *T, const double *r, double *output_b, double *output_dbdp, double *output_dbdr) const
Definition: PvtLiveOil.hpp:39