version 3.8
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// 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_CUBIC
14#define DUMUX_MATERIAL_FLUIDMATRIX_TWO_P_INTERFACIAL_AREA_EXPONENTIAL_CUBIC
15
16#include <cmath>
17#include <dune/common/exceptions.hh>
18#include <dune/common/float_cmp.hh>
20
21namespace Dumux::FluidMatrix {
22
29{
30public:
31
32 template<class Scalar>
33 struct Params
34 {
35 Params(Scalar swr = 0, Scalar a1 = 0, Scalar a2 = 0, Scalar a3 = 0)
36 : swr_(swr), a1_(a1), a2_(a2), a3_(a3) {}
37
38 Scalar swr() const { return swr_; }
39 void setSwr(Scalar swr) { swr_ = swr; }
40
41 Scalar a1() const { return a1_; }
42 void setA1(Scalar a1) { a1_ = a1; }
43
44 Scalar a2() const { return a2_; }
45 void setA2(Scalar a2) { a2_ = a2; }
46
47 Scalar a3() const { return a3_; }
48 void setA3(Scalar a3) { a3_ = a3; }
49
50 bool operator== (const Params& p) const
51 {
52 return Dune::FloatCmp::eq(swr(), p.swr(), 1e-6)
53 && Dune::FloatCmp::eq(a1(), p.a1(), 1e-6)
54 && Dune::FloatCmp::eq(a2(), p.a2(), 1e-6)
55 && Dune::FloatCmp::eq(a3(), p.a3(), 1e-6);
56 }
57
58 private:
59 Scalar swr_, a1_, a2_, a3_;
60 };
61
66 template<class Scalar = double>
67 static Params<Scalar> makeParams(const std::string& paramGroup)
68 {
69 const auto swr = getParamFromGroup<Scalar>(paramGroup, "Swr", 0.0);
70 const auto a1 = getParamFromGroup<Scalar>(paramGroup, "A1");
71 const auto a2 = getParamFromGroup<Scalar>(paramGroup, "A2");
72 const auto a3 = getParamFromGroup<Scalar>(paramGroup, "A3");
73 return {swr, a1, a2, a3};
74 }
75
86 template<class Scalar>
87 static Scalar area(const Scalar sw, const Scalar pc, const Params<Scalar>& params)
88 {
89 const Scalar a1 = params.a1();
90 const Scalar a2 = params.a2();
91 const Scalar a3 = params.a3();
92 const Scalar swr = params.swr();
93 const Scalar factor = (swr-sw)*(1.0-sw) ;
94 return a1*factor + a2*factor*pc + a3*factor*pc*pc;
95 }
96
104 template<class Scalar>
105 static Scalar darea_dpc(const Scalar sw, const Scalar pc, const Params<Scalar>& params)
106 {
107 const Scalar swr = params.swr();
108 const Scalar a1 = params.a1();
109 const Scalar a2 = params.a2();
110 const Scalar a3 = params.a3();
111 return a2*(swr-sw)*(1-sw) + a3*(swr-sw)*(1-sw)*pc;
112 }
113
121 template<class Scalar>
122 static Scalar darea_dsw(const Scalar sw, const Scalar pc, const Params<Scalar>& params)
123 {
124 const Scalar swr = params.swr();
125 const Scalar a1 = params.a1();
126 const Scalar a2 = params.a2();
127 const Scalar a3 = params.a3();
128 const Scalar derivativeFactor = (sw-1.0)+(sw-swr);
129 return a1 * derivativeFactor + a2 * derivativeFactor * pc + a3 * derivativeFactor * pc*pc ;
130 }
131};
132
133} // namespace Dumux::FluidMatrix
134
135#endif
Implementation of the polynomial of second order relating specific interfacial area to wetting phase ...
Definition: polynomialedgezero2ndorder.hh:29
static Params< Scalar > makeParams(const std::string &paramGroup)
Construct from a subgroup from the global parameter tree.
Definition: polynomialedgezero2ndorder.hh:67
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:105
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:122
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:87
Definition: brookscorey.hh:23
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
void setA2(Scalar a2)
Definition: polynomialedgezero2ndorder.hh:45
Scalar a1() const
Definition: polynomialedgezero2ndorder.hh:41
void setSwr(Scalar swr)
Definition: polynomialedgezero2ndorder.hh:39
Scalar a3() const
Definition: polynomialedgezero2ndorder.hh:47
Params(Scalar swr=0, Scalar a1=0, Scalar a2=0, Scalar a3=0)
Definition: polynomialedgezero2ndorder.hh:35
bool operator==(const Params &p) const
Definition: polynomialedgezero2ndorder.hh:50
void setA3(Scalar a3)
Definition: polynomialedgezero2ndorder.hh:48
void setA1(Scalar a1)
Definition: polynomialedgezero2ndorder.hh:42
Scalar swr() const
Definition: polynomialedgezero2ndorder.hh:38
Scalar a2() const
Definition: polynomialedgezero2ndorder.hh:44