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<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");
177 {
return volumeFlux_[phaseIdx]; }
189 assert(phaseIdx < numPhases);
191 return upIdx_[phaseIdx];
203 assert(phaseIdx < numPhases);
205 return dnIdx_[phaseIdx];
208 void updateSolvent(
const ElementContext& elemCtx,
unsigned scvfIdx,
unsigned timeIdx)
209 { asImp_().updateVolumeFluxTrans(elemCtx, scvfIdx, timeIdx); }
211 void updatePolymer(
const ElementContext& elemCtx,
unsigned scvfIdx,
unsigned timeIdx)
212 { asImp_().updateShearMultipliers(elemCtx, scvfIdx, timeIdx); }
219 Valgrind::SetUndefined(*
this);
222 static_assert(std::is_same_v<Discretization, EcfvDiscretization<TypeTag>>);
224 const auto& stencil = elemCtx.stencil(timeIdx);
225 const auto& scvf = stencil.interiorFace(scvfIdx);
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];
241 const auto& intQuantsIn = elemCtx.intensiveQuantities(interiorDofIdx_, timeIdx);
242 const auto& intQuantsEx = elemCtx.intensiveQuantities(exteriorDofIdx_, timeIdx);
244 const Scalar zIn = dofCenterDepth_(elemCtx, interiorDofIdx_, timeIdx);
245 const Scalar zEx = dofCenterDepth_(elemCtx, exteriorDofIdx_, timeIdx);
248 const Scalar distZ = zIn - zEx;
250 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
251 if (!FluidSystem::phaseIsActive(phaseIdx)) {
257 if (intQuantsIn.mobility(phaseIdx) <= 0.0 &&
258 intQuantsEx.mobility(phaseIdx) <= 0.0)
260 upIdx_[phaseIdx] = interiorDofIdx_;
261 dnIdx_[phaseIdx] = exteriorDofIdx_;
262 pressureDifference_[phaseIdx] = 0.0;
263 volumeFlux_[phaseIdx] = 0.0;
269 const Evaluation& rhoIn = intQuantsIn.fluidState().density(phaseIdx);
270 const Scalar rhoEx = Toolbox::value(intQuantsEx.fluidState().density(phaseIdx));
271 const Evaluation rhoAvg = (rhoIn + rhoEx) / 2;
273 const Evaluation& pressureInterior = intQuantsIn.fluidState().pressure(phaseIdx);
274 Evaluation pressureExterior = Toolbox::value(intQuantsEx.fluidState().pressure(phaseIdx));
276 pressureExterior += rhoAvg * (distZ * g);
278 pressureDifference_[phaseIdx] = pressureExterior - pressureInterior;
284 if (pressureDifference_[phaseIdx] > 0.0) {
285 upIdx_[phaseIdx] = exteriorDofIdx_;
286 dnIdx_[phaseIdx] = interiorDofIdx_;
288 else if (pressureDifference_[phaseIdx] < 0.0) {
289 upIdx_[phaseIdx] = interiorDofIdx_;
290 dnIdx_[phaseIdx] = exteriorDofIdx_;
295 const Scalar Vin = elemCtx.dofVolume(interiorDofIdx_, 0);
296 const Scalar Vex = elemCtx.dofVolume(exteriorDofIdx_, 0);
298 upIdx_[phaseIdx] = interiorDofIdx_;
299 dnIdx_[phaseIdx] = exteriorDofIdx_;
301 else if (Vin < Vex) {
302 upIdx_[phaseIdx] = exteriorDofIdx_;
303 dnIdx_[phaseIdx] = interiorDofIdx_;
310 upIdx_[phaseIdx] = interiorDofIdx_;
311 dnIdx_[phaseIdx] = exteriorDofIdx_;
314 upIdx_[phaseIdx] = exteriorDofIdx_;
315 dnIdx_[phaseIdx] = interiorDofIdx_;
324 const auto& up = elemCtx.intensiveQuantities(upstreamIdx, timeIdx);
326 if (upstreamIdx == interiorDofIdx_) {
327 volumeFlux_[phaseIdx] =
328 pressureDifference_[phaseIdx] * up.mobility(phaseIdx) * (-trans);
331 volumeFlux_[phaseIdx] =
332 pressureDifference_[phaseIdx] * (Toolbox::value(up.mobility(phaseIdx)) * (-trans));
340 template <
class Flu
idState>
344 const FluidState& exFluidState)
346 const auto& stencil = elemCtx.stencil(timeIdx);
347 const auto& scvf = stencil.boundaryFace(scvfIdx);
349 interiorDofIdx_ = scvf.interiorIndex();
351 const Scalar trans = transmissibilityBoundary_(elemCtx, scvfIdx, timeIdx);
356 const Scalar g = elemCtx.problem().gravity()[dimWorld - 1];
358 const auto& intQuantsIn = elemCtx.intensiveQuantities(interiorDofIdx_, timeIdx);
365 const Scalar zIn = dofCenterDepth_(elemCtx, interiorDofIdx_, timeIdx);
366 const Scalar zEx = scvf.integrationPos()[dimWorld - 1];
370 const Scalar distZ = zIn - zEx;
372 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
373 if (!FluidSystem::phaseIsActive(phaseIdx)) {
379 const Evaluation& rhoIn = intQuantsIn.fluidState().density(phaseIdx);
380 const auto& rhoEx = exFluidState.density(phaseIdx);
381 const Evaluation rhoAvg = (rhoIn + rhoEx) / 2;
383 const Evaluation& pressureInterior = intQuantsIn.fluidState().pressure(phaseIdx);
384 Evaluation pressureExterior = exFluidState.pressure(phaseIdx);
385 pressureExterior += rhoAvg * (distZ * g);
387 pressureDifference_[phaseIdx] = pressureExterior - pressureInterior;
393 if (pressureDifference_[phaseIdx] > 0.0) {
394 upIdx_[phaseIdx] = -1;
395 dnIdx_[phaseIdx] = interiorDofIdx_;
398 upIdx_[phaseIdx] = interiorDofIdx_;
399 dnIdx_[phaseIdx] = -1;
403 if (upstreamIdx == interiorDofIdx_) {
407 const auto& up = elemCtx.intensiveQuantities(upstreamIdx, timeIdx);
409 volumeFlux_[phaseIdx] =
410 pressureDifference_[phaseIdx] * up.mobility(phaseIdx) * (-trans);
415 const auto& matParams =
416 elemCtx.problem().materialLawParams(elemCtx,
419 std::array<typename FluidState::Scalar,numPhases> kr;
420 MaterialLaw::relativePermeabilities(kr, matParams, exFluidState);
422 const auto& mob = kr[phaseIdx] / exFluidState.viscosity(phaseIdx);
423 volumeFlux_[phaseIdx] =
424 pressureDifference_[phaseIdx] * mob * (-trans);
439 Scalar transmissibility_(
const ElementContext& elemCtx,
unsigned scvfIdx,
unsigned timeIdx)
const
441 const auto& stencil = elemCtx.stencil(timeIdx);
442 const auto& face = stencil.interiorFace(scvfIdx);
443 const auto& interiorPos = stencil.subControlVolume(face.interiorIndex()).globalPos();
444 const auto& exteriorPos = stencil.subControlVolume(face.exteriorIndex()).globalPos();
445 const auto distVec0 = face.integrationPos() - interiorPos;
446 const auto distVec1 = face.integrationPos() - exteriorPos;
447 const Scalar ndotDistIn = std::abs(face.normal() * distVec0);
448 const Scalar ndotDistExt = std::abs(face.normal() * distVec1);
450 const Scalar distSquaredIn = distVec0 * distVec0;
451 const Scalar distSquaredExt = distVec1 * distVec1;
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]);
466 const Scalar& K0 = K0mat[idx][idx];
467 const Scalar& K1 = K1mat[idx][idx];
468 const Scalar T0 = K0 * ndotDistIn / distSquaredIn;
469 const Scalar T1 = K1 * ndotDistExt / distSquaredExt;
470 return T0 * T1 / (T0 + T1);
472 Scalar transmissibilityBoundary_(
const ElementContext& elemCtx,
unsigned scvfIdx,
unsigned timeIdx)
const
474 const auto& stencil = elemCtx.stencil(timeIdx);
475 const auto& face = stencil.interiorFace(scvfIdx);
476 const auto& interiorPos = stencil.subControlVolume(face.interiorIndex()).globalPos();
477 const auto distVec0 = face.integrationPos() - interiorPos;
478 const Scalar ndotDistIn = face.normal() * distVec0;
479 const Scalar distSquaredIn = distVec0 * distVec0;
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]);
493 const Scalar& K0 = K0mat[idx][idx];
494 const Scalar T0 = K0 * ndotDistIn / distSquaredIn;
498 template <
class Context>
499 Scalar dofCenterDepth_(
const Context& context,
unsigned spaceIdx,
unsigned timeIdx)
const
501 const auto& pos = context.pos(spaceIdx, timeIdx);
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
void updateSolvent(const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Definition: transfluxmodule.hh:208
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
void calculateBoundaryFluxes_(const ElementContext &, unsigned, unsigned)
Definition: transfluxmodule.hh:435
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
void updatePolymer(const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Definition: transfluxmodule.hh:211
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: blackoilboundaryratevector.hh:39
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