checkFluidSystem.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*/
27#ifndef OPM_CHECK_FLUIDSYSTEM_HPP
28#define OPM_CHECK_FLUIDSYSTEM_HPP
29
30// include all fluid systems in opm-material
39
40// include all fluid states
47
48#include <dune/common/classname.hh>
49
50#include <iostream>
51#include <string>
52
57template <class ScalarT,
58 class FluidSystem,
61 : protected BaseFluidState
62{
63public:
64 enum { numPhases = FluidSystem::numPhases };
65 enum { numComponents = FluidSystem::numComponents };
66
67 typedef ScalarT Scalar;
68
70 {
71 // initially, do not allow anything
72 allowTemperature(false);
73 allowPressure(false);
74 allowComposition(false);
75 allowDensity(false);
76
77 // do not allow accessing any phase
78 restrictToPhase(1000);
79 }
80
81 void allowTemperature(bool yesno)
82 { allowTemperature_ = yesno; }
83
84 void allowPressure(bool yesno)
85 { allowPressure_ = yesno; }
86
87 void allowComposition(bool yesno)
88 { allowComposition_ = yesno; }
89
90 void allowDensity(bool yesno)
91 { allowDensity_ = yesno; }
92
93 void restrictToPhase(int phaseIdx)
94 { restrictPhaseIdx_ = phaseIdx; }
95
96 BaseFluidState& base()
97 { return *static_cast<BaseFluidState*>(this); }
98
99 const BaseFluidState& base() const
100 { return *static_cast<const BaseFluidState*>(this); }
101
102 auto temperature(unsigned phaseIdx) const
103 -> decltype(this->base().temperature(phaseIdx))
104 {
105 assert(allowTemperature_);
106 assert(restrictPhaseIdx_ < 0 || restrictPhaseIdx_ == static_cast<int>(phaseIdx));
107 return this->base().temperature(phaseIdx);
108 }
109
110 auto pressure(unsigned phaseIdx) const
111 -> decltype(this->base().pressure(phaseIdx))
112 {
113 assert(allowPressure_);
114 assert(restrictPhaseIdx_ < 0 || restrictPhaseIdx_ == static_cast<int>(phaseIdx));
115 return this->base().pressure(phaseIdx);
116 }
117
118 auto moleFraction(unsigned phaseIdx, unsigned compIdx) const
119 -> decltype(this->base().moleFraction(phaseIdx, compIdx))
120 {
121 assert(allowComposition_);
122 assert(restrictPhaseIdx_ < 0 || restrictPhaseIdx_ == static_cast<int>(phaseIdx));
123 return this->base().moleFraction(phaseIdx, compIdx);
124 }
125
126 auto massFraction(unsigned phaseIdx, unsigned compIdx) const
127 -> decltype(this->base().massFraction(phaseIdx, compIdx))
128 {
129 assert(allowComposition_);
130 assert(restrictPhaseIdx_ < 0 || restrictPhaseIdx_ == static_cast<int>(phaseIdx));
131 return this->base().massFraction(phaseIdx, compIdx);
132 }
133
134 auto averageMolarMass(unsigned phaseIdx) const
135 -> decltype(this->base().averageMolarMass(phaseIdx))
136 {
137 assert(allowComposition_);
138 assert(restrictPhaseIdx_ < 0 || restrictPhaseIdx_ == static_cast<int>(phaseIdx));
139 return this->base().averageMolarMass(phaseIdx);
140 }
141
142 auto molarity(unsigned phaseIdx, unsigned compIdx) const
143 -> decltype(this->base().molarity(phaseIdx, compIdx))
144 {
145 assert(allowDensity_ && allowComposition_);
146 assert(restrictPhaseIdx_ < 0 || restrictPhaseIdx_ == static_cast<int>(phaseIdx));
147 return this->base().molarity(phaseIdx, compIdx);
148 }
149
150 auto molarDensity(unsigned phaseIdx) const
151 -> decltype(this->base().molarDensity(phaseIdx))
152 {
153 assert(allowDensity_);
154 assert(restrictPhaseIdx_ < 0 || restrictPhaseIdx_ == static_cast<int>(phaseIdx));
155 return this->base().molarDensity(phaseIdx);
156 }
157
158 auto molarVolume(unsigned phaseIdx) const
159 -> decltype(this->base().molarVolume(phaseIdx))
160 {
161 assert(allowDensity_);
162 assert(restrictPhaseIdx_ < 0 || restrictPhaseIdx_ == static_cast<int>(phaseIdx));
163 return this->base().molarVolume(phaseIdx);
164 }
165
166 auto density(unsigned phaseIdx) const
167 -> decltype(this->base().density(phaseIdx))
168 {
169 assert(allowDensity_);
170 assert(restrictPhaseIdx_ < 0 || restrictPhaseIdx_ == static_cast<int>(phaseIdx));
171 return this->base().density(phaseIdx);
172 }
173
174 auto saturation(unsigned phaseIdx) const
175 -> decltype(this->base().saturation(phaseIdx))
176 {
177 assert(false);
178 return this->base().saturation(phaseIdx);
179 }
180
181 auto fugacity(unsigned phaseIdx, unsigned compIdx) const
182 -> decltype(this->base().fugacity(phaseIdx, compIdx))
183 {
184 assert(false);
185 return this->base().fugacity(phaseIdx, compIdx);
186 }
187
188 auto fugacityCoefficient(unsigned phaseIdx, unsigned compIdx) const
189 -> decltype(this->base().fugacityCoefficient(phaseIdx, compIdx))
190 {
191 assert(false);
192 return this->base().fugacityCoefficient(phaseIdx, compIdx);
193 }
194
195 auto enthalpy(unsigned phaseIdx) const
196 -> decltype(this->base().enthalpy(phaseIdx))
197 {
198 assert(false);
199 return this->base().enthalpy(phaseIdx);
200 }
201
202 auto internalEnergy(unsigned phaseIdx) const
203 -> decltype(this->base().internalEnergy(phaseIdx))
204 {
205 assert(false);
206 return this->base().internalEnergy(phaseIdx);
207 }
208
209 auto viscosity(unsigned phaseIdx) const
210 -> decltype(this->base().viscosity(phaseIdx))
211 {
212 assert(false);
213 return this->base().viscosity(phaseIdx);
214 }
215
216private:
217 bool allowSaturation_;
218 bool allowTemperature_;
219 bool allowPressure_;
220 bool allowComposition_;
221 bool allowDensity_;
222 int restrictPhaseIdx_;
223};
224
225template <class Scalar, class BaseFluidState>
226void checkFluidState(const BaseFluidState& fs)
227{
228 // fluid states must be copy-able
229 [[maybe_unused]] BaseFluidState tmpFs(fs);
230 tmpFs = fs;
231
232 // a fluid state must provide a checkDefined() method
233 fs.checkDefined();
234
235 // fluid states must export the types which they use as Scalars
236 typedef typename BaseFluidState::Scalar FsScalar;
237 static_assert(std::is_same<FsScalar, Scalar>::value,
238 "Fluid states must export the type they are given as scalar in an unmodified way");
239
240 // make sure the fluid state provides all mandatory methods
241 while (false) {
242 Scalar val = 1.0;
243
244 val = 2*val; // get rid of GCC warning (only occurs with paranoid warning flags)
245
246 val = fs.temperature(/*phaseIdx=*/0);
247 val = fs.pressure(/*phaseIdx=*/0);
248 val = fs.moleFraction(/*phaseIdx=*/0, /*compIdx=*/0);
249 val = fs.massFraction(/*phaseIdx=*/0, /*compIdx=*/0);
250 val = fs.averageMolarMass(/*phaseIdx=*/0);
251 val = fs.molarity(/*phaseIdx=*/0, /*compIdx=*/0);
252 val = fs.molarDensity(/*phaseIdx=*/0);
253 val = fs.molarVolume(/*phaseIdx=*/0);
254 val = fs.density(/*phaseIdx=*/0);
255 val = fs.saturation(/*phaseIdx=*/0);
256 val = fs.fugacity(/*phaseIdx=*/0, /*compIdx=*/0);
257 val = fs.fugacityCoefficient(/*phaseIdx=*/0, /*compIdx=*/0);
258 val = fs.enthalpy(/*phaseIdx=*/0);
259 val = fs.internalEnergy(/*phaseIdx=*/0);
260 val = fs.viscosity(/*phaseIdx=*/0);
261 };
262}
263
267template <class Scalar, class FluidSystem, class RhsEval, class LhsEval>
269{
270 std::cout << "Testing fluid system '" << Dune::className<FluidSystem>() << "'\n";
271
272 // make sure the fluid system provides the number of phases and
273 // the number of components
274 enum { numPhases = FluidSystem::numPhases };
275 enum { numComponents = FluidSystem::numComponents };
276
278 FluidState fs;
279 fs.allowTemperature(true);
280 fs.allowPressure(true);
281 fs.allowComposition(true);
282 fs.restrictToPhase(-1);
283
284 // initialize memory the fluid state
285 fs.base().setTemperature(273.15 + 20.0);
286 for (int phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
287 fs.base().setPressure(phaseIdx, 1e5);
288 fs.base().setSaturation(phaseIdx, 1.0/numPhases);
289 for (int compIdx = 0; compIdx < numComponents; ++ compIdx) {
290 fs.base().setMoleFraction(phaseIdx, compIdx, 1.0/numComponents);
291 }
292 }
293
294 static_assert(std::is_same<typename FluidSystem::Scalar, Scalar>::value,
295 "The type used for floating point used by the fluid system must be the same"
296 " as the one passed to the checkFluidSystem() function");
297
298 // check whether the parameter cache adheres to the API
299 typedef typename FluidSystem::template ParameterCache<LhsEval> ParameterCache;
300
301 ParameterCache paramCache;
302 try { paramCache.updateAll(fs); } catch (...) {};
303 try { paramCache.updateAll(fs, /*except=*/ParameterCache::None); } catch (...) {};
304 try { paramCache.updateAll(fs, /*except=*/ParameterCache::Temperature | ParameterCache::Pressure | ParameterCache::Composition); } catch (...) {};
305 try { paramCache.updateAllPressures(fs); } catch (...) {};
306
307 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
308 fs.restrictToPhase(static_cast<int>(phaseIdx));
309 try { paramCache.updatePhase(fs, phaseIdx); } catch (...) {};
310 try { paramCache.updatePhase(fs, phaseIdx, /*except=*/ParameterCache::None); } catch (...) {};
311 try { paramCache.updatePhase(fs, phaseIdx, /*except=*/ParameterCache::Temperature | ParameterCache::Pressure | ParameterCache::Composition); } catch (...) {};
312 try { paramCache.updateTemperature(fs, phaseIdx); } catch (...) {};
313 try { paramCache.updatePressure(fs, phaseIdx); } catch (...) {};
314 try { paramCache.updateComposition(fs, phaseIdx); } catch (...) {};
315 try { paramCache.updateSingleMoleFraction(fs, phaseIdx, /*compIdx=*/0); } catch (...) {};
316 }
317
318 // some value to make sure the return values of the fluid system
319 // are convertible to scalars
320 LhsEval val = 0.0;
321 Scalar scalarVal = 0.0;
322
323 scalarVal = 2*scalarVal; // get rid of GCC warning (only occurs with paranoid warning flags)
324 val = 2*val; // get rid of GCC warning (only occurs with paranoid warning flags)
325
326 // actually check the fluid system API
327 try { FluidSystem::init(); } catch (...) {};
328 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
329 fs.restrictToPhase(static_cast<int>(phaseIdx));
330 fs.allowPressure(FluidSystem::isCompressible(phaseIdx));
331 fs.allowComposition(true);
332 fs.allowDensity(false);
333 try { [[maybe_unused]] auto tmpVal = FluidSystem::density(fs, paramCache, phaseIdx); static_assert(std::is_same<decltype(tmpVal), RhsEval>::value, "The default return value must be the scalar used by the fluid state!"); } catch (...) {};
334 try { val = FluidSystem::template density<FluidState, LhsEval>(fs, paramCache, phaseIdx); } catch (...) {};
335 try { scalarVal = FluidSystem::template density<FluidState, Scalar>(fs, paramCache, phaseIdx); } catch (...) {};
336
337 fs.allowPressure(true);
338 fs.allowDensity(true);
339 try { [[maybe_unused]] auto tmpVal = FluidSystem::viscosity(fs, paramCache, phaseIdx); static_assert(std::is_same<decltype(tmpVal), RhsEval>::value, "The default return value must be the scalar used by the fluid state!"); } catch (...) {};
340 try { [[maybe_unused]] auto tmpVal = FluidSystem::enthalpy(fs, paramCache, phaseIdx); static_assert(std::is_same<decltype(tmpVal), RhsEval>::value, "The default return value must be the scalar used by the fluid state!"); } catch (...) {};
341 try { [[maybe_unused]] auto tmpVal = FluidSystem::heatCapacity(fs, paramCache, phaseIdx); static_assert(std::is_same<decltype(tmpVal), RhsEval>::value, "The default return value must be the scalar used by the fluid state!"); } catch (...) {};
342 try { [[maybe_unused]] auto tmpVal = FluidSystem::thermalConductivity(fs, paramCache, phaseIdx); static_assert(std::is_same<decltype(tmpVal), RhsEval>::value, "The default return value must be the scalar used by the fluid state!"); } catch (...) {};
343 try { val = FluidSystem::template viscosity<FluidState, LhsEval>(fs, paramCache, phaseIdx); } catch (...) {};
344 try { val = FluidSystem::template enthalpy<FluidState, LhsEval>(fs, paramCache, phaseIdx); } catch (...) {};
345 try { val = FluidSystem::template heatCapacity<FluidState, LhsEval>(fs, paramCache, phaseIdx); } catch (...) {};
346 try { val = FluidSystem::template thermalConductivity<FluidState, LhsEval>(fs, paramCache, phaseIdx); } catch (...) {};
347 try { scalarVal = FluidSystem::template viscosity<FluidState, Scalar>(fs, paramCache, phaseIdx); } catch (...) {};
348 try { scalarVal = FluidSystem::template enthalpy<FluidState, Scalar>(fs, paramCache, phaseIdx); } catch (...) {};
349 try { scalarVal = FluidSystem::template heatCapacity<FluidState, Scalar>(fs, paramCache, phaseIdx); } catch (...) {};
350 try { scalarVal = FluidSystem::template thermalConductivity<FluidState, Scalar>(fs, paramCache, phaseIdx); } catch (...) {};
351
352 for (unsigned compIdx = 0; compIdx < numComponents; ++ compIdx) {
353 fs.allowComposition(!FluidSystem::isIdealMixture(phaseIdx));
354 try { [[maybe_unused]] auto tmpVal = FluidSystem::fugacityCoefficient(fs, paramCache, phaseIdx, compIdx); static_assert(std::is_same<decltype(tmpVal), RhsEval>::value, "The default return value must be the scalar used by the fluid state!"); } catch (...) {};
355 try { val = FluidSystem::template fugacityCoefficient<FluidState, LhsEval>(fs, paramCache, phaseIdx, compIdx); } catch (...) {};
356 try { scalarVal = FluidSystem::template fugacityCoefficient<FluidState, Scalar>(fs, paramCache, phaseIdx, compIdx); } catch (...) {};
357 fs.allowComposition(true);
358 try { [[maybe_unused]] auto tmpVal = FluidSystem::diffusionCoefficient(fs, paramCache, phaseIdx, compIdx); static_assert(std::is_same<decltype(tmpVal), RhsEval>::value, "The default return value must be the scalar used by the fluid state!"); } catch (...) {};
359 try { val = FluidSystem::template diffusionCoefficient<FluidState, LhsEval>(fs, paramCache, phaseIdx, compIdx); } catch (...) {};
360 try { scalarVal = FluidSystem::template fugacityCoefficient<FluidState, Scalar>(fs, paramCache, phaseIdx, compIdx); } catch (...) {};
361 }
362 }
363
364 // test for phaseName(), isLiquid() and isIdealGas()
365 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
366 [[maybe_unused]] std::string name = FluidSystem::phaseName(phaseIdx);
367 bool bVal = FluidSystem::isLiquid(phaseIdx);
368 bVal = FluidSystem::isIdealGas(phaseIdx);
369 bVal = !bVal; // get rid of GCC warning (only occurs with paranoid warning flags)
370 }
371
372 // test for molarMass() and componentName()
373 for (unsigned compIdx = 0; compIdx < numComponents; ++ compIdx) {
374 val = FluidSystem::molarMass(compIdx);
375 std::string name = FluidSystem::componentName(compIdx);
376 }
377
378 std::cout << "----------------------------------\n";
379}
380
381#endif
Represents all relevant thermodynamic quantities of a multi-phase, multi-component fluid system assum...
Represents all relevant thermodynamic quantities of a multi-phase, multi-component fluid system assum...
Represents all relevant thermodynamic quantities of a multi-phase, multi-component fluid system not a...
void checkFluidState(const BaseFluidState &fs)
Definition: checkFluidSystem.hpp:226
void checkFluidSystem()
Checks whether a fluid system adheres to the specification.
Definition: checkFluidSystem.hpp:268
This is a fluid state which makes sure that only the quantities allowed are accessed.
Definition: checkFluidSystem.hpp:62
auto fugacity(unsigned phaseIdx, unsigned compIdx) const -> decltype(this->base().fugacity(phaseIdx, compIdx))
Definition: checkFluidSystem.hpp:181
void restrictToPhase(int phaseIdx)
Definition: checkFluidSystem.hpp:93
auto molarity(unsigned phaseIdx, unsigned compIdx) const -> decltype(this->base().molarity(phaseIdx, compIdx))
Definition: checkFluidSystem.hpp:142
auto massFraction(unsigned phaseIdx, unsigned compIdx) const -> decltype(this->base().massFraction(phaseIdx, compIdx))
Definition: checkFluidSystem.hpp:126
void allowDensity(bool yesno)
Definition: checkFluidSystem.hpp:90
HairSplittingFluidState()
Definition: checkFluidSystem.hpp:69
auto moleFraction(unsigned phaseIdx, unsigned compIdx) const -> decltype(this->base().moleFraction(phaseIdx, compIdx))
Definition: checkFluidSystem.hpp:118
auto averageMolarMass(unsigned phaseIdx) const -> decltype(this->base().averageMolarMass(phaseIdx))
Definition: checkFluidSystem.hpp:134
@ numComponents
Definition: checkFluidSystem.hpp:65
auto fugacityCoefficient(unsigned phaseIdx, unsigned compIdx) const -> decltype(this->base().fugacityCoefficient(phaseIdx, compIdx))
Definition: checkFluidSystem.hpp:188
auto temperature(unsigned phaseIdx) const -> decltype(this->base().temperature(phaseIdx))
Definition: checkFluidSystem.hpp:102
auto saturation(unsigned phaseIdx) const -> decltype(this->base().saturation(phaseIdx))
Definition: checkFluidSystem.hpp:174
ScalarT Scalar
Definition: checkFluidSystem.hpp:67
void allowComposition(bool yesno)
Definition: checkFluidSystem.hpp:87
const BaseFluidState & base() const
Definition: checkFluidSystem.hpp:99
auto density(unsigned phaseIdx) const -> decltype(this->base().density(phaseIdx))
Definition: checkFluidSystem.hpp:166
auto viscosity(unsigned phaseIdx) const -> decltype(this->base().viscosity(phaseIdx))
Definition: checkFluidSystem.hpp:209
BaseFluidState & base()
Definition: checkFluidSystem.hpp:96
@ numPhases
Definition: checkFluidSystem.hpp:64
auto pressure(unsigned phaseIdx) const -> decltype(this->base().pressure(phaseIdx))
Definition: checkFluidSystem.hpp:110
auto internalEnergy(unsigned phaseIdx) const -> decltype(this->base().internalEnergy(phaseIdx))
Definition: checkFluidSystem.hpp:202
void allowTemperature(bool yesno)
Definition: checkFluidSystem.hpp:81
auto molarDensity(unsigned phaseIdx) const -> decltype(this->base().molarDensity(phaseIdx))
Definition: checkFluidSystem.hpp:150
auto molarVolume(unsigned phaseIdx) const -> decltype(this->base().molarVolume(phaseIdx))
Definition: checkFluidSystem.hpp:158
void allowPressure(bool yesno)
Definition: checkFluidSystem.hpp:84
auto enthalpy(unsigned phaseIdx) const -> decltype(this->base().enthalpy(phaseIdx))
Definition: checkFluidSystem.hpp:195