version 3.10-dev
exponentialcubic.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//
16#ifndef DUMUX_MATERIAL_FLUIDMATRIX_TWO_P_INTERFACIAL_AREA_EXPONENTIAL_CUBIC
17#define DUMUX_MATERIAL_FLUIDMATRIX_TWO_P_INTERFACIAL_AREA_EXPONENTIAL_CUBIC
18
19#include <cmath>
20#include <dune/common/exceptions.hh>
21#include <dune/common/float_cmp.hh>
23
24namespace Dumux::FluidMatrix {
25
32{
33public:
34
35 template<class Scalar>
36 struct Params
37 {
38 Params(Scalar a1 = 0, Scalar a2 = 0, Scalar a3 = 0)
39 : a1_(a1), a2_(a2), a3_(a3) {}
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(a1(), p.a1(), 1e-6)
53 && Dune::FloatCmp::eq(a2(), p.a2(), 1e-6)
54 && Dune::FloatCmp::eq(a3(), p.a3(), 1e-6);
55 }
56
57 private:
58 Scalar a1_, a2_, a3_;
59 };
60
65 template<class Scalar = double>
66 static Params<Scalar> makeParams(const std::string& paramGroup)
67 {
68 const auto a1 = getParamFromGroup<Scalar>(paramGroup, "A1");
69 const auto a2 = getParamFromGroup<Scalar>(paramGroup, "A2");
70 const auto a3 = getParamFromGroup<Scalar>(paramGroup, "A3");
71 return {a1, a2, a3};
72 }
73
85 template<class Scalar>
86 static Scalar area(const Scalar swe, const Scalar pc, const Params<Scalar>& params)
87 {
88 // TODO think about awn surface for relative saturation
89 const Scalar a1 = params.a1();
90 const Scalar a2 = params.a2();
91 const Scalar a3 = params.a3();
92
93 using std::exp;
94 return a1 * exp( a2 * swe) + a3 * pc * pc * pc ;
95 }
96};
97
98} // namespace Dumux::FluidMatrix
99
100#endif
Implementation of a exponential function relating specific interfacial area to wetting phase saturati...
Definition: exponentialcubic.hh:32
static Params< Scalar > makeParams(const std::string &paramGroup)
Construct from a subgroup from the global parameter tree.
Definition: exponentialcubic.hh:66
static Scalar area(const Scalar swe, const Scalar pc, const Params< Scalar > &params)
The interfacial area.
Definition: exponentialcubic.hh:86
Definition: brookscorey.hh:23
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
void setA2(Scalar a2)
Definition: exponentialcubic.hh:45
Params(Scalar a1=0, Scalar a2=0, Scalar a3=0)
Definition: exponentialcubic.hh:38
Scalar a1() const
Definition: exponentialcubic.hh:41
void setA1(Scalar a1)
Definition: exponentialcubic.hh:42
bool operator==(const Params &p) const
Definition: exponentialcubic.hh:50
Scalar a2() const
Definition: exponentialcubic.hh:44
Scalar a3() const
Definition: exponentialcubic.hh:47
void setA3(Scalar a3)
Definition: exponentialcubic.hh:48