version 3.10-dev
pcmax.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_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>
23
24namespace Dumux::FluidMatrix {
25
32{
33public:
34
35 template<class Scalar>
36 struct Params
37 {
38 Params(Scalar pcMax = 0, Scalar a1 = 0, Scalar a2 = 0, Scalar a3 = 0)
39 : pcMax_(pcMax), a1_(a1), a2_(a2), a3_(a3) {}
40
41 Scalar pcMax() const { return pcMax_; }
42 void setPcMax(Scalar pcMax) { pcMax_ = pcMax; }
43
44 Scalar a1() const { return a1_; }
45 void setA1(Scalar a1) { a1_ = a1; }
46
47 Scalar a2() const { return a2_; }
48 void setA2(Scalar a2) { a2_ = a2; }
49
50 Scalar a3() const { return a3_; }
51 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
69 template<class Scalar = double>
70 static Params<Scalar> makeParams(const std::string& paramGroup)
71 {
72 const auto pcMax = getParamFromGroup<Scalar>(paramGroup, "PcMax");
73 const auto a1 = getParamFromGroup<Scalar>(paramGroup, "A1");
74 const auto a2 = getParamFromGroup<Scalar>(paramGroup, "A2");
75 const auto a3 = getParamFromGroup<Scalar>(paramGroup, "A3");
76 return {pcMax, a1, a2, a3};
77 }
78
90 template<class Scalar>
91 static Scalar area(const Scalar swe, const Scalar pc, const Params<Scalar>& params)
92 {
93 const Scalar a1 = params.a1();
94 const Scalar a2 = params.a2();
95 const Scalar a3 = params.a3();
96 const Scalar pcMax = params.pcMax();
97 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
Implementation of the polynomial of second order relating specific interfacial area to wetting phase ...
Definition: pcmax.hh:32
static Scalar area(const Scalar swe, const Scalar pc, const Params< Scalar > &params)
The interfacial area.
Definition: pcmax.hh:91
static Params< Scalar > makeParams(const std::string &paramGroup)
Construct from a subgroup from the global parameter tree.
Definition: pcmax.hh:70
Definition: brookscorey.hh:23
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Scalar a2() const
Definition: pcmax.hh:47
void setA1(Scalar a1)
Definition: pcmax.hh:45
Params(Scalar pcMax=0, Scalar a1=0, Scalar a2=0, Scalar a3=0)
Definition: pcmax.hh:38
void setA3(Scalar a3)
Definition: pcmax.hh:51
Scalar a3() const
Definition: pcmax.hh:50
bool operator==(const Params &p) const
Definition: pcmax.hh:53
Scalar a1() const
Definition: pcmax.hh:44
void setA2(Scalar a2)
Definition: pcmax.hh:48
void setPcMax(Scalar pcMax)
Definition: pcmax.hh:42
Scalar pcMax() const
Definition: pcmax.hh:41