24#ifndef DUMUX_FVVELOCITY2P_HH
25#define DUMUX_FVVELOCITY2P_HH
27#include <dune/common/float_cmp.hh>
28#include <dune/grid/common/gridenums.hh>
61template<
class TypeTag>
77 using PrimaryVariables =
typename SolutionTypes::PrimaryVariables;
80 using Element =
typename GridView::Traits::template Codim<0>::Entity;
81 using Intersection =
typename GridView::Intersection;
83 using Geometry =
typename Element::Geometry;
84 using JacobianTransposed =
typename Geometry::JacobianTransposed;
88 dim = GridView::dimension, dimWorld = GridView::dimensionworld
93 pw = Indices::pressureW,
94 pn = Indices::pressureNw,
95 pGlobal = Indices::pressureGlobal,
96 vw = Indices::velocityW,
97 vn = Indices::velocityNw,
98 vt = Indices::velocityTotal,
99 sw = Indices::saturationW,
100 sn = Indices::saturationNw,
101 pressureIdx = Indices::pressureIdx,
102 saturationIdx = Indices::saturationIdx,
103 eqIdxPress = Indices::pressureEqIdx,
104 eqIdxSat = Indices::satEqIdx
108 wPhaseIdx = Indices::wPhaseIdx, nPhaseIdx = Indices::nPhaseIdx, numPhases = getPropValue<TypeTag, Properties::NumPhases>()
111 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
112 using GravityVector = Dune::FieldVector<Scalar, dimWorld>;
113 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
121 problem_(problem), gravity_(problem.gravity())
123 if (getPropValue<TypeTag, Properties::EnableCompressibility>() && velocityType_ == vt)
125 DUNE_THROW(Dune::NotImplemented,
126 "Total velocity - global pressure - model cannot be used with compressible fluids!");
128 if (velocityType_ != vw && velocityType_ != vn && velocityType_ != vt)
130 DUNE_THROW(Dune::NotImplemented,
"Velocity type not supported!");
133 density_[wPhaseIdx] = 0.;
134 density_[nPhaseIdx] = 0.;
135 viscosity_[wPhaseIdx] = 0.;
136 viscosity_[nPhaseIdx] = 0.;
138 vtkOutputLevel_ = getParam<int>(
"Vtk.OutputLevel");
144 if (!compressibility_)
146 const auto element = *problem_.gridView().template begin<0> ();
147 FluidState fluidState;
148 fluidState.setPressure(wPhaseIdx, problem_.referencePressure(element));
149 fluidState.setPressure(nPhaseIdx, problem_.referencePressure(element));
150 fluidState.setTemperature(problem_.temperature(element));
151 fluidState.setSaturation(wPhaseIdx, 1.);
152 fluidState.setSaturation(nPhaseIdx, 0.);
184 template<
class MultiWriter>
187 if (vtkOutputLevel_ > 0)
189 Dune::BlockVector < Dune::FieldVector<Scalar, dim> > &velocity = *(writer.template allocateManagedBuffer<Scalar,
190 dim>(problem_.gridView().size(0)));
191 Dune::BlockVector < Dune::FieldVector<Scalar, dim> > &velocitySecondPhase =
192 *(writer.template allocateManagedBuffer<Scalar, dim>(problem_.gridView().size(0)));
195 for (
const auto& element : elements(problem_.gridView()))
198 int eIdxGlobal = problem_.variables().index(element);
200 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
202 const typename Element::Geometry& geometry = element.geometry();
205 const auto refElement = referenceElement(geometry);
207 const int numberOfFaces=refElement.size(1);
208 std::vector<Scalar> fluxW(numberOfFaces,0);
209 std::vector<Scalar> fluxNw(numberOfFaces,0);
212 for (
const auto& intersection : intersections(problem_.gridView(), element))
214 int isIndex = intersection.indexInInside();
216 fluxW[isIndex] += intersection.geometry().volume()
217 * (intersection.centerUnitOuterNormal() * cellData.fluxData().velocity(wPhaseIdx, isIndex));
218 fluxNw[isIndex] += intersection.geometry().volume()
219 * (intersection.centerUnitOuterNormal() * cellData.fluxData().velocity(nPhaseIdx, isIndex));
224 Dune::FieldVector<Scalar, dim> refVelocity;
226 if (refElement.type().isSimplex()) {
227 for (
int dimIdx = 0; dimIdx < dim; dimIdx++)
229 refVelocity[dimIdx] = -fluxW[dim - 1 - dimIdx];
230 for (
int fIdx = 0; fIdx < dim + 1; fIdx++)
232 refVelocity[dimIdx] += fluxW[fIdx]/(dim + 1);
237 else if (refElement.type().isCube()){
238 for (
int i = 0; i < dim; i++)
239 refVelocity[i] = 0.5 * (fluxW[2*i + 1] - fluxW[2*i]);
243 DUNE_THROW(Dune::NotImplemented,
"velocity output for prism/pyramid not implemented");
246 const Dune::FieldVector<Scalar, dim> localPos =
247 refElement.position(0, 0);
250 const JacobianTransposed jacobianT =
251 geometry.jacobianTransposed(localPos);
254 Dune::FieldVector < Scalar, dim > elementVelocity(0);
255 jacobianT.umtv(refVelocity, elementVelocity);
256 elementVelocity /= geometry.integrationElement(localPos);
258 velocity[eIdxGlobal] = elementVelocity;
263 if (refElement.type().isSimplex()) {
264 for (
int dimIdx = 0; dimIdx < dim; dimIdx++)
266 refVelocity[dimIdx] = -fluxNw[dim - 1 - dimIdx];
267 for (
int fIdx = 0; fIdx < dim + 1; fIdx++)
269 refVelocity[dimIdx] += fluxNw[fIdx]/(dim + 1);
274 else if (refElement.type().isCube()){
275 for (
int i = 0; i < dim; i++)
276 refVelocity[i] = 0.5 * (fluxNw[2*i + 1] - fluxNw[2*i]);
280 DUNE_THROW(Dune::NotImplemented,
"velocity output for prism/pyramid not implemented");
285 jacobianT.umtv(refVelocity, elementVelocity);
286 elementVelocity /= geometry.integrationElement(localPos);
288 velocitySecondPhase[eIdxGlobal] = elementVelocity;
292 if (velocityType_ == vt)
294 writer.attachCellData(velocity,
"total velocity", dim);
298 writer.attachCellData(velocity,
"wetting-velocity", dim);
299 writer.attachCellData(velocitySecondPhase,
"nonwetting-velocity", dim);
308 const GravityVector& gravity_;
309 Scalar density_[numPhases];
310 Scalar viscosity_[numPhases];
315 static const int velocityType_ = getPropValue<TypeTag, Properties::VelocityFormulation>();
316 static const bool compressibility_ = getPropValue<TypeTag, Properties::EnableCompressibility>();
318 static const int pressureType_ = getPropValue<TypeTag, Properties::PressureFormulation>();
320 static const int saturationType_ = getPropValue<TypeTag, Properties::SaturationFormulation>();
331template<
class TypeTag>
334 auto elementI = intersection.inside();
335 auto elementJ = intersection.outside();
337 int eIdxGlobalJ = problem_.variables().index(elementJ);
339 CellData& cellDataJ = problem_.variables().cellData(eIdxGlobalJ);
342 const GlobalPosition& globalPosI = (elementI).geometry().center();
343 const GlobalPosition& globalPosJ = (elementJ).geometry().center();
346 Scalar lambdaWI = cellData.mobility(wPhaseIdx);
347 Scalar lambdaNwI = cellData.mobility(nPhaseIdx);
348 Scalar lambdaWJ = cellDataJ.mobility(wPhaseIdx);
349 Scalar lambdaNwJ = cellDataJ.mobility(nPhaseIdx);
352 Scalar pcI = cellData.capillaryPressure();
353 Scalar pcJ = cellDataJ.capillaryPressure();
356 int isIndexI = intersection.indexInInside();
357 int isIndexJ = intersection.indexInOutside();
360 const Dune::FieldVector<Scalar, dim>& unitOuterNormal = intersection.centerUnitOuterNormal();
363 GlobalPosition distVec = globalPosJ - globalPosI;
366 Scalar dist = distVec.two_norm();
369 DimMatrix meanPermeability(0);
371 problem_.spatialParams().meanK(meanPermeability, problem_.spatialParams().intrinsicPermeability(elementI),
372 problem_.spatialParams().intrinsicPermeability(elementJ));
378 Scalar potentialDiffW = cellData.potential(wPhaseIdx) - cellDataJ.potential(wPhaseIdx);
379 Scalar potentialDiffNw = cellData.potential(nPhaseIdx) - cellDataJ.potential(nPhaseIdx);
381 if (compressibility_)
383 density_[wPhaseIdx] = (potentialDiffW > 0.) ? cellData.density(wPhaseIdx) : cellDataJ.density(wPhaseIdx);
384 density_[nPhaseIdx] = (potentialDiffNw > 0.) ? cellData.density(nPhaseIdx) : cellDataJ.density(nPhaseIdx);
386 density_[wPhaseIdx] =
387 (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffW, 0.0, 1.0e-30)) ? 0.5 * (cellData.density(wPhaseIdx) + cellDataJ.density(wPhaseIdx)) :
389 density_[nPhaseIdx] =
390 (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffNw, 0.0, 1.0e-30)) ? 0.5 * (cellData.density(nPhaseIdx) + cellDataJ.density(nPhaseIdx)) :
393 potentialDiffW = (cellData.pressure(wPhaseIdx) - cellDataJ.pressure(wPhaseIdx));
394 potentialDiffNw = (cellData.pressure(nPhaseIdx) - cellDataJ.pressure(nPhaseIdx));
396 potentialDiffW += density_[wPhaseIdx] * (distVec * gravity_);
397 potentialDiffNw += density_[nPhaseIdx] * (distVec * gravity_);
401 cellData.fluxData().setUpwindPotential(wPhaseIdx, isIndexI, potentialDiffW);
402 cellData.fluxData().setUpwindPotential(nPhaseIdx, isIndexI, potentialDiffNw);
404 cellDataJ.fluxData().setUpwindPotential(wPhaseIdx, isIndexJ, -potentialDiffW);
405 cellDataJ.fluxData().setUpwindPotential(nPhaseIdx, isIndexJ, -potentialDiffNw);
408 Scalar lambdaW = (potentialDiffW > 0.) ? lambdaWI : lambdaWJ;
409 lambdaW = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffW, 0.0, 1.0e-30)) ? 0.5 * (lambdaWI + lambdaWJ) : lambdaW;
410 Scalar lambdaNw = (potentialDiffNw > 0.) ? lambdaNwI : lambdaNwJ;
411 lambdaNw = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffNw, 0.0, 1.0e-30)) ? 0.5 * (lambdaNwI + lambdaNwJ) : lambdaNw;
413 if (compressibility_)
415 density_[wPhaseIdx] = (potentialDiffW > 0.) ? cellData.density(wPhaseIdx) : cellDataJ.density(wPhaseIdx);
416 density_[nPhaseIdx] = (potentialDiffNw > 0.) ? cellData.density(nPhaseIdx) : cellDataJ.density(nPhaseIdx);
418 density_[wPhaseIdx] =
419 (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffW, 0.0, 1.0e-30)) ? 0.5 * (cellData.density(wPhaseIdx) + cellDataJ.density(wPhaseIdx)) :
421 density_[nPhaseIdx] =
422 (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffNw, 0.0, 1.0e-30)) ? 0.5 * (cellData.density(nPhaseIdx) + cellDataJ.density(nPhaseIdx)) :
429 Dune::FieldVector<Scalar, dimWorld> velocityW(unitOuterNormal);
430 Dune::FieldVector<Scalar, dimWorld> velocityNw(unitOuterNormal);
434 Scalar areaScaling = (unitOuterNormal * distVec);
437 Scalar gravityTermW = (gravity_ * distVec) * density_[wPhaseIdx] * areaScaling;
438 Scalar gravityTermNw = (gravity_ * distVec) * density_[nPhaseIdx] * areaScaling;
441 switch (pressureType_)
445 velocityW *= lambdaW * scalarPerm
446 * ((cellData.pressure(wPhaseIdx) - cellDataJ.pressure(wPhaseIdx)) / dist + gravityTermW);
447 velocityNw *= lambdaNw * scalarPerm
448 * ((cellData.pressure(wPhaseIdx) - cellDataJ.pressure(wPhaseIdx)) / dist + gravityTermNw)
449 + 0.5 * (lambdaNwI + lambdaNwJ) * scalarPerm * (pcI - pcJ) / dist;
454 velocityW *= lambdaW * scalarPerm
455 * ((cellData.pressure(nPhaseIdx) - cellDataJ.pressure(nPhaseIdx)) / dist + gravityTermW)
456 - 0.5 * (lambdaWI + lambdaWJ) * scalarPerm * (pcI - pcJ) / dist;
457 velocityNw *= lambdaNw * scalarPerm
458 * ((cellData.pressure(nPhaseIdx) - cellDataJ.pressure(nPhaseIdx)) / dist + gravityTermNw);
463 velocityW *= (lambdaW + lambdaNw) * scalarPerm * (cellData.globalPressure() - cellDataJ.globalPressure()) / dist
464 + scalarPerm * (lambdaW * gravityTermW + lambdaNw * gravityTermNw);
471 cellData.fluxData().setVelocity(wPhaseIdx, isIndexI, velocityW);
472 cellData.fluxData().setVelocity(nPhaseIdx, isIndexI, velocityNw);
473 cellData.fluxData().setVelocityMarker(isIndexI);
475 cellDataJ.fluxData().setVelocity(wPhaseIdx, isIndexJ, velocityW);
476 cellDataJ.fluxData().setVelocity(nPhaseIdx, isIndexJ, velocityNw);
477 cellDataJ.fluxData().setVelocityMarker(isIndexJ);
491template<
class TypeTag>
494 auto element = intersection.inside();
497 int isIndex = intersection.indexInInside();
500 const Dune::FieldVector<Scalar, dim>& unitOuterNormal = intersection.centerUnitOuterNormal();
502 BoundaryTypes bcType;
504 problem_.boundaryTypes(bcType, intersection);
505 PrimaryVariables boundValues(0.0);
507 if (bcType.isDirichlet(eqIdxPress))
509 problem_.dirichlet(boundValues, intersection);
512 const GlobalPosition& globalPosI = element.geometry().center();
515 const GlobalPosition& globalPosJ = intersection.geometry().center();
518 Scalar lambdaWI = cellData.mobility(wPhaseIdx);
519 Scalar lambdaNwI = cellData.mobility(nPhaseIdx);
520 Scalar fractionalWI = cellData.fracFlowFunc(wPhaseIdx);
521 Scalar fractionalNwI = cellData.fracFlowFunc(nPhaseIdx);
524 Scalar pcI = cellData.capillaryPressure();
527 GlobalPosition distVec = globalPosJ - globalPosI;
530 Scalar dist = distVec.two_norm();
534 DimMatrix meanPermeability(0);
536 problem_.spatialParams().meanK(meanPermeability, problem_.spatialParams().intrinsicPermeability(element));
544 if (bcType.isDirichlet(eqIdxSat))
546 switch (saturationType_)
550 satW = boundValues[saturationIdx];
551 satNw = 1 - boundValues[saturationIdx];
556 satW = 1 - boundValues[saturationIdx];
557 satNw = boundValues[saturationIdx];
564 satW = cellData.saturation(wPhaseIdx);
565 satNw = cellData.saturation(nPhaseIdx);
568 const Scalar pressBound = boundValues[pressureIdx];
573 const auto fluidMatrixInteraction = Deprecated::makePcKrSw(Scalar{}, problem_.spatialParams(), element);
575 const Scalar pcBound = fluidMatrixInteraction.pc(satW);
578 Scalar pressWBound = 0;
579 Scalar pressNwBound = 0;
580 if (pressureType_ == pw)
582 pressWBound = pressBound;
583 pressNwBound = pressBound + pcBound;
585 else if (pressureType_ == pn)
587 pressWBound = pressBound - pcBound;
588 pressNwBound = pressBound;
592 const Scalar
temperature = problem_.temperature(element);
594 Scalar densityWBound = density_[wPhaseIdx];
595 Scalar densityNwBound = density_[nPhaseIdx];
596 Scalar viscosityWBound = viscosity_[wPhaseIdx];
597 Scalar viscosityNwBound = viscosity_[nPhaseIdx];
599 if (compressibility_)
601 FluidState fluidState;
602 fluidState.setSaturation(wPhaseIdx, satW);
603 fluidState.setSaturation(nPhaseIdx, satNw);
605 fluidState.setPressure(wPhaseIdx, pressWBound);
606 fluidState.setPressure(nPhaseIdx, pressNwBound);
614 Scalar lambdaWBound = fluidMatrixInteraction.krw(satW)
616 Scalar lambdaNwBound = fluidMatrixInteraction.krn(satW)
619 Scalar potentialDiffW = cellData.fluxData().upwindPotential(wPhaseIdx, isIndex);
620 Scalar potentialDiffNw = cellData.fluxData().upwindPotential(nPhaseIdx, isIndex);
622 if (compressibility_)
624 density_[wPhaseIdx] = (potentialDiffW > 0.) ? cellData.density(wPhaseIdx) : densityWBound;
625 density_[wPhaseIdx] =
626 (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffW, 0.0, 1.0e-30)) ? 0.5 * (cellData.density(wPhaseIdx) + densityWBound) : density_[wPhaseIdx];
627 density_[nPhaseIdx] = (potentialDiffNw > 0.) ? cellData.density(nPhaseIdx) : densityNwBound;
628 density_[nPhaseIdx] =
629 (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffNw, 0.0, 1.0e-30)) ? 0.5 * (cellData.density(nPhaseIdx) + densityNwBound) : density_[nPhaseIdx];
633 if (pressureType_ == pGlobal)
635 potentialDiffW = (cellData.globalPressure() - pressBound - fractionalNwI * (pcI - pcBound));
636 potentialDiffNw = (cellData.globalPressure() - pressBound + fractionalWI * (pcI - pcBound));
640 potentialDiffW = (cellData.pressure(wPhaseIdx) - pressWBound);
641 potentialDiffNw = (cellData.pressure(nPhaseIdx) - pressNwBound);
644 potentialDiffW += density_[wPhaseIdx] * (distVec * gravity_);
645 potentialDiffNw += density_[nPhaseIdx] * (distVec * gravity_);
648 cellData.fluxData().setUpwindPotential(wPhaseIdx, isIndex, potentialDiffW);
649 cellData.fluxData().setUpwindPotential(nPhaseIdx, isIndex, potentialDiffNw);
652 Scalar lambdaW = (potentialDiffW > 0.) ? lambdaWI : lambdaWBound;
653 lambdaW = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffW, 0.0, 1.0e-30)) ? 0.5 * (lambdaWI + lambdaWBound) : lambdaW;
654 Scalar lambdaNw = (potentialDiffNw > 0.) ? lambdaNwI : lambdaNwBound;
655 lambdaNw = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffNw, 0.0, 1.0e-30)) ? 0.5 * (lambdaNwI + lambdaNwBound) : lambdaNw;
657 if (compressibility_)
659 density_[wPhaseIdx] = (potentialDiffW > 0.) ? cellData.density(wPhaseIdx) : densityWBound;
660 density_[wPhaseIdx] =
661 (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffW, 0.0, 1.0e-30)) ? 0.5 * (cellData.density(wPhaseIdx) + densityWBound) : density_[wPhaseIdx];
662 density_[nPhaseIdx] = (potentialDiffNw > 0.) ? cellData.density(nPhaseIdx) : densityNwBound;
663 density_[nPhaseIdx] =
664 (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffNw, 0.0, 1.0e-30)) ? 0.5 * (cellData.density(nPhaseIdx) + densityNwBound) : density_[nPhaseIdx];
670 Dune::FieldVector<Scalar, dimWorld> velocityW(unitOuterNormal);
671 Dune::FieldVector<Scalar, dimWorld> velocityNw(unitOuterNormal);
675 Scalar areaScaling = (unitOuterNormal * distVec);
678 Scalar gravityTermW = (gravity_ * distVec) * density_[wPhaseIdx] * areaScaling;
679 Scalar gravityTermNw = (gravity_ * distVec) * density_[nPhaseIdx] * areaScaling;
682 switch (pressureType_)
686 velocityW *= lambdaW * scalarPerm * ((cellData.pressure(wPhaseIdx) - pressBound) / dist + gravityTermW);
687 velocityNw *= lambdaNw * scalarPerm * ((cellData.pressure(wPhaseIdx) - pressBound) / dist + gravityTermNw)
688 + 0.5 * (lambdaNwI + lambdaNwBound) * scalarPerm * (pcI - pcBound) / dist;
693 velocityW *= lambdaW * scalarPerm * ((cellData.pressure(nPhaseIdx) - pressBound) / dist + gravityTermW)
694 - 0.5 * (lambdaWI + lambdaWBound) * scalarPerm * (pcI - pcBound) / dist;
695 velocityNw *= lambdaNw * scalarPerm * ((cellData.pressure(nPhaseIdx) - pressBound) / dist + gravityTermNw);
700 velocityW *= (lambdaW + lambdaNw) * scalarPerm * (cellData.globalPressure() - pressBound) / dist
701 + scalarPerm * (lambdaW * gravityTermW + lambdaNw * gravityTermNw);
708 cellData.fluxData().setVelocity(wPhaseIdx, isIndex, velocityW);
709 cellData.fluxData().setVelocity(nPhaseIdx, isIndex, velocityNw);
710 cellData.fluxData().setVelocityMarker(isIndex);
714 else if (bcType.isNeumann(eqIdxPress))
716 problem_.neumann(boundValues, intersection);
718 Dune::FieldVector<Scalar, dimWorld> velocityW(unitOuterNormal);
719 Dune::FieldVector<Scalar, dimWorld> velocityNw(unitOuterNormal);
721 velocityW *= boundValues[wPhaseIdx];
722 velocityNw *= boundValues[nPhaseIdx];
724 if (!compressibility_)
726 velocityW /= density_[wPhaseIdx];
727 velocityNw /= density_[nPhaseIdx];
731 cellData.fluxData().setUpwindPotential(wPhaseIdx, isIndex, boundValues[wPhaseIdx]);
732 cellData.fluxData().setUpwindPotential(nPhaseIdx, isIndex, boundValues[nPhaseIdx]);
734 cellData.fluxData().setVelocity(wPhaseIdx, isIndex, velocityW);
735 cellData.fluxData().setVelocity(nPhaseIdx, isIndex, velocityNw);
736 cellData.fluxData().setVelocityMarker(isIndex);
740 DUNE_THROW(Dune::NotImplemented,
"No valid boundary condition type defined for pressure equation!");
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
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 addOutputVtkFields(MultiWriter &writer)
Adds velocity output to the output file.
Definition: 2p/sequential/diffusion/cellcentered/velocity.hh:185
void calculateVelocity(const Intersection &, CellData &)
Calculates the velocity at a cell-cell interface.
Definition: 2p/sequential/diffusion/cellcentered/velocity.hh:332
bool calculateVelocityInTransport()
Indicates if velocity is reconstructed in the pressure step or in the transport step.
Definition: 2p/sequential/diffusion/cellcentered/velocity.hh:171
FVVelocity2P(Problem &problem)
Constructs a FVVelocity2P object.
Definition: 2p/sequential/diffusion/cellcentered/velocity.hh:120
void initialize()
For initialization.
Definition: 2p/sequential/diffusion/cellcentered/velocity.hh:142
void calculateVelocityOnBoundary(const Intersection &, CellData &)
Calculates the velocity at a boundary.
Definition: 2p/sequential/diffusion/cellcentered/velocity.hh:492
Specifies the properties for immiscible 2p diffusion/pressure models.