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