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