25#ifndef DUMUX_2P1C_SPURIOUS_FLUX_BLOCKING_DARCYS_LAW_HH
26#define DUMUX_2P1C_SPURIOUS_FLUX_BLOCKING_DARCYS_LAW_HH
37 template<
class TypeTag,
class MyTypeTag>
47template <
class TypeTag>
53 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
54 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
56 using ElemFluxVarCache =
typename GridFluxVariablesCache::LocalView;
63 using Element =
typename GridView::template Codim<0>::Entity;
64 using CoordScalar =
typename GridView::ctype;
66 enum { dim = GridView::dimension};
67 enum { dimWorld = GridView::dimensionworld};
72 liquidPhaseIdx = FluidSystem::liquidPhaseIdx,
73 gasPhaseIdx = FluidSystem::gasPhaseIdx,
76 using DimWorldMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
77 using GlobalPosition = Dune::FieldVector<Scalar, dimWorld>;
82 static Scalar
flux(
const Problem& problem,
83 const Element& element,
84 const FVElementGeometry& fvGeometry,
85 const ElementVolumeVariables& elemVolVars,
86 const SubControlVolumeFace& scvf,
88 const ElemFluxVarCache& elemFluxVarCache)
90 const Scalar
flux = ParentType::flux(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, elemFluxVarCache);
93 if((!getPropValue<TypeTag, Properties::UseBlockingOfSpuriousFlow>()) || phaseIdx != liquidPhaseIdx)
96 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
97 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
99 const Scalar factor = std::signbit(
flux) ? factor_(outsideVolVars, insideVolVars, phaseIdx) :
100 factor_(insideVolVars, outsideVolVars, phaseIdx);
102 return flux * factor;
113 static Scalar factor_(
const VolumeVariables &up,
const VolumeVariables &dn,
const int phaseIdx)
117 const Scalar tDn = dn.temperature();
118 const Scalar tUp = up.temperature();
120 const Scalar sgDn = dn.saturation(gasPhaseIdx);
121 const Scalar sgUp = up.saturation(gasPhaseIdx);
123 bool upIsNotSteam =
false;
124 bool downIsSteam =
false;
125 bool spuriousFlow =
false;
133 if(upIsNotSteam && downIsSteam && tDn > tUp && phaseIdx == liquidPhaseIdx)
138 Scalar deltaT = tDn - tUp;
144 factor = 1-(deltaT/15);
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Define some often used mathematical functions.
The available discretization methods in Dumux.
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
a tag to mark properties as undefined
Definition: propertysystem.hh:35
forward declaration of the method-specific implementation
Definition: flux/box/darcyslaw.hh:44
Determines whether blocking of spurious flow is used or not.
Definition: porousmediumflow/2p1c/darcyslaw.hh:38
Specialization of Darcy's Law for the two-phase one-component model, including a the possibility to...
Definition: porousmediumflow/2p1c/darcyslaw.hh:49
static Scalar flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const int phaseIdx, const ElemFluxVarCache &elemFluxVarCache)
Compute the Darcy flux and use a blocking factor, if specified.
Definition: porousmediumflow/2p1c/darcyslaw.hh:82
Declares all properties used in Dumux.
Advective fluxes according to Darcy's law.