3.3.0
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> >
46class [[deprecated("Use FluidMatrix::MPLinearMaterial. Will be removed after 3.3")]] MpLinearMaterial
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::clamp;
94 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
95 values[phaseIdx] = clamp(state.saturation(phaseIdx), 0.0, 1.0);
96 }
97};
98} // end namespace Dumux
99
100#include <dune/common/fvector.hh>
102
103namespace Dumux::FluidMatrix {
104
113template<class S, int numFluidPhases>
115: public Adapter<MPLinearMaterial<S, numFluidPhases>, MultiPhasePcKrSw>
116{
117 struct Params
118 {
119 Params(const std::array<S, numFluidPhases>& pcMaxSat,
120 const std::array<S, numFluidPhases>& pcMinSat)
121 : pcMaxSat_(pcMaxSat), pcMinSat_(pcMinSat) {}
122
123 S pcMaxSat(int phaseIdx) const { return pcMaxSat_[phaseIdx]; }
124 S pcMinSat(int phaseIdx) const { return pcMinSat_[phaseIdx]; }
125 private:
126 std::array<S, numFluidPhases> pcMaxSat_;
127 std::array<S, numFluidPhases> pcMinSat_;
128 };
129
130public:
131 using BasicParams = Params;
132 using Scalar = S;
133
134 MPLinearMaterial(const BasicParams& basicParams)
135 : basicParams_(basicParams)
136 {}
137
149 template <class FluidState>
150 auto capillaryPressures(const FluidState& state, int wPhaseIdx = 0) const
151 {
152 static_assert(FluidState::numPhases == numFluidPhases, "FluidState doesn't match the number of fluid phases!");
153 Dune::FieldVector<typename FluidState::Scalar, numFluidPhases> values;
154 for (int phaseIdx = 0; phaseIdx < numFluidPhases; ++phaseIdx)
155 {
156 const Scalar saturation = state.saturation(phaseIdx);
157 values[phaseIdx] =
158 saturation*basicParams_.pcMaxSat(phaseIdx) +
159 (1 - saturation)*basicParams_.pcMinSat(phaseIdx);
160 }
161 return values;
162 }
163
169 template <class FluidState>
170 auto relativePermeabilities(const FluidState& state, int wPhaseIdx = 0) const
171 {
172 static_assert(FluidState::numPhases == numFluidPhases, "FluidState doesn't match the number of fluid phases!");
173 using std::clamp;
174 Dune::FieldVector<typename FluidState::Scalar, FluidState::numPhases> values;
175 for (int phaseIdx = 0; phaseIdx < FluidState::numPhases; ++phaseIdx)
176 values[phaseIdx] = clamp(state.saturation(phaseIdx), 0.0, 1.0);
177 return values;
178 }
179private:
180 BasicParams basicParams_;
181};
182
183} // end namespace Dumux::FluidMatrix
184
185#endif
Reference implementation of parameters for the M-phase linear material material.
Wrapper type to combine an arbitrary number of different laws for fluid-matrix interaction (e....
Definition: adapt.hh:29
std::string saturation(int phaseIdx) noexcept
I/O name of saturation for multiphase systems.
Definition: name.hh:43
Definition: brookscorey.hh:286
Adapter to inherit from, allowing the inheriting class to be wrapped by the makeFluidMatrixInteractio...
Definition: fluidmatrixinteraction.hh:67
Implements a linear saturation-capillary pressure relation.
Definition: mplinearmaterial.hh:47
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
Implements a linear saturation-capillary pressure relation.
Definition: mplinearmaterial.hh:116
auto relativePermeabilities(const FluidState &state, int wPhaseIdx=0) const
The relative permeability of all phases.
Definition: mplinearmaterial.hh:170
S Scalar
Definition: mplinearmaterial.hh:132
Params BasicParams
Definition: mplinearmaterial.hh:131
auto capillaryPressures(const FluidState &state, int wPhaseIdx=0) const
The linear capillary pressure-saturation curve.
Definition: mplinearmaterial.hh:150
MPLinearMaterial(const BasicParams &basicParams)
Definition: mplinearmaterial.hh:134