2011-12-16 13:24:29 +00:00
|
|
|
// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
|
|
|
|
|
// vi: set et ts=4 sw=4 sts=4:
|
2011-12-16 13:24:20 +00:00
|
|
|
/*****************************************************************************
|
2012-08-31 14:38:58 +02:00
|
|
|
* Copyright (C) 2011-2012 by Andreas Lauser *
|
2011-12-16 13:24:20 +00:00
|
|
|
* *
|
|
|
|
|
* This program is free software: you can redistribute it and/or modify *
|
|
|
|
|
* it under the terms of the GNU General Public License as published by *
|
|
|
|
|
* the Free Software Foundation, either version 2 of the License, or *
|
|
|
|
|
* (at your option) any later version. *
|
|
|
|
|
* *
|
|
|
|
|
* This program is distributed in the hope that it will be useful, *
|
|
|
|
|
* but WITHOUT ANY WARRANTY; without even the implied warranty of *
|
|
|
|
|
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
|
|
|
|
|
* GNU General Public License for more details. *
|
|
|
|
|
* *
|
|
|
|
|
* You should have received a copy of the GNU General Public License *
|
|
|
|
|
* along with this program. If not, see <http://www.gnu.org/licenses/>. *
|
|
|
|
|
*****************************************************************************/
|
|
|
|
|
/*!
|
|
|
|
|
* \file
|
2012-11-02 13:45:42 +01:00
|
|
|
* \copydoc checkFluidSystem
|
2011-12-16 13:24:20 +00:00
|
|
|
*/
|
2012-11-18 16:42:57 +01:00
|
|
|
#ifndef EWOMS_CHECK_FLUIDSYSTEM_HH
|
|
|
|
|
#define EWOMS_CHECK_FLUIDSYSTEM_HH
|
2011-12-16 13:24:20 +00:00
|
|
|
|
2012-11-18 16:42:57 +01:00
|
|
|
// include all fluid systems in ewoms-stable
|
|
|
|
|
#include <ewoms/material/fluidsystems/1pfluidsystem.hh>
|
|
|
|
|
#include <ewoms/material/fluidsystems/2pimmisciblefluidsystem.hh>
|
|
|
|
|
#include <ewoms/material/fluidsystems/h2on2fluidsystem.hh>
|
|
|
|
|
#include <ewoms/material/fluidsystems/h2on2liquidphasefluidsystem.hh>
|
|
|
|
|
#include <ewoms/material/fluidsystems/h2oairfluidsystem.hh>
|
|
|
|
|
#include <ewoms/material/fluidsystems/h2oairmesitylenefluidsystem.hh>
|
|
|
|
|
#include <ewoms/material/fluidsystems/h2oairxylenefluidsystem.hh>
|
|
|
|
|
#include <ewoms/material/fluidsystems/spe5fluidsystem.hh>
|
2011-12-16 13:24:20 +00:00
|
|
|
|
2011-12-16 13:24:44 +00:00
|
|
|
// include all fluid states
|
2012-11-18 16:42:57 +01:00
|
|
|
#include <ewoms/material/fluidstates/pressureoverlayfluidstate.hh>
|
|
|
|
|
#include <ewoms/material/fluidstates/saturationoverlayfluidstate.hh>
|
|
|
|
|
#include <ewoms/material/fluidstates/temperatureoverlayfluidstate.hh>
|
|
|
|
|
#include <ewoms/material/fluidstates/compositionalfluidstate.hh>
|
|
|
|
|
#include <ewoms/material/fluidstates/nonequilibriumfluidstate.hh>
|
|
|
|
|
#include <ewoms/material/fluidstates/immisciblefluidstate.hh>
|
2011-12-16 13:24:20 +00:00
|
|
|
|
|
|
|
|
#include <dune/common/classname.hh>
|
|
|
|
|
|
2012-03-20 11:13:22 +01:00
|
|
|
#include <iostream>
|
|
|
|
|
#include <string>
|
|
|
|
|
|
2012-11-02 13:45:42 +01:00
|
|
|
/*!
|
|
|
|
|
* \brief This is a fluid state which makes sure that only the quantities
|
|
|
|
|
* allowed are accessed.
|
|
|
|
|
*/
|
2012-02-03 12:34:25 +00:00
|
|
|
template <class Scalar,
|
2011-12-16 13:24:56 +00:00
|
|
|
class FluidSystem,
|
2012-11-18 16:42:57 +01:00
|
|
|
class BaseFluidState = Ewoms::CompositionalFluidState<Scalar, FluidSystem> >
|
2012-02-03 12:34:25 +00:00
|
|
|
class HairSplittingFluidState
|
2011-12-16 13:24:56 +00:00
|
|
|
: protected BaseFluidState
|
2011-12-16 13:24:20 +00:00
|
|
|
{
|
|
|
|
|
public:
|
|
|
|
|
enum { numPhases = FluidSystem::numPhases };
|
|
|
|
|
enum { numComponents = FluidSystem::numComponents };
|
2012-02-03 12:34:25 +00:00
|
|
|
|
2011-12-16 13:24:20 +00:00
|
|
|
HairSplittingFluidState()
|
|
|
|
|
{
|
|
|
|
|
// initially, do not allow anything
|
|
|
|
|
allowTemperature(false);
|
|
|
|
|
allowPressure(false);
|
|
|
|
|
allowComposition(false);
|
|
|
|
|
allowDensity(false);
|
2012-02-03 12:34:25 +00:00
|
|
|
|
2011-12-16 13:24:23 +00:00
|
|
|
// do not allow accessing any phase
|
|
|
|
|
restrictToPhase(1000);
|
2011-12-16 13:24:20 +00:00
|
|
|
}
|
|
|
|
|
|
2011-12-16 13:24:23 +00:00
|
|
|
void allowTemperature(bool yesno)
|
|
|
|
|
{ allowTemperature_ = yesno; }
|
2011-12-16 13:24:20 +00:00
|
|
|
|
|
|
|
|
void allowPressure(bool yesno)
|
|
|
|
|
{ allowPressure_ = yesno; }
|
|
|
|
|
|
|
|
|
|
void allowComposition(bool yesno)
|
|
|
|
|
{ allowComposition_ = yesno; }
|
2012-02-03 12:34:25 +00:00
|
|
|
|
|
|
|
|
void allowDensity(bool yesno)
|
2011-12-16 13:24:20 +00:00
|
|
|
{ allowDensity_ = yesno; }
|
|
|
|
|
|
2012-02-03 12:34:25 +00:00
|
|
|
void restrictToPhase(int phaseIdx)
|
2011-12-16 13:24:23 +00:00
|
|
|
{ restrictPhaseIdx_ = phaseIdx; }
|
2011-12-16 13:24:20 +00:00
|
|
|
|
|
|
|
|
Scalar temperature(int phaseIdx) const
|
2011-12-16 13:24:23 +00:00
|
|
|
{
|
|
|
|
|
assert(allowTemperature_);
|
|
|
|
|
assert(restrictPhaseIdx_ < 0 || restrictPhaseIdx_ == phaseIdx);
|
2012-09-05 18:08:13 +02:00
|
|
|
DUNE_UNUSED Scalar tmp = BaseFluidState::temperature(phaseIdx);
|
|
|
|
|
return 1e100;
|
2011-12-16 13:24:23 +00:00
|
|
|
}
|
2012-02-03 12:34:25 +00:00
|
|
|
|
2011-12-16 13:24:20 +00:00
|
|
|
Scalar pressure(int phaseIdx) const
|
2012-02-03 12:34:25 +00:00
|
|
|
{
|
2011-12-16 13:24:23 +00:00
|
|
|
assert(allowPressure_);
|
|
|
|
|
assert(restrictPhaseIdx_ < 0 || restrictPhaseIdx_ == phaseIdx);
|
2012-09-05 18:08:13 +02:00
|
|
|
DUNE_UNUSED Scalar tmp = BaseFluidState::pressure(phaseIdx);
|
|
|
|
|
return 1e100;
|
2011-12-16 13:24:23 +00:00
|
|
|
}
|
2011-12-16 13:24:20 +00:00
|
|
|
|
|
|
|
|
Scalar moleFraction(int phaseIdx, int compIdx) const
|
2012-02-03 12:34:25 +00:00
|
|
|
{
|
|
|
|
|
assert(allowComposition_);
|
2011-12-16 13:24:23 +00:00
|
|
|
assert(restrictPhaseIdx_ < 0 || restrictPhaseIdx_ == phaseIdx);
|
2012-09-05 18:08:13 +02:00
|
|
|
DUNE_UNUSED Scalar tmp = BaseFluidState::moleFraction(phaseIdx, compIdx);
|
|
|
|
|
return 1e100;
|
2011-12-16 13:24:23 +00:00
|
|
|
}
|
2011-12-16 13:24:20 +00:00
|
|
|
|
|
|
|
|
Scalar massFraction(int phaseIdx, int compIdx) const
|
2012-02-03 12:34:25 +00:00
|
|
|
{
|
|
|
|
|
assert(allowComposition_);
|
2011-12-16 13:24:23 +00:00
|
|
|
assert(restrictPhaseIdx_ < 0 || restrictPhaseIdx_ == phaseIdx);
|
2012-09-05 18:08:13 +02:00
|
|
|
DUNE_UNUSED Scalar tmp = BaseFluidState::massFraction(phaseIdx, compIdx);
|
|
|
|
|
return 1e100;
|
2011-12-16 13:24:23 +00:00
|
|
|
}
|
2011-12-16 13:24:20 +00:00
|
|
|
|
|
|
|
|
Scalar averageMolarMass(int phaseIdx) const
|
2011-12-16 13:24:23 +00:00
|
|
|
{
|
2012-02-03 12:34:25 +00:00
|
|
|
assert(allowComposition_);
|
2011-12-16 13:24:23 +00:00
|
|
|
assert(restrictPhaseIdx_ < 0 || restrictPhaseIdx_ == phaseIdx);
|
2012-09-05 18:08:13 +02:00
|
|
|
DUNE_UNUSED Scalar tmp = BaseFluidState::averageMolarMass(phaseIdx);
|
|
|
|
|
return 1e100;
|
2011-12-16 13:24:23 +00:00
|
|
|
}
|
2011-12-16 13:24:20 +00:00
|
|
|
|
|
|
|
|
Scalar molarity(int phaseIdx, int compIdx) const
|
2011-12-16 13:24:23 +00:00
|
|
|
{
|
|
|
|
|
assert(allowDensity_ && allowComposition_);
|
|
|
|
|
assert(restrictPhaseIdx_ < 0 || restrictPhaseIdx_ == phaseIdx);
|
2012-09-05 18:08:13 +02:00
|
|
|
DUNE_UNUSED Scalar tmp = BaseFluidState::molarity(phaseIdx, compIdx);
|
|
|
|
|
return 1e100;
|
2011-12-16 13:24:23 +00:00
|
|
|
}
|
2011-12-16 13:24:20 +00:00
|
|
|
|
|
|
|
|
Scalar molarDensity(int phaseIdx) const
|
2011-12-16 13:24:23 +00:00
|
|
|
{
|
|
|
|
|
assert(allowDensity_);
|
|
|
|
|
assert(restrictPhaseIdx_ < 0 || restrictPhaseIdx_ == phaseIdx);
|
2012-09-05 18:08:13 +02:00
|
|
|
DUNE_UNUSED Scalar tmp = BaseFluidState::molarDensity(phaseIdx);
|
|
|
|
|
return 1e100;
|
2011-12-16 13:24:23 +00:00
|
|
|
}
|
2011-12-16 13:24:20 +00:00
|
|
|
|
|
|
|
|
Scalar molarVolume(int phaseIdx) const
|
2011-12-16 13:24:23 +00:00
|
|
|
{
|
2012-02-03 12:34:25 +00:00
|
|
|
assert(allowDensity_);
|
2011-12-16 13:24:23 +00:00
|
|
|
assert(restrictPhaseIdx_ < 0 || restrictPhaseIdx_ == phaseIdx);
|
2012-09-05 18:08:13 +02:00
|
|
|
DUNE_UNUSED Scalar tmp = BaseFluidState::molarVolume(phaseIdx);
|
|
|
|
|
return 1e100;
|
2011-12-16 13:24:23 +00:00
|
|
|
}
|
2011-12-16 13:24:20 +00:00
|
|
|
|
|
|
|
|
Scalar density(int phaseIdx) const
|
2011-12-16 13:24:23 +00:00
|
|
|
{
|
|
|
|
|
assert(allowDensity_);
|
|
|
|
|
assert(restrictPhaseIdx_ < 0 || restrictPhaseIdx_ == phaseIdx);
|
2012-09-05 18:08:13 +02:00
|
|
|
DUNE_UNUSED Scalar tmp = BaseFluidState::density(phaseIdx);
|
|
|
|
|
return 1e100;
|
2011-12-16 13:24:56 +00:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
Scalar saturation(int phaseIdx) const
|
2012-02-03 12:34:25 +00:00
|
|
|
{
|
2011-12-16 13:24:56 +00:00
|
|
|
assert(false);
|
2012-09-05 18:08:13 +02:00
|
|
|
DUNE_UNUSED Scalar tmp = BaseFluidState::saturation(phaseIdx);
|
|
|
|
|
return 1e100;
|
2011-12-16 13:24:56 +00:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
Scalar fugacity(int phaseIdx, int compIdx) const
|
|
|
|
|
{
|
|
|
|
|
assert(false);
|
2012-09-05 18:08:13 +02:00
|
|
|
DUNE_UNUSED Scalar tmp = BaseFluidState::fugacity(phaseIdx, compIdx);
|
|
|
|
|
return 1e100;
|
2011-12-16 13:24:56 +00:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
Scalar fugacityCoefficient(int phaseIdx, int compIdx) const
|
|
|
|
|
{
|
|
|
|
|
assert(false);
|
2012-09-05 18:08:13 +02:00
|
|
|
DUNE_UNUSED Scalar tmp = BaseFluidState::fugacityCoefficient(phaseIdx, compIdx);
|
|
|
|
|
return 1e100;
|
2011-12-16 13:24:56 +00:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
Scalar enthalpy(int phaseIdx) const
|
|
|
|
|
{
|
|
|
|
|
assert(false);
|
2012-09-05 18:08:13 +02:00
|
|
|
DUNE_UNUSED Scalar tmp = BaseFluidState::enthalpy(phaseIdx);
|
|
|
|
|
return 1e100;
|
2011-12-16 13:24:56 +00:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
Scalar internalEnergy(int phaseIdx) const
|
|
|
|
|
{
|
|
|
|
|
assert(false);
|
2012-09-05 18:08:13 +02:00
|
|
|
DUNE_UNUSED Scalar tmp = BaseFluidState::internalEnergy(phaseIdx);
|
|
|
|
|
return 1e100;
|
2011-12-16 13:24:56 +00:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
Scalar viscosity(int phaseIdx) const
|
2012-02-03 12:34:25 +00:00
|
|
|
{
|
2011-12-16 13:24:56 +00:00
|
|
|
assert(false);
|
2012-09-05 18:08:13 +02:00
|
|
|
DUNE_UNUSED Scalar tmp = BaseFluidState::viscosity(phaseIdx);
|
|
|
|
|
return 1e100;
|
2011-12-16 13:24:23 +00:00
|
|
|
}
|
2011-12-16 13:24:20 +00:00
|
|
|
|
|
|
|
|
private:
|
|
|
|
|
bool allowSaturation_;
|
|
|
|
|
bool allowTemperature_;
|
|
|
|
|
bool allowPressure_;
|
|
|
|
|
bool allowComposition_;
|
|
|
|
|
bool allowDensity_;
|
2011-12-16 13:24:23 +00:00
|
|
|
int restrictPhaseIdx_;
|
2011-12-16 13:24:20 +00:00
|
|
|
};
|
|
|
|
|
|
2011-12-16 13:24:56 +00:00
|
|
|
template <class Scalar, class BaseFluidState>
|
|
|
|
|
void checkFluidState(const BaseFluidState &fs)
|
|
|
|
|
{
|
|
|
|
|
// fluid states must be copy-able
|
|
|
|
|
BaseFluidState tmpFs(fs);
|
|
|
|
|
tmpFs = fs;
|
|
|
|
|
|
|
|
|
|
// a fluid state must provide a checkDefined() method
|
|
|
|
|
fs.checkDefined();
|
2012-02-03 12:34:25 +00:00
|
|
|
|
2011-12-16 13:24:56 +00:00
|
|
|
// make sure the fluid state provides all mandatory methods
|
|
|
|
|
while (false) {
|
2012-05-23 15:52:18 +02:00
|
|
|
Scalar DUNE_UNUSED val;
|
2011-12-16 13:24:56 +00:00
|
|
|
|
|
|
|
|
val = fs.temperature(/*phaseIdx=*/0);
|
|
|
|
|
val = fs.pressure(/*phaseIdx=*/0);
|
|
|
|
|
val = fs.moleFraction(/*phaseIdx=*/0, /*compIdx=*/0);
|
|
|
|
|
val = fs.massFraction(/*phaseIdx=*/0, /*compIdx=*/0);
|
|
|
|
|
val = fs.averageMolarMass(/*phaseIdx=*/0);
|
|
|
|
|
val = fs.molarity(/*phaseIdx=*/0, /*compIdx=*/0);
|
|
|
|
|
val = fs.molarDensity(/*phaseIdx=*/0);
|
|
|
|
|
val = fs.molarVolume(/*phaseIdx=*/0);
|
|
|
|
|
val = fs.density(/*phaseIdx=*/0);
|
|
|
|
|
val = fs.saturation(/*phaseIdx=*/0);
|
|
|
|
|
val = fs.fugacity(/*phaseIdx=*/0, /*compIdx=*/0);
|
|
|
|
|
val = fs.fugacityCoefficient(/*phaseIdx=*/0, /*compIdx=*/0);
|
|
|
|
|
val = fs.enthalpy(/*phaseIdx=*/0);
|
|
|
|
|
val = fs.internalEnergy(/*phaseIdx=*/0);
|
|
|
|
|
val = fs.viscosity(/*phaseIdx=*/0);
|
|
|
|
|
};
|
|
|
|
|
}
|
|
|
|
|
|
2012-11-02 13:45:42 +01:00
|
|
|
/*!
|
|
|
|
|
* \brief Checks whether a fluid system adheres to the specification.
|
|
|
|
|
*/
|
2011-12-16 13:24:20 +00:00
|
|
|
template <class Scalar, class FluidSystem>
|
|
|
|
|
void checkFluidSystem()
|
|
|
|
|
{
|
|
|
|
|
std::cout << "Testing fluid system '" << Dune::className<FluidSystem>() << "'\n";
|
|
|
|
|
|
2011-12-16 13:24:23 +00:00
|
|
|
// make sure the fluid system provides the number of phases and
|
|
|
|
|
// the number of components
|
|
|
|
|
enum { numPhases = FluidSystem::numPhases };
|
|
|
|
|
enum { numComponents = FluidSystem::numComponents };
|
2011-12-16 13:24:20 +00:00
|
|
|
|
2011-12-16 13:24:23 +00:00
|
|
|
HairSplittingFluidState<Scalar, FluidSystem> fs;
|
2011-12-16 13:24:20 +00:00
|
|
|
fs.allowTemperature(true);
|
|
|
|
|
fs.allowPressure(true);
|
|
|
|
|
fs.allowComposition(true);
|
2011-12-16 13:24:23 +00:00
|
|
|
fs.restrictToPhase(-1);
|
2011-12-16 13:24:20 +00:00
|
|
|
|
|
|
|
|
// check whether the parameter cache adheres to the API
|
|
|
|
|
typedef typename FluidSystem::ParameterCache PC;
|
|
|
|
|
PC paramCache;
|
|
|
|
|
try { paramCache.updateAll(fs); } catch (...) {};
|
|
|
|
|
try { paramCache.updateAll(fs, /*except=*/PC::None); } catch (...) {};
|
|
|
|
|
try { paramCache.updateAll(fs, /*except=*/PC::Temperature | PC::Pressure | PC::Composition); } catch (...) {};
|
|
|
|
|
try { paramCache.updateAllPressures(fs); } catch (...) {};
|
2011-12-16 13:24:23 +00:00
|
|
|
|
|
|
|
|
for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
|
|
|
|
|
fs.restrictToPhase(phaseIdx);
|
|
|
|
|
try { paramCache.updatePhase(fs, phaseIdx); } catch (...) {};
|
|
|
|
|
try { paramCache.updatePhase(fs, phaseIdx, /*except=*/PC::None); } catch (...) {};
|
|
|
|
|
try { paramCache.updatePhase(fs, phaseIdx, /*except=*/PC::Temperature | PC::Pressure | PC::Composition); } catch (...) {};
|
|
|
|
|
try { paramCache.updateTemperature(fs, phaseIdx); } catch (...) {};
|
2012-02-06 17:47:48 +00:00
|
|
|
try { paramCache.updatePressure(fs, phaseIdx); } catch (...) {};
|
2011-12-16 13:24:23 +00:00
|
|
|
try { paramCache.updateComposition(fs, phaseIdx); } catch (...) {};
|
|
|
|
|
try { paramCache.updateSingleMoleFraction(fs, phaseIdx, /*compIdx=*/0); } catch (...) {};
|
|
|
|
|
}
|
2011-12-16 13:24:20 +00:00
|
|
|
|
|
|
|
|
// some value to make sure the return values of the fluid system
|
|
|
|
|
// are convertible to scalars
|
2012-05-23 15:52:18 +02:00
|
|
|
Scalar DUNE_UNUSED val;
|
2011-12-16 13:24:20 +00:00
|
|
|
|
|
|
|
|
// actually check the fluid system API
|
2012-06-13 16:14:06 +02:00
|
|
|
try { FluidSystem::init(); } catch (...) {};
|
2011-12-16 13:24:20 +00:00
|
|
|
for (int phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
|
2011-12-16 13:24:23 +00:00
|
|
|
fs.restrictToPhase(phaseIdx);
|
2011-12-16 13:24:20 +00:00
|
|
|
fs.allowPressure(FluidSystem::isCompressible(phaseIdx));
|
|
|
|
|
fs.allowComposition(true);
|
|
|
|
|
fs.allowDensity(false);
|
|
|
|
|
try { val = FluidSystem::density(fs, paramCache, phaseIdx); } catch (...) {};
|
2012-10-28 21:45:08 +01:00
|
|
|
|
2011-12-16 13:24:20 +00:00
|
|
|
fs.allowPressure(true);
|
|
|
|
|
fs.allowDensity(true);
|
|
|
|
|
try { val = FluidSystem::viscosity(fs, paramCache, phaseIdx); } catch (...) {};
|
|
|
|
|
try { val = FluidSystem::enthalpy(fs, paramCache, phaseIdx); } catch (...) {};
|
|
|
|
|
try { val = FluidSystem::heatCapacity(fs, paramCache, phaseIdx); } catch (...) {};
|
|
|
|
|
try { val = FluidSystem::thermalConductivity(fs, paramCache, phaseIdx); } catch (...) {};
|
2012-10-28 21:45:08 +01:00
|
|
|
|
2011-12-16 13:24:20 +00:00
|
|
|
for (int compIdx = 0; compIdx < numComponents; ++ compIdx) {
|
2012-01-27 16:17:41 +00:00
|
|
|
fs.allowComposition(!FluidSystem::isIdealMixture(phaseIdx));
|
2011-12-16 13:24:20 +00:00
|
|
|
try { val = FluidSystem::fugacityCoefficient(fs, paramCache, phaseIdx, compIdx); } catch (...) {};
|
2012-01-27 16:17:41 +00:00
|
|
|
fs.allowComposition(true);
|
2011-12-16 13:24:20 +00:00
|
|
|
try { val = FluidSystem::diffusionCoefficient(fs, paramCache, phaseIdx, compIdx); } catch (...) {};
|
|
|
|
|
}
|
2012-10-28 21:45:08 +01:00
|
|
|
|
2011-12-16 13:24:20 +00:00
|
|
|
}
|
|
|
|
|
|
2012-01-27 16:17:41 +00:00
|
|
|
// test for phaseName(), isLiquid() and isIdealGas()
|
2011-12-16 13:24:20 +00:00
|
|
|
for (int phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
|
2012-05-23 15:52:18 +02:00
|
|
|
std::string DUNE_UNUSED name = FluidSystem::phaseName(phaseIdx);
|
|
|
|
|
bool DUNE_UNUSED bVal = FluidSystem::isLiquid(phaseIdx);
|
2012-01-27 16:17:41 +00:00
|
|
|
bVal = FluidSystem::isIdealGas(phaseIdx);
|
2011-12-16 13:24:20 +00:00
|
|
|
}
|
2012-02-03 12:34:25 +00:00
|
|
|
|
2011-12-16 13:24:23 +00:00
|
|
|
// test for componentName()
|
2011-12-16 13:24:20 +00:00
|
|
|
for (int compIdx = 0; compIdx < numComponents; ++ compIdx) {
|
|
|
|
|
val = FluidSystem::molarMass(compIdx);
|
2012-05-23 15:52:18 +02:00
|
|
|
std::string DUNE_UNUSED name = FluidSystem::componentName(compIdx);
|
2011-12-16 13:24:20 +00:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
std::cout << "----------------------------------\n";
|
2012-01-03 12:57:12 +00:00
|
|
|
}
|
2011-12-16 13:24:20 +00:00
|
|
|
|
|
|
|
|
#endif
|