3.3.0
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// remove from here after release 3.3 /////////////
31
32#include <algorithm>
33#include <cassert>
35
36namespace Dumux {
37
38
43template <class MaterialLaw, int numFluidPhases>
45{
46 static_assert(AlwaysFalse<MaterialLaw>::value, "Adapter not implemented for the specified number of phases");
47};
48
49template <class MaterialLaw>
50class [[deprecated("Use new material laws and FluidMatrix::MPAdapter instead!")]] MPAdapter<MaterialLaw, 2 /*numFluidPhases*/>
51{
52public:
53 using Params = typename MaterialLaw::Params;
61 template <class ContainerT, class FluidState>
62 static void capillaryPressures(ContainerT &values,
63 const Params &params,
64 const FluidState &state,
65 int wPhaseIdx)
66 {
67 assert(values.size() == 2);
68 const int nPhaseIdx = 1 - wPhaseIdx;
69 // nonwetting phase gets the capillary pressure added
70 values[nPhaseIdx] = 0;
71 // wetting phase does not get anything added
72 values[wPhaseIdx] = - MaterialLaw::pc(params, state.saturation(wPhaseIdx));
73 }
74
82 template <class ContainerT, class FluidState>
83 static void relativePermeabilities(ContainerT &values,
84 const Params &params,
85 const FluidState &state,
86 int wPhaseIdx)
87 {
88 assert(values.size() == 2);
89 const int nPhaseIdx = 1 - wPhaseIdx;
90 values[wPhaseIdx] = MaterialLaw::krw(params, state.saturation(wPhaseIdx));
91 values[nPhaseIdx] = MaterialLaw::krn(params, state.saturation(wPhaseIdx));
92 }
93};
94
95
96} // end namespace Dumux
97
98// remove until here after release 3.3 /////////////
99
100#include <dune/common/fvector.hh>
103
104namespace Dumux::FluidMatrix {
105
110template <class MaterialLaw, int numFluidPhases = MaterialLaw::numFluidPhases()>
112{
113 static_assert(AlwaysFalse<MaterialLaw>::value, "Adapter not implemented for the specified number of phases");
114};
115
116
117template<class MaterialLaw>
118class MPAdapter<MaterialLaw, 2>
119: public Adapter<MPAdapter<MaterialLaw, 2>, MultiPhasePcKrSw>
120{
121public:
122 using Scalar = typename MaterialLaw::Scalar;
123
124 MPAdapter(const MaterialLaw& pcKrS)
125 : pcKrS_(pcKrS)
126 {}
127
133 template <class FluidState>
134 auto capillaryPressures(const FluidState& state, int wPhaseIdx) const
135 {
136 Dune::FieldVector<typename FluidState::Scalar, 2> values;
137
138 const int nPhaseIdx = 1 - wPhaseIdx;
139 // non-wetting phase gets the capillary pressure added
140 values[nPhaseIdx] = 0;
141 // wetting phase does not get anything added
142 values[wPhaseIdx] = - pcKrS_.pc(state.saturation(wPhaseIdx));
143
144 return values;
145 }
146
152 template <class FluidState>
153 auto relativePermeabilities(const FluidState& state, int wPhaseIdx) const
154 {
155 Dune::FieldVector<typename FluidState::Scalar, 2> values;
156
157 const int nPhaseIdx = 1 - wPhaseIdx;
158 values[wPhaseIdx] = pcKrS_.krw(state.saturation(wPhaseIdx));
159 values[nPhaseIdx] = pcKrS_.krn(state.saturation(wPhaseIdx));
160
161 return values;
162 }
163private:
164 const MaterialLaw& pcKrS_;
165};
166
167} // end namespace Dumux::FluidMatrix
168
169
170#endif
Wrapper type to combine an arbitrary number of different laws for fluid-matrix interaction (e....
Definition: adapt.hh:29
Definition: brookscorey.hh:286
Template which always yields a false value.
Definition: typetraits.hh:36
Adapter to inherit from, allowing the inheriting class to be wrapped by the makeFluidMatrixInteractio...
Definition: fluidmatrixinteraction.hh:67
An adapter for mpnc to use the capillary pressure-saturation relationships.
Definition: mpadapter.hh:45
static void capillaryPressures(ContainerT &values, const Params &params, const FluidState &state, int wPhaseIdx)
The capillary pressure-saturation curve.
Definition: mpadapter.hh:62
static void relativePermeabilities(ContainerT &values, const Params &params, const FluidState &state, int wPhaseIdx)
The relative permeability of all phases.
Definition: mpadapter.hh:83
typename MaterialLaw::Params Params
Definition: mpadapter.hh:53
An adapter for mpnc to use the capillary pressure-saturation relationships.
Definition: mpadapter.hh:112
auto capillaryPressures(const FluidState &state, int wPhaseIdx) const
The capillary pressure-saturation curve.
Definition: mpadapter.hh:134
typename MaterialLaw::Scalar Scalar
Definition: mpadapter.hh:122
auto relativePermeabilities(const FluidState &state, int wPhaseIdx) const
The relative permeability of all phases.
Definition: mpadapter.hh:153
MPAdapter(const MaterialLaw &pcKrS)
Definition: mpadapter.hh:124