3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
heatpipelaw.hh
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program 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 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program 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 this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
25#ifndef HEATPIPELAW_HH
26#define HEATPIPELAW_HH
27
28#include "heatpipelawparams.hh"
29
31
32#include <algorithm>
33
34#include <math.h>
35#include <assert.h>
36
37namespace Dumux {
38
47template <class ScalarT, class ParamsT = HeatPipeLawParams<ScalarT> >
49{
50public:
51 using Params = ParamsT;
52 using Scalar = typename Params::Scalar;
53
60 static Scalar pc(const Params &params, Scalar Sw)
61 {
62 Scalar Sn = 1 - Sw;
63 Scalar p0Gamma = params.p0()*params.gamma();
64
65 // regularization
66 if (Sn >= 1.0) {
67 Scalar y = p0Gamma*( (1.263*1.0 - 2.120)*1.0 + 1.417)*1.0;
68 Scalar m = p0Gamma*((3*1.263*1.0 - 2*2.120)*1.0 + 1.417);
69 return (Sn - 1)*m + y;
70 }
71 else if (Sn <= 0.0) {
72 Scalar y = 0.0;
73 Scalar m = p0Gamma*1.417;
74 return Sn*m + y;
75 }
76
77 return p0Gamma*((1.263*Sn - 2.120)*Sn + 1.417) * Sn;
78 }
79
87 static Scalar Sw(const Params &params, Scalar pC)
88 {
89 DUNE_THROW(Dune::NotImplemented, "HeatPipeLaw::Sw");
90 }
91
98 static Scalar dpC_dSw(const Params &params, Scalar Sw)
99 {
100 Scalar Sn = 1 - Sw;
101 Scalar p0Gamma = params.p0()*params.gamma();
102 if (Sn > 1.0)
103 Sn = 1.0;
104 else if (Sn <= 0.0) {
105 Scalar m = -p0Gamma*1.417;
106 return m;
107 }
108
109 Scalar m = - p0Gamma*((3*1.263*Sn - 2*2.120)*Sn + 1.417);
110 return m;
111 }
112
119 static Scalar dSw_dpC(const Params &params, Scalar pC)
120 {
121 DUNE_THROW(Dune::NotImplemented, "HeatPipeLaw::dSw_dpC");
122 }
123
130 static Scalar krw(const Params &params, Scalar Sw)
131 {
132 return kr_(Sw);
133 }
134
141 static Scalar krn(const Params &params, Scalar Sw)
142 {
143 Scalar Sn = 1 - Sw;
144 return kr_(Sn);
145 }
146
147private:
148 static Scalar kr_(Scalar S)
149 {
150 const Scalar eps = 0.95;
151 if (S >= 1)
152 return 1;
153 else if (S <= 0)
154 return 0;
155 else if (S > eps) {
156 // regularize
158 Spline sp(eps, 1.0, // x1, x2
159 eps*eps*eps, 1, // y1, y2
160 3*eps*eps, 0); // m1, m2
161 return sp.eval(S);
162 }
163
164 return S*S*S;
165 }
166};
167
168} // end namespace Dumux
169
170#endif
Provides 3rd order polynomial splines.
Specification of the material params for the heat pipe's capillary pressure model.
Definition: adapt.hh:29
A 3rd order polynomial spline.
Definition: spline.hh:55
Implementation of the capillary pressure <-> saturation relation for the heatpipe problem.
Definition: heatpipelaw.hh:49
static Scalar dpC_dSw(const Params &params, Scalar Sw)
Returns the partial derivative of the capillary pressure to the effective saturation.
Definition: heatpipelaw.hh:98
static Scalar krn(const Params &params, Scalar Sw)
The relative permeability for the non-wetting phase.
Definition: heatpipelaw.hh:141
typename Params::Scalar Scalar
Definition: heatpipelaw.hh:52
static Scalar krw(const Params &params, Scalar Sw)
The relative permeability for the wetting phase.
Definition: heatpipelaw.hh:130
ParamsT Params
Definition: heatpipelaw.hh:51
static Scalar pc(const Params &params, Scalar Sw)
The capillary pressure-saturation curve.
Definition: heatpipelaw.hh:60
static Scalar Sw(const Params &params, Scalar pC)
The saturation-capillary pressure curve.
Definition: heatpipelaw.hh:87
static Scalar dSw_dpC(const Params &params, Scalar pC)
Returns the partial derivative of the effective saturation to the capillary pressure.
Definition: heatpipelaw.hh:119