31#ifndef EWOMS_TRANS_FLUX_MODULE_HH
32#define EWOMS_TRANS_FLUX_MODULE_HH
36#include <opm/material/common/Valgrind.hpp>
37#include <dune/common/fvector.hh>
38#include <dune/common/fmatrix.hh>
42template <
class TypeTag>
43class TransIntensiveQuantities;
45template <
class TypeTag>
46class TransExtensiveQuantities;
48template <
class TypeTag>
49class TransBaseProblem;
54template <
class TypeTag>
71template <
class TypeTag>
78template <
class TypeTag>
83 void update_(
const ElementContext&,
unsigned,
unsigned)
90template <
class TypeTag>
103 enum { dimWorld = GridView::dimensionworld };
104 enum { numPhases = FluidSystem::numPhases };
106 typedef MathToolbox<Evaluation> Toolbox;
107 typedef Dune::FieldVector<Scalar, dimWorld> DimVector;
108 typedef Dune::FieldVector<Evaluation, dimWorld> EvalDimVector;
109 typedef Dune::FieldMatrix<Scalar, dimWorld, dimWorld> DimMatrix;
117 throw std::logic_error(
"The ECL transmissibility module does not provide an explicit intrinsic permeability");
128 throw std::logic_error(
"The ECL transmissibility module does not provide explicit potential gradients");
138 {
return pressureDifference_[phaseIdx]; }
148 throw std::logic_error(
"The ECL transmissibility module does not provide explicit filter velocities");
161 {
return volumeFlux_[phaseIdx]; }
173 assert(phaseIdx < numPhases);
175 return upIdx_[phaseIdx];
187 assert(phaseIdx < numPhases);
189 return dnIdx_[phaseIdx];
192 void updateSolvent(
const ElementContext& elemCtx,
unsigned scvfIdx,
unsigned timeIdx)
193 { asImp_().updateVolumeFluxTrans(elemCtx, scvfIdx, timeIdx); }
195 void updatePolymer(
const ElementContext& elemCtx,
unsigned scvfIdx,
unsigned timeIdx)
196 { asImp_().updateShearMultipliers(elemCtx, scvfIdx, timeIdx); }
203 Valgrind::SetUndefined(*
this);
206 static const bool isEcfv = std::is_same<Discretization, EcfvDiscretization<TypeTag> >::value;
208 static_assert(isEcfv);
210 const auto& stencil = elemCtx.stencil(timeIdx);
211 const auto& scvf = stencil.interiorFace(scvfIdx);
213 interiorDofIdx_ = scvf.interiorIndex();
214 exteriorDofIdx_ = scvf.exteriorIndex();
215 assert(interiorDofIdx_ != exteriorDofIdx_);
217 unsigned I = stencil.globalSpaceIndex(interiorDofIdx_);
218 unsigned J = stencil.globalSpaceIndex(exteriorDofIdx_);
220 Scalar trans = transmissibility_(elemCtx, scvfIdx, timeIdx);
225 Scalar g = elemCtx.problem().gravity()[dimWorld - 1];
227 const auto& intQuantsIn = elemCtx.intensiveQuantities(interiorDofIdx_, timeIdx);
228 const auto& intQuantsEx = elemCtx.intensiveQuantities(exteriorDofIdx_, timeIdx);
230 Scalar zIn = dofCenterDepth_(elemCtx, interiorDofIdx_, timeIdx);
231 Scalar zEx = dofCenterDepth_(elemCtx, exteriorDofIdx_, timeIdx);
234 Scalar distZ = zIn - zEx;
236 for (
unsigned phaseIdx=0; phaseIdx < numPhases; phaseIdx++) {
237 if (!FluidSystem::phaseIsActive(phaseIdx))
242 if (intQuantsIn.mobility(phaseIdx) <= 0.0 &&
243 intQuantsEx.mobility(phaseIdx) <= 0.0)
245 upIdx_[phaseIdx] = interiorDofIdx_;
246 dnIdx_[phaseIdx] = exteriorDofIdx_;
247 pressureDifference_[phaseIdx] = 0.0;
248 volumeFlux_[phaseIdx] = 0.0;
254 const Evaluation& rhoIn = intQuantsIn.fluidState().density(phaseIdx);
255 Scalar rhoEx = Toolbox::value(intQuantsEx.fluidState().density(phaseIdx));
256 Evaluation rhoAvg = (rhoIn + rhoEx)/2;
258 const Evaluation& pressureInterior = intQuantsIn.fluidState().pressure(phaseIdx);
259 Evaluation pressureExterior = Toolbox::value(intQuantsEx.fluidState().pressure(phaseIdx));
261 pressureExterior += rhoAvg*(distZ*g);
263 pressureDifference_[phaseIdx] = pressureExterior - pressureInterior;
269 if (pressureDifference_[phaseIdx] > 0.0) {
270 upIdx_[phaseIdx] = exteriorDofIdx_;
271 dnIdx_[phaseIdx] = interiorDofIdx_;
273 else if (pressureDifference_[phaseIdx] < 0.0) {
274 upIdx_[phaseIdx] = interiorDofIdx_;
275 dnIdx_[phaseIdx] = exteriorDofIdx_;
280 Scalar Vin = elemCtx.dofVolume(interiorDofIdx_, 0);
281 Scalar Vex = elemCtx.dofVolume(exteriorDofIdx_, 0);
283 upIdx_[phaseIdx] = interiorDofIdx_;
284 dnIdx_[phaseIdx] = exteriorDofIdx_;
286 else if (Vin < Vex) {
287 upIdx_[phaseIdx] = exteriorDofIdx_;
288 dnIdx_[phaseIdx] = interiorDofIdx_;
295 upIdx_[phaseIdx] = interiorDofIdx_;
296 dnIdx_[phaseIdx] = exteriorDofIdx_;
299 upIdx_[phaseIdx] = exteriorDofIdx_;
300 dnIdx_[phaseIdx] = interiorDofIdx_;
309 const auto& up = elemCtx.intensiveQuantities(upstreamIdx, timeIdx);
311 if (upstreamIdx == interiorDofIdx_)
312 volumeFlux_[phaseIdx] =
313 pressureDifference_[phaseIdx]*up.mobility(phaseIdx)*(-trans);
315 volumeFlux_[phaseIdx] =
316 pressureDifference_[phaseIdx]*(Toolbox::value(up.mobility(phaseIdx))*(-trans));
323 template <
class Flu
idState>
327 const FluidState& exFluidState)
329 const auto& stencil = elemCtx.stencil(timeIdx);
330 const auto& scvf = stencil.boundaryFace(scvfIdx);
332 interiorDofIdx_ = scvf.interiorIndex();
334 Scalar trans = transmissibilityBoundary_(elemCtx, scvfIdx, timeIdx);
339 Scalar g = elemCtx.problem().gravity()[dimWorld - 1];
341 const auto& intQuantsIn = elemCtx.intensiveQuantities(interiorDofIdx_, timeIdx);
348 Scalar zIn = dofCenterDepth_(elemCtx, interiorDofIdx_, timeIdx);
349 Scalar zEx = scvf.integrationPos()[dimWorld - 1];
353 Scalar distZ = zIn - zEx;
355 for (
unsigned phaseIdx=0; phaseIdx < numPhases; phaseIdx++) {
356 if (!FluidSystem::phaseIsActive(phaseIdx))
361 const Evaluation& rhoIn = intQuantsIn.fluidState().density(phaseIdx);
362 const auto& rhoEx = exFluidState.density(phaseIdx);
363 Evaluation rhoAvg = (rhoIn + rhoEx)/2;
365 const Evaluation& pressureInterior = intQuantsIn.fluidState().pressure(phaseIdx);
366 Evaluation pressureExterior = exFluidState.pressure(phaseIdx);
367 pressureExterior += rhoAvg*(distZ*g);
369 pressureDifference_[phaseIdx] = pressureExterior - pressureInterior;
375 if (pressureDifference_[phaseIdx] > 0.0) {
376 upIdx_[phaseIdx] = -1;
377 dnIdx_[phaseIdx] = interiorDofIdx_;
380 upIdx_[phaseIdx] = interiorDofIdx_;
381 dnIdx_[phaseIdx] = -1;
385 if (upstreamIdx == interiorDofIdx_) {
390 const auto& up = elemCtx.intensiveQuantities(upstreamIdx, timeIdx);
392 volumeFlux_[phaseIdx] =
393 pressureDifference_[phaseIdx]*up.mobility(phaseIdx)*(-trans);
399 const auto& matParams =
400 elemCtx.problem().materialLawParams(elemCtx,
403 std::array<typename FluidState::Scalar,numPhases> kr;
404 MaterialLaw::relativePermeabilities(kr, matParams, exFluidState);
406 const auto& mob = kr[phaseIdx]/exFluidState.viscosity(phaseIdx);
407 volumeFlux_[phaseIdx] =
408 pressureDifference_[phaseIdx]*mob*(-trans);
425 Scalar transmissibility_(
const ElementContext& elemCtx,
unsigned scvfIdx,
unsigned timeIdx)
const
427 const auto& stencil = elemCtx.stencil(timeIdx);
428 const auto& face = stencil.interiorFace(scvfIdx);
429 const auto& interiorPos = stencil.subControlVolume(face.interiorIndex()).globalPos();
430 const auto& exteriorPos = stencil.subControlVolume(face.exteriorIndex()).globalPos();
431 auto distVec0 = face.integrationPos() - interiorPos;
432 auto distVec1 = face.integrationPos() - exteriorPos;
433 Scalar ndotDistIn = std::abs(face.normal() * distVec0);
434 Scalar ndotDistExt = std::abs(face.normal() * distVec1);
436 Scalar distSquaredIn = distVec0 * distVec0;
437 Scalar distSquaredExt = distVec1 * distVec1;
438 const auto& K0mat = elemCtx.problem().intrinsicPermeability(elemCtx, face.interiorIndex(), timeIdx);
439 const auto& K1mat = elemCtx.problem().intrinsicPermeability(elemCtx, face.exteriorIndex(), timeIdx);
445 for (
unsigned i = 0; i < dimWorld; ++ i){
446 if (std::abs(face.normal()[i]) > val) {
447 val = std::abs(face.normal()[i]);
451 const Scalar& K0 = K0mat[idx][idx];
452 const Scalar& K1 = K1mat[idx][idx];
453 const Scalar T0 = K0 * ndotDistIn / distSquaredIn;
454 const Scalar T1 = K1 * ndotDistExt / distSquaredExt;
455 return T0 * T1 / (T0 + T1);
457 Scalar transmissibilityBoundary_(
const ElementContext& elemCtx,
unsigned scvfIdx,
unsigned timeIdx)
const
459 const auto& stencil = elemCtx.stencil(timeIdx);
460 const auto& face = stencil.interiorFace(scvfIdx);
461 const auto& interiorPos = stencil.subControlVolume(face.interiorIndex()).globalPos();
462 auto distVec0 = face.integrationPos() - interiorPos;
463 Scalar ndotDistIn = face.normal() * distVec0;
464 Scalar distSquaredIn = distVec0 * distVec0;
465 const auto& K0mat = elemCtx.problem().intrinsicPermeability(elemCtx, face.interiorIndex(), timeIdx);
471 for (
unsigned i = 0; i < dimWorld; ++ i){
472 if (std::abs(face.normal()[i]) > val) {
473 val = std::abs(face.normal()[i]);
477 const Scalar& K0 = K0mat[idx][idx];
478 const Scalar T0 = K0 * ndotDistIn / distSquaredIn;
482 template <
class Context>
483 Scalar dofCenterDepth_(
const Context& context,
unsigned spaceIdx,
unsigned timeIdx)
const
485 const auto& pos = context.pos(spaceIdx, timeIdx);
486 return pos[dimWorld-1];
489 Implementation& asImp_()
490 {
return *
static_cast<Implementation*
>(
this); }
492 const Implementation& asImp_()
const
493 {
return *
static_cast<const Implementation*
>(
this); }
496 Evaluation volumeFlux_[numPhases];
500 Evaluation pressureDifference_[numPhases];
503 unsigned short interiorDofIdx_;
504 unsigned short exteriorDofIdx_;
505 short upIdx_[numPhases];
506 short dnIdx_[numPhases];
Provides the defaults for the parameters required by the transmissibility based volume flux calculati...
Definition: transfluxmodule.hh:73
Provides the transmissibility based flux module.
Definition: transfluxmodule.hh:92
const Evaluation & volumeFlux(unsigned phaseIdx) const
Return the volume flux of a fluid phase at the face's integration point .
Definition: transfluxmodule.hh:160
unsigned downstreamIndex_(unsigned phaseIdx) const
Returns the local index of the degree of freedom in which is in downstream direction.
Definition: transfluxmodule.hh:185
const EvalDimVector & filterVelocity(unsigned) const
Return the filter velocity of a fluid phase at the face's integration point [m/s].
Definition: transfluxmodule.hh:146
void calculateGradients_(const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Update the required gradients for interior faces.
Definition: transfluxmodule.hh:201
void updateSolvent(const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Definition: transfluxmodule.hh:192
const Evaluation & pressureDifference(unsigned phaseIdx) const
Return the gravity corrected pressure difference between the interior and the exterior of a face.
Definition: transfluxmodule.hh:137
void calculateFluxes_(const ElementContext &, unsigned, unsigned)
Update the volumetric fluxes for all fluid phases on the interior faces of the context.
Definition: transfluxmodule.hh:417
void calculateBoundaryFluxes_(const ElementContext &, unsigned, unsigned)
Definition: transfluxmodule.hh:420
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:126
void calculateBoundaryGradients_(const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx, const FluidState &exFluidState)
Update the required gradients for boundary faces.
Definition: transfluxmodule.hh:324
const DimMatrix & intrinsicPermeability() const
Return the intrinsic permeability tensor at a face [m^2].
Definition: transfluxmodule.hh:115
unsigned upstreamIndex_(unsigned phaseIdx) const
Returns the local index of the degree of freedom in which is in upstream direction.
Definition: transfluxmodule.hh:171
void updatePolymer(const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Definition: transfluxmodule.hh:195
Provides the intensive quantities for the transmissibility based flux module.
Definition: transfluxmodule.hh:80
void update_(const ElementContext &, unsigned, unsigned)
Definition: transfluxmodule.hh:83
Defines the common properties required by the porous medium multi-phase models.
Definition: blackoilboundaryratevector.hh:37
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:235
Specifies a flux module which uses transmissibilities.
Definition: transfluxmodule.hh:56
static void registerParameters()
Register all run-time parameters for the flux module.
Definition: transfluxmodule.hh:63