blackoilprimaryvariables.hh
Go to the documentation of this file.
1 // -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2 // vi: set et ts=4 sw=4 sts=4:
3 /*
4  Copyright (C) 2011-2014 by Andreas Lauser
5 
6  This file is part of the Open Porous Media project (OPM).
7 
8  OPM is free software: you can redistribute it and/or modify
9  it under the terms of the GNU General Public License as published by
10  the Free Software Foundation, either version 2 of the License, or
11  (at your option) any later version.
12 
13  OPM is distributed in the hope that it will be useful,
14  but WITHOUT ANY WARRANTY; without even the implied warranty of
15  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16  GNU General Public License for more details.
17 
18  You should have received a copy of the GNU General Public License
19  along with OPM. If not, see <http://www.gnu.org/licenses/>.
20 */
26 #ifndef EWOMS_BLACK_OIL_PRIMARY_VARIABLES_HH
27 #define EWOMS_BLACK_OIL_PRIMARY_VARIABLES_HH
28 
29 #include "blackoilproperties.hh"
30 
32 
33 #include <dune/common/fvector.hh>
34 
35 #include <opm/material/constraintsolvers/NcpFlash.hpp>
36 #include <opm/material/fluidstates/CompositionalFluidState.hpp>
37 
38 namespace Ewoms {
44 template <class TypeTag>
46 {
47  typedef FvBasePrimaryVariables<TypeTag> ParentType;
48 
49  typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
50  typedef typename GET_PROP_TYPE(TypeTag, Evaluation) Evaluation;
51  typedef typename GET_PROP_TYPE(TypeTag, Indices) Indices;
52  typedef typename GET_PROP_TYPE(TypeTag, FluidSystem) FluidSystem;
53  typedef typename GET_PROP_TYPE(TypeTag, MaterialLaw) MaterialLaw;
54  typedef typename GET_PROP_TYPE(TypeTag, MaterialLawParams) MaterialLawParams;
55 
56  // number of equations
57  enum { numEq = GET_PROP_VALUE(TypeTag, NumEq) };
58 
59  // primary variable indices
60  enum { gasPressureIdx = Indices::gasPressureIdx };
61  enum { waterSaturationIdx = Indices::waterSaturationIdx };
62  enum { switchIdx = Indices::switchIdx };
63 
64  // phase indices from the fluid system
65  enum { numPhases = GET_PROP_VALUE(TypeTag, NumPhases) };
66  enum { gasPhaseIdx = FluidSystem::gasPhaseIdx };
67  enum { waterPhaseIdx = FluidSystem::waterPhaseIdx };
68  enum { oilPhaseIdx = FluidSystem::oilPhaseIdx };
69 
70  // component indices from the fluid system
71  enum { numComponents = GET_PROP_VALUE(TypeTag, NumComponents) };
72  enum { gasCompIdx = FluidSystem::gasCompIdx };
73  enum { waterCompIdx = FluidSystem::waterCompIdx };
74  enum { oilCompIdx = FluidSystem::oilCompIdx };
75 
76  typedef typename Opm::MathToolbox<Evaluation> Toolbox;
77  typedef Dune::FieldVector<Scalar, numComponents> ComponentVector;
78 
79  static_assert(numPhases == 3, "The black-oil model assumes three phases!");
80  static_assert(numComponents == 3, "The black-oil model assumes three components!");
81 
82 public:
83  enum SwitchingVarMeaning {
84  GasSaturation,
85  OilMoleFractionInGas,
86  GasMoleFractionInOil
87  };
88 
90  : ParentType()
91  {
92  Valgrind::SetUndefined(*this);
93  pvtRegionIdx_ = 0;
94  temperature_ = FluidSystem::surfaceTemperature;
95  }
96 
100  BlackOilPrimaryVariables(Scalar value)
101  : ParentType(value)
102  {
103  Valgrind::SetUndefined(switchingVarMeaning_);
104  pvtRegionIdx_ = 0;
105  temperature_ = FluidSystem::surfaceTemperature;
106  }
107 
112  : ParentType(value)
113  , temperature_(FluidSystem::surfaceTemperature)
114  , switchingVarMeaning_(value.switchingVarMeaning_)
115  , pvtRegionIdx_(value.pvtRegionIdx_)
116  { }
117 
118  void setPvtRegionIndex(int value)
119  { pvtRegionIdx_ = value; }
120 
121  unsigned pvtRegionIndex() const
122  { return pvtRegionIdx_; }
123 
124  SwitchingVarMeaning switchingVarMeaning() const
125  { return switchingVarMeaning_; }
126 
127  void setSwitchingVarMeaning(SwitchingVarMeaning newMeaning)
128  { switchingVarMeaning_ = newMeaning; }
129 
133  template <class FluidState>
134  void assignMassConservative(const FluidState &fluidState,
135  const MaterialLawParams &matParams,
136  bool isInEquilibrium = false)
137  {
138  typedef typename std::remove_reference<typename FluidState::Scalar>::type ConstEvaluation;
139  typedef typename std::remove_const<ConstEvaluation>::type FsEvaluation;
140  typedef typename Opm::MathToolbox<FsEvaluation> FsToolbox;
141 
142 #ifndef NDEBUG
143  // make sure the temperature is the same in all fluid phases
144  for (int phaseIdx = 1; phaseIdx < numPhases; ++phaseIdx) {
145  Valgrind::CheckDefined(fluidState.temperature(0));
146  Valgrind::CheckDefined(fluidState.temperature(phaseIdx));
147 
148  assert(fluidState.temperature(0) == fluidState.temperature(phaseIdx));
149  }
150 #endif // NDEBUG
151 
152  // for the equilibrium case, we don't need complicated
153  // computations.
154  if (isInEquilibrium) {
155  assignNaive(fluidState);
156  return;
157  }
158 
159  // If your compiler bails out here, you're probably not using a suitable black
160  // oil fluid system.
161  typename FluidSystem::ParameterCache paramCache;
162  paramCache.setRegionIndex(pvtRegionIdx_);
163 
164  // create a mutable fluid state with well defined densities based on the input
165  typedef Opm::NcpFlash<Scalar, FluidSystem> NcpFlash;
166  typedef Opm::CompositionalFluidState<Scalar, FluidSystem> FlashFluidState;
167  FlashFluidState fsFlash;
168  fsFlash.setTemperature(FsToolbox::value(fluidState.temperature(/*phaseIdx=*/0)));
169  for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
170  fsFlash.setPressure(phaseIdx, FsToolbox::value(fluidState.pressure(phaseIdx)));
171  fsFlash.setSaturation(phaseIdx, FsToolbox::value(fluidState.saturation(phaseIdx)));
172  for (int compIdx = 0; compIdx < numComponents; ++compIdx)
173  fsFlash.setMoleFraction(phaseIdx, compIdx, FsToolbox::value(fluidState.moleFraction(phaseIdx, compIdx)));
174  }
175 
176  for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
177  Scalar rho = FluidSystem::template density<FlashFluidState, Scalar>(fsFlash, paramCache, phaseIdx);
178  fsFlash.setDensity(phaseIdx, rho);
179  }
180 
181  // calculate the "global molarities"
182  ComponentVector globalMolarities(0.0);
183  for (int compIdx = 0; compIdx < numComponents; ++compIdx) {
184  for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
185  globalMolarities[compIdx] +=
186  fsFlash.saturation(phaseIdx) * fsFlash.molarity(phaseIdx, compIdx);
187  }
188  }
189 
190  // use a flash calculation to calculate a fluid state in
191  // thermodynamic equilibrium
192 
193  // run the flash calculation
194  NcpFlash::template solve<MaterialLaw>(fsFlash, paramCache, matParams, globalMolarities);
195 
196  // use the result to assign the primary variables
197  assignNaive(fsFlash);
198  }
199 
203  template <class FluidState>
204  void assignNaive(const FluidState &fluidState)
205  {
206  typedef typename std::remove_reference<typename FluidState::Scalar>::type ConstEvaluation;
207  typedef typename std::remove_const<ConstEvaluation>::type FsEvaluation;
208  typedef typename Opm::MathToolbox<FsEvaluation> FsToolbox;
209 
210  temperature_ = FsToolbox::value(fluidState.temperature(gasPhaseIdx));
211 
212  // the pressure of the first phase and the saturation of water
213  // are always primary variables
214  (*this)[gasPressureIdx] = FsToolbox::value(fluidState.pressure(gasPhaseIdx));
215  (*this)[waterSaturationIdx] = FsToolbox::value(fluidState.saturation(waterPhaseIdx));
216 
217  bool gasPresent = (fluidState.saturation(gasPhaseIdx) > 0.0);
218  bool oilPresent = (fluidState.saturation(oilPhaseIdx) > 0.0);
219 
220  if ((gasPresent && oilPresent) || (!gasPresent && !oilPresent))
221  // gas and oil
222  switchingVarMeaning_ = GasSaturation;
223  else if (oilPresent) {
224  // only oil
225  if (FluidSystem::enableDissolvedGas())
226  switchingVarMeaning_ = GasMoleFractionInOil;
227  else
228  switchingVarMeaning_ = GasSaturation;
229  }
230  else {
231  assert(gasPresent);
232  // only gas
233  if (FluidSystem::enableVaporizedOil())
234  switchingVarMeaning_ = OilMoleFractionInGas;
235  else
236  switchingVarMeaning_ = GasSaturation;
237  }
238 
239  // depending on the phases present, we use either Sg, x_o^G or x_g^O as the
240  // switching primary variable
241  if (switchingVarMeaning() == GasSaturation)
242  (*this)[switchIdx] = FsToolbox::value(fluidState.saturation(gasPhaseIdx));
243  else if (switchingVarMeaning() == GasMoleFractionInOil)
244  (*this)[switchIdx] = FsToolbox::value(fluidState.moleFraction(oilPhaseIdx, gasCompIdx));
245  else {
246  assert(switchingVarMeaning() == OilMoleFractionInGas);
247  (*this)[switchIdx] = FsToolbox::value(fluidState.moleFraction(gasPhaseIdx, oilCompIdx));
248  }
249  }
250 
257  bool adaptSwitchingVariable()
258  {
259  Scalar pg = (*this)[Indices::gasPressureIdx];
260 
261  if (switchingVarMeaning() == GasSaturation) {
262  // both hydrocarbon phases
263  Scalar Sw = (*this)[Indices::waterSaturationIdx];
264  Scalar Sg = (*this)[Indices::switchIdx];
265  Scalar So = 1 - Sw - Sg;
266 
267  if (Sg < 0.0 && So > 0.0 && FluidSystem::enableDissolvedGas()) {
268  // we switch to the gas mole fraction in the oil phase if some oil phase
269  // is present and dissolved gas is enabled
270  Scalar xoGsat = FluidSystem::saturatedOilGasMoleFraction(temperature_,
271  pg,
272  pvtRegionIdx_);
273  setSwitchingVarMeaning(GasMoleFractionInOil);
274  (*this)[Indices::switchIdx] = xoGsat;
275  return true;
276  }
277 
278  if (So < 0.0 && Sg > 0.0 && FluidSystem::enableVaporizedOil()) {
279  // we switch to the oil mole fraction in the gas phase if some gas phase
280  // is present and vaporized oil is enabled. TODO (?): use oil instead of
281  // gas pressure!
282  Scalar xgOsat = FluidSystem::saturatedGasOilMoleFraction(temperature_,
283  pg,
284  pvtRegionIdx_);
285  setSwitchingVarMeaning(OilMoleFractionInGas);
286  (*this)[Indices::switchIdx] = xgOsat;
287  return true;
288  }
289  return false;
290  }
291  else if (switchingVarMeaning() == OilMoleFractionInGas) {
292  // only gas. oil appears as soon as more oil is present as can be dissolved
293  Scalar Sw = (*this)[Indices::waterSaturationIdx];
294  Scalar xgO = (*this)[Indices::switchIdx];
295  Scalar xgOsat = FluidSystem::saturatedGasOilMoleFraction(temperature_,
296  pg,
297  pvtRegionIdx_);
298 
299  if (xgO > xgOsat) {
300  setSwitchingVarMeaning(GasSaturation);
301  (*this)[Indices::switchIdx] = 1 - Sw;
302  return true;
303  }
304 
305  return false;
306  }
307  else {
308  assert(switchingVarMeaning() == GasMoleFractionInOil);
309 
310  // only oil. gas appears as soon as more gas is present as can be dissolved
311  Scalar xoG = (*this)[Indices::switchIdx];
312  Scalar xoGsat = FluidSystem::saturatedOilGasMoleFraction(temperature_,
313  pg,
314  pvtRegionIdx_);
315 
316  if (xoG > xoGsat) {
317  setSwitchingVarMeaning(GasSaturation);
318  (*this)[Indices::switchIdx] = 0.0;
319  return true;
320  }
321 
322  return false;
323  }
324 
325  assert(false);
326  return false;
327  }
328 
329  BlackOilPrimaryVariables& operator=(const BlackOilPrimaryVariables& other)
330  {
331  ParentType::operator=(other);
332  temperature_ = other.temperature_;
333  switchingVarMeaning_ = other.switchingVarMeaning_;
334  pvtRegionIdx_ = other.pvtRegionIdx_;
335  return *this;
336  }
337 
338  BlackOilPrimaryVariables& operator=(Scalar value)
339  {
340  for (int i = 0; i < numEq; ++i)
341  (*this)[i] = value;
342 
343  return *this;
344  }
345 
346 private:
347  Scalar temperature_; // so far this is just a pseudo-primary variable
348  SwitchingVarMeaning switchingVarMeaning_;
349  unsigned char pvtRegionIdx_;
350 };
351 
352 
353 } // namespace Ewoms
354 
355 #endif
Represents the primary variables used by the a model.
#define GET_PROP_VALUE(TypeTag, PropTagName)
Access the value attribute of a property for a type tag.
Definition: propertysystem.hh:468
Represents the primary variables used by the black-oil model.
Definition: blackoilprimaryvariables.hh:45
Represents the primary variables used by the a model.
Definition: fvbaseprimaryvariables.hh:41
Declares the properties required by the black oil model.
Definition: baseauxiliarymodule.hh:35
void assignNaive(const FluidState &fluidState)
Assign the primary variables "somehow" from a fluid state.
Definition: fvbaseprimaryvariables.hh:98