3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
polynomialedgezero2ndorder.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 DUMUX_MATERIAL_FLUIDMATRIX_TWO_P_INTERFACIAL_AREA_EXPONENTIAL_CUBIC
26#define DUMUX_MATERIAL_FLUIDMATRIX_TWO_P_INTERFACIAL_AREA_EXPONENTIAL_CUBIC
27
28#include <cmath>
29#include <dune/common/exceptions.hh>
30#include <dune/common/float_cmp.hh>
32
33namespace Dumux::FluidMatrix {
34
41{
42public:
43
44 template<class Scalar>
45 struct Params
46 {
47 Params(Scalar swr = 0, Scalar a1 = 0, Scalar a2 = 0, Scalar a3 = 0)
48 : swr_(swr), a1_(a1), a2_(a2), a3_(a3) {}
49
50 Scalar swr() const { return swr_; }
51 void setSwr(Scalar swr) { swr_ = swr; }
52
53 Scalar a1() const { return a1_; }
54 void setA1(Scalar a1) { a1_ = a1; }
55
56 Scalar a2() const { return a2_; }
57 void setA2(Scalar a2) { a2_ = a2; }
58
59 Scalar a3() const { return a3_; }
60 void setA3(Scalar a3) { a3_ = a3; }
61
62 bool operator== (const Params& p) const
63 {
64 return Dune::FloatCmp::eq(swr(), p.swr(), 1e-6)
65 && Dune::FloatCmp::eq(a1(), p.a1(), 1e-6)
66 && Dune::FloatCmp::eq(a2(), p.a2(), 1e-6)
67 && Dune::FloatCmp::eq(a3(), p.a3(), 1e-6);
68 }
69
70 private:
71 Scalar swr_, a1_, a2_, a3_;
72 };
73
78 template<class Scalar = double>
79 static Params<Scalar> makeParams(const std::string& paramGroup)
80 {
81 const auto swr = getParamFromGroup<Scalar>(paramGroup, "Swr", 0.0);
82 const auto a1 = getParamFromGroup<Scalar>(paramGroup, "A1");
83 const auto a2 = getParamFromGroup<Scalar>(paramGroup, "A2");
84 const auto a3 = getParamFromGroup<Scalar>(paramGroup, "A3");
85 return {swr, a1, a2, a3};
86 }
87
98 template<class Scalar>
99 static Scalar area(const Scalar sw, const Scalar pc, const Params<Scalar>& params)
100 {
101 const Scalar a1 = params.a1();
102 const Scalar a2 = params.a2();
103 const Scalar a3 = params.a3();
104 const Scalar swr = params.swr();
105 const Scalar factor = (swr-sw)*(1.0-sw) ;
106 return a1*factor + a2*factor*pc + a3*factor*pc*pc;
107 }
108
116 template<class Scalar>
117 static Scalar darea_dpc(const Scalar sw, const Scalar pc, const Params<Scalar>& params)
118 {
119 const Scalar swr = params.swr();
120 const Scalar a1 = params.a1();
121 const Scalar a2 = params.a2();
122 const Scalar a3 = params.a3();
123 return a2*(swr-sw)*(1-sw) + a3*(swr-sw)*(1-sw)*pc;
124 }
125
133 template<class Scalar>
134 static Scalar darea_dsw(const Scalar sw, const Scalar pc, const Params<Scalar>& params)
135 {
136 const Scalar swr = params.swr();
137 const Scalar a1 = params.a1();
138 const Scalar a2 = params.a2();
139 const Scalar a3 = params.a3();
140 const Scalar derivativeFactor = (sw-1.0)+(sw-swr);
141 return a1 * derivativeFactor + a2 * derivativeFactor * pc + a3 * derivativeFactor * pc*pc ;
142 }
143};
144
145} // namespace Dumux::FluidMatrix
146
147#endif
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Definition: brookscorey.hh:35
Implementation of the polynomial of second order relating specific interfacial area to wetting phase ...
Definition: polynomialedgezero2ndorder.hh:41
static Params< Scalar > makeParams(const std::string &paramGroup)
Construct from a subgroup from the global parameter tree.
Definition: polynomialedgezero2ndorder.hh:79
static Scalar darea_dpc(const Scalar sw, const Scalar pc, const Params< Scalar > &params)
the derivative of specific interfacial area function w.r.t. capillary pressure
Definition: polynomialedgezero2ndorder.hh:117
static Scalar darea_dsw(const Scalar sw, const Scalar pc, const Params< Scalar > &params)
the derivative of specific interfacial area function w.r.t. saturation
Definition: polynomialedgezero2ndorder.hh:134
static Scalar area(const Scalar sw, const Scalar pc, const Params< Scalar > &params)
The interfacial area the suggested (as estimated from pore network models) awn surface: ].
Definition: polynomialedgezero2ndorder.hh:99
void setA2(Scalar a2)
Definition: polynomialedgezero2ndorder.hh:57
Scalar a1() const
Definition: polynomialedgezero2ndorder.hh:53
void setSwr(Scalar swr)
Definition: polynomialedgezero2ndorder.hh:51
Scalar a3() const
Definition: polynomialedgezero2ndorder.hh:59
Params(Scalar swr=0, Scalar a1=0, Scalar a2=0, Scalar a3=0)
Definition: polynomialedgezero2ndorder.hh:47
bool operator==(const Params &p) const
Definition: polynomialedgezero2ndorder.hh:62
void setA3(Scalar a3)
Definition: polynomialedgezero2ndorder.hh:60
void setA1(Scalar a1)
Definition: polynomialedgezero2ndorder.hh:54
Scalar swr() const
Definition: polynomialedgezero2ndorder.hh:50
Scalar a2() const
Definition: polynomialedgezero2ndorder.hh:56