version 3.8
spe5parametercache.hh
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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef SPE5_PARAMETER_CACHE_HH
13#define SPE5_PARAMETER_CACHE_HH
14
15#include <cassert>
16
19
22
23namespace Dumux {
24
30template <class Scalar, class FluidSystem>
32 : public ParameterCacheBase<Spe5ParameterCache<Scalar, FluidSystem> >
33{
36
38
39 enum { numPhases = FluidSystem::numPhases };
40
41 enum { wPhaseIdx = FluidSystem::wPhaseIdx };
42 enum { oPhaseIdx = FluidSystem::oPhaseIdx };
43 enum { gPhaseIdx = FluidSystem::gPhaseIdx };
44
45public:
46 // types of the parameter objects for each phase
47 using OilPhaseParams = PengRobinsonParamsMixture<Scalar, FluidSystem, oPhaseIdx, /*useSpe5=*/true>;
48 using GasPhaseParams = PengRobinsonParamsMixture<Scalar, FluidSystem, gPhaseIdx, /*useSpe5=*/true>;
49
54 {
55 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
56 VmUpToDate_[phaseIdx] = false;
57 }
58
66 template <class FluidState>
67 void updatePhase(const FluidState &fs,
68 int phaseIdx,
69 int except = ParentType::None)
70 {
71 updateEosParams(fs, phaseIdx, except);
72
73 // if we don't need to recalculate the molar volume, we exit
74 // here
75 if (VmUpToDate_[phaseIdx])
76 return;
77
78 // update the phase's molar volume
79 updateMolarVolume_(fs, phaseIdx);
80 }
81
95 template <class FluidState>
96 void updateSingleMoleFraction(const FluidState &fs,
97 int phaseIdx,
98 int compIdx)
99 {
100 if (phaseIdx == oPhaseIdx)
102 else if (phaseIdx == gPhaseIdx)
104
105 // update the phase's molar volume
106 updateMolarVolume_(fs, phaseIdx);
107 }
108
113 Scalar a(int phaseIdx) const
114 {
115 switch (phaseIdx)
116 {
117 case oPhaseIdx: return oilPhaseParams_.a();
118 case gPhaseIdx: return gasPhaseParams_.a();
119 default:
120 DUNE_THROW(Dune::InvalidStateException,
121 "The a() parameter is only defined for "
122 "oil and gas phases");
123 }
124 }
125
130 Scalar b(int phaseIdx) const
131 {
132 switch (phaseIdx)
133 {
134 case oPhaseIdx: return oilPhaseParams_.b();
135 case gPhaseIdx: return gasPhaseParams_.b();
136 default:
137 DUNE_THROW(Dune::InvalidStateException,
138 "The b() parameter is only defined for "
139 "oil and gas phases");
140 }
141 }
142
150 Scalar aPure(int phaseIdx, int compIdx) const
151 {
152 switch (phaseIdx)
153 {
154 case oPhaseIdx: return oilPhaseParams_.pureParams(compIdx).a();
155 case gPhaseIdx: return gasPhaseParams_.pureParams(compIdx).a();
156 default:
157 DUNE_THROW(Dune::InvalidStateException,
158 "The a() parameter is only defined for "
159 "oil and gas phases");
160 }
161 }
162
169 Scalar bPure(int phaseIdx, int compIdx) const
170 {
171 switch (phaseIdx)
172 {
173 case oPhaseIdx: return oilPhaseParams_.pureParams(compIdx).b();
174 case gPhaseIdx: return gasPhaseParams_.pureParams(compIdx).b();
175 default:
176 DUNE_THROW(Dune::InvalidStateException,
177 "The b() parameter is only defined for "
178 "oil and gas phases");
179 }
180 }
181
186 Scalar molarVolume(int phaseIdx) const
187 {
188 assert(VmUpToDate_[phaseIdx]);
189 return Vm_[phaseIdx];
190 }
191
192
198 { return oilPhaseParams_; }
199
205 { return gasPhaseParams_; }
206
211 template <class FluidState>
212 void updateEosParams(const FluidState &fs,
213 int phaseIdx,
214 int exceptQuantities = ParentType::None)
215 {
216 if (!(exceptQuantities & ParentType::Temperature))
217 {
218 updatePure_(fs, phaseIdx);
219 updateMix_(fs, phaseIdx);
220 VmUpToDate_[phaseIdx] = false;
221 }
222 else if (!(exceptQuantities & ParentType::Composition))
223 {
224 updateMix_(fs, phaseIdx);
225 VmUpToDate_[phaseIdx] = false;
226 }
227 else if (!(exceptQuantities & ParentType::Pressure)) {
228 VmUpToDate_[phaseIdx] = false;
229 }
230 }
231
232protected:
239 template <class FluidState>
240 void updatePure_(const FluidState &fs, int phaseIdx)
241 {
242 Scalar T = fs.temperature(phaseIdx);
243 Scalar p = fs.pressure(phaseIdx);
244
245 switch (phaseIdx)
246 {
247 case oPhaseIdx: oilPhaseParams_.updatePure(T, p); break;
248 case gPhaseIdx: gasPhaseParams_.updatePure(T, p); break;
249 }
250 }
251
259 template <class FluidState>
260 void updateMix_(const FluidState &fs, int phaseIdx)
261 {
262 switch (phaseIdx)
263 {
264 case oPhaseIdx:
266 break;
267 case gPhaseIdx:
269 break;
270 case wPhaseIdx:
271 break;
272 }
273 }
274
275 template <class FluidState>
276 void updateMolarVolume_(const FluidState &fs,
277 int phaseIdx)
278 {
279 VmUpToDate_[phaseIdx] = true;
280
281 // calculate molar volume of the phase (we will need this for the
282 // fugacity coefficients and the density anyway)
283 switch (phaseIdx) {
284 case gPhaseIdx: {
285 // calculate molar volumes for the given composition. although
286 // this isn't a Peng-Robinson parameter strictly speaking, the
287 // molar volume appears in basically every quantity the fluid
288 // system can get queried, so it is okay to calculate it
289 // here...
290 Vm_[gPhaseIdx] =
292 *this,
293 phaseIdx,
294 /*isGasPhase=*/true);
295 break;
296 }
297 case oPhaseIdx: {
298 // calculate molar volumes for the given composition. although
299 // this isn't a Peng-Robinson parameter strictly speaking, the
300 // molar volume appears in basically every quantity the fluid
301 // system can get queried, so it is okay to calculate it
302 // here...
303 Vm_[oPhaseIdx] =
305 *this,
306 phaseIdx,
307 /*isGasPhase=*/false);
308 break;
309 }
310 case wPhaseIdx: {
311 // Density of water in the stock tank (i.e. atmospheric
312 // pressure) is specified as 62.4 lb/ft^3 by the SPE-5
313 // paper. Also 1 lb = 0.4535923 and 1 ft = 0.3048 m.
314 const Scalar stockTankWaterDensity = 62.4 * 0.45359237 / 0.028316847;
315 // Water compressibility is specified as 3.3e-6 per psi
316 // overpressure, where 1 psi = 6894.7573 Pa
317 Scalar overPressure = fs.pressure(wPhaseIdx) - 1.013e5; // [Pa]
318 Scalar waterDensity =
319 stockTankWaterDensity * (1 + 3.3e-6*overPressure/6894.7573);
320
321 // convert water density [kg/m^3] to molar volume [m^3/mol]
322 Vm_[wPhaseIdx] = fs.averageMolarMass(wPhaseIdx)/waterDensity;
323 break;
324 }
325 default:
326 DUNE_THROW(Dune::InvalidStateException, "invalid phaseIdx " << phaseIdx);
327 }
328 }
329
330 bool VmUpToDate_[numPhases];
331 Scalar Vm_[numPhases];
332
335};
336
337} // end namespace Dumux
338
339#endif
The base class of the parameter cache classes for fluid systems.
Definition: parametercachebase.hh:23
@ Composition
Definition: parametercachebase.hh:29
@ Pressure
Definition: parametercachebase.hh:28
@ Temperature
Definition: parametercachebase.hh:27
@ None
Definition: parametercachebase.hh:26
Implements the Peng-Robinson equation of state for liquids and gases.
Definition: pengrobinson.hh:48
static Scalar computeMolarVolume(const FluidState &fs, Params &params, int phaseIdx, bool isGasPhase, bool handleUnphysicalPhase=true)
Computes molar volumes where the Peng-Robinson EOS is true.
Definition: pengrobinson.hh:144
Scalar a() const
Returns the attractive parameter 'a' of the Peng-Robinson fluid.
Definition: pengrobinsonparams.hh:38
Scalar b() const
Returns the repulsive parameter 'b' of the Peng-Robinson fluid.
Definition: pengrobinsonparams.hh:46
void updateMix(const FluidState &fs)
Calculates the "a" and "b" Peng-Robinson parameters for the mixture.
Definition: pengrobinsonparamsmixture.hh:124
const PureParams & pureParams(int compIdx) const
Return the Peng-Robinson parameters of a pure substance,.
Definition: pengrobinsonparamsmixture.hh:180
void updateSingleMoleFraction(const FluidState &fs, int compIdx)
Calculates the "a" and "b" Peng-Robinson parameters for the mixture provided that only a single mole ...
Definition: pengrobinsonparamsmixture.hh:170
void updatePure(const FluidState &fluidState)
Update Peng-Robinson parameters for the pure components.
Definition: pengrobinsonparamsmixture.hh:64
Specifies the parameters required by the SPE5 problem which are despondent on the thermodynamic state...
Definition: spe5parametercache.hh:33
bool VmUpToDate_[numPhases]
Definition: spe5parametercache.hh:330
Scalar b(int phaseIdx) const
The Peng-Robinson co-volume for a phase.
Definition: spe5parametercache.hh:130
const GasPhaseParams & gasPhaseParams() const
Returns the Peng-Robinson mixture parameters for the gas phase.
Definition: spe5parametercache.hh:204
Scalar aPure(int phaseIdx, int compIdx) const
The Peng-Robinson attractive parameter for a pure component given the same temperature and pressure o...
Definition: spe5parametercache.hh:150
void updateMolarVolume_(const FluidState &fs, int phaseIdx)
Definition: spe5parametercache.hh:276
void updatePhase(const FluidState &fs, int phaseIdx, int except=ParentType::None)
Update all parameters required by the fluid system to calculate some quantities for the phase.
Definition: spe5parametercache.hh:67
Scalar Vm_[numPhases]
Definition: spe5parametercache.hh:331
Scalar a(int phaseIdx) const
The Peng-Robinson attractive parameter for a phase.
Definition: spe5parametercache.hh:113
void updateEosParams(const FluidState &fs, int phaseIdx, int exceptQuantities=ParentType::None)
Update all parameters required by the equation of state to calculate some quantities for the phase.
Definition: spe5parametercache.hh:212
void updatePure_(const FluidState &fs, int phaseIdx)
Update all parameters of a phase which only depend on temperature and/or pressure.
Definition: spe5parametercache.hh:240
const OilPhaseParams & oilPhaseParams() const
Returns the Peng-Robinson mixture parameters for the oil phase.
Definition: spe5parametercache.hh:197
Spe5ParameterCache()
The constructor.
Definition: spe5parametercache.hh:53
GasPhaseParams gasPhaseParams_
Definition: spe5parametercache.hh:334
OilPhaseParams oilPhaseParams_
Definition: spe5parametercache.hh:333
void updateMix_(const FluidState &fs, int phaseIdx)
Update all parameters of a phase which depend on the fluid composition. It is assumed that updatePure...
Definition: spe5parametercache.hh:260
void updateSingleMoleFraction(const FluidState &fs, int phaseIdx, int compIdx)
Update all cached parameters of a specific fluid phase which depend on the mole fraction of a single ...
Definition: spe5parametercache.hh:96
Scalar molarVolume(int phaseIdx) const
Returns the molar volume of a phase .
Definition: spe5parametercache.hh:186
Scalar bPure(int phaseIdx, int compIdx) const
The Peng-Robinson co-volume for a pure component given the same temperature and pressure of the phase...
Definition: spe5parametercache.hh:169
Material properties of pure water .
Definition: adapt.hh:17
The base class of the parameter cache classes for fluid systems.
Implements the Peng-Robinson equation of state for liquids and gases.
The mixing rule for the oil and the gas phases of the SPE5 problem.