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