opm-common
Loading...
Searching...
No Matches
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 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 2 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 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
28#ifndef OPM_FLUID_STATE_COMPOSITION_MODULES_HPP
29#define OPM_FLUID_STATE_COMPOSITION_MODULES_HPP
30
33
34#include <algorithm>
35#include <array>
36#include <cmath>
37
38namespace Opm {
39
44template <class Scalar,
45 class FluidSystem,
46 class Implementation>
47class FluidStateExplicitCompositionModule
48{
49 enum { numPhases = FluidSystem::numPhases };
50 enum { numComponents = FluidSystem::numComponents };
51
52public:
53 FluidStateExplicitCompositionModule()
54 {
55 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
56 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
57 moleFraction_[phaseIdx][compIdx] = 0.0;
58 // TODO: possilby we should begin with totalMoleFractions_ here
59 // the current implementation assuming we have the moleFractions_ for each
60 // individual phases first.
61
62 Valgrind::SetDefined(moleFraction_);
63 Valgrind::SetUndefined(averageMolarMass_);
64 Valgrind::SetUndefined(sumMoleFractions_);
65 }
66
70 const Scalar& moleFraction(unsigned phaseIdx, unsigned compIdx) const
71 { return moleFraction_[phaseIdx][compIdx]; }
72
76 const Scalar& moleFraction(unsigned compIdx) const
77 { return totalModelFractions_[compIdx]; }
78
82 Scalar massFraction(unsigned phaseIdx, unsigned compIdx) const
83 {
84 return
85 abs(sumMoleFractions_[phaseIdx])
86 *moleFraction_[phaseIdx][compIdx]
87 *FluidSystem::molarMass(compIdx)
88 / max(1e-40, abs(averageMolarMass_[phaseIdx]));
89 }
90
99 const Scalar& averageMolarMass(unsigned phaseIdx) const
100 { return averageMolarMass_[phaseIdx]; }
101
111 Scalar molarity(unsigned phaseIdx, unsigned compIdx) const
112 { return asImp_().molarDensity(phaseIdx)*moleFraction(phaseIdx, compIdx); }
113
119 void setMoleFraction(unsigned phaseIdx, unsigned compIdx, const Scalar& value)
120 {
122 Valgrind::SetUndefined(sumMoleFractions_[phaseIdx]);
123 Valgrind::SetUndefined(averageMolarMass_[phaseIdx]);
124 Valgrind::SetUndefined(moleFraction_[phaseIdx][compIdx]);
125
126 moleFraction_[phaseIdx][compIdx] = value;
127
128 // re-calculate the mean molar mass
129 sumMoleFractions_[phaseIdx] = 0.0;
130 averageMolarMass_[phaseIdx] = 0.0;
131 for (unsigned compJIdx = 0; compJIdx < numComponents; ++compJIdx) {
132 sumMoleFractions_[phaseIdx] += moleFraction_[phaseIdx][compJIdx];
133 averageMolarMass_[phaseIdx] += moleFraction_[phaseIdx][compJIdx]*FluidSystem::molarMass(compJIdx);
134 }
135 }
136
140 void setMoleFraction(unsigned compIdx, const Scalar& value)
141 {
143 totalModelFractions_[compIdx] = value;
144 }
145
146 void setCompressFactor(unsigned phaseIdx, const Scalar& value) {
148 Z_[phaseIdx] = value;
149 }
150
151 Scalar compressFactor(unsigned phaseIdx) const {
152 return Z_[phaseIdx];
153 }
154
159 template <class FluidState>
160 void assign(const FluidState& fs)
161 {
162 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
163 averageMolarMass_[phaseIdx] = 0;
164 sumMoleFractions_[phaseIdx] = 0;
165 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
166 moleFraction_[phaseIdx][compIdx] =
167 decay<Scalar>(fs.moleFraction(phaseIdx, compIdx));
168
169 averageMolarMass_[phaseIdx] += moleFraction_[phaseIdx][compIdx]*FluidSystem::molarMass(compIdx);
170 sumMoleFractions_[phaseIdx] += moleFraction_[phaseIdx][compIdx];
171 }
172 }
173 }
174
183 void checkDefined() const
184 {
185 Valgrind::CheckDefined(moleFraction_);
186 Valgrind::CheckDefined(averageMolarMass_);
187 Valgrind::CheckDefined(sumMoleFractions_);
190 }
191
192 const Scalar& K(unsigned compIdx) const
193 {
194 return K_[compIdx];
195 }
196
200 void setKvalue(unsigned compIdx, const Scalar& value)
201 {
202 K_[compIdx] = value;
203 }
204
208 const Scalar& L() const
209 {
210 return L_;
211 }
212
216 void setLvalue(const Scalar& value)
217 {
218 L_ = value;
219 }
220
225 Scalar wilsonK_(unsigned compIdx) const
226 {
227 const auto& acf = FluidSystem::acentricFactor(compIdx);
228 const auto& T_crit = FluidSystem::criticalTemperature(compIdx);
229 const auto& T = asImp_().temperature(0);
230 const auto& p_crit = FluidSystem::criticalPressure(compIdx);
231 const auto& p = asImp_().pressure(0); //for now assume no capillary pressure
232
233 const auto tmp = exp(5.37 * (1+acf) * (1-T_crit/T)) * (p_crit/p);
234 return tmp;
235 }
236
237protected:
238 const Implementation& asImp_() const
239 {
240 return *static_cast<const Implementation*>(this);
241 }
242
243 std::array<std::array<Scalar,numComponents>,numPhases> moleFraction_{};
244 std::array<Scalar,numPhases> averageMolarMass_{};
245 std::array<Scalar,numPhases> sumMoleFractions_{}; // per phase
246 std::array<Scalar, numComponents> totalModelFractions_{}; // total mole fractions for each component
247 std::array<Scalar,numPhases> Z_{};
248 std::array<Scalar,numComponents> K_{};
249 Scalar L_{};
250};
251
256template <class Scalar,
257 class FluidSystem,
258 class Implementation>
259class FluidStateImmiscibleCompositionModule
260{
261 enum { numPhases = FluidSystem::numPhases };
262
263public:
264 enum { numComponents = FluidSystem::numComponents };
265 static_assert(static_cast<int>(numPhases) == static_cast<int>(numComponents),
266 "The number of phases must be the same as the number of (pseudo-) components if you assume immiscibility");
267
268 FluidStateImmiscibleCompositionModule() = default;
269
273 Scalar moleFraction(unsigned phaseIdx, unsigned compIdx) const
274 { return (phaseIdx == compIdx)?1.0:0.0; }
275
279 Scalar massFraction(unsigned phaseIdx, unsigned compIdx) const
280 { return (phaseIdx == compIdx)?1.0:0.0; }
281
290 Scalar averageMolarMass(unsigned phaseIdx) const
291 { return FluidSystem::molarMass(/*compIdx=*/phaseIdx); }
292
302 Scalar molarity(unsigned phaseIdx, unsigned compIdx) const
303 { return asImp_().molarDensity(phaseIdx)*moleFraction(phaseIdx, compIdx); }
304
309 template <class FluidState>
310 void assign(const FluidState& /* fs */)
311 { }
312
321 void checkDefined() const
322 { }
323
324protected:
325 const Implementation& asImp_() const
326 { return *static_cast<const Implementation*>(this); }
327};
328
333template <class Scalar>
334class FluidStateNullCompositionModule
335{
336public:
337 enum { numComponents = 0 };
338
339 FluidStateNullCompositionModule() = default;
340
344 Scalar moleFraction(unsigned /* phaseIdx */, unsigned /* compIdx */) const
345 { throw std::logic_error("Mole fractions are not provided by this fluid state"); }
346
350 Scalar massFraction(unsigned /* phaseIdx */, unsigned /* compIdx */) const
351 { throw std::logic_error("Mass fractions are not provided by this fluid state"); }
352
361 Scalar averageMolarMass(unsigned /* phaseIdx */) const
362 { throw std::logic_error("Mean molar masses are not provided by this fluid state"); }
363
373 Scalar molarity(unsigned /* phaseIdx */, unsigned /* compIdx */) const
374 { throw std::logic_error("Molarities are not provided by this fluid state"); }
375
384 void checkDefined() const
385 { }
386};
387
388
389} // namespace Opm
390
391#endif
A traits class which provides basic mathematical functions for arbitrary scalar floating point values...
Some templates to wrap the valgrind client request macros.
OPM_HOST_DEVICE void SetDefined(const T &value)
Make the memory on which an object resides defined.
Definition Valgrind.hpp:225
OPM_HOST_DEVICE void SetUndefined(const T &value)
Make the memory on which an object resides undefined in valgrind runs.
Definition Valgrind.hpp:174
OPM_HOST_DEVICE bool CheckDefined(const T &value)
Make valgrind complain if any of the memory occupied by an object is undefined.
Definition Valgrind.hpp:76
void setLvalue(const Scalar &value)
Set the L value [-].
Definition FluidStateCompositionModules.hpp:216
Scalar molarity(unsigned phaseIdx, unsigned compIdx) const
The concentration of a component in a phase [mol/m^3].
Definition FluidStateCompositionModules.hpp:111
const Scalar & L() const
The L value of a composition [-].
Definition FluidStateCompositionModules.hpp:208
Scalar wilsonK_(unsigned compIdx) const
Wilson formula to calculate K.
Definition FluidStateCompositionModules.hpp:225
void assign(const FluidState &fs)
Retrieve all parameters from an arbitrary fluid state.
Definition FluidStateCompositionModules.hpp:160
const Scalar & averageMolarMass(unsigned phaseIdx) const
The mean molar mass of a fluid phase [kg/mol].
Definition FluidStateCompositionModules.hpp:99
const Scalar & moleFraction(unsigned phaseIdx, unsigned compIdx) const
The mole fraction of a component in a phase [].
Definition FluidStateCompositionModules.hpp:70
const Scalar & moleFraction(unsigned compIdx) const
The total mole fraction of a component [].
Definition FluidStateCompositionModules.hpp:76
void setMoleFraction(unsigned compIdx, const Scalar &value)
Set the total mole fraction of a component.
Definition FluidStateCompositionModules.hpp:140
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:119
void setKvalue(unsigned compIdx, const Scalar &value)
Set the K value of a component [-].
Definition FluidStateCompositionModules.hpp:200
Scalar massFraction(unsigned phaseIdx, unsigned compIdx) const
The mass fraction of a component in a phase [].
Definition FluidStateCompositionModules.hpp:82
void checkDefined() const
Make sure that all attributes are defined.
Definition FluidStateCompositionModules.hpp:183
Scalar massFraction(unsigned phaseIdx, unsigned compIdx) const
The mass fraction of a component in a phase [].
Definition FluidStateCompositionModules.hpp:279
Scalar molarity(unsigned phaseIdx, unsigned compIdx) const
The concentration of a component in a phase [mol/m^3].
Definition FluidStateCompositionModules.hpp:302
void checkDefined() const
Make sure that all attributes are defined.
Definition FluidStateCompositionModules.hpp:321
void assign(const FluidState &)
Retrieve all parameters from an arbitrary fluid state.
Definition FluidStateCompositionModules.hpp:310
Scalar moleFraction(unsigned phaseIdx, unsigned compIdx) const
The mole fraction of a component in a phase [].
Definition FluidStateCompositionModules.hpp:273
Scalar averageMolarMass(unsigned phaseIdx) const
The mean molar mass of a fluid phase [kg/mol].
Definition FluidStateCompositionModules.hpp:290
Scalar moleFraction(unsigned, unsigned) const
The mole fraction of a component in a phase [].
Definition FluidStateCompositionModules.hpp:344
Scalar averageMolarMass(unsigned) const
The mean molar mass of a fluid phase [kg/mol].
Definition FluidStateCompositionModules.hpp:361
Scalar molarity(unsigned, unsigned) const
The concentration of a component in a phase [mol/m^3].
Definition FluidStateCompositionModules.hpp:373
void checkDefined() const
Make sure that all attributes are defined.
Definition FluidStateCompositionModules.hpp:384
Scalar massFraction(unsigned, unsigned) const
The mass fraction of a component in a phase [].
Definition FluidStateCompositionModules.hpp:350
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30