ThermalWaterPvtWrapper.hpp
Go to the documentation of this file.
1 /*
2  Copyright 2015 Andreas Lauser
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_THERMAL_WATER_PVT_WRAPPER_HPP
21 #define OPM_THERMAL_WATER_PVT_WRAPPER_HPP
22 
24 #include <opm/common/ErrorMacros.hpp>
25 
26 #include <opm/parser/eclipse/EclipseState/EclipseState.hpp>
27 
28 #include <vector>
29 
30 namespace Opm
31 {
35  {
36  public:
38  {}
39 
40 
42  void initFromDeck(std::shared_ptr<const PvtInterface> isothermalPvt,
43  Opm::DeckConstPtr deck,
44  Opm::EclipseStateConstPtr eclipseState)
45  {
46  isothermalPvt_ = isothermalPvt;
47  watvisctTables_ = 0;
48 
49  // stuff which we need to get from the PVTW keyword
50  Opm::DeckKeywordConstPtr pvtwKeyword = deck->getKeyword("PVTW");
51  int numRegions = pvtwKeyword->size();
52  pvtwRefPress_.resize(numRegions);
53  pvtwRefB_.resize(numRegions);
54  pvtwCompressibility_.resize(numRegions);
55  pvtwViscosity_.resize(numRegions);
56  pvtwViscosibility_.resize(numRegions);
57  for (int regionIdx = 0; regionIdx < numRegions; ++ regionIdx) {
58  Opm::DeckRecordConstPtr pvtwRecord = pvtwKeyword->getRecord(regionIdx);
59  pvtwRefPress_[regionIdx] = pvtwRecord->getItem("P_REF")->getSIDouble(0);
60  pvtwRefB_[regionIdx] = pvtwRecord->getItem("WATER_VOL_FACTOR")->getSIDouble(0);
61  pvtwViscosity_[regionIdx] = pvtwRecord->getItem("WATER_VISCOSITY")->getSIDouble(0);
62  pvtwViscosibility_[regionIdx] = pvtwRecord->getItem("WATER_VISCOSIBILITY")->getSIDouble(0);
63  }
64 
65  // quantities required for the temperature dependence of the viscosity
66  // (basically we expect well-behaved VISCREF and WATVISCT keywords.)
67  if (deck->hasKeyword("VISCREF")) {
68  auto tables = eclipseState->getTableManager();
69  watvisctTables_ = &tables->getWatvisctTables();
70  Opm::DeckKeywordConstPtr viscrefKeyword = deck->getKeyword("VISCREF");
71 
72  assert(int(watvisctTables_->size()) == numRegions);
73  assert(int(viscrefKeyword->size()) == numRegions);
74 
75  viscrefPress_.resize(numRegions);
76  for (int regionIdx = 0; regionIdx < numRegions; ++ regionIdx) {
77  Opm::DeckRecordConstPtr viscrefRecord = viscrefKeyword->getRecord(regionIdx);
78 
79  viscrefPress_[regionIdx] = viscrefRecord->getItem("REFERENCE_PRESSURE")->getSIDouble(0);
80  }
81  }
82 
83  // quantities required for the temperature dependence of the density
84  if (deck->hasKeyword("WATDENT")) {
85  DeckKeywordConstPtr watdentKeyword = deck->getKeyword("WATDENT");
86 
87  assert(int(watdentKeyword->size()) == numRegions);
88 
89  watdentRefTemp_.resize(numRegions);
90  watdentCT1_.resize(numRegions);
91  watdentCT2_.resize(numRegions);
92  for (int regionIdx = 0; regionIdx < numRegions; ++regionIdx) {
93  Opm::DeckRecordConstPtr watdentRecord = watdentKeyword->getRecord(regionIdx);
94 
95  watdentRefTemp_[regionIdx] = watdentRecord->getItem("REFERENCE_TEMPERATURE")->getSIDouble(0);
96  watdentCT1_[regionIdx] = watdentRecord->getItem("EXPANSION_COEFF_LINEAR")->getSIDouble(0);
97  watdentCT2_[regionIdx] = watdentRecord->getItem("EXPANSION_COEFF_QUADRATIC")->getSIDouble(0);
98  }
99  }
100  }
101 
102  virtual void mu(const int n,
103  const int* pvtRegionIdx,
104  const double* p,
105  const double* T,
106  const double* z,
107  double* output_mu) const
108  {
109  if (watvisctTables_)
110  // TODO: temperature dependence for viscosity depending on z
111  OPM_THROW(std::runtime_error,
112  "temperature dependent viscosity as a function of z "
113  "is not yet implemented!");
114 
115  // compute the isothermal viscosity
116  isothermalPvt_->mu(n, pvtRegionIdx, p, T, z, output_mu);
117  }
118 
119  virtual void mu(const int n,
120  const int* pvtRegionIdx,
121  const double* p,
122  const double* T,
123  const double* r,
124  double* output_mu,
125  double* output_dmudp,
126  double* output_dmudr) const
127  {
128  // compute the isothermal viscosity and its derivatives
129  isothermalPvt_->mu(n, pvtRegionIdx, p, T, r, output_mu, output_dmudp, output_dmudr);
130 
131  if (!watvisctTables_)
132  // isothermal case
133  return;
134 
135  // temperature dependence
136  for (int i = 0; i < n; ++i) {
137  int tableIdx = getTableIndex_(pvtRegionIdx, i);
138 
139  // calculate the viscosity of the isothermal keyword for the reference
140  // pressure given by the VISCREF keyword.
141  double x = -pvtwViscosibility_[tableIdx]*(viscrefPress_[tableIdx] - pvtwRefPress_[tableIdx]);
142  double muRef = pvtwViscosity_[tableIdx]/(1.0 + x + 0.5*x*x);
143 
144  // compute the viscosity deviation due to temperature
145  double alpha;
146  {
147  const WatvisctTable& watVisctTable = watvisctTables_->getTable<WatvisctTable>(tableIdx);
148  double muWatvisct = watVisctTable.evaluate("Viscosity", T[i]);
149  alpha = muWatvisct/muRef;
150  }
151 
152  output_mu[i] *= alpha;
153  output_dmudp[i] *= alpha;
154  output_dmudr[i] *= alpha;
155  // TODO (?): derivative of viscosity w.r.t. temperature.
156  }
157  }
158 
159  virtual void mu(const int n,
160  const int* pvtRegionIdx,
161  const double* p,
162  const double* T,
163  const double* r,
164  const PhasePresence* cond,
165  double* output_mu,
166  double* output_dmudp,
167  double* output_dmudr) const
168  {
169  // compute the isothermal viscosity and its derivatives
170  isothermalPvt_->mu(n, pvtRegionIdx, p, T, r, cond, output_mu, output_dmudp, output_dmudr);
171 
172  if (!watvisctTables_)
173  // isothermal case
174  return;
175 
176  // temperature dependence
177  for (int i = 0; i < n; ++i) {
178  int tableIdx = getTableIndex_(pvtRegionIdx, i);
179 
180  // calculate the viscosity of the isothermal keyword for the reference
181  // pressure given by the VISCREF keyword.
182  double x = -pvtwViscosibility_[tableIdx]*(viscrefPress_[tableIdx] - pvtwRefPress_[tableIdx]);
183  double muRef = pvtwViscosity_[tableIdx]/(1.0 + x + 0.5*x*x);
184 
185  // compute the viscosity deviation due to temperature
186  double alpha;
187  {
188  const WatvisctTable& watVisctTable = watvisctTables_->getTable<WatvisctTable>(tableIdx);
189  double muWatvisct = watVisctTable.evaluate("Viscosity", T[i]);
190  alpha = muWatvisct/muRef;
191  }
192  output_mu[i] *= alpha;
193  output_dmudp[i] *= alpha;
194  output_dmudr[i] *= alpha;
195  // TODO (?): derivative of viscosity w.r.t. temperature.
196  }
197  }
198 
199  virtual void B(const int n,
200  const int* pvtRegionIdx,
201  const double* p,
202  const double* T,
203  const double* z,
204  double* output_B) const
205  {
206  if (watdentRefTemp_.empty()) {
207  // isothermal case
208  isothermalPvt_->B(n, pvtRegionIdx, p, T, z, output_B);
209  return;
210  }
211 
212  // This changes how the water density depends on pressure compared to what's
213  // used for the PVTW keyword, but it seems to be what Eclipse does. For
214  // details, see the documentation for the WATDENT keyword in the Eclipse RM.
215  for (int i = 0; i < n; ++i) {
216  int tableIdx = getTableIndex_(pvtRegionIdx, i);
217  double BwRef = pvtwRefB_[tableIdx];
218  double TRef = watdentRefTemp_[tableIdx];
219  double X = pvtwCompressibility_[tableIdx]*(p[i] - pvtwRefPress_[tableIdx]);
220  double cT1 = watdentCT1_[tableIdx];
221  double cT2 = watdentCT2_[tableIdx];
222  double Y = T[i] - TRef;
223  double Bw = BwRef*(1 - X)*(1 + cT1*Y + cT2*Y*Y);
224  output_B[i] = Bw;
225  }
226  }
227 
228  virtual void dBdp(const int n,
229  const int* pvtRegionIdx,
230  const double* p,
231  const double* T,
232  const double* z,
233  double* output_B,
234  double* output_dBdp) const
235  {
236  if (watdentRefTemp_.empty()) {
237  // isothermal case
238  isothermalPvt_->dBdp(n, pvtRegionIdx, p, T, z, output_B, output_dBdp);
239  return;
240  }
241 
242  // This changes how the water density depends on pressure. This is awkward,
243  // but it seems to be what Eclipse does. See the documentation for the
244  // WATDENT keyword in the Eclipse RM
245  for (int i = 0; i < n; ++i) {
246  int tableIdx = getTableIndex_(pvtRegionIdx, i);
247  double BwRef = pvtwRefB_[tableIdx];
248  double TRef = watdentRefTemp_[tableIdx];
249  double X = pvtwCompressibility_[tableIdx]*(p[i] - pvtwRefPress_[tableIdx]);
250  double cT1 = watdentCT1_[tableIdx];
251  double cT2 = watdentCT2_[tableIdx];
252  double Y = T[i] - TRef;
253  double Bw = BwRef*(1 - X)*(1 + cT1*Y + cT2*Y*Y);
254  output_B[i] = Bw;
255  }
256 
257  std::fill(output_dBdp, output_dBdp + n, 0.0);
258  }
259 
260  virtual void b(const int n,
261  const int* pvtRegionIdx,
262  const double* p,
263  const double* T,
264  const double* r,
265  double* output_b,
266  double* output_dbdp,
267  double* output_dbdr) const
268  {
269  if (watdentRefTemp_.empty()) {
270  // isothermal case
271  isothermalPvt_->b(n, pvtRegionIdx, p, T, r, output_b, output_dbdp, output_dbdr);
272  return;
273  }
274 
275  // This changes how the water density depends on pressure. This is awkward,
276  // but it seems to be what Eclipse does. See the documentation for the
277  // WATDENT keyword in the Eclipse RM
278  for (int i = 0; i < n; ++i) {
279  int tableIdx = getTableIndex_(pvtRegionIdx, i);
280  double BwRef = pvtwRefB_[tableIdx];
281  double TRef = watdentRefTemp_[tableIdx];
282  double X = pvtwCompressibility_[tableIdx]*(p[i] - pvtwRefPress_[tableIdx]);
283  double cT1 = watdentCT1_[tableIdx];
284  double cT2 = watdentCT2_[tableIdx];
285  double Y = T[i] - TRef;
286  double Bw = BwRef*(1 - X)*(1 + cT1*Y + cT2*Y*Y);
287  output_b[i] = 1.0/Bw;
288  }
289 
290  std::fill(output_dbdp, output_dbdp + n, 0.0);
291  std::fill(output_dbdr, output_dbdr + n, 0.0);
292  }
293 
294  virtual void b(const int n,
295  const int* pvtRegionIdx,
296  const double* p,
297  const double* T,
298  const double* r,
299  const PhasePresence* cond,
300  double* output_b,
301  double* output_dbdp,
302  double* output_dbdr) const
303  {
304  if (watdentRefTemp_.empty()) {
305  // isothermal case
306  isothermalPvt_->b(n, pvtRegionIdx, p, T, r, cond, output_b, output_dbdp, output_dbdr);
307  return;
308  }
309 
310  // This changes pressure dependence of the water density, but it seems to be
311  // what Eclipse does. See the documentation for the WATDENT keyword in the
312  // Eclipse RM
313  for (int i = 0; i < n; ++i) {
314  int tableIdx = getTableIndex_(pvtRegionIdx, i);
315  double BwRef = pvtwRefB_[tableIdx];
316  double TRef = watdentRefTemp_[tableIdx];
317  double X = pvtwCompressibility_[tableIdx]*(p[i] - pvtwRefPress_[tableIdx]);
318  double cT1 = watdentCT1_[tableIdx];
319  double cT2 = watdentCT2_[tableIdx];
320  double Y = T[i] - TRef;
321  double Bw = BwRef*(1 - X)*(1 + cT1*Y + cT2*Y*Y);
322  output_b[i] = 1.0/Bw;
323  }
324 
325  std::fill(output_dbdp, output_dbdp + n, 0.0);
326  std::fill(output_dbdr, output_dbdr + n, 0.0);
327 
328  }
329 
330  virtual void rsSat(const int n,
331  const int* pvtRegionIdx,
332  const double* p,
333  double* output_rsSat,
334  double* output_drsSatdp) const
335  {
336  isothermalPvt_->rsSat(n, pvtRegionIdx, p, output_rsSat, output_drsSatdp);
337  }
338 
339  virtual void rvSat(const int n,
340  const int* pvtRegionIdx,
341  const double* p,
342  double* output_rvSat,
343  double* output_drvSatdp) const
344  {
345  isothermalPvt_->rvSat(n, pvtRegionIdx, p, output_rvSat, output_drvSatdp);
346  }
347 
348  virtual void R(const int n,
349  const int* pvtRegionIdx,
350  const double* p,
351  const double* z,
352  double* output_R) const
353  {
354  isothermalPvt_->R(n, pvtRegionIdx, p, z, output_R);
355  }
356 
357  virtual void dRdp(const int n,
358  const int* pvtRegionIdx,
359  const double* p,
360  const double* z,
361  double* output_R,
362  double* output_dRdp) const
363  {
364  isothermalPvt_->dRdp(n, pvtRegionIdx, p, z, output_R, output_dRdp);
365  }
366 
367  private:
368  int getTableIndex_(const int* pvtTableIdx, int cellIdx) const
369  {
370  if (!pvtTableIdx)
371  return 0;
372  return pvtTableIdx[cellIdx];
373  }
374 
375  // the PVT propertied for the isothermal case
376  std::shared_ptr<const PvtInterface> isothermalPvt_;
377 
378  // The PVT properties needed for temperature dependence. We need to store one
379  // value per PVT region.
380  std::vector<double> viscrefPress_;
381 
382  std::vector<double> watdentRefTemp_;
383  std::vector<double> watdentCT1_;
384  std::vector<double> watdentCT2_;
385 
386  std::vector<double> pvtwRefPress_;
387  std::vector<double> pvtwRefB_;
388  std::vector<double> pvtwCompressibility_;
389  std::vector<double> pvtwViscosity_;
390  std::vector<double> pvtwViscosibility_;
391 
392  const TableContainer* watvisctTables_;
393  };
394 
395 }
396 
397 #endif
398 
Definition: BlackoilPhases.hpp:48
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
Definition: ThermalWaterPvtWrapper.hpp:260
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.
Definition: ThermalWaterPvtWrapper.hpp:199
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 and z.
Definition: ThermalWaterPvtWrapper.hpp:228
virtual void mu(const int n, const int *pvtRegionIdx, const double *p, const double *T, const double *r, const PhasePresence *cond, double *output_mu, double *output_dmudp, double *output_dmudr) const
Definition: ThermalWaterPvtWrapper.hpp:159
void initFromDeck(std::shared_ptr< const PvtInterface > isothermalPvt, Opm::DeckConstPtr deck, Opm::EclipseStateConstPtr eclipseState)
set the tables which specify the temperature dependence of the water viscosity
Definition: ThermalWaterPvtWrapper.hpp:42
virtual void b(const int n, const int *pvtRegionIdx, const double *p, const double *T, const double *r, const PhasePresence *cond, double *output_b, double *output_dbdp, double *output_dbdr) const
Definition: ThermalWaterPvtWrapper.hpp:294
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: ThermalWaterPvtWrapper.hpp:348
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...
Definition: ThermalWaterPvtWrapper.hpp:330
virtual void mu(const int n, const int *pvtRegionIdx, const double *p, const double *T, const double *r, double *output_mu, double *output_dmudp, double *output_dmudr) const
Definition: ThermalWaterPvtWrapper.hpp:119
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.
Definition: ThermalWaterPvtWrapper.hpp:102
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.
Definition: ThermalWaterPvtWrapper.hpp:339
ThermalWaterPvtWrapper()
Definition: ThermalWaterPvtWrapper.hpp:37
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.
Definition: ThermalWaterPvtWrapper.hpp:357
Definition: PvtInterface.hpp:30
Definition: ThermalWaterPvtWrapper.hpp:34