MiscibleMultiPhaseComposition.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  Copyright (C) 2011-2012 by Klaus Mosthaf
6  Copyright (C) 2012 by Bernd Flemisch
7 
8  This file is part of the Open Porous Media project (OPM).
9 
10  OPM is free software: you can redistribute it and/or modify
11  it under the terms of the GNU General Public License as published by
12  the Free Software Foundation, either version 2 of the License, or
13  (at your option) any later version.
14 
15  OPM is distributed in the hope that it will be useful,
16  but WITHOUT ANY WARRANTY; without even the implied warranty of
17  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18  GNU General Public License for more details.
19 
20  You should have received a copy of the GNU General Public License
21  along with OPM. If not, see <http://www.gnu.org/licenses/>.
22 */
27 #ifndef OPM_MISCIBLE_MULTIPHASE_COMPOSITION_HPP
28 #define OPM_MISCIBLE_MULTIPHASE_COMPOSITION_HPP
29 
31 
32 #include <dune/common/fvector.hh>
33 #include <dune/common/fmatrix.hh>
34 
35 #include <opm/common/Exceptions.hpp>
36 #include <opm/common/ErrorMacros.hpp>
38 
39 namespace Opm {
40 
48 template <class Scalar>
50 {
51 public:
53  {}
54 
55  MMPCAuxConstraint(unsigned phaseIdx, unsigned compIdx, Scalar value)
56  : phaseIdx_(phaseIdx)
57  , compIdx_(compIdx)
58  , value_(value)
59  {}
60 
64  void set(unsigned phaseIdx, unsigned compIdx, Scalar value)
65  {
66  phaseIdx_ = phaseIdx;
67  compIdx_ = compIdx;
68  value_ = value;
69  }
70 
75  unsigned phaseIdx() const
76  { return phaseIdx_; }
77 
82  unsigned compIdx() const
83  { return compIdx_; }
84 
89  Scalar value() const
90  { return value_; }
91 
92 private:
93  unsigned phaseIdx_;
94  unsigned compIdx_;
95  Scalar value_;
96 };
97 
121 template <class Scalar, class FluidSystem, class Evaluation = Scalar>
123 {
124  static const int numPhases = FluidSystem::numPhases;
125  static const int numComponents = FluidSystem::numComponents;
126 
128 
129  static_assert(numPhases <= numComponents,
130  "This solver requires that the number fluid phases is smaller or equal "
131  "to the number of components");
132 
133 
134 public:
158  template <class FluidState, class ParameterCache>
159  static void solve(FluidState &fluidState,
160  ParameterCache &paramCache,
161  int phasePresence,
162  const MMPCAuxConstraint<Evaluation> *auxConstraints,
163  unsigned numAuxConstraints,
164  bool setViscosity,
165  bool setInternalEnergy)
166  {
168  static_assert(std::is_same<typename FluidState::Scalar, Evaluation>::value,
169  "The scalar type of the fluid state must be 'Evaluation'");
170 
171 #ifndef NDEBUG
172  // currently this solver can only handle fluid systems which
173  // assume ideal mixtures of all fluids. TODO: relax this
174  // (requires solving a non-linear system of equations, i.e. using
175  // newton method.)
176  for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
177  assert(FluidSystem::isIdealMixture(phaseIdx));
178  }
179 #endif
180 
181  // compute all fugacity coefficients
182  for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
183  paramCache.updatePhase(fluidState, phaseIdx);
184 
185  // since we assume ideal mixtures, the fugacity
186  // coefficients of the components cannot depend on
187  // composition, i.e. the parameters in the cache are valid
188  for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
189  Evaluation fugCoeff = FsToolbox::template toLhs<Evaluation>(
190  FluidSystem::fugacityCoefficient(fluidState, paramCache, phaseIdx, compIdx));
191  fluidState.setFugacityCoefficient(phaseIdx, compIdx, fugCoeff);
192  }
193  }
194 
195  // create the linear system of equations which defines the
196  // mole fractions
197  static const int numEq = numComponents*numPhases;
198  Dune::FieldMatrix<Evaluation, numEq, numEq> M(Toolbox::createConstant(0.0));
199  Dune::FieldVector<Evaluation, numEq> x(Toolbox::createConstant(0.0));
200  Dune::FieldVector<Evaluation, numEq> b(Toolbox::createConstant(0.0));
201 
202  // assemble the equations expressing the fact that the
203  // fugacities of each component are equal in all phases
204  for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
205  const Evaluation& entryCol1 =
206  fluidState.fugacityCoefficient(/*phaseIdx=*/0, compIdx)
207  *fluidState.pressure(/*phaseIdx=*/0);
208  unsigned col1Idx = compIdx;
209 
210  for (unsigned phaseIdx = 1; phaseIdx < numPhases; ++phaseIdx) {
211  unsigned rowIdx = (phaseIdx - 1)*numComponents + compIdx;
212  unsigned col2Idx = phaseIdx*numComponents + compIdx;
213 
214  const Evaluation& entryCol2 =
215  fluidState.fugacityCoefficient(phaseIdx, compIdx)
216  *fluidState.pressure(phaseIdx);
217 
218  M[rowIdx][col1Idx] = entryCol1;
219  M[rowIdx][col2Idx] = -entryCol2;
220  }
221  }
222 
223  // assemble the equations expressing the assumption that the
224  // sum of all mole fractions in each phase must be 1 for the
225  // phases present.
226  unsigned presentPhases = 0;
227  for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
228  if (!(phasePresence & (1 << phaseIdx)))
229  continue;
230 
231  unsigned rowIdx = numComponents*(numPhases - 1) + presentPhases;
232  presentPhases += 1;
233 
234  b[rowIdx] = Toolbox::createConstant(1.0);
235  for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
236  unsigned colIdx = phaseIdx*numComponents + compIdx;
237 
238  M[rowIdx][colIdx] = Toolbox::createConstant(1.0);
239  }
240  }
241 
242  assert(presentPhases + numAuxConstraints == numComponents);
243 
244  // incorperate the auxiliary equations, i.e., the explicitly given mole fractions
245  for (unsigned auxEqIdx = 0; auxEqIdx < numAuxConstraints; ++auxEqIdx) {
246  unsigned rowIdx = numComponents*(numPhases - 1) + presentPhases + auxEqIdx;
247  b[rowIdx] = auxConstraints[auxEqIdx].value();
248 
249  unsigned colIdx = auxConstraints[auxEqIdx].phaseIdx()*numComponents + auxConstraints[auxEqIdx].compIdx();
250  M[rowIdx][colIdx] = 1.0;
251  }
252 
253  // solve for all mole fractions
254  try {
255  Dune::FMatrixPrecision<Scalar>::set_singular_limit(1e-50);
256  M.solve(x, b);
257  }
258  catch (const Dune::FMatrixError &e) {
259  OPM_THROW(NumericalProblem,
260  "Numerical problem in MiscibleMultiPhaseComposition::solve(): " << e.what() << "; M="<<M);
261  }
262  catch (...) {
263  throw;
264  }
265 
266 
267  // set all mole fractions and the additional quantities in
268  // the fluid state
269  for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
270  for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
271  unsigned rowIdx = phaseIdx*numComponents + compIdx;
272  fluidState.setMoleFraction(phaseIdx, compIdx, x[rowIdx]);
273  }
274  paramCache.updateComposition(fluidState, phaseIdx);
275 
276  const Evaluation& rho = FluidSystem::density(fluidState, paramCache, phaseIdx);
277  fluidState.setDensity(phaseIdx, rho);
278 
279  if (setViscosity) {
280  const Evaluation& mu = FluidSystem::viscosity(fluidState, paramCache, phaseIdx);
281  fluidState.setViscosity(phaseIdx, mu);
282  }
283 
284  if (setInternalEnergy) {
285  const Evaluation& h = FluidSystem::enthalpy(fluidState, paramCache, phaseIdx);
286  fluidState.setEnthalpy(phaseIdx, h);
287  }
288  }
289  }
290 
298  template <class FluidState, class ParameterCache>
299  static void solve(FluidState &fluidState,
300  ParameterCache &paramCache,
301  bool setViscosity,
302  bool setInternalEnergy)
303  {
304  solve(fluidState,
305  paramCache,
306  /*phasePresence=*/0xffffff,
307  /*numAuxConstraints=*/0,
308  /*auxConstraints=*/0,
309  setViscosity,
310  setInternalEnergy);
311  }
312 };
313 
314 } // namespace Opm
315 
316 #endif
Definition: MathToolbox.hpp:39
Definition: Air_Mesitylene.hpp:31
Some templates to wrap the valgrind client request macros.
void set(unsigned phaseIdx, unsigned compIdx, Scalar value)
Specify the auxiliary constraint.
Definition: MiscibleMultiPhaseComposition.hpp:64
MMPCAuxConstraint()
Definition: MiscibleMultiPhaseComposition.hpp:52
Specifies an auxiliary constraint for the MiscibleMultiPhaseComposition constraint solver...
Definition: MiscibleMultiPhaseComposition.hpp:49
unsigned phaseIdx() const
Returns the index of the fluid phase for which the auxiliary constraint is specified.
Definition: MiscibleMultiPhaseComposition.hpp:75
Computes the composition of all phases of a N-phase, N-component fluid system assuming that all N pha...
Definition: MiscibleMultiPhaseComposition.hpp:122
unsigned compIdx() const
Returns the index of the component for which the auxiliary constraint is specified.
Definition: MiscibleMultiPhaseComposition.hpp:82
Scalar value() const
Returns value of the mole fraction of the auxiliary constraint.
Definition: MiscibleMultiPhaseComposition.hpp:89
MMPCAuxConstraint(unsigned phaseIdx, unsigned compIdx, Scalar value)
Definition: MiscibleMultiPhaseComposition.hpp:55
A traits class which provides basic mathematical functions for arbitrary scalar floating point values...