31#ifndef EWOMS_TRANS_FLUX_MODULE_HH
32#define EWOMS_TRANS_FLUX_MODULE_HH
34#include <dune/common/fmatrix.hh>
35#include <dune/common/fvector.hh>
37#include <opm/material/common/MathToolbox.hpp>
38#include <opm/material/common/Valgrind.hpp>
53template <
class TypeTag>
54class TransIntensiveQuantities;
56template <
class TypeTag>
57class TransExtensiveQuantities;
59template <
class TypeTag>
60class TransBaseProblem;
65template <
class TypeTag>
83template <
class TypeTag>
90template <
class TypeTag>
96 void update_(
const ElementContext&,
unsigned,
unsigned)
103template <
class TypeTag>
116 enum { dimWorld = GridView::dimensionworld };
117 enum { numPhases = FluidSystem::numPhases };
119 typedef MathToolbox<Evaluation> Toolbox;
120 typedef Dune::FieldVector<Evaluation, dimWorld> EvalDimVector;
121 typedef Dune::FieldMatrix<Scalar, dimWorld, dimWorld> DimMatrix;
129 throw std::logic_error(
"The ECL transmissibility module does not "
130 "provide an explicit intrinsic permeability");
141 throw std::logic_error(
"The ECL transmissibility module does not "
142 "provide explicit potential gradients");
152 {
return pressureDifference_[phaseIdx]; }
162 throw std::logic_error(
"The ECL transmissibility module does not "
163 "provide explicit filter velocities");
176 {
return volumeFlux_[phaseIdx]; }
188 assert(phaseIdx < numPhases);
190 return upIdx_[phaseIdx];
202 assert(phaseIdx < numPhases);
204 return dnIdx_[phaseIdx];
207 void updateSolvent(
const ElementContext& elemCtx,
unsigned scvfIdx,
unsigned timeIdx)
208 { asImp_().updateVolumeFluxTrans(elemCtx, scvfIdx, timeIdx); }
210 void updatePolymer(
const ElementContext& elemCtx,
unsigned scvfIdx,
unsigned timeIdx)
211 { asImp_().updateShearMultipliers(elemCtx, scvfIdx, timeIdx); }
218 Valgrind::SetUndefined(*
this);
221 static_assert(std::is_same_v<Discretization, EcfvDiscretization<TypeTag>>);
223 const auto& stencil = elemCtx.stencil(timeIdx);
224 const auto& scvf = stencil.interiorFace(scvfIdx);
226 interiorDofIdx_ = scvf.interiorIndex();
227 exteriorDofIdx_ = scvf.exteriorIndex();
228 assert(interiorDofIdx_ != exteriorDofIdx_);
230 const unsigned I = stencil.globalSpaceIndex(interiorDofIdx_);
231 const unsigned J = stencil.globalSpaceIndex(exteriorDofIdx_);
233 const Scalar trans = transmissibility_(elemCtx, scvfIdx, timeIdx);
238 const Scalar g = elemCtx.problem().gravity()[dimWorld - 1];
240 const auto& intQuantsIn = elemCtx.intensiveQuantities(interiorDofIdx_, timeIdx);
241 const auto& intQuantsEx = elemCtx.intensiveQuantities(exteriorDofIdx_, timeIdx);
243 const Scalar zIn = dofCenterDepth_(elemCtx, interiorDofIdx_, timeIdx);
244 const Scalar zEx = dofCenterDepth_(elemCtx, exteriorDofIdx_, timeIdx);
247 const Scalar distZ = zIn - zEx;
249 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
250 if (!FluidSystem::phaseIsActive(phaseIdx)) {
256 if (intQuantsIn.mobility(phaseIdx) <= 0.0 &&
257 intQuantsEx.mobility(phaseIdx) <= 0.0)
259 upIdx_[phaseIdx] = interiorDofIdx_;
260 dnIdx_[phaseIdx] = exteriorDofIdx_;
261 pressureDifference_[phaseIdx] = 0.0;
262 volumeFlux_[phaseIdx] = 0.0;
268 const Evaluation& rhoIn = intQuantsIn.fluidState().density(phaseIdx);
269 const Scalar rhoEx = Toolbox::value(intQuantsEx.fluidState().density(phaseIdx));
270 const Evaluation rhoAvg = (rhoIn + rhoEx) / 2;
272 const Evaluation& pressureInterior = intQuantsIn.fluidState().pressure(phaseIdx);
273 Evaluation pressureExterior = Toolbox::value(intQuantsEx.fluidState().pressure(phaseIdx));
275 pressureExterior += rhoAvg * (distZ * g);
277 pressureDifference_[phaseIdx] = pressureExterior - pressureInterior;
283 if (pressureDifference_[phaseIdx] > 0.0) {
284 upIdx_[phaseIdx] = exteriorDofIdx_;
285 dnIdx_[phaseIdx] = interiorDofIdx_;
287 else if (pressureDifference_[phaseIdx] < 0.0) {
288 upIdx_[phaseIdx] = interiorDofIdx_;
289 dnIdx_[phaseIdx] = exteriorDofIdx_;
294 const Scalar Vin = elemCtx.dofVolume(interiorDofIdx_, 0);
295 const Scalar Vex = elemCtx.dofVolume(exteriorDofIdx_, 0);
297 upIdx_[phaseIdx] = interiorDofIdx_;
298 dnIdx_[phaseIdx] = exteriorDofIdx_;
300 else if (Vin < Vex) {
301 upIdx_[phaseIdx] = exteriorDofIdx_;
302 dnIdx_[phaseIdx] = interiorDofIdx_;
309 upIdx_[phaseIdx] = interiorDofIdx_;
310 dnIdx_[phaseIdx] = exteriorDofIdx_;
313 upIdx_[phaseIdx] = exteriorDofIdx_;
314 dnIdx_[phaseIdx] = interiorDofIdx_;
323 const auto& up = elemCtx.intensiveQuantities(upstreamIdx, timeIdx);
325 if (upstreamIdx == interiorDofIdx_) {
326 volumeFlux_[phaseIdx] =
327 pressureDifference_[phaseIdx] * up.mobility(phaseIdx) * (-trans);
330 volumeFlux_[phaseIdx] =
331 pressureDifference_[phaseIdx] * (Toolbox::value(up.mobility(phaseIdx)) * (-trans));
339 template <
class Flu
idState>
343 const FluidState& exFluidState)
345 const auto& stencil = elemCtx.stencil(timeIdx);
346 const auto& scvf = stencil.boundaryFace(scvfIdx);
348 interiorDofIdx_ = scvf.interiorIndex();
350 const Scalar trans = transmissibilityBoundary_(elemCtx, scvfIdx, timeIdx);
355 const Scalar g = elemCtx.problem().gravity()[dimWorld - 1];
357 const auto& intQuantsIn = elemCtx.intensiveQuantities(interiorDofIdx_, timeIdx);
364 const Scalar zIn = dofCenterDepth_(elemCtx, interiorDofIdx_, timeIdx);
365 const Scalar zEx = scvf.integrationPos()[dimWorld - 1];
369 const Scalar distZ = zIn - zEx;
371 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
372 if (!FluidSystem::phaseIsActive(phaseIdx)) {
378 const Evaluation& rhoIn = intQuantsIn.fluidState().density(phaseIdx);
379 const auto& rhoEx = exFluidState.density(phaseIdx);
380 const Evaluation rhoAvg = (rhoIn + rhoEx) / 2;
382 const Evaluation& pressureInterior = intQuantsIn.fluidState().pressure(phaseIdx);
383 Evaluation pressureExterior = exFluidState.pressure(phaseIdx);
384 pressureExterior += rhoAvg * (distZ * g);
386 pressureDifference_[phaseIdx] = pressureExterior - pressureInterior;
392 if (pressureDifference_[phaseIdx] > 0.0) {
393 upIdx_[phaseIdx] = -1;
394 dnIdx_[phaseIdx] = interiorDofIdx_;
397 upIdx_[phaseIdx] = interiorDofIdx_;
398 dnIdx_[phaseIdx] = -1;
402 if (upstreamIdx == interiorDofIdx_) {
406 const auto& up = elemCtx.intensiveQuantities(upstreamIdx, timeIdx);
408 volumeFlux_[phaseIdx] =
409 pressureDifference_[phaseIdx] * up.mobility(phaseIdx) * (-trans);
414 const auto& matParams =
415 elemCtx.problem().materialLawParams(elemCtx,
418 std::array<typename FluidState::Scalar,numPhases> kr;
419 MaterialLaw::relativePermeabilities(kr, matParams, exFluidState);
421 const auto& mob = kr[phaseIdx] / exFluidState.viscosity(phaseIdx);
422 volumeFlux_[phaseIdx] =
423 pressureDifference_[phaseIdx] * mob * (-trans);
438 Scalar transmissibility_(
const ElementContext& elemCtx,
unsigned scvfIdx,
unsigned timeIdx)
const
440 const auto& stencil = elemCtx.stencil(timeIdx);
441 const auto& face = stencil.interiorFace(scvfIdx);
442 const auto& interiorPos = stencil.subControlVolume(face.interiorIndex()).globalPos();
443 const auto& exteriorPos = stencil.subControlVolume(face.exteriorIndex()).globalPos();
444 const auto distVec0 = face.integrationPos() - interiorPos;
445 const auto distVec1 = face.integrationPos() - exteriorPos;
446 const Scalar ndotDistIn = std::abs(face.normal() * distVec0);
447 const Scalar ndotDistExt = std::abs(face.normal() * distVec1);
449 const Scalar distSquaredIn = distVec0 * distVec0;
450 const Scalar distSquaredExt = distVec1 * distVec1;
451 const auto& K0mat = elemCtx.problem().intrinsicPermeability(elemCtx, face.interiorIndex(), timeIdx);
452 const auto& K1mat = elemCtx.problem().intrinsicPermeability(elemCtx, face.exteriorIndex(), timeIdx);
459 for (
unsigned i = 0; i < dimWorld; ++i) {
460 if (std::abs(face.normal()[i]) > val) {
461 val = std::abs(face.normal()[i]);
465 const Scalar& K0 = K0mat[idx][idx];
466 const Scalar& K1 = K1mat[idx][idx];
467 const Scalar T0 = K0 * ndotDistIn / distSquaredIn;
468 const Scalar T1 = K1 * ndotDistExt / distSquaredExt;
469 return T0 * T1 / (T0 + T1);
471 Scalar transmissibilityBoundary_(
const ElementContext& elemCtx,
unsigned scvfIdx,
unsigned timeIdx)
const
473 const auto& stencil = elemCtx.stencil(timeIdx);
474 const auto& face = stencil.interiorFace(scvfIdx);
475 const auto& interiorPos = stencil.subControlVolume(face.interiorIndex()).globalPos();
476 const auto distVec0 = face.integrationPos() - interiorPos;
477 const Scalar ndotDistIn = face.normal() * distVec0;
478 const Scalar distSquaredIn = distVec0 * distVec0;
479 const auto& K0mat = elemCtx.problem().intrinsicPermeability(elemCtx, face.interiorIndex(), timeIdx);
486 for (
unsigned i = 0; i < dimWorld; ++i) {
487 if (std::abs(face.normal()[i]) > val) {
488 val = std::abs(face.normal()[i]);
492 const Scalar& K0 = K0mat[idx][idx];
493 const Scalar T0 = K0 * ndotDistIn / distSquaredIn;
497 template <
class Context>
498 Scalar dofCenterDepth_(
const Context& context,
unsigned spaceIdx,
unsigned timeIdx)
const
500 const auto& pos = context.pos(spaceIdx, timeIdx);
501 return pos[dimWorld-1];
504 Implementation& asImp_()
505 {
return *
static_cast<Implementation*
>(
this); }
507 const Implementation& asImp_()
const
508 {
return *
static_cast<const Implementation*
>(
this); }
511 std::array<Evaluation, numPhases> volumeFlux_;
515 std::array<Evaluation, numPhases> pressureDifference_;
518 unsigned short interiorDofIdx_{};
519 unsigned short exteriorDofIdx_{};
520 std::array<short, numPhases> upIdx_{};
521 std::array<short, numPhases> dnIdx_{};
Provides the defaults for the parameters required by the transmissibility based volume flux calculati...
Definition: transfluxmodule.hh:85
Provides the transmissibility based flux module.
Definition: transfluxmodule.hh:105
const Evaluation & volumeFlux(unsigned phaseIdx) const
Return the volume flux of a fluid phase at the face's integration point .
Definition: transfluxmodule.hh:175
unsigned downstreamIndex_(unsigned phaseIdx) const
Returns the local index of the degree of freedom in which is in downstream direction.
Definition: transfluxmodule.hh:200
const EvalDimVector & filterVelocity(unsigned) const
Return the filter velocity of a fluid phase at the face's integration point [m/s].
Definition: transfluxmodule.hh:160
void calculateGradients_(const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Update the required gradients for interior faces.
Definition: transfluxmodule.hh:216
void updateSolvent(const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Definition: transfluxmodule.hh:207
const Evaluation & pressureDifference(unsigned phaseIdx) const
Return the gravity corrected pressure difference between the interior and the exterior of a face.
Definition: transfluxmodule.hh:151
void calculateFluxes_(const ElementContext &, unsigned, unsigned)
Update the volumetric fluxes for all fluid phases on the interior faces of the context.
Definition: transfluxmodule.hh:431
void calculateBoundaryFluxes_(const ElementContext &, unsigned, unsigned)
Definition: transfluxmodule.hh:434
const EvalDimVector & potentialGrad(unsigned) const
Return the pressure potential gradient of a fluid phase at the face's integration point [Pa/m].
Definition: transfluxmodule.hh:139
void calculateBoundaryGradients_(const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx, const FluidState &exFluidState)
Update the required gradients for boundary faces.
Definition: transfluxmodule.hh:340
const DimMatrix & intrinsicPermeability() const
Return the intrinsic permeability tensor at a face [m^2].
Definition: transfluxmodule.hh:127
unsigned upstreamIndex_(unsigned phaseIdx) const
Returns the local index of the degree of freedom in which is in upstream direction.
Definition: transfluxmodule.hh:186
void updatePolymer(const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Definition: transfluxmodule.hh:210
Provides the intensive quantities for the transmissibility based flux module.
Definition: transfluxmodule.hh:92
void update_(const ElementContext &, unsigned, unsigned)
Definition: transfluxmodule.hh:96
Declare the properties used by the infrastructure code of the finite volume discretizations.
Defines the common properties required by the porous medium multi-phase models.
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
Specifies a flux module which uses transmissibilities.
Definition: transfluxmodule.hh:67
static void registerParameters()
Register all run-time parameters for the flux module.
Definition: transfluxmodule.hh:75