PvtDeadSpline.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_PVTDEADSPLINE_HEADER_INCLUDED
21 #define OPM_PVTDEADSPLINE_HEADER_INCLUDED
22 
25 
26 #include <opm/parser/eclipse/EclipseState/EclipseState.hpp>
27 #include <opm/parser/eclipse/EclipseState/Tables/TableContainer.hpp>
28 
29 #include <vector>
30 
31 namespace Opm
32 {
33 
41  class PvtDeadSpline : public PvtInterface
42  {
43  public:
44  PvtDeadSpline();
45 
46  void initFromOil(const TableContainer& pvdoTables,
47  int numSamples);
48  void initFromGas(const TableContainer& pvdgTables,
49  int numSamples);
50 
51  virtual ~PvtDeadSpline();
52 
54  virtual void mu(const int n,
55  const int* pvtTableIdx,
56  const double* p,
57  const double* T,
58  const double* z,
59  double* output_mu) const;
60 
63  virtual void mu(const int n,
64  const int* pvtTableIdx,
65  const double* p,
66  const double* T,
67  const double* r,
68  double* output_mu,
69  double* output_dmudp,
70  double* output_dmudr) const;
71 
74  virtual void mu(const int n,
75  const int* pvtTableIdx,
76  const double* p,
77  const double* T,
78  const double* r,
79  const PhasePresence* cond,
80  double* output_mu,
81  double* output_dmudp,
82  double* output_dmudr) const;
83 
85  virtual void B(const int n,
86  const int* pvtTableIdx,
87  const double* p,
88  const double* T,
89  const double* z,
90  double* output_B) const;
91 
93  virtual void dBdp(const int n,
94  const int* pvtTableIdx,
95  const double* p,
96  const double* T,
97  const double* z,
98  double* output_B,
99  double* output_dBdp) const;
100 
103  virtual void b(const int n,
104  const int* pvtTableIdx,
105  const double* p,
106  const double* T,
107  const double* r,
108  double* output_b,
109  double* output_dbdp,
110  double* output_dbdr) const;
111 
114  virtual void b(const int n,
115  const int* pvtTableIdx,
116  const double* p,
117  const double* T,
118  const double* r,
119  const PhasePresence* cond,
120  double* output_b,
121  double* output_dbdp,
122  double* output_dbdr) const;
123 
125  virtual void rsSat(const int n,
126  const int* pvtTableIdx,
127  const double* p,
128  double* output_rsSat,
129  double* output_drsSatdp) const;
130 
132  virtual void rvSat(const int n,
133  const int* pvtTableIdx,
134  const double* p,
135  double* output_rvSat,
136  double* output_drvSatdp) const;
137 
139  virtual void R(const int n,
140  const int* pvtTableIdx,
141  const double* p,
142  const double* z,
143  double* output_R) const;
144 
146  virtual void dRdp(const int n,
147  const int* pvtTableIdx,
148  const double* p,
149  const double* z,
150  double* output_R,
151  double* output_dRdp) const;
152  private:
153  int getTableIndex_(const int* pvtTableIdx, int cellIdx) const
154  {
155  if (!pvtTableIdx)
156  return 0;
157  return pvtTableIdx[cellIdx];
158  }
159 
160  // PVT properties of dry gas or dead oil. We need to store one
161  // table per PVT region.
162  std::vector<UniformTableLinear<double> > b_;
163  std::vector<UniformTableLinear<double> > viscosity_;
164  };
165 
166 }
167 
168 #endif // OPM_PVTDEADSPLINE_HEADER_INCLUDED
169 
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...
Definition: BlackoilPhases.hpp:48
void initFromGas(const TableContainer &pvdgTables, int numSamples)
Definition: AnisotropicEikonal.hpp:43
virtual ~PvtDeadSpline()
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 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
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.
void initFromOil(const TableContainer &pvdoTables, int numSamples)
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 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.
Definition: PvtDeadSpline.hpp:41
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 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.