version 3.8
exponential.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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
13#ifndef DUMUX_MATERIAL_FLUIDMATRIX_TWO_P_INTERFACIAL_AREA_EXPONENTIAL
14#define DUMUX_MATERIAL_FLUIDMATRIX_TWO_P_INTERFACIAL_AREA_EXPONENTIAL
15
16#include <cmath>
17#include <dune/common/exceptions.hh>
18#include <dune/common/float_cmp.hh>
20
21namespace Dumux::FluidMatrix {
22
30{
31public:
32
33 template<class Scalar>
34 struct Params
35 {
36 Params(Scalar swr = 0, Scalar a1 = 0, Scalar a2 = 0, Scalar a3 = 0)
37 : swr_(swr), a1_(a1), a2_(a2), a3_(a3) {}
38
39 Scalar swr() const { return swr_; }
40 void setSwr(Scalar swr) { swr_ = swr; }
41
42 Scalar a1() const { return a1_; }
43 void setA1(Scalar a1) { a1_ = a1; }
44
45 Scalar a2() const { return a2_; }
46 void setA2(Scalar a2) { a2_ = a2; }
47
48 Scalar a3() const { return a3_; }
49 void setA3(Scalar a3) { a3_ = a3; }
50
51 bool operator== (const Params& p) const
52 {
53 return Dune::FloatCmp::eq(swr(), p.swr(), 1e-6)
54 && Dune::FloatCmp::eq(a1(), p.a1(), 1e-6)
55 && Dune::FloatCmp::eq(a2(), p.a2(), 1e-6)
56 && Dune::FloatCmp::eq(a3(), p.a3(), 1e-6);
57 }
58
59 private:
60 Scalar swr_, a1_, a2_, a3_;
61 };
62
67 template<class Scalar = double>
68 static Params<Scalar> makeParams(const std::string& paramGroup)
69 {
70 const auto swr = getParamFromGroup<Scalar>(paramGroup, "Swr", 0.0);
71 const auto a1 = getParamFromGroup<Scalar>(paramGroup, "A1");
72 const auto a2 = getParamFromGroup<Scalar>(paramGroup, "A2");
73 const auto a3 = getParamFromGroup<Scalar>(paramGroup, "A3");
74 return {swr, a1, a2, a3};
75 }
76
77
89 template<class Scalar>
90 static Scalar area(const Scalar swe, const Scalar pc, const Params<Scalar>& params)
91 {
92 const Scalar a1 = params.a1();
93 const Scalar a2 = params.a2();
94 const Scalar a3 = params.a3();
95 const Scalar swr = params.swr();
96 using std::exp;
97 return a1 * (swr-swe) * (1-swe) + a2 * (swr-swe) * (1-swe) * exp( a3 * pc) ;
98 }
99
106 template<class Scalar>
107 static Scalar darea_dpc(const Scalar swe, const Scalar pc, const Params<Scalar>& params)
108 {
109 const Scalar a2 = params.a2();
110 const Scalar a3 = params.a3();
111 const Scalar swr = params.swr();
112 using std::exp;
113 return a2 * a3 * (swr-swe) * (1-swe) * exp(a3*pc);
114 }
115
122 template<class Scalar>
123 static Scalar darea_dsw(const Scalar swe, const Scalar pc, const Params<Scalar>& params)
124 {
125 Scalar a1 = params.a1();
126 Scalar a2 = params.a2();
127 Scalar a3 = params.a3();
128 Scalar swr = params.swr();
129 using std::exp;
130 return - a1 *(1+swr-2*swe) - a2 * exp(a3*pc) * (1+swr-2*swe);
131 }
132
133};
134
135} // namespace Dumux::FluidMatrix
136
137#endif
Implementation of the exponential function relating specific interfacial area to wetting phase satura...
Definition: exponential.hh:30
static Scalar darea_dsw(const Scalar swe, const Scalar pc, const Params< Scalar > &params)
the derivative of specific interfacial area function w.r.t. saturation
Definition: exponential.hh:123
static Params< Scalar > makeParams(const std::string &paramGroup)
Construct from a subgroup from the global parameter tree.
Definition: exponential.hh:68
static Scalar darea_dpc(const Scalar swe, const Scalar pc, const Params< Scalar > &params)
the derivative of specific interfacial area function w.r.t. capillary pressure
Definition: exponential.hh:107
static Scalar area(const Scalar swe, const Scalar pc, const Params< Scalar > &params)
The interfacial area.
Definition: exponential.hh:90
Definition: brookscorey.hh:23
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
void setSwr(Scalar swr)
Definition: exponential.hh:40
Scalar a3() const
Definition: exponential.hh:48
Params(Scalar swr=0, Scalar a1=0, Scalar a2=0, Scalar a3=0)
Definition: exponential.hh:36
Scalar a2() const
Definition: exponential.hh:45
void setA1(Scalar a1)
Definition: exponential.hh:43
void setA2(Scalar a2)
Definition: exponential.hh:46
Scalar a1() const
Definition: exponential.hh:42
Scalar swr() const
Definition: exponential.hh:39
bool operator==(const Params &p) const
Definition: exponential.hh:51
void setA3(Scalar a3)
Definition: exponential.hh:49