My Project
Loading...
Searching...
No Matches
EclMultiplexerMaterial.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_MULTIPLEXER_MATERIAL_HPP
28#define OPM_ECL_MULTIPLEXER_MATERIAL_HPP
29
32#include "EclStone1Material.hpp"
33#include "EclStone2Material.hpp"
35
36#include <opm/common/TimingMacros.hpp>
37
38#include <algorithm>
39#include <stdexcept>
40
41namespace Opm {
42
49template <class TraitsT,
50 class GasOilMaterialLawT,
51 class OilWaterMaterialLawT,
52 class GasWaterMaterialLawT,
53 class ParamsT = EclMultiplexerMaterialParams<TraitsT,
54 GasOilMaterialLawT,
55 OilWaterMaterialLawT,
56 GasWaterMaterialLawT> >
57class EclMultiplexerMaterial : public TraitsT
58{
59public:
60 using GasOilMaterialLaw = GasOilMaterialLawT;
61 using OilWaterMaterialLaw = OilWaterMaterialLawT;
62 using GasWaterMaterialLaw = GasWaterMaterialLawT;
63
68
69 // some safety checks
70 static_assert(TraitsT::numPhases == 3,
71 "The number of phases considered by this capillary pressure "
72 "law is always three!");
73 static_assert(GasOilMaterialLaw::numPhases == 2,
74 "The number of phases considered by the gas-oil capillary "
75 "pressure law must be two!");
76 static_assert(OilWaterMaterialLaw::numPhases == 2,
77 "The number of phases considered by the oil-water capillary "
78 "pressure law must be two!");
79 static_assert(GasWaterMaterialLaw::numPhases == 2,
80 "The number of phases considered by the gas-water capillary "
81 "pressure law must be two!");
82 static_assert(std::is_same<typename GasOilMaterialLaw::Scalar,
83 typename OilWaterMaterialLaw::Scalar>::value,
84 "The two two-phase capillary pressure laws must use the same "
85 "type of floating point values.");
86
87 using Traits = TraitsT;
88 using Params = ParamsT;
89 using Scalar = typename Traits::Scalar;
90
91 static constexpr int numPhases = 3;
92 static constexpr int waterPhaseIdx = Traits::wettingPhaseIdx;
93 static constexpr int oilPhaseIdx = Traits::nonWettingPhaseIdx;
94 static constexpr int gasPhaseIdx = Traits::gasPhaseIdx;
95
98 static constexpr bool implementsTwoPhaseApi = false;
99
102 static constexpr bool implementsTwoPhaseSatApi = false;
103
106 static constexpr bool isSaturationDependent = true;
107
110 static constexpr bool isPressureDependent = false;
111
114 static constexpr bool isTemperatureDependent = false;
115
118 static constexpr bool isCompositionDependent = false;
119
134 template <class ContainerT, class FluidState>
135 static void capillaryPressures(ContainerT& values,
136 const Params& params,
137 const FluidState& fluidState)
138 {
139 OPM_TIMEFUNCTION_LOCAL();
140 switch (params.approach()) {
141 case EclMultiplexerApproach::Stone1:
143 params.template getRealParams<EclMultiplexerApproach::Stone1>(),
144 fluidState);
145 break;
146
147 case EclMultiplexerApproach::Stone2:
149 params.template getRealParams<EclMultiplexerApproach::Stone2>(),
150 fluidState);
151 break;
152
153 case EclMultiplexerApproach::Default:
155 params.template getRealParams<EclMultiplexerApproach::Default>(),
156 fluidState);
157 break;
158
159 case EclMultiplexerApproach::TwoPhase:
161 params.template getRealParams<EclMultiplexerApproach::TwoPhase>(),
162 fluidState);
163 break;
164
165 case EclMultiplexerApproach::OnePhase:
166 values[0] = 0.0;
167 break;
168 }
169 }
170
171 /*
172 * Hysteresis parameters for oil-water
173 * @see EclHysteresisTwoPhaseLawParams::pcSwMdc(...)
174 * @see EclHysteresisTwoPhaseLawParams::krnSwMdc(...)
175 * \param params Parameters
176 */
177 static void oilWaterHysteresisParams(Scalar& pcSwMdc,
178 Scalar& krnSwMdc,
179 const Params& params)
180 {
181 OPM_TIMEFUNCTION_LOCAL();
182 switch (params.approach()) {
183 case EclMultiplexerApproach::Stone1:
184 Stone1Material::oilWaterHysteresisParams(pcSwMdc, krnSwMdc,
185 params.template getRealParams<EclMultiplexerApproach::Stone1>());
186 break;
187
188 case EclMultiplexerApproach::Stone2:
189 Stone2Material::oilWaterHysteresisParams(pcSwMdc, krnSwMdc,
190 params.template getRealParams<EclMultiplexerApproach::Stone2>());
191 break;
192
193 case EclMultiplexerApproach::Default:
194 DefaultMaterial::oilWaterHysteresisParams(pcSwMdc, krnSwMdc,
195 params.template getRealParams<EclMultiplexerApproach::Default>());
196 break;
197
198 case EclMultiplexerApproach::TwoPhase:
199 TwoPhaseMaterial::oilWaterHysteresisParams(pcSwMdc, krnSwMdc,
200 params.template getRealParams<EclMultiplexerApproach::TwoPhase>());
201 break;
202
203 case EclMultiplexerApproach::OnePhase:
204 // Do nothing.
205 break;
206 }
207 }
208
209 /*
210 * Hysteresis parameters for oil-water
211 * @see EclHysteresisTwoPhaseLawParams::pcSwMdc(...)
212 * @see EclHysteresisTwoPhaseLawParams::krnSwMdc(...)
213 * \param params Parameters
214 */
215 static void setOilWaterHysteresisParams(const Scalar& pcSwMdc,
216 const Scalar& krnSwMdc,
217 Params& params)
218 {
219 OPM_TIMEFUNCTION_LOCAL();
220 switch (params.approach()) {
221 case EclMultiplexerApproach::Stone1:
222 Stone1Material::setOilWaterHysteresisParams(pcSwMdc, krnSwMdc,
223 params.template getRealParams<EclMultiplexerApproach::Stone1>());
224 break;
225
226 case EclMultiplexerApproach::Stone2:
227 Stone2Material::setOilWaterHysteresisParams(pcSwMdc, krnSwMdc,
228 params.template getRealParams<EclMultiplexerApproach::Stone2>());
229 break;
230
231 case EclMultiplexerApproach::Default:
232 DefaultMaterial::setOilWaterHysteresisParams(pcSwMdc, krnSwMdc,
233 params.template getRealParams<EclMultiplexerApproach::Default>());
234 break;
235
236 case EclMultiplexerApproach::TwoPhase:
237 TwoPhaseMaterial::setOilWaterHysteresisParams(pcSwMdc, krnSwMdc,
238 params.template getRealParams<EclMultiplexerApproach::TwoPhase>());
239 break;
240
241 case EclMultiplexerApproach::OnePhase:
242 // Do nothing.
243 break;
244 }
245 }
246
247 /*
248 * Hysteresis parameters for gas-oil
249 * @see EclHysteresisTwoPhaseLawParams::pcSwMdc(...)
250 * @see EclHysteresisTwoPhaseLawParams::krnSwMdc(...)
251 * \param params Parameters
252 */
253 static void gasOilHysteresisParams(Scalar& pcSwMdc,
254 Scalar& krnSwMdc,
255 const Params& params)
256 {
257 OPM_TIMEFUNCTION_LOCAL();
258 switch (params.approach()) {
259 case EclMultiplexerApproach::Stone1:
260 Stone1Material::gasOilHysteresisParams(pcSwMdc, krnSwMdc,
261 params.template getRealParams<EclMultiplexerApproach::Stone1>());
262 break;
263
264 case EclMultiplexerApproach::Stone2:
265 Stone2Material::gasOilHysteresisParams(pcSwMdc, krnSwMdc,
266 params.template getRealParams<EclMultiplexerApproach::Stone2>());
267 break;
268
269 case EclMultiplexerApproach::Default:
270 DefaultMaterial::gasOilHysteresisParams(pcSwMdc, krnSwMdc,
271 params.template getRealParams<EclMultiplexerApproach::Default>());
272 break;
273
274 case EclMultiplexerApproach::TwoPhase:
275 TwoPhaseMaterial::gasOilHysteresisParams(pcSwMdc, krnSwMdc,
276 params.template getRealParams<EclMultiplexerApproach::TwoPhase>());
277 break;
278
279 case EclMultiplexerApproach::OnePhase:
280 // Do nothing.
281 break;
282 }
283 }
284
285 static Scalar trappedGasSaturation(const Params& params)
286 {
287 OPM_TIMEFUNCTION_LOCAL();
288 switch (params.approach()) {
289 case EclMultiplexerApproach::Stone1:
290 return params.template getRealParams<EclMultiplexerApproach::Stone1>().gasOilParams().SnTrapped();
291 case EclMultiplexerApproach::Stone2:
292 return params.template getRealParams<EclMultiplexerApproach::Stone2>().gasOilParams().SnTrapped();
293 case EclMultiplexerApproach::Default:
294 return params.template getRealParams<EclMultiplexerApproach::Default>().gasOilParams().SnTrapped();
295 case EclMultiplexerApproach::TwoPhase:
296 if(params.template getRealParams<EclMultiplexerApproach::TwoPhase>().approach() == EclTwoPhaseApproach::GasOil)
297 return params.template getRealParams<EclMultiplexerApproach::TwoPhase>().gasOilParams().SnTrapped();
298 if(params.template getRealParams<EclMultiplexerApproach::TwoPhase>().approach() == EclTwoPhaseApproach::GasWater)
299 return params.template getRealParams<EclMultiplexerApproach::TwoPhase>().gasWaterParams().SnTrapped();
300 return 0.0; // oil-water case
301 case EclMultiplexerApproach::OnePhase:
302 return 0.0;
303 }
304 return 0.0;
305 }
306
307 /*
308 * Hysteresis parameters for gas-oil
309 * @see EclHysteresisTwoPhaseLawParams::pcSwMdc(...)
310 * @see EclHysteresisTwoPhaseLawParams::krnSwMdc(...)
311 * \param params Parameters
312 */
313 static void setGasOilHysteresisParams(const Scalar& pcSwMdc,
314 const Scalar& krnSwMdc,
315 Params& params)
316 {
317 OPM_TIMEFUNCTION_LOCAL();
318 switch (params.approach()) {
319 case EclMultiplexerApproach::Stone1:
320 Stone1Material::setGasOilHysteresisParams(pcSwMdc, krnSwMdc,
321 params.template getRealParams<EclMultiplexerApproach::Stone1>());
322 break;
323
324 case EclMultiplexerApproach::Stone2:
325 Stone2Material::setGasOilHysteresisParams(pcSwMdc, krnSwMdc,
326 params.template getRealParams<EclMultiplexerApproach::Stone2>());
327 break;
328
329 case EclMultiplexerApproach::Default:
330 DefaultMaterial::setGasOilHysteresisParams(pcSwMdc, krnSwMdc,
331 params.template getRealParams<EclMultiplexerApproach::Default>());
332 break;
333
334 case EclMultiplexerApproach::TwoPhase:
335 TwoPhaseMaterial::setGasOilHysteresisParams(pcSwMdc, krnSwMdc,
336 params.template getRealParams<EclMultiplexerApproach::TwoPhase>());
337 break;
338
339 case EclMultiplexerApproach::OnePhase:
340 // Do nothing.
341 break;
342 }
343 }
344
354 template <class FluidState, class Evaluation = typename FluidState::Scalar>
355 static Evaluation pcgn(const Params& /* params */,
356 const FluidState& /* fs */)
357 {
358 throw std::logic_error("Not implemented: pcgn()");
359 }
360
370 template <class FluidState, class Evaluation = typename FluidState::Scalar>
371 static Evaluation pcnw(const Params& /* params */,
372 const FluidState& /* fs */)
373 {
374 throw std::logic_error("Not implemented: pcnw()");
375 }
376
380 template <class ContainerT, class FluidState>
381 static void saturations(ContainerT& /* values */,
382 const Params& /* params */,
383 const FluidState& /* fs */)
384 {
385 throw std::logic_error("Not implemented: saturations()");
386 }
387
391 template <class FluidState, class Evaluation = typename FluidState::Scalar>
392 static Evaluation Sg(const Params& /* params */,
393 const FluidState& /* fluidState */)
394 {
395 throw std::logic_error("Not implemented: Sg()");
396 }
397
401 template <class FluidState, class Evaluation = typename FluidState::Scalar>
402 static Evaluation Sn(const Params& /* params */,
403 const FluidState& /* fluidState */)
404 {
405 throw std::logic_error("Not implemented: Sn()");
406 }
407
411 template <class FluidState, class Evaluation = typename FluidState::Scalar>
412 static Evaluation Sw(const Params& /* params */,
413 const FluidState& /* fluidState */)
414 {
415 throw std::logic_error("Not implemented: Sw()");
416 }
417
433 template <class ContainerT, class FluidState>
434 static void relativePermeabilities(ContainerT& values,
435 const Params& params,
436 const FluidState& fluidState)
437 {
438 OPM_TIMEFUNCTION_LOCAL();
439 switch (params.approach()) {
440 case EclMultiplexerApproach::Stone1:
442 params.template getRealParams<EclMultiplexerApproach::Stone1>(),
443 fluidState);
444 break;
445
446 case EclMultiplexerApproach::Stone2:
448 params.template getRealParams<EclMultiplexerApproach::Stone2>(),
449 fluidState);
450 break;
451
452 case EclMultiplexerApproach::Default:
454 params.template getRealParams<EclMultiplexerApproach::Default>(),
455 fluidState);
456 break;
457
458 case EclMultiplexerApproach::TwoPhase:
460 params.template getRealParams<EclMultiplexerApproach::TwoPhase>(),
461 fluidState);
462 break;
463
464 case EclMultiplexerApproach::OnePhase:
465 values[0] = 1.0;
466 break;
467
468 default:
469 throw std::logic_error("Not implemented: relativePermeabilities() option for unknown EclMultiplexerApproach (="
470 + std::to_string(static_cast<int>(params.approach())) + ")");
471 }
472 }
473
477 template <class Evaluation, class FluidState>
478 static Evaluation relpermOilInOilGasSystem(const Params& params,
479 const FluidState& fluidState)
480 {
481 OPM_TIMEFUNCTION_LOCAL();
482 switch (params.approach()) {
483 case EclMultiplexerApproach::Stone1:
484 return Stone1Material::template relpermOilInOilGasSystem<Evaluation>
485 (params.template getRealParams<EclMultiplexerApproach::Stone1>(),
486 fluidState);
487
488 case EclMultiplexerApproach::Stone2:
489 return Stone2Material::template relpermOilInOilGasSystem<Evaluation>
490 (params.template getRealParams<EclMultiplexerApproach::Stone2>(),
491 fluidState);
492
493 case EclMultiplexerApproach::Default:
494 return DefaultMaterial::template relpermOilInOilGasSystem<Evaluation>
495 (params.template getRealParams<EclMultiplexerApproach::Default>(),
496 fluidState);
497
498 default:
499 throw std::logic_error {
500 "relpermOilInOilGasSystem() is specific to three phases"
501 };
502 }
503 }
504
508 template <class Evaluation, class FluidState>
509 static Evaluation relpermOilInOilWaterSystem(const Params& params,
510 const FluidState& fluidState)
511 {
512 OPM_TIMEFUNCTION_LOCAL();
513 switch (params.approach()) {
514 case EclMultiplexerApproach::Stone1:
515 return Stone1Material::template relpermOilInOilWaterSystem<Evaluation>
516 (params.template getRealParams<EclMultiplexerApproach::Stone1>(),
517 fluidState);
518
519 case EclMultiplexerApproach::Stone2:
520 return Stone2Material::template relpermOilInOilWaterSystem<Evaluation>
521 (params.template getRealParams<EclMultiplexerApproach::Stone2>(),
522 fluidState);
523
524 case EclMultiplexerApproach::Default:
525 return DefaultMaterial::template relpermOilInOilWaterSystem<Evaluation>
526 (params.template getRealParams<EclMultiplexerApproach::Default>(),
527 fluidState);
528
529 default:
530 throw std::logic_error {
531 "relpermOilInOilWaterSystem() is specific to three phases"
532 };
533 }
534 }
535
539 template <class FluidState, class Evaluation = typename FluidState::Scalar>
540 static Evaluation krg(const Params& /* params */,
541 const FluidState& /* fluidState */)
542 {
543 throw std::logic_error("Not implemented: krg()");
544 }
545
549 template <class FluidState, class Evaluation = typename FluidState::Scalar>
550 static Evaluation krw(const Params& /* params */,
551 const FluidState& /* fluidState */)
552 {
553 throw std::logic_error("Not implemented: krw()");
554 }
555
559 template <class FluidState, class Evaluation = typename FluidState::Scalar>
560 static Evaluation krn(const Params& /* params */,
561 const FluidState& /* fluidState */)
562 {
563 throw std::logic_error("Not implemented: krn()");
564 }
565
566
574 template <class FluidState>
575 static bool updateHysteresis(Params& params, const FluidState& fluidState)
576 {
577 OPM_TIMEFUNCTION_LOCAL();
578 switch (params.approach()) {
579 case EclMultiplexerApproach::Stone1:
580 return Stone1Material::updateHysteresis(params.template getRealParams<EclMultiplexerApproach::Stone1>(),
581 fluidState);
582 break;
583
584 case EclMultiplexerApproach::Stone2:
585 return Stone2Material::updateHysteresis(params.template getRealParams<EclMultiplexerApproach::Stone2>(),
586 fluidState);
587 break;
588
589 case EclMultiplexerApproach::Default:
590 return DefaultMaterial::updateHysteresis(params.template getRealParams<EclMultiplexerApproach::Default>(),
591 fluidState);
592 break;
593
594 case EclMultiplexerApproach::TwoPhase:
595 return TwoPhaseMaterial::updateHysteresis(params.template getRealParams<EclMultiplexerApproach::TwoPhase>(),
596 fluidState);
597 break;
598 case EclMultiplexerApproach::OnePhase:
599 return false;
600 break;
601 }
602 return false;
603 }
604};
605
606} // namespace Opm
607
608#endif
Implements the default three phase capillary pressure law used by the ECLipse simulator.
Multiplexer implementation for the parameters required by the multiplexed three-phase material law.
Implements the second phase capillary pressure/relperm law suggested by Stone as used by the ECLipse ...
Implements the second phase capillary pressure/relperm law suggested by Stone as used by the ECLipse ...
Implements a multiplexer class that provides ECL saturation functions for twophase simulations.
Implements the default three phase capillary pressure law used by the ECLipse simulator.
Definition EclDefaultMaterial.hpp:61
static void capillaryPressures(ContainerT &values, const Params &params, const FluidState &state)
Implements the default three phase capillary pressure law used by the ECLipse simulator.
Definition EclDefaultMaterial.hpp:136
static bool updateHysteresis(Params &params, const FluidState &fluidState)
Update the hysteresis parameters after a time step.
Definition EclDefaultMaterial.hpp:439
static void relativePermeabilities(ContainerT &values, const Params &params, const FluidState &fluidState)
The relative permeability of all phases.
Definition EclDefaultMaterial.hpp:319
Implements a multiplexer class that provides all three phase capillary pressure laws used by the ECLi...
Definition EclMultiplexerMaterial.hpp:58
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 EclMultiplexerMaterial.hpp:135
static Evaluation relpermOilInOilWaterSystem(const Params &params, const FluidState &fluidState)
The relative permeability of oil in oil/water system.
Definition EclMultiplexerMaterial.hpp:509
static Evaluation krw(const Params &, const FluidState &)
The relative permeability of the wetting phase.
Definition EclMultiplexerMaterial.hpp:550
static Evaluation Sg(const Params &, const FluidState &)
The saturation of the gas phase.
Definition EclMultiplexerMaterial.hpp:392
static Evaluation relpermOilInOilGasSystem(const Params &params, const FluidState &fluidState)
The relative permeability of oil in oil/gas system.
Definition EclMultiplexerMaterial.hpp:478
static Evaluation pcgn(const Params &, const FluidState &)
Capillary pressure between the gas and the non-wetting liquid (i.e., oil) phase.
Definition EclMultiplexerMaterial.hpp:355
static Evaluation pcnw(const Params &, const FluidState &)
Capillary pressure between the non-wetting liquid (i.e., oil) and the wetting liquid (i....
Definition EclMultiplexerMaterial.hpp:371
static constexpr bool isPressureDependent
Specify whether the quantities defined by this material law are dependent on the absolute pressure.
Definition EclMultiplexerMaterial.hpp:110
static constexpr bool implementsTwoPhaseSatApi
Specify whether this material law implements the two-phase convenience API which only depends on the ...
Definition EclMultiplexerMaterial.hpp:102
static constexpr bool isSaturationDependent
Specify whether the quantities defined by this material law are saturation dependent.
Definition EclMultiplexerMaterial.hpp:106
static constexpr bool isTemperatureDependent
Specify whether the quantities defined by this material law are temperature dependent.
Definition EclMultiplexerMaterial.hpp:114
static constexpr bool implementsTwoPhaseApi
Specify whether this material law implements the two-phase convenience API.
Definition EclMultiplexerMaterial.hpp:98
static void relativePermeabilities(ContainerT &values, const Params &params, const FluidState &fluidState)
The relative permeability of all phases.
Definition EclMultiplexerMaterial.hpp:434
static Evaluation Sn(const Params &, const FluidState &)
The saturation of the non-wetting (i.e., oil) phase.
Definition EclMultiplexerMaterial.hpp:402
static Evaluation krg(const Params &, const FluidState &)
The relative permeability of the gas phase.
Definition EclMultiplexerMaterial.hpp:540
static constexpr bool isCompositionDependent
Specify whether the quantities defined by this material law are dependent on the phase composition.
Definition EclMultiplexerMaterial.hpp:118
static void saturations(ContainerT &, const Params &, const FluidState &)
The inverse of the capillary pressure.
Definition EclMultiplexerMaterial.hpp:381
static Evaluation krn(const Params &, const FluidState &)
The relative permeability of the non-wetting (i.e., oil) phase.
Definition EclMultiplexerMaterial.hpp:560
static Evaluation Sw(const Params &, const FluidState &)
The saturation of the wetting (i.e., water) phase.
Definition EclMultiplexerMaterial.hpp:412
static bool updateHysteresis(Params &params, const FluidState &fluidState)
Update the hysteresis parameters after a time step.
Definition EclMultiplexerMaterial.hpp:575
Implements the second phase capillary pressure/relperm law suggested by Stone as used by the ECLipse ...
Definition EclStone1Material.hpp:60
static void relativePermeabilities(ContainerT &values, const Params &params, const FluidState &fluidState)
The relative permeability of all phases.
Definition EclStone1Material.hpp:313
static void capillaryPressures(ContainerT &values, const Params &params, const FluidState &state)
Implements the default three phase capillary pressure law used by the ECLipse simulator.
Definition EclStone1Material.hpp:135
static bool updateHysteresis(Params &params, const FluidState &fluidState)
Update the hysteresis parameters after a time step.
Definition EclStone1Material.hpp:420
Implements the second phase capillary pressure/relperm law suggested by Stone as used by the ECLipse ...
Definition EclStone2Material.hpp:61
static bool updateHysteresis(Params &params, const FluidState &fluidState)
Update the hysteresis parameters after a time step.
Definition EclStone2Material.hpp:404
static void capillaryPressures(ContainerT &values, const Params &params, const FluidState &state)
Implements the default three phase capillary pressure law used by the ECLipse simulator.
Definition EclStone2Material.hpp:136
static void relativePermeabilities(ContainerT &values, const Params &params, const FluidState &fluidState)
The relative permeability of all phases.
Definition EclStone2Material.hpp:314
Implements a multiplexer class that provides ECL saturation functions for twophase simulations.
Definition EclTwoPhaseMaterial.hpp:57
static void relativePermeabilities(ContainerT &values, const Params &params, const FluidState &fluidState)
The relative permeability of all phases.
Definition EclTwoPhaseMaterial.hpp:343
static bool updateHysteresis(Params &params, const FluidState &fluidState)
Update the hysteresis parameters after a time step.
Definition EclTwoPhaseMaterial.hpp:420
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
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30