SomertonThermalConductionLaw.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_SOMERTON_THERMAL_CONDUCTION_LAW_HPP
28#define OPM_SOMERTON_THERMAL_CONDUCTION_LAW_HPP
29
31
33
36
37#include <algorithm>
38
39namespace Opm {
40
60template <class FluidSystem,
61 class ScalarT,
62 class ParamsT = SomertonThermalConductionLawParams<FluidSystem::numPhases, ScalarT> >
64{
65 enum { numPhases = FluidSystem::numPhases };
66
67public:
68 typedef ParamsT Params;
69 typedef typename Params::Scalar Scalar;
70
89 template <class FluidState, class Evaluation = typename FluidState::Scalar>
90 static Evaluation thermalConductivity(const Params& params,
91 const FluidState& fluidState)
92 {
93 Valgrind::CheckDefined(params.vacuumLambda());
94
95 Evaluation lambda = 0;
96 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
97 Valgrind::CheckDefined(params.fullySaturatedLambda(phaseIdx));
98
99 if (FluidSystem::isLiquid(phaseIdx)) {
100 const auto& sat = decay<Evaluation>(fluidState.saturation(phaseIdx));
101 lambda +=
102 regularizedSqrt_(max(0.0, min(1.0, sat)))
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
115protected:
116 template <class Evaluation>
117 static Evaluation regularizedSqrt_(const Evaluation& x)
118 {
119 typedef Spline<Scalar> Spline;
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 static const Spline sqrtRegSpline(0, xMin, // x0, x1
126 0, sqrtXMin, // y0, y1
127 fPrime0, fPrimeXMin); // m0, m1
128
129 if (x > xMin)
130 return sqrt(x);
131 else if (x <= 0)
132 return fPrime0 * x;
133 else
134 return sqrtRegSpline.eval(x);
135 }
136};
137} // namespace Opm
138
139#endif
A traits class which provides basic mathematical functions for arbitrary scalar floating point values...
Some templates to wrap the valgrind client request macros.
Implements the Somerton law of thermal conductivity in a porous medium.
Definition: SomertonThermalConductionLaw.hpp:64
Params::Scalar Scalar
Definition: SomertonThermalConductionLaw.hpp:69
static Evaluation regularizedSqrt_(const Evaluation &x)
Definition: SomertonThermalConductionLaw.hpp:117
static Evaluation thermalConductivity(const Params &params, const FluidState &fluidState)
Given a fluid state, return the effective thermal conductivity [W/m^2 / (K/m)] of the porous medium.
Definition: SomertonThermalConductionLaw.hpp:90
ParamsT Params
Definition: SomertonThermalConductionLaw.hpp:68
Class implementing cubic splines.
Definition: Spline.hpp:91
Evaluation eval(const Evaluation &x, bool extrapolate=false) const
Evaluate the spline at a given position.
Definition: Spline.hpp:797
bool CheckDefined(const T &value)
Make valgrind complain if any of the memory occupied by an object is undefined.
Definition: Valgrind.hpp:74
Definition: Air_Mesitylene.hpp:34
Evaluation sqrt(const Evaluation &value)
Definition: MathToolbox.hpp:399
ReturnEval_< Evaluation1, Evaluation2 >::type min(const Evaluation1 &arg1, const Evaluation2 &arg2)
Definition: MathToolbox.hpp:346
ReturnEval_< Evaluation1, Evaluation2 >::type max(const Evaluation1 &arg1, const Evaluation2 &arg2)
Definition: MathToolbox.hpp:341