3.4
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
30#include <algorithm>
31#include <dune/common/fvector.hh>
33
34namespace Dumux::FluidMatrix {
35
44template<class S, int numFluidPhases>
46: public Adapter<MPLinearMaterial<S, numFluidPhases>, MultiPhasePcKrSw>
47{
48 struct Params
49 {
50 Params(const std::array<S, numFluidPhases>& pcMaxSat,
51 const std::array<S, numFluidPhases>& pcMinSat)
52 : pcMaxSat_(pcMaxSat), pcMinSat_(pcMinSat) {}
53
54 S pcMaxSat(int phaseIdx) const { return pcMaxSat_[phaseIdx]; }
55 S pcMinSat(int phaseIdx) const { return pcMinSat_[phaseIdx]; }
56 private:
57 std::array<S, numFluidPhases> pcMaxSat_;
58 std::array<S, numFluidPhases> pcMinSat_;
59 };
60
61public:
62 using BasicParams = Params;
63 using Scalar = S;
64
65 MPLinearMaterial(const BasicParams& basicParams)
66 : basicParams_(basicParams)
67 {}
68
80 template <class FluidState>
81 auto capillaryPressures(const FluidState& state, int wPhaseIdx = 0) const
82 {
83 static_assert(FluidState::numPhases == numFluidPhases, "FluidState doesn't match the number of fluid phases!");
84 Dune::FieldVector<typename FluidState::Scalar, numFluidPhases> values;
85 for (int phaseIdx = 0; phaseIdx < numFluidPhases; ++phaseIdx)
86 {
87 const Scalar saturation = state.saturation(phaseIdx);
88 values[phaseIdx] =
89 saturation*basicParams_.pcMaxSat(phaseIdx) +
90 (1 - saturation)*basicParams_.pcMinSat(phaseIdx);
91 }
92 return values;
93 }
94
100 template <class FluidState>
101 auto relativePermeabilities(const FluidState& state, int wPhaseIdx = 0) const
102 {
103 static_assert(FluidState::numPhases == numFluidPhases, "FluidState doesn't match the number of fluid phases!");
104 using std::clamp;
105 Dune::FieldVector<typename FluidState::Scalar, FluidState::numPhases> values;
106 for (int phaseIdx = 0; phaseIdx < FluidState::numPhases; ++phaseIdx)
107 values[phaseIdx] = clamp(state.saturation(phaseIdx), 0.0, 1.0);
108 return values;
109 }
110private:
111 BasicParams basicParams_;
112};
113
114} // end namespace Dumux::FluidMatrix
115
116#endif
Wrapper type to combine an arbitrary number of different laws for fluid-matrix interaction (e....
std::string saturation(int phaseIdx) noexcept
I/O name of saturation for multiphase systems.
Definition: name.hh:43
Definition: brookscorey.hh:35
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
auto relativePermeabilities(const FluidState &state, int wPhaseIdx=0) const
The relative permeability of all phases.
Definition: mplinearmaterial.hh:101
S Scalar
Definition: mplinearmaterial.hh:63
Params BasicParams
Definition: mplinearmaterial.hh:62
auto capillaryPressures(const FluidState &state, int wPhaseIdx=0) const
The linear capillary pressure-saturation curve.
Definition: mplinearmaterial.hh:81
MPLinearMaterial(const BasicParams &basicParams)
Definition: mplinearmaterial.hh:65