24#ifndef DUMUX_CAPILLARYDIFFUSION_HH
25#define DUMUX_CAPILLARYDIFFUSION_HH
47template<
class TypeTag>
57 using MaterialLaw =
typename SpatialParams::MaterialLaw;
66 using PrimaryVariables =
typename SolutionTypes::PrimaryVariables;
70 dim = GridView::dimension, dimWorld = GridView::dimensionworld
74 wPhaseIdx = Indices::wPhaseIdx, nPhaseIdx = Indices::nPhaseIdx,
75 pressEqIdx = Indices::pressureEqIdx
78 using Intersection =
typename GridView::Intersection;
79 using DimVector = Dune::FieldVector<Scalar, dim>;
80 using GlobalPosition = Dune::FieldVector<Scalar, dimWorld>;
81 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
94 void getFlux (DimVector& flux,
const Intersection& intersection, Scalar satI, Scalar satJ,
95 const DimVector& pcGradient)
const
97 auto element = intersection.inside();
99 const GlobalPosition& globalPos = element.geometry().center();
101 int globalIdxI = problem_.variables().index(element);
102 CellData& CellDataI = problem_.variables().cellData(globalIdxI);
107 Scalar
temperature = problem_.temperature(element);
108 Scalar referencePressure = problem_.referencePressure(element);
112 Scalar mobilityWI = 0;
113 Scalar mobilityNwI = 0;
117 mobilityWI = CellDataI.mobility(wPhaseIdx);
118 mobilityNwI = CellDataI.mobility(nPhaseIdx);
122 FluidState fluidState;
123 fluidState.setPressure(wPhaseIdx, referencePressure);
124 fluidState.setPressure(nPhaseIdx, referencePressure);
126 mobilityWI = MaterialLaw::krw(problem_.spatialParams().materialLawParams(element), satI);
128 mobilityNwI = MaterialLaw::krn(problem_.spatialParams().materialLawParams(element), satI);
132 DimMatrix meanPermeability(0);
134 if (intersection.neighbor())
137 auto neighbor = intersection.outside();
139 int globalIdxJ = problem_.variables().index(neighbor);
140 CellData& cellDataJ = problem_.variables().cellData(globalIdxJ);
143 const GlobalPosition& globalPosNeighbor = neighbor.geometry().center();
146 DimVector distVec = globalPosNeighbor - globalPos;
149 Scalar dist = distVec.two_norm();
151 DimVector unitDistVec(distVec);
155 problem_.spatialParams().meanK(meanPermeability,
156 problem_.spatialParams().intrinsicPermeability(element),
157 problem_.spatialParams().intrinsicPermeability(neighbor));
160 Scalar mobilityWJ = 0;
161 Scalar mobilityNwJ = 0;
165 mobilityWJ = cellDataJ.mobility(wPhaseIdx);
166 mobilityNwJ = cellDataJ.mobility(nPhaseIdx);
170 FluidState fluidState;
171 fluidState.setPressure(wPhaseIdx, referencePressure);
172 fluidState.setPressure(nPhaseIdx, referencePressure);
175 mobilityWJ = MaterialLaw::krw(problem_.spatialParams().materialLawParams(neighbor), satJ);
177 mobilityNwJ = MaterialLaw::krn(problem_.spatialParams().materialLawParams(neighbor), satJ);
180 Scalar mobilityWMean = 0.5*(mobilityWI + mobilityWJ);
181 Scalar mobilityNwMean = 0.5*(mobilityNwI + mobilityNwJ);
182 mobBar = mobilityWMean*mobilityNwMean/(mobilityWMean+mobilityNwMean);
186 BoundaryTypes bcTypes;
187 problem_.boundaryTypes(bcTypes, intersection);
188 if (bcTypes.isNeumann(pressEqIdx))
190 PrimaryVariables priVars;
191 problem_.neumann(priVars, intersection);
192 if (priVars[wPhaseIdx] == 0)
199 problem_.spatialParams().meanK(meanPermeability,
200 problem_.spatialParams().intrinsicPermeability(element));
202 Scalar mobilityWJ = 0;
203 Scalar mobilityNwJ = 0;
206 FluidState fluidState;
207 fluidState.setPressure(wPhaseIdx, referencePressure);
208 fluidState.setPressure(nPhaseIdx, referencePressure);
210 mobilityWJ = MaterialLaw::krw(problem_.spatialParams().materialLawParams(element), satJ);
212 mobilityNwJ = MaterialLaw::krn(problem_.spatialParams().materialLawParams(element), satJ);
215 Scalar mobWMean = 0.5 * (mobilityWI + mobilityWJ);
216 Scalar mobNwMean = 0.5 * (mobilityNwI + mobilityNwJ);
218 mobBar = mobWMean * mobNwMean / (mobWMean + mobNwMean);
222 meanPermeability.mv(pcGradient, flux);
234 :
DiffusivePart<TypeTag>(problem), problem_(problem), preComput_(
getPropValue<TypeTag, Properties::PrecomputedConstRels>())
239 const bool preComput_;
Base class for defining a diffusive 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 GetProp
get the type of a property (equivalent to old macro GET_PROP(...))
Definition: propertysystem.hh:140
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 temperature() noexcept
I/O name of temperature for equilibrium models.
Definition: name.hh:51
std::string viscosity(int phaseIdx) noexcept
I/O name of viscosity for multiphase systems.
Definition: name.hh:74
Class for defining the diffusive capillary pressure term of a saturation equation.
Definition: capillarydiffusion.hh:49
void getFlux(DimVector &flux, const Intersection &intersection, Scalar satI, Scalar satJ, const DimVector &pcGradient) const
Returns capillary diffusion term for current element face.
Definition: capillarydiffusion.hh:94
CapillaryDiffusion(Problem &problem)
Constructs a CapillaryDiffusion object.
Definition: capillarydiffusion.hh:233
Base class for defining the diffusive part of the saturation transport equation.
Definition: diffusivepart.hh:38
Base file for properties related to sequential models.