opm-common
Loading...
Searching...
No Matches
EclTwoPhaseMaterial.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_ECL_TWO_PHASE_MATERIAL_HPP
28#define OPM_ECL_TWO_PHASE_MATERIAL_HPP
29
31
32#include <opm/common/TimingMacros.hpp>
33
36
37namespace Opm {
38
48template <class TraitsT,
49 class GasOilMaterialLawT,
50 class OilWaterMaterialLawT,
51 class GasWaterMaterialLawT,
52 class ParamsT = EclTwoPhaseMaterialParams<TraitsT,
53 typename GasOilMaterialLawT::Params,
54 typename OilWaterMaterialLawT::Params,
55 typename GasWaterMaterialLawT::Params> >
56class EclTwoPhaseMaterial : public TraitsT
57{
58public:
59 using GasOilMaterialLaw = GasOilMaterialLawT;
60 using OilWaterMaterialLaw = OilWaterMaterialLawT;
61 using GasWaterMaterialLaw = GasWaterMaterialLawT;
62
63 // some safety checks
64 static_assert(TraitsT::numPhases == 3,
65 "The number of phases considered by this capillary pressure "
66 "law is always three!");
67 static_assert(GasOilMaterialLaw::numPhases == 2,
68 "The number of phases considered by the gas-oil capillary "
69 "pressure law must be two!");
70 static_assert(OilWaterMaterialLaw::numPhases == 2,
71 "The number of phases considered by the oil-water capillary "
72 "pressure law must be two!");
73 static_assert(GasWaterMaterialLaw::numPhases == 2,
74 "The number of phases considered by the gas-water capillary "
75 "pressure law must be two!");
76 static_assert(std::is_same<typename GasOilMaterialLaw::Scalar,
77 typename OilWaterMaterialLaw::Scalar>::value,
78 "The two two-phase capillary pressure laws must use the same "
79 "type of floating point values.");
80
81 using Traits = TraitsT;
82 using Params = ParamsT;
83 using Scalar = typename Traits::Scalar;
84
85 static constexpr int numPhases = 3;
86 static constexpr int waterPhaseIdx = Traits::wettingPhaseIdx;
87 static constexpr int oilPhaseIdx = Traits::nonWettingPhaseIdx;
88 static constexpr int gasPhaseIdx = Traits::gasPhaseIdx;
89
92 static constexpr bool implementsTwoPhaseApi = false;
93
96 static constexpr bool implementsTwoPhaseSatApi = false;
97
100 static constexpr bool isSaturationDependent = true;
101
104 static constexpr bool isPressureDependent = false;
105
108 static constexpr bool isTemperatureDependent = false;
109
112 static constexpr bool isCompositionDependent = false;
113
114 template <class ContainerT, class FluidState>
115 static Scalar relpermOilInOilGasSystem(const Params& /*params*/,
116 const FluidState& /*fluidState*/) {
117 throw std::logic_error {
118 "relpermOilInOilGasSystem() is specific to three phases"
119 };
120 }
121 template <class ContainerT, class FluidState>
122 static Scalar relpermOilInOilWaterSystem(const Params& /*params*/,
123 const FluidState& /*fluidState*/) {
124 throw std::logic_error {
125 "relpermOilInOilWaterSystem() is specific to three phases"
126 };
127 }
128
143
144 template <class ContainerT, class FluidState, class ...Args>
145 static void capillaryPressures(ContainerT& values,
146 const Params& params,
147 const FluidState& fluidState)
148 {
149 OPM_TIMEFUNCTION_LOCAL(Subsystem::SatProps);
150 using Evaluation = typename std::remove_reference<decltype(values[0])>::type;
151
152 switch (params.approach()) {
153 case EclTwoPhaseApproach::GasOil: {
154 const Evaluation& So =
155 decay<Evaluation>(fluidState.saturation(oilPhaseIdx));
156
157 values[oilPhaseIdx] = 0.0;
158 values[gasPhaseIdx] = GasOilMaterialLaw::template twoPhaseSatPcnw<Evaluation, Args...>(params.gasOilParams(), So);
159 break;
160 }
161
162 case EclTwoPhaseApproach::OilWater: {
163 const Evaluation& Sw =
164 decay<Evaluation>(fluidState.saturation(waterPhaseIdx));
165
166 values[waterPhaseIdx] = 0.0;
167 values[oilPhaseIdx] = OilWaterMaterialLaw::template twoPhaseSatPcnw<Evaluation, Args...>(params.oilWaterParams(), Sw);
168 break;
169 }
170
171 case EclTwoPhaseApproach::GasWater: {
172 const Evaluation& Sw =
173 decay<Evaluation>(fluidState.saturation(waterPhaseIdx));
174
175 values[waterPhaseIdx] = 0.0;
176 values[gasPhaseIdx] = GasWaterMaterialLaw::template twoPhaseSatPcnw<Evaluation, Args...>(params.gasWaterParams(), Sw);
177 break;
178 }
179
180 }
181 }
182
183 /*
184 * Hysteresis parameters for oil-water
185 * @see EclHysteresisTwoPhaseLawParams::soMax(...)
186 * @see EclHysteresisTwoPhaseLawParams::swMax(...)
187 * @see EclHysteresisTwoPhaseLawParams::swMin(...)
188 * \param params Parameters
189 */
190 static void oilWaterHysteresisParams(Scalar& soMax,
191 Scalar& swMax,
192 Scalar& swMin,
193 const Params& params)
194 {
195 if constexpr (Traits::enableHysteresis) {
196 soMax = 1.0 - params.oilWaterParams().krnSwMdc();
197 swMax = params.oilWaterParams().krwSwMdc();
198 swMin = params.oilWaterParams().pcSwMdc();
202 }
203 }
204
205 /*
206 * Hysteresis parameters for oil-water
207 * @see EclHysteresisTwoPhaseLawParams::soMax(...)
208 * @see EclHysteresisTwoPhaseLawParams::swMax(...)
209 * @see EclHysteresisTwoPhaseLawParams::swMin(...)
210 * \param params Parameters
211 */
212 static void setOilWaterHysteresisParams(const Scalar& soMax,
213 const Scalar& swMax,
214 const Scalar& swMin,
215 Params& params)
216 {
217 if constexpr (Traits::enableHysteresis) {
218 params.oilWaterParams().update(swMin, swMax, 1.0 - soMax);
219 }
220 }
221
222
223 /*
224 * Hysteresis parameters for gas-oil
225 * @see EclHysteresisTwoPhaseLawParams::sgmax(...)
226 * @see EclHysteresisTwoPhaseLawParams::shmax(...)
227 * @see EclHysteresisTwoPhaseLawParams::somin(...)
228 * \param params Parameters
229 */
230 static void gasOilHysteresisParams(Scalar& sgmax,
231 Scalar& shmax,
232 Scalar& somin,
233 const Params& params)
234 {
235 if constexpr (Traits::enableHysteresis) {
236 sgmax = 1.0 - params.gasOilParams().krnSwMdc();
237 shmax = params.gasOilParams().krwSwMdc();
238 somin = params.gasOilParams().pcSwMdc();
239 Valgrind::CheckDefined(sgmax);
240 Valgrind::CheckDefined(shmax);
241 Valgrind::CheckDefined(somin);
242 }
243 }
244
245 /*
246 * Hysteresis parameters for gas-oil
247 * @see EclHysteresisTwoPhaseLawParams::sgmax(...)
248 * @see EclHysteresisTwoPhaseLawParams::shmax(...)
249 * \param params Parameters
250 */
251 static void setGasOilHysteresisParams(const Scalar& sgmax,
252 const Scalar& shmax,
253 const Scalar& somin,
254 Params& params)
255 {
256 if constexpr (Traits::enableHysteresis) {
257 params.gasOilParams().update(somin , shmax, 1.0 - sgmax);
258 }
259 }
260
261 static Scalar trappedGasSaturation(const Params& params, bool maximumTrapping){
262 if(params.approach() == EclTwoPhaseApproach::GasOil)
263 return params.gasOilParams().SnTrapped(maximumTrapping);
264 if(params.approach() == EclTwoPhaseApproach::GasWater)
265 return params.gasWaterParams().SnTrapped(maximumTrapping);
266 return 0.0; // oil-water case
267 }
268
269 static Scalar strandedGasSaturation(const Params& params, Scalar Sg, Scalar Kg){
270 if(params.approach() == EclTwoPhaseApproach::GasOil)
271 return params.gasOilParams().SnStranded(Sg, Kg);
272 if(params.approach() == EclTwoPhaseApproach::GasWater)
273 return params.gasWaterParams().SnStranded(Sg, Kg);
274 return 0.0; // oil-water case
275 }
276
277 static Scalar trappedOilSaturation(const Params& params, bool maximumTrapping){
278 if(params.approach() == EclTwoPhaseApproach::GasOil)
279 return params.gasOilParams().SwTrapped();
280 if(params.approach() == EclTwoPhaseApproach::OilWater)
281 return params.oilWaterParams().SnTrapped(maximumTrapping);
282 return 0.0; // gas-water case
283 }
284
285 static Scalar trappedWaterSaturation(const Params& params){
286 if(params.approach() == EclTwoPhaseApproach::GasWater)
287 return params.gasWaterParams().SwTrapped();
288 if(params.approach() == EclTwoPhaseApproach::OilWater)
289 return params.oilWaterParams().SwTrapped();
290 return 0.0; // gas-oil case
291 }
292
293
303 template <class FluidState, class Evaluation = typename FluidState::Scalar>
304 static Evaluation pcgn(const Params& /* params */,
305 const FluidState& /* fs */)
306 {
307 throw std::logic_error("Not implemented: pcgn()");
308 }
309
319 template <class FluidState, class Evaluation = typename FluidState::Scalar>
320 static Evaluation pcnw(const Params& /* params */,
321 const FluidState& /* fs */)
322 {
323 throw std::logic_error("Not implemented: pcnw()");
324 }
325
329 template <class ContainerT, class FluidState>
330 static void saturations(ContainerT& /* values */,
331 const Params& /* params */,
332 const FluidState& /* fs */)
333 {
334 throw std::logic_error("Not implemented: saturations()");
335 }
336
340 template <class FluidState, class Evaluation = typename FluidState::Scalar>
341 static Evaluation Sg(const Params& /* params */,
342 const FluidState& /* fluidState */)
343 {
344 throw std::logic_error("Not implemented: Sg()");
345 }
346
350 template <class FluidState, class Evaluation = typename FluidState::Scalar>
351 static Evaluation Sn(const Params& /* params */,
352 const FluidState& /* fluidState */)
353 {
354 throw std::logic_error("Not implemented: Sn()");
355 }
356
360 template <class FluidState, class Evaluation = typename FluidState::Scalar>
361 static Evaluation Sw(const Params& /* params */,
362 const FluidState& /* fluidState */)
363 {
364 throw std::logic_error("Not implemented: Sw()");
365 }
366
382 template <class ContainerT, class FluidState, class ...Args>
383 static void relativePermeabilities(ContainerT& values,
384 const Params& params,
385 const FluidState& fluidState)
386 {
387 OPM_TIMEFUNCTION_LOCAL(Subsystem::SatProps);
388 using Evaluation = typename std::remove_reference<decltype(values[0])>::type;
389
390 switch (params.approach()) {
391 case EclTwoPhaseApproach::GasOil: {
392 const Evaluation& So =
393 decay<Evaluation>(fluidState.saturation(oilPhaseIdx));
394
395 values[oilPhaseIdx] = GasOilMaterialLaw::template twoPhaseSatKrw<Evaluation, Args...>(params.gasOilParams(), So);
396 values[gasPhaseIdx] = GasOilMaterialLaw::template twoPhaseSatKrn<Evaluation, Args...>(params.gasOilParams(), So);
397 break;
398 }
399
400 case EclTwoPhaseApproach::OilWater: {
401 const Evaluation& sw =
402 decay<Evaluation>(fluidState.saturation(waterPhaseIdx));
403
404 values[waterPhaseIdx] = OilWaterMaterialLaw::template twoPhaseSatKrw<Evaluation, Args...>(params.oilWaterParams(), sw);
405 values[oilPhaseIdx] = OilWaterMaterialLaw::template twoPhaseSatKrn<Evaluation, Args...>(params.oilWaterParams(), sw);
406 break;
407 }
408
409 case EclTwoPhaseApproach::GasWater: {
410 const Evaluation& sw =
411 decay<Evaluation>(fluidState.saturation(waterPhaseIdx));
412
413 values[waterPhaseIdx] = GasWaterMaterialLaw::template twoPhaseSatKrw<Evaluation, Args...>(params.gasWaterParams(), sw);
414 values[gasPhaseIdx] = GasWaterMaterialLaw::template twoPhaseSatKrn<Evaluation, Args...>(params.gasWaterParams(), sw);
415
416 break;
417 }
418 }
419 }
420
424 template <class FluidState, class Evaluation = typename FluidState::Scalar>
425 static Evaluation krg(const Params& /* params */,
426 const FluidState& /* fluidState */)
427 {
428 throw std::logic_error("Not implemented: krg()");
429 }
430
434 template <class FluidState, class Evaluation = typename FluidState::Scalar>
435 static Evaluation krw(const Params& /* params */,
436 const FluidState& /* fluidState */)
437 {
438 throw std::logic_error("Not implemented: krw()");
439 }
440
444 template <class FluidState, class Evaluation = typename FluidState::Scalar>
445 static Evaluation krn(const Params& /* params */,
446 const FluidState& /* fluidState */)
447 {
448 throw std::logic_error("Not implemented: krn()");
449 }
450
451
459 template <class FluidState>
460 static bool updateHysteresis(Params& params, const FluidState& fluidState)
461 {
462 if constexpr (Traits::enableHysteresis) {
463 OPM_TIMEFUNCTION_LOCAL(Subsystem::SatProps);
464 switch (params.approach()) {
465 case EclTwoPhaseApproach::GasOil: {
466 Scalar So = scalarValue(fluidState.saturation(oilPhaseIdx));
467 return params.gasOilParams().update(/*pcSw=*/So, /*krwSw=*/So, /*krnSw=*/So);
468 }
469 case EclTwoPhaseApproach::OilWater: {
470 Scalar sw = scalarValue(fluidState.saturation(waterPhaseIdx));
471 return params.oilWaterParams().update(/*pcSw=*/sw, /*krwSw=*/sw, /*krnSw=*/sw);
472 }
473 case EclTwoPhaseApproach::GasWater: {
474 Scalar sw = scalarValue(fluidState.saturation(waterPhaseIdx));
475 return params.gasWaterParams().update(/*pcSw=*/sw, /*krwSw=*/sw, /*krnSw=*/sw);
476 }
477 }
478 return false;
479 } else {
480 return false;
481 }
482 }
483};
484
485} // namespace Opm
486
487#endif
Implementation for the parameters required by the material law for two-phase simulations.
A traits class which provides basic mathematical functions for arbitrary scalar floating point values...
Some templates to wrap the valgrind client request macros.
OPM_HOST_DEVICE bool CheckDefined(const T &value)
Make valgrind complain if any of the memory occupied by an object is undefined.
Definition Valgrind.hpp:76
Implementation for the parameters required by the material law for two-phase simulations.
Definition EclTwoPhaseMaterialParams.hpp:51
const GasWaterParams & gasWaterParams() const
The parameter object for the gas-water twophase law.
Definition EclTwoPhaseMaterialParams.hpp:113
const GasOilParams & gasOilParams() const
The parameter object for the gas-oil twophase law.
Definition EclTwoPhaseMaterialParams.hpp:77
const OilWaterParams & oilWaterParams() const
The parameter object for the oil-water twophase law.
Definition EclTwoPhaseMaterialParams.hpp:95
Implements a multiplexer class that provides ECL saturation functions for twophase simulations.
Definition EclTwoPhaseMaterial.hpp:57
static constexpr bool isTemperatureDependent
Definition EclTwoPhaseMaterial.hpp:108
static void capillaryPressures(ContainerT &values, const Params &params, const FluidState &fluidState)
Implements the multiplexer three phase capillary pressure law used by the ECLipse simulator.
Definition EclTwoPhaseMaterial.hpp:145
static constexpr bool implementsTwoPhaseApi
Definition EclTwoPhaseMaterial.hpp:92
static Evaluation Sw(const Params &, const FluidState &)
Definition EclTwoPhaseMaterial.hpp:361
static Evaluation pcgn(const Params &, const FluidState &)
Capillary pressure between the gas and the non-wetting liquid (i.e., oil) phase.
Definition EclTwoPhaseMaterial.hpp:304
static void saturations(ContainerT &, const Params &, const FluidState &)
The inverse of the capillary pressure.
Definition EclTwoPhaseMaterial.hpp:330
static Evaluation krn(const Params &, const FluidState &)
The relative permeability of the non-wetting (i.e., oil) phase.
Definition EclTwoPhaseMaterial.hpp:445
static Evaluation Sn(const Params &, const FluidState &)
The saturation of the non-wetting (i.e., oil) phase.
Definition EclTwoPhaseMaterial.hpp:351
static constexpr bool isSaturationDependent
Definition EclTwoPhaseMaterial.hpp:100
static constexpr bool isPressureDependent
Definition EclTwoPhaseMaterial.hpp:104
static Evaluation krg(const Params &, const FluidState &)
The relative permeability of the gas phase.
Definition EclTwoPhaseMaterial.hpp:425
static Evaluation Sg(const Params &, const FluidState &)
The saturation of the gas phase.
Definition EclTwoPhaseMaterial.hpp:341
static bool updateHysteresis(Params &params, const FluidState &fluidState)
Update the hysteresis parameters after a time step.
Definition EclTwoPhaseMaterial.hpp:460
static constexpr bool isCompositionDependent
Definition EclTwoPhaseMaterial.hpp:112
static Evaluation pcnw(const Params &, const FluidState &)
Capillary pressure between the non-wetting liquid (i.e., oil) and the wetting liquid (i....
Definition EclTwoPhaseMaterial.hpp:320
static void relativePermeabilities(ContainerT &values, const Params &params, const FluidState &fluidState)
The relative permeability of all phases.
Definition EclTwoPhaseMaterial.hpp:383
static Evaluation krw(const Params &, const FluidState &)
The relative permeability of the wetting phase.
Definition EclTwoPhaseMaterial.hpp:435
static constexpr bool implementsTwoPhaseSatApi
Definition EclTwoPhaseMaterial.hpp:96
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30