24#ifndef DUMUX_FVTRANSPORT2P2C_HH
25#define DUMUX_FVTRANSPORT2P2C_HH
28#include <unordered_map>
30#include <dune/grid/common/gridenums.hh>
31#include <dune/common/float_cmp.hh>
58template<
class TypeTag>
66 using SpatialParams =
typename GET_PROP_TYPE(TypeTag, SpatialParams);
67 using MaterialLaw =
typename SpatialParams::MaterialLaw;
69 using Indices =
typename GET_PROP_TYPE(TypeTag, ModelTraits)::Indices;
70 using BoundaryTypes =
typename GET_PROP_TYPE(TypeTag, BoundaryTypes);
72 using FluidSystem =
typename GET_PROP_TYPE(TypeTag, FluidSystem);
73 using FluidState =
typename GET_PROP_TYPE(TypeTag, FluidState);
77 using TransportSolutionType =
typename GET_PROP_TYPE(TypeTag, TransportSolutionType);
81 dim = GridView::dimension, dimWorld = GridView::dimensionworld
85 pw = Indices::pressureW,
86 pn = Indices::pressureN
90 wPhaseIdx = Indices::wPhaseIdx, nPhaseIdx = Indices::nPhaseIdx,
91 wCompIdx = Indices::wPhaseIdx, nCompIdx = Indices::nPhaseIdx,
92 contiWEqIdx=Indices::contiWEqIdx, contiNEqIdx=Indices::contiNEqIdx,
97 using Element =
typename GridView::Traits::template Codim<0>::Entity;
98 using Intersection =
typename GridView::Intersection;
100 using GlobalPosition = Dune::FieldVector<Scalar, dimWorld>;
101 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
102 using PhaseVector = Dune::FieldVector<Scalar, NumPhases>;
103 using ComponentVector = Dune::FieldVector<Scalar, NumComponents>;
104 using PrimaryVariables =
typename GET_PROP_TYPE(TypeTag, PrimaryVariables);
125 void innerUpdate(TransportSolutionType& updateVec);
128 virtual void update(
const Scalar t, Scalar& dt, TransportSolutionType& updateVec,
bool impet =
false);
138 const Intersection& intersection, CellData& cellDataI);
142 const Intersection& intersection,
const CellData& cellDataI);
145 const Intersection& intersection,
146 FluidState& BCfluidState,
147 PhaseVector& pressBound);
159 for (
int eqNumber = 0; eqNumber < transportedQuantities; eqNumber++)
170 template<
class MultiWriter>
173 if(
problem().vtkOutputLevel()>3)
176 int size =
problem_.gridView().size(0);
177 ScalarSolutionType *totalC1PV = writer.allocateManagedBuffer(size);
178 ScalarSolutionType *totalC2PV = writer.allocateManagedBuffer(size);
181 writer.attachCellData(*totalC1PV,
"total Concentration w-Comp - vector");
182 writer.attachCellData(*totalC2PV,
"total Concentration n-Comp - vector");
189 int eIdxGlobal =
problem().variables().index(element);
196 int eIdxGlobal =
problem().variables().index(element);
197 CellData& cellData =
problem().variables().cellData(eIdxGlobal);
245 template<
class DataEntry>
249 for(
int compIdx = 0; compIdx < numComp; compIdx++)
251 if (entry[compIdx] < -1.0e-6)
283 Scalar cFLFactor = getParam<Scalar>(
"Impet.CFLFactor");
285 subCFLFactor_ = min(getParam<Scalar>(
"Impet.SubCFLFactor"), cFLFactor);
286 verbosity_ = getParam<int>(
"TimeManager.SubTimestepVerbosity");
291 std::cout<<
"max CFL-Number of "<<cFLFactor<<
", max Sub-CFL-Number of "
292 <<
subCFLFactor_<<
": Enable local time-stepping!" << std::endl;
336 Implementation &asImp_()
337 {
return *
static_cast<Implementation *
>(
this); }
340 const Implementation &asImp_()
const
341 {
return *
static_cast<const Implementation *
>(
this); }
362template<
class TypeTag>
364 TransportSolutionType& updateVec,
bool impet)
369 unsigned int size = problem_.gridView().size(0);
370 if (localTimeStepping_)
372 if (this->timeStepData_.size() != size)
373 this->timeStepData_.resize(size);
381 updateVec[wCompIdx].resize(problem_.gridView().size(0));
382 updateVec[nCompIdx].resize(problem_.gridView().size(0));
383 updateVec[wCompIdx] = 0;
384 updateVec[nCompIdx] = 0;
387 int restrictingCell = -1;
389 ComponentVector entries(0.);
392 for (
const auto& element : elements(problem().gridView()))
395 int eIdxGlobalI = problem().variables().index(element);
396 CellData& cellDataI = problem().variables().cellData(eIdxGlobalI);
399 double sumfactorin = 0;
400 double sumfactorout = 0;
403 for (
const auto& intersection :
intersections(problem().gridView(), element))
405 int indexInInside = intersection.indexInInside();
408 if (intersection.neighbor())
409 getFlux(entries, timestepFlux, intersection, cellDataI);
412 if (intersection.boundary())
413 getFluxOnBoundary(entries, timestepFlux, intersection, cellDataI);
416 if (localTimeStepping_)
419 if (localData.
faceTargetDt[indexInInside] < accumulatedDt_ + dtThreshold_)
421 localData.
faceFluxes[indexInInside] = entries;
427 updateVec[wCompIdx][eIdxGlobalI] += entries[wCompIdx];
428 updateVec[nCompIdx][eIdxGlobalI] += entries[nCompIdx];
432 sumfactorin += timestepFlux[0];
433 sumfactorout += timestepFlux[1];
437 if (localTimeStepping_)
440 for (
int i=0; i < 2*dim; i++)
442 updateVec[wCompIdx][eIdxGlobalI] += localData.
faceFluxes[i][wCompIdx];
443 updateVec[nCompIdx][eIdxGlobalI] += localData.
faceFluxes[i][nCompIdx];
448 PrimaryVariables q(NAN);
449 problem().source(q, element);
450 updateVec[wCompIdx][eIdxGlobalI] += q[contiWEqIdx];
451 updateVec[nCompIdx][eIdxGlobalI] += q[contiNEqIdx];
455 sumfactorin = max(sumfactorin,sumfactorout)
456 / problem().spatialParams().porosity(element);
459 if (localTimeStepping_)
461 timeStepData_[eIdxGlobalI].dt = 1./sumfactorin;
462 if ( 1./sumfactorin < dt)
465 restrictingCell= eIdxGlobalI;
470 if ( 1./sumfactorin < dt)
473 restrictingCell= eIdxGlobalI;
481 using ElementMapper =
typename SolutionTypes::ElementMapper;
483 for (
int i = 0; i < updateVec.size(); i++)
485 DataHandle dataHandle(problem_.variables().elementMapper(), updateVec[i]);
486 problem_.gridView().template communicate<DataHandle>(dataHandle,
487 Dune::InteriorBorder_All_Interface,
488 Dune::ForwardCommunication);
491 if (localTimeStepping_)
495 TimeDataHandle timeDataHandle(problem_.elementMapper(), timeStepData_);
496 problem_.gridView().template communicate<TimeDataHandle>(timeDataHandle,
497 Dune::InteriorBorder_All_Interface,
498 Dune::ForwardCommunication);
501 dt = problem_.gridView().comm().min(dt);
506 Dune::dinfo <<
"Timestep restricted by CellIdx " << restrictingCell <<
" leads to dt = "
507 <<dt * getParam<Scalar>(
"Impet.CFLFactor")<< std::endl;
508 if (averagedFaces_ != 0)
509 Dune::dinfo <<
" Averageing done for " << averagedFaces_ <<
" faces. "<< std::endl;
518template<
class TypeTag>
521 if (this->enableLocalTimeStepping())
522 this->innerUpdate(updateVector);
524 updateConcentrations(updateVector, problem().timeManager().timeStepSize());
534template<
class TypeTag>
537 updateConcentrations(updateVector, dt);
547template<
class TypeTag>
551 for (
int i = 0; i< problem().gridView().size(0); i++)
553 CellData& cellDataI = problem().variables().cellData(i);
556 totalConcentration_[compIdx][i] += (updateVector[compIdx][i]*=dt);
557 cellDataI.setMassConcentration(compIdx, totalConcentration_[compIdx][i]);
578template<
class TypeTag>
580 Dune::FieldVector<Scalar, 2>& timestepFlux,
581 const Intersection& intersection,
587 auto elementI = intersection.inside();
588 int eIdxGlobalI = problem().variables().index(elementI);
591 const GlobalPosition globalPos = elementI.geometry().center();
592 const GlobalPosition& gravity_ = problem().gravity();
594 Scalar volume = elementI.geometry().volume();
597 Scalar pressI = problem().pressureModel().pressure(eIdxGlobalI);
598 Scalar pcI = cellDataI.capillaryPressure();
599 DimMatrix K_I(problem().spatialParams().intrinsicPermeability(elementI));
601 PhaseVector SmobI(0.);
603 SmobI[wPhaseIdx] = max((cellDataI.saturation(wPhaseIdx)
604 - problem().spatialParams().materialLawParams(elementI).swr())
606 SmobI[nPhaseIdx] = max((cellDataI.saturation(nPhaseIdx)
607 - problem().spatialParams().materialLawParams(elementI).snr())
610 Scalar densityWI (0.), densityNWI(0.);
611 densityWI= cellDataI.density(wPhaseIdx);
612 densityNWI = cellDataI.density(nPhaseIdx);
615 GlobalPosition unitOuterNormal = intersection.centerUnitOuterNormal();
617 unitOuterNormal *= -1.0;
618 Scalar faceArea = intersection.geometry().volume();
625 Dune::FieldVector<Scalar, 2> factor (0.);
626 Dune::FieldVector<Scalar, 2> updFactor (0.);
628 PhaseVector potential(0.);
631 auto neighbor = intersection.outside();
632 int eIdxGlobalJ = problem().variables().index(neighbor);
633 CellData& cellDataJ = problem().variables().cellData(eIdxGlobalJ);
636 const GlobalPosition& globalPosNeighbor = neighbor.geometry().center();
639 GlobalPosition distVec = globalPosNeighbor - globalPos;
641 Scalar dist = distVec.two_norm();
643 GlobalPosition unitDistVec(distVec);
647 Scalar densityWJ (0.), densityNWJ(0.);
648 densityWJ = cellDataJ.density(wPhaseIdx);
649 densityNWJ = cellDataJ.density(nPhaseIdx);
652 double densityW_mean = (densityWI + densityWJ) * 0.5;
653 double densityNW_mean = (densityNWI + densityNWJ) * 0.5;
655 double pressJ = problem().pressureModel().pressure(eIdxGlobalJ);
656 Scalar pcJ = cellDataJ.capillaryPressure();
659 DimMatrix meanK_(0.);
662 problem().spatialParams().intrinsicPermeability(neighbor));
663 Dune::FieldVector<Scalar,dim> K(0);
664 meanK_.umv(unitDistVec,K);
667 switch (pressureType)
671 potential[wPhaseIdx] = (pressI - pressJ) / (dist);
672 potential[nPhaseIdx] = (pressI - pressJ + pcI - pcJ) / (dist);
677 potential[wPhaseIdx] = (pressI - pressJ - pcI + pcJ) / (dist);
678 potential[nPhaseIdx] = (pressI - pressJ) / (dist);
683 potential[nPhaseIdx] += (gravity_ * unitDistVec) * densityNW_mean;
684 potential[wPhaseIdx] += (gravity_ * unitDistVec) * densityW_mean;
686 potential[wPhaseIdx] *= fabs(K * unitOuterNormal);
687 potential[nPhaseIdx] *= fabs(K * unitOuterNormal);
690 Dune::FieldVector<bool, NumPhases> doUpwinding(
true);
691 PhaseVector lambda(0.);
692 for(
int phaseIdx = 0; phaseIdx <
NumPhases; phaseIdx++)
695 if(phaseIdx == wPhaseIdx)
696 contiEqIdx = contiWEqIdx;
698 contiEqIdx = contiNEqIdx;
700 if(!impet_ || restrictFluxInTransport_==0)
702 if (potential[phaseIdx] > 0.)
704 lambda[phaseIdx] = cellDataI.mobility(phaseIdx);
705 cellDataI.setUpwindCell(intersection.indexInInside(), contiEqIdx,
true);
708 else if(potential[phaseIdx] < 0.)
710 lambda[phaseIdx] = cellDataJ.mobility(phaseIdx);
711 cellDataI.setUpwindCell(intersection.indexInInside(), contiEqIdx,
false);
715 doUpwinding[phaseIdx] =
false;
716 cellDataI.setUpwindCell(intersection.indexInInside(), contiEqIdx,
false);
717 cellDataJ.setUpwindCell(intersection.indexInOutside(), contiEqIdx,
false);
722 bool wasUpwindCell = cellDataI.isUpwindCell(intersection.indexInInside(), contiEqIdx);
724 if (potential[phaseIdx] > 0. && wasUpwindCell)
725 lambda[phaseIdx] = cellDataI.mobility(phaseIdx);
726 else if (potential[phaseIdx] < 0. && !wasUpwindCell)
727 lambda[phaseIdx] = cellDataJ.mobility(phaseIdx);
729 else if(restrictFluxInTransport_ == 2)
730 doUpwinding[phaseIdx] =
false;
734 if (potential[phaseIdx] > 0. && (Dune::FloatCmp::ne<Scalar, Dune::FloatCmp::absolute>(cellDataJ.mobility(phaseIdx), 0.0, 1.0e-30)
735 || (cellDataI.wasRefined() && cellDataJ.wasRefined() && elementI.father() == neighbor.father())))
736 lambda[phaseIdx] = cellDataI.mobility(phaseIdx);
737 else if (potential[phaseIdx] < 0. && (Dune::FloatCmp::ne<Scalar, Dune::FloatCmp::absolute>(cellDataI.mobility(phaseIdx), 0.0, 1.0e-30)
738 || (cellDataI.wasRefined() && cellDataJ.wasRefined() && elementI.father() == neighbor.father())))
739 lambda[phaseIdx] = cellDataJ.mobility(phaseIdx);
741 doUpwinding[phaseIdx] =
false;
747 if(!doUpwinding[phaseIdx])
750 if(cellDataI.mobility(phaseIdx)+cellDataJ.mobility(phaseIdx)==0.)
752 potential[phaseIdx] = 0;
757 fluxEntries[wCompIdx] -= potential[phaseIdx] * faceArea / volume
758 *
harmonicMean(cellDataI.massFraction(phaseIdx, wCompIdx) * cellDataI.mobility(phaseIdx) * cellDataI.density(phaseIdx),
759 cellDataJ.massFraction(phaseIdx, wCompIdx) * cellDataJ.mobility(phaseIdx) * cellDataJ.density(phaseIdx));
760 fluxEntries[nCompIdx] -= potential[phaseIdx] * faceArea / volume
761 *
harmonicMean(cellDataI.massFraction(phaseIdx, nCompIdx) * cellDataI.mobility(phaseIdx) * cellDataI.density(phaseIdx),
762 cellDataJ.massFraction(phaseIdx, nCompIdx) * cellDataJ.mobility(phaseIdx) * cellDataJ.density(phaseIdx));
766 timestepFlux[0] += max(0.,
767 - potential[phaseIdx] * faceArea / volume
768 *
harmonicMean(cellDataI.mobility(phaseIdx),cellDataJ.mobility(phaseIdx)));
770 timestepFlux[1] += max(0.,
771 potential[phaseIdx] * faceArea / volume
772 *
harmonicMean(cellDataI.mobility(phaseIdx),cellDataJ.mobility(phaseIdx))/SmobI[phaseIdx]);
775 if(!(cellDataI.wasRefined() && cellDataJ.wasRefined() && elementI.father() == neighbor.father())
776 && eIdxGlobalI > eIdxGlobalJ)
779 #if DUNE_MINIMAL_DEBUG_LEVEL < 3
781 if(eIdxGlobalI > eIdxGlobalJ)
782 Dune::dinfo <<
"harmonicMean flux of phase" << phaseIdx <<
" used from cell" << eIdxGlobalI<<
" into " << eIdxGlobalJ
783 <<
" ; TE upwind I = "<< cellDataI.isUpwindCell(intersection.indexInInside(), contiEqIdx)
784 <<
" but pot = "<< potential[phaseIdx] << std::endl;
789 potential[phaseIdx] = 0;
795 double velocityJIw = max((-lambda[wPhaseIdx] * potential[wPhaseIdx]) * faceArea / volume, 0.0);
796 double velocityIJw = max(( lambda[wPhaseIdx] * potential[wPhaseIdx]) * faceArea / volume, 0.0);
797 double velocityJIn = max((-lambda[nPhaseIdx] * potential[nPhaseIdx]) * faceArea / volume, 0.0);
798 double velocityIJn = max(( lambda[nPhaseIdx] * potential[nPhaseIdx]) * faceArea / volume, 0.0);
801 timestepFlux[0] += velocityJIw + velocityJIn;
803 double foutw = velocityIJw/SmobI[wPhaseIdx];
804 double foutn = velocityIJn/SmobI[nPhaseIdx];
807 if (isnan(foutw) || isinf(foutw) || foutw < 0) foutw = 0;
808 if (isnan(foutn) || isinf(foutn) || foutn < 0) foutn = 0;
809 timestepFlux[1] += foutw + foutn;
811 fluxEntries[wCompIdx] +=
812 velocityJIw * cellDataJ.massFraction(wPhaseIdx, wCompIdx) * densityWJ
813 - velocityIJw * cellDataI.massFraction(wPhaseIdx, wCompIdx) * densityWI
814 + velocityJIn * cellDataJ.massFraction(nPhaseIdx, wCompIdx) * densityNWJ
815 - velocityIJn * cellDataI.massFraction(nPhaseIdx, wCompIdx) * densityNWI;
816 fluxEntries[nCompIdx] +=
817 velocityJIw * cellDataJ.massFraction(wPhaseIdx, nCompIdx) * densityWJ
818 - velocityIJw * cellDataI.massFraction(wPhaseIdx, nCompIdx) * densityWI
819 + velocityJIn * cellDataJ.massFraction(nPhaseIdx, nCompIdx) * densityNWJ
820 - velocityIJn * cellDataI.massFraction(nPhaseIdx, nCompIdx) * densityNWI;
883template<
class TypeTag>
885 Dune::FieldVector<Scalar, 2>& timestepFlux,
886 const Intersection& intersection,
887 const CellData& cellDataI)
891 auto elementI = intersection.inside();
892 int eIdxGlobalI = problem().variables().index(elementI);
895 const GlobalPosition globalPos = elementI.geometry().center();
898 Scalar volume = elementI.geometry().volume();
899 const GlobalPosition& gravity_ = problem().gravity();
901 Scalar pressI = problem().pressureModel().pressure(eIdxGlobalI);
902 Scalar pcI = cellDataI.capillaryPressure();
903 DimMatrix K_I(problem().spatialParams().intrinsicPermeability(elementI));
905 if(regulateBoundaryPermeability)
907 int axis = intersection.indexInInside() / 2;
908 if(K_I[axis][axis] < minimalBoundaryPermeability)
909 K_I[axis][axis] = minimalBoundaryPermeability;
912 Scalar SwmobI = max((cellDataI.saturation(wPhaseIdx)
913 - problem().spatialParams().materialLawParams(elementI).swr())
915 Scalar SnmobI = max((cellDataI.saturation(nPhaseIdx)
916 - problem().spatialParams().materialLawParams(elementI).snr())
919 Scalar densityWI (0.), densityNWI(0.);
920 densityWI= cellDataI.density(wPhaseIdx);
921 densityNWI = cellDataI.density(nPhaseIdx);
924 GlobalPosition unitOuterNormal = intersection.centerUnitOuterNormal();
926 unitOuterNormal *= -1.0;
927 Scalar faceArea = intersection.geometry().volume();
930 Dune::FieldVector<Scalar, 2> factor (0.);
931 Dune::FieldVector<Scalar, 2> updFactor (0.);
933 PhaseVector potential(0.);
935 const GlobalPosition& globalPosFace = intersection.geometry().center();
938 GlobalPosition distVec = globalPosFace - globalPos;
940 Scalar dist = distVec.two_norm();
942 GlobalPosition unitDistVec(distVec);
946 FluidState BCfluidState;
949 BoundaryTypes bcTypes;
950 problem().boundaryTypes(bcTypes, intersection);
953 if (bcTypes.isDirichlet(contiWEqIdx))
956 PhaseVector pressBound(0.);
960 this->evalBoundary(globalPosFace,
966 Scalar densityWBound = BCfluidState.density(wPhaseIdx);
967 Scalar densityNWBound = BCfluidState.density(nPhaseIdx);
971 pcBound = (BCfluidState.pressure(nPhaseIdx) - BCfluidState.pressure(wPhaseIdx));
973 double densityW_mean = (densityWI + densityWBound) / 2;
974 double densityNW_mean = (densityNWI + densityNWBound) / 2;
977 Dune::FieldVector<Scalar,dim> K(0);
978 K_I.umv(unitDistVec,K);
981 switch (pressureType)
985 potential[wPhaseIdx] = (pressI - pressBound[wPhaseIdx]) / dist;
986 potential[nPhaseIdx] = (pressI + pcI - pressBound[wPhaseIdx] - pcBound)
992 potential[wPhaseIdx] = (pressI - pcI - pressBound[nPhaseIdx] + pcBound)
994 potential[nPhaseIdx] = (pressI - pressBound[nPhaseIdx]) / dist;
998 potential[wPhaseIdx] += (gravity_ * unitDistVec) * densityW_mean;
999 potential[nPhaseIdx] += (gravity_ * unitDistVec) * densityNW_mean;
1001 potential[wPhaseIdx] *= fabs(K * unitOuterNormal);
1002 potential[nPhaseIdx] *= fabs(K * unitOuterNormal);
1005 PhaseVector lambda(0.);
1006 if (potential[wPhaseIdx] >= 0.)
1007 lambda[wPhaseIdx] = cellDataI.mobility(wPhaseIdx);
1011 lambda[wPhaseIdx] = BCfluidState.saturation(wPhaseIdx) / viscosityWBound;
1013 lambda[wPhaseIdx] = MaterialLaw::krw(
1014 problem().spatialParams().materialLawParams(elementI), BCfluidState.saturation(wPhaseIdx))
1017 if (potential[nPhaseIdx] >= 0.)
1018 lambda[nPhaseIdx] = cellDataI.mobility(nPhaseIdx);
1022 lambda[nPhaseIdx] = BCfluidState.saturation(nPhaseIdx) / viscosityNWBound;
1024 lambda[nPhaseIdx] = MaterialLaw::krn(
1025 problem().spatialParams().materialLawParams(elementI), BCfluidState.saturation(wPhaseIdx))
1030 double velocityJIw = max((-lambda[wPhaseIdx] * potential[wPhaseIdx]) * faceArea / volume, 0.0);
1031 double velocityIJw = max(( lambda[wPhaseIdx] * potential[wPhaseIdx]) * faceArea / volume, 0.0);
1032 double velocityJIn = max((-lambda[nPhaseIdx] * potential[nPhaseIdx]) * faceArea / volume, 0.0);
1033 double velocityIJn = max(( lambda[nPhaseIdx] * potential[nPhaseIdx]) * faceArea / volume, 0.0);
1036 timestepFlux[0] = velocityJIw + velocityJIn;
1038 double foutw = velocityIJw/SwmobI;
1039 double foutn = velocityIJn/SnmobI;
1042 if (isnan(foutw) || isinf(foutw) || foutw < 0) foutw = 0;
1043 if (isnan(foutn) || isinf(foutn) || foutn < 0) foutn = 0;
1044 timestepFlux[1] = foutw + foutn;
1046 fluxEntries[wCompIdx] =
1047 + velocityJIw * BCfluidState.massFraction(wPhaseIdx, wCompIdx) * densityWBound
1048 - velocityIJw * cellDataI.massFraction(wPhaseIdx, wCompIdx) * densityWI
1049 + velocityJIn * BCfluidState.massFraction(nPhaseIdx, wCompIdx) * densityNWBound
1050 - velocityIJn * cellDataI.massFraction(nPhaseIdx, wCompIdx) * densityNWI ;
1051 fluxEntries[nCompIdx] =
1052 velocityJIw * BCfluidState.massFraction(wPhaseIdx, nCompIdx) * densityWBound
1053 - velocityIJw * cellDataI.massFraction(wPhaseIdx, nCompIdx) * densityWI
1054 + velocityJIn * BCfluidState.massFraction(nPhaseIdx, nCompIdx) * densityNWBound
1055 - velocityIJn * cellDataI.massFraction(nPhaseIdx, nCompIdx) * densityNWI ;
1057 else if (bcTypes.isNeumann(contiWEqIdx))
1060 PrimaryVariables J(NAN);
1061 problem().neumann(J, intersection);
1062 fluxEntries[wCompIdx] = - J[contiWEqIdx] * faceArea / volume;
1063 fluxEntries[nCompIdx] = - J[contiNEqIdx] * faceArea / volume;
1066 #define cflIgnoresNeumann
1067 #ifdef cflIgnoresNeumann
1068 timestepFlux[0] = 0;
1069 timestepFlux[1] = 0;
1071 double inflow = updFactor[wCompIdx] / densityW + updFactor[nCompIdx] / densityNW;
1074 timestepFlux[0] = updFactor[wCompIdx] / densityW
1075 + updFactor[nCompIdx] / densityNW;
1076 timestepFlux[1] = -(updFactor[wCompIdx] / densityW /SwmobI
1077 + updFactor[nCompIdx] / densityNW / SnmobI);
1081 timestepFlux[0] = -(updFactor[wCompIdx] / densityW
1082 + updFactor[nCompIdx] / densityNW);
1083 timestepFlux[1] = updFactor[wCompIdx] / densityW /SwmobI
1084 + updFactor[nCompIdx] / densityNW / SnmobI;
1106template<
class TypeTag>
1108 const Intersection& intersection,
1109 FluidState &BCfluidState,
1110 PhaseVector &pressBound)
1115 auto element = intersection.inside();
1117 PrimaryVariables primaryVariablesOnBoundary(0.);
1118 problem().dirichlet(primaryVariablesOnBoundary, intersection);
1121 typename Indices::BoundaryFormulation bcType;
1122 problem().boundaryFormulation(bcType, intersection);
1125 Scalar satBound = primaryVariablesOnBoundary[contiWEqIdx];
1128 Scalar pcBound = MaterialLaw::pc(problem().spatialParams().materialLawParams(element),
1130 switch (pressureType)
1134 pressBound[wPhaseIdx] = primaryVariablesOnBoundary[Indices::pressureEqIdx];
1135 pressBound[nPhaseIdx] = primaryVariablesOnBoundary[Indices::pressureEqIdx] + pcBound;
1140 pressBound[wPhaseIdx] = primaryVariablesOnBoundary[Indices::pressureEqIdx] - pcBound;
1141 pressBound[nPhaseIdx] = primaryVariablesOnBoundary[Indices::pressureEqIdx];
1147 pressBound[wPhaseIdx] = pressBound[nPhaseIdx] = primaryVariablesOnBoundary[Indices::pressureEqIdx];
1150 problem().spatialParams().
porosity(element), problem().temperatureAtPos(globalPosFace));
1152 else if (bcType == Indices::concentration)
1155 pressBound[wPhaseIdx] = pressBound[nPhaseIdx] = primaryVariablesOnBoundary[Indices::pressureEqIdx];
1156 Scalar Z0Bound = primaryVariablesOnBoundary[contiWEqIdx];
1158 problem().spatialParams().
porosity(element), problem().temperatureAtPos(globalPosFace));
1162 Scalar pcBound = MaterialLaw::pc(problem().spatialParams().materialLawParams(element),
1163 BCfluidState.saturation(wPhaseIdx));
1166 for(
int iter=0; iter < maxiter; iter++)
1169 switch (pressureType)
1173 pressBound[wPhaseIdx] = primaryVariablesOnBoundary[Indices::pressureEqIdx];
1174 pressBound[nPhaseIdx] = primaryVariablesOnBoundary[Indices::pressureEqIdx]
1180 pressBound[wPhaseIdx] = primaryVariablesOnBoundary[Indices::pressureEqIdx]
1182 pressBound[nPhaseIdx] = primaryVariablesOnBoundary[Indices::pressureEqIdx];
1188 Scalar oldPc = pcBound;
1191 problem().spatialParams().
porosity(element), problem().temperatureAtPos(globalPosFace));
1192 pcBound = MaterialLaw::pc(problem().spatialParams().materialLawParams(element),
1193 BCfluidState.saturation(wPhaseIdx));
1197 if (abs(oldPc - pcBound) < 10.0)
1203 DUNE_THROW(Dune::NotImplemented,
"Boundary Formulation neither Concentration nor Saturation??");
1206template<
class TypeTag>
1209 dt = std::numeric_limits<Scalar>::max();
1212 for (
const auto& element : elements(problem_.gridView()))
1215 if (element.partitionType() != Dune::InteriorEntity)
1222 int eIdxGlobalI = problem_.variables().index(element);
1227 using FaceDt = std::unordered_map<int, Scalar>;
1231 for (
const auto& intersection :
intersections(problem_.gridView(), element))
1233 int indexInInside = intersection.indexInInside();
1235 if (intersection.neighbor())
1237 auto neighbor = intersection.outside();
1238 int eIdxGlobalJ = problem_.variables().index(neighbor);
1240 int levelI = element.level();
1241 int levelJ = neighbor.level();
1244 if (eIdxGlobalI < eIdxGlobalJ && levelI <= levelJ)
1248 int indexInOutside = intersection.indexInOutside();
1250 if (localDataI.
faceTargetDt[indexInInside] < accumulatedDt_ + dtThreshold_
1251 || localDataJ.
faceTargetDt[indexInOutside] < accumulatedDt_ + dtThreshold_)
1253 Scalar timeStep = min(localDataI.
dt, localDataJ.
dt);
1255 if (levelI < levelJ)
1257 typename FaceDt::iterator it = faceDt.find(indexInInside);
1258 if (it != faceDt.end())
1260 it->second = min(it->second, timeStep);
1264 faceDt.insert(std::make_pair(indexInInside, timeStep));
1269 localDataI.
faceTargetDt[indexInInside] += subCFLFactor_ * timeStep;
1270 localDataJ.
faceTargetDt[indexInOutside] += subCFLFactor_ * timeStep;
1273 dt = min(dt, timeStep);
1277 else if (intersection.boundary())
1279 if (localDataI.
faceTargetDt[indexInInside] < accumulatedDt_ + dtThreshold_)
1281 localDataI.
faceTargetDt[indexInInside] += subCFLFactor_ * localDataI.
dt;
1282 dt = min(dt, subCFLFactor_ * localDataI.
dt);
1286 if (faceDt.size() > 0)
1288 typename FaceDt::iterator it = faceDt.begin();
1289 for (;it != faceDt.end();++it)
1291 localDataI.
faceTargetDt[it->first] += subCFLFactor_ * it->second;
1294 for (
const auto& intersection :
intersections(problem_.gridView(), element))
1296 if (intersection.neighbor())
1298 int indexInInside = intersection.indexInInside();
1300 it = faceDt.find(indexInInside);
1301 if (it != faceDt.end())
1303 auto neighbor = intersection.outside();
1304 int eIdxGlobalJ = problem_.variables().index(neighbor);
1308 int indexInOutside = intersection.indexInOutside();
1310 localDataJ.
faceTargetDt[indexInOutside] += subCFLFactor_ * it->second;
1320 using ElementMapper =
typename SolutionTypes::ElementMapper;
1323 TimeDataHandle timeDataHandle(problem_.elementMapper(), timeStepData_);
1324 problem_.gridView().template communicate<TimeDataHandle>(timeDataHandle,
1325 Dune::InteriorBorder_All_Interface,
1326 Dune::ForwardCommunication);
1328 dt = problem_.gridView().comm().min(dt);
1332template<
class TypeTag>
1335 if (localTimeStepping_)
1337 Scalar realDt = problem_.timeManager().timeStepSize();
1339 Scalar subDt = realDt;
1341 updatedTargetDt_(subDt);
1343 Scalar accumulatedDtOld = accumulatedDt_;
1344 accumulatedDt_ += subDt;
1346 Scalar t = problem_.timeManager().time();
1348 if (accumulatedDt_ < realDt)
1353 Scalar dtCorrection = min(0.0, realDt - accumulatedDt_);
1354 subDt += dtCorrection;
1357 std::cout<<
" Sub-time-step size: "<<subDt<< std::endl;
1359 bool stopTimeStep =
false;
1360 int size = problem_.gridView().size(0);
1361 for (
int i = 0; i < size; i++)
1365 for (
int eqNumber = 0; eqNumber < transportedQuantities; eqNumber++)
1367 newVal[eqNumber] = totalConcentration_[eqNumber][i];
1368 newVal[eqNumber] += updateVec[eqNumber][i] * subDt;
1370 if (!asImp_().inPhysicalRange(newVal))
1372 stopTimeStep =
true;
1381 rank = problem_.gridView().comm().rank();
1383 rank = problem_.gridView().comm().max(rank);
1384 problem_.gridView().comm().broadcast(&stopTimeStep,1,rank);
1388 if (stopTimeStep && accumulatedDtOld > dtThreshold_)
1390 problem_.timeManager().setTimeStepSize(accumulatedDtOld);
1395 asImp_().updateTransportedQuantity(updateVec, subDt);
1399 if (accumulatedDt_ >= realDt)
1404 problem_.pressureModel().updateMaterialLaws();
1405 problem_.model().updateTransport(t, subDt, updateVec);
1407 updatedTargetDt_(subDt);
1409 accumulatedDtOld = accumulatedDt_;
1410 accumulatedDt_ += subDt;
1415 asImp_().updateTransportedQuantity(updateVec, realDt);
1418 resetTimeStepData_();
#define GET_PROP_VALUE(TypeTag, PropTagName)
Definition: propertysystemmacros.hh:282
#define GET_PROP(TypeTag, PropTagName)
Definition: propertysystemmacros.hh:281
#define GET_PROP_TYPE(TypeTag, PropTagName)
Definition: propertysystemmacros.hh:283
Define some often used mathematical functions.
Contains a class to exchange entries of a vector.
Determines the pressures and saturations of all fluid phases given the total mass of all components.
constexpr Scalar harmonicMean(Scalar x, Scalar y, Scalar wx=1.0, Scalar wy=1.0) noexcept
Calculate the (weighted) harmonic mean of two scalar values.
Definition: math.hh:68
void harmonicMeanMatrix(Dune::FieldMatrix< Scalar, m, n > &K, const Dune::FieldMatrix< Scalar, m, n > &Ki, const Dune::FieldMatrix< Scalar, m, n > &Kj)
Calculate the harmonic mean of a fixed-size matrix.
Definition: math.hh:108
T getParam(Args &&... args)
A free function to get a parameter from the parameter tree singleton.
Definition: parameters.hh:428
Dune::IteratorRange< typename MultiDomainGlue< DomainGridView, TargetGridView, DomainMapper, TargetMapper >::Intersections::const_iterator > intersections(const MultiDomainGlue< DomainGridView, TargetGridView, DomainMapper, TargetMapper > &glue)
Range generator to iterate with range-based for loops over all intersections as follows: for (const a...
Definition: glue.hh:62
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
Property tag EnableCapillarity
Returns whether capillarity is regarded.
Definition: porousmediumflow/2p2c/sequential/properties.hh:60
Property tag BoundaryMobility
Definition: porousmediumflow/2p2c/sequential/properties.hh:61
Property tag TransportModel
The type of the discretization of a transport model.
Definition: porousmediumflow/sequential/properties.hh:66
Property tag NumComponents
Number of components in the system.
Definition: porousmediumflow/sequential/properties.hh:70
Property tag NumEq
Number of equations in the system of PDEs.
Definition: porousmediumflow/sequential/properties.hh:68
Property tag RegulateBoundaryPermeability
A minimum permeability can be assigned via the runtime-Parameter SpatialParams.minBoundaryPermeabilit...
Definition: porousmediumflow/2p2c/sequential/properties.hh:63
Property tag NumPhases
Number of phases in the system.
Definition: porousmediumflow/sequential/properties.hh:69
Property tag SolutionTypes
Definition: porousmediumflow/sequential/properties.hh:58
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 porosity() noexcept
I/O name of porosity.
Definition: name.hh:139
void resize(const Vector &v, std::size_t size)
Definition: test_isvalid.cc:26
A data handle class to exchange entries of a vector.
Definition: vectorexchange.hh:40
Flash calculation routines for compositional sequential models.
Definition: compositionalflash.hh:51
static void saturationFlash2p2c(FluidState &fluidState, const Scalar &saturation, const PhaseVector &phasePressure, const Scalar &porosity, const Scalar &temperature)
2p2c flash for constant p & T if the saturation instead of concentration (feed mass fraction) is know...
Definition: compositionalflash.hh:224
static void concentrationFlash2p2c(FluidState &fluidState, const Scalar &Z0, const PhaseVector &phasePressure, const Scalar &porosity, const Scalar &temperature)
2p2c Flash for constant p & T if concentration (feed mass fraction) is given.
Definition: compositionalflash.hh:92
Compositional transport step in a Finite Volume discretization.
Definition: fvtransport.hh:60
virtual void update(const Scalar t, Scalar &dt, TransportSolutionType &updateVec, bool impet=false)
Calculate the update vector and determine timestep size This method calculates the update vector of ...
Definition: fvtransport.hh:363
void evalBoundary(GlobalPosition globalPosFace, const Intersection &intersection, FluidState &BCfluidState, PhaseVector &pressBound)
Evaluate the boundary conditions.
Definition: fvtransport.hh:1107
void deserializeEntity(std::istream &instream, const Element &element)
Function needed for restart option of the transport model: Read in.
Definition: fvtransport.hh:194
bool regulateBoundaryPermeability
Enables regulation of permeability in the direction of a Dirichlet Boundary Condition.
Definition: fvtransport.hh:310
void updateConcentrations(TransportSolutionType &updateVector, Scalar dt)
Updates the concentrations once an update is calculated.
Definition: fvtransport.hh:548
bool localTimeStepping_
Definition: fvtransport.hh:330
Scalar subCFLFactor_
Definition: fvtransport.hh:329
int restrictFluxInTransport_
Restriction of flux on new pressure field if direction reverses from the pressure equation.
Definition: fvtransport.hh:308
FVTransport2P2C(Problem &problem)
Constructs a FVTransport2P2C object Currently, the compositional transport scheme can not be applied ...
Definition: fvtransport.hh:273
const Scalar dtThreshold_
Threshold for sub-time-stepping routine.
Definition: fvtransport.hh:317
void initialize()
Set the initial values before the first pressure equation This method is called before first pressure...
Definition: fvtransport.hh:154
void addOutputVtkFields(MultiWriter &writer)
Write transport variables into the output files.
Definition: fvtransport.hh:171
void getSource(Scalar &update, const Element &element, CellData &cellDataI)
Definition: fvtransport.hh:236
Scalar & totalConcentration(int compIdx, int eIdxGlobal)
Return the the total concentration stored in the transport vector To get real cell values,...
Definition: fvtransport.hh:231
Problem & problem_
Definition: fvtransport.hh:301
void innerUpdate(TransportSolutionType &updateVec)
Definition: fvtransport.hh:1333
int averagedFaces_
number of faces were flux was restricted
Definition: fvtransport.hh:303
void updatedTargetDt_(Scalar &dt)
Definition: fvtransport.hh:1207
void resetTimeStepData_()
Definition: fvtransport.hh:323
Dune::FieldVector< Scalar, 2 > EntryType
Definition: fvtransport.hh:108
bool impet_
indicating if we are in an estimate (false) or real impet (true) step.
Definition: fvtransport.hh:302
void serializeEntity(std::ostream &outstream, const Element &element)
Function needed for restart option of the transport model: Write out.
Definition: fvtransport.hh:187
void getTransportedQuantity(TransportSolutionType &transportedQuantity)
Return the vector of the transported quantity For an immiscible IMPES scheme, this is the saturation....
Definition: fvtransport.hh:212
void getFlux(ComponentVector &fluxEntries, EntryType ×tepFlux, const Intersection &intersection, CellData &cellDataI)
Get flux at an interface between two cells The flux through is calculated according to the underlyin...
Definition: fvtransport.hh:579
Problem & problem()
Acess function for the current problem.
Definition: fvtransport.hh:122
std::vector< LocalTimesteppingData > timeStepData_
Stores data for sub-time-stepping.
Definition: fvtransport.hh:319
int verbosity_
Definition: fvtransport.hh:331
bool inPhysicalRange(DataEntry &entry)
Function to control the abort of the transport-sub-time-stepping depending on a physical parameter ra...
Definition: fvtransport.hh:246
TransportSolutionType totalConcentration_
private vector of transported primary variables
Definition: fvtransport.hh:300
Scalar accumulatedDt_
Current time-interval in sub-time-stepping routine.
Definition: fvtransport.hh:315
virtual ~FVTransport2P2C()
Definition: fvtransport.hh:295
void updateTransportedQuantity(TransportSolutionType &updateVector)
Updates the transported quantity once an update is calculated.
Definition: fvtransport.hh:519
void getFluxOnBoundary(ComponentVector &fluxEntries, EntryType ×tepFlux, const Intersection &intersection, const CellData &cellDataI)
Get flux on Boundary.
Definition: fvtransport.hh:884
static const int pressureType
gives kind of pressure used ( , , )
Definition: fvtransport.hh:306
bool switchNormals
Definition: fvtransport.hh:313
bool enableLocalTimeStepping()
Function to check if local time stepping is activated.
Definition: fvtransport.hh:260
Dune::FieldVector< Scalar, 2 > TimeStepFluxType
Definition: fvtransport.hh:109
Scalar minimalBoundaryPermeability
Minimal limit for the boundary permeability.
Definition: fvtransport.hh:312
Definition: fvtransport.hh:113
Dune::FieldVector< EntryType, 2 *dim > faceFluxes
Definition: fvtransport.hh:114
Dune::FieldVector< Scalar, 2 *dim > faceTargetDt
Definition: fvtransport.hh:115
LocalTimesteppingData()
Definition: fvtransport.hh:117
Scalar dt
Definition: fvtransport.hh:116
Defines the properties required for the sequential 2p2c models.