3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
mplinearmaterial.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 *****************************************************************************/
27#ifndef DUMUX_MP_LINEAR_MATERIAL_HH
28#define DUMUX_MP_LINEAR_MATERIAL_HH
29
31
32#include <algorithm>
33
34namespace Dumux {
35
45template <int numPhasesV, class ScalarT, class ParamsT = MpLinearMaterialParams<numPhasesV, ScalarT> >
47{
48public:
49 using Params = ParamsT;
50 using Scalar = typename Params::Scalar;
51 enum { numPhases = numPhasesV };
52
66 template <class ContainerT, class FluidState>
67 static void capillaryPressures(ContainerT &values,
68 const Params &params,
69 const FluidState &state,
70 int wPhaseIdx = 0)
71 {
72 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
73 Scalar S = state.saturation(phaseIdx);
74 values[phaseIdx] =
75 S*params.pcMaxSat(phaseIdx) +
76 (1 - S)*params.pcMinSat(phaseIdx);
77 }
78 }
79
87 template <class ContainerT, class FluidState>
88 static void relativePermeabilities(ContainerT &values,
89 const Params &params,
90 const FluidState &state,
91 int wPhaseIdx = 0)
92 {
93 using std::max;
94 using std::min;
95 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
96 values[phaseIdx] = max(min(state.saturation(phaseIdx),1.0),0.0);
97 }
98};
99} // end namespace Dumux
100
101#endif
Reference implementation of parameters for the M-phase linear material material.
Definition: adapt.hh:29
Implements a linear saturation-capillary pressure relation.
Definition: mplinearmaterial.hh:47
@ numPhases
Definition: mplinearmaterial.hh:51
typename Params::Scalar Scalar
Definition: mplinearmaterial.hh:50
ParamsT Params
Definition: mplinearmaterial.hh:49
static void capillaryPressures(ContainerT &values, const Params &params, const FluidState &state, int wPhaseIdx=0)
The linear capillary pressure-saturation curve.
Definition: mplinearmaterial.hh:67
static void relativePermeabilities(ContainerT &values, const Params &params, const FluidState &state, int wPhaseIdx=0)
The relative permeability of all phases.
Definition: mplinearmaterial.hh:88