24#ifndef DUMUX_FVVELOCITY2P_ADAPTIVE_HH
25#define DUMUX_FVVELOCITY2P_ADAPTIVE_HH
27#include <dune/common/float_cmp.hh>
38template<
class TypeTag>
55 using Intersection =
typename GridView::Intersection;
59 dim = GridView::dimension, dimWorld = GridView::dimensionworld
64 pw = Indices::pressureW,
65 pn = Indices::pressureNw,
66 pGlobal = Indices::pressureGlobal,
67 vw = Indices::velocityW,
68 vn = Indices::velocityNw,
69 vt = Indices::velocityTotal
73 wPhaseIdx = Indices::wPhaseIdx, nPhaseIdx = Indices::nPhaseIdx, numPhases = getPropValue<TypeTag, Properties::NumPhases>()
76 using Element =
typename GridView::template Codim<0>::Entity;
78 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
79 using GravityVector = Dune::FieldVector<Scalar, dimWorld>;
80 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
88 :
ParentType(problem), problem_(problem), gravity_(problem.gravity())
90 if (getPropValue<TypeTag, Properties::EnableCompressibility>() && velocityType_ == vt)
92 DUNE_THROW(Dune::NotImplemented,
"Total velocity - global pressure - model cannot be used with compressible fluids!");
94 if (velocityType_ != vw && velocityType_ != vn && velocityType_ != vt)
96 DUNE_THROW(Dune::NotImplemented,
"Velocity type not supported!");
99 density_[wPhaseIdx] = 0.;
100 density_[nPhaseIdx] = 0.;
101 viscosity_[wPhaseIdx] = 0.;
102 viscosity_[nPhaseIdx] = 0.;
110 if (!compressibility_)
112 const auto element = *problem_.gridView().template begin<0>();
113 FluidState fluidState;
114 fluidState.setPressure(wPhaseIdx, problem_.referencePressure(element));
115 fluidState.setPressure(nPhaseIdx, problem_.referencePressure(element));
116 fluidState.setTemperature(problem_.temperature(element));
117 fluidState.setSaturation(wPhaseIdx, 1.);
118 fluidState.setSaturation(nPhaseIdx, 0.);
142 const GravityVector& gravity_;
143 Scalar density_[numPhases];
144 Scalar viscosity_[numPhases];
147 static const int velocityType_ = getPropValue<TypeTag, Properties::VelocityFormulation>();
148 static const bool compressibility_ = getPropValue<TypeTag, Properties::EnableCompressibility>();
150 static const int pressureType_ = getPropValue<TypeTag, Properties::PressureFormulation>();
152 static const int saturationType_ = getPropValue<TypeTag, Properties::SaturationFormulation>();
162template<
class TypeTag>
165 auto elementI = intersection.inside();
166 auto elementJ = intersection.outside();
168 if (elementI.level() == elementJ.level())
170 ParentType::calculateVelocity(intersection, cellData);
172 else if (elementJ.level() == elementI.level() + 1 && dim == 2)
174 CellData& cellDataJ = problem_.variables().cellData(problem_.variables().index(elementJ));
177 const GlobalPosition& globalPosI = elementI.geometry().center();
178 const GlobalPosition& globalPosJ = elementJ.geometry().center();
181 Scalar lambdaWI = cellData.mobility(wPhaseIdx);
182 Scalar lambdaNwI = cellData.mobility(nPhaseIdx);
183 Scalar fractionalWI = cellData.fracFlowFunc(wPhaseIdx);
184 Scalar fractionalNwI = cellData.fracFlowFunc(nPhaseIdx);
185 Scalar lambdaWJ = cellDataJ.mobility(wPhaseIdx);
186 Scalar lambdaNwJ = cellDataJ.mobility(nPhaseIdx);
187 Scalar fractionalWJ = cellDataJ.fracFlowFunc(wPhaseIdx);
188 Scalar fractionalNwJ = cellDataJ.fracFlowFunc(nPhaseIdx);
191 Scalar pcI = cellData.capillaryPressure();
192 Scalar pcJ = cellDataJ.capillaryPressure();
195 int isIndexI = intersection.indexInInside();
197 Scalar faceArea = intersection.geometry().volume();
198 Scalar faceAreaSum = faceArea;
201 const Dune::FieldVector<Scalar, dim>& unitOuterNormal = intersection.centerUnitOuterNormal();
205 int isIndexJ = intersection.indexInOutside();
208 auto elementK = intersection.outside();
215 for (
const auto& intersectionI : intersections(problem_.gridView(), elementI))
217 if (intersectionI.neighbor())
219 auto neighbor2 = intersectionI.outside();
223 if (neighbor2 != elementJ && intersectionI.indexInInside() == isIndexI)
225 globalIdxK = problem_.variables().index(neighbor2);
226 elementK = neighbor2;
227 faceAreaSum += intersectionI.geometry().volume();
234 CellData& cellDataK = problem_.variables().cellData(globalIdxK);
237 const GlobalPosition& globalPosInterface = intersection.geometry().center();
239 Dune::FieldVector < Scalar, dimWorld > distVec = globalPosInterface - globalPosI;
240 Scalar lI = distVec * unitOuterNormal;
241 distVec = globalPosJ - globalPosInterface;
242 Scalar lJ = distVec * unitOuterNormal;
245 DimMatrix permeabilityI(0);
246 DimMatrix permeabilityJ(0);
247 DimMatrix permeabilityK(0);
249 problem_.spatialParams().meanK(permeabilityI,
250 problem_.spatialParams().intrinsicPermeability(elementI));
251 problem_.spatialParams().meanK(permeabilityJ,
252 problem_.spatialParams().intrinsicPermeability(elementJ));
253 problem_.spatialParams().meanK(permeabilityK,
254 problem_.spatialParams().intrinsicPermeability(elementK));
257 Scalar kI, kJ, kK, ng, kMean;
258 Dune::FieldVector < Scalar, dim > permI(0);
259 Dune::FieldVector < Scalar, dim > permJ(0);
260 Dune::FieldVector < Scalar, dim > permK(0);
262 permeabilityI.mv(unitOuterNormal, permI);
263 permeabilityJ.mv(unitOuterNormal, permJ);
264 permeabilityK.mv(unitOuterNormal, permK);
267 kI = unitOuterNormal * permI;
268 kJ = unitOuterNormal * permJ;
269 kK = unitOuterNormal * permK;
270 kMean = kI * kJ * kK * l / (kJ * kK * lI + kI * (kJ + kK) / 2 * lJ);
272 ng = gravity_ * unitOuterNormal;
274 Scalar fractionalWK = cellDataK.fracFlowFunc(wPhaseIdx);
275 Scalar fractionalNwK = cellDataK.fracFlowFunc(nPhaseIdx);
277 Scalar pcK = cellDataK.capillaryPressure();
278 Scalar pcJK = (pcJ + pcK) / 2;
282 Scalar potentialDiffWIJ = cellData.fluxData().upwindPotential(wPhaseIdx, isIndexI);
283 Scalar potentialDiffNwIJ = cellData.fluxData().upwindPotential(nPhaseIdx, isIndexI);
284 Scalar potentialDiffWIK = potentialDiffWIJ;
285 Scalar potentialDiffNwIK = potentialDiffNwIJ;
289 Scalar rhoMeanWIJ = density_[wPhaseIdx];
290 Scalar rhoMeanNwIJ = density_[nPhaseIdx];
291 Scalar rhoMeanWIK = density_[wPhaseIdx];
292 Scalar rhoMeanNwIK = density_[nPhaseIdx];
294 if (compressibility_)
296 rhoMeanWIJ = (lJ * cellData.density(wPhaseIdx) + lI * cellDataJ.density(wPhaseIdx)) / l;
297 rhoMeanNwIJ = (lJ * cellData.density(nPhaseIdx) + lI * cellDataJ.density(nPhaseIdx)) / l;
298 rhoMeanWIK = (lJ * cellData.density(wPhaseIdx) + lI * cellDataK.density(wPhaseIdx)) / l;
299 rhoMeanNwIK = (lJ * cellData.density(nPhaseIdx) + lI * cellDataK.density(nPhaseIdx)) / l;
302 Scalar fMeanWIJ = (lJ * fractionalWI + lI * fractionalWJ) / l;
303 Scalar fMeanNwIJ = (lJ * fractionalNwI + lI * fractionalNwJ) / l;
304 Scalar fMeanWIK = (lJ * fractionalWI + lI * fractionalWK) / l;
305 Scalar fMeanNwIK = (lJ * fractionalNwI + lI * fractionalNwK) / l;
307 Scalar densityWIJ = density_[wPhaseIdx];
308 Scalar densityNwIJ = density_[nPhaseIdx];
309 Scalar densityWIK = density_[wPhaseIdx];
310 Scalar densityNwIK = density_[nPhaseIdx];
312 if (compressibility_)
315 densityWIJ = (potentialDiffWIJ > 0.) ? cellData.density(wPhaseIdx) : cellDataJ.density(wPhaseIdx);
316 densityNwIJ = (potentialDiffNwIJ > 0.) ? cellData.density(nPhaseIdx) : cellDataJ.density(nPhaseIdx);
317 densityWIK = (potentialDiffWIK > 0.) ? cellData.density(wPhaseIdx) : cellDataK.density(wPhaseIdx);
318 densityNwIK = (potentialDiffNwIK > 0.) ? cellData.density(nPhaseIdx) : cellDataK.density(nPhaseIdx);
320 densityWIJ = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffWIJ, 0.0, 1.0e-30)) ? rhoMeanWIJ : densityWIJ;
321 densityNwIJ = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffNwIJ, 0.0, 1.0e-30)) ? rhoMeanNwIJ : densityNwIJ;
322 densityWIK = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffWIK, 0.0, 1.0e-30)) ? rhoMeanWIK : densityWIK;
323 densityNwIK = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffNwIK, 0.0, 1.0e-30)) ? rhoMeanNwIK : densityNwIK;
326 Scalar fractionalWIJ = (potentialDiffWIJ > 0.) ? fractionalWI : fractionalWJ;
327 Scalar fractionalNwIJ = (potentialDiffNwIJ > 0.) ? fractionalNwI : fractionalNwJ;
328 Scalar fractionalWIK = (potentialDiffWIK > 0.) ? fractionalWI : fractionalWK;
329 Scalar fractionalNwIK = (potentialDiffNwIK > 0.) ? fractionalNwI : fractionalNwK;
331 fractionalWIJ = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffWIJ, 0.0, 1.0e-30)) ? fMeanWIJ : fractionalWIJ;
332 fractionalNwIJ = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffNwIJ, 0.0, 1.0e-30)) ? fMeanNwIJ : fractionalNwIJ;
333 fractionalWIK = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffWIK, 0.0, 1.0e-30)) ? fMeanWIK : fractionalWIK;
334 fractionalNwIK = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffNwIK, 0.0, 1.0e-30)) ? fMeanNwIK : fractionalNwIK;
336 switch (pressureType_)
340 Scalar pressJK = (cellDataJ.globalPressure() + cellDataK.globalPressure()) / 2;
342 potentialDiffWIJ = (cellData.globalPressure() - fMeanNwIJ * pcI - (pressJK - fMeanNwIJ * pcJK)) / l
343 + (densityWIJ - lJ / l * (kI + kK) / kI * (densityWIK - densityWIJ) / 2) * ng;
344 potentialDiffNwIJ = (cellData.globalPressure() + fMeanWIJ * pcI - (pressJK + fMeanWIJ * pcJK)) / l
345 + (densityNwIJ - lJ / l * (kI + kK) / kI * (densityNwIK - densityNwIJ) / 2) * ng;
346 potentialDiffWIK = (cellData.globalPressure() - fMeanNwIK * pcI - (pressJK - fMeanNwIK * pcJK)) / l
347 + (densityWIK - lJ / l * (kI + kK) / kI * (densityWIJ - densityWIK) / 2) * ng;
348 potentialDiffNwIK = (cellData.globalPressure() + fMeanWIK * pcI - (pressJK + fMeanWIK * pcJK)) / l
349 + (densityNwIK - lJ / l * (kI + kK) / kI * (densityNwIJ - densityNwIK) / 2) * ng;
354 potentialDiffWIJ = (cellData.pressure(wPhaseIdx)
355 - 0.5 * (cellDataJ.pressure(wPhaseIdx) + cellDataK.pressure(wPhaseIdx))) / l
356 + (densityWIJ - lJ / l * (kI + kK) / kI * (densityWIK - densityWIJ) / 2) * ng;
357 potentialDiffNwIJ = (cellData.pressure(nPhaseIdx)
358 - (0.5 * (cellDataJ.pressure(nPhaseIdx) + cellDataK.pressure(nPhaseIdx)))) / l
359 + (densityNwIJ - lJ / l * (kI + kK) / kI * (densityNwIK - densityNwIJ) / 2) * ng;
360 potentialDiffWIK = (cellData.pressure(wPhaseIdx)
361 - 0.5 * (cellDataJ.pressure(wPhaseIdx) + cellDataK.pressure(wPhaseIdx))) / l
362 + (densityWIK - lJ / l * (kI + kK) / kI * (densityWIJ - densityWIK) / 2) * ng;
363 potentialDiffNwIK = (cellData.pressure(nPhaseIdx)
364 - (0.5 * (cellDataJ.pressure(nPhaseIdx) + cellDataK.pressure(nPhaseIdx)))) / l
365 + (densityNwIK - lJ / l * (kI + kK) / kI * (densityNwIJ - densityNwIK) / 2) * ng;
373 cellData.fluxData().addUpwindPotential(wPhaseIdx, isIndexI, potentialDiffWIJ);
374 cellData.fluxData().addUpwindPotential(nPhaseIdx, isIndexI, potentialDiffNwIJ);
375 cellDataJ.fluxData().setUpwindPotential(wPhaseIdx, isIndexJ, -potentialDiffWIJ);
376 cellDataJ.fluxData().setUpwindPotential(nPhaseIdx, isIndexJ, -potentialDiffNwIJ);
379 Scalar lambdaWIJ = (potentialDiffWIJ > 0.) ? lambdaWI : lambdaWJ;
380 lambdaWIJ = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffWIJ, 0.0, 1.0e-30)) ? 0.5 * (lambdaWI + lambdaWJ) : lambdaWIJ;
381 Scalar lambdaNwIJ = (potentialDiffNwIJ > 0.) ? lambdaNwI : lambdaNwJ;
382 lambdaNwIJ = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffNwIJ, 0.0, 1.0e-30)) ? 0.5 * (lambdaNwI + lambdaNwJ) : lambdaNwIJ;
384 if (compressibility_)
386 densityWIJ = (potentialDiffWIJ > 0.) ? cellData.density(wPhaseIdx) : cellDataJ.density(wPhaseIdx);
387 densityNwIJ = (potentialDiffNwIJ > 0.) ? cellData.density(nPhaseIdx) : cellDataJ.density(nPhaseIdx);
388 densityWIJ = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffWIJ, 0.0, 1.0e-30)) ? rhoMeanWIJ : densityWIJ;
389 densityNwIJ = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffNwIJ, 0.0, 1.0e-30)) ? rhoMeanNwIJ : densityNwIJ;
390 densityWIK = (potentialDiffWIK > 0.) ? cellData.density(wPhaseIdx) : cellDataK.density(wPhaseIdx);
391 densityNwIK = (potentialDiffNwIK > 0.) ? cellData.density(nPhaseIdx) : cellDataK.density(nPhaseIdx);
392 densityWIK = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffWIK, 0.0, 1.0e-30)) ? rhoMeanWIK : densityWIK;
393 densityNwIK = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffNwIK, 0.0, 1.0e-30)) ? rhoMeanNwIK : densityNwIK;
397 Dune::FieldVector < Scalar, dimWorld > velocityW(unitOuterNormal);
398 Dune::FieldVector < Scalar, dimWorld > velocityNw(unitOuterNormal);
399 Dune::FieldVector < Scalar, dimWorld > gravityTermW(unitOuterNormal);
400 Dune::FieldVector < Scalar, dimWorld > gravityTermNw(unitOuterNormal);
402 gravityTermW *= lambdaWIJ * kMean * ng;
403 gravityTermW *= densityWIJ - (lJ / l) * (kI + kK) / kI * (densityWIK - densityWIJ) / 2;
404 gravityTermNw *= lambdaNwIJ * kMean * ng;
405 gravityTermNw *= densityNwIJ - (lJ / l) * (kI + kK) / kI * (densityNwIK - densityNwIJ) / 2;
407 switch (pressureType_)
411 velocityW *= lambdaWIJ * kMean * (cellData.pressure(wPhaseIdx) -
412 (cellDataJ.pressure(wPhaseIdx) + cellDataK.pressure(wPhaseIdx)) / 2.0) / l;
413 velocityNw *= lambdaNwIJ * kMean * (cellData.pressure(nPhaseIdx) -
414 (cellDataJ.pressure(nPhaseIdx) + cellDataK.pressure(nPhaseIdx)) / 2.0) / l;
416 velocityW += gravityTermW;
417 velocityNw += gravityTermNw;
422 velocityNw *= lambdaNwIJ * kMean * (cellData.pressure(nPhaseIdx) -
423 (cellDataJ.pressure(nPhaseIdx) + cellDataK.pressure(nPhaseIdx)) / 2.0) / l;
424 velocityW *= lambdaWIJ * kMean * (cellData.pressure(wPhaseIdx) -
425 (cellDataJ.pressure(wPhaseIdx) + cellDataK.pressure(wPhaseIdx)) / 2.0) / l;
427 velocityW += gravityTermW;
428 velocityNw += gravityTermNw;
433 Scalar pressJK = (cellDataJ.globalPressure() + cellDataK.globalPressure()) / 2;
435 velocityW *= lambdaWIJ * kMean * (cellData.globalPressure() - pressJK) / l;
436 velocityW += gravityTermW;
437 velocityW += gravityTermNw;
443 cellDataJ.fluxData().setVelocity(wPhaseIdx, isIndexJ, velocityW);
444 cellDataJ.fluxData().setVelocity(nPhaseIdx, isIndexJ, velocityNw);
445 cellDataJ.fluxData().setVelocityMarker(isIndexJ);
448 velocityW *= faceArea/faceAreaSum;
449 velocityNw *= faceArea/faceAreaSum;
450 cellData.fluxData().addVelocity(wPhaseIdx, isIndexI, velocityW);
451 cellData.fluxData().addVelocity(nPhaseIdx, isIndexI, velocityNw);
453 else if (elementI.level() > elementJ.level() && dim == 3)
455 int globalIdxJ = problem_.variables().index(elementJ);
457 CellData& cellDataJ = problem_.variables().cellData(globalIdxJ);
460 const GlobalPosition& globalPosI = elementI.geometry().center();
461 const GlobalPosition& globalPosJ = elementJ.geometry().center();
464 Scalar lambdaWI = cellData.mobility(wPhaseIdx);
465 Scalar lambdaNwI = cellData.mobility(nPhaseIdx);
466 Scalar lambdaWJ = cellDataJ.mobility(wPhaseIdx);
467 Scalar lambdaNwJ = cellDataJ.mobility(nPhaseIdx);
470 Scalar pcI = cellData.capillaryPressure();
471 Scalar pcJ = cellDataJ.capillaryPressure();
474 int isIndexI = intersection.indexInInside();
475 int isIndexJ = intersection.indexInOutside();
478 const Dune::FieldVector<Scalar, dim>& unitOuterNormal = intersection.centerUnitOuterNormal();
481 GlobalPosition distVec = globalPosJ - globalPosI;
484 Scalar dist = distVec.two_norm();
487 DimMatrix meanPermeability(0);
489 problem_.spatialParams().meanK(meanPermeability, problem_.spatialParams().intrinsicPermeability(elementI),
490 problem_.spatialParams().intrinsicPermeability(elementJ));
496 Scalar potentialDiffW = cellData.potential(wPhaseIdx) - cellDataJ.potential(wPhaseIdx);
497 Scalar potentialDiffNw = cellData.potential(nPhaseIdx) - cellDataJ.potential(nPhaseIdx);
499 if (compressibility_)
501 density_[wPhaseIdx] = (potentialDiffW > 0.) ? cellData.density(wPhaseIdx) : cellDataJ.density(wPhaseIdx);
502 density_[nPhaseIdx] = (potentialDiffNw > 0.) ? cellData.density(nPhaseIdx) : cellDataJ.density(nPhaseIdx);
504 density_[wPhaseIdx] =
505 (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffW, 0.0, 1.0e-30)) ? 0.5 * (cellData.density(wPhaseIdx) + cellDataJ.density(wPhaseIdx)) :
507 density_[nPhaseIdx] =
508 (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffNw, 0.0, 1.0e-30)) ? 0.5 * (cellData.density(nPhaseIdx) + cellDataJ.density(nPhaseIdx)) :
511 potentialDiffW = (cellData.pressure(wPhaseIdx) - cellDataJ.pressure(wPhaseIdx));
512 potentialDiffNw = (cellData.pressure(nPhaseIdx) - cellDataJ.pressure(nPhaseIdx));
515 potentialDiffW += density_[wPhaseIdx] * (distVec * gravity_);
516 potentialDiffNw += density_[nPhaseIdx] * (distVec * gravity_);
520 cellData.fluxData().setUpwindPotential(wPhaseIdx, isIndexI, potentialDiffW);
521 cellData.fluxData().setUpwindPotential(nPhaseIdx, isIndexI, potentialDiffNw);
523 cellDataJ.fluxData().setUpwindPotential(wPhaseIdx, isIndexJ, -potentialDiffW);
524 cellDataJ.fluxData().setUpwindPotential(nPhaseIdx, isIndexJ, -potentialDiffNw);
527 Scalar lambdaW = (potentialDiffW > 0.) ? lambdaWI : lambdaWJ;
528 lambdaW = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffW, 0.0, 1.0e-30 )) ? 0.5 * (lambdaWI + lambdaWJ) : lambdaW;
529 Scalar lambdaNw = (potentialDiffNw > 0.) ? lambdaNwI : lambdaNwJ;
530 lambdaNw = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffNw, 0.0, 1.0e-30)) ? 0.5 * (lambdaNwI + lambdaNwJ) : lambdaNw;
532 if (compressibility_)
534 density_[wPhaseIdx] = (potentialDiffW > 0.) ? cellData.density(wPhaseIdx) : cellDataJ.density(wPhaseIdx);
535 density_[nPhaseIdx] = (potentialDiffNw > 0.) ? cellData.density(nPhaseIdx) : cellDataJ.density(nPhaseIdx);
537 density_[wPhaseIdx] =
538 (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffW, 0.0, 1.0e-30)) ? 0.5 * (cellData.density(wPhaseIdx) + cellDataJ.density(wPhaseIdx)) :
540 density_[nPhaseIdx] =
541 (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffNw, 0.0, 1.0e-30)) ? 0.5 * (cellData.density(nPhaseIdx) + cellDataJ.density(nPhaseIdx)) :
548 Dune::FieldVector<Scalar, dimWorld> velocityW(unitOuterNormal);
549 Dune::FieldVector<Scalar, dimWorld> velocityNw(unitOuterNormal);
553 Scalar areaScaling = (unitOuterNormal * distVec);
556 Scalar gravityTermW = (gravity_ * distVec) * density_[wPhaseIdx] * areaScaling;
557 Scalar gravityTermNw = (gravity_ * distVec) * density_[nPhaseIdx] * areaScaling;
560 switch (pressureType_)
564 velocityW *= lambdaW * scalarPerm
565 * ((cellData.pressure(wPhaseIdx) - cellDataJ.pressure(wPhaseIdx)) / dist + gravityTermW);
566 velocityNw *= lambdaNw * scalarPerm
567 * ((cellData.pressure(wPhaseIdx) - cellDataJ.pressure(wPhaseIdx)) / dist + gravityTermNw)
568 + 0.5 * (lambdaNwI + lambdaNwJ) * scalarPerm * (pcI - pcJ) / dist;
573 velocityW *= lambdaW * scalarPerm
574 * ((cellData.pressure(nPhaseIdx) - cellDataJ.pressure(nPhaseIdx)) / dist + gravityTermW)
575 - 0.5 * (lambdaWI + lambdaWJ) * scalarPerm * (pcI - pcJ) / dist;
576 velocityNw *= lambdaNw * scalarPerm
577 * ((cellData.pressure(nPhaseIdx) - cellDataJ.pressure(nPhaseIdx)) / dist + gravityTermNw);
582 velocityW *= (lambdaW + lambdaNw) * scalarPerm * (cellData.globalPressure() - cellDataJ.globalPressure()) / dist
583 + scalarPerm * (lambdaW * gravityTermW + lambdaNw * gravityTermNw);
590 cellData.fluxData().setVelocity(wPhaseIdx, isIndexI, velocityW);
591 cellData.fluxData().setVelocity(nPhaseIdx, isIndexI, velocityNw);
592 cellData.fluxData().setVelocityMarker(isIndexI);
595 Scalar weightingFactor = pow(0.5, (dim - 1)*(elementI.level() - elementJ.level()));
597 velocityW *= weightingFactor;
598 velocityNw *= weightingFactor;
600 cellDataJ.fluxData().addVelocity(wPhaseIdx, isIndexJ, velocityW);
601 cellDataJ.fluxData().addVelocity(nPhaseIdx, isIndexJ, velocityNw);
602 cellDataJ.fluxData().setVelocityMarker(isIndexJ);
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
Determines the velocity from a finite volume solution of the pressure equation of a sequential model ...
Definition: 2p/sequential/diffusion/cellcentered/velocity.hh:63
void initialize()
For initialization.
Definition: 2p/sequential/diffusion/cellcentered/velocity.hh:142
Determines the velocity from a finite volume solution of the pressure equation of a sequential model ...
Definition: velocityadaptive.hh:40
void initialize()
For initialization.
Definition: velocityadaptive.hh:106
bool calculateVelocityInTransport()
Indicates if velocity is reconstructed in the pressure step or in the transport step.
Definition: velocityadaptive.hh:135
FVVelocity2PAdaptive(Problem &problem)
Constructs a FVVelocity2PAdaptive object.
Definition: velocityadaptive.hh:87
void calculateVelocity(const Intersection &intersection, CellData &cellData)
Calculates the velocity at a cell-cell interface.
Definition: velocityadaptive.hh:163
Velocity field from a finite volume solution of a pressure equation.