27#ifndef DUMUX_MP_LINEAR_MATERIAL_HH
28#define DUMUX_MP_LINEAR_MATERIAL_HH
31#include <dune/common/fvector.hh>
44template<
class S,
int numFlu
idPhases>
46:
public Adapter<MPLinearMaterial<S, numFluidPhases>, MultiPhasePcKrSw>
50 Params(
const std::array<S, numFluidPhases>& pcMaxSat,
51 const std::array<S, numFluidPhases>& pcMinSat)
52 : pcMaxSat_(pcMaxSat), pcMinSat_(pcMinSat) {}
54 S pcMaxSat(
int phaseIdx)
const {
return pcMaxSat_[phaseIdx]; }
55 S pcMinSat(
int phaseIdx)
const {
return pcMinSat_[phaseIdx]; }
57 std::array<S, numFluidPhases> pcMaxSat_;
58 std::array<S, numFluidPhases> pcMinSat_;
66 : basicParams_(basicParams)
80 template <
class Flu
idState>
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)
87 const Scalar saturation = state.saturation(phaseIdx);
89 saturation*basicParams_.pcMaxSat(phaseIdx) +
90 (1 - saturation)*basicParams_.pcMinSat(phaseIdx);
100 template <
class Flu
idState>
103 static_assert(FluidState::numPhases == numFluidPhases,
"FluidState doesn't match the number of fluid phases!");
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);
Wrapper type to combine an arbitrary number of different laws for fluid-matrix interaction (e....
Definition brookscorey.hh:35
Adapter to inherit from, allowing the inheriting class to be wrapped by the makeFluidMatrixInteractio...
Definition fluidmatrixinteraction.hh:67
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