24#ifndef DUMUX_CAPILLARYDIFFUSION_HH
25#define DUMUX_CAPILLARYDIFFUSION_HH
47template<
class TypeTag>
65 using PrimaryVariables =
typename SolutionTypes::PrimaryVariables;
69 dim = GridView::dimension, dimWorld = GridView::dimensionworld
73 wPhaseIdx = Indices::wPhaseIdx, nPhaseIdx = Indices::nPhaseIdx,
74 pressEqIdx = Indices::pressureEqIdx
77 using Intersection =
typename GridView::Intersection;
78 using DimVector = Dune::FieldVector<Scalar, dim>;
79 using GlobalPosition = Dune::FieldVector<Scalar, dimWorld>;
80 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
93 void getFlux (DimVector& flux,
const Intersection& intersection, Scalar satI, Scalar satJ,
94 const DimVector& pcGradient)
const
96 auto element = intersection.inside();
98 const GlobalPosition& globalPos = element.geometry().center();
100 int globalIdxI = problem_.variables().index(element);
101 CellData& CellDataI = problem_.variables().cellData(globalIdxI);
106 Scalar
temperature = problem_.temperature(element);
107 Scalar referencePressure = problem_.referencePressure(element);
111 Scalar mobilityWI = 0;
112 Scalar mobilityNwI = 0;
114 const auto fluidMatrixInteraction = problem_.spatialParams().fluidMatrixInteractionAtPos(element.geometry().center());
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 = fluidMatrixInteraction.krw(satI);
128 mobilityNwI = fluidMatrixInteraction.krn(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 const auto fluidMatrixInteractionNeighbor = problem_.spatialParams().fluidMatrixInteractionAtPos(neighbor.geometry().center());
176 mobilityWJ = fluidMatrixInteractionNeighbor.krw(satJ);
178 mobilityNwJ = fluidMatrixInteractionNeighbor.krn(satJ);
181 Scalar mobilityWMean = 0.5*(mobilityWI + mobilityWJ);
182 Scalar mobilityNwMean = 0.5*(mobilityNwI + mobilityNwJ);
183 mobBar = mobilityWMean*mobilityNwMean/(mobilityWMean+mobilityNwMean);
187 BoundaryTypes bcTypes;
188 problem_.boundaryTypes(bcTypes, intersection);
189 if (bcTypes.isNeumann(pressEqIdx))
191 PrimaryVariables priVars;
192 problem_.neumann(priVars, intersection);
193 if (priVars[wPhaseIdx] == 0)
200 problem_.spatialParams().meanK(meanPermeability,
201 problem_.spatialParams().intrinsicPermeability(element));
203 Scalar mobilityWJ = 0;
204 Scalar mobilityNwJ = 0;
207 FluidState fluidState;
208 fluidState.setPressure(wPhaseIdx, referencePressure);
209 fluidState.setPressure(nPhaseIdx, referencePressure);
211 mobilityWJ = fluidMatrixInteraction.krw(satJ);
213 mobilityNwJ = fluidMatrixInteraction.krn(satJ);
216 Scalar mobWMean = 0.5 * (mobilityWI + mobilityWJ);
217 Scalar mobNwMean = 0.5 * (mobilityNwI + mobilityNwJ);
219 mobBar = mobWMean * mobNwMean / (mobWMean + mobNwMean);
223 meanPermeability.mv(pcGradient, flux);
235 :
DiffusivePart<TypeTag>(problem), problem_(problem), preComput_(
getPropValue<TypeTag, Properties::PrecomputedConstRels>())
240 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:154
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type GetProp
get the type of a property
Definition: propertysystem.hh:141
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:150
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:93
CapillaryDiffusion(Problem &problem)
Constructs a CapillaryDiffusion object.
Definition: capillarydiffusion.hh:234
Base class for defining the diffusive part of the saturation transport equation.
Definition: diffusivepart.hh:38
Base file for properties related to sequential models.