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 };
120 typedef Dune::FieldVector<Scalar, dimWorld> DimVector;
121 typedef Dune::FieldVector<Evaluation, dimWorld> EvalDimVector;
122 typedef Dune::FieldMatrix<Scalar, dimWorld, dimWorld> DimMatrix;
130 throw std::logic_error(
"The ECL transmissibility module does not "
131 "provide an explicit intrinsic permeability");
142 throw std::logic_error(
"The ECL transmissibility module does not "
143 "provide explicit potential gradients");
153 {
return pressureDifference_[
phaseIdx]; }
163 throw std::logic_error(
"The ECL transmissibility module does not "
164 "provide explicit filter velocities");
208 void updateSolvent(
const ElementContext& elemCtx,
unsigned scvfIdx,
unsigned timeIdx)
211 void updatePolymer(
const ElementContext& elemCtx,
unsigned scvfIdx,
unsigned timeIdx)
219 Valgrind::SetUndefined(*
this);
222 static_assert(std::is_same_v<Discretization, EcfvDiscretization<TypeTag>>);
224 const auto& stencil = elemCtx.stencil(
timeIdx);
227 interiorDofIdx_ =
scvf.interiorIndex();
228 exteriorDofIdx_ =
scvf.exteriorIndex();
229 assert(interiorDofIdx_ != exteriorDofIdx_);
231 const unsigned I = stencil.globalSpaceIndex(interiorDofIdx_);
232 const unsigned J = stencil.globalSpaceIndex(exteriorDofIdx_);
234 const Scalar trans = transmissibility_(elemCtx,
scvfIdx,
timeIdx);
239 const Scalar
g = elemCtx.problem().gravity()[dimWorld - 1];
244 const Scalar
zIn = dofCenterDepth_(elemCtx, interiorDofIdx_,
timeIdx);
245 const Scalar
zEx = dofCenterDepth_(elemCtx, exteriorDofIdx_,
timeIdx);
251 if (!FluidSystem::phaseIsActive(
phaseIdx)) {
262 pressureDifference_[
phaseIdx] = 0.0;
284 if (pressureDifference_[
phaseIdx] > 0.0) {
288 else if (pressureDifference_[
phaseIdx] < 0.0) {
295 const Scalar Vin = elemCtx.dofVolume(interiorDofIdx_, 0);
296 const Scalar Vex = elemCtx.dofVolume(exteriorDofIdx_, 0);
301 else if (Vin < Vex) {
340 template <
class Flu
idState>
344 const FluidState& exFluidState)
346 const auto& stencil = elemCtx.stencil(
timeIdx);
349 interiorDofIdx_ =
scvf.interiorIndex();
351 const Scalar trans = transmissibilityBoundary_(elemCtx,
scvfIdx,
timeIdx);
356 const Scalar
g = elemCtx.problem().gravity()[dimWorld - 1];
365 const Scalar
zIn = dofCenterDepth_(elemCtx, interiorDofIdx_,
timeIdx);
366 const Scalar
zEx =
scvf.integrationPos()[dimWorld - 1];
373 if (!FluidSystem::phaseIsActive(
phaseIdx)) {
393 if (pressureDifference_[
phaseIdx] > 0.0) {
416 elemCtx.problem().materialLawParams(elemCtx,
419 std::array<typename FluidState::Scalar,numPhases>
kr;
420 MaterialLaw::relativePermeabilities(
kr,
matParams, exFluidState);
435 void calculateBoundaryFluxes_(
const ElementContext&,
unsigned,
unsigned)
439 Scalar transmissibility_(
const ElementContext& elemCtx,
unsigned scvfIdx,
unsigned timeIdx)
const
441 const auto& stencil = elemCtx.stencil(
timeIdx);
443 const auto&
interiorPos = stencil.subControlVolume(
face.interiorIndex()).globalPos();
444 const auto&
exteriorPos = stencil.subControlVolume(
face.exteriorIndex()).globalPos();
452 const auto&
K0mat = elemCtx.problem().intrinsicPermeability(elemCtx,
face.interiorIndex(),
timeIdx);
453 const auto&
K1mat = elemCtx.problem().intrinsicPermeability(elemCtx,
face.exteriorIndex(),
timeIdx);
460 for (
unsigned i = 0; i < dimWorld; ++i) {
461 if (std::abs(
face.normal()[i]) >
val) {
462 val = std::abs(
face.normal()[i]);
472 Scalar transmissibilityBoundary_(
const ElementContext& elemCtx,
unsigned scvfIdx,
unsigned timeIdx)
const
474 const auto& stencil = elemCtx.stencil(
timeIdx);
476 const auto&
interiorPos = stencil.subControlVolume(
face.interiorIndex()).globalPos();
480 const auto&
K0mat = elemCtx.problem().intrinsicPermeability(elemCtx,
face.interiorIndex(),
timeIdx);
487 for (
unsigned i = 0; i < dimWorld; ++i) {
488 if (std::abs(
face.normal()[i]) >
val) {
489 val = std::abs(
face.normal()[i]);
498 template <
class Context>
499 Scalar dofCenterDepth_(
const Context& context,
unsigned spaceIdx,
unsigned timeIdx)
const
502 return pos[dimWorld-1];
505 Implementation& asImp_()
506 {
return *
static_cast<Implementation*
>(
this); }
508 const Implementation& asImp_()
const
509 {
return *
static_cast<const Implementation*
>(
this); }
512 std::array<Evaluation, numPhases> volumeFlux_;
516 std::array<Evaluation, numPhases> pressureDifference_;
519 unsigned short interiorDofIdx_{};
520 unsigned short exteriorDofIdx_{};
521 std::array<short, numPhases> upIdx_{};
522 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:176
unsigned downstreamIndex_(unsigned phaseIdx) const
Returns the local index of the degree of freedom in which is in downstream direction.
Definition transfluxmodule.hh:201
const EvalDimVector & filterVelocity(unsigned) const
Return the filter velocity of a fluid phase at the face's integration point [m/s].
Definition transfluxmodule.hh:161
void calculateGradients_(const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Update the required gradients for interior faces.
Definition transfluxmodule.hh:217
const Evaluation & pressureDifference(unsigned phaseIdx) const
Return the gravity corrected pressure difference between the interior and the exterior of a face.
Definition transfluxmodule.hh:152
void calculateFluxes_(const ElementContext &, unsigned, unsigned)
Update the volumetric fluxes for all fluid phases on the interior faces of the context.
Definition transfluxmodule.hh:432
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:140
void calculateBoundaryGradients_(const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx, const FluidState &exFluidState)
Update the required gradients for boundary faces.
Definition transfluxmodule.hh:341
const DimMatrix & intrinsicPermeability() const
Return the intrinsic permeability tensor at a face [m^2].
Definition transfluxmodule.hh:128
unsigned upstreamIndex_(unsigned phaseIdx) const
Returns the local index of the degree of freedom in which is in upstream direction.
Definition transfluxmodule.hh:187
Provides the intensive quantities for the transmissibility based flux module.
Definition transfluxmodule.hh:92
The base class for the element-centered finite-volume discretization scheme.
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.
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:43
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:240
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