FluidStateCompositionModules.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 */
26 #ifndef OPM_FLUID_STATE_COMPOSITION_MODULES_HPP
27 #define OPM_FLUID_STATE_COMPOSITION_MODULES_HPP
28 
31 
32 #include <opm/common/ErrorMacros.hpp>
33 #include <opm/common/Exceptions.hpp>
34 
35 #include <algorithm>
36 #include <cmath>
37 
38 namespace Opm {
39 
44 template <class Scalar,
45  class FluidSystem,
46  class Implementation>
48 {
49  enum { numPhases = FluidSystem::numPhases };
50  enum { numComponents = FluidSystem::numComponents };
51 
52 public:
54  {
58  }
59 
63  const Scalar& moleFraction(unsigned phaseIdx, unsigned compIdx) const
64  { return moleFraction_[phaseIdx][compIdx]; }
65 
69  Scalar massFraction(unsigned phaseIdx, unsigned compIdx) const
70  {
71  typedef Opm::MathToolbox<Scalar> Toolbox;
72 
73  return
75  *moleFraction_[phaseIdx][compIdx]
76  *FluidSystem::molarMass(compIdx)
77  / Toolbox::max(1e-40, Toolbox::abs(averageMolarMass_[phaseIdx]));
78  }
79 
88  const Scalar& averageMolarMass(unsigned phaseIdx) const
89  { return averageMolarMass_[phaseIdx]; }
90 
100  Scalar molarity(unsigned phaseIdx, unsigned compIdx) const
101  { return asImp_().molarDensity(phaseIdx)*moleFraction(phaseIdx, compIdx); }
102 
108  void setMoleFraction(unsigned phaseIdx, unsigned compIdx, const Scalar& value)
109  {
110  Valgrind::CheckDefined(value);
113  Valgrind::SetUndefined(moleFraction_[phaseIdx][compIdx]);
114 
115  moleFraction_[phaseIdx][compIdx] = value;
116 
117  // re-calculate the mean molar mass
118  sumMoleFractions_[phaseIdx] = 0.0;
119  averageMolarMass_[phaseIdx] = 0.0;
120  for (unsigned compJIdx = 0; compJIdx < numComponents; ++compJIdx) {
121  sumMoleFractions_[phaseIdx] += moleFraction_[phaseIdx][compJIdx];
122  averageMolarMass_[phaseIdx] += moleFraction_[phaseIdx][compJIdx]*FluidSystem::molarMass(compJIdx);
123  }
124  }
125 
130  template <class FluidState>
131  void assign(const FluidState& fs)
132  {
133  typedef typename FluidState::Scalar FsScalar;
134  typedef Opm::MathToolbox<FsScalar> FsToolbox;
135 
136  for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
137  averageMolarMass_[phaseIdx] = 0;
138  sumMoleFractions_[phaseIdx] = 0;
139  for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
140  moleFraction_[phaseIdx][compIdx] =
141  FsToolbox::template toLhs<Scalar>(fs.moleFraction(phaseIdx, compIdx));
142 
143  averageMolarMass_[phaseIdx] += moleFraction_[phaseIdx][compIdx]*FluidSystem::molarMass(compIdx);
144  sumMoleFractions_[phaseIdx] += moleFraction_[phaseIdx][compIdx];
145  }
146  }
147  }
148 
157  void checkDefined() const
158  {
162  }
163 
164 protected:
165  const Implementation &asImp_() const
166  { return *static_cast<const Implementation*>(this); }
167 
168  Scalar moleFraction_[numPhases][numComponents];
169  Scalar averageMolarMass_[numPhases];
170  Scalar sumMoleFractions_[numPhases];
171 };
172 
177 template <class Scalar,
178  class FluidSystem,
179  class Implementation>
181 {
182  enum { numPhases = FluidSystem::numPhases };
183 
184 public:
185  enum { numComponents = FluidSystem::numComponents };
186  static_assert((int) numPhases == (int) numComponents,
187  "The number of phases must be the same as the number of (pseudo-) components if you assume immiscibility");
188 
190  { }
191 
195  Scalar moleFraction(unsigned phaseIdx, unsigned compIdx) const
196  { return (phaseIdx == compIdx)?1.0:0.0; }
197 
201  Scalar massFraction(unsigned phaseIdx, unsigned compIdx) const
202  { return (phaseIdx == compIdx)?1.0:0.0; }
203 
212  Scalar averageMolarMass(unsigned phaseIdx) const
213  { return FluidSystem::molarMass(/*compIdx=*/phaseIdx); }
214 
224  Scalar molarity(unsigned phaseIdx, unsigned compIdx) const
225  { return asImp_().molarDensity(phaseIdx)*moleFraction(phaseIdx, compIdx); }
226 
231  template <class FluidState>
232  void assign(const FluidState& /* fs */)
233  { }
234 
243  void checkDefined() const
244  { }
245 
246 protected:
247  const Implementation &asImp_() const
248  { return *static_cast<const Implementation*>(this); }
249 };
250 
255 template <class Scalar>
257 {
258 public:
259  enum { numComponents = 0 };
260 
262  { }
263 
267  Scalar moleFraction(int /* phaseIdx */, int /* compIdx */) const
268  { OPM_THROW(std::logic_error, "Mole fractions are not provided by this fluid state"); }
269 
273  Scalar massFraction(int /* phaseIdx */, int /* compIdx */) const
274  { OPM_THROW(std::logic_error, "Mass fractions are not provided by this fluid state"); }
275 
284  Scalar averageMolarMass(int /* phaseIdx */) const
285  { OPM_THROW(std::logic_error, "Mean molar masses are not provided by this fluid state"); }
286 
296  Scalar molarity(int /* phaseIdx */, int /* compIdx */) const
297  { OPM_THROW(std::logic_error, "Molarities are not provided by this fluid state"); }
298 
307  void checkDefined() const
308  { }
309 };
310 
311 
312 } // namespace Opm
313 
314 #endif
void SetUndefined(const T &value OPM_UNUSED)
Make the memory on which an object resides undefined in valgrind runs.
Definition: Valgrind.hpp:138
Definition: FluidStateCompositionModules.hpp:185
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: MathToolbox.hpp:39
Definition: Air_Mesitylene.hpp:31
Scalar moleFraction(int, int) const
The mole fraction of a component in a phase [].
Definition: FluidStateCompositionModules.hpp:267
Scalar moleFraction_[numPhases][numComponents]
Definition: FluidStateCompositionModules.hpp:168
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
Scalar molarity(unsigned phaseIdx, unsigned compIdx) const
The concentration of a component in a phase [mol/m^3].
Definition: FluidStateCompositionModules.hpp:100
Scalar massFraction(unsigned phaseIdx, unsigned compIdx) const
The mass fraction of a component in a phase [].
Definition: FluidStateCompositionModules.hpp:69
const Scalar & averageMolarMass(unsigned phaseIdx) const
The mean molar mass of a fluid phase [kg/mol].
Definition: FluidStateCompositionModules.hpp:88
Module for the modular fluid state which provides the phase compositions assuming immiscibility...
Definition: FluidStateCompositionModules.hpp:180
Module for the modular fluid state which stores the phase compositions explicitly in terms of mole fr...
Definition: FluidStateCompositionModules.hpp:47
Scalar averageMolarMass(int) const
The mean molar mass of a fluid phase [kg/mol].
Definition: FluidStateCompositionModules.hpp:284
const Implementation & asImp_() const
Definition: FluidStateCompositionModules.hpp:165
Scalar massFraction(int, int) const
The mass fraction of a component in a phase [].
Definition: FluidStateCompositionModules.hpp:273
Scalar averageMolarMass_[numPhases]
Definition: FluidStateCompositionModules.hpp:169
Evaluation< Scalar, VarSetTag, numVars > abs(const Evaluation< Scalar, VarSetTag, numVars > &)
Definition: Math.hpp:41
void checkDefined() const
Make sure that all attributes are defined.
Definition: FluidStateCompositionModules.hpp:157
FluidStateExplicitCompositionModule()
Definition: FluidStateCompositionModules.hpp:53
void checkDefined() const
Make sure that all attributes are defined.
Definition: FluidStateCompositionModules.hpp:307
const Scalar & moleFraction(unsigned phaseIdx, unsigned compIdx) const
The mole fraction of a component in a phase [].
Definition: FluidStateCompositionModules.hpp:63
FluidStateNullCompositionModule()
Definition: FluidStateCompositionModules.hpp:261
Scalar sumMoleFractions_[numPhases]
Definition: FluidStateCompositionModules.hpp:170
void assign(const FluidState &fs)
Retrieve all parameters from an arbitrary fluid state.
Definition: FluidStateCompositionModules.hpp:131
Scalar molarity(int, int) const
The concentration of a component in a phase [mol/m^3].
Definition: FluidStateCompositionModules.hpp:296
void SetDefined(const T &value OPM_UNUSED)
Make the memory on which an object resides defined.
Definition: Valgrind.hpp:188
void setMoleFraction(unsigned phaseIdx, unsigned compIdx, const Scalar &value)
Set the mole fraction of a component in a phase [] and update the average molar mass [kg/mol] accordi...
Definition: FluidStateCompositionModules.hpp:108
Module for the modular fluid state which does not store the compositions but throws std::logic_error ...
Definition: FluidStateCompositionModules.hpp:256
Definition: FluidStateCompositionModules.hpp:259
A traits class which provides basic mathematical functions for arbitrary scalar floating point values...