BlackoilCo2PVT.hpp
Go to the documentation of this file.
1 /*
2  Copyright 2010 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_BLACKOILCO2PVT_HEADER_INCLUDED
21 #define OPM_BLACKOILCO2PVT_HEADER_INCLUDED
22 #define OPM_BLACKOILPVT_HEADER_INCLUDED
23 
24 #define OPM_DEPRECATED __attribute__((deprecated))
25 #define OPM_DEPRECATED_MSG(msg) __attribute__((deprecated))
26 
27 #include <opm/material/fluidsystems/BrineCO2FluidSystem.hpp>
28 #include <opm/material/fluidstates/CompositionalFluidState.hpp>
29 #include <opm/material/constraintsolvers/MiscibleMultiPhaseComposition.hpp>
30 #include <opm/material/constraintsolvers/ComputeFromReferencePhase.hpp>
31 
32 #include <opm/common/Exceptions.hpp>
33 #include <opm/common/ErrorMacros.hpp>
37 
38 #include <opm/parser/eclipse/Deck/Deck.hpp>
39 
40 #include <string>
41 #include <iostream>
42 #include <fstream>
43 
44 
45 namespace Opm
46 {
48 {
49 public:
50 
51  typedef Opm::FluidSystems::BrineCO2</*Scalar=*/double, Opm::Benchmark3::CO2Tables> FluidSystem;
52  typedef Opm::CompositionalFluidState<double, FluidSystem> CompositionalFluidState;
53 
54  void init(Opm::DeckConstPtr deck);
55 
56  void generateBlackOilTables(double temperature);
57 
58  double getViscosity(double press,
59  const CompVec& surfvol,
60  PhaseIndex phase) const;
61  CompVec surfaceDensities() const;
62  double getSaturation(double press,
63  const CompVec& surfvol,
64  PhaseIndex phase) const;
65  double B (double press,
66  const CompVec& surfvol,
67  PhaseIndex phase) const;
68  double dBdp(double press,
69  const CompVec& surfvol,
70  PhaseIndex phase) const;
71  double R (double press,
72  const CompVec& surfvol,
73  PhaseIndex phase) const;
74  double dRdp(double press,
75  const CompVec& surfvol,
76  PhaseIndex phase) const;
77 
78  void getViscosity(const std::vector<PhaseVec>& pressures,
79  const std::vector<CompVec>& surfvol,
80  std::vector<PhaseVec>& output) const;
81  void B(const std::vector<PhaseVec>& pressures,
82  const std::vector<CompVec>& surfvol,
83  std::vector<PhaseVec>& output) const;
84  void dBdp(const std::vector<PhaseVec>& pressures,
85  const std::vector<CompVec>& surfvol,
86  std::vector<PhaseVec>& output_B,
87  std::vector<PhaseVec>& output_dBdp) const;
88  void R(const std::vector<PhaseVec>& pressures,
89  const std::vector<CompVec>& surfvol,
90  std::vector<PhaseVec>& output) const;
91  void dRdp(const std::vector<PhaseVec>& pressures,
92  const std::vector<CompVec>& surfvol,
93  std::vector<PhaseVec>& output_R,
94  std::vector<PhaseVec>& output_dRdp) const;
95 
96 private:
97  CompVec surfaceDensities_;
98  FluidSystem brineCo2_;
99  double temperature_;
100  struct SubState
101  {
102  Dune::FieldVector<double,2> density;
103  Dune::FieldMatrix<double,2,2> massfrac;
104  Dune::FieldVector<double,2> phaseVolume;
105  Dune::FieldVector<double,2> phaseViscosity;
106  double saturation;
107  };
108  void computeState(SubState& ss, double zBrine, double zCO2, double pressure) const;
109  enum {
110  wPhase = FluidSystem::liquidPhaseIdx,
111  nPhase = FluidSystem::gasPhaseIdx,
112 
113  wComp = FluidSystem::BrineIdx,
114  nComp = FluidSystem::CO2Idx
115  };
116 
117 }; // class BlackoilCo2PVT
118 
119 // ------------ Method implementations --------------
120 
121 void BlackoilCo2PVT::init(Opm::DeckConstPtr deck)
122 {
123  surfaceDensities_[Water] = 1000.;
124  surfaceDensities_[Gas] = 2.0;
125  surfaceDensities_[Oil] = 1000.;
126 
127  temperature_ = 300.;
128 
129  brineCo2_.init();
130 }
131 
133 {
134  std::cout << "\n Generating pvt tables for the eclipse black oil formulation\n using the oil component as brine and the gas component as co_2." << std::endl;
135  if (std::fabs(temperature-400.) > 100.0) {
136  std::cout << "T=" << temperature << " is outside the allowed range [300,500] Kelvin" << std::endl;
137  exit(-1);
138  }
139 
140  temperature_ = temperature;
141 
142  CompVec z;
143  z[Water] = 0.0;
144  z[Oil] = 1.0;
145  z[Gas] = 1.0e6;
146 
147  std::ofstream black("blackoil_pvt");
148  black.precision(8);
149  black << std::fixed << std::showpoint;
150 
151  const double pMin=150.e5;
152  const double pMax=400.e5;
153  const unsigned int np=11;
154  std::vector<double> pValue(np+1,0.0);
155  std::vector<double> rs(np+1,0.0);
156 
157  pValue[0] = 101330.0;
158  rs[0] = 0.0;
159 
160  // Buble points
161  z[Gas] = 1.0e4;
162  for (unsigned int i=0; i<np; ++i) {
163  pValue[i+1] = pMin + i*(pMax-pMin)/(np-1);
164  rs[i+1] = R(pValue[i+1], z, Liquid);
165  }
166 
167  const unsigned int np_fill_in = 10;
168  const double dr = (rs[1] - rs[0])/(np_fill_in+1);
169  const double dp = (pValue[1] - pValue[0])/(np_fill_in+1);
170  rs.insert(rs.begin()+1,np_fill_in,0.0);
171  pValue.insert(pValue.begin()+1,np_fill_in,0.0);
172  for (unsigned int i=1; i<=np_fill_in; ++i) {
173  rs[i] = rs[i-1] + dr;
174  pValue[i] = pValue[i-1] + dp;
175  }
176 
177  // Brine with dissolved co_2 ("live oil")
178  black << "PVTO\n";
179  black << "-- Rs Pbub Bo Vo\n";
180  black << "-- sm3/sm3 barsa rm3/sm3 cP\n";
181 
182  // Undersaturated
183  for (unsigned int i=0; i<np+np_fill_in; ++i) {
184  z[Gas] = rs[i];
185  black << std::setw(14) << rs[i];
186  for (unsigned int j=i; j<np+1+np_fill_in; ++j) {
187  if (j<=np_fill_in) {
188  if (j==i) black << std::setw(14) << pValue[j]*1.e-5 << std::setw(14) << 1.0-j*0.001 << std::setw(14) << 1.06499;
189  continue;
190  }
191  if (j>i) black << std::endl << std::setw(14) << ' ';
192  black << std::setw(14) << pValue[j]*1.e-5
193  << std::setw(14) << B(pValue[j], z, Liquid)
194  << std::setw(14) << getViscosity(pValue[j], z, Liquid)*1.e3;
195  }
196  black << " /" << std::endl;
197  }
198  black << "/ " << std::endl;
199 
200  // We provide tables for co_2 both with and without dissolved water:
201 
202  // Co_2 neglecting dissolved water ("dry gas")
203  black << "\nPVDG\n";
204  black << "-- Pg Bg Vg\n";
205  black << "-- barsa rm3/sm3 cP\n";
206 
207  for (unsigned int i=0; i<np; ++i) {
208  z[Oil] = 0.0;
209  z[Gas] = 1.0;
210  black << std::setw(14) << pValue[i+np_fill_in+1]*1.e-5
211  << std::setw(14) << B(pValue[i+np_fill_in+1], z, Vapour)
212  << std::setw(14) << getViscosity(pValue[i+np_fill_in+1], z, Vapour)*1.e3
213  << std::endl;
214  }
215  black << "/ " << std::endl;
216 
217  // Co_2 with dissolved water ("wet gas")
218  black << "\nPVTG\n";
219  black << "-- Pg Rv Bg Vg\n";
220  black << "-- barsa sm3/sm3 rm3/sm3 cP\n";
221  for (unsigned int i=0; i<np; ++i) {
222  z[Oil] = 1000.0;
223  z[Gas] = 1.0;
224  black << std::setw(14) << pValue[i+np_fill_in+1]*1.e-5
225  << std::setw(14) << R(pValue[i+np_fill_in+1], z, Vapour)
226  << std::setw(14) << B(pValue[i+np_fill_in+1], z, Vapour)
227  << std::setw(14) << getViscosity(pValue[i+np_fill_in+1], z, Vapour)*1.e3
228  << std::endl;
229  z[Oil] = 0.0;
230  black << std::setw(14) << ' '
231  << std::setw(14) << R(pValue[i+np_fill_in+1], z, Vapour)
232  << std::setw(14) << B(pValue[i+np_fill_in+1], z, Vapour)
233  << std::setw(14) << getViscosity(pValue[i+np_fill_in+1], z, Vapour)*1.e3
234  << " /" << std::endl;
235  }
236  black << "/ " << std::endl;
237  black << std::endl;
238  std::cout << " Pvt tables for temperature=" << temperature << " Kelvin is written to file blackoil_pvt. " << std::endl;
239  std::cout << " NOTE that the file contains tables for both PVDG (dry gas) and PVTG (wet gas)." << std::endl;
240 }
241 
243 {
244  return surfaceDensities_;
245 }
246 
247 double BlackoilCo2PVT::getViscosity(double press, const CompVec& surfvol, PhaseIndex phase) const
248 {
249  SubState ss;
250  computeState(ss, surfvol[Oil], surfvol[Gas], press);
251  switch(phase) {
252  case Aqua: return 1.0e-10;
253  case Liquid: return ss.phaseViscosity[wPhase];
254  case Vapour: return ss.phaseViscosity[nPhase];
255  };
256  return 1.0e-10;
257 }
258 
259 
260 double BlackoilCo2PVT::getSaturation(double press, const CompVec& surfvol, PhaseIndex phase) const
261 {
262  SubState ss;
263  computeState(ss, surfvol[Oil], surfvol[Gas], press);
264  switch(phase) {
265  case Aqua: return 0.0;
266  case Liquid: return ss.saturation;
267  case Vapour: return 1.0 - ss.saturation;
268  };
269  return 0.0;
270 }
271 
272 
273 double BlackoilCo2PVT::B(double press, const CompVec& surfvol, PhaseIndex phase) const
274 {
275  SubState ss;
276  computeState(ss, surfvol[Oil], surfvol[Gas], press);
277  switch(phase) {
278  case Aqua: return 1.0;
279  case Liquid: return surfaceDensities_[Oil]/(ss.massfrac[wPhase][wComp]*ss.density[wPhase]+1.0e-10);
280  case Vapour: return surfaceDensities_[Gas]/(ss.massfrac[nPhase][nComp]*ss.density[nPhase]+1.0e-10);
281  };
282  return 1.0;
283 }
284 
285 double BlackoilCo2PVT::dBdp(double press, const CompVec& surfvol, PhaseIndex phase) const
286 {
287  const double dp = 100.;
288  return (B(press+dp,surfvol,phase)-B(press,surfvol,phase))/dp;
289 }
290 
291 double BlackoilCo2PVT::R(double press, const CompVec& surfvol, PhaseIndex phase) const
292 {
293  SubState ss;
294  computeState(ss, surfvol[Oil], surfvol[Gas], press);
295  switch(phase) {
296  case Aqua: return 0.0;
297  case Liquid: return (ss.massfrac[wPhase][nComp]*surfaceDensities_[Oil])/(ss.massfrac[wPhase][wComp]*surfaceDensities_[Gas]+1.0e-10);
298  case Vapour: return (ss.massfrac[nPhase][wComp]*surfaceDensities_[Gas])/(ss.massfrac[nPhase][nComp]*surfaceDensities_[Oil]+1.0e-10);
299  };
300  return 0.0;
301 }
302 
303 double BlackoilCo2PVT::dRdp(double press, const CompVec& surfvol, PhaseIndex phase) const
304 {
305  const double dp = 100.;
306  return (R(press+dp,surfvol,phase)-R(press,surfvol,phase))/dp;
307 }
308 
309 void BlackoilCo2PVT::getViscosity(const std::vector<PhaseVec>& pressures,
310  const std::vector<CompVec>& surfvol,
311  std::vector<PhaseVec>& output) const
312 {
313  int num = pressures.size();
314  output.resize(num);
315  SubState ss;
316  for (int i = 0; i < num; ++i) {
317  computeState(ss, surfvol[i][Oil], surfvol[i][Gas], pressures[i][Liquid]);
318  output[i][Aqua] = 1.0e-10;
319  output[i][Liquid] = ss.phaseViscosity[wPhase];
320  output[i][Vapour] = ss.phaseViscosity[nPhase];
321  }
322 }
323 
324 void BlackoilCo2PVT::B(const std::vector<PhaseVec>& pressures,
325  const std::vector<CompVec>& surfvol,
326  std::vector<PhaseVec>& output) const
327 {
328  int num = pressures.size();
329  output.resize(num);
330  SubState ss;
331  for (int i = 0; i < num; ++i) {
332  computeState(ss, surfvol[i][Oil], surfvol[i][Gas], pressures[i][Liquid]);
333  output[i][Aqua] = 1.0;
334  output[i][Liquid] = surfaceDensities_[Oil]/(ss.massfrac[wPhase][wComp]*ss.density[wPhase]+1.0e-10);
335  output[i][Vapour] = surfaceDensities_[Gas]/(ss.massfrac[nPhase][nComp]*ss.density[nPhase]+1.0e-10);
336  }
337 }
338 
339 void BlackoilCo2PVT::dBdp(const std::vector<PhaseVec>& pressures,
340  const std::vector<CompVec>& surfvol,
341  std::vector<PhaseVec>& output_B,
342  std::vector<PhaseVec>& output_dBdp) const
343 {
344  int num = pressures.size();
345  output_B.resize(num);
346  output_dBdp.resize(num);
347  SubState ss;
348  const double dp = 100.;
349  for (int i = 0; i < num; ++i) {
350  computeState(ss, surfvol[i][Oil], surfvol[i][Gas], pressures[i][Liquid]);
351  output_B[i][Aqua] = 1.0;
352  output_B[i][Liquid] = surfaceDensities_[Oil]/(ss.massfrac[wPhase][wComp]*ss.density[wPhase]+1.0e-10);
353  output_B[i][Vapour] = surfaceDensities_[Gas]/(ss.massfrac[nPhase][nComp]*ss.density[nPhase]+1.0e-10);
354  computeState(ss, surfvol[i][Oil], surfvol[i][Gas], pressures[i][Liquid]+dp);
355  output_dBdp[i][Aqua] = 0.0;
356  output_dBdp[i][Liquid] = (surfaceDensities_[Oil]/(ss.massfrac[wPhase][wComp]*ss.density[wPhase]+1.0e-10) - output_B[i][Liquid])/dp;
357  output_dBdp[i][Vapour] = (surfaceDensities_[Gas]/(ss.massfrac[nPhase][nComp]*ss.density[nPhase]+1.0e-10) - output_B[i][Vapour])/dp;
358  }
359 }
360 
361 void BlackoilCo2PVT::R(const std::vector<PhaseVec>& pressures,
362  const std::vector<CompVec>& surfvol,
363  std::vector<PhaseVec>& output) const
364 {
365  int num = pressures.size();
366  output.resize(num);
367  SubState ss;
368  for (int i = 0; i < num; ++i) {
369  computeState(ss, surfvol[i][Oil], surfvol[i][Gas], pressures[i][Liquid]);
370  output[i][Aqua] = 0.0;
371  output[i][Liquid] = (ss.massfrac[wPhase][nComp]*surfaceDensities_[Oil])/(ss.massfrac[wPhase][wComp]*surfaceDensities_[Gas]+1.0e-10);
372  output[i][Vapour] = (ss.massfrac[nPhase][wComp]*surfaceDensities_[Gas])/(ss.massfrac[nPhase][nComp]*surfaceDensities_[Oil]+1.0e-10);
373  }
374 }
375 
376 void BlackoilCo2PVT::dRdp(const std::vector<PhaseVec>& pressures,
377  const std::vector<CompVec>& surfvol,
378  std::vector<PhaseVec>& output_R,
379  std::vector<PhaseVec>& output_dRdp) const
380 {
381  int num = pressures.size();
382  output_R.resize(num);
383  output_dRdp.resize(num);
384  SubState ss;
385  const double dp = 100.;
386  for (int i = 0; i < num; ++i) {
387  computeState(ss, surfvol[i][Oil], surfvol[i][Gas], pressures[i][Liquid]);
388  output_R[i][Aqua] = 0.0;
389  output_R[i][Liquid] = (ss.massfrac[wPhase][nComp]*surfaceDensities_[Oil])/(ss.massfrac[wPhase][wComp]*surfaceDensities_[Gas]+1.0e-10);
390  output_R[i][Vapour] = (ss.massfrac[nPhase][wComp]*surfaceDensities_[Gas])/(ss.massfrac[nPhase][nComp]*surfaceDensities_[Oil]+1.0e-10);
391  computeState(ss, surfvol[i][Oil], surfvol[i][Gas], pressures[i][Liquid]+dp);
392  output_dRdp[i][Aqua] = 0.0;
393  output_dRdp[i][Liquid] = ((ss.massfrac[wPhase][nComp]*surfaceDensities_[Oil])/(ss.massfrac[wPhase][wComp]*surfaceDensities_[Gas]+1.0e-10) - output_R[i][Liquid])/dp;
394  output_dRdp[i][Vapour] = ((ss.massfrac[nPhase][wComp]*surfaceDensities_[Gas])/(ss.massfrac[nPhase][nComp]*surfaceDensities_[Oil]+1.0e-10) - output_R[i][Vapour])/dp;
395  }
396 }
397 
398 void BlackoilCo2PVT::computeState(BlackoilCo2PVT::SubState& ss, double zBrine, double zCO2, double pressure) const
399 {
400 
401  CompositionalFluidState fluidState;
402  fluidState.setTemperature(temperature_);
403  fluidState.setPressure(wPhase, pressure);
404  fluidState.setPressure(nPhase, pressure);
405 
406  double massH20 = surfaceDensities_[Oil]*zBrine;
407  double massCO2 = surfaceDensities_[Gas]*zCO2;
408 
409  // A priori, assume presence of both phases
410  typename FluidSystem::ParameterCache paramCache;
411  typedef Opm::MiscibleMultiPhaseComposition</*Scalar=*/double, FluidSystem> MMPC;
412  MMPC::solve(fluidState, paramCache, /*setViscosity=*/false, /*setEnthalpy=*/false);
413  ss.density[wPhase] = fluidState.density(wPhase);
414  ss.density[nPhase] = fluidState.density(nPhase);
415  ss.massfrac[wPhase][nComp] = fluidState.massFraction(wPhase, nComp);
416  ss.massfrac[nPhase][wComp] = fluidState.massFraction(nPhase, wComp);
417  ss.massfrac[wPhase][wComp] = 1.0 - ss.massfrac[wPhase][nComp];
418  ss.massfrac[nPhase][nComp] = 1.0 - ss.massfrac[nPhase][wComp];
419 
420  double detX = ss.massfrac[wPhase][wComp]*ss.massfrac[nPhase][nComp]-ss.massfrac[wPhase][nComp]*ss.massfrac[nPhase][wComp];
421  ss.phaseVolume[wPhase] = (massH20*ss.massfrac[nPhase][nComp] - massCO2*ss.massfrac[nPhase][wComp])/(ss.density[wPhase]*detX);
422  ss.phaseVolume[nPhase] = (massCO2*ss.massfrac[wPhase][wComp] - massH20*ss.massfrac[wPhase][nComp])/(ss.density[nPhase]*detX);
423 
424  // Determine number of phase
425  if (ss.phaseVolume[wPhase] > 0.0 && ss.phaseVolume[nPhase] > 0.0) { // Both phases
426  ss.saturation = ss.phaseVolume[wPhase]/(ss.phaseVolume[wPhase]+ss.phaseVolume[nPhase]);
427  fluidState.setSaturation(wPhase, ss.saturation);
428  fluidState.setSaturation(nPhase, 1.0 - ss.saturation);
429  }
430  else if (ss.phaseVolume[wPhase] <= 0.0) { // Wetting phase only
431  ss.saturation = 0.0;
432  // Gas phase:
433  ss.massfrac[nPhase][nComp] = massCO2/(massCO2+massH20);
434  ss.massfrac[nPhase][wComp] = 1.0 - ss.massfrac[nPhase][nComp];
435  double M1 = FluidSystem::molarMass(wComp);
436  double M2 = FluidSystem::molarMass(nComp);
437  double avgMolarMass = M1*M2/(M2 + ss.massfrac[nPhase][nComp]*(M1 - M2));
438  fluidState.setMoleFraction(nPhase, nComp, ss.massfrac[nPhase][nComp]*avgMolarMass/M2);
439  fluidState.setMoleFraction(nPhase, wComp, ss.massfrac[nPhase][wComp]*avgMolarMass/M1);
440  ss.density[nPhase] = brineCo2_.density(fluidState, paramCache, nPhase);
441  fluidState.setDensity(nPhase, ss.density[nPhase]);
442  ss.phaseVolume[nPhase] = (massH20+massCO2)/ss.density[nPhase];
443  fluidState.setSaturation(nPhase, 1.0 - ss.saturation);
444  // Virtual properties of non-existing liquid phase:
445  paramCache.updatePhase(fluidState, /*phaseIdx=*/nPhase);
446  typedef Opm::ComputeFromReferencePhase</*Scalar=*/double, FluidSystem> CFRP;
447  CFRP::solve(fluidState,
448  paramCache,
449  /*refPhaseIdx=*/nPhase,
450  /*setViscosity=*/false,
451  /*setEnthalpy=*/false);
452  ss.massfrac[wPhase][wComp] = fluidState.massFraction(wPhase, wComp);
453  ss.massfrac[wPhase][nComp] = fluidState.massFraction(wPhase, nComp);
454  ss.density[wPhase] = fluidState.density(wPhase);
455  ss.phaseVolume[wPhase] = 0.0;
456  fluidState.setSaturation(wPhase, ss.saturation);
457  }
458  else if (ss.phaseVolume[nPhase] <= 0.0) { // Non-wetting phase only
459  ss.saturation = 1.0;
460  // Liquid phase:
461  ss.massfrac[wPhase][wComp] = massH20/(massCO2+massH20);
462  ss.massfrac[wPhase][nComp] = 1.0 - ss.massfrac[wPhase][wComp];
463  double M1 = FluidSystem::molarMass(wComp);
464  double M2 = FluidSystem::molarMass(nComp);
465  double avgMolarMass = M1*M2/(M2 + ss.massfrac[wPhase][nComp]*(M1 - M2));
466  fluidState.setMoleFraction(wPhase, nComp, ss.massfrac[wPhase][nComp]*avgMolarMass/M2);
467  fluidState.setMoleFraction(wPhase, wComp, ss.massfrac[wPhase][wComp]*avgMolarMass/M1);
468  ss.density[wPhase] = brineCo2_.density(fluidState, paramCache, wPhase);
469  fluidState.setDensity(wPhase, ss.density[wPhase]);
470  ss.phaseVolume[wPhase] = (massH20+massCO2)/ss.density[wPhase];
471  fluidState.setSaturation(wPhase, ss.saturation);
472  // Virtual properties of non-existing gas phase:
473  paramCache.updatePhase(fluidState, /*phaseIdx=*/nPhase);
474  typedef ComputeFromReferencePhase</*Scalar=*/double, FluidSystem> CFRP;
475  CFRP::solve(fluidState,
476  paramCache,
477  /*refPhaseIdx=*/wPhase,
478  /*setViscosity=*/false,
479  /*setEnthalpy=*/false);
480  ss.massfrac[nPhase][nComp] = fluidState.massFraction(nPhase, nComp);
481  ss.massfrac[nPhase][wComp] = fluidState.massFraction(nPhase, wComp);
482  ss.density[nPhase] = fluidState.density(nPhase);
483  ss.phaseVolume[nPhase] = 0.0;
484  fluidState.setSaturation(nPhase, 1.0 - ss.saturation);
485  }
486 
487  ss.phaseViscosity[wPhase] = brineCo2_.viscosity(fluidState, paramCache, wPhase);
488  ss.phaseViscosity[nPhase] = brineCo2_.viscosity(fluidState, paramCache, nPhase);
489 
490 }
491 
492 
493 
495 
496 } // Opm
497 
498 
499 #endif // OPM_BLACKOILCO2PVT_HEADER_INCLUDED
double dRdp(double press, const CompVec &surfvol, PhaseIndex phase) const
Definition: BlackoilCo2PVT.hpp:303
Definition: BlackoilCo2PVT.hpp:47
Definition: BlackoilFluid.hpp:31
Dune::FieldVector< Scalar, numComponents > CompVec
Definition: BlackoilDefs.hpp:40
double B(double press, const CompVec &surfvol, PhaseIndex phase) const
Definition: BlackoilCo2PVT.hpp:273
Definition: BlackoilDefs.hpp:36
Definition: benchmark3co2tables.hh:25261
double dBdp(double press, const CompVec &surfvol, PhaseIndex phase) const
Definition: BlackoilCo2PVT.hpp:285
double R(double press, const CompVec &surfvol, PhaseIndex phase) const
Definition: BlackoilCo2PVT.hpp:291
Definition: BlackoilDefs.hpp:37
Definition: BlackoilDefs.hpp:36
Provides the class with the tabulated values of CO2 for the benchmark3 problem.
Opm::CompositionalFluidState< double, FluidSystem > CompositionalFluidState
Definition: BlackoilCo2PVT.hpp:52
void generateBlackOilTables(double temperature)
Definition: BlackoilCo2PVT.hpp:132
Definition: BlackoilDefs.hpp:30
BlackoilCo2PVT BlackoilPVT
Definition: BlackoilCo2PVT.hpp:494
void init(Opm::DeckConstPtr deck)
Definition: BlackoilCo2PVT.hpp:121
double getSaturation(double press, const CompVec &surfvol, PhaseIndex phase) const
Definition: BlackoilCo2PVT.hpp:260
Definition: BlackoilDefs.hpp:37
double getViscosity(double press, const CompVec &surfvol, PhaseIndex phase) const
Definition: BlackoilCo2PVT.hpp:247
Definition: BlackoilDefs.hpp:36
Definition: BlackoilDefs.hpp:37
CompVec surfaceDensities() const
Definition: BlackoilCo2PVT.hpp:242
Opm::FluidSystems::BrineCO2< double, Opm::Benchmark3::CO2Tables > FluidSystem
Definition: BlackoilCo2PVT.hpp:51
PhaseIndex
Definition: BlackoilDefs.hpp:37