3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
mplinearmaterialparams.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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
25#ifndef MP_LINEAR_MATERIAL_PARAMS_HH
26#define MP_LINEAR_MATERIAL_PARAMS_HH
27
28namespace Dumux {
29
35template<int numPhasesV, class ScalarT>
37{
38public:
39 using Scalar = ScalarT;
40 enum { numPhases = numPhasesV };
41
42 [[deprecated("Use FluidMatrix::MPLinearMaterial. Class will be removed after 3.3.")]]
44 {
45 for (int i = 0; i < numPhases; ++i) {
46 setPcMinSat(i, 0.0);
47 setPcMaxSat(i, 0.0);
48 }
49 }
50
58 Scalar sReg(int phaseIdx) const
59 { return 0.10; }
60
65 Scalar pcMinSat(int phaseIdx) const
66 { return pcMinSat_[phaseIdx]; }
67
73 void setPcMinSat(int phaseIdx, Scalar val)
74 { pcMinSat_[phaseIdx] = val; }
75
80 Scalar pcMaxSat(int phaseIdx) const
81 { return pcMaxSat_[phaseIdx]; }
82
88 void setPcMaxSat(int phaseIdx, Scalar val)
89 { pcMaxSat_[phaseIdx] = val; }
90
99 Scalar krLowS(int phaseIdx) const
100 { return 0.05; }
101
110 Scalar krHighS(int phaseIdx) const
111 { return 0.95; }
112
113private:
114 Scalar pcMaxSat_[numPhases];
115 Scalar pcMinSat_[numPhases];
116};
117} // namespace Dumux
118
119#endif
Definition: adapt.hh:29
Reference implementation of params for the linear M-phase material material.
Definition: mplinearmaterialparams.hh:37
Scalar krHighS(int phaseIdx) const
Return the threshold saturation of the respective phase above which the relative permeability gets re...
Definition: mplinearmaterialparams.hh:110
Scalar pcMinSat(int phaseIdx) const
Return the capillary pressure in for a phase at .
Definition: mplinearmaterialparams.hh:65
Scalar krLowS(int phaseIdx) const
Return the threshold saturation respective phase below which the relative permeability gets regulariz...
Definition: mplinearmaterialparams.hh:99
Scalar sReg(int phaseIdx) const
Return the threshold saturation at which the relative permeability starts to get regularized.
Definition: mplinearmaterialparams.hh:58
ScalarT Scalar
Definition: mplinearmaterialparams.hh:39
void setPcMaxSat(int phaseIdx, Scalar val)
Set the capillary pressure in for a phase at .
Definition: mplinearmaterialparams.hh:88
Scalar pcMaxSat(int phaseIdx) const
Return the capillary pressure in for a phase at .
Definition: mplinearmaterialparams.hh:80
void setPcMinSat(int phaseIdx, Scalar val)
Set the capillary pressure in for a phase at .
Definition: mplinearmaterialparams.hh:73
@ numPhases
Definition: mplinearmaterialparams.hh:40
MpLinearMaterialParams()
Definition: mplinearmaterialparams.hh:43