27#ifndef OPM_LINEAR_MATERIAL_HPP
28#define OPM_LINEAR_MATERIAL_HPP
49template <
class TraitsT,
class ParamsT = LinearMaterialParams<TraitsT> >
53 typedef TraitsT Traits;
54 typedef ParamsT Params;
55 typedef typename Traits::Scalar Scalar;
58 using Traits::numPhases;
96 template <
class ContainerT,
class Flu
idState>
99 const FluidState& state)
101 typedef typename std::remove_reference<
decltype(values[0])>::type Evaluation;
103 for (
unsigned phaseIdx = 0; phaseIdx < Traits::numPhases; ++phaseIdx) {
104 const Evaluation& S =
105 decay<Evaluation>(state.saturation(phaseIdx));
106 Valgrind::CheckDefined(S);
109 S*params.pcMaxSat(phaseIdx) +
110 (1.0 - S)*params.pcMinSat(phaseIdx);
117 template <
class ContainerT,
class Flu
idState>
122 throw std::runtime_error(
"Not implemented: LinearMaterial::saturations()");
128 template <
class ContainerT,
class Flu
idState>
131 const FluidState& state)
133 typedef typename std::remove_reference<
decltype(values[0])>::type Evaluation;
135 for (
unsigned phaseIdx = 0; phaseIdx < Traits::numPhases; ++phaseIdx) {
136 const Evaluation& S =
137 decay<Evaluation>(state.saturation(phaseIdx));
138 Valgrind::CheckDefined(S);
140 values[phaseIdx] = max(min(S,1.0),0.0);
147 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
148 static Evaluation
pcnw(
const Params& params,
const FluidState& fs)
150 const Evaluation&
Sw =
151 decay<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
152 Valgrind::CheckDefined(
Sw);
154 const Evaluation& wPhasePressure =
155 Sw*params.pcMaxSat(Traits::wettingPhaseIdx) +
156 (1.0 -
Sw)*params.pcMinSat(Traits::wettingPhaseIdx);
158 const Evaluation&
Sn =
159 decay<Evaluation>(fs.saturation(Traits::nonWettingPhaseIdx));
160 Valgrind::CheckDefined(
Sn);
162 const Evaluation& nPhasePressure =
163 Sn*params.pcMaxSat(Traits::nonWettingPhaseIdx) +
164 (1.0 -
Sn)*params.pcMinSat(Traits::nonWettingPhaseIdx);
166 return nPhasePressure - wPhasePressure;
170 template <
class Evaluation = Scalar>
171 static typename std::enable_if<Traits::numPhases == 2, Evaluation>::type
172 twoPhaseSatPcnw(
const Params& params,
const Evaluation&
Sw)
174 const Evaluation& wPhasePressure =
175 Sw*params.pcMaxSat(Traits::wettingPhaseIdx) +
176 (1.0 -
Sw)*params.pcMinSat(Traits::wettingPhaseIdx);
178 const Evaluation& nPhasePressure =
179 (1.0 -
Sw)*params.pcMaxSat(Traits::nonWettingPhaseIdx) +
180 Sw*params.pcMinSat(Traits::nonWettingPhaseIdx);
182 return nPhasePressure - wPhasePressure;
189 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
190 static Evaluation
Sw(
const Params& ,
const FluidState& )
191 {
throw std::runtime_error(
"Not implemented: Sw()"); }
193 template <
class Evaluation = Scalar>
194 static typename std::enable_if<Traits::numPhases == 2, Evaluation>::type
195 twoPhaseSatSw(
const Params& ,
const Evaluation& )
196 {
throw std::runtime_error(
"Not implemented: twoPhaseSatSw()"); }
202 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
203 static Evaluation
Sn(
const Params& ,
const FluidState& )
204 {
throw std::runtime_error(
"Not implemented: Sn()"); }
206 template <
class Evaluation = Scalar>
207 static typename std::enable_if<Traits::numPhases == 2, Evaluation>::type
208 twoPhaseSatSn(
const Params& ,
const Evaluation& )
209 {
throw std::runtime_error(
"Not implemented: twoPhaseSatSn()"); }
217 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
218 static typename std::enable_if<Traits::numPhases == 3, Evaluation>::type
219 Sg(
const Params& ,
const FluidState& )
220 {
throw std::runtime_error(
"Not implemented: Sg()"); }
225 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
226 static Evaluation
krw(
const Params& ,
const FluidState& fs)
228 const Evaluation& sw =
229 decay<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
230 return max(0.0, min(1.0, sw));
233 template <
class Evaluation = Scalar>
234 static typename std::enable_if<Traits::numPhases == 2, Evaluation>::type
235 twoPhaseSatKrw(
const Params& ,
const Evaluation&
Sw)
236 {
return max(0.0, min(1.0,
Sw)); }
241 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
242 static Evaluation
krn(
const Params& ,
const FluidState& fs)
244 const Evaluation& sn =
245 decay<Evaluation>(fs.saturation(Traits::nonWettingPhaseIdx));
246 return max(0.0, min(1.0, sn));
249 template <
class Evaluation = Scalar>
250 static typename std::enable_if<Traits::numPhases == 2, Evaluation>::type
251 twoPhaseSatKrn(
const Params& ,
const Evaluation&
Sw)
253 return max(0.0, min(1.0,
Sw));
261 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
262 static typename std::enable_if<Traits::numPhases == 3, Evaluation>::type
263 krg(
const Params& ,
const FluidState& fs)
265 const Evaluation& sg =
266 decay<Evaluation>(fs.saturation(Traits::gasPhaseIdx));
267 return max(0.0, min(1.0, sg));
275 template <
class Flu
idState,
class Evaluation = Scalar>
276 static typename std::enable_if<Traits::numPhases == 3, Evaluation>::type
277 pcgn(
const Params& params,
const FluidState& fs)
279 const Evaluation& sn =
280 decay<Evaluation>(fs.saturation(Traits::nonWettingPhaseIdx));
281 Valgrind::CheckDefined(sn);
283 const Evaluation& nPhasePressure =
284 Sn*params.pcMaxSat(Traits::nonWettingPhaseIdx) +
285 (1.0 -
Sn)*params.pcMinSat(Traits::nonWettingPhaseIdx);
287 const Evaluation& sg =
288 decay<Evaluation>(fs.saturation(Traits::gasPhaseIdx));
289 Valgrind::CheckDefined(sg);
291 const Evaluation& gPhasePressure =
292 Sg*params.pcMaxSat(Traits::gasPhaseIdx) +
293 (1.0 -
Sg)*params.pcMinSat(Traits::gasPhaseIdx);
295 return gPhasePressure - nPhasePressure;
Reference implementation of params for the linear M-phase material material.
Some templates to wrap the valgrind client request macros.
Implements a linear saturation-capillary pressure relation.
Definition LinearMaterial.hpp:51
static Evaluation pcnw(const Params ¶ms, const FluidState &fs)
The difference between the pressures of the non-wetting and wetting phase.
Definition LinearMaterial.hpp:148
static Evaluation krn(const Params &, const FluidState &fs)
The relative permability of the liquid non-wetting phase.
Definition LinearMaterial.hpp:242
static void relativePermeabilities(ContainerT &values, const Params &, const FluidState &state)
The relative permeability of all phases.
Definition LinearMaterial.hpp:129
static const bool implementsTwoPhaseSatApi
Specify whether this material law implements the two-phase convenience API which only depends on the ...
Definition LinearMaterial.hpp:66
static const bool isPressureDependent
Specify whether the quantities defined by this material law are dependent on the absolute pressure.
Definition LinearMaterial.hpp:74
static std::enable_if< Traits::numPhases==3, Evaluation >::type pcgn(const Params ¶ms, const FluidState &fs)
The difference between the pressures of the gas and the non-wetting phase.
Definition LinearMaterial.hpp:277
static const bool isSaturationDependent
Specify whether the quantities defined by this material law are saturation dependent.
Definition LinearMaterial.hpp:70
static std::enable_if< Traits::numPhases==3, Evaluation >::type Sg(const Params &, const FluidState &)
Calculate gas phase saturation given that the rest of the fluid state has been initialized.
Definition LinearMaterial.hpp:219
static Evaluation Sw(const Params &, const FluidState &)
Calculate wetting phase saturation given that the rest of the fluid state has been initialized.
Definition LinearMaterial.hpp:190
static Evaluation Sn(const Params &, const FluidState &)
Calculate non-wetting liquid phase saturation given that the rest of the fluid state has been initial...
Definition LinearMaterial.hpp:203
static const bool isTemperatureDependent
Specify whether the quantities defined by this material law are temperature dependent.
Definition LinearMaterial.hpp:78
static void saturations(ContainerT &, const Params &, const FluidState &)
The inverse of the capillary pressure.
Definition LinearMaterial.hpp:118
static const bool isCompositionDependent
Specify whether the quantities defined by this material law are dependent on the phase composition.
Definition LinearMaterial.hpp:82
static const bool implementsTwoPhaseApi
Specify whether this material law implements the two-phase convenience API.
Definition LinearMaterial.hpp:62
static void capillaryPressures(ContainerT &values, const Params ¶ms, const FluidState &state)
The linear capillary pressure-saturation curve.
Definition LinearMaterial.hpp:97
static std::enable_if< Traits::numPhases==3, Evaluation >::type krg(const Params &, const FluidState &fs)
The relative permability of the gas phase.
Definition LinearMaterial.hpp:263
static Evaluation krw(const Params &, const FluidState &fs)
The relative permability of the wetting phase.
Definition LinearMaterial.hpp:226
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30