version 3.8
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// 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_LINEAR_MATERIAL_HH
16#define DUMUX_MP_LINEAR_MATERIAL_HH
17
18#include <algorithm>
19#include <dune/common/fvector.hh>
21
22namespace Dumux::FluidMatrix {
23
32template<class S, int numFluidPhases>
34: public Adapter<MPLinearMaterial<S, numFluidPhases>, MultiPhasePcKrSw>
35{
36 struct Params
37 {
38 Params(const std::array<S, numFluidPhases>& pcMaxSat,
39 const std::array<S, numFluidPhases>& pcMinSat)
40 : pcMaxSat_(pcMaxSat), pcMinSat_(pcMinSat) {}
41
42 S pcMaxSat(int phaseIdx) const { return pcMaxSat_[phaseIdx]; }
43 S pcMinSat(int phaseIdx) const { return pcMinSat_[phaseIdx]; }
44 private:
45 std::array<S, numFluidPhases> pcMaxSat_;
46 std::array<S, numFluidPhases> pcMinSat_;
47 };
48
49public:
50 using BasicParams = Params;
51 using Scalar = S;
52
53 MPLinearMaterial(const BasicParams& basicParams)
54 : basicParams_(basicParams)
55 {}
56
68 template <class FluidState>
69 auto capillaryPressures(const FluidState& state, int wPhaseIdx = 0) const
70 {
71 static_assert(FluidState::numPhases == numFluidPhases, "FluidState doesn't match the number of fluid phases!");
72 Dune::FieldVector<typename FluidState::Scalar, numFluidPhases> values;
73 for (int phaseIdx = 0; phaseIdx < numFluidPhases; ++phaseIdx)
74 {
75 const Scalar saturation = state.saturation(phaseIdx);
76 values[phaseIdx] =
77 saturation*basicParams_.pcMaxSat(phaseIdx) +
78 (1 - saturation)*basicParams_.pcMinSat(phaseIdx);
79 }
80 return values;
81 }
82
88 template <class FluidState>
89 auto relativePermeabilities(const FluidState& state, int wPhaseIdx = 0) const
90 {
91 static_assert(FluidState::numPhases == numFluidPhases, "FluidState doesn't match the number of fluid phases!");
92 using std::clamp;
93 Dune::FieldVector<typename FluidState::Scalar, FluidState::numPhases> values;
94 for (int phaseIdx = 0; phaseIdx < FluidState::numPhases; ++phaseIdx)
95 values[phaseIdx] = clamp(state.saturation(phaseIdx), 0.0, 1.0);
96 return values;
97 }
98private:
99 BasicParams basicParams_;
100};
101
102} // end namespace Dumux::FluidMatrix
103
104#endif
Implements a linear saturation-capillary pressure relation.
Definition: mplinearmaterial.hh:35
auto relativePermeabilities(const FluidState &state, int wPhaseIdx=0) const
The relative permeability of all phases.
Definition: mplinearmaterial.hh:89
S Scalar
Definition: mplinearmaterial.hh:51
Params BasicParams
Definition: mplinearmaterial.hh:50
auto capillaryPressures(const FluidState &state, int wPhaseIdx=0) const
The linear capillary pressure-saturation curve.
Definition: mplinearmaterial.hh:69
MPLinearMaterial(const BasicParams &basicParams)
Definition: mplinearmaterial.hh:53
Wrapper type to combine an arbitrary number of different laws for fluid-matrix interaction (e....
Definition: brookscorey.hh:23
std::string saturation(int phaseIdx) noexcept
I/O name of saturation for multiphase systems.
Definition: name.hh:31
Adapter to inherit from, allowing the inheriting class to be wrapped by the makeFluidMatrixInteractio...
Definition: fluidmatrixinteraction.hh:55