3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
mpadapter.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_ADAPTER_HH
28#define DUMUX_MP_ADAPTER_HH
29
30#include <algorithm>
31#include <cassert>
33
34namespace Dumux {
35
40template <class MaterialLaw, int numPhases>
42{
43 static_assert(AlwaysFalse<MaterialLaw>::value, "Adapter not implemented for the specified number of phases");
44};
45
46template <class MaterialLaw>
47class MPAdapter<MaterialLaw, 2 /*numPhases*/>
48{
49public:
50 using Params = typename MaterialLaw::Params;
58 template <class ContainerT, class FluidState>
59 static void capillaryPressures(ContainerT &values,
60 const Params &params,
61 const FluidState &state,
62 int wPhaseIdx)
63 {
64 assert(values.size() == 2);
65 const int nPhaseIdx = 1 - wPhaseIdx;
66 // non-wetting phase gets the capillary pressure added
67 values[nPhaseIdx] = 0;
68 // wetting phase does not get anything added
69 values[wPhaseIdx] = - MaterialLaw::pc(params, state.saturation(wPhaseIdx));
70 }
71
79 template <class ContainerT, class FluidState>
80 static void relativePermeabilities(ContainerT &values,
81 const Params &params,
82 const FluidState &state,
83 int wPhaseIdx)
84 {
85 assert(values.size() == 2);
86 const int nPhaseIdx = 1 - wPhaseIdx;
87 values[wPhaseIdx] = MaterialLaw::krw(params, state.saturation(wPhaseIdx));
88 values[nPhaseIdx] = MaterialLaw::krn(params, state.saturation(wPhaseIdx));
89 }
90};
91} // end namespace Dumux
92
93#endif
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
Template which always yields a false value.
Definition: typetraits.hh:37
An adapter for mpnc to use the capillary pressure-saturation relationships.
Definition: mpadapter.hh:42
static void capillaryPressures(ContainerT &values, const Params &params, const FluidState &state, int wPhaseIdx)
The capillary pressure-saturation curve.
Definition: mpadapter.hh:59
static void relativePermeabilities(ContainerT &values, const Params &params, const FluidState &state, int wPhaseIdx)
The relative permeability of all phases.
Definition: mpadapter.hh:80
typename MaterialLaw::Params Params
Definition: mpadapter.hh:50