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