version 3.10-dev
interfacialarea.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//
12#ifndef DUMUX_MATERIAL_FLUIDMATRIX_TWO_P_INTERFACIAL_AREA_HH
13#define DUMUX_MATERIAL_FLUIDMATRIX_TWO_P_INTERFACIAL_AREA_HH
14
19
20namespace Dumux::FluidMatrix {
21
38template<class ScalarType,
39 class BaseLaw,
40 template<class> class InterfaceType,
41 class Regularization = NoRegularization,
42 class EffToAbsPolicy = TwoPEffToAbsDefaultPolicy>
43class InterfacialArea : public Adapter<InterfacialArea<ScalarType, BaseLaw, InterfaceType, Regularization, EffToAbsPolicy>, InterfaceType>
44{
45public:
46
47 using Scalar = ScalarType;
48
49 using BasicParams = typename BaseLaw::template Params<Scalar>;
50 using EffToAbsParams = typename EffToAbsPolicy::template Params<Scalar>;
51 using RegularizationParams = typename Regularization::template Params<Scalar>;
52
53 using EffToAbs = EffToAbsPolicy;
54
58 static constexpr bool isRegularized()
59 { return !std::is_same<Regularization, NoRegularization>::value; }
60
65 InterfacialArea() = delete;
66
71 explicit InterfacialArea(const std::string& paramGroup)
72 : basicParams_(makeBasicParams(paramGroup))
73 , effToAbsParams_(makeEffToAbsParams(paramGroup))
74 {
75 if constexpr (isRegularized())
76 regularization_.init(this, paramGroup);
77 }
78
83 InterfacialArea(const BasicParams& baseParams,
85 const RegularizationParams& regParams = {})
86 : basicParams_(baseParams)
87 , effToAbsParams_(effToAbsParams)
88 {
89 if constexpr (isRegularized())
90 regularization_.init(this, baseParams, effToAbsParams, regParams);
91 }
92
96 template<bool enableRegularization = isRegularized()>
97 Scalar area(const Scalar sw, const Scalar pc) const
98 {
99 const auto swe = EffToAbs::swToSwe(sw, effToAbsParams_);
100 if constexpr (enableRegularization)
101 {
102 const auto regularized = regularization_.area(swe, pc);
103 if (regularized)
104 return regularized.value();
105 }
106
107 return BaseLaw::area(swe, pc, basicParams_);
108 }
109
113 template<bool enableRegularization = isRegularized()>
114 Scalar darea_dpc(const Scalar sw, const Scalar pc) const
115 {
116 if constexpr (enableRegularization)
117 {
118 const auto regularized = regularization_.darea_dpc(sw, pc);
119 if (regularized)
120 return regularized.value();
121 }
122
123 return BaseLaw::darea_dpc(sw, pc, basicParams_);
124 }
125
129 template<bool enableRegularization = isRegularized()>
130 Scalar darea_dsw(const Scalar sw, const Scalar pc) const
131 {
132 if constexpr (enableRegularization)
133 {
134 const auto regularized = regularization_.darea_dsw(sw, pc);
135 if (regularized)
136 return regularized.value();
137 }
138
139 return BaseLaw::darea_dsw(sw, pc, basicParams_);
140 }
141
145 bool operator== (const InterfacialArea& o) const
146 {
147 return basicParams_ == o.basicParams_
148 && effToAbsParams_ == o.effToAbsParams_
149 && regularization_ == o.regularization_;
150 }
151
156 static BasicParams makeBasicParams(const std::string& paramGroup)
157 {
158 return BaseLaw::template makeParams<Scalar>(paramGroup);
159 }
160
165 { return basicParams_; }
166
171 static EffToAbsParams makeEffToAbsParams(const std::string& paramGroup)
172 {
173 return EffToAbs::template makeParams<Scalar>(paramGroup);
174 }
175
180 { return effToAbsParams_; }
181
182private:
183 BasicParams basicParams_;
184 EffToAbsParams effToAbsParams_;
185 Regularization regularization_;
186};
187
188} // end namespace Dumux::FluidMatrix
189
190#endif
Wrapper class to implement regularized laws (pc-sw-a) with a conversion policy between absolution and...
Definition: interfacialarea.hh:44
EffToAbsPolicy EffToAbs
Definition: interfacialarea.hh:53
const BasicParams & basicParams() const
Return the base law's parameters.
Definition: interfacialarea.hh:164
bool operator==(const InterfacialArea &o) const
Equality comparison with another instance.
Definition: interfacialarea.hh:145
ScalarType Scalar
Definition: interfacialarea.hh:47
Scalar darea_dpc(const Scalar sw, const Scalar pc) const
The partial derivative of the capillary pressure w.r.t. the saturation.
Definition: interfacialarea.hh:114
const EffToAbsParams & effToAbsParams() const
Return the parameters of the EffToAbs policy.
Definition: interfacialarea.hh:179
typename EffToAbsPolicy::template Params< Scalar > EffToAbsParams
Definition: interfacialarea.hh:50
InterfacialArea()=delete
Deleted default constructor (so we are never in an undefined state)
Scalar area(const Scalar sw, const Scalar pc) const
The capillary pressure-saturation curve.
Definition: interfacialarea.hh:97
Scalar darea_dsw(const Scalar sw, const Scalar pc) const
The partial derivative of the saturation to the capillary pressure.
Definition: interfacialarea.hh:130
static EffToAbsParams makeEffToAbsParams(const std::string &paramGroup)
Create the parameters of the EffToAbs policy using input file parameters.
Definition: interfacialarea.hh:171
InterfacialArea(const BasicParams &baseParams, const EffToAbsParams &effToAbsParams={}, const RegularizationParams &regParams={})
Construct from parameter structs.
Definition: interfacialarea.hh:83
static constexpr bool isRegularized()
Return whether this law is regularized.
Definition: interfacialarea.hh:58
static BasicParams makeBasicParams(const std::string &paramGroup)
Create the base law's parameters using input file parameters.
Definition: interfacialarea.hh:156
InterfacialArea(const std::string &paramGroup)
Construct from a subgroup from the global parameter tree.
Definition: interfacialarea.hh:71
typename Regularization::template Params< Scalar > RegularizationParams
Definition: interfacialarea.hh:51
typename BaseLaw::template Params< Scalar > BasicParams
Definition: interfacialarea.hh:49
This is a policy for 2p material laws how to convert absolute to relative saturations and vice versa.
Wrapper type to combine an arbitrary number of different laws for fluid-matrix interaction (e....
Definition: brookscorey.hh:23
A tag to turn off regularization and it's overhead.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Adapter to inherit from, allowing the inheriting class to be wrapped by the makeFluidMatrixInteractio...
Definition: fluidmatrixinteraction.hh:55