27#ifndef DUMUX_MP_LINEAR_MATERIAL_HH
28#define DUMUX_MP_LINEAR_MATERIAL_HH
45template <
int numPhasesV,
class ScalarT,
class ParamsT = MpLinearMaterialParams<numPhasesV, ScalarT> >
46class [[deprecated("Use FluidMatrix::MPLinearMaterial. Will be removed after 3.3")]]
MpLinearMaterial
50 using Scalar =
typename Params::Scalar;
51 enum { numPhases = numPhasesV };
66 template <
class ContainerT,
class Flu
idState>
69 const FluidState &state,
72 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
73 Scalar S = state.saturation(phaseIdx);
75 S*params.pcMaxSat(phaseIdx) +
76 (1 - S)*params.pcMinSat(phaseIdx);
87 template <
class ContainerT,
class Flu
idState>
90 const FluidState &state,
94 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
95 values[phaseIdx] = clamp(state.saturation(phaseIdx), 0.0, 1.0);
100#include <dune/common/fvector.hh>
113template<
class S,
int numFlu
idPhases>
115:
public Adapter<MPLinearMaterial<S, numFluidPhases>, MultiPhasePcKrSw>
119 Params(
const std::array<S, numFluidPhases>& pcMaxSat,
120 const std::array<S, numFluidPhases>& pcMinSat)
121 : pcMaxSat_(pcMaxSat), pcMinSat_(pcMinSat) {}
123 S pcMaxSat(
int phaseIdx)
const {
return pcMaxSat_[phaseIdx]; }
124 S pcMinSat(
int phaseIdx)
const {
return pcMinSat_[phaseIdx]; }
126 std::array<S, numFluidPhases> pcMaxSat_;
127 std::array<S, numFluidPhases> pcMinSat_;
135 : basicParams_(basicParams)
149 template <
class Flu
idState>
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)
159 (1 -
saturation)*basicParams_.pcMinSat(phaseIdx);
169 template <
class Flu
idState>
172 static_assert(FluidState::numPhases == numFluidPhases,
"FluidState doesn't match the number of fluid phases!");
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);
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....
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 ¶ms, const FluidState &state, int wPhaseIdx=0)
The linear capillary pressure-saturation curve.
Definition: mplinearmaterial.hh:67
static void relativePermeabilities(ContainerT &values, const Params ¶ms, 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