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 uses a maximum value of capillary pressure. | ||
12 | * | ||
13 | * \note It is used for calculating the interfacial area between the nonwetting and wetting phase | ||
14 | * by Nuske 2014 (https://elib.uni-stuttgart.de/handle/11682/614, page 60) \cite nuske2014. | ||
15 | */ | ||
16 | #ifndef DUMUX_MATERIAL_FLUIDMATRIX_TWO_P_INTERFACIAL_AREA_PC_MAX | ||
17 | #define DUMUX_MATERIAL_FLUIDMATRIX_TWO_P_INTERFACIAL_AREA_PC_MAX | ||
18 | |||
19 | #include <cmath> | ||
20 | #include <dune/common/exceptions.hh> | ||
21 | #include <dune/common/float_cmp.hh> | ||
22 | #include <dumux/common/parameters.hh> | ||
23 | |||
24 | namespace Dumux::FluidMatrix { | ||
25 | |||
26 | /*! | ||
27 | * \ingroup Fluidmatrixinteractions | ||
28 | * \brief Implementation of the polynomial of second order relating | ||
29 | * specific interfacial area to wetting phase saturation and capillary pressure. | ||
30 | */ | ||
31 | class InterfacialAreaPcMax | ||
32 | { | ||
33 | public: | ||
34 | |||
35 | template<class Scalar> | ||
36 | struct Params | ||
37 | { | ||
38 | 3 | Params(Scalar pcMax = 0, Scalar a1 = 0, Scalar a2 = 0, Scalar a3 = 0) | |
39 |
1/2✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
|
3 | : pcMax_(pcMax), a1_(a1), a2_(a2), a3_(a3) {} |
40 | |||
41 | ✗ | Scalar pcMax() const { return pcMax_; } | |
42 |
1/2✓ Branch 1 taken 3 times.
✗ Branch 2 not taken.
|
3 | void setPcMax(Scalar pcMax) { pcMax_ = pcMax; } |
43 | |||
44 | ✗ | Scalar a1() const { return a1_; } | |
45 |
1/2✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
|
1 | void setA1(Scalar a1) { a1_ = a1; } |
46 | |||
47 | ✗ | Scalar a2() const { return a2_; } | |
48 |
1/2✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
|
1 | void setA2(Scalar a2) { a2_ = a2; } |
49 | |||
50 | ✗ | Scalar a3() const { return a3_; } | |
51 |
1/2✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
|
1 | void setA3(Scalar a3) { a3_ = a3; } |
52 | |||
53 | bool operator== (const Params& p) const | ||
54 | { | ||
55 | return Dune::FloatCmp::eq(pcMax(), p.pcMax(), 1e-6) | ||
56 | && Dune::FloatCmp::eq(a1(), p.a1(), 1e-6) | ||
57 | && Dune::FloatCmp::eq(a2(), p.a2(), 1e-6) | ||
58 | && Dune::FloatCmp::eq(a3(), p.a3(), 1e-6); | ||
59 | } | ||
60 | |||
61 | private: | ||
62 | Scalar pcMax_, a1_, a2_, a3_; | ||
63 | }; | ||
64 | |||
65 | /*! | ||
66 | * \brief Construct from a subgroup from the global parameter tree | ||
67 | * \note This will give you nice error messages if a mandatory parameter is missing | ||
68 | */ | ||
69 | template<class Scalar = double> | ||
70 | 2 | static Params<Scalar> makeParams(const std::string& paramGroup) | |
71 | { | ||
72 | 2 | const auto pcMax = getParamFromGroup<Scalar>(paramGroup, "PcMax"); | |
73 | 2 | const auto a1 = getParamFromGroup<Scalar>(paramGroup, "A1"); | |
74 | 2 | const auto a2 = getParamFromGroup<Scalar>(paramGroup, "A2"); | |
75 | 2 | const auto a3 = getParamFromGroup<Scalar>(paramGroup, "A3"); | |
76 | 4 | return {pcMax, a1, a2, a3}; | |
77 | } | ||
78 | |||
79 | /*! | ||
80 | * \brief The interfacial area | ||
81 | * | ||
82 | * the suggested (as estimated from pore network models) interfacial area between the nonwetting and wetting phase: | ||
83 | * \f$\mathrm{ | ||
84 | a_{wn} = a_{1} (p_{c { \sf max} } - p_{c}) (1-S_{w}) + a_{2} (p_{c {\sf max} } -p_{c})^2 (1-S_{w}) + a_{3} (p_{c {\sf max} }- p_{c}) (1-S_{w})^2 | ||
85 | }\f$ | ||
86 | * \param swe Effective saturation of the wetting phase | ||
87 | * \param pc Capillary pressure in \f$\mathrm{[Pa]}\f$ | ||
88 | * \param params parameter container for the coefficients of the surface | ||
89 | */ | ||
90 | template<class Scalar> | ||
91 | static Scalar area(const Scalar swe, const Scalar pc, const Params<Scalar>& params) | ||
92 | { | ||
93 | 1180508 | const Scalar a1 = params.a1(); | |
94 | 1180508 | const Scalar a2 = params.a2(); | |
95 | 1180508 | const Scalar a3 = params.a3(); | |
96 | 1180508 | const Scalar pcMax = params.pcMax(); | |
97 | 1180508 | return a1 * (pcMax-pc) * (1.0-swe) + a2*(pcMax-pc)*(pcMax-pc) * (1.-swe) + a3 * (pcMax-pc)*(1-swe)*(1-swe); | |
98 | } | ||
99 | }; | ||
100 | |||
101 | } // namespace Dumux::FluidMatrix | ||
102 | |||
103 | #endif | ||
104 |