version 3.8
singleshapelocalrules.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_PNM_2P_SINGLE_SHAPE_LOCAL_RULES_HH
14#define DUMUX_PNM_2P_SINGLE_SHAPE_LOCAL_RULES_HH
15
20
22
28template<class ScalarType,
29 class BaseLaw,
30 class Regularization = Dumux::FluidMatrix::NoRegularization>
31class SingleShapeTwoPLocalRules : public Dumux::FluidMatrix::Adapter<SingleShapeTwoPLocalRules<ScalarType, BaseLaw, Regularization>, Dumux::FluidMatrix::PcKrSw>
32{
33public:
34
35 using Scalar = ScalarType;
36
37 using BasicParams = typename BaseLaw::template Params<Scalar>;
38 using RegularizationParams = typename Regularization::template Params<Scalar>;
39
40 static constexpr bool supportsMultipleGeometries()
41 { return false; }
42
43 template<class SpatialParams, class Element, class SubControlVolume, class ElemSol>
44 void updateParams(const SpatialParams& spatialParams,
45 const Element& element,
46 const SubControlVolume& scv,
47 const ElemSol& elemSol)
48 {
49 basicParams_.update(spatialParams, element, scv, elemSol);
50 regularization_.updateParams(spatialParams, element, scv, elemSol);
51 }
52
56 static constexpr int numFluidPhases()
57 { return 2; }
58
62 static constexpr bool isRegularized()
63 { return !std::is_same<Regularization, Dumux::FluidMatrix::NoRegularization>::value; }
64
69 SingleShapeTwoPLocalRules(const BasicParams& baseParams = {},
70 const RegularizationParams& regParams = {},
71 const std::string& paramGroup = "")
72 : basicParams_(baseParams)
73
74 {
75 if constexpr (isRegularized())
76 regularization_.init(this, regParams, paramGroup);
77 }
78
82 template<bool enableRegularization = isRegularized()>
83 Scalar pc(const Scalar sw) const
84 {
85 if constexpr (enableRegularization)
86 {
87 const auto regularized = regularization_.pc(sw);
88 if (regularized)
89 return regularized.value();
90 }
91
92 return BaseLaw::pc(sw, basicParams_);
93 }
94
98 template<bool enableRegularization = isRegularized()>
99 Scalar dpc_dsw(const Scalar sw) const
100 {
101 if constexpr (enableRegularization)
102 {
103 const auto regularized = regularization_.dpc_dsw(sw);
104 if (regularized)
105 return regularized.value();
106 }
107
108 return BaseLaw::dpc_dsw(sw, basicParams_);
109 }
110
114 template<bool enableRegularization = isRegularized()>
115 Scalar sw(const Scalar pc) const
116 {
117 if constexpr (enableRegularization)
118 {
119 const auto regularized = regularization_.sw(pc);
120 if (regularized)
121 return regularized.value();
122 }
123
124 return BaseLaw::sw(pc, basicParams_);
125 }
126
130 template<bool enableRegularization = isRegularized()>
131 Scalar dsw_dpc(const Scalar pc) const
132 {
133 if constexpr (enableRegularization)
134 {
135 const auto regularized = regularization_.dsw_dpc(pc);
136 if (regularized)
137 return regularized.value();
138 }
139
140 return BaseLaw::dsw_dpc(pc, basicParams_);
141 }
142
147 template<bool enableRegularization = isRegularized()>
148 Scalar krw(const Scalar sw) const
149 { return 1.0; }
150
155 template<bool enableRegularization = isRegularized()>
157 { return 0; }
158
163 template<bool enableRegularization = isRegularized()>
164 Scalar krn(const Scalar sw) const
165 { return 1.0; }
166
171 template<bool enableRegularization = isRegularized()>
173 { return 0.0; }
174
179 {
180 return basicParams_ == o.basicParams_
181 && regularization_ == o.regularization_;
182 }
183
188 { return basicParams_; }
189
190
191private:
192 BasicParams basicParams_;
193 Regularization regularization_;
194};
195
196} // end namespace Dumux::PoreNetwork::FluidMatrix
197
198#endif
Base class for all standard pore-local pc-Sw curves.
Definition: singleshapelocalrules.hh:32
void updateParams(const SpatialParams &spatialParams, const Element &element, const SubControlVolume &scv, const ElemSol &elemSol)
Definition: singleshapelocalrules.hh:44
Scalar dpc_dsw(const Scalar sw) const
The partial derivative of the capillary pressure w.r.t. the saturation.
Definition: singleshapelocalrules.hh:99
Scalar krw(const Scalar sw) const
The relative permeability for the wetting phase.
Definition: singleshapelocalrules.hh:148
static constexpr bool supportsMultipleGeometries()
Definition: singleshapelocalrules.hh:40
ScalarType Scalar
Definition: singleshapelocalrules.hh:35
Scalar dsw_dpc(const Scalar pc) const
The partial derivative of the saturation to the capillary pressure.
Definition: singleshapelocalrules.hh:131
Scalar sw(const Scalar pc) const
The saturation-capillary-pressure curve.
Definition: singleshapelocalrules.hh:115
Scalar dkrn_dsw(const Scalar sw) const
The derivative of the relative permeability for the non-wetting phase w.r.t. saturation.
Definition: singleshapelocalrules.hh:172
Scalar pc(const Scalar sw) const
The capillary pressure-saturation curve.
Definition: singleshapelocalrules.hh:83
Scalar dkrw_dsw(const Scalar sw) const
The derivative of the relative permeability for the wetting phase w.r.t. saturation.
Definition: singleshapelocalrules.hh:156
typename BaseLaw::template Params< Scalar > BasicParams
Definition: singleshapelocalrules.hh:37
typename Regularization::template Params< Scalar > RegularizationParams
Definition: singleshapelocalrules.hh:38
static constexpr bool isRegularized()
Return whether this law is regularized.
Definition: singleshapelocalrules.hh:62
const BasicParams & basicParams() const
Return the base law's parameters.
Definition: singleshapelocalrules.hh:187
bool operator==(const SingleShapeTwoPLocalRules &o) const
Equality comparison with another instance.
Definition: singleshapelocalrules.hh:178
SingleShapeTwoPLocalRules(const BasicParams &baseParams={}, const RegularizationParams &regParams={}, const std::string &paramGroup="")
Construct from parameter structs.
Definition: singleshapelocalrules.hh:69
Scalar krn(const Scalar sw) const
The relative permeability for the non-wetting phase.
Definition: singleshapelocalrules.hh:164
static constexpr int numFluidPhases()
Return the number of fluid phases.
Definition: singleshapelocalrules.hh:56
The base class for spatial parameters for pore-network models.
Definition: porenetwork/common/spatialparams.hh:33
Wrapper type to combine an arbitrary number of different laws for fluid-matrix interaction (e....
Definition: localrulesforplatonicbody.hh:28
A tag to turn off regularization and it's overhead.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
This file contains functions related to calculate pore-body properties.
Adapter to inherit from, allowing the inheriting class to be wrapped by the makeFluidMatrixInteractio...
Definition: fluidmatrixinteraction.hh:55
A tag to turn off regularization and it's overhead.
Definition: noregularization.hh:22