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 <array>
35#include <cmath>
36
37namespace Opm {
38
43template <class ValueType,
44 class FluidSystem,
45 class Implementation>
46class FluidStateExplicitCompositionModule
47{
48 enum { numPhases = FluidSystem::numPhases };
49 enum { numComponents = FluidSystem::numComponents };
50
51public:
52 FluidStateExplicitCompositionModule()
53 {
54 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
55 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
56 moleFraction_[phaseIdx][compIdx] = 0.0;
57 // TODO: possilby we should begin with totalMoleFractions_ here
58 // the current implementation assuming we have the moleFractions_ for each
59 // individual phases first.
60
61 Valgrind::SetDefined(moleFraction_);
62 Valgrind::SetUndefined(averageMolarMass_);
63 Valgrind::SetUndefined(sumMoleFractions_);
64 }
65
69 const ValueType& moleFraction(unsigned phaseIdx, unsigned compIdx) const
70 { return moleFraction_[phaseIdx][compIdx]; }
71
75 const ValueType& moleFraction(unsigned compIdx) const
76 { return totalModelFractions_[compIdx]; }
77
81 ValueType massFraction(unsigned phaseIdx, unsigned compIdx) const
82 {
83 return
84 abs(sumMoleFractions_[phaseIdx])
85 *moleFraction_[phaseIdx][compIdx]
86 *FluidSystem::molarMass(compIdx)
87 / max(1e-40, abs(averageMolarMass_[phaseIdx]));
88 }
89
98 const ValueType& averageMolarMass(unsigned phaseIdx) const
99 { return averageMolarMass_[phaseIdx]; }
100
110 ValueType molarity(unsigned phaseIdx, unsigned compIdx) const
111 { return asImp_().molarDensity(phaseIdx)*moleFraction(phaseIdx, compIdx); }
112
118 void setMoleFraction(unsigned phaseIdx, unsigned compIdx, const ValueType& value)
119 {
121 Valgrind::SetUndefined(sumMoleFractions_[phaseIdx]);
122 Valgrind::SetUndefined(averageMolarMass_[phaseIdx]);
123 Valgrind::SetUndefined(moleFraction_[phaseIdx][compIdx]);
124
125 moleFraction_[phaseIdx][compIdx] = value;
126
127 // re-calculate the mean molar mass
128 sumMoleFractions_[phaseIdx] = 0.0;
129 averageMolarMass_[phaseIdx] = 0.0;
130 for (unsigned compJIdx = 0; compJIdx < numComponents; ++compJIdx) {
131 sumMoleFractions_[phaseIdx] += moleFraction_[phaseIdx][compJIdx];
132 averageMolarMass_[phaseIdx] += moleFraction_[phaseIdx][compJIdx]*FluidSystem::molarMass(compJIdx);
133 }
134 }
135
139 void setMoleFraction(unsigned compIdx, const ValueType& value)
140 {
142 totalModelFractions_[compIdx] = value;
143 }
144
145 void setCompressFactor(unsigned phaseIdx, const ValueType& value) {
147 Z_[phaseIdx] = value;
148 }
149
150 ValueType compressFactor(unsigned phaseIdx) const {
151 return Z_[phaseIdx];
152 }
153
158 template <class FluidState>
159 void assign(const FluidState& fs)
160 {
161 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
162 averageMolarMass_[phaseIdx] = 0;
163 sumMoleFractions_[phaseIdx] = 0;
164 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
165 moleFraction_[phaseIdx][compIdx] =
166 decay<ValueType>(fs.moleFraction(phaseIdx, compIdx));
167
168 averageMolarMass_[phaseIdx] += moleFraction_[phaseIdx][compIdx]*FluidSystem::molarMass(compIdx);
169 sumMoleFractions_[phaseIdx] += moleFraction_[phaseIdx][compIdx];
170 }
171 }
172 }
173
182 void checkDefined() const
183 {
184 Valgrind::CheckDefined(moleFraction_);
185 Valgrind::CheckDefined(averageMolarMass_);
186 Valgrind::CheckDefined(sumMoleFractions_);
189 }
190
191 const ValueType& K(unsigned compIdx) const
192 {
193 return K_[compIdx];
194 }
195
199 void setKvalue(unsigned compIdx, const ValueType& value)
200 {
201 K_[compIdx] = value;
202 }
203
207 const ValueType& L() const
208 {
209 return L_;
210 }
211
215 void setLvalue(const ValueType& value)
216 {
217 L_ = value;
218 }
219
224 ValueType wilsonK_(unsigned compIdx) const
225 {
226 const auto& acf = FluidSystem::acentricFactor(compIdx);
227 const auto& T_crit = FluidSystem::criticalTemperature(compIdx);
228 const auto& T = asImp_().temperature(0);
229 const auto& p_crit = FluidSystem::criticalPressure(compIdx);
230 const auto& p = asImp_().pressure(0); //for now assume no capillary pressure
231
232 const auto tmp = exp(5.37 * (1+acf) * (1-T_crit/T)) * (p_crit/p);
233 return tmp;
234 }
235
236protected:
237 const Implementation& asImp_() const
238 {
239 return *static_cast<const Implementation*>(this);
240 }
241
242 std::array<std::array<ValueType,numComponents>,numPhases> moleFraction_{};
243 std::array<ValueType,numPhases> averageMolarMass_{};
244 std::array<ValueType,numPhases> sumMoleFractions_{}; // per phase
245 std::array<ValueType, numComponents> totalModelFractions_{}; // total mole fractions for each component
246 std::array<ValueType,numPhases> Z_{};
247 std::array<ValueType,numComponents> K_{};
248 ValueType L_{};
249};
250
255template <class ValueType,
256 class FluidSystem,
257 class Implementation>
258class FluidStateImmiscibleCompositionModule
259{
260 enum { numPhases = FluidSystem::numPhases };
261
262public:
263 enum { numComponents = FluidSystem::numComponents };
264 static_assert(static_cast<int>(numPhases) == static_cast<int>(numComponents),
265 "The number of phases must be the same as the number of (pseudo-) components if you assume immiscibility");
266
267 FluidStateImmiscibleCompositionModule() = default;
268
272 ValueType moleFraction(unsigned phaseIdx, unsigned compIdx) const
273 { return (phaseIdx == compIdx)?1.0:0.0; }
274
278 ValueType massFraction(unsigned phaseIdx, unsigned compIdx) const
279 { return (phaseIdx == compIdx)?1.0:0.0; }
280
289 ValueType averageMolarMass(unsigned phaseIdx) const
290 { return FluidSystem::molarMass(/*compIdx=*/phaseIdx); }
291
301 ValueType molarity(unsigned phaseIdx, unsigned compIdx) const
302 { return asImp_().molarDensity(phaseIdx)*moleFraction(phaseIdx, compIdx); }
303
308 template <class FluidState>
309 void assign(const FluidState& /* fs */)
310 { }
311
320 void checkDefined() const
321 { }
322
323protected:
324 const Implementation& asImp_() const
325 { return *static_cast<const Implementation*>(this); }
326};
327
332template <class ValueT>
333class FluidStateNullCompositionModule
334{
335public:
336 enum { numComponents = 0 };
337
338 FluidStateNullCompositionModule() = default;
339
343 ValueT moleFraction(unsigned /* phaseIdx */, unsigned /* compIdx */) const
344 { throw std::logic_error("Mole fractions are not provided by this fluid state"); }
345
349 ValueT massFraction(unsigned /* phaseIdx */, unsigned /* compIdx */) const
350 { throw std::logic_error("Mass fractions are not provided by this fluid state"); }
351
360 ValueT averageMolarMass(unsigned /* phaseIdx */) const
361 { throw std::logic_error("Mean molar masses are not provided by this fluid state"); }
362
372 ValueT molarity(unsigned /* phaseIdx */, unsigned /* compIdx */) const
373 { throw std::logic_error("Molarities are not provided by this fluid state"); }
374
383 void checkDefined() const
384 { }
385};
386
387
388} // namespace Opm
389
390#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 setMoleFraction(unsigned phaseIdx, unsigned compIdx, const ValueType &value)
Set the mole fraction of a component in a phase [] and update the average molar mass [kg/mol] accordi...
Definition FluidStateCompositionModules.hpp:118
ValueType molarity(unsigned phaseIdx, unsigned compIdx) const
The concentration of a component in a phase [mol/m^3].
Definition FluidStateCompositionModules.hpp:110
const ValueType & averageMolarMass(unsigned phaseIdx) const
The mean molar mass of a fluid phase [kg/mol].
Definition FluidStateCompositionModules.hpp:98
void setMoleFraction(unsigned compIdx, const ValueType &value)
Set the total mole fraction of a component.
Definition FluidStateCompositionModules.hpp:139
void assign(const FluidState &fs)
Retrieve all parameters from an arbitrary fluid state.
Definition FluidStateCompositionModules.hpp:159
ValueType wilsonK_(unsigned compIdx) const
Wilson formula to calculate K.
Definition FluidStateCompositionModules.hpp:224
const ValueType & moleFraction(unsigned phaseIdx, unsigned compIdx) const
The mole fraction of a component in a phase [].
Definition FluidStateCompositionModules.hpp:69
const ValueType & L() const
The L value of a composition [-].
Definition FluidStateCompositionModules.hpp:207
const ValueType & moleFraction(unsigned compIdx) const
The total mole fraction of a component [].
Definition FluidStateCompositionModules.hpp:75
void setLvalue(const ValueType &value)
Set the L value [-].
Definition FluidStateCompositionModules.hpp:215
ValueType massFraction(unsigned phaseIdx, unsigned compIdx) const
The mass fraction of a component in a phase [].
Definition FluidStateCompositionModules.hpp:81
void setKvalue(unsigned compIdx, const ValueType &value)
Set the K value of a component [-].
Definition FluidStateCompositionModules.hpp:199
void checkDefined() const
Make sure that all attributes are defined.
Definition FluidStateCompositionModules.hpp:182
ValueType molarity(unsigned phaseIdx, unsigned compIdx) const
The concentration of a component in a phase [mol/m^3].
Definition FluidStateCompositionModules.hpp:301
void assign(const FluidState &)
Retrieve all parameters from an arbitrary fluid state.
Definition FluidStateCompositionModules.hpp:309
ValueType averageMolarMass(unsigned phaseIdx) const
The mean molar mass of a fluid phase [kg/mol].
Definition FluidStateCompositionModules.hpp:289
ValueType moleFraction(unsigned phaseIdx, unsigned compIdx) const
The mole fraction of a component in a phase [].
Definition FluidStateCompositionModules.hpp:272
void checkDefined() const
Make sure that all attributes are defined.
Definition FluidStateCompositionModules.hpp:320
ValueType massFraction(unsigned phaseIdx, unsigned compIdx) const
The mass fraction of a component in a phase [].
Definition FluidStateCompositionModules.hpp:278
ValueT massFraction(unsigned, unsigned) const
The mass fraction of a component in a phase [].
Definition FluidStateCompositionModules.hpp:349
ValueT averageMolarMass(unsigned) const
The mean molar mass of a fluid phase [kg/mol].
Definition FluidStateCompositionModules.hpp:360
ValueT moleFraction(unsigned, unsigned) const
The mole fraction of a component in a phase [].
Definition FluidStateCompositionModules.hpp:343
void checkDefined() const
Make sure that all attributes are defined.
Definition FluidStateCompositionModules.hpp:383
ValueT molarity(unsigned, unsigned) const
The concentration of a component in a phase [mol/m^3].
Definition FluidStateCompositionModules.hpp:372
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30