3.5-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 <dune/common/fvector.hh>
33
34namespace Dumux::FluidMatrix {
35
40template <class MaterialLaw, int numFluidPhases = std::decay_t<MaterialLaw>::numFluidPhases()>
42{
43 static_assert(AlwaysFalse<MaterialLaw>::value, "Adapter not implemented for the specified number of phases");
44};
45
46
47template<class MaterialLaw>
48class MPAdapter<MaterialLaw, 2>
49: public Adapter<MPAdapter<MaterialLaw, 2>, MultiPhasePcKrSw>
50{
51public:
52 using Scalar = typename std::decay_t<MaterialLaw>::Scalar;
53
54 MPAdapter(MaterialLaw&& pcKrS)
55 : pcKrS_(std::forward<MaterialLaw>(pcKrS))
56 {}
57
63 template <class FluidState>
64 auto capillaryPressures(const FluidState& state, int wPhaseIdx) const
65 {
66 Dune::FieldVector<typename FluidState::Scalar, 2> values;
67
68 const int nPhaseIdx = 1 - wPhaseIdx;
69 // non-wetting phase gets the capillary pressure added
70 values[nPhaseIdx] = 0;
71 // wetting phase does not get anything added
72 values[wPhaseIdx] = - pcKrS_.pc(state.saturation(wPhaseIdx));
73
74 return values;
75 }
76
82 template <class FluidState>
83 auto relativePermeabilities(const FluidState& state, int wPhaseIdx) const
84 {
85 Dune::FieldVector<typename FluidState::Scalar, 2> values;
86
87 const int nPhaseIdx = 1 - wPhaseIdx;
88 values[wPhaseIdx] = pcKrS_.krw(state.saturation(wPhaseIdx));
89 values[nPhaseIdx] = pcKrS_.krn(state.saturation(wPhaseIdx));
90
91 return values;
92 }
93private:
94 MaterialLaw pcKrS_;
95};
96
103template<typename T>
105
106} // end namespace Dumux::FluidMatrix
107
108
109#endif
Type traits.
Wrapper type to combine an arbitrary number of different laws for fluid-matrix interaction (e....
MPAdapter(T &&) -> MPAdapter< T >
Deduction guide for the MPAdapter class. Makes sure that MPAdapter stores a copy of T if the construc...
Definition: brookscorey.hh:35
Template which always yields a false value.
Definition: typetraits.hh:35
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:42
auto capillaryPressures(const FluidState &state, int wPhaseIdx) const
The capillary pressure-saturation curve.
Definition: mpadapter.hh:64
auto relativePermeabilities(const FluidState &state, int wPhaseIdx) const
The relative permeability of all phases.
Definition: mpadapter.hh:83
MPAdapter(MaterialLaw &&pcKrS)
Definition: mpadapter.hh:54
typename std::decay_t< MaterialLaw >::Scalar Scalar
Definition: mpadapter.hh:52