24#ifndef DUMUX_GRAVITYPART_HH
25#define DUMUX_GRAVITYPART_HH
49template<
class TypeTag>
67 dim = GridView::dimension, dimWorld = GridView::dimensionworld
71 wPhaseIdx = Indices::wPhaseIdx, nPhaseIdx = Indices::nPhaseIdx, numPhases = getPropValue<TypeTag, Properties::NumPhases>()
74 using Intersection =
typename GridView::Intersection;
75 using DimVector = Dune::FieldVector<Scalar, dim>;
76 using GlobalPosition = Dune::FieldVector<Scalar, dimWorld>;
77 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
89 void getFlux(DimVector& flux,
const Intersection& intersection,
const Scalar satI,
const Scalar satJ)
const
91 auto element = intersection.inside();
93 int globalIdxI = problem_.variables().index(element);
94 CellData& cellDataI = problem_.variables().cellData(globalIdxI);
96 int indexInInside = intersection.indexInInside();
100 Scalar lambdaNwI = 0;
102 Scalar lambdaNwJ = 0;
107 const auto fluidMatrixInteraction = Deprecated::makePcKrSw(Scalar{}, problem_.spatialParams(), element);
111 lambdaWI=cellDataI.mobility(wPhaseIdx);
112 lambdaNwI=cellDataI.mobility(nPhaseIdx);
116 lambdaWI = fluidMatrixInteraction.krw(satI);
117 lambdaWI /= viscosity_[wPhaseIdx];
118 lambdaNwI = fluidMatrixInteraction.krn(satI);
119 lambdaNwI /= viscosity_[nPhaseIdx];
122 Scalar potentialDiffW = cellDataI.fluxData().upwindPotential(wPhaseIdx, indexInInside);
123 Scalar potentialDiffNw = cellDataI.fluxData().upwindPotential(nPhaseIdx, indexInInside);
125 DimMatrix meanPermeability(0);
126 GlobalPosition distVec(0);
130 if (intersection.neighbor())
133 auto neighbor = intersection.outside();
135 int globalIdxJ = problem_.variables().index(neighbor);
136 CellData& cellDataJ = problem_.variables().cellData(globalIdxJ);
138 distVec = neighbor.geometry().center() - element.geometry().center();
141 problem_.spatialParams().meanK(meanPermeability,
142 problem_.spatialParams().intrinsicPermeability(element),
143 problem_.spatialParams().intrinsicPermeability(neighbor));
148 lambdaWJ=cellDataJ.mobility(wPhaseIdx);
149 lambdaNwJ=cellDataJ.mobility(nPhaseIdx);
156 const auto fluidMatrixInteractionNeighbor = Deprecated::makePcKrSw(Scalar{}, problem_.spatialParams(), neighbor);
158 lambdaWJ = fluidMatrixInteractionNeighbor.krw(satJ);
159 lambdaWJ /= viscosity_[wPhaseIdx];
160 lambdaNwJ = fluidMatrixInteractionNeighbor.krn(satJ);
161 lambdaNwJ /= viscosity_[nPhaseIdx];
164 lambdaW = (potentialDiffW >= 0) ? lambdaWI : lambdaWJ;
165 lambdaW = (potentialDiffW == 0) ? 0.5 * (lambdaWI + lambdaWJ) : lambdaW;
166 lambdaNw = (potentialDiffNw >= 0) ? lambdaNwI : lambdaNwJ;
167 lambdaNw = (potentialDiffNw == 0) ? 0.5 * (lambdaNwI + lambdaNwJ) : lambdaNw;
172 problem_.spatialParams().meanK(meanPermeability,
173 problem_.spatialParams().intrinsicPermeability(element));
175 distVec = intersection.geometry().center() - element.geometry().center();
178 lambdaWJ = fluidMatrixInteraction.krw(satJ);
179 lambdaWJ /= viscosity_[wPhaseIdx];
180 lambdaNwJ = fluidMatrixInteraction.krn(satJ);
181 lambdaNwJ /= viscosity_[nPhaseIdx];
184 lambdaW = (potentialDiffW > 0) ? lambdaWI : lambdaWJ;
185 lambdaNw = (potentialDiffNw > 0) ? lambdaNwI : lambdaNwJ;
189 const Dune::FieldVector<Scalar, dim>& unitOuterNormal = intersection.centerUnitOuterNormal();
190 Scalar dist = distVec.two_norm();
193 Scalar areaScaling = (unitOuterNormal * distVec);
200 Scalar scalarGravity = problem_.gravity() * distVec;
202 flux = unitOuterNormal;
205 flux *= lambdaW*lambdaNw/(lambdaW+lambdaNw) * scalarPerm * (density_[wPhaseIdx] - density_[nPhaseIdx]) * scalarGravity * areaScaling;
220 auto element = *problem_.gridView().template begin<0> ();
221 FluidState fluidState;
222 fluidState.setPressure(wPhaseIdx, problem_.referencePressure(element));
223 fluidState.setPressure(nPhaseIdx, problem_.referencePressure(element));
224 fluidState.setTemperature(problem_.temperature(element));
225 fluidState.setSaturation(wPhaseIdx, 1.);
226 fluidState.setSaturation(nPhaseIdx, 0.);
235 const bool preComput_;
237 Scalar density_[numPhases];
238 Scalar viscosity_[numPhases];
Base class for defining a convective part of the saturation transport equation.
constexpr auto getPropValue()
get the value data member of a property
Definition: propertysystem.hh:153
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:149
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:51
void initialize()
For initialization.
Definition: gravitypart.hh:218
GravityPart(Problem &problem)
Constructs a GravityPart object.
Definition: gravitypart.hh:213
void getFlux(DimVector &flux, const Intersection &intersection, const Scalar satI, const Scalar satJ) const
Returns convective term for current element face.
Definition: gravitypart.hh:89
Base file for properties related to sequential models.