24#ifndef DUMUX_GRAVITYPART_HH
25#define DUMUX_GRAVITYPART_HH
47template<
class TypeTag>
54 using Indices =
typename GET_PROP_TYPE(TypeTag, ModelTraits)::Indices;
56 using SpatialParams =
typename GET_PROP_TYPE(TypeTag, SpatialParams);
57 using MaterialLaw =
typename SpatialParams::MaterialLaw;
59 using FluidSystem =
typename GET_PROP_TYPE(TypeTag, FluidSystem);
60 using FluidState =
typename GET_PROP_TYPE(TypeTag, FluidState);
66 dim = GridView::dimension, dimWorld = GridView::dimensionworld
73 using Intersection =
typename GridView::Intersection;
74 using DimVector = Dune::FieldVector<Scalar, dim>;
75 using GlobalPosition = Dune::FieldVector<Scalar, dimWorld>;
76 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
88 void getFlux(DimVector& flux,
const Intersection& intersection,
const Scalar satI,
const Scalar satJ)
const
90 auto element = intersection.inside();
92 int globalIdxI = problem_.variables().index(element);
93 CellData& cellDataI = problem_.variables().cellData(globalIdxI);
95 int indexInInside = intersection.indexInInside();
101 Scalar lambdaNwJ = 0;
105 lambdaWI=cellDataI.mobility(wPhaseIdx);
106 lambdaNwI=cellDataI.mobility(nPhaseIdx);
110 lambdaWI = MaterialLaw::krw(problem_.spatialParams().materialLawParams(element), satI);
111 lambdaWI /= viscosity_[wPhaseIdx];
112 lambdaNwI = MaterialLaw::krn(problem_.spatialParams().materialLawParams(element), satI);
113 lambdaNwI /= viscosity_[nPhaseIdx];
116 Scalar potentialDiffW = cellDataI.fluxData().upwindPotential(wPhaseIdx, indexInInside);
117 Scalar potentialDiffNw = cellDataI.fluxData().upwindPotential(nPhaseIdx, indexInInside);
119 DimMatrix meanPermeability(0);
120 GlobalPosition distVec(0);
124 if (intersection.neighbor())
127 auto neighbor = intersection.outside();
129 int globalIdxJ = problem_.variables().index(neighbor);
130 CellData& cellDataJ = problem_.variables().cellData(globalIdxJ);
132 distVec = neighbor.geometry().center() - element.geometry().center();
135 problem_.spatialParams().meanK(meanPermeability,
136 problem_.spatialParams().intrinsicPermeability(element),
137 problem_.spatialParams().intrinsicPermeability(neighbor));
142 lambdaWJ=cellDataJ.mobility(wPhaseIdx);
143 lambdaNwJ=cellDataJ.mobility(nPhaseIdx);
147 lambdaWJ = MaterialLaw::krw(problem_.spatialParams().materialLawParams(neighbor), satJ);
148 lambdaWJ /= viscosity_[wPhaseIdx];
149 lambdaNwJ = MaterialLaw::krn(problem_.spatialParams().materialLawParams(neighbor), satJ);
150 lambdaNwJ /= viscosity_[nPhaseIdx];
153 lambdaW = (potentialDiffW >= 0) ? lambdaWI : lambdaWJ;
154 lambdaW = (potentialDiffW == 0) ? 0.5 * (lambdaWI + lambdaWJ) : lambdaW;
155 lambdaNw = (potentialDiffNw >= 0) ? lambdaNwI : lambdaNwJ;
156 lambdaNw = (potentialDiffNw == 0) ? 0.5 * (lambdaNwI + lambdaNwJ) : lambdaNw;
161 problem_.spatialParams().meanK(meanPermeability,
162 problem_.spatialParams().intrinsicPermeability(element));
164 distVec = intersection.geometry().center() - element.geometry().center();
167 lambdaWJ = MaterialLaw::krw(problem_.spatialParams().materialLawParams(element), satJ);
168 lambdaWJ /= viscosity_[wPhaseIdx];
169 lambdaNwJ = MaterialLaw::krn(problem_.spatialParams().materialLawParams(element), satJ);
170 lambdaNwJ /= viscosity_[nPhaseIdx];
173 lambdaW = (potentialDiffW > 0) ? lambdaWI : lambdaWJ;
174 lambdaNw = (potentialDiffNw > 0) ? lambdaNwI : lambdaNwJ;
178 const Dune::FieldVector<Scalar, dim>& unitOuterNormal = intersection.centerUnitOuterNormal();
179 Scalar dist = distVec.two_norm();
182 Scalar areaScaling = (unitOuterNormal * distVec);
189 Scalar scalarGravity = problem_.gravity() * distVec;
191 flux = unitOuterNormal;
194 flux *= lambdaW*lambdaNw/(lambdaW+lambdaNw) * scalarPerm * (density_[wPhaseIdx] - density_[nPhaseIdx]) * scalarGravity * areaScaling;
209 auto element = *problem_.gridView().template begin<0> ();
210 FluidState fluidState;
211 fluidState.setPressure(wPhaseIdx, problem_.referencePressure(element));
212 fluidState.setPressure(nPhaseIdx, problem_.referencePressure(element));
213 fluidState.setTemperature(problem_.temperature(element));
214 fluidState.setSaturation(wPhaseIdx, 1.);
215 fluidState.setSaturation(nPhaseIdx, 0.);
224 const bool preComput_;
226 Scalar density_[numPhases];
227 Scalar viscosity_[numPhases];
#define GET_PROP_VALUE(TypeTag, PropTagName)
Definition: propertysystemmacros.hh:282
#define GET_PROP_TYPE(TypeTag, PropTagName)
Definition: propertysystemmacros.hh:283
Base class for defining a convective part of the saturation transport equation.
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
Property tag PrecomputedConstRels
Definition: porousmediumflow/2p/sequential/transport/cellcentered/properties.hh:46
Property tag NumPhases
Number of phases in the system.
Definition: porousmediumflow/sequential/properties.hh:69
std::string viscosity(int phaseIdx) noexcept
I/O name of viscosity for multiphase systems.
Definition: name.hh:74
std::string permeability() noexcept
I/O name of permeability.
Definition: name.hh:143
std::string density(int phaseIdx) noexcept
I/O name of density for multiphase systems.
Definition: name.hh:65
Base class for defining a convective part of the saturation transport equation.
Definition: convectivepart.hh:39
Class for defining the gravity term of a two-phase flow saturation equation.
Definition: gravitypart.hh:49
void initialize()
For initialization.
Definition: gravitypart.hh:207
GravityPart(Problem &problem)
Constructs a GravityPart object.
Definition: gravitypart.hh:202
void getFlux(DimVector &flux, const Intersection &intersection, const Scalar satI, const Scalar satJ) const
Returns convective term for current element face.
Definition: gravitypart.hh:88
Base file for properties related to sequential models.