Line | Branch | Exec | Source |
---|---|---|---|
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 | // | ||
7 | /*! | ||
8 | * \file | ||
9 | * \ingroup Fluidmatrixinteractions | ||
10 | * \brief Specification of a function relating volume specific interfacial area to capillary pressure and saturation. | ||
11 | * This parametrization is a second order polynomial. | ||
12 | */ | ||
13 | #ifndef DUMUX_MATERIAL_FLUIDMATRIX_TWO_P_INTERFACIAL_AREA_POLYNOMIAL_SECOND_ORDER_HH | ||
14 | #define DUMUX_MATERIAL_FLUIDMATRIX_TWO_P_INTERFACIAL_AREA_POLYNOMIAL_SECOND_ORDER_HH | ||
15 | |||
16 | #include <cmath> | ||
17 | #include <dune/common/exceptions.hh> | ||
18 | #include <dune/common/float_cmp.hh> | ||
19 | #include <dumux/common/parameters.hh> | ||
20 | |||
21 | namespace Dumux::FluidMatrix { | ||
22 | |||
23 | /*! | ||
24 | * \ingroup Fluidmatrixinteractions | ||
25 | * \brief Implementation of the polynomial of second order relating | ||
26 | * specific interfacial area to wetting phase saturation and capillary pressure as suggested by Joekar-Niasar(2008) \cite joekar2008 . | ||
27 | */ | ||
28 | class InterfacialAreaPolynomialSecondOrder | ||
29 | { | ||
30 | public: | ||
31 | |||
32 | template<class Scalar> | ||
33 | struct Params | ||
34 | { | ||
35 | 1 | Params(Scalar a00 = 0, Scalar a01 = 0, Scalar a02 = 0, Scalar a10 = 0, Scalar a11 = 0, Scalar a20 = 0) | |
36 |
1/2✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
|
1 | : a00_(a00), a01_(a01), a02_(a02), a10_(a10), a11_(a11), a20_(a20) {} |
37 | |||
38 | ✗ | Scalar a00() const { return a00_; } | |
39 | void setA00(Scalar a00) { a00_ = a00; } | ||
40 | |||
41 | ✗ | Scalar a01() const { return a01_; } | |
42 | void setA01(Scalar a01) { a01_ = a01; } | ||
43 | |||
44 | ✗ | Scalar a02() const { return a02_; } | |
45 | void setA02(Scalar a02) { a02_ = a02; } | ||
46 | |||
47 | ✗ | Scalar a10() const { return a10_; } | |
48 | void setA10(Scalar a10) { a10_ = a10; } | ||
49 | |||
50 | ✗ | Scalar a11() const { return a11_; } | |
51 | void setA11(Scalar a11) { a11_ = a11; } | ||
52 | |||
53 | ✗ | Scalar a20() const { return a20_; } | |
54 | void setA20(Scalar a20) { a20_ = a20; } | ||
55 | |||
56 | bool operator== (const Params& p) const | ||
57 | { | ||
58 | return Dune::FloatCmp::eq(a00(), p.a00(), 1e-6) | ||
59 | && Dune::FloatCmp::eq(a01(), p.a01(), 1e-6) | ||
60 | && Dune::FloatCmp::eq(a02(), p.a02(), 1e-6) | ||
61 | && Dune::FloatCmp::eq(a10(), p.a10(), 1e-6) | ||
62 | && Dune::FloatCmp::eq(a11(), p.a11(), 1e-6) | ||
63 | && Dune::FloatCmp::eq(a20(), p.a20(), 1e-6); | ||
64 | } | ||
65 | |||
66 | private: | ||
67 | Scalar a00_, a01_, a02_, a10_, a11_, a20_; | ||
68 | }; | ||
69 | |||
70 | /*! | ||
71 | * \brief Construct from a subgroup from the global parameter tree | ||
72 | * \note This will give you nice error messages if a mandatory parameter is missing | ||
73 | */ | ||
74 | template<class Scalar = double> | ||
75 | static Params<Scalar> makeParams(const std::string& paramGroup) | ||
76 | { | ||
77 | const auto a00 = getParamFromGroup<Scalar>(paramGroup, "A00", 0.0); | ||
78 | const auto a01 = getParamFromGroup<Scalar>(paramGroup, "A01"); | ||
79 | const auto a02 = getParamFromGroup<Scalar>(paramGroup, "A02"); | ||
80 | const auto a10 = getParamFromGroup<Scalar>(paramGroup, "A10"); | ||
81 | const auto a11 = getParamFromGroup<Scalar>(paramGroup, "A11"); | ||
82 | const auto a20 = getParamFromGroup<Scalar>(paramGroup, "A20"); | ||
83 | return {a00, a01, a02, a10, a11, a20}; | ||
84 | } | ||
85 | |||
86 | /*! | ||
87 | * \brief The interfacial area | ||
88 | * | ||
89 | * the suggested (as estimated from pore network models) awn surface: | ||
90 | * \f$[ | ||
91 | a_{wn} = a_{00} + a_{10}S_{w} + a_{20}S_{w}^2 + a_{11} S_{w} p_{c} + a_{01} p_{c} + a_{02}p_{c}^2 | ||
92 | \f$] | ||
93 | * \param params parameter container for the coefficients of the surface | ||
94 | * \param swe Effective saturation of the wetting phase | ||
95 | * \param pc Capillary pressure in \f$\mathrm{[Pa]}\f$ | ||
96 | */ | ||
97 | template<class Scalar> | ||
98 | static Scalar area(const Scalar swe, const Scalar pc, const Params<Scalar>& params) | ||
99 | { | ||
100 | ✗ | const Scalar a00 = params.a00(); | |
101 | ✗ | const Scalar a10 = params.a10(); | |
102 | ✗ | const Scalar a20 = params.a20(); | |
103 | ✗ | const Scalar a11 = params.a11(); | |
104 | ✗ | const Scalar a01 = params.a01(); | |
105 | ✗ | const Scalar a02 = params.a02(); | |
106 | |||
107 | ✗ | return a00 + a10 * swe + a20 * swe*swe + a11*swe*pc + a01*pc + a02*pc*pc; | |
108 | } | ||
109 | |||
110 | /*! | ||
111 | * \brief the derivative of specific interfacial area function w.r.t. capillary pressure | ||
112 | * | ||
113 | * \param params parameter container for the coefficients of the surface | ||
114 | * \param swe Effective saturation of the wetting phase | ||
115 | * \param pc Capillary pressure in \f$\mathrm{[Pa]}\f$ | ||
116 | */ | ||
117 | template<class Scalar> | ||
118 | static Scalar darea_dpc(const Scalar swe, const Scalar pc, const Params<Scalar>& params) | ||
119 | { | ||
120 | return params.a11()*swe + params.a01() + 2.0*params.a02() * pc; | ||
121 | } | ||
122 | |||
123 | /*! | ||
124 | * \brief the derivative of specific interfacial area function w.r.t. saturation | ||
125 | * | ||
126 | * \param params parameter container for the coefficients of the surface | ||
127 | * \param swe Effective saturation of the wetting phase | ||
128 | * \param pc Capillary pressure in \f$\mathrm{[Pa]}\f$ | ||
129 | */ | ||
130 | template<class Scalar> | ||
131 | static Scalar darea_dsw(const Scalar swe, const Scalar pc, const Params<Scalar>& params) | ||
132 | { | ||
133 | return params.a11()*pc + params.a10() + 2.0*params.a20()*swe; | ||
134 | } | ||
135 | }; | ||
136 | |||
137 | } // end namespace Dumux::FluidMatrix | ||
138 | |||
139 | #endif | ||
140 |