PvtLiveGas.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_PVTLIVEGAS_HEADER_INCLUDED
21 #define OPM_PVTLIVEGAS_HEADER_INCLUDED
22 
24 
25 #include <opm/parser/eclipse/EclipseState/EclipseState.hpp>
26 
27 #include <vector>
28 
29 namespace Opm
30 {
38  class PvtLiveGas : public PvtInterface
39  {
40  public:
41  PvtLiveGas(const std::vector<Opm::PvtgTable>& pvtgTables);
42  virtual ~PvtLiveGas();
43 
45  virtual void mu(const int n,
46  const int* pvtRegionIdx,
47  const double* p,
48  const double* T,
49  const double* z,
50  double* output_mu) const;
51 
54  virtual void mu(const int n,
55  const int* pvtRegionIdx,
56  const double* p,
57  const double* T,
58  const double* r,
59  double* output_mu,
60  double* output_dmudp,
61  double* output_dmudr) const;
62 
65  virtual void mu(const int n,
66  const int* pvtRegionIdx,
67  const double* p,
68  const double* T,
69  const double* r,
70  const PhasePresence* cond,
71  double* output_mu,
72  double* output_dmudp,
73  double* output_dmudr) const;
74 
76  virtual void B(const int n,
77  const int* pvtRegionIdx,
78  const double* p,
79  const double* T,
80  const double* z,
81  double* output_B) const;
82 
84  virtual void dBdp(const int n,
85  const int* pvtRegionIdx,
86  const double* p,
87  const double* T,
88  const double* z,
89  double* output_B,
90  double* output_dBdp) const;
91 
94  virtual void b(const int n,
95  const int* pvtRegionIdx,
96  const double* p,
97  const double* T,
98  const double* r,
99  double* output_b,
100  double* output_dbdp,
101  double* output_dbdr) const;
102 
105  virtual void b(const int n,
106  const int* pvtRegionIdx,
107  const double* p,
108  const double* T,
109  const double* r,
110  const PhasePresence* cond,
111  double* output_b,
112  double* output_dbdp,
113  double* output_dbdr) const;
114 
115 
116 
118  virtual void rsSat(const int n,
119  const int* pvtRegionIdx,
120  const double* p,
121  double* output_rsSat,
122  double* output_drsSatdp) const;
123 
125  virtual void rvSat(const int n,
126  const int* pvtRegionIdx,
127  const double* p,
128  double* output_rvSat,
129  double* output_drvSatdp) const;
130 
132  virtual void R(const int n,
133  const int* pvtRegionIdx,
134  const double* p,
135  const double* z,
136  double* output_R) const;
137 
139  virtual void dRdp(const int n,
140  const int* pvtRegionIdx,
141  const double* p,
142  const double* z,
143  double* output_R,
144  double* output_dRdp) const;
145 
146  protected:
147  int getTableIndex_(const int* pvtTableIdx, int cellIdx) const
148  {
149  if (!pvtTableIdx)
150  return 0;
151  return pvtTableIdx[cellIdx];
152  }
153 
154  double evalB(double press, const double* surfvol, int pvtTableIdx) const;
155  void evalBDeriv(double press, const double* surfvol, int pvtTableIdx, double& B, double& dBdp) const;
156  double evalR(double press, const double* surfvol, int pvtTableIdx) const;
157  void evalRDeriv(double press, const double* surfvol, int pvtTableIdx, double& R, double& dRdp) const;
158 
159  // item: 1=>B 2=>mu;
160  double miscible_gas(const double press,
161  const double* surfvol,
162  const int pvtTableIdx,
163  const int item,
164  const bool deriv = false) const;
165  double miscible_gas(const double press,
166  const double r,
167  const PhasePresence& cond,
168  const int pvtTableIdx,
169  const int item,
170  const int deriv = 0) const;
171  // PVT properties of wet gas (with vaporised oil). We need to
172  // store one table per PVT region.
173  std::vector< std::vector<std::vector<double> > > saturated_gas_table_;
174  std::vector< std::vector<std::vector<std::vector<double> > > > undersat_gas_tables_;
175  };
176 
177 }
178 
179 #endif // OPM_PVTLIVEGAS_HEADER_INCLUDED
180 
Definition: BlackoilPhases.hpp:48
Definition: AnisotropicEikonal.hpp:43
virtual void B(const int n, const int *pvtRegionIdx, const double *p, const double *T, const double *z, double *output_B) const
Formation volume factor as a function of p, T and z.
PvtLiveGas(const std::vector< Opm::PvtgTable > &pvtgTables)
virtual ~PvtLiveGas()
int getTableIndex_(const int *pvtTableIdx, int cellIdx) const
Definition: PvtLiveGas.hpp:147
virtual void b(const int n, const int *pvtRegionIdx, const double *p, const double *T, const double *r, double *output_b, double *output_dbdp, double *output_dbdr) const
double evalB(double press, const double *surfvol, int pvtTableIdx) const
double evalR(double press, const double *surfvol, int pvtTableIdx) const
virtual void rsSat(const int n, const int *pvtRegionIdx, 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...
void evalRDeriv(double press, const double *surfvol, int pvtTableIdx, double &R, double &dRdp) const
virtual void mu(const int n, const int *pvtRegionIdx, const double *p, const double *T, const double *z, double *output_mu) const
Viscosity as a function of p, T and z.
virtual void R(const int n, const int *pvtRegionIdx, const double *p, const double *z, double *output_R) const
Solution factor as a function of p and z.
Definition: PvtLiveGas.hpp:38
double miscible_gas(const double press, const double *surfvol, const int pvtTableIdx, const int item, const bool deriv=false) const
std::vector< std::vector< std::vector< std::vector< double > > > > undersat_gas_tables_
Definition: PvtLiveGas.hpp:174
void evalBDeriv(double press, const double *surfvol, int pvtTableIdx, double &B, double &dBdp) const
virtual void rvSat(const int n, const int *pvtRegionIdx, 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.
std::vector< std::vector< std::vector< double > > > saturated_gas_table_
Definition: PvtLiveGas.hpp:173
Definition: PvtInterface.hpp:30
virtual void dBdp(const int n, const int *pvtRegionIdx, 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.
virtual void dRdp(const int n, const int *pvtRegionIdx, const double *p, const double *z, double *output_R, double *output_dRdp) const
Solution factor and p-derivative as functions of p and z.