NewTranFluxModule.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*/
31#ifndef OPM_NEWTRAN_FLUX_MODULE_HPP
32#define OPM_NEWTRAN_FLUX_MODULE_HPP
33
34#include <dune/common/fvector.hh>
35#include <dune/common/fmatrix.hh>
36
37#include <opm/common/OpmLog/OpmLog.hpp>
38
39#include <opm/input/eclipse/EclipseState/Grid/FaceDir.hpp>
40
41#include <opm/material/common/MathToolbox.hpp>
42#include <opm/material/common/Valgrind.hpp>
43
44#include <opm/models/discretization/common/fvbaseproperties.hh>
45#include <opm/models/blackoil/blackoilproperties.hh>
46#include <opm/models/utils/signum.hh>
47
48#include <array>
49
50namespace Opm {
51
52template <class TypeTag>
53class NewTranIntensiveQuantities;
54
55template <class TypeTag>
56class NewTranExtensiveQuantities;
57
58template <class TypeTag>
59class NewTranBaseProblem;
60
65template <class TypeTag>
67{
71
75 static void registerParameters()
76 { }
77};
78
84template <class TypeTag>
86{ };
87
92template <class TypeTag>
94{
95 using ElementContext = GetPropType<TypeTag, Properties::ElementContext>;
96protected:
97 void update_(const ElementContext&, unsigned, unsigned)
98 { }
99};
100
105template <class TypeTag>
107{
108 using Implementation = GetPropType<TypeTag, Properties::ExtensiveQuantities>;
109
110 using IntensiveQuantities = GetPropType<TypeTag, Properties::IntensiveQuantities>;
111 using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
112 using ElementContext = GetPropType<TypeTag, Properties::ElementContext>;
113 using Scalar = GetPropType<TypeTag, Properties::Scalar>;
114 using Evaluation = GetPropType<TypeTag, Properties::Evaluation>;
115 using GridView = GetPropType<TypeTag, Properties::GridView>;
116 using MaterialLaw = GetPropType<TypeTag, Properties::MaterialLaw>;
117
118 enum { dimWorld = GridView::dimensionworld };
119 enum { gasPhaseIdx = FluidSystem::gasPhaseIdx };
120 enum { numPhases = FluidSystem::numPhases };
121 enum { enableSolvent = getPropValue<TypeTag, Properties::EnableSolvent>() };
122 enum { enableExtbo = getPropValue<TypeTag, Properties::EnableExtbo>() };
123 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
124
125 using Toolbox = MathToolbox<Evaluation>;
126 using DimVector = Dune::FieldVector<Scalar, dimWorld>;
127 using EvalDimVector = Dune::FieldVector<Evaluation, dimWorld>;
128 using DimMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
129
130public:
134 const DimMatrix& intrinsicPermeability() const
135 {
136 throw std::invalid_argument("The ECL transmissibility module does not provide an explicit intrinsic permeability");
137 }
138
145 const EvalDimVector& potentialGrad(unsigned) const
146 {
147 throw std::invalid_argument("The ECL transmissibility module does not provide explicit potential gradients");
148 }
149
156 const Evaluation& pressureDifference(unsigned phaseIdx) const
157 { return pressureDifference_[phaseIdx]; }
158
165 const EvalDimVector& filterVelocity(unsigned) const
166 {
167 throw std::invalid_argument("The ECL transmissibility module does not provide explicit filter velocities");
168 }
169
179 const Evaluation& volumeFlux(unsigned phaseIdx) const
180 { return volumeFlux_[phaseIdx]; }
181
182protected:
190 unsigned upstreamIndex_(unsigned phaseIdx) const
191 {
192 assert(phaseIdx < numPhases);
193
194 return upIdx_[phaseIdx];
195 }
196
204 unsigned downstreamIndex_(unsigned phaseIdx) const
205 {
206 assert(phaseIdx < numPhases);
207
208 return dnIdx_[phaseIdx];
209 }
210
211 void updateSolvent(const ElementContext& elemCtx, unsigned scvfIdx, unsigned timeIdx)
212 { asImp_().updateVolumeFluxTrans(elemCtx, scvfIdx, timeIdx); }
213
214 void updatePolymer(const ElementContext& elemCtx, unsigned scvfIdx, unsigned timeIdx)
215 { asImp_().updateShearMultipliers(elemCtx, scvfIdx, timeIdx); }
216
217public:
218
219 static void volumeAndPhasePressureDifferences(std::array<short, numPhases>& upIdx,
220 std::array<short, numPhases>& dnIdx,
221 Evaluation (&volumeFlux)[numPhases],
222 Evaluation (&pressureDifferences)[numPhases],
223 const ElementContext& elemCtx,
224 unsigned scvfIdx,
225 unsigned timeIdx)
226 {
227 const auto& problem = elemCtx.problem();
228 const auto& stencil = elemCtx.stencil(timeIdx);
229 const auto& scvf = stencil.interiorFace(scvfIdx);
230 unsigned interiorDofIdx = scvf.interiorIndex();
231 unsigned exteriorDofIdx = scvf.exteriorIndex();
232
233 assert(interiorDofIdx != exteriorDofIdx);
234
235 unsigned I = stencil.globalSpaceIndex(interiorDofIdx);
236 unsigned J = stencil.globalSpaceIndex(exteriorDofIdx);
237 Scalar trans = problem.transmissibility(elemCtx, interiorDofIdx, exteriorDofIdx);
238 Scalar faceArea = scvf.area();
239 Scalar thpres = problem.thresholdPressure(I, J);
240
241 // estimate the gravity correction: for performance reasons we use a simplified
242 // approach for this flux module that assumes that gravity is constant and always
243 // acts into the downwards direction. (i.e., no centrifuge experiments, sorry.)
244 Scalar g = elemCtx.problem().gravity()[dimWorld - 1];
245
246 const auto& intQuantsIn = elemCtx.intensiveQuantities(interiorDofIdx, timeIdx);
247 const auto& intQuantsEx = elemCtx.intensiveQuantities(exteriorDofIdx, timeIdx);
248
249 // this is quite hacky because the dune grid interface does not provide a
250 // cellCenterDepth() method (so we ask the problem to provide it). The "good"
251 // solution would be to take the Z coordinate of the element centroids, but since
252 // ECL seems to like to be inconsistent on that front, it needs to be done like
253 // here...
254 Scalar zIn = problem.dofCenterDepth(elemCtx, interiorDofIdx, timeIdx);
255 Scalar zEx = problem.dofCenterDepth(elemCtx, exteriorDofIdx, timeIdx);
256
257 // the distances from the DOF's depths. (i.e., the additional depth of the
258 // exterior DOF)
259 Scalar distZ = zIn - zEx;
260
261 Scalar Vin = elemCtx.dofVolume(interiorDofIdx, /*timeIdx=*/0);
262 Scalar Vex = elemCtx.dofVolume(exteriorDofIdx, /*timeIdx=*/0);
263
264 for (unsigned phaseIdx=0; phaseIdx < numPhases; phaseIdx++) {
265 if (!FluidSystem::phaseIsActive(phaseIdx))
266 continue;
267 calculatePhasePressureDiff_(upIdx[phaseIdx],
268 dnIdx[phaseIdx],
269 pressureDifferences[phaseIdx],
270 intQuantsIn,
271 intQuantsEx,
272 phaseIdx,//input
273 interiorDofIdx,//input
274 exteriorDofIdx,//input
275 Vin,
276 Vex,
277 I,
278 J,
279 distZ*g,
280 thpres);
281 if (pressureDifferences[phaseIdx] == 0) {
282 volumeFlux[phaseIdx] = 0.0;
283 continue;
284 }
285
286 const bool upwindIsInterior = (static_cast<unsigned>(upIdx[phaseIdx]) == interiorDofIdx);
287 const IntensiveQuantities& up = upwindIsInterior ? intQuantsIn : intQuantsEx;
288 // Use arithmetic average (more accurate with harmonic, but that requires recomputing the transmissbility)
289 const Evaluation transMult = (intQuantsIn.rockCompTransMultiplier() + Toolbox::value(intQuantsEx.rockCompTransMultiplier()))/2;
290
291 const auto& materialLawManager = problem.materialLawManager();
292 FaceDir::DirEnum facedir = FaceDir::DirEnum::Unknown;
293 if (materialLawManager->hasDirectionalRelperms()) {
294 facedir = scvf.faceDirFromDirId(); // direction (X, Y, or Z) of the face
295 }
296 if (upwindIsInterior)
297 volumeFlux[phaseIdx] =
298 pressureDifferences[phaseIdx]*up.mobility(phaseIdx, facedir)*transMult*(-trans/faceArea);
299 else
300 volumeFlux[phaseIdx] =
301 pressureDifferences[phaseIdx]*
302 (Toolbox::value(up.mobility(phaseIdx, facedir))*transMult*(-trans/faceArea));
303 }
304 }
305
306 template<class EvalType>
307 static void calculatePhasePressureDiff_(short& upIdx,
308 short& dnIdx,
309 EvalType& pressureDifference,
310 const IntensiveQuantities& intQuantsIn,
311 const IntensiveQuantities& intQuantsEx,
312 const unsigned phaseIdx,
313 const unsigned interiorDofIdx,
314 const unsigned exteriorDofIdx,
315 const Scalar Vin,
316 const Scalar Vex,
317 const unsigned globalIndexIn,
318 const unsigned globalIndexEx,
319 const Scalar distZg,
320 const Scalar thpres
321 )
322 {
323
324 // check shortcut: if the mobility of the phase is zero in the interior as
325 // well as the exterior DOF, we can skip looking at the phase.
326 if (intQuantsIn.mobility(phaseIdx) <= 0.0 &&
327 intQuantsEx.mobility(phaseIdx) <= 0.0)
328 {
329 upIdx = interiorDofIdx;
330 dnIdx = exteriorDofIdx;
331 pressureDifference = 0.0;
332 return;
333 }
334
335 // do the gravity correction: compute the hydrostatic pressure for the
336 // external at the depth of the internal one
337 const Evaluation& rhoIn = intQuantsIn.fluidState().density(phaseIdx);
338 Scalar rhoEx = Toolbox::value(intQuantsEx.fluidState().density(phaseIdx));
339 Evaluation rhoAvg = (rhoIn + rhoEx)/2;
340
341 const Evaluation& pressureInterior = intQuantsIn.fluidState().pressure(phaseIdx);
342 Evaluation pressureExterior = Toolbox::value(intQuantsEx.fluidState().pressure(phaseIdx));
343 if (enableExtbo) // added stability; particulary useful for solvent migrating in pure water
344 // where the solvent fraction displays a 0/1 behaviour ...
345 pressureExterior += Toolbox::value(rhoAvg)*(distZg);
346 else
347 pressureExterior += rhoAvg*(distZg);
348
349 pressureDifference = pressureExterior - pressureInterior;
350
351 // decide the upstream index for the phase. for this we make sure that the
352 // degree of freedom which is regarded upstream if both pressures are equal
353 // is always the same: if the pressure is equal, the DOF with the lower
354 // global index is regarded to be the upstream one.
355 if (pressureDifference > 0.0) {
356 upIdx = exteriorDofIdx;
357 dnIdx = interiorDofIdx;
358 }
359 else if (pressureDifference < 0.0) {
360 upIdx = interiorDofIdx;
361 dnIdx = exteriorDofIdx;
362 }
363 else {
364 // if the pressure difference is zero, we chose the DOF which has the
365 // larger volume associated to it as upstream DOF
366 if (Vin > Vex) {
367 upIdx = interiorDofIdx;
368 dnIdx = exteriorDofIdx;
369 }
370 else if (Vin < Vex) {
371 upIdx = exteriorDofIdx;
372 dnIdx = interiorDofIdx;
373 }
374 else {
375 assert(Vin == Vex);
376 // if the volumes are also equal, we pick the DOF which exhibits the
377 // smaller global index
378 if (globalIndexIn < globalIndexEx) {
379 upIdx = interiorDofIdx;
380 dnIdx = exteriorDofIdx;
381 }
382 else {
383 upIdx = exteriorDofIdx;
384 dnIdx = interiorDofIdx;
385 }
386 }
387 }
388
389 // apply the threshold pressure for the intersection. note that the concept
390 // of threshold pressure is a quite big hack that only makes sense for ECL
391 // datasets. (and even there, its physical justification is quite
392 // questionable IMO.)
393 if (thpres > 0.0) {
394 if (std::abs(Toolbox::value(pressureDifference)) > thpres) {
395 if (pressureDifference < 0.0)
396 pressureDifference += thpres;
397 else
398 pressureDifference -= thpres;
399 }
400 else {
401 pressureDifference = 0.0;
402 }
403 }
404 }
405
406protected:
410 void calculateGradients_(const ElementContext& elemCtx, unsigned scvfIdx, unsigned timeIdx)
411 {
412 Valgrind::SetUndefined(*this);
413
414 volumeAndPhasePressureDifferences(upIdx_ , dnIdx_, volumeFlux_, pressureDifference_, elemCtx, scvfIdx, timeIdx);
415 }
416
420 template <class FluidState>
421 void calculateBoundaryGradients_(const ElementContext& elemCtx,
422 unsigned scvfIdx,
423 unsigned timeIdx,
424 const FluidState& exFluidState)
425 {
426 const auto& scvf = elemCtx.stencil(timeIdx).boundaryFace(scvfIdx);
427 const Scalar faceArea = scvf.area();
428 const Scalar zEx = scvf.integrationPos()[dimWorld - 1];
429 const auto& problem = elemCtx.problem();
430 const unsigned globalSpaceIdx = elemCtx.globalSpaceIndex(0, timeIdx);
431 const auto& intQuantsIn = elemCtx.intensiveQuantities(0, timeIdx);
432
434 globalSpaceIdx,
435 intQuantsIn,
436 scvfIdx,
437 faceArea,
438 zEx,
439 exFluidState,
440 upIdx_,
441 dnIdx_,
442 volumeFlux_,
443 pressureDifference_);
444
445 // Treating solvent here and not in the static method, since that would require more
446 // extensive refactoring. It means that the TpfaLinearizer will not support bcs for solvent until this is
447 // addressed.
448 if constexpr (enableSolvent) {
449 if (upIdx_[gasPhaseIdx] == 0) {
450 const Scalar trans = problem.transmissibilityBoundary(globalSpaceIdx, scvfIdx);
451 const Scalar transModified = trans * Toolbox::value(intQuantsIn.rockCompTransMultiplier());
452 const auto solventFlux = pressureDifference_[gasPhaseIdx] * intQuantsIn.mobility(gasPhaseIdx) * (-transModified/faceArea);
453 asImp_().setSolventVolumeFlux(solventFlux);
454 } else {
455 asImp_().setSolventVolumeFlux(0.0);
456 }
457 }
458 }
459
460public:
464 template <class Problem, class FluidState, class EvaluationContainer>
465 static void calculateBoundaryGradients_(const Problem& problem,
466 const unsigned globalSpaceIdx,
467 const IntensiveQuantities& intQuantsIn,
468 const unsigned bfIdx,
469 const double faceArea,
470 const double zEx,
471 const FluidState& exFluidState,
472 std::array<short, numPhases>& upIdx,
473 std::array<short, numPhases>& dnIdx,
474 EvaluationContainer& volumeFlux,
475 EvaluationContainer& pressureDifference)
476 {
477
478 bool enableBoundaryMassFlux = problem.nonTrivialBoundaryConditions();
479 if (!enableBoundaryMassFlux)
480 return;
481
482 Scalar trans = problem.transmissibilityBoundary(globalSpaceIdx, bfIdx);
483
484 // estimate the gravity correction: for performance reasons we use a simplified
485 // approach for this flux module that assumes that gravity is constant and always
486 // acts into the downwards direction. (i.e., no centrifuge experiments, sorry.)
487 Scalar g = problem.gravity()[dimWorld - 1];
488
489 // this is quite hacky because the dune grid interface does not provide a
490 // cellCenterDepth() method (so we ask the problem to provide it). The "good"
491 // solution would be to take the Z coordinate of the element centroids, but since
492 // ECL seems to like to be inconsistent on that front, it needs to be done like
493 // here...
494 Scalar zIn = problem.dofCenterDepth(globalSpaceIdx);
495
496 // the distances from the DOF's depths. (i.e., the additional depth of the
497 // exterior DOF)
498 Scalar distZ = zIn - zEx;
499
500 for (unsigned phaseIdx=0; phaseIdx < numPhases; phaseIdx++) {
501 if (!FluidSystem::phaseIsActive(phaseIdx))
502 continue;
503
504 // do the gravity correction: compute the hydrostatic pressure for the
505 // integration position
506 const Evaluation& rhoIn = intQuantsIn.fluidState().density(phaseIdx);
507 const auto& rhoEx = exFluidState.density(phaseIdx);
508 Evaluation rhoAvg = (rhoIn + rhoEx)/2;
509
510 const Evaluation& pressureInterior = intQuantsIn.fluidState().pressure(phaseIdx);
511 Evaluation pressureExterior = exFluidState.pressure(phaseIdx);
512 pressureExterior += rhoAvg*(distZ*g);
513
514 pressureDifference[phaseIdx] = pressureExterior - pressureInterior;
515
516 // decide the upstream index for the phase. for this we make sure that the
517 // degree of freedom which is regarded upstream if both pressures are equal
518 // is always the same: if the pressure is equal, the DOF with the lower
519 // global index is regarded to be the upstream one.
520 const unsigned interiorDofIdx = 0; // Valid only for cell-centered FV.
521 if (pressureDifference[phaseIdx] > 0.0) {
522 upIdx[phaseIdx] = -1;
523 dnIdx[phaseIdx] = interiorDofIdx;
524 }
525 else {
526 upIdx[phaseIdx] = interiorDofIdx;
527 dnIdx[phaseIdx] = -1;
528 }
529
530 Evaluation transModified = trans;
531
532 if (upIdx[phaseIdx] == interiorDofIdx) {
533
534 // this is slightly hacky because in the automatic differentiation case, it
535 // only works for the element centered finite volume method.
536 const auto& up = intQuantsIn;
537
538 // deal with water induced rock compaction
539 const Scalar transMult = Toolbox::value(up.rockCompTransMultiplier());
540 transModified *= transMult;
541
542 volumeFlux[phaseIdx] =
543 pressureDifference[phaseIdx]*up.mobility(phaseIdx)*(-transModified/faceArea);
544 }
545 else {
546 // compute the phase mobility using the material law parameters of the
547 // interior element. TODO: this could probably be done more efficiently
548 const auto& matParams = problem.materialLawParams(globalSpaceIdx);
549 std::array<typename FluidState::Scalar,numPhases> kr;
550 MaterialLaw::relativePermeabilities(kr, matParams, exFluidState);
551
552 const auto& mob = kr[phaseIdx]/exFluidState.viscosity(phaseIdx);
553 volumeFlux[phaseIdx] =
554 pressureDifference[phaseIdx]*mob*(-transModified/faceArea);
555 }
556 }
557 }
558
559protected:
560
564 void calculateFluxes_(const ElementContext&, unsigned, unsigned)
565 { }
566
567 void calculateBoundaryFluxes_(const ElementContext&, unsigned, unsigned)
568 {}
569
570private:
571 Implementation& asImp_()
572 { return *static_cast<Implementation*>(this); }
573
574 const Implementation& asImp_() const
575 { return *static_cast<const Implementation*>(this); }
576
577 // the volumetric flux of all phases [m^3/s]
578 Evaluation volumeFlux_[numPhases];
579
580 // the difference in effective pressure between the exterior and the interior degree
581 // of freedom [Pa]
582 Evaluation pressureDifference_[numPhases];
583
584 // the local indices of the interior and exterior degrees of freedom
585 std::array<short, numPhases> upIdx_;
586 std::array<short, numPhases> dnIdx_;
587 };
588
589} // namespace Opm
590
591#endif // OPM_NEWTRAN_FLUX_MODULE_HPP
Provides the defaults for the parameters required by the transmissibility based volume flux calculati...
Definition: NewTranFluxModule.hpp:86
Provides the ECL flux module.
Definition: NewTranFluxModule.hpp:107
const Evaluation & volumeFlux(unsigned phaseIdx) const
Return the volume flux of a fluid phase at the face's integration point .
Definition: NewTranFluxModule.hpp:179
void updateSolvent(const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Definition: NewTranFluxModule.hpp:211
void calculateBoundaryFluxes_(const ElementContext &, unsigned, unsigned)
Definition: NewTranFluxModule.hpp:567
const DimMatrix & intrinsicPermeability() const
Return the intrinsic permeability tensor at a face [m^2].
Definition: NewTranFluxModule.hpp:134
static void calculatePhasePressureDiff_(short &upIdx, short &dnIdx, EvalType &pressureDifference, const IntensiveQuantities &intQuantsIn, const IntensiveQuantities &intQuantsEx, const unsigned phaseIdx, const unsigned interiorDofIdx, const unsigned exteriorDofIdx, const Scalar Vin, const Scalar Vex, const unsigned globalIndexIn, const unsigned globalIndexEx, const Scalar distZg, const Scalar thpres)
Definition: NewTranFluxModule.hpp:307
void calculateBoundaryGradients_(const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx, const FluidState &exFluidState)
Update the required gradients for boundary faces.
Definition: NewTranFluxModule.hpp:421
static void calculateBoundaryGradients_(const Problem &problem, const unsigned globalSpaceIdx, const IntensiveQuantities &intQuantsIn, const unsigned bfIdx, const double faceArea, const double zEx, const FluidState &exFluidState, std::array< short, numPhases > &upIdx, std::array< short, numPhases > &dnIdx, EvaluationContainer &volumeFlux, EvaluationContainer &pressureDifference)
Update the required gradients for boundary faces.
Definition: NewTranFluxModule.hpp:465
unsigned upstreamIndex_(unsigned phaseIdx) const
Returns the local index of the degree of freedom in which is in upstream direction.
Definition: NewTranFluxModule.hpp:190
void updatePolymer(const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Definition: NewTranFluxModule.hpp:214
void calculateFluxes_(const ElementContext &, unsigned, unsigned)
Update the volumetric fluxes for all fluid phases on the interior faces of the context.
Definition: NewTranFluxModule.hpp:564
const EvalDimVector & potentialGrad(unsigned) const
Return the pressure potential gradient of a fluid phase at the face's integration point [Pa/m].
Definition: NewTranFluxModule.hpp:145
static void volumeAndPhasePressureDifferences(std::array< short, numPhases > &upIdx, std::array< short, numPhases > &dnIdx, Evaluation(&volumeFlux)[numPhases], Evaluation(&pressureDifferences)[numPhases], const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Definition: NewTranFluxModule.hpp:219
unsigned downstreamIndex_(unsigned phaseIdx) const
Returns the local index of the degree of freedom in which is in downstream direction.
Definition: NewTranFluxModule.hpp:204
const EvalDimVector & filterVelocity(unsigned) const
Return the filter velocity of a fluid phase at the face's integration point [m/s].
Definition: NewTranFluxModule.hpp:165
const Evaluation & pressureDifference(unsigned phaseIdx) const
Return the gravity corrected pressure difference between the interior and the exterior of a face.
Definition: NewTranFluxModule.hpp:156
void calculateGradients_(const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Update the required gradients for interior faces.
Definition: NewTranFluxModule.hpp:410
Provides the intensive quantities for the ECL flux module.
Definition: NewTranFluxModule.hpp:94
void update_(const ElementContext &, unsigned, unsigned)
Definition: NewTranFluxModule.hpp:97
Definition: BlackoilPhases.hpp:27
Specifies a flux module which uses ECL transmissibilities.
Definition: NewTranFluxModule.hpp:67
static void registerParameters()
Register all run-time parameters for the flux module.
Definition: NewTranFluxModule.hpp:75