24#ifndef DUMUX_FVVELOCITY2P_HH
25#define DUMUX_FVVELOCITY2P_HH
27#include <dune/common/float_cmp.hh>
28#include <dune/grid/common/gridenums.hh>
59template<
class TypeTag>
67 using MaterialLaw =
typename SpatialParams::MaterialLaw;
76 using PrimaryVariables =
typename SolutionTypes::PrimaryVariables;
79 using Element =
typename GridView::Traits::template Codim<0>::Entity;
80 using Intersection =
typename GridView::Intersection;
82 using Geometry =
typename Element::Geometry;
83 using JacobianTransposed =
typename Geometry::JacobianTransposed;
87 dim = GridView::dimension, dimWorld = GridView::dimensionworld
90 using ReferenceElements = Dune::ReferenceElements<Scalar, dim>;
94 pw = Indices::pressureW,
95 pn = Indices::pressureNw,
96 pGlobal = Indices::pressureGlobal,
97 vw = Indices::velocityW,
98 vn = Indices::velocityNw,
99 vt = Indices::velocityTotal,
100 sw = Indices::saturationW,
101 sn = Indices::saturationNw,
102 pressureIdx = Indices::pressureIdx,
103 saturationIdx = Indices::saturationIdx,
104 eqIdxPress = Indices::pressureEqIdx,
105 eqIdxSat = Indices::satEqIdx
109 wPhaseIdx = Indices::wPhaseIdx, nPhaseIdx = Indices::nPhaseIdx, numPhases = getPropValue<TypeTag, Properties::NumPhases>()
112 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
113 using GravityVector = Dune::FieldVector<Scalar, dimWorld>;
114 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
122 problem_(problem), gravity_(problem.gravity())
124 if (getPropValue<TypeTag, Properties::EnableCompressibility>() && velocityType_ == vt)
126 DUNE_THROW(Dune::NotImplemented,
127 "Total velocity - global pressure - model cannot be used with compressible fluids!");
129 if (velocityType_ != vw && velocityType_ != vn && velocityType_ != vt)
131 DUNE_THROW(Dune::NotImplemented,
"Velocity type not supported!");
134 density_[wPhaseIdx] = 0.;
135 density_[nPhaseIdx] = 0.;
136 viscosity_[wPhaseIdx] = 0.;
137 viscosity_[nPhaseIdx] = 0.;
139 vtkOutputLevel_ = getParam<int>(
"Vtk.OutputLevel");
145 if (!compressibility_)
147 const auto element = *problem_.gridView().template begin<0> ();
148 FluidState fluidState;
149 fluidState.setPressure(wPhaseIdx, problem_.referencePressure(element));
150 fluidState.setPressure(nPhaseIdx, problem_.referencePressure(element));
151 fluidState.setTemperature(problem_.temperature(element));
152 fluidState.setSaturation(wPhaseIdx, 1.);
153 fluidState.setSaturation(nPhaseIdx, 0.);
185 template<
class MultiWriter>
188 if (vtkOutputLevel_ > 0)
190 Dune::BlockVector < Dune::FieldVector<Scalar, dim> > &velocity = *(writer.template allocateManagedBuffer<Scalar,
191 dim>(problem_.gridView().size(0)));
192 Dune::BlockVector < Dune::FieldVector<Scalar, dim> > &velocitySecondPhase =
193 *(writer.template allocateManagedBuffer<Scalar, dim>(problem_.gridView().size(0)));
196 for (
const auto& element : elements(problem_.gridView()))
199 int eIdxGlobal = problem_.variables().index(element);
201 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
203 const typename Element::Geometry& geometry = element.geometry();
206 using ReferenceElements = Dune::ReferenceElements<Scalar, dim>;
207 const auto refElement = ReferenceElements::general(geometry.type());
209 const int numberOfFaces=refElement.size(1);
210 std::vector<Scalar> fluxW(numberOfFaces,0);
211 std::vector<Scalar> fluxNw(numberOfFaces,0);
214 for (
const auto& intersection : intersections(problem_.gridView(), element))
216 int isIndex = intersection.indexInInside();
218 fluxW[isIndex] += intersection.geometry().volume()
219 * (intersection.centerUnitOuterNormal() * cellData.fluxData().velocity(wPhaseIdx, isIndex));
220 fluxNw[isIndex] += intersection.geometry().volume()
221 * (intersection.centerUnitOuterNormal() * cellData.fluxData().velocity(nPhaseIdx, isIndex));
226 Dune::FieldVector<Scalar, dim> refVelocity;
228 if (refElement.type().isSimplex()) {
229 for (
int dimIdx = 0; dimIdx < dim; dimIdx++)
231 refVelocity[dimIdx] = -fluxW[dim - 1 - dimIdx];
232 for (
int fIdx = 0; fIdx < dim + 1; fIdx++)
234 refVelocity[dimIdx] += fluxW[fIdx]/(dim + 1);
239 else if (refElement.type().isCube()){
240 for (
int i = 0; i < dim; i++)
241 refVelocity[i] = 0.5 * (fluxW[2*i + 1] - fluxW[2*i]);
245 DUNE_THROW(Dune::NotImplemented,
"velocity output for prism/pyramid not implemented");
248 const Dune::FieldVector<Scalar, dim> localPos =
249 refElement.position(0, 0);
252 const JacobianTransposed jacobianT =
253 geometry.jacobianTransposed(localPos);
256 Dune::FieldVector < Scalar, dim > elementVelocity(0);
257 jacobianT.umtv(refVelocity, elementVelocity);
258 elementVelocity /= geometry.integrationElement(localPos);
260 velocity[eIdxGlobal] = elementVelocity;
265 if (refElement.type().isSimplex()) {
266 for (
int dimIdx = 0; dimIdx < dim; dimIdx++)
268 refVelocity[dimIdx] = -fluxNw[dim - 1 - dimIdx];
269 for (
int fIdx = 0; fIdx < dim + 1; fIdx++)
271 refVelocity[dimIdx] += fluxNw[fIdx]/(dim + 1);
276 else if (refElement.type().isCube()){
277 for (
int i = 0; i < dim; i++)
278 refVelocity[i] = 0.5 * (fluxNw[2*i + 1] - fluxNw[2*i]);
282 DUNE_THROW(Dune::NotImplemented,
"velocity output for prism/pyramid not implemented");
287 jacobianT.umtv(refVelocity, elementVelocity);
288 elementVelocity /= geometry.integrationElement(localPos);
290 velocitySecondPhase[eIdxGlobal] = elementVelocity;
294 if (velocityType_ == vt)
296 writer.attachCellData(velocity,
"total velocity", dim);
300 writer.attachCellData(velocity,
"wetting-velocity", dim);
301 writer.attachCellData(velocitySecondPhase,
"non-wetting-velocity", dim);
310 const GravityVector& gravity_;
311 Scalar density_[numPhases];
312 Scalar viscosity_[numPhases];
317 static const int velocityType_ = getPropValue<TypeTag, Properties::VelocityFormulation>();
318 static const bool compressibility_ = getPropValue<TypeTag, Properties::EnableCompressibility>();
320 static const int pressureType_ = getPropValue<TypeTag, Properties::PressureFormulation>();
322 static const int saturationType_ = getPropValue<TypeTag, Properties::SaturationFormulation>();
333template<
class TypeTag>
336 auto elementI = intersection.inside();
337 auto elementJ = intersection.outside();
339 int eIdxGlobalJ = problem_.variables().index(elementJ);
341 CellData& cellDataJ = problem_.variables().cellData(eIdxGlobalJ);
344 const GlobalPosition& globalPosI = (elementI).geometry().center();
345 const GlobalPosition& globalPosJ = (elementJ).geometry().center();
348 Scalar lambdaWI = cellData.mobility(wPhaseIdx);
349 Scalar lambdaNwI = cellData.mobility(nPhaseIdx);
350 Scalar lambdaWJ = cellDataJ.mobility(wPhaseIdx);
351 Scalar lambdaNwJ = cellDataJ.mobility(nPhaseIdx);
354 Scalar pcI = cellData.capillaryPressure();
355 Scalar pcJ = cellDataJ.capillaryPressure();
358 int isIndexI = intersection.indexInInside();
359 int isIndexJ = intersection.indexInOutside();
362 const Dune::FieldVector<Scalar, dim>& unitOuterNormal = intersection.centerUnitOuterNormal();
365 GlobalPosition distVec = globalPosJ - globalPosI;
368 Scalar dist = distVec.two_norm();
371 DimMatrix meanPermeability(0);
373 problem_.spatialParams().meanK(meanPermeability, problem_.spatialParams().intrinsicPermeability(elementI),
374 problem_.spatialParams().intrinsicPermeability(elementJ));
380 Scalar potentialDiffW = cellData.potential(wPhaseIdx) - cellDataJ.potential(wPhaseIdx);
381 Scalar potentialDiffNw = cellData.potential(nPhaseIdx) - cellDataJ.potential(nPhaseIdx);
383 if (compressibility_)
385 density_[wPhaseIdx] = (potentialDiffW > 0.) ? cellData.density(wPhaseIdx) : cellDataJ.density(wPhaseIdx);
386 density_[nPhaseIdx] = (potentialDiffNw > 0.) ? cellData.density(nPhaseIdx) : cellDataJ.density(nPhaseIdx);
388 density_[wPhaseIdx] =
389 (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffW, 0.0, 1.0e-30)) ? 0.5 * (cellData.density(wPhaseIdx) + cellDataJ.density(wPhaseIdx)) :
391 density_[nPhaseIdx] =
392 (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffNw, 0.0, 1.0e-30)) ? 0.5 * (cellData.density(nPhaseIdx) + cellDataJ.density(nPhaseIdx)) :
395 potentialDiffW = (cellData.pressure(wPhaseIdx) - cellDataJ.pressure(wPhaseIdx));
396 potentialDiffNw = (cellData.pressure(nPhaseIdx) - cellDataJ.pressure(nPhaseIdx));
398 potentialDiffW += density_[wPhaseIdx] * (distVec * gravity_);
399 potentialDiffNw += density_[nPhaseIdx] * (distVec * gravity_);
403 cellData.fluxData().setUpwindPotential(wPhaseIdx, isIndexI, potentialDiffW);
404 cellData.fluxData().setUpwindPotential(nPhaseIdx, isIndexI, potentialDiffNw);
406 cellDataJ.fluxData().setUpwindPotential(wPhaseIdx, isIndexJ, -potentialDiffW);
407 cellDataJ.fluxData().setUpwindPotential(nPhaseIdx, isIndexJ, -potentialDiffNw);
410 Scalar lambdaW = (potentialDiffW > 0.) ? lambdaWI : lambdaWJ;
411 lambdaW = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffW, 0.0, 1.0e-30)) ? 0.5 * (lambdaWI + lambdaWJ) : lambdaW;
412 Scalar lambdaNw = (potentialDiffNw > 0.) ? lambdaNwI : lambdaNwJ;
413 lambdaNw = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffNw, 0.0, 1.0e-30)) ? 0.5 * (lambdaNwI + lambdaNwJ) : lambdaNw;
415 if (compressibility_)
417 density_[wPhaseIdx] = (potentialDiffW > 0.) ? cellData.density(wPhaseIdx) : cellDataJ.density(wPhaseIdx);
418 density_[nPhaseIdx] = (potentialDiffNw > 0.) ? cellData.density(nPhaseIdx) : cellDataJ.density(nPhaseIdx);
420 density_[wPhaseIdx] =
421 (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffW, 0.0, 1.0e-30)) ? 0.5 * (cellData.density(wPhaseIdx) + cellDataJ.density(wPhaseIdx)) :
423 density_[nPhaseIdx] =
424 (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffNw, 0.0, 1.0e-30)) ? 0.5 * (cellData.density(nPhaseIdx) + cellDataJ.density(nPhaseIdx)) :
431 Dune::FieldVector<Scalar, dimWorld> velocityW(unitOuterNormal);
432 Dune::FieldVector<Scalar, dimWorld> velocityNw(unitOuterNormal);
436 Scalar areaScaling = (unitOuterNormal * distVec);
439 Scalar gravityTermW = (gravity_ * distVec) * density_[wPhaseIdx] * areaScaling;
440 Scalar gravityTermNw = (gravity_ * distVec) * density_[nPhaseIdx] * areaScaling;
443 switch (pressureType_)
447 velocityW *= lambdaW * scalarPerm
448 * ((cellData.pressure(wPhaseIdx) - cellDataJ.pressure(wPhaseIdx)) / dist + gravityTermW);
449 velocityNw *= lambdaNw * scalarPerm
450 * ((cellData.pressure(wPhaseIdx) - cellDataJ.pressure(wPhaseIdx)) / dist + gravityTermNw)
451 + 0.5 * (lambdaNwI + lambdaNwJ) * scalarPerm * (pcI - pcJ) / dist;
456 velocityW *= lambdaW * scalarPerm
457 * ((cellData.pressure(nPhaseIdx) - cellDataJ.pressure(nPhaseIdx)) / dist + gravityTermW)
458 - 0.5 * (lambdaWI + lambdaWJ) * scalarPerm * (pcI - pcJ) / dist;
459 velocityNw *= lambdaNw * scalarPerm
460 * ((cellData.pressure(nPhaseIdx) - cellDataJ.pressure(nPhaseIdx)) / dist + gravityTermNw);
465 velocityW *= (lambdaW + lambdaNw) * scalarPerm * (cellData.globalPressure() - cellDataJ.globalPressure()) / dist
466 + scalarPerm * (lambdaW * gravityTermW + lambdaNw * gravityTermNw);
473 cellData.fluxData().setVelocity(wPhaseIdx, isIndexI, velocityW);
474 cellData.fluxData().setVelocity(nPhaseIdx, isIndexI, velocityNw);
475 cellData.fluxData().setVelocityMarker(isIndexI);
477 cellDataJ.fluxData().setVelocity(wPhaseIdx, isIndexJ, velocityW);
478 cellDataJ.fluxData().setVelocity(nPhaseIdx, isIndexJ, velocityNw);
479 cellDataJ.fluxData().setVelocityMarker(isIndexJ);
493template<
class TypeTag>
496 auto element = intersection.inside();
499 int isIndex = intersection.indexInInside();
502 const Dune::FieldVector<Scalar, dim>& unitOuterNormal = intersection.centerUnitOuterNormal();
504 BoundaryTypes bcType;
506 problem_.boundaryTypes(bcType, intersection);
507 PrimaryVariables boundValues(0.0);
509 if (bcType.isDirichlet(eqIdxPress))
511 problem_.dirichlet(boundValues, intersection);
514 const GlobalPosition& globalPosI = element.geometry().center();
517 const GlobalPosition& globalPosJ = intersection.geometry().center();
520 Scalar lambdaWI = cellData.mobility(wPhaseIdx);
521 Scalar lambdaNwI = cellData.mobility(nPhaseIdx);
522 Scalar fractionalWI = cellData.fracFlowFunc(wPhaseIdx);
523 Scalar fractionalNwI = cellData.fracFlowFunc(nPhaseIdx);
526 Scalar pcI = cellData.capillaryPressure();
529 GlobalPosition distVec = globalPosJ - globalPosI;
532 Scalar dist = distVec.two_norm();
536 DimMatrix meanPermeability(0);
538 problem_.spatialParams().meanK(meanPermeability, problem_.spatialParams().intrinsicPermeability(element));
546 if (bcType.isDirichlet(eqIdxSat))
548 switch (saturationType_)
552 satW = boundValues[saturationIdx];
553 satNw = 1 - boundValues[saturationIdx];
558 satW = 1 - boundValues[saturationIdx];
559 satNw = boundValues[saturationIdx];
566 satW = cellData.saturation(wPhaseIdx);
567 satNw = cellData.saturation(nPhaseIdx);
570 Scalar pressBound = boundValues[pressureIdx];
571 Scalar pcBound = MaterialLaw::pc(problem_.spatialParams().materialLawParams(element), satW);
574 Scalar pressWBound = 0;
575 Scalar pressNwBound = 0;
576 if (pressureType_ == pw)
578 pressWBound = pressBound;
579 pressNwBound = pressBound + pcBound;
581 else if (pressureType_ == pn)
583 pressWBound = pressBound - pcBound;
584 pressNwBound = pressBound;
588 Scalar
temperature = problem_.temperature(element);
590 Scalar densityWBound = density_[wPhaseIdx];
591 Scalar densityNwBound = density_[nPhaseIdx];
592 Scalar viscosityWBound = viscosity_[wPhaseIdx];
593 Scalar viscosityNwBound = viscosity_[nPhaseIdx];
595 if (compressibility_)
597 FluidState fluidState;
598 fluidState.setSaturation(wPhaseIdx, satW);
599 fluidState.setSaturation(nPhaseIdx, satNw);
601 fluidState.setPressure(wPhaseIdx, pressWBound);
602 fluidState.setPressure(nPhaseIdx, pressNwBound);
610 Scalar lambdaWBound = MaterialLaw::krw(problem_.spatialParams().materialLawParams(element), satW)
612 Scalar lambdaNwBound = MaterialLaw::krn(problem_.spatialParams().materialLawParams(element), satW)
615 Scalar potentialDiffW = cellData.fluxData().upwindPotential(wPhaseIdx, isIndex);
616 Scalar potentialDiffNw = cellData.fluxData().upwindPotential(nPhaseIdx, isIndex);
618 if (compressibility_)
620 density_[wPhaseIdx] = (potentialDiffW > 0.) ? cellData.density(wPhaseIdx) : densityWBound;
621 density_[wPhaseIdx] =
622 (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffW, 0.0, 1.0e-30)) ? 0.5 * (cellData.density(wPhaseIdx) + densityWBound) : density_[wPhaseIdx];
623 density_[nPhaseIdx] = (potentialDiffNw > 0.) ? cellData.density(nPhaseIdx) : densityNwBound;
624 density_[nPhaseIdx] =
625 (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffNw, 0.0, 1.0e-30)) ? 0.5 * (cellData.density(nPhaseIdx) + densityNwBound) : density_[nPhaseIdx];
629 if (pressureType_ == pGlobal)
631 potentialDiffW = (cellData.globalPressure() - pressBound - fractionalNwI * (pcI - pcBound));
632 potentialDiffNw = (cellData.globalPressure() - pressBound + fractionalWI * (pcI - pcBound));
636 potentialDiffW = (cellData.pressure(wPhaseIdx) - pressWBound);
637 potentialDiffNw = (cellData.pressure(nPhaseIdx) - pressNwBound);
640 potentialDiffW += density_[wPhaseIdx] * (distVec * gravity_);
641 potentialDiffNw += density_[nPhaseIdx] * (distVec * gravity_);
644 cellData.fluxData().setUpwindPotential(wPhaseIdx, isIndex, potentialDiffW);
645 cellData.fluxData().setUpwindPotential(nPhaseIdx, isIndex, potentialDiffNw);
648 Scalar lambdaW = (potentialDiffW > 0.) ? lambdaWI : lambdaWBound;
649 lambdaW = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffW, 0.0, 1.0e-30)) ? 0.5 * (lambdaWI + lambdaWBound) : lambdaW;
650 Scalar lambdaNw = (potentialDiffNw > 0.) ? lambdaNwI : lambdaNwBound;
651 lambdaNw = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffNw, 0.0, 1.0e-30)) ? 0.5 * (lambdaNwI + lambdaNwBound) : lambdaNw;
653 if (compressibility_)
655 density_[wPhaseIdx] = (potentialDiffW > 0.) ? cellData.density(wPhaseIdx) : densityWBound;
656 density_[wPhaseIdx] =
657 (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffW, 0.0, 1.0e-30)) ? 0.5 * (cellData.density(wPhaseIdx) + densityWBound) : density_[wPhaseIdx];
658 density_[nPhaseIdx] = (potentialDiffNw > 0.) ? cellData.density(nPhaseIdx) : densityNwBound;
659 density_[nPhaseIdx] =
660 (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffNw, 0.0, 1.0e-30)) ? 0.5 * (cellData.density(nPhaseIdx) + densityNwBound) : density_[nPhaseIdx];
666 Dune::FieldVector<Scalar, dimWorld> velocityW(unitOuterNormal);
667 Dune::FieldVector<Scalar, dimWorld> velocityNw(unitOuterNormal);
671 Scalar areaScaling = (unitOuterNormal * distVec);
674 Scalar gravityTermW = (gravity_ * distVec) * density_[wPhaseIdx] * areaScaling;
675 Scalar gravityTermNw = (gravity_ * distVec) * density_[nPhaseIdx] * areaScaling;
678 switch (pressureType_)
682 velocityW *= lambdaW * scalarPerm * ((cellData.pressure(wPhaseIdx) - pressBound) / dist + gravityTermW);
683 velocityNw *= lambdaNw * scalarPerm * ((cellData.pressure(wPhaseIdx) - pressBound) / dist + gravityTermNw)
684 + 0.5 * (lambdaNwI + lambdaNwBound) * scalarPerm * (pcI - pcBound) / dist;
689 velocityW *= lambdaW * scalarPerm * ((cellData.pressure(nPhaseIdx) - pressBound) / dist + gravityTermW)
690 - 0.5 * (lambdaWI + lambdaWBound) * scalarPerm * (pcI - pcBound) / dist;
691 velocityNw *= lambdaNw * scalarPerm * ((cellData.pressure(nPhaseIdx) - pressBound) / dist + gravityTermNw);
696 velocityW *= (lambdaW + lambdaNw) * scalarPerm * (cellData.globalPressure() - pressBound) / dist
697 + scalarPerm * (lambdaW * gravityTermW + lambdaNw * gravityTermNw);
704 cellData.fluxData().setVelocity(wPhaseIdx, isIndex, velocityW);
705 cellData.fluxData().setVelocity(nPhaseIdx, isIndex, velocityNw);
706 cellData.fluxData().setVelocityMarker(isIndex);
710 else if (bcType.isNeumann(eqIdxPress))
712 problem_.neumann(boundValues, intersection);
714 Dune::FieldVector<Scalar, dimWorld> velocityW(unitOuterNormal);
715 Dune::FieldVector<Scalar, dimWorld> velocityNw(unitOuterNormal);
717 velocityW *= boundValues[wPhaseIdx];
718 velocityNw *= boundValues[nPhaseIdx];
720 if (!compressibility_)
722 velocityW /= density_[wPhaseIdx];
723 velocityNw /= density_[nPhaseIdx];
727 cellData.fluxData().setUpwindPotential(wPhaseIdx, isIndex, boundValues[wPhaseIdx]);
728 cellData.fluxData().setUpwindPotential(nPhaseIdx, isIndex, boundValues[nPhaseIdx]);
730 cellData.fluxData().setVelocity(wPhaseIdx, isIndex, velocityW);
731 cellData.fluxData().setVelocity(nPhaseIdx, isIndex, velocityNw);
732 cellData.fluxData().setVelocityMarker(isIndex);
736 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:61
void addOutputVtkFields(MultiWriter &writer)
Adds velocity output to the output file.
Definition: 2p/sequential/diffusion/cellcentered/velocity.hh:186
void calculateVelocity(const Intersection &, CellData &)
Calculates the velocity at a cell-cell interface.
Definition: 2p/sequential/diffusion/cellcentered/velocity.hh:334
bool calculateVelocityInTransport()
Indicates if velocity is reconstructed in the pressure step or in the transport step.
Definition: 2p/sequential/diffusion/cellcentered/velocity.hh:172
FVVelocity2P(Problem &problem)
Constructs a FVVelocity2P object.
Definition: 2p/sequential/diffusion/cellcentered/velocity.hh:121
void initialize()
For initialization.
Definition: 2p/sequential/diffusion/cellcentered/velocity.hh:143
void calculateVelocityOnBoundary(const Intersection &, CellData &)
Calculates the velocity at a boundary.
Definition: 2p/sequential/diffusion/cellcentered/velocity.hh:494
Specifies the properties for immiscible 2p diffusion/pressure models.