24#ifndef DUMUX_CAPILLARYDIFFUSION_HH
25#define DUMUX_CAPILLARYDIFFUSION_HH
49template<
class TypeTag>
67 using PrimaryVariables =
typename SolutionTypes::PrimaryVariables;
71 dim = GridView::dimension, dimWorld = GridView::dimensionworld
75 wPhaseIdx = Indices::wPhaseIdx, nPhaseIdx = Indices::nPhaseIdx,
76 pressEqIdx = Indices::pressureEqIdx
79 using Intersection =
typename GridView::Intersection;
80 using DimVector = Dune::FieldVector<Scalar, dim>;
81 using GlobalPosition = Dune::FieldVector<Scalar, dimWorld>;
82 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
95 void getFlux (DimVector& flux,
const Intersection& intersection, Scalar satI, Scalar satJ,
96 const DimVector& pcGradient)
const
98 auto element = intersection.inside();
100 const GlobalPosition& globalPos = element.geometry().center();
102 int globalIdxI = problem_.variables().index(element);
103 CellData& CellDataI = problem_.variables().cellData(globalIdxI);
108 Scalar
temperature = problem_.temperature(element);
109 Scalar referencePressure = problem_.referencePressure(element);
113 Scalar mobilityWI = 0;
114 Scalar mobilityNwI = 0;
119 const auto fluidMatrixInteraction = Deprecated::makePcKrSw(Scalar{}, problem_.spatialParams(), element);
123 mobilityWI = CellDataI.mobility(wPhaseIdx);
124 mobilityNwI = CellDataI.mobility(nPhaseIdx);
128 FluidState fluidState;
129 fluidState.setPressure(wPhaseIdx, referencePressure);
130 fluidState.setPressure(nPhaseIdx, referencePressure);
132 mobilityWI = fluidMatrixInteraction.krw(satI);
134 mobilityNwI = fluidMatrixInteraction.krn(satI);
138 DimMatrix meanPermeability(0);
140 if (intersection.neighbor())
143 auto neighbor = intersection.outside();
145 int globalIdxJ = problem_.variables().index(neighbor);
146 CellData& cellDataJ = problem_.variables().cellData(globalIdxJ);
149 const GlobalPosition& globalPosNeighbor = neighbor.geometry().center();
152 DimVector distVec = globalPosNeighbor - globalPos;
155 Scalar dist = distVec.two_norm();
157 DimVector unitDistVec(distVec);
161 problem_.spatialParams().meanK(meanPermeability,
162 problem_.spatialParams().intrinsicPermeability(element),
163 problem_.spatialParams().intrinsicPermeability(neighbor));
166 Scalar mobilityWJ = 0;
167 Scalar mobilityNwJ = 0;
171 mobilityWJ = cellDataJ.mobility(wPhaseIdx);
172 mobilityNwJ = cellDataJ.mobility(nPhaseIdx);
176 FluidState fluidState;
177 fluidState.setPressure(wPhaseIdx, referencePressure);
178 fluidState.setPressure(nPhaseIdx, referencePressure);
184 const auto fluidMatrixInteractionNeighbor = Deprecated::makePcKrSw(Scalar{}, problem_.spatialParams(), neighbor);
186 mobilityWJ = fluidMatrixInteractionNeighbor.krw(satJ);
188 mobilityNwJ = fluidMatrixInteractionNeighbor.krn(satJ);
191 Scalar mobilityWMean = 0.5*(mobilityWI + mobilityWJ);
192 Scalar mobilityNwMean = 0.5*(mobilityNwI + mobilityNwJ);
193 mobBar = mobilityWMean*mobilityNwMean/(mobilityWMean+mobilityNwMean);
197 BoundaryTypes bcTypes;
198 problem_.boundaryTypes(bcTypes, intersection);
199 if (bcTypes.isNeumann(pressEqIdx))
201 PrimaryVariables priVars;
202 problem_.neumann(priVars, intersection);
203 if (priVars[wPhaseIdx] == 0)
210 problem_.spatialParams().meanK(meanPermeability,
211 problem_.spatialParams().intrinsicPermeability(element));
213 Scalar mobilityWJ = 0;
214 Scalar mobilityNwJ = 0;
217 FluidState fluidState;
218 fluidState.setPressure(wPhaseIdx, referencePressure);
219 fluidState.setPressure(nPhaseIdx, referencePressure);
221 mobilityWJ = fluidMatrixInteraction.krw(satJ);
223 mobilityNwJ = fluidMatrixInteraction.krn(satJ);
226 Scalar mobWMean = 0.5 * (mobilityWI + mobilityWJ);
227 Scalar mobNwMean = 0.5 * (mobilityNwI + mobilityNwJ);
229 mobBar = mobWMean * mobNwMean / (mobWMean + mobNwMean);
233 meanPermeability.mv(pcGradient, flux);
245 :
DiffusivePart<TypeTag>(problem), problem_(problem), preComput_(
getPropValue<TypeTag, Properties::PrecomputedConstRels>())
250 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:51
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:95
CapillaryDiffusion(Problem &problem)
Constructs a CapillaryDiffusion object.
Definition: capillarydiffusion.hh:244
Base class for defining the diffusive part of the saturation transport equation.
Definition: diffusivepart.hh:38
Base file for properties related to sequential models.