28#ifndef OPM_FLUID_STATE_COMPOSITION_MODULES_HPP
29#define OPM_FLUID_STATE_COMPOSITION_MODULES_HPP
43template <
class ValueType,
46class FluidStateExplicitCompositionModule
48 enum { numPhases = FluidSystem::numPhases };
49 enum { numComponents = FluidSystem::numComponents };
52 FluidStateExplicitCompositionModule()
54 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
55 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
56 moleFraction_[phaseIdx][compIdx] = 0.0;
69 const ValueType&
moleFraction(
unsigned phaseIdx,
unsigned compIdx)
const
70 {
return moleFraction_[phaseIdx][compIdx]; }
76 {
return totalModelFractions_[compIdx]; }
84 abs(sumMoleFractions_[phaseIdx])
85 *moleFraction_[phaseIdx][compIdx]
86 *FluidSystem::molarMass(compIdx)
87 / max(1e-40, abs(averageMolarMass_[phaseIdx]));
99 {
return averageMolarMass_[phaseIdx]; }
110 ValueType
molarity(
unsigned phaseIdx,
unsigned compIdx)
const
111 {
return asImp_().molarDensity(phaseIdx)*
moleFraction(phaseIdx, compIdx); }
125 moleFraction_[phaseIdx][compIdx] = value;
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);
142 totalModelFractions_[compIdx] = value;
145 void setCompressFactor(
unsigned phaseIdx,
const ValueType& value) {
147 Z_[phaseIdx] = value;
150 ValueType compressFactor(
unsigned phaseIdx)
const {
158 template <
class Flu
idState>
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));
168 averageMolarMass_[phaseIdx] += moleFraction_[phaseIdx][compIdx]*FluidSystem::molarMass(compIdx);
169 sumMoleFractions_[phaseIdx] += moleFraction_[phaseIdx][compIdx];
191 const ValueType& K(
unsigned compIdx)
const
199 void setKvalue(
unsigned compIdx,
const ValueType& value)
207 const ValueType&
L()
const
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);
232 const auto tmp = exp(5.37 * (1+acf) * (1-T_crit/T)) * (p_crit/p);
237 const Implementation& asImp_()
const
239 return *
static_cast<const Implementation*
>(
this);
242 std::array<std::array<ValueType,numComponents>,numPhases> moleFraction_{};
243 std::array<ValueType,numPhases> averageMolarMass_{};
244 std::array<ValueType,numPhases> sumMoleFractions_{};
245 std::array<ValueType, numComponents> totalModelFractions_{};
246 std::array<ValueType,numPhases> Z_{};
247 std::array<ValueType,numComponents> K_{};
255template <
class ValueType,
257 class Implementation>
258class FluidStateImmiscibleCompositionModule
260 enum { numPhases = FluidSystem::numPhases };
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");
267 FluidStateImmiscibleCompositionModule() =
default;
273 {
return (phaseIdx == compIdx)?1.0:0.0; }
279 {
return (phaseIdx == compIdx)?1.0:0.0; }
290 {
return FluidSystem::molarMass(phaseIdx); }
301 ValueType
molarity(
unsigned phaseIdx,
unsigned compIdx)
const
302 {
return asImp_().molarDensity(phaseIdx)*
moleFraction(phaseIdx, compIdx); }
308 template <
class Flu
idState>
324 const Implementation& asImp_()
const
325 {
return *
static_cast<const Implementation*
>(
this); }
332template <
class ValueT>
333class FluidStateNullCompositionModule
336 enum { numComponents = 0 };
338 FluidStateNullCompositionModule() =
default;
344 {
throw std::logic_error(
"Mole fractions are not provided by this fluid state"); }
350 {
throw std::logic_error(
"Mass fractions are not provided by this fluid state"); }
361 {
throw std::logic_error(
"Mean molar masses are not provided by this fluid state"); }
373 {
throw std::logic_error(
"Molarities are not provided by this fluid state"); }
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