ImmiscibleFlash.hpp
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-2013 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 */
25 #ifndef OPM_IMMISCIBLE_FLASH_HPP
26 #define OPM_IMMISCIBLE_FLASH_HPP
27 
28 #include <opm/common/ErrorMacros.hpp>
29 #include <opm/common/Exceptions.hpp>
31 
32 #include <dune/common/fvector.hh>
33 #include <dune/common/fmatrix.hh>
34 
35 #include <limits>
36 #include <iostream>
37 
38 namespace Opm {
39 
66 template <class Scalar, class FluidSystem>
68 {
69  static const int numPhases = FluidSystem::numPhases;
70  static const int numComponents = FluidSystem::numComponents;
71  static_assert(numPhases == numComponents,
72  "Immiscibility assumes that the number of phases"
73  " is equal to the number of components");
74 
75  typedef typename FluidSystem::ParameterCache ParameterCache;
76 
77  static const int numEq = numPhases;
78 
79  typedef Dune::FieldMatrix<Scalar, numEq, numEq> Matrix;
80  typedef Dune::FieldVector<Scalar, numEq> Vector;
81 
82 public:
83  typedef Dune::FieldVector<Scalar, numComponents> ComponentVector;
84 
88  template <class FluidState>
89  static void guessInitial(FluidState &fluidState,
90  ParameterCache &/*paramCache*/,
91  const ComponentVector &globalMolarities)
92  {
93  // the sum of all molarities
94  Scalar sumMoles = 0;
95  for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx)
96  sumMoles += globalMolarities[compIdx];
97 
98  for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
99  // pressure. use atmospheric pressure as initial guess
100  fluidState.setPressure(phaseIdx, 2e5);
101 
102  // saturation. assume all fluids to be equally distributed
103  fluidState.setSaturation(phaseIdx, 1.0/numPhases);
104  }
105  }
106 
113  template <class MaterialLaw, class FluidState>
114  static void solve(FluidState &fluidState,
115  ParameterCache &paramCache,
116  const typename MaterialLaw::Params &matParams,
117  const ComponentVector &globalMolarities)
118  {
119  Dune::FMatrixPrecision<Scalar>::set_singular_limit(1e-25);
120 
122  // Check if all fluid phases are incompressible
124  bool allIncompressible = true;
125  for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
126  if (FluidSystem::isCompressible(phaseIdx)) {
127  allIncompressible = false;
128  break;
129  }
130  };
131 
132  if (allIncompressible) {
133  paramCache.updateAll(fluidState);
134  for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
135  Scalar rho = FluidSystem::density(fluidState, paramCache, phaseIdx);
136  fluidState.setDensity(phaseIdx, rho);
137 
138  Scalar saturation =
139  globalMolarities[/*compIdx=*/phaseIdx]
140  / fluidState.molarDensity(phaseIdx);
141  fluidState.setSaturation(phaseIdx, saturation);
142  }
143  };
144 
146  // Newton method
148 
149  // Jacobian matrix
150  Matrix J;
151  // solution, i.e. phase composition
152  Vector deltaX;
153  // right hand side
154  Vector b;
155 
157  Valgrind::SetUndefined(deltaX);
159 
160  completeFluidState_<MaterialLaw>(fluidState, paramCache, matParams);
161 
162  const int nMax = 50; // <- maximum number of newton iterations
163  for (int nIdx = 0; nIdx < nMax; ++nIdx) {
164  // calculate Jacobian matrix and right hand side
165  linearize_<MaterialLaw>(J, b, fluidState, paramCache, matParams, globalMolarities);
168 
169  // Solve J*x = b
170  deltaX = 0;
171 
172  try { J.solve(deltaX, b); }
173  catch (Dune::FMatrixError e)
174  {
175  /*
176  std::cout << "error: " << e << "\n";
177  std::cout << "b: " << b << "\n";
178  */
179 
180  throw Opm::NumericalProblem(e.what());
181  }
182  Valgrind::CheckDefined(deltaX);
183 
184  /*
185  printFluidState_(fluidState);
186  std::cout << "J:\n";
187  for (int i = 0; i < numEq; ++i) {
188  for (int j = 0; j < numEq; ++j) {
189  std::ostringstream os;
190  os << J[i][j];
191 
192  std::string s(os.str());
193  do {
194  s += " ";
195  } while (s.size() < 20);
196  std::cout << s;
197  }
198  std::cout << "\n";
199  };
200  std::cout << "residual: " << b << "\n";
201  std::cout << "deltaX: " << deltaX << "\n";
202  std::cout << "---------------\n";
203  */
204 
205  // update the fluid quantities.
206  Scalar relError = update_<MaterialLaw>(fluidState, paramCache, matParams, deltaX);
207 
208  if (relError < 1e-9)
209  return;
210  }
211 
212  /*
213  printFluidState_(fluidState);
214  std::cout << "globalMolarities: ";
215  for (unsigned compIdx = 0; compIdx < numComponents; ++ compIdx)
216  std::cout << globalMolarities[compIdx] << " ";
217  std::cout << "\n";
218  */
219 
220  OPM_THROW(Opm::NumericalProblem,
221  "Flash calculation failed."
222  " {c_alpha^kappa} = {" << globalMolarities << "}, T = "
223  << fluidState.temperature(/*phaseIdx=*/0));
224  }
225 
226 
227 protected:
228  template <class FluidState>
229  static void printFluidState_(const FluidState &fs)
230  {
231  std::cout << "saturations: ";
232  for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
233  std::cout << fs.saturation(phaseIdx) << " ";
234  std::cout << "\n";
235 
236  std::cout << "pressures: ";
237  for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
238  std::cout << fs.pressure(phaseIdx) << " ";
239  std::cout << "\n";
240 
241  std::cout << "densities: ";
242  for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
243  std::cout << fs.density(phaseIdx) << " ";
244  std::cout << "\n";
245 
246  std::cout << "global component molarities: ";
247  for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
248  Scalar sum = 0;
249  for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
250  sum += fs.saturation(phaseIdx)*fs.molarity(phaseIdx, compIdx);
251  }
252  std::cout << sum << " ";
253  }
254  std::cout << "\n";
255  }
256 
257  template <class MaterialLaw, class FluidState>
258  static void linearize_(Matrix &J,
259  Vector &b,
260  FluidState &fluidState,
261  ParameterCache &paramCache,
262  const typename MaterialLaw::Params &matParams,
263  const ComponentVector &globalMolarities)
264  {
265  FluidState origFluidState(fluidState);
266  ParameterCache origParamCache(paramCache);
267 
268  Vector tmp;
269 
270  // reset jacobian
271  J = 0;
272 
274  calculateDefect_(b, fluidState, fluidState, globalMolarities);
276 
277  // assemble jacobian matrix
278  for (unsigned pvIdx = 0; pvIdx < numEq; ++ pvIdx) {
280  // approximately calculate partial derivatives of the
281  // fugacity defect of all components in regard to the mole
282  // fraction of the i-th component. This is done via
283  // forward differences
284 
285  // deviate the mole fraction of the i-th component
286  Scalar xI = getQuantity_(fluidState, pvIdx);
287  const Scalar eps = 1e-10/quantityWeight_(fluidState, pvIdx);
288  setQuantity_<MaterialLaw>(fluidState, paramCache, matParams, pvIdx, xI + eps);
289  assert(std::abs(getQuantity_(fluidState, pvIdx) - (xI + eps))
290  <= std::max(1.0, std::abs(xI))*std::numeric_limits<Scalar>::epsilon()*100);
291 
292  // compute derivative of the defect
293  calculateDefect_(tmp, origFluidState, fluidState, globalMolarities);
294  tmp -= b;
295  tmp /= eps;
296  // store derivative in jacobian matrix
297  for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx)
298  J[eqIdx][pvIdx] = tmp[eqIdx];
299 
300  // fluid state and parameter cache to their original values
301  fluidState = origFluidState;
302  paramCache = origParamCache;
303 
304  // end forward differences
306  }
307  }
308 
309  template <class FluidState>
310  static void calculateDefect_(Vector &b,
311  const FluidState &/*fluidStateEval*/,
312  const FluidState &fluidState,
313  const ComponentVector &globalMolarities)
314  {
315  // global molarities are given
316  for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
317  b[phaseIdx] =
318  fluidState.saturation(phaseIdx)
319  * fluidState.molarity(phaseIdx, /*compIdx=*/phaseIdx);
320 
321  b[phaseIdx] -= globalMolarities[/*compIdx=*/phaseIdx];
322  }
323  }
324 
325  template <class MaterialLaw, class FluidState>
326  static Scalar update_(FluidState &fluidState,
327  ParameterCache &paramCache,
328  const typename MaterialLaw::Params &matParams,
329  const Vector &deltaX)
330  {
331  Scalar relError = 0;
332  for (unsigned pvIdx = 0; pvIdx < numEq; ++ pvIdx) {
333  Scalar tmp = getQuantity_(fluidState, pvIdx);
334  Scalar delta = deltaX[pvIdx];
335 
336  relError = std::max(relError, std::abs(delta)*quantityWeight_(fluidState, pvIdx));
337 
338  if (isSaturationIdx_(pvIdx)) {
339  // dampen to at most 20% change in saturation per
340  // iteration
341  delta = std::min(0.2, std::max(-0.2, delta));
342  }
343  else if (isPressureIdx_(pvIdx)) {
344  // dampen to at most 30% change in pressure per
345  // iteration
346  delta = std::min(0.30*fluidState.pressure(0), std::max(-0.30*fluidState.pressure(0), delta));
347  };
348 
349  setQuantityRaw_(fluidState, pvIdx, tmp - delta);
350  }
351 
352  completeFluidState_<MaterialLaw>(fluidState, paramCache, matParams);
353 
354  return relError;
355  }
356 
357  template <class MaterialLaw, class FluidState>
358  static void completeFluidState_(FluidState &fluidState,
359  ParameterCache &paramCache,
360  const typename MaterialLaw::Params &matParams)
361  {
362  // calculate the saturation of the last phase as a function of
363  // the other saturations
364  Scalar sumSat = 0.0;
365  for (unsigned phaseIdx = 0; phaseIdx < numPhases - 1; ++phaseIdx)
366  sumSat += fluidState.saturation(phaseIdx);
367 
368  if (sumSat > 1.0) {
369  // make sure that the last saturation does not become
370  // negative
371  for (unsigned phaseIdx = 0; phaseIdx < numPhases - 1; ++phaseIdx)
372  {
373  Scalar S = fluidState.saturation(phaseIdx);
374  fluidState.setSaturation(phaseIdx, S/sumSat);
375  }
376  sumSat = 1;
377  }
378 
379  // set the last saturation
380  fluidState.setSaturation(/*phaseIdx=*/numPhases - 1, 1.0 - sumSat);
381 
382  // update the pressures using the material law (saturations
383  // and first pressure are already set because it is implicitly
384  // solved for.)
385  ComponentVector pC;
386  MaterialLaw::capillaryPressures(pC, matParams, fluidState);
387  for (unsigned phaseIdx = 1; phaseIdx < numPhases; ++phaseIdx)
388  fluidState.setPressure(phaseIdx,
389  fluidState.pressure(0)
390  + (pC[phaseIdx] - pC[0]));
391 
392  // update the parameter cache
393  paramCache.updateAll(fluidState, /*except=*/ParameterCache::Temperature|ParameterCache::Composition);
394 
395  // update all densities
396  for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
397  Scalar rho = FluidSystem::density(fluidState, paramCache, phaseIdx);
398  fluidState.setDensity(phaseIdx, rho);
399  }
400  }
401 
402  static bool isPressureIdx_(unsigned pvIdx)
403  { return pvIdx == 0; }
404 
405  static bool isSaturationIdx_(unsigned pvIdx)
406  { return 1 <= pvIdx && pvIdx < numPhases; }
407 
408  // retrieves a quantity from the fluid state
409  template <class FluidState>
410  static Scalar getQuantity_(const FluidState &fs, unsigned pvIdx)
411  {
412  assert(pvIdx < numEq);
413 
414  // first pressure
415  if (pvIdx < 1) {
416  unsigned phaseIdx = 0;
417  return fs.pressure(phaseIdx);
418  }
419  // saturations
420  else {
421  assert(pvIdx < numPhases);
422  unsigned phaseIdx = pvIdx - 1;
423  return fs.saturation(phaseIdx);
424  }
425  }
426 
427  // set a quantity in the fluid state
428  template <class MaterialLaw, class FluidState>
429  static void setQuantity_(FluidState &fs,
430  ParameterCache &paramCache,
431  const typename MaterialLaw::Params &matParams,
432  unsigned pvIdx,
433  Scalar value)
434  {
435  assert(0 <= pvIdx && pvIdx < numEq);
436 
437  if (pvIdx < 1) {
438  // -> first pressure
439  Scalar delta = value - fs.pressure(0);
440 
441  // set all pressures. here we assume that the capillary
442  // pressure does not depend on absolute pressure.
443  for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
444  fs.setPressure(phaseIdx, fs.pressure(phaseIdx) + delta);
445  paramCache.updateAllPressures(fs);
446 
447  // update all densities
448  for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
449  Scalar rho = FluidSystem::density(fs, paramCache, phaseIdx);
450  fs.setDensity(phaseIdx, rho);
451  }
452  }
453  else {
454  assert(pvIdx < numPhases);
455 
456  // -> first M-1 saturations
457  Scalar delta = value - fs.saturation(/*phaseIdx=*/pvIdx - 1);
458  fs.setSaturation(/*phaseIdx=*/pvIdx - 1, value);
459 
460  // set last saturation (-> minus the change of the
461  // saturation of the other phases)
462  fs.setSaturation(/*phaseIdx=*/numPhases - 1,
463  fs.saturation(numPhases - 1) - delta);
464 
465  // update all fluid pressures using the capillary pressure
466  // law
467  ComponentVector pC;
468  MaterialLaw::capillaryPressures(pC, matParams, fs);
469  for (unsigned phaseIdx = 1; phaseIdx < numPhases; ++phaseIdx)
470  fs.setPressure(phaseIdx,
471  fs.pressure(0)
472  + (pC[phaseIdx] - pC[0]));
473  paramCache.updateAllPressures(fs);
474 
475  // update all densities
476  for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
477  Scalar rho = FluidSystem::density(fs, paramCache, phaseIdx);
478  fs.setDensity(phaseIdx, rho);
479  }
480  }
481  }
482 
483  // set a quantity in the fluid state
484  template <class FluidState>
485  static void setQuantityRaw_(FluidState &fs, unsigned pvIdx, Scalar value)
486  {
487  assert(pvIdx < numEq);
488 
489  // first pressure
490  if (pvIdx < 1) {
491  unsigned phaseIdx = 0;
492  fs.setPressure(phaseIdx, value);
493  }
494  // saturations
495  else {
496  assert(pvIdx < numPhases);
497  unsigned phaseIdx = pvIdx - 1;
498 
499  // make sure that the first M-1 saturations does not get
500  // negative
501  value = std::max(0.0, value);
502  fs.setSaturation(phaseIdx, value);
503  }
504  }
505 
506  template <class FluidState>
507  static Scalar quantityWeight_(const FluidState &/*fs*/, unsigned pvIdx)
508  {
509  // first pressure
510  if (pvIdx < 1)
511  return 1.0/1e5;
512  // first M - 1 saturations
513  else {
514  assert(pvIdx < numPhases);
515  return 1.0;
516  }
517  }
518 };
519 
520 } // namespace Opm
521 
522 #endif
void SetUndefined(const T &value OPM_UNUSED)
Make the memory on which an object resides undefined in valgrind runs.
Definition: Valgrind.hpp:138
bool CheckDefined(const T &value OPM_UNUSED)
Make valgrind complain if any of the memory occupied by an object is undefined.
Definition: Valgrind.hpp:74
Definition: Air_Mesitylene.hpp:31
Some templates to wrap the valgrind client request macros.
Evaluation< Scalar, VarSetTag, numVars > max(const Evaluation< Scalar, VarSetTag, numVars > &x1, const Evaluation< Scalar, VarSetTag, numVars > &x2)
Definition: Math.hpp:114
Determines the pressures and saturations of all fluid phases given the total mass of all components...
Definition: ImmiscibleFlash.hpp:67
Evaluation< Scalar, VarSetTag, numVars > abs(const Evaluation< Scalar, VarSetTag, numVars > &)
Definition: Math.hpp:41
Evaluation< Scalar, VarSetTag, numVars > min(const Evaluation< Scalar, VarSetTag, numVars > &x1, const Evaluation< Scalar, VarSetTag, numVars > &x2)
Definition: Math.hpp:61