PvtConstCompr.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_PVTCONSTCOMPR_HEADER_INCLUDED
21 #define OPM_PVTCONSTCOMPR_HEADER_INCLUDED
22 
24 #include <opm/common/ErrorMacros.hpp>
25 
26 #include <vector>
27 #include <algorithm>
28 
29 
30 namespace Opm
31 {
32 
44  class PvtConstCompr : public PvtInterface
45  {
46  public:
48  {}
49 
50  void initFromWater(Opm::DeckKeywordConstPtr pvtwKeyword)
51  {
52  int numRegions = pvtwKeyword->size();
53 
54  ref_press_.resize(numRegions);
55  ref_B_.resize(numRegions);
56  comp_.resize(numRegions);
57  viscosity_.resize(numRegions);
58  visc_comp_.resize(numRegions);
59 
60  for (int regionIdx = 0; regionIdx < numRegions; ++ regionIdx) {
61  Opm::DeckRecordConstPtr pvtwRecord = pvtwKeyword->getRecord(regionIdx);
62 
63  ref_press_[regionIdx] = pvtwRecord->getItem("P_REF")->getSIDouble(0);
64  ref_B_[regionIdx] = pvtwRecord->getItem("WATER_VOL_FACTOR")->getSIDouble(0);
65  comp_[regionIdx] = pvtwRecord->getItem("WATER_COMPRESSIBILITY")->getSIDouble(0);
66  viscosity_[regionIdx] = pvtwRecord->getItem("WATER_VISCOSITY")->getSIDouble(0);
67  visc_comp_[regionIdx] = pvtwRecord->getItem("WATER_VISCOSIBILITY")->getSIDouble(0);
68  }
69  }
70 
71  void initFromOil(Opm::DeckKeywordConstPtr pvcdoKeyword)
72  {
73  int numRegions = pvcdoKeyword->size();
74 
75  ref_press_.resize(numRegions);
76  ref_B_.resize(numRegions);
77  comp_.resize(numRegions);
78  viscosity_.resize(numRegions);
79  visc_comp_.resize(numRegions);
80 
81  for (int regionIdx = 0; regionIdx < numRegions; ++ regionIdx) {
82  Opm::DeckRecordConstPtr pvcdoRecord = pvcdoKeyword->getRecord(regionIdx);
83 
84  ref_press_[regionIdx] = pvcdoRecord->getItem("P_REF")->getSIDouble(0);
85  ref_B_[regionIdx] = pvcdoRecord->getItem("OIL_VOL_FACTOR")->getSIDouble(0);
86  comp_[regionIdx] = pvcdoRecord->getItem("OIL_COMPRESSIBILITY")->getSIDouble(0);
87  viscosity_[regionIdx] = pvcdoRecord->getItem("OIL_VISCOSITY")->getSIDouble(0);
88  visc_comp_[regionIdx] = pvcdoRecord->getItem("OIL_VISCOSIBILITY")->getSIDouble(0);
89  }
90  }
91 
96  explicit PvtConstCompr(double visc)
97  : ref_press_(1, 0.0),
98  ref_B_(1, 1.0),
99  comp_(1, 0.0),
100  viscosity_(1, visc),
101  visc_comp_(1, 0.0)
102  {
103  }
104 
105  virtual ~PvtConstCompr()
106  {
107  }
108 
109  virtual void mu(const int n,
110  const int* pvtRegionIdx,
111  const double* p,
112  const double* /*T*/,
113  const double* /*z*/,
114  double* output_mu) const
115  {
116 // #pragma omp parallel for
117  for (int i = 0; i < n; ++i) {
118  // Computing a polynomial approximation to the exponential.
119  int tableIdx = getTableIndex_(pvtRegionIdx, i);
120  double x = -visc_comp_[tableIdx]*(p[i] - ref_press_[tableIdx]);
121  output_mu[i] = viscosity_[tableIdx]/(1.0 + x + 0.5*x*x);
122  }
123  }
124 
125  virtual void mu(const int n,
126  const int* pvtRegionIdx,
127  const double* p,
128  const double* /*T*/,
129  const double* /*r*/,
130  double* output_mu,
131  double* output_dmudp,
132  double* output_dmudr) const
133  {
134  // #pragma omp parallel for
135  for (int i = 0; i < n; ++i) {
136  // Computing a polynomial approximation to the exponential.
137  int tableIdx = getTableIndex_(pvtRegionIdx, i);
138  double x = -visc_comp_[tableIdx]*(p[i] - ref_press_[tableIdx]);
139  double d = (1.0 + x + 0.5*x*x);
140  output_mu[i] = viscosity_[tableIdx]/d;
141  output_dmudp[i] = (viscosity_[tableIdx]/(d*d))*(1+x) * visc_comp_[tableIdx];
142  }
143  std::fill(output_dmudr, output_dmudr + n, 0.0);
144  }
145 
146  virtual void mu(const int n,
147  const int* pvtRegionIdx,
148  const double* p,
149  const double* /*T*/,
150  const double* /*r*/,
151  const PhasePresence* /*cond*/,
152  double* output_mu,
153  double* output_dmudp,
154  double* output_dmudr) const
155  {
156  // #pragma omp parallel for
157  for (int i = 0; i < n; ++i) {
158  // Computing a polynomial approximation to the exponential.
159  int tableIdx = getTableIndex_(pvtRegionIdx, i);
160  double x = -visc_comp_[tableIdx]*(p[i] - ref_press_[tableIdx]);
161  double d = (1.0 + x + 0.5*x*x);
162  output_mu[i] = viscosity_[tableIdx]/d;
163  output_dmudp[i] = (viscosity_[tableIdx]/(d*d))*(1+x) * visc_comp_[tableIdx];
164  }
165  std::fill(output_dmudr, output_dmudr + n, 0.0);
166  }
167 
168  virtual void B(const int n,
169  const int* pvtRegionIdx,
170  const double* p,
171  const double* /*T*/,
172  const double* /*z*/,
173  double* output_B) const
174  {
175 // #pragma omp parallel for
176  for (int i = 0; i < n; ++i) {
177  // Computing a polynomial approximation to the exponential.
178  int tableIdx = getTableIndex_(pvtRegionIdx, i);
179  double x = comp_[tableIdx]*(p[i] - ref_press_[tableIdx]);
180  output_B[i] = ref_B_[tableIdx]/(1.0 + x + 0.5*x*x);
181  }
182  }
183 
184  virtual void dBdp(const int n,
185  const int* pvtRegionIdx,
186  const double* p,
187  const double* /*T*/,
188  const double* /*z*/,
189  double* output_B,
190  double* output_dBdp) const
191  {
192 // #pragma omp parallel for
193  for (int i = 0; i < n; ++i) {
194  int tableIdx = getTableIndex_(pvtRegionIdx, i);
195  double x = comp_[tableIdx]*(p[i] - ref_press_[tableIdx]);
196  double d = (1.0 + x + 0.5*x*x);
197  output_B[i] = ref_B_[tableIdx]/d;
198  output_dBdp[i] = (-ref_B_[tableIdx]/(d*d))*(1 + x) * comp_[tableIdx];
199  }
200  }
201 
202  virtual void b(const int n,
203  const int* pvtRegionIdx,
204  const double* p,
205  const double* /*T*/,
206  const double* /*r*/,
207  double* output_b,
208  double* output_dbdp,
209  double* output_dbdr) const
210  {
211  // #pragma omp parallel for
212  for (int i = 0; i < n; ++i) {
213  // Computing a polynomial approximation to the exponential.
214  int tableIdx = getTableIndex_(pvtRegionIdx, i);
215  double x = comp_[tableIdx]*(p[i] - ref_press_[tableIdx]);
216  double d = (1.0 + x + 0.5*x*x);
217 
218  // b = 1/B = d/ref_B_B;
219  output_b[i] = d/ref_B_[tableIdx];
220  output_dbdp[i] = (1 + x) * comp_[tableIdx]/ref_B_[tableIdx];
221  }
222 
223  std::fill(output_dbdr, output_dbdr + n, 0.0);
224  }
225 
226  virtual void b(const int n,
227  const int* pvtRegionIdx,
228  const double* p,
229  const double* /*T*/,
230  const double* /*r*/,
231  const PhasePresence* /*cond*/,
232  double* output_b,
233  double* output_dbdp,
234  double* output_dbdr) const
235  {
236  // #pragma omp parallel for
237  for (int i = 0; i < n; ++i) {
238  // Computing a polynomial approximation to the exponential.
239  int tableIdx = getTableIndex_(pvtRegionIdx, i);
240  double x = comp_[tableIdx]*(p[i] - ref_press_[tableIdx]);
241  double d = (1.0 + x + 0.5*x*x);
242 
243  // b = 1/B = d/ref_B_[tableIdx]B;
244  output_b[i] = d/ref_B_[tableIdx];
245  output_dbdp[i] = (1 + x) * comp_[tableIdx]/ref_B_[tableIdx];
246  }
247  std::fill(output_dbdr, output_dbdr + n, 0.0);
248  }
249 
250  virtual void rsSat(const int n,
251  const int* /*pvtRegionIdx*/,
252  const double* /*p*/,
253  double* output_rsSat,
254  double* output_drsSatdp) const
255  {
256  std::fill(output_rsSat, output_rsSat + n, 0.0);
257  std::fill(output_drsSatdp, output_drsSatdp + n, 0.0);
258  }
259 
260  virtual void rvSat(const int n,
261  const int* /*pvtRegionIdx*/,
262  const double* /*p*/,
263  double* output_rvSat,
264  double* output_drvSatdp) const
265  {
266  std::fill(output_rvSat, output_rvSat + n, 0.0);
267  std::fill(output_drvSatdp, output_drvSatdp + n, 0.0);
268  }
269 
270  virtual void R(const int n,
271  const int* /*pvtRegionIdx*/,
272  const double* /*p*/,
273  const double* /*z*/,
274  double* output_R) const
275  {
276  std::fill(output_R, output_R + n, 0.0);
277  }
278 
279  virtual void dRdp(const int n,
280  const int* /*pvtRegionIdx*/,
281  const double* /*p*/,
282  const double* /*z*/,
283  double* output_R,
284  double* output_dRdp) const
285  {
286  std::fill(output_R, output_R + n, 0.0);
287  std::fill(output_dRdp, output_dRdp + n, 0.0);
288  }
289 
290  private:
291  int getTableIndex_(const int* pvtTableIdx, int cellIdx) const
292  {
293  if (!pvtTableIdx)
294  return 0;
295  return pvtTableIdx[cellIdx];
296  }
297 
298  // The PVT properties. We need to store one value per PVT
299  // region.
300  std::vector<double> ref_press_;
301  std::vector<double> ref_B_;
302  std::vector<double> comp_;
303  std::vector<double> viscosity_;
304  std::vector<double> visc_comp_;
305  };
306 
307 }
308 
309 
310 #endif // OPM_PVTCONSTCOMPR_HEADER_INCLUDED
311 
Definition: BlackoilPhases.hpp:48
Definition: AnisotropicEikonal.hpp:43
virtual void mu(const int n, const int *pvtRegionIdx, const double *p, const double *, const double *, const PhasePresence *, double *output_mu, double *output_dmudp, double *output_dmudr) const
Definition: PvtConstCompr.hpp:146
virtual void dRdp(const int n, const int *, const double *, const double *, double *output_R, double *output_dRdp) const
Solution factor and p-derivative as functions of p and z.
Definition: PvtConstCompr.hpp:279
virtual void R(const int n, const int *, const double *, const double *, double *output_R) const
Solution factor as a function of p and z.
Definition: PvtConstCompr.hpp:270
PvtConstCompr()
Definition: PvtConstCompr.hpp:47
virtual void b(const int n, const int *pvtRegionIdx, const double *p, const double *, const double *, double *output_b, double *output_dbdp, double *output_dbdr) const
Definition: PvtConstCompr.hpp:202
virtual void dBdp(const int n, const int *pvtRegionIdx, const double *p, const double *, const double *, double *output_B, double *output_dBdp) const
Formation volume factor and p-derivative as functions of p and z.
Definition: PvtConstCompr.hpp:184
void initFromOil(Opm::DeckKeywordConstPtr pvcdoKeyword)
Definition: PvtConstCompr.hpp:71
virtual void B(const int n, const int *pvtRegionIdx, const double *p, const double *, const double *, double *output_B) const
Formation volume factor as a function of p, T and z.
Definition: PvtConstCompr.hpp:168
virtual void b(const int n, const int *pvtRegionIdx, const double *p, const double *, const double *, const PhasePresence *, double *output_b, double *output_dbdp, double *output_dbdr) const
Definition: PvtConstCompr.hpp:226
virtual void mu(const int n, const int *pvtRegionIdx, const double *p, const double *, const double *, double *output_mu) const
Viscosity as a function of p, T and z.
Definition: PvtConstCompr.hpp:109
Definition: PvtConstCompr.hpp:44
void initFromWater(Opm::DeckKeywordConstPtr pvtwKeyword)
Definition: PvtConstCompr.hpp:50
Definition: PvtInterface.hpp:30
virtual void rsSat(const int n, const int *, const double *, double *output_rsSat, double *output_drsSatdp) const
Solution gas/oil ratio and its derivatives at saturated conditions as a function of p...
Definition: PvtConstCompr.hpp:250
virtual void rvSat(const int n, const int *, const double *, double *output_rvSat, double *output_drvSatdp) const
Vapor oil/gas ratio and its derivatives at saturated conditions as a function of p.
Definition: PvtConstCompr.hpp:260
PvtConstCompr(double visc)
Create a PVT object with a given viscosity that assumes all fluid phases to be incompressible.
Definition: PvtConstCompr.hpp:96
const double visc
Definition: Units.hpp:152
virtual void mu(const int n, const int *pvtRegionIdx, const double *p, const double *, const double *, double *output_mu, double *output_dmudp, double *output_dmudr) const
Definition: PvtConstCompr.hpp:125
virtual ~PvtConstCompr()
Definition: PvtConstCompr.hpp:105