My Project
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 <stdexcept>
31 
32 #include "NullParameterCache.hpp"
33 
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  throw std::runtime_error("Not implemented: The fluid system '"+Dune::className<Implementation>()+"' does not provide a phaseName() method!");
77  }
78 
84  static bool isLiquid(unsigned /*phaseIdx*/)
85  {
86  throw std::runtime_error("Not implemented: The fluid system '"+Dune::className<Implementation>()+"' does not provide a isLiquid() method!");
87  }
88 
103  static bool isIdealMixture(unsigned /*phaseIdx*/)
104  {
105  throw std::runtime_error("Not implemented: The fluid system '"+Dune::className<Implementation>()+"' does not provide a isIdealMixture() method!");
106  }
107 
117  static bool isCompressible(unsigned /*phaseIdx*/)
118  {
119  throw std::runtime_error("Not implemented: The fluid system '"+Dune::className<Implementation>()+"' does not provide a isCompressible() method!");
120  }
121 
128  static bool isIdealGas(unsigned /*phaseIdx*/)
129  {
130  throw std::runtime_error("Not implemented: The fluid system '"+Dune::className<Implementation>()+"' does not provide a isIdealGas() method!");
131  }
132 
138  static const char* componentName(unsigned /*compIdx*/)
139  {
140  throw std::runtime_error("Not implemented: The fluid system '"+Dune::className<Implementation>()+"' does not provide a componentName() method!");
141  }
142 
148  static Scalar molarMass(unsigned /*compIdx*/)
149  {
150  throw std::runtime_error("Not implemented: The fluid system '"+Dune::className<Implementation>()+"' does not provide a molarMass() method!");
151  }
152 
158  static Scalar acentricFactor(unsigned /*compIdx*/)
159  {
160  throw std::runtime_error("Not implemented: The fluid system '"+Dune::className<Implementation>()+"' does not provide a acentricFactor() method!");
161  }
162 
166  static void init()
167  { }
168 
175  template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCache>
176  static LhsEval density(const FluidState& /*fluidState*/,
177  const ParamCache& /*paramCache*/,
178  unsigned /*phaseIdx*/)
179  {
180  throw std::runtime_error("Not implemented: The fluid system '"+Dune::className<Implementation>()+"' does not provide a density() method!");
181  }
182 
197  template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCache>
198  static LhsEval fugacityCoefficient(const FluidState& /*fluidState*/,
199  ParamCache& /*paramCache*/,
200  unsigned /*phaseIdx*/,
201  unsigned /*compIdx*/)
202  {
203  throw std::runtime_error("Not implemented: The fluid system '"+Dune::className<Implementation>()+"' does not provide a fugacityCoefficient() method!");
204  }
205 
212  template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCache>
213  static LhsEval viscosity(const FluidState& /*fluidState*/,
214  ParamCache& /*paramCache*/,
215  unsigned /*phaseIdx*/)
216  {
217  throw std::runtime_error("Not implemented: The fluid system '"+Dune::className<Implementation>()+"' does not provide a viscosity() method!");
218  }
219 
237  template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCache>
238  static LhsEval diffusionCoefficient(const FluidState& /*fluidState*/,
239  ParamCache& /*paramCache*/,
240  unsigned /*phaseIdx*/,
241  unsigned /*compIdx*/)
242  {
243  throw std::runtime_error("Not implemented: The fluid system '"+Dune::className<Implementation>()+"' does not provide a diffusionCoefficient() method!");
244  }
245 
253  template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCache>
254  static LhsEval enthalpy(const FluidState& /*fluidState*/,
255  ParamCache& /*paramCache*/,
256  unsigned /*phaseIdx*/)
257  {
258  throw std::runtime_error("Not implemented: The fluid system '"+Dune::className<Implementation>()+"' does not provide an enthalpy() method!");
259  }
260 
267  template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCache>
268  static LhsEval thermalConductivity(const FluidState& /*fluidState*/,
269  ParamCache& /*paramCache*/,
270  unsigned /*phaseIdx*/)
271  {
272  throw std::runtime_error("Not implemented: The fluid system '"+Dune::className<Implementation>()+"' does not provide a thermalConductivity() method!");
273  }
274 
281  template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCache>
282  static LhsEval heatCapacity(const FluidState& /*fluidState*/,
283  ParamCache& /*paramCache*/,
284  unsigned /*phaseIdx*/)
285  {
286  throw std::runtime_error("Not implemented: The fluid system '"+Dune::className<Implementation>()+"' does not provide a heatCapacity() method!");
287  }
288 
289 
291  static unsigned phaseIsActive(unsigned /*phaseIdx*/)
292  {
293  return true;
294  }
295 };
296 
297 } // namespace Opm
298 
299 #endif
A parameter cache which does nothing.
The base class for all fluid systems.
Definition: BaseFluidSystem.hpp:44
static LhsEval heatCapacity(const FluidState &, ParamCache &, unsigned)
Specific isobaric heat capacity of a fluid phase [J/kg].
Definition: BaseFluidSystem.hpp:282
static unsigned phaseIsActive(unsigned)
Returns whether a fluid phase is active.
Definition: BaseFluidSystem.hpp:291
static LhsEval enthalpy(const FluidState &, ParamCache &, unsigned)
Given a phase's composition, temperature, pressure and density, calculate its specific enthalpy [J/kg...
Definition: BaseFluidSystem.hpp:254
static const int numPhases
Number of fluid phases in the fluid system.
Definition: BaseFluidSystem.hpp:67
static bool isIdealGas(unsigned)
Returns true if and only if a fluid phase is assumed to be an ideal gas.
Definition: BaseFluidSystem.hpp:128
static bool isCompressible(unsigned)
Returns true if and only if a fluid phase is assumed to be compressible.
Definition: BaseFluidSystem.hpp:117
static Scalar molarMass(unsigned)
Return the molar mass of a component in [kg/mol].
Definition: BaseFluidSystem.hpp:148
static bool isIdealMixture(unsigned)
Returns true if and only if a fluid phase is assumed to be an ideal mixture.
Definition: BaseFluidSystem.hpp:103
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:238
static bool isLiquid(unsigned)
Return whether a phase is liquid.
Definition: BaseFluidSystem.hpp:84
static void init()
Initialize the fluid system's static parameters.
Definition: BaseFluidSystem.hpp:166
ScalarT Scalar
The type used for scalar quantities.
Definition: BaseFluidSystem.hpp:49
static const int numComponents
Number of chemical species in the fluid system.
Definition: BaseFluidSystem.hpp:64
static LhsEval density(const FluidState &, const ParamCache &, unsigned)
Calculate the density [kg/m^3] of a fluid phase.
Definition: BaseFluidSystem.hpp:176
static const char * componentName(unsigned)
Return the human readable name of a component.
Definition: BaseFluidSystem.hpp:138
static Scalar acentricFactor(unsigned)
Return the acetntric factor of a component.
Definition: BaseFluidSystem.hpp:158
static char * phaseName(unsigned)
Return the human readable name of a fluid phase.
Definition: BaseFluidSystem.hpp:74
static LhsEval fugacityCoefficient(const FluidState &, ParamCache &, unsigned, unsigned)
Calculate the fugacity coefficient [Pa] of an individual component in a fluid phase.
Definition: BaseFluidSystem.hpp:198
static LhsEval viscosity(const FluidState &, ParamCache &, unsigned)
Calculate the dynamic viscosity of a fluid phase [Pa*s].
Definition: BaseFluidSystem.hpp:213
static LhsEval thermalConductivity(const FluidState &, ParamCache &, unsigned)
Thermal conductivity of a fluid phase [W/(m K)].
Definition: BaseFluidSystem.hpp:268
The type of the fluid system's parameter cache.
Definition: BaseFluidSystem.hpp:59