version 3.8
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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
15#ifndef DUMUX_MP_ADAPTER_HH
16#define DUMUX_MP_ADAPTER_HH
17
18#include <dune/common/fvector.hh>
21
22namespace Dumux::FluidMatrix {
23
28template <class MaterialLaw, int numFluidPhases = std::decay_t<MaterialLaw>::numFluidPhases()>
30{
31 static_assert(AlwaysFalse<MaterialLaw>::value, "Adapter not implemented for the specified number of phases");
32};
33
34
35template<class MaterialLaw>
36class MPAdapter<MaterialLaw, 2>
37: public Adapter<MPAdapter<MaterialLaw, 2>, MultiPhasePcKrSw>
38{
39public:
40 using Scalar = typename std::decay_t<MaterialLaw>::Scalar;
41
42 MPAdapter(MaterialLaw&& pcKrS)
43 : pcKrS_(std::forward<MaterialLaw>(pcKrS))
44 {}
45
51 template <class FluidState>
52 auto capillaryPressures(const FluidState& state, int wPhaseIdx) const
53 {
54 Dune::FieldVector<typename FluidState::Scalar, 2> values;
55
56 const int nPhaseIdx = 1 - wPhaseIdx;
57 // non-wetting phase gets the capillary pressure added
58 values[nPhaseIdx] = 0;
59 // wetting phase does not get anything added
60 values[wPhaseIdx] = - pcKrS_.pc(state.saturation(wPhaseIdx));
61
62 return values;
63 }
64
70 template <class FluidState>
71 auto relativePermeabilities(const FluidState& state, int wPhaseIdx) const
72 {
73 Dune::FieldVector<typename FluidState::Scalar, 2> values;
74
75 const int nPhaseIdx = 1 - wPhaseIdx;
76 values[wPhaseIdx] = pcKrS_.krw(state.saturation(wPhaseIdx));
77 values[nPhaseIdx] = pcKrS_.krn(state.saturation(wPhaseIdx));
78
79 return values;
80 }
81private:
82 MaterialLaw pcKrS_;
83};
84
91template<typename T>
93
94} // end namespace Dumux::FluidMatrix
95
96
97#endif
auto capillaryPressures(const FluidState &state, int wPhaseIdx) const
The capillary pressure-saturation curve.
Definition: mpadapter.hh:52
auto relativePermeabilities(const FluidState &state, int wPhaseIdx) const
The relative permeability of all phases.
Definition: mpadapter.hh:71
MPAdapter(MaterialLaw &&pcKrS)
Definition: mpadapter.hh:42
typename std::decay_t< MaterialLaw >::Scalar Scalar
Definition: mpadapter.hh:40
An adapter for mpnc to use the capillary pressure-saturation relationships.
Definition: mpadapter.hh:30
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:23
Template which always yields a false value.
Definition: common/typetraits/typetraits.hh:24
Adapter to inherit from, allowing the inheriting class to be wrapped by the makeFluidMatrixInteractio...
Definition: fluidmatrixinteraction.hh:55