00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00027 #ifndef OPM_BASE_FLUID_SYSTEM_HPP
00028 #define OPM_BASE_FLUID_SYSTEM_HPP
00029
00030 #include "NullParameterCache.hpp"
00031
00032 #include <opm/common/Exceptions.hpp>
00033 #include <opm/common/ErrorMacros.hpp>
00034 #include <dune/common/classname.hh>
00035
00036 namespace Opm {
00037
00042 template <class ScalarT, class Implementation>
00043 class BaseFluidSystem
00044 {
00045 public:
00049 typedef ScalarT Scalar;
00050
00058 template <class Evaluation>
00059 struct ParameterCache {
00060 ParameterCache() = delete;
00061 };
00062
00064 static const int numComponents = -1000;
00065
00067 static const int numPhases = -2000;
00068
00074 static char* phaseName(unsigned )
00075 {
00076 OPM_THROW(std::runtime_error,
00077 "Not implemented: The fluid system '" << Dune::className<Implementation>() << "' does not provide a phaseName() method!");
00078 }
00079
00085 static bool isLiquid(unsigned )
00086 {
00087 OPM_THROW(std::runtime_error,
00088 "Not implemented: The fluid system '" << Dune::className<Implementation>() << "' does not provide a isLiquid() method!");
00089 }
00090
00105 static bool isIdealMixture(unsigned )
00106 {
00107 OPM_THROW(std::runtime_error,
00108 "Not implemented: The fluid system '" << Dune::className<Implementation>() << "' does not provide a isIdealMixture() method!");
00109 }
00110
00120 static bool isCompressible(unsigned )
00121 {
00122 OPM_THROW(std::runtime_error,
00123 "Not implemented: The fluid system '" << Dune::className<Implementation>() << "' does not provide a isCompressible() method!");
00124 }
00125
00132 static bool isIdealGas(unsigned )
00133 {
00134 OPM_THROW(std::runtime_error,
00135 "Not implemented: The fluid system '" << Dune::className<Implementation>() << "' does not provide a isIdealGas() method!");
00136 }
00137
00143 static const char* componentName(unsigned )
00144 {
00145 OPM_THROW(std::runtime_error,
00146 "Not implemented: The fluid system '" << Dune::className<Implementation>() << "' does not provide a componentName() method!");
00147 }
00148
00154 static Scalar molarMass(unsigned )
00155 {
00156 OPM_THROW(std::runtime_error,
00157 "Not implemented: The fluid system '" << Dune::className<Implementation>() << "' does not provide a molarMass() method!");
00158 }
00159
00163 static void init()
00164 { }
00165
00172 template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCache>
00173 static LhsEval density(const FluidState& ,
00174 const ParamCache& ,
00175 unsigned )
00176 {
00177 OPM_THROW(std::runtime_error,
00178 "Not implemented: The fluid system '" << Dune::className<Implementation>() << "' does not provide a density() method!");
00179 }
00180
00195 template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCache>
00196 static LhsEval fugacityCoefficient(const FluidState& ,
00197 ParamCache& ,
00198 unsigned ,
00199 unsigned )
00200 {
00201 OPM_THROW(std::runtime_error, "Not implemented: The fluid system '" << Dune::className<Implementation>() << "' does not provide a fugacityCoefficient() method!");
00202 }
00203
00210 template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCache>
00211 static LhsEval viscosity(const FluidState& ,
00212 ParamCache& ,
00213 unsigned )
00214 {
00215 OPM_THROW(std::runtime_error, "Not implemented: The fluid system '" << Dune::className<Implementation>() << "' does not provide a viscosity() method!");
00216 }
00217
00235 template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCache>
00236 static LhsEval diffusionCoefficient(const FluidState& ,
00237 ParamCache& ,
00238 unsigned ,
00239 unsigned )
00240 {
00241 OPM_THROW(std::runtime_error, "Not implemented: The fluid system '" << Dune::className<Implementation>() << "' does not provide a diffusionCoefficient() method!");
00242 }
00243
00251 template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCache>
00252 static LhsEval enthalpy(const FluidState& ,
00253 ParamCache& ,
00254 unsigned )
00255 {
00256 OPM_THROW(std::runtime_error, "Not implemented: The fluid system '" << Dune::className<Implementation>() << "' does not provide an enthalpy() method!");
00257 }
00258
00265 template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCache>
00266 static LhsEval thermalConductivity(const FluidState& ,
00267 ParamCache& ,
00268 unsigned )
00269 {
00270 OPM_THROW(std::runtime_error, "Not implemented: The fluid system '" << Dune::className<Implementation>() << "' does not provide a thermalConductivity() method!");
00271 }
00272
00279 template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCache>
00280 static LhsEval heatCapacity(const FluidState& ,
00281 ParamCache& ,
00282 unsigned )
00283 {
00284 OPM_THROW(std::runtime_error, "Not implemented: The fluid system '" << Dune::className<Implementation>() << "' does not provide a heatCapacity() method!");
00285 }
00286 };
00287
00288 }
00289
00290 #endif