opm-common
Loading...
Searching...
No Matches
EffToAbsLaw.hpp
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 This file is part of the Open Porous Media project (OPM).
5
6 OPM 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 2 of the License, or
9 (at your option) any later version.
10
11 OPM 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 OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
27#ifndef OPM_EFF_TO_ABS_LAW_HPP
28#define OPM_EFF_TO_ABS_LAW_HPP
29
30#include "EffToAbsLawParams.hpp"
31
33
34namespace Opm {
69template <class EffLawT, class ParamsT = EffToAbsLawParams<typename EffLawT::Params, EffLawT::numPhases> >
70class EffToAbsLaw : public EffLawT::Traits
71{
72 typedef EffLawT EffLaw;
73
74public:
75 typedef typename EffLaw::Traits Traits;
76 typedef ParamsT Params;
77 typedef typename EffLaw::Scalar Scalar;
78
80 static const int numPhases = EffLaw::numPhases;
81
84 static const bool implementsTwoPhaseApi = EffLaw::implementsTwoPhaseApi;
85
88 static const bool implementsTwoPhaseSatApi = EffLaw::implementsTwoPhaseSatApi;
89
92 static const bool isSaturationDependent = EffLaw::isSaturationDependent;
93
96 static const bool isPressureDependent = EffLaw::isPressureDependent;
97
100 static const bool isTemperatureDependent = EffLaw::isTemperatureDependent;
101
104 static const bool isCompositionDependent = EffLaw::isCompositionDependent;
105
116 template <class Container, class FluidState>
117 static void capillaryPressures(Container& values, const Params& params, const FluidState& fs)
118 {
119 typedef SaturationOverlayFluidState<FluidState> OverlayFluidState;
120
121 OverlayFluidState overlayFs(fs);
122 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
123 overlayFs.setSaturation(phaseIdx,
124 effectiveSaturation(params,
125 fs.saturation(phaseIdx),
126 phaseIdx));
127 }
128
129 EffLaw::template capillaryPressures<Container, OverlayFluidState>(values, params, overlayFs);
130 }
131
142 template <class Container, class FluidState>
143 static void relativePermeabilities(Container& values, const Params& params, const FluidState& fs)
144 {
145 typedef SaturationOverlayFluidState<FluidState> OverlayFluidState;
146
147 OverlayFluidState overlayFs(fs);
148 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
149 overlayFs.setSaturation(phaseIdx,
150 effectiveSaturation(params,
151 fs.saturation(phaseIdx),
152 phaseIdx));
153 }
154
155 EffLaw::template relativePermeabilities<Container, OverlayFluidState>(values, params, overlayFs);
156 }
157
170 template <class FluidState, class Evaluation = typename FluidState::ValueType>
171 static Evaluation pcnw(const Params& params, const FluidState& fs)
172 {
173 typedef SaturationOverlayFluidState<FluidState> OverlayFluidState;
174
175 static_assert(FluidState::numPhases == numPhases,
176 "The fluid state and the material law must exhibit the same "
177 "number of phases!");
178
179 OverlayFluidState overlayFs(fs);
180 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
181 overlayFs.setSaturation(phaseIdx,
182 effectiveSaturation(params,
183 fs.saturation(phaseIdx),
184 phaseIdx));
185 }
186
187 return EffLaw::template pcnw<OverlayFluidState, Evaluation>(params, overlayFs);
188 }
189
190 template <class Evaluation>
191 static typename std::enable_if<implementsTwoPhaseSatApi, Evaluation>::type
192 twoPhaseSatPcnw(const Params& params, const Evaluation& SwAbs)
193 {
194 const Evaluation& SwEff = effectiveSaturation(params, SwAbs, Traits::wettingPhaseIdx);
195
196 return EffLaw::twoPhaseSatPcnw(params, SwEff);
197 }
198
202 template <class Container, class FluidState>
203 static void saturations(Container& values, const Params& params, const FluidState& fs)
204 {
205 EffLaw::template saturations<Container, FluidState>(values, params, fs);
206 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
207 values[phaseIdx] = absoluteSaturation(params, values[phaseIdx], phaseIdx);
208 }
209 }
210
215 template <class FluidState, class Evaluation = typename FluidState::ValueType>
216 static Evaluation Sw(const Params& params, const FluidState& fs)
217 {
218 return absoluteSaturation(params,
219 EffLaw::template Sw<FluidState, Evaluation>(params, fs),
220 Traits::wettingPhaseIdx);
221 }
222
223 template <class Evaluation>
224 static typename std::enable_if<implementsTwoPhaseSatApi, Evaluation>::type
225 twoPhaseSatSw(const Params& params, const Evaluation& Sw)
226 { return absoluteSaturation(params,
227 EffLaw::twoPhaseSatSw(params, Sw),
228 Traits::wettingPhaseIdx); }
229
234 template <class FluidState, class Evaluation = typename FluidState::ValueType>
235 static Evaluation Sn(const Params& params, const FluidState& fs)
236 {
237 return absoluteSaturation(params,
238 EffLaw::template Sn<FluidState, Evaluation>(params, fs),
239 Traits::nonWettingPhaseIdx);
240 }
241
242 template <class Evaluation>
243 static typename std::enable_if<implementsTwoPhaseSatApi, Evaluation>::type
244 twoPhaseSatSn(const Params& params, const Evaluation& Sw)
245 {
246 return absoluteSaturation(params,
247 EffLaw::twoPhaseSatSn(params, Sw),
248 Traits::nonWettingPhaseIdx);
249 }
250
257 template <class FluidState, class Evaluation = typename FluidState::ValueType>
258 static typename std::enable_if< (Traits::numPhases > 2), Evaluation>::type
259 Sg(const Params& params, const FluidState& fs)
260 {
261 return absoluteSaturation(params,
262 EffLaw::template Sg<FluidState, Evaluation>(params, fs),
263 Traits::gasPhaseIdx);
264 }
265
276 template <class FluidState, class Evaluation = typename FluidState::ValueType>
277 static Evaluation krw(const Params& params, const FluidState& fs)
278 {
279 typedef SaturationOverlayFluidState<FluidState> OverlayFluidState;
280
281 static_assert(FluidState::numPhases == numPhases,
282 "The fluid state and the material law must exhibit the same "
283 "number of phases!");
284
285 OverlayFluidState overlayFs(fs);
286 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
287 overlayFs.setSaturation(phaseIdx,
288 effectiveSaturation(params,
289 fs.saturation(phaseIdx),
290 phaseIdx));
291 }
292
293 return EffLaw::template krw<OverlayFluidState, Evaluation>(params, overlayFs);
294 }
295
296 template <class Evaluation>
297 static typename std::enable_if<implementsTwoPhaseSatApi, Evaluation>::type
298 twoPhaseSatKrw(const Params& params, const Evaluation& Sw)
299 { return EffLaw::twoPhaseSatKrw(params, effectiveSaturation(params, Sw, Traits::wettingPhaseIdx)); }
300
304 template <class FluidState, class Evaluation = typename FluidState::ValueType>
305 static Evaluation krn(const Params& params, const FluidState& fs)
306 {
307 typedef SaturationOverlayFluidState<FluidState> OverlayFluidState;
308
309 static_assert(FluidState::numPhases == numPhases,
310 "The fluid state and the material law must exhibit the same "
311 "number of phases!");
312
313 OverlayFluidState overlayFs(fs);
314 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
315 overlayFs.setSaturation(phaseIdx,
316 effectiveSaturation(params,
317 fs.saturation(phaseIdx),
318 phaseIdx));
319 }
320
321 return EffLaw::template krn<OverlayFluidState, Evaluation>(params, overlayFs);
322 }
323
324 template <class Evaluation>
325 static typename std::enable_if<implementsTwoPhaseSatApi, Evaluation>::type
326 twoPhaseSatKrn(const Params& params, const Evaluation& Sw)
327 { return EffLaw::twoPhaseSatKrn(params, effectiveSaturation(params, Sw, Traits::wettingPhaseIdx)); }
328
334 template <class FluidState, class Evaluation = typename FluidState::ValueType>
335 static typename std::enable_if< (Traits::numPhases > 2), Evaluation>::type
336 krg(const Params& params, const FluidState& fs)
337 {
338 typedef SaturationOverlayFluidState<FluidState> OverlayFluidState;
339
340 static_assert(FluidState::numPhases == numPhases,
341 "The fluid state and the material law must exhibit the same "
342 "number of phases!");
343
344 OverlayFluidState overlayFs(fs);
345 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
346 overlayFs.setSaturation(phaseIdx,
347 effectiveSaturation(params,
348 fs.saturation(phaseIdx),
349 phaseIdx));
350 }
351
352 return EffLaw::template krg<OverlayFluidState, Evaluation>(params, overlayFs);
353 }
354
358 template <class Evaluation>
359 static Evaluation effectiveSaturation(const Params& params, const Evaluation& S, unsigned phaseIdx)
360 { return (S - params.residualSaturation(phaseIdx))/(1.0 - params.sumResidualSaturations()); }
361
365 template <class Evaluation>
366 static Evaluation absoluteSaturation(const Params& params, const Evaluation& S, unsigned phaseIdx)
367 { return S*(1.0 - params.sumResidualSaturations()) + params.residualSaturation(phaseIdx); }
368
369private:
378 static Scalar dSeff_dSabs_(const Params& params, int /*phaseIdx*/)
379 { return 1.0/(1 - params.sumResidualSaturations()); }
380
389 static Scalar dSabs_dSeff_(const Params& params, int /*phaseIdx*/)
390 { return 1 - params.sumResidualSaturations(); }
391};
392} // namespace Opm
393
394#endif
A default implementation of the parameters for the adapter class to convert material laws from effect...
This is a fluid state which allows to set the fluid saturations and takes all other quantities from a...
This material law takes a material law defined for effective saturations and converts it to a materia...
Definition EffToAbsLaw.hpp:71
static const bool implementsTwoPhaseApi
Specify whether this material law implements the two-phase convenience API.
Definition EffToAbsLaw.hpp:84
static const bool isPressureDependent
Specify whether the quantities defined by this material law are dependent on the absolute pressure.
Definition EffToAbsLaw.hpp:96
static const bool isTemperatureDependent
Specify whether the quantities defined by this material law are temperature dependent.
Definition EffToAbsLaw.hpp:100
static void capillaryPressures(Container &values, const Params &params, const FluidState &fs)
The capillary pressure-saturation curves depending on absolute saturations.
Definition EffToAbsLaw.hpp:117
static const bool isCompositionDependent
Specify whether the quantities defined by this material law are dependent on the phase composition.
Definition EffToAbsLaw.hpp:104
static Evaluation pcnw(const Params &params, const FluidState &fs)
The capillary pressure-saturation curve.
Definition EffToAbsLaw.hpp:171
static Evaluation krn(const Params &params, const FluidState &fs)
The relative permeability of the non-wetting phase.
Definition EffToAbsLaw.hpp:305
static Evaluation Sn(const Params &params, const FluidState &fs)
Calculate non-wetting liquid phase saturation given that the rest of the fluid state has been initial...
Definition EffToAbsLaw.hpp:235
static const bool isSaturationDependent
Specify whether the quantities defined by this material law are saturation dependent.
Definition EffToAbsLaw.hpp:92
static Evaluation krw(const Params &params, const FluidState &fs)
The relative permeability for the wetting phase.
Definition EffToAbsLaw.hpp:277
static std::enable_if<(Traits::numPhases >2), Evaluation >::type Sg(const Params &params, const FluidState &fs)
Calculate gas phase saturation given that the rest of the fluid state has been initialized.
Definition EffToAbsLaw.hpp:259
static const bool implementsTwoPhaseSatApi
Specify whether this material law implements the two-phase convenience API which only depends on the ...
Definition EffToAbsLaw.hpp:88
static Evaluation absoluteSaturation(const Params &params, const Evaluation &S, unsigned phaseIdx)
Convert an effective saturation to an absolute one.
Definition EffToAbsLaw.hpp:366
static const int numPhases
The number of fluid phases.
Definition EffToAbsLaw.hpp:80
static std::enable_if<(Traits::numPhases >2), Evaluation >::type krg(const Params &params, const FluidState &fs)
The relative permability of the gas phase.
Definition EffToAbsLaw.hpp:336
static Evaluation Sw(const Params &params, const FluidState &fs)
Calculate wetting liquid phase saturation given that the rest of the fluid state has been initialized...
Definition EffToAbsLaw.hpp:216
static Evaluation effectiveSaturation(const Params &params, const Evaluation &S, unsigned phaseIdx)
Convert an absolute saturation to an effective one.
Definition EffToAbsLaw.hpp:359
static void relativePermeabilities(Container &values, const Params &params, const FluidState &fs)
The relative permeability-saturation curves depending on absolute saturations.
Definition EffToAbsLaw.hpp:143
static void saturations(Container &values, const Params &params, const FluidState &fs)
The saturation-capillary pressure curves.
Definition EffToAbsLaw.hpp:203
This is a fluid state which allows to set the fluid saturations and takes all other quantities from a...
Definition SaturationOverlayFluidState.hpp:44
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30