31#ifndef OPM_NEWTRAN_FLUX_MODULE_HPP
32#define OPM_NEWTRAN_FLUX_MODULE_HPP
34#include <dune/common/fvector.hh>
35#include <dune/common/fmatrix.hh>
37#include <opm/common/ErrorMacros.hpp>
38#include <opm/common/OpmLog/OpmLog.hpp>
40#include <opm/input/eclipse/EclipseState/Grid/FaceDir.hpp>
42#include <opm/material/common/MathToolbox.hpp>
43#include <opm/material/common/Valgrind.hpp>
44#include <opm/material/thermal/EnergyModuleType.hpp>
53#include <opm/common/utility/gpuDecorators.hpp>
59template <
class TypeTag>
60class NewTranIntensiveQuantities;
62template <
class TypeTag>
63class NewTranExtensiveQuantities;
65template <
class TypeTag>
66class NewTranBaseProblem;
72template <
class TypeTag>
91template <
class TypeTag>
99template <
class TypeTag>
104 void update_(
const ElementContext&,
unsigned,
unsigned)
112template <
class TypeTag>
125 enum { dimWorld = GridView::dimensionworld };
126 enum { gasPhaseIdx = FluidSystem::gasPhaseIdx };
127 enum { numPhases = FluidSystem::numPhases };
128 static constexpr bool enableSolvent = getPropValue<TypeTag, Properties::EnableSolvent>();
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>();
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>;
141 using ConvectiveMixingModule = BlackOilConvectiveMixingModule<TypeTag, enableConvectiveMixing>;
150 OPM_THROW(std::invalid_argument,
"The ECL transmissibility module does not provide an explicit intrinsic permeability");
161 OPM_THROW(std::invalid_argument,
"The ECL transmissibility module does not provide explicit potential gradients");
171 {
return pressureDifference_[phaseIdx]; }
181 OPM_THROW(std::invalid_argument,
"The ECL transmissibility module does not provide explicit filter velocities");
193 OPM_HOST_DEVICE
const Evaluation&
volumeFlux(
unsigned phaseIdx)
const
194 {
return volumeFlux_[phaseIdx]; }
206 assert(phaseIdx < numPhases);
208 return upIdx_[phaseIdx];
220 assert(phaseIdx < numPhases);
222 return dnIdx_[phaseIdx];
225 OPM_HOST_DEVICE
void updateSolvent(
const ElementContext& elemCtx,
unsigned scvfIdx,
unsigned timeIdx)
227 if constexpr (enableSolvent) {
228 asImp_().updateVolumeFluxTrans(elemCtx, scvfIdx, timeIdx);
232 OPM_HOST_DEVICE
void updatePolymer(
const ElementContext& elemCtx,
unsigned scvfIdx,
unsigned timeIdx)
234 if constexpr (enablePolymer) {
235 asImp_().updateShearMultipliers(elemCtx, scvfIdx, timeIdx);
242 std::array<short, numPhases>& dnIdx,
244 Evaluation (&pressureDifferences)[numPhases],
245 const ElementContext& elemCtx,
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();
255 assert(interiorDofIdx != exteriorDofIdx);
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);
267 Scalar g = elemCtx.problem().gravity()[dimWorld - 1];
269 const auto& intQuantsIn = elemCtx.intensiveQuantities(interiorDofIdx, timeIdx);
270 const auto& intQuantsEx = elemCtx.intensiveQuantities(exteriorDofIdx, timeIdx);
277 Scalar zIn = problem.dofCenterDepth(elemCtx, interiorDofIdx, timeIdx);
278 Scalar zEx = problem.dofCenterDepth(elemCtx, exteriorDofIdx, timeIdx);
282 Scalar distZ = zIn - zEx;
284 Scalar Vin = elemCtx.dofVolume(interiorDofIdx, 0);
285 Scalar Vex = elemCtx.dofVolume(exteriorDofIdx, 0);
287 for (
unsigned phaseIdx=0; phaseIdx < numPhases; phaseIdx++) {
288 if (!FluidSystem::phaseIsActive(phaseIdx))
292 pressureDifferences[phaseIdx],
305 problem.moduleParams());
307 const bool upwindIsInterior = (
static_cast<unsigned>(upIdx[phaseIdx]) == interiorDofIdx);
308 const IntensiveQuantities& up = upwindIsInterior ? intQuantsIn : intQuantsEx;
310 const Evaluation transMult = (intQuantsIn.rockCompTransMultiplier() + Toolbox::value(intQuantsEx.rockCompTransMultiplier()))/2;
312 const auto& materialLawManager = problem.materialLawManager();
313 FaceDir::DirEnum facedir = FaceDir::DirEnum::Unknown;
314 if (materialLawManager->hasDirectionalRelperms()) {
315 facedir = scvf.faceDirFromDirId();
317 if (upwindIsInterior)
319 pressureDifferences[phaseIdx]*up.mobility(phaseIdx, facedir)*transMult*(-trans/faceArea);
322 pressureDifferences[phaseIdx]*
323 (Toolbox::value(up.mobility(phaseIdx, facedir))*transMult*(-trans/faceArea));
327 template<
class EvalType>
331 const IntensiveQuantities& intQuantsIn,
332 const IntensiveQuantities& intQuantsEx,
333 const unsigned phaseIdx,
334 const unsigned interiorDofIdx,
335 const unsigned exteriorDofIdx,
338 const unsigned globalIndexIn,
339 const unsigned globalIndexEx,
341 const Scalar thpresInToEx,
342 const Scalar thpresExToIn,
348 if (intQuantsIn.mobility(phaseIdx) <= 0.0 &&
349 intQuantsEx.mobility(phaseIdx) <= 0.0)
351 upIdx = interiorDofIdx;
352 dnIdx = exteriorDofIdx;
359 const Evaluation& rhoIn = intQuantsIn.fluidState().density(phaseIdx);
360 Scalar rhoEx = Toolbox::value(intQuantsEx.fluidState().density(phaseIdx));
361 Evaluation rhoAvg = (rhoIn + rhoEx)/2;
363 if constexpr (enableConvectiveMixing) {
367 const Evaluation& pressureInterior = intQuantsIn.fluidState().pressure(phaseIdx);
368 Evaluation pressureExterior = Toolbox::value(intQuantsEx.fluidState().pressure(phaseIdx));
371 pressureExterior += Toolbox::value(rhoAvg)*(distZg);
373 pressureExterior += rhoAvg*(distZg);
382 upIdx = exteriorDofIdx;
383 dnIdx = interiorDofIdx;
386 upIdx = interiorDofIdx;
387 dnIdx = exteriorDofIdx;
393 upIdx = interiorDofIdx;
394 dnIdx = exteriorDofIdx;
396 else if (Vin < Vex) {
397 upIdx = exteriorDofIdx;
398 dnIdx = interiorDofIdx;
404 if (globalIndexIn < globalIndexEx) {
405 upIdx = interiorDofIdx;
406 dnIdx = exteriorDofIdx;
409 upIdx = exteriorDofIdx;
410 dnIdx = interiorDofIdx;
437 OPM_HOST_DEVICE
void calculateGradients_(
const ElementContext& elemCtx,
unsigned scvfIdx,
unsigned timeIdx)
439 Valgrind::SetUndefined(*
this);
447 template <
class Flu
idState>
451 const FluidState& exFluidState)
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);
470 pressureDifference_);
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);
482 asImp_().setSolventVolumeFlux(0.0);
491 template <
class Problem,
class Flu
idState,
class EvaluationContainer>
493 const unsigned globalSpaceIdx,
494 const IntensiveQuantities& intQuantsIn,
495 const unsigned bfIdx,
496 const double faceArea,
498 const FluidState& exFluidState,
499 std::array<short, numPhases>& upIdx,
500 std::array<short, numPhases>& dnIdx,
505 bool enableBoundaryMassFlux = problem.nonTrivialBoundaryConditions();
506 if (!enableBoundaryMassFlux)
509 Scalar trans = problem.transmissibilityBoundary(globalSpaceIdx, bfIdx);
514 Scalar g = problem.gravity()[dimWorld - 1];
521 Scalar zIn = problem.dofCenterDepth(globalSpaceIdx);
525 Scalar distZ = zIn - zEx;
527 for (
unsigned phaseIdx=0; phaseIdx < numPhases; phaseIdx++) {
528 if (!FluidSystem::phaseIsActive(phaseIdx))
533 const Evaluation& rhoIn = intQuantsIn.fluidState().density(phaseIdx);
534 const auto& rhoEx = exFluidState.density(phaseIdx);
535 Evaluation rhoAvg = (rhoIn + rhoEx)/2;
537 const Evaluation& pressureInterior = intQuantsIn.fluidState().pressure(phaseIdx);
538 Evaluation pressureExterior = exFluidState.pressure(phaseIdx);
539 pressureExterior += rhoAvg*(distZ*g);
547 const unsigned interiorDofIdx = 0;
549 upIdx[phaseIdx] = -1;
550 dnIdx[phaseIdx] = interiorDofIdx;
553 upIdx[phaseIdx] = interiorDofIdx;
554 dnIdx[phaseIdx] = -1;
557 Evaluation transModified = trans;
559 if (upIdx[phaseIdx] == interiorDofIdx) {
563 const auto& up = intQuantsIn;
566 const Scalar transMult = Toolbox::value(up.rockCompTransMultiplier());
567 transModified *= transMult;
575 const auto& matParams = problem.materialLawParams(globalSpaceIdx);
576 std::array<typename FluidState::ValueType,numPhases> kr;
577 MaterialLaw::relativePermeabilities(kr, matParams, exFluidState);
579 const auto& mob = kr[phaseIdx]/exFluidState.viscosity(phaseIdx);
598 Implementation& asImp_()
599 {
return *
static_cast<Implementation*
>(
this); }
601 const Implementation& asImp_()
const
602 {
return *
static_cast<const Implementation*
>(
this); }
605 Evaluation volumeFlux_[numPhases];
609 Evaluation pressureDifference_[numPhases];
612 std::array<short, numPhases> upIdx_;
613 std::array<short, numPhases> dnIdx_;
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