24#ifndef DUMUX_FVSATURATION2P_HH
25#define DUMUX_FVSATURATION2P_HH
76template<
class TypeTag>
86 dim = GridView::dimension, dimWorld = GridView::dimensionworld
106 using PrimaryVariables =
typename SolutionTypes::PrimaryVariables;
112 pw = Indices::pressureW,
113 pn = Indices::pressureNw,
114 pGlobal = Indices::pressureGlobal,
115 vw = Indices::velocityW,
116 vn = Indices::velocityNw,
117 vt = Indices::velocityTotal,
118 sw = Indices::saturationW,
119 sn = Indices::saturationNw
123 wPhaseIdx = Indices::wPhaseIdx,
124 nPhaseIdx = Indices::nPhaseIdx,
125 saturationIdx = Indices::saturationIdx,
126 satEqIdx = Indices::satEqIdx,
127 numPhases = getPropValue<TypeTag, Properties::NumPhases>()
132 using Element =
typename GridView::Traits::template Codim<0>::Entity;
133 using Intersection =
typename GridView::Intersection;
135 using GlobalPosition = Dune::FieldVector<Scalar, dimWorld>;
136 using DimVector = Dune::FieldVector<Scalar, dim>;
141 return *capillaryFlux_;
146 return *capillaryFlux_;
151 return *gravityFlux_;
156 return *gravityFlux_;
171 void getFlux(Scalar&
update,
const Intersection& intersection, CellData& cellDataI);
177 void getSource(Scalar&
update,
const Element& element, CellData& cellDataI);
194 int size = problem_.gridView().size(0);
195 transportedQuantity.resize(size);
196 for (
int i = 0; i < size; i++)
198 switch (saturationType_)
202 transportedQuantity[i] = problem_.variables().cellData(i).saturation(wPhaseIdx);
207 transportedQuantity[i] = problem_.variables().cellData(i).saturation(nPhaseIdx);
222 int size = problem_.gridView().size(0);
223 for (
int i = 0; i < size; i++)
225 CellData& cellData = problem_.variables().cellData(i);
226 switch (saturationType_)
230 Scalar sat = transportedQuantity[i];
232 cellData.setSaturation(wPhaseIdx, sat);
233 cellData.setSaturation(nPhaseIdx, 1 - sat);
238 Scalar sat = transportedQuantity[i];
240 cellData.setSaturation(wPhaseIdx,1 -sat);
241 cellData.setSaturation(nPhaseIdx, sat);
254 template<
class DataEntry>
257 return (entry > -1e-6 && entry < 1.0 + 1e-6);
270 asImp_().updateSaturationSolution(updateVec);
281 asImp_().updateSaturationSolution(updateVec, dt);
291 Scalar dt = problem_.timeManager().timeStepSize();
292 int size = problem_.gridView().size(0);
293 for (
int i = 0; i < size; i++)
295 asImp_().updateSaturationSolution(i, updateVec[i][0], dt);
307 int size = problem_.gridView().size(0);
308 for (
int i = 0; i < size; i++)
310 asImp_().updateSaturationSolution(i, updateVec[i][0], dt);
325 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
327 switch (saturationType_)
331 Scalar sat = cellData.saturation(wPhaseIdx) + dt*
update;
333 cellData.setSaturation(wPhaseIdx, sat);
334 cellData.setSaturation(nPhaseIdx, 1 - sat);
339 Scalar sat = cellData.saturation(nPhaseIdx) + dt*
update;
341 cellData.setSaturation(wPhaseIdx,1 -sat);
342 cellData.setSaturation(nPhaseIdx, sat);
360 template<
class MultiWriter>
363 int size = problem_.gridView().size(0);
364 TransportSolutionType *
saturation = writer.allocateManagedBuffer(size);
365 TransportSolutionType *saturationSecond = 0;
367 if (vtkOutputLevel_ > 0)
369 saturationSecond = writer.allocateManagedBuffer(size);
372 for (
int i = 0; i < size; i++)
374 CellData& cellData = problem_.variables().cellData(i);
376 if (saturationType_ == sw)
378 (*saturation)[i] = cellData.saturation(wPhaseIdx);
379 if (vtkOutputLevel_ > 0)
381 (*saturationSecond)[i] = cellData.saturation(nPhaseIdx);
384 else if (saturationType_ == sn)
386 (*saturation)[i] = cellData.saturation(nPhaseIdx);
387 if (vtkOutputLevel_ > 0)
389 (*saturationSecond)[i] = cellData.saturation(wPhaseIdx);
394 if (saturationType_ == sw)
396 writer.attachCellData(*
saturation,
"wetting saturation");
397 if (vtkOutputLevel_ > 0)
399 writer.attachCellData(*saturationSecond,
"nonwetting saturation");
402 else if (saturationType_ == sn)
404 writer.attachCellData(*
saturation,
"nonwetting saturation");
405 if (vtkOutputLevel_ > 0)
407 writer.attachCellData(*saturationSecond,
"wetting saturation");
411 if (
velocity().calculateVelocityInTransport())
412 velocity().addOutputVtkFields(writer);
424 int eIdxGlobal = problem_.variables().index(element);
427 switch (saturationType_)
430 sat = problem_.variables().cellData(eIdxGlobal).saturation(wPhaseIdx);
433 sat = problem_.variables().cellData(eIdxGlobal).saturation(nPhaseIdx);
447 int eIdxGlobal = problem_.variables().index(element);
452 switch (saturationType_)
455 problem_.variables().cellData(eIdxGlobal).setSaturation(wPhaseIdx, sat);
456 problem_.variables().cellData(eIdxGlobal).setSaturation(nPhaseIdx, 1-sat);
459 problem_.variables().cellData(eIdxGlobal).setSaturation(nPhaseIdx, sat);
460 problem_.variables().cellData(eIdxGlobal).setSaturation(wPhaseIdx, 1-sat);
472 ParentType(problem), problem_(problem), threshold_(1e-6),
473 switchNormals_(
getParam<bool>(
"Impet.SwitchNormals"))
475 if (compressibility_ && velocityType_ == vt)
477 DUNE_THROW(Dune::NotImplemented,
478 "Total velocity - global pressure - model cannot be used with compressible fluids!");
480 if (saturationType_ != sw && saturationType_ != sn)
482 DUNE_THROW(Dune::NotImplemented,
"Saturation type not supported!");
484 if (pressureType_ != pw && pressureType_ != pn && pressureType_ != pGlobal)
486 DUNE_THROW(Dune::NotImplemented,
"Pressure type not supported!");
488 if (velocityType_ != vw && velocityType_ != vn && velocityType_ != vt)
490 DUNE_THROW(Dune::NotImplemented,
"Velocity type not supported!");
493 capillaryFlux_ = std::make_shared<CapillaryFlux>(problem);
494 gravityFlux_ = std::make_shared<GravityFlux>(problem);
495 velocity_ = std::make_shared<Velocity>(problem);
497 vtkOutputLevel_ = getParam<int>(
"Vtk.OutputLevel");
498 porosityThreshold_ = getParam<Scalar>(
"Impet.PorosityThreshold");
503 Implementation &asImp_()
504 {
return *
static_cast<Implementation *
>(
this); }
507 const Implementation &asImp_()
const
508 {
return *
static_cast<const Implementation *
>(
this); }
511 std::shared_ptr<Velocity> velocity_;
512 std::shared_ptr<CapillaryFlux> capillaryFlux_;
513 std::shared_ptr<GravityFlux> gravityFlux_;
516 Scalar porosityThreshold_;
518 static const bool compressibility_ = getPropValue<TypeTag, Properties::EnableCompressibility>();
519 static const int saturationType_ = getPropValue<TypeTag, Properties::SaturationFormulation>();
520 static const int velocityType_ = getPropValue<TypeTag, Properties::VelocityFormulation>();
521 static const int pressureType_ = getPropValue<TypeTag, Properties::PressureFormulation>();
523 Scalar density_[numPhases];
524 Scalar viscosity_[numPhases];
526 const Scalar threshold_;
527 const bool switchNormals_;
541template<
class TypeTag>
544 auto elementI = intersection.inside();
545 auto elementJ = intersection.outside();
547 const CellData& cellDataJ = problem_.variables().cellData(problem_.variables().index(elementJ));
550 const GlobalPosition& globalPosI = elementI.geometry().center();
551 const GlobalPosition& globalPosJ = elementJ.geometry().center();
554 Scalar volume = elementI.geometry().volume();
556 Scalar
porosity = max(problem_.spatialParams().porosity(elementI), porosityThreshold_);
558 if (compressibility_)
560 viscosity_[wPhaseIdx] = cellDataI.viscosity(wPhaseIdx);
561 viscosity_[nPhaseIdx] = cellDataI.viscosity(nPhaseIdx);
565 int isIndex = intersection.indexInInside();
567 GlobalPosition unitOuterNormal = intersection.centerUnitOuterNormal();
569 unitOuterNormal *= -1.0;
571 Scalar faceArea = intersection.geometry().volume();
573 if (velocity().calculateVelocityInTransport() && !cellDataI.fluxData().haveVelocity(isIndex))
574 velocity().calculateVelocity(intersection, cellDataI);
577 Scalar factorW = (cellDataI.fluxData().velocity(wPhaseIdx, isIndex) * unitOuterNormal) * faceArea;
578 Scalar factorNw = (cellDataI.fluxData().velocity(nPhaseIdx, isIndex) * unitOuterNormal) * faceArea;
579 Scalar factorTotal = factorW + factorNw;
582 GlobalPosition distVec = globalPosJ - globalPosI;
584 Scalar dist = distVec.two_norm();
586 bool takeNeighbor = (elementI.level() < elementJ.level());
589 (takeNeighbor) ? !cellDataJ.fluxData().isUpwindCell(wPhaseIdx, intersection.indexInOutside()) :
590 cellDataI.fluxData().isUpwindCell(wPhaseIdx, isIndex);
592 (takeNeighbor) ? !cellDataJ.fluxData().isUpwindCell(nPhaseIdx, intersection.indexInOutside()) :
593 cellDataI.fluxData().isUpwindCell(nPhaseIdx, isIndex);
601 lambdaW = cellDataI.mobility(wPhaseIdx);
602 if (compressibility_)
604 lambdaW /= cellDataI.density(wPhaseIdx);
609 lambdaW = cellDataJ.mobility(wPhaseIdx);
610 if (compressibility_)
612 lambdaW /= cellDataJ.density(wPhaseIdx);
618 lambdaNw = cellDataI.mobility(nPhaseIdx);
619 if (compressibility_)
621 lambdaNw /= cellDataI.density(nPhaseIdx);
626 lambdaNw = cellDataJ.mobility(nPhaseIdx);
627 if (compressibility_)
629 lambdaNw /= cellDataJ.density(nPhaseIdx);
633 switch (velocityType_)
638 this->evalCflFluxFunction().addFlux(lambdaW, lambdaNw, viscosity_[wPhaseIdx],
639 viscosity_[nPhaseIdx], factorTotal, intersection);
642 Scalar satWI = cellDataI.saturation(wPhaseIdx);
643 Scalar satWJ = cellDataJ.saturation(wPhaseIdx);
645 Scalar pcI = cellDataI.capillaryPressure();
646 Scalar pcJ = cellDataJ.capillaryPressure();
649 GlobalPosition pcGradient = unitOuterNormal;
650 pcGradient *= (pcI - pcJ) / dist;
654 capillaryFlux().getFlux(flux, intersection, satWI, satWJ, pcGradient);
655 Scalar capillaryFlux = (flux * unitOuterNormal * faceArea);
658 gravityFlux().getFlux(flux, intersection, satWI, satWJ);
659 Scalar gravityFlux = (flux * unitOuterNormal * faceArea);
661 switch (saturationType_)
666 factorTotal *= lambdaW / (lambdaW + lambdaNw);
672 factorTotal *= lambdaNw / (lambdaW + lambdaNw);
678 factorTotal -= capillaryFlux;
679 factorTotal += gravityFlux;
682 this->evalCflFluxFunction().addFlux(lambdaW, lambdaNw, viscosity_[wPhaseIdx],
683 viscosity_[nPhaseIdx], 10 * capillaryFlux, intersection);
684 this->evalCflFluxFunction().addFlux(lambdaW, lambdaNw, viscosity_[wPhaseIdx],
685 viscosity_[nPhaseIdx], 10 * gravityFlux, intersection);
691 if (compressibility_)
693 factorW /= cellDataI.density(wPhaseIdx);
694 factorNw /= cellDataI.density(nPhaseIdx);
698 this->evalCflFluxFunction().addFlux(lambdaW, lambdaNw, viscosity_[wPhaseIdx],
699 viscosity_[nPhaseIdx], factorW, intersection, wPhaseIdx);
700 this->evalCflFluxFunction().addFlux(lambdaW, lambdaNw, viscosity_[wPhaseIdx],
701 viscosity_[nPhaseIdx], factorNw, intersection, nPhaseIdx);
707 switch (velocityType_)
710 update -= factorTotal / (volume *
porosity);
713 switch (saturationType_)
716 update -= factorW / (volume *
porosity);
719 update -= factorNw / (volume *
porosity);
734template<
class TypeTag>
737 auto elementI = intersection.inside();
740 const GlobalPosition& globalPosI = elementI.geometry().center();
743 const GlobalPosition& globalPosJ = intersection.geometry().center();
746 Scalar volume = elementI.geometry().volume();
748 Scalar
porosity = max(problem_.spatialParams().porosity(elementI), porosityThreshold_);
750 if (compressibility_)
752 viscosity_[wPhaseIdx] = cellDataI.viscosity(wPhaseIdx);
753 viscosity_[nPhaseIdx] = cellDataI.viscosity(nPhaseIdx);
757 int isIndex = intersection.indexInInside();
759 GlobalPosition unitOuterNormal = intersection.centerUnitOuterNormal();
761 unitOuterNormal *= -1.0;
763 Scalar faceArea = intersection.geometry().volume();
765 if (velocity().calculateVelocityInTransport())
766 velocity().calculateVelocityOnBoundary(intersection, cellDataI);
769 Scalar factorW = (cellDataI.fluxData().velocity(wPhaseIdx, isIndex) * unitOuterNormal) * faceArea;
770 Scalar factorNw = (cellDataI.fluxData().velocity(nPhaseIdx, isIndex) * unitOuterNormal) * faceArea;
771 Scalar factorTotal = factorW + factorNw;
774 GlobalPosition distVec = globalPosJ - globalPosI;
777 Scalar dist = distVec.two_norm();
780 BoundaryTypes bcType;
781 problem_.boundaryTypes(bcType, intersection);
782 PrimaryVariables boundValues(0.0);
784 if (bcType.isDirichlet(satEqIdx))
786 problem_.dirichlet(boundValues, intersection);
788 Scalar satBound = boundValues[saturationIdx];
791 Scalar satWI = cellDataI.saturation(wPhaseIdx);
792 Scalar satWBound = 0;
793 switch (saturationType_)
797 satWBound = satBound;
802 satWBound = 1 - satBound;
810 const auto fluidMatrixInteraction = Deprecated::makePcKrSw(Scalar{}, problem_.spatialParams(), elementI);
812 const Scalar pcBound = fluidMatrixInteraction.pc(satWBound);
818 if (cellDataI.fluxData().isUpwindCell(wPhaseIdx, isIndex))
820 lambdaW = cellDataI.mobility(wPhaseIdx);
821 if (compressibility_)
823 lambdaW /= cellDataI.density(wPhaseIdx);
828 if (compressibility_)
830 lambdaW = fluidMatrixInteraction.krw(satWBound)
835 lambdaW = fluidMatrixInteraction.krw(satWBound)
836 / viscosity_[wPhaseIdx];
840 if (cellDataI.fluxData().isUpwindCell(nPhaseIdx, isIndex))
842 lambdaNw = cellDataI.mobility(nPhaseIdx);
843 if (compressibility_)
845 lambdaNw /= cellDataI.density(nPhaseIdx);
850 if (compressibility_)
852 lambdaNw = fluidMatrixInteraction.krn(satWBound)
857 lambdaNw = fluidMatrixInteraction.krn(satWBound)
858 / viscosity_[nPhaseIdx];
862 switch (velocityType_)
866 this->evalCflFluxFunction().addFlux(lambdaW, lambdaNw, viscosity_[wPhaseIdx],
867 viscosity_[nPhaseIdx], factorTotal, intersection);
869 Scalar pcI = cellDataI.capillaryPressure();
872 GlobalPosition pcGradient = unitOuterNormal;
873 pcGradient *= (pcI - pcBound) / dist;
877 capillaryFlux().getFlux(flux, intersection, satWI, satWBound, pcGradient);
878 Scalar capillaryFlux = flux * unitOuterNormal * faceArea;
881 gravityFlux().getFlux(flux, intersection, satWI, satWBound);
882 Scalar gravityFlux = flux * unitOuterNormal * faceArea;
884 switch (saturationType_)
889 factorTotal *= lambdaW / (lambdaW + lambdaNw);
895 factorTotal *= lambdaNw / (lambdaW + lambdaNw);
897 this->evalCflFluxFunction().addFlux(lambdaW, lambdaNw, viscosity_[wPhaseIdx],
898 viscosity_[nPhaseIdx], factorW, intersection, wPhaseIdx);
899 this->evalCflFluxFunction().addFlux(lambdaW, lambdaNw, viscosity_[wPhaseIdx],
900 viscosity_[nPhaseIdx], factorNw, intersection, nPhaseIdx);
906 factorTotal -= capillaryFlux;
907 factorTotal += gravityFlux;
910 this->evalCflFluxFunction().addFlux(lambdaW, lambdaNw, viscosity_[wPhaseIdx],
911 viscosity_[nPhaseIdx], 10 * capillaryFlux, intersection);
912 this->evalCflFluxFunction().addFlux(lambdaW, lambdaNw, viscosity_[wPhaseIdx],
913 viscosity_[nPhaseIdx], 10 * gravityFlux, intersection);
919 if (compressibility_)
921 factorW /= cellDataI.density(wPhaseIdx);
922 factorNw /= cellDataI.density(nPhaseIdx);
926 this->evalCflFluxFunction().addFlux(lambdaW, lambdaNw, viscosity_[wPhaseIdx],
927 viscosity_[nPhaseIdx], factorW, intersection, wPhaseIdx);
928 this->evalCflFluxFunction().addFlux(lambdaW, lambdaNw, viscosity_[wPhaseIdx],
929 viscosity_[nPhaseIdx], factorNw, intersection, nPhaseIdx);
936 if (bcType.isNeumann(satEqIdx))
938 problem_.neumann(boundValues, intersection);
939 factorW = boundValues[wPhaseIdx];
940 factorNw = boundValues[nPhaseIdx];
942 factorNw *= faceArea;
945 Scalar lambdaW, lambdaNw;
947 lambdaW = cellDataI.mobility(wPhaseIdx);
948 lambdaNw = cellDataI.mobility(nPhaseIdx);
949 if (compressibility_)
951 lambdaW /= cellDataI.density(wPhaseIdx);
952 lambdaNw /= cellDataI.density(nPhaseIdx);
953 factorW /= cellDataI.density(wPhaseIdx);
954 factorNw /= cellDataI.density(nPhaseIdx);
958 factorW /= density_[wPhaseIdx];
959 factorNw /= density_[nPhaseIdx];
962 switch (velocityType_)
967 this->evalCflFluxFunction().addFlux(lambdaW, lambdaNw, viscosity_[wPhaseIdx],
968 viscosity_[nPhaseIdx], factorW + factorNw, intersection);
973 this->evalCflFluxFunction().addFlux(lambdaW, lambdaNw, viscosity_[wPhaseIdx],
974 viscosity_[nPhaseIdx], factorW, intersection, wPhaseIdx);
975 this->evalCflFluxFunction().addFlux(lambdaW, lambdaNw, viscosity_[wPhaseIdx],
976 viscosity_[nPhaseIdx], factorNw, intersection, nPhaseIdx);
983 if (bcType.isOutflow(satEqIdx))
986 Scalar lambdaW = cellDataI.mobility(wPhaseIdx);
987 Scalar lambdaNw = cellDataI.mobility(nPhaseIdx);
988 if (compressibility_)
990 lambdaW /= cellDataI.density(wPhaseIdx);
991 lambdaNw /= cellDataI.density(nPhaseIdx);
994 if (velocityType_ == vt)
996 switch (saturationType_)
1001 factorTotal *= lambdaW / (lambdaW + lambdaNw);
1002 this->evalCflFluxFunction().addFlux(lambdaW, lambdaNw, viscosity_[wPhaseIdx],
1003 viscosity_[nPhaseIdx], factorTotal, intersection);
1009 factorTotal *= lambdaNw / (lambdaW + lambdaNw);
1010 this->evalCflFluxFunction().addFlux(lambdaW, lambdaNw, viscosity_[wPhaseIdx],
1011 viscosity_[nPhaseIdx], factorTotal, intersection);
1018 this->evalCflFluxFunction().addFlux(lambdaW, lambdaNw, viscosity_[wPhaseIdx],
1019 viscosity_[nPhaseIdx], factorW, intersection, wPhaseIdx);
1020 this->evalCflFluxFunction().addFlux(lambdaW, lambdaNw, viscosity_[wPhaseIdx],
1021 viscosity_[nPhaseIdx], factorNw, intersection, nPhaseIdx);
1024 switch (velocityType_)
1027 update -= factorTotal / (volume *
porosity);
1030 switch (saturationType_)
1033 update -= factorW / (volume *
porosity);
1036 update -= factorNw / (volume *
porosity);
1050template<
class TypeTag>
1054 Scalar volume = element.geometry().volume();
1057 Scalar
porosity = max(problem_.spatialParams().porosity(element), porosityThreshold_);
1059 if (compressibility_)
1061 viscosity_[wPhaseIdx] = cellDataI.viscosity(wPhaseIdx);
1062 viscosity_[nPhaseIdx] = cellDataI.viscosity(nPhaseIdx);
1065 PrimaryVariables sourceVec(0.0);
1066 problem_.source(sourceVec, element);
1068 if (compressibility_)
1070 sourceVec[wPhaseIdx] /= cellDataI.density(wPhaseIdx);
1071 sourceVec[nPhaseIdx] /= cellDataI.density(nPhaseIdx);
1075 sourceVec[wPhaseIdx] /= density_[wPhaseIdx];
1076 sourceVec[nPhaseIdx] /= density_[nPhaseIdx];
1080 Scalar lambdaW = cellDataI.mobility(wPhaseIdx);
1081 Scalar lambdaNw = cellDataI.mobility(nPhaseIdx);
1082 if (compressibility_)
1084 lambdaW /= cellDataI.density(wPhaseIdx);
1085 lambdaNw /= cellDataI.density(nPhaseIdx);
1088 switch (saturationType_)
1092 if (sourceVec[wPhaseIdx] < 0 && cellDataI.saturation(wPhaseIdx) < threshold_)
1093 sourceVec[wPhaseIdx] = 0.0;
1095 update += sourceVec[wPhaseIdx] /
porosity;
1100 if (sourceVec[nPhaseIdx] < 0 && cellDataI.saturation(nPhaseIdx) < threshold_)
1101 sourceVec[nPhaseIdx] = 0.0;
1103 update += sourceVec[nPhaseIdx] /
porosity;
1108 switch (velocityType_)
1113 this->evalCflFluxFunction().addFlux(lambdaW, lambdaNw, viscosity_[wPhaseIdx], viscosity_[nPhaseIdx],
1114 (sourceVec[wPhaseIdx] + sourceVec[nPhaseIdx]) * -1 * volume, element);
1120 this->evalCflFluxFunction().addFlux(lambdaW, lambdaNw, viscosity_[wPhaseIdx],
1121 viscosity_[nPhaseIdx], sourceVec[wPhaseIdx] * -1 * volume, element,
1123 this->evalCflFluxFunction().addFlux(lambdaW, lambdaNw, viscosity_[wPhaseIdx],
1124 viscosity_[nPhaseIdx], sourceVec[nPhaseIdx] * -1 * volume, element,
1132template<
class TypeTag>
1135 ParentType::initialize();
1137 if (!compressibility_)
1139 const auto element = *problem_.gridView().template begin<0>();
1140 FluidState fluidState;
1141 fluidState.setPressure(wPhaseIdx, problem_.referencePressure(element));
1142 fluidState.setPressure(nPhaseIdx, problem_.referencePressure(element));
1143 fluidState.setTemperature(problem_.temperature(element));
1144 fluidState.setSaturation(wPhaseIdx, 1.);
1145 fluidState.setSaturation(nPhaseIdx, 0.);
1153 for (
const auto& element : elements(problem_.gridView()))
1155 PrimaryVariables initSol(0.0);
1156 problem_.initial(initSol, element);
1158 int eIdxGlobal = problem_.variables().index(element);
1160 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
1162 switch (saturationType_)
1166 cellData.setSaturation(wPhaseIdx, initSol[saturationIdx]);
1167 cellData.setSaturation(nPhaseIdx, 1 - initSol[saturationIdx]);
1172 cellData.setSaturation(wPhaseIdx,1 -initSol[saturationIdx]);
1173 cellData.setSaturation(nPhaseIdx, initSol[saturationIdx]);
1180 velocity_->initialize();
1181 capillaryFlux_->initialize();
1182 gravityFlux_->initialize();
1192template<
class TypeTag>
1196 for (
const auto& element : elements(problem_.gridView()))
1198 int eIdxGlobal = problem_.variables().index(element);
1200 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
1203 Scalar satW = cellData.saturation(wPhaseIdx);
1208 const auto fluidMatrixInteraction = Deprecated::makePcKrSw(Scalar{}, problem_.spatialParams(), element);
1210 const Scalar pc = fluidMatrixInteraction.pc(satW);
1212 cellData.setCapillaryPressure(pc);
1215 const Scalar mobilityW = fluidMatrixInteraction.krw(satW) / viscosity_[wPhaseIdx];
1216 const Scalar mobilityNw = fluidMatrixInteraction.krn(satW) / viscosity_[nPhaseIdx];
1219 cellData.setMobility(wPhaseIdx, mobilityW);
1220 cellData.setMobility(nPhaseIdx, mobilityNw);
1223 cellData.setFracFlowFunc(wPhaseIdx, mobilityW / (mobilityW + mobilityNw));
1224 cellData.setFracFlowFunc(nPhaseIdx, mobilityNw / (mobilityW + mobilityNw));
Finite Volume discretization of a transport equation.
T getParam(Args &&... args)
A free function to get a parameter from the parameter tree singleton.
Definition: parameters.hh:364
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 saturation(int phaseIdx) noexcept
I/O name of saturation for multiphase systems.
Definition: name.hh:43
std::string viscosity(int phaseIdx) noexcept
I/O name of viscosity for multiphase systems.
Definition: name.hh:74
std::string density(int phaseIdx) noexcept
I/O name of density for multiphase systems.
Definition: name.hh:65
std::string porosity() noexcept
I/O name of porosity.
Definition: name.hh:139
The finite volume discretization of a saturation transport equation.
Definition: saturation.hh:78
const CapillaryFlux & capillaryFlux() const
Definition: saturation.hh:144
void getFlux(Scalar &update, const Intersection &intersection, CellData &cellDataI)
Function which calculates the flux update.
Definition: saturation.hh:542
Velocity & velocity()
Definition: saturation.hh:160
GravityFlux & gravityFlux()
Definition: saturation.hh:149
void updateTransportedQuantity(TransportSolutionType &updateVec)
Updates the primary transport variable.
Definition: saturation.hh:265
void updateSaturationSolution(int eIdxGlobal, Scalar update, Scalar dt)
Updates the saturation solution of a cell.
Definition: saturation.hh:323
void getFluxOnBoundary(Scalar &update, const Intersection &intersection, CellData &cellDataI)
Function which calculates the boundary flux update.
Definition: saturation.hh:735
void updateSaturationSolution(TransportSolutionType &updateVec)
Globally updates the saturation solution.
Definition: saturation.hh:289
CapillaryFlux & capillaryFlux()
Definition: saturation.hh:139
Velocity & velocity() const
Definition: saturation.hh:165
void updateSaturationSolution(TransportSolutionType &updateVec, Scalar dt)
Globally updates the saturation solution.
Definition: saturation.hh:305
void initialize()
Sets the initial solution.
Definition: saturation.hh:1133
void getTransportedQuantity(TransportSolutionType &transportedQuantity)
Writes the current values of the primary transport variable into the transportedQuantity-vector (come...
Definition: saturation.hh:192
const GravityFlux & gravityFlux() const
Definition: saturation.hh:154
void serializeEntity(std::ostream &outstream, const Element &element)
Function for serialization of the primary transport variable.
Definition: saturation.hh:422
void deserializeEntity(std::istream &instream, const Element &element)
Function for deserialization of the primary transport variable.
Definition: saturation.hh:445
void updateMaterialLaws()
Update the values of the material laws and constitutive relations.
Definition: saturation.hh:1193
void addOutputVtkFields(MultiWriter &writer)
Adds saturation output to the output file.
Definition: saturation.hh:361
FVSaturation2P(Problem &problem)
Constructs a FVSaturation2P object.
Definition: saturation.hh:471
void getSource(Scalar &update, const Element &element, CellData &cellDataI)
Function which calculates the source update.
Definition: saturation.hh:1051
void setTransportedQuantity(TransportSolutionType &transportedQuantity)
Writes the current values of the primary transport variable into the variable container.
Definition: saturation.hh:220
void updateTransportedQuantity(TransportSolutionType &updateVec, Scalar dt)
Updates the primary transport variable.
Definition: saturation.hh:279
bool inPhysicalRange(DataEntry &entry)
Check if saturation is in physical range.
Definition: saturation.hh:255
The finite volume discretization of a transport equation.
Definition: transport.hh:52
bool enableLocalTimeStepping()
Definition: transport.hh:202
void update(const Scalar t, Scalar &dt, TransportSolutionType &updateVec, bool impet=false)
Calculate the update vector.
Definition: transport.hh:270
void innerUpdate(TransportSolutionType &updateVec)
Definition: transport.hh:564
Specifies the properties for immiscible 2p transport.