BaseFluidSystem.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_BASE_FLUID_SYSTEM_HPP
28 #define OPM_BASE_FLUID_SYSTEM_HPP
29 
30 #include "NullParameterCache.hpp"
31 
32 #include <opm/common/Exceptions.hpp>
33 #include <opm/common/ErrorMacros.hpp>
34 #include <dune/common/classname.hh>
35 
36 namespace Opm {
37 
42 template <class ScalarT, class Implementation>
44 {
45 public:
49  typedef ScalarT Scalar;
50 
58  template <class Evaluation>
59  struct ParameterCache {
60  ParameterCache() = delete; // derived fluid systems must specify this class!
61  };
62 
64  static const int numComponents = -1000;
65 
67  static const int numPhases = -2000;
68 
74  static char* phaseName(unsigned /*phaseIdx*/)
75  {
76  OPM_THROW(std::runtime_error,
77  "Not implemented: The fluid system '" << Dune::className<Implementation>() << "' does not provide a phaseName() method!");
78  }
79 
85  static bool isLiquid(unsigned /*phaseIdx*/)
86  {
87  OPM_THROW(std::runtime_error,
88  "Not implemented: The fluid system '" << Dune::className<Implementation>() << "' does not provide a isLiquid() method!");
89  }
90 
105  static bool isIdealMixture(unsigned /*phaseIdx*/)
106  {
107  OPM_THROW(std::runtime_error,
108  "Not implemented: The fluid system '" << Dune::className<Implementation>() << "' does not provide a isIdealMixture() method!");
109  }
110 
120  static bool isCompressible(unsigned /*phaseIdx*/)
121  {
122  OPM_THROW(std::runtime_error,
123  "Not implemented: The fluid system '" << Dune::className<Implementation>() << "' does not provide a isCompressible() method!");
124  }
125 
132  static bool isIdealGas(unsigned /*phaseIdx*/)
133  {
134  OPM_THROW(std::runtime_error,
135  "Not implemented: The fluid system '" << Dune::className<Implementation>() << "' does not provide a isIdealGas() method!");
136  }
137 
143  static const char* componentName(unsigned /*compIdx*/)
144  {
145  OPM_THROW(std::runtime_error,
146  "Not implemented: The fluid system '" << Dune::className<Implementation>() << "' does not provide a componentName() method!");
147  }
148 
154  static Scalar molarMass(unsigned /*compIdx*/)
155  {
156  OPM_THROW(std::runtime_error,
157  "Not implemented: The fluid system '" << Dune::className<Implementation>() << "' does not provide a molarMass() method!");
158  }
159 
163  static void init()
164  { }
165 
172  template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCache>
173  static LhsEval density(const FluidState& /*fluidState*/,
174  const ParamCache& /*paramCache*/,
175  unsigned /*phaseIdx*/)
176  {
177  OPM_THROW(std::runtime_error,
178  "Not implemented: The fluid system '" << Dune::className<Implementation>() << "' does not provide a density() method!");
179  }
180 
195  template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCache>
196  static LhsEval fugacityCoefficient(const FluidState& /*fluidState*/,
197  ParamCache& /*paramCache*/,
198  unsigned /*phaseIdx*/,
199  unsigned /*compIdx*/)
200  {
201  OPM_THROW(std::runtime_error, "Not implemented: The fluid system '" << Dune::className<Implementation>() << "' does not provide a fugacityCoefficient() method!");
202  }
203 
210  template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCache>
211  static LhsEval viscosity(const FluidState& /*fluidState*/,
212  ParamCache& /*paramCache*/,
213  unsigned /*phaseIdx*/)
214  {
215  OPM_THROW(std::runtime_error, "Not implemented: The fluid system '" << Dune::className<Implementation>() << "' does not provide a viscosity() method!");
216  }
217 
235  template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCache>
236  static LhsEval diffusionCoefficient(const FluidState& /*fluidState*/,
237  ParamCache& /*paramCache*/,
238  unsigned /*phaseIdx*/,
239  unsigned /*compIdx*/)
240  {
241  OPM_THROW(std::runtime_error, "Not implemented: The fluid system '" << Dune::className<Implementation>() << "' does not provide a diffusionCoefficient() method!");
242  }
243 
251  template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCache>
252  static LhsEval enthalpy(const FluidState& /*fluidState*/,
253  ParamCache& /*paramCache*/,
254  unsigned /*phaseIdx*/)
255  {
256  OPM_THROW(std::runtime_error, "Not implemented: The fluid system '" << Dune::className<Implementation>() << "' does not provide an enthalpy() method!");
257  }
258 
265  template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCache>
266  static LhsEval thermalConductivity(const FluidState& /*fluidState*/,
267  ParamCache& /*paramCache*/,
268  unsigned /*phaseIdx*/)
269  {
270  OPM_THROW(std::runtime_error, "Not implemented: The fluid system '" << Dune::className<Implementation>() << "' does not provide a thermalConductivity() method!");
271  }
272 
279  template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCache>
280  static LhsEval heatCapacity(const FluidState& /*fluidState*/,
281  ParamCache& /*paramCache*/,
282  unsigned /*phaseIdx*/)
283  {
284  OPM_THROW(std::runtime_error, "Not implemented: The fluid system '" << Dune::className<Implementation>() << "' does not provide a heatCapacity() method!");
285  }
286 };
287 
288 } // namespace Opm
289 
290 #endif
static LhsEval diffusionCoefficient(const FluidState &, ParamCache &, unsigned, unsigned)
Calculate the binary molecular diffusion coefficient for a component in a fluid phase [mol^2 * s / (k...
Definition: BaseFluidSystem.hpp:236
static bool isIdealMixture(unsigned)
Returns true if and only if a fluid phase is assumed to be an ideal mixture.
Definition: BaseFluidSystem.hpp:105
static const int numComponents
Number of chemical species in the fluid system.
Definition: BaseFluidSystem.hpp:64
static const int numPhases
Number of fluid phases in the fluid system.
Definition: BaseFluidSystem.hpp:67
static void init()
Initialize the fluid system&#39;s static parameters.
Definition: BaseFluidSystem.hpp:163
Definition: Air_Mesitylene.hpp:33
static bool isIdealGas(unsigned)
Returns true if and only if a fluid phase is assumed to be an ideal gas.
Definition: BaseFluidSystem.hpp:132
static LhsEval thermalConductivity(const FluidState &, ParamCache &, unsigned)
Thermal conductivity of a fluid phase [W/(m K)].
Definition: BaseFluidSystem.hpp:266
A parameter cache which does nothing.
The type of the fluid system&#39;s parameter cache.
Definition: BaseFluidSystem.hpp:59
static LhsEval enthalpy(const FluidState &, ParamCache &, unsigned)
Given a phase&#39;s composition, temperature, pressure and density, calculate its specific enthalpy [J/kg...
Definition: BaseFluidSystem.hpp:252
static LhsEval viscosity(const FluidState &, ParamCache &, unsigned)
Calculate the dynamic viscosity of a fluid phase [Pa*s].
Definition: BaseFluidSystem.hpp:211
static LhsEval density(const FluidState &, const ParamCache &, unsigned)
Calculate the density [kg/m^3] of a fluid phase.
Definition: BaseFluidSystem.hpp:173
The base class for all fluid systems.
Definition: BaseFluidSystem.hpp:43
static Scalar molarMass(unsigned)
Return the molar mass of a component in [kg/mol].
Definition: BaseFluidSystem.hpp:154
static char * phaseName(unsigned)
Return the human readable name of a fluid phase.
Definition: BaseFluidSystem.hpp:74
ScalarT Scalar
The type used for scalar quantities.
Definition: BaseFluidSystem.hpp:49
static const char * componentName(unsigned)
Return the human readable name of a component.
Definition: BaseFluidSystem.hpp:143
static bool isLiquid(unsigned)
Return whether a phase is liquid.
Definition: BaseFluidSystem.hpp:85
static LhsEval heatCapacity(const FluidState &, ParamCache &, unsigned)
Specific isobaric heat capacity of a fluid phase [J/kg].
Definition: BaseFluidSystem.hpp:280
static LhsEval fugacityCoefficient(const FluidState &, ParamCache &, unsigned, unsigned)
Calculate the fugacity coefficient [Pa] of an individual component in a fluid phase.
Definition: BaseFluidSystem.hpp:196
static bool isCompressible(unsigned)
Returns true if and only if a fluid phase is assumed to be compressible.
Definition: BaseFluidSystem.hpp:120