Somerton.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  Copyright (C) 2008-2013 by Andreas Lauser
5 
6  This file is part of the Open Porous Media project (OPM).
7 
8  OPM is free software: you can redistribute it and/or modify
9  it under the terms of the GNU General Public License as published by
10  the Free Software Foundation, either version 2 of the License, or
11  (at your option) any later version.
12 
13  OPM is distributed in the hope that it will be useful,
14  but WITHOUT ANY WARRANTY; without even the implied warranty of
15  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16  GNU General Public License for more details.
17 
18  You should have received a copy of the GNU General Public License
19  along with OPM. If not, see <http://www.gnu.org/licenses/>.
20 */
25 #ifndef OPM_SOMERTON_HPP
26 #define OPM_SOMERTON_HPP
27 
28 #include "SomertonParams.hpp"
29 
31 
34 
35 #include <algorithm>
36 
37 namespace Opm
38 {
58 template <class FluidSystem,
59  class ScalarT,
60  class ParamsT = SomertonParams<FluidSystem::numPhases, ScalarT> >
61 class Somerton
62 {
63  enum { numPhases = FluidSystem::numPhases };
64 
65 public:
66  typedef ParamsT Params;
67  typedef typename Params::Scalar Scalar;
68 
87  template <class FluidState, class Evaluation = Scalar>
88  static Evaluation heatConductivity(const Params &params,
89  const FluidState &fluidState)
90  {
91  typedef Opm::MathToolbox<Evaluation> Toolbox;
92 
93  Valgrind::CheckDefined(params.vacuumLambda());
94 
95  Evaluation lambda = 0;
96  for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
97  Valgrind::CheckDefined(params.fullySaturatedLambda(phaseIdx));
98 
99  if (FluidSystem::isLiquid(phaseIdx)) {
100  const auto& sat = Toolbox::template toLhs<Evaluation>(fluidState.saturation(phaseIdx));
101  lambda +=
103  * (params.fullySaturatedLambda(phaseIdx) - params.vacuumLambda());
104  }
105  else { // gas phase
106  lambda += params.fullySaturatedLambda(phaseIdx) - params.vacuumLambda();
107  }
108  };
109 
110  lambda += params.vacuumLambda();
111  assert(lambda >= 0);
112  return lambda;
113  }
114 
115 protected:
116  template <class Evaluation>
117  static Evaluation regularizedSqrt_(const Evaluation& x)
118  {
119  typedef Opm::MathToolbox<Evaluation> Toolbox;
120 
121  static const Scalar xMin = 1e-2;
122  static const Scalar sqrtXMin = std::sqrt(xMin);
123  static const Scalar fPrimeXMin = 1.0/(2*std::sqrt(xMin));
124  static const Scalar fPrime0 = 2*fPrimeXMin;
125  typedef Opm::Spline<Scalar> Spline;
126  static const Spline sqrtRegSpline(0, xMin, // x0, x1
127  0, sqrtXMin, // y0, y1
128  fPrime0, fPrimeXMin); // m0, m1
129 
130  if (x > xMin)
131  return Toolbox::sqrt(x);
132  else if (x <= 0)
133  return fPrime0 * x;
134  else
135  return sqrtRegSpline.eval(x);
136  }
137 };
138 } // namespace Opm
139 
140 #endif
bool CheckDefined(const T &value OPM_UNUSED)
Make valgrind complain if any of the memory occupied by an object is undefined.
Definition: Valgrind.hpp:74
Definition: MathToolbox.hpp:39
Definition: Air_Mesitylene.hpp:31
Params::Scalar Scalar
Definition: Somerton.hpp:67
Some templates to wrap the valgrind client request macros.
Evaluation< Scalar, VarSetTag, numVars > max(const Evaluation< Scalar, VarSetTag, numVars > &x1, const Evaluation< Scalar, VarSetTag, numVars > &x2)
Definition: Math.hpp:114
Evaluation< Scalar, VarSetTag, numVars > sqrt(const Evaluation< Scalar, VarSetTag, numVars > &x)
Definition: Math.hpp:278
Class implementing cubic splines.
ParamsT Params
Definition: Somerton.hpp:66
Class implementing cubic splines.
Definition: Spline.hpp:89
Evaluation< Scalar, VarSetTag, numVars > min(const Evaluation< Scalar, VarSetTag, numVars > &x1, const Evaluation< Scalar, VarSetTag, numVars > &x2)
Definition: Math.hpp:61
Implements the Somerton law of heat conductivity in a porous medium.
Definition: Somerton.hpp:61
Scalar eval(Scalar x, bool extrapolate=false) const
Evaluate the spline at a given position.
Definition: Spline.hpp:809
The default implementation of a parameter object for the Somerton heatconduction law.
static Evaluation heatConductivity(const Params &params, const FluidState &fluidState)
Given a fluid state, return the effective heat conductivity [W/m^2 / (K/m)] of the porous medium...
Definition: Somerton.hpp:88
static Evaluation regularizedSqrt_(const Evaluation &x)
Definition: Somerton.hpp:117
A traits class which provides basic mathematical functions for arbitrary scalar floating point values...