24#ifndef DUMUX_FVPRESSURE2P2C_HH
25#define DUMUX_FVPRESSURE2P2C_HH
28#include <dune/istl/bvector.hh>
29#include <dune/istl/operators.hh>
30#include <dune/istl/solvers.hh>
31#include <dune/istl/preconditioners.hh>
32#include <dune/common/float_cmp.hh>
92 dim = GridView::dimension, dimWorld = GridView::dimensionworld
96 pw = Indices::pressureW,
97 pn = Indices::pressureN,
98 pGlobal = Indices::pressureGlobal
102 wPhaseIdx = Indices::wPhaseIdx, nPhaseIdx = Indices::nPhaseIdx,
103 wCompIdx = Indices::wPhaseIdx, nCompIdx = Indices::nPhaseIdx,
104 contiWEqIdx = Indices::contiWEqIdx, contiNEqIdx = Indices::contiNEqIdx
121 using Element =
typename GridView::Traits::template Codim<0>::Entity;
122 using Intersection =
typename GridView::Intersection;
125 using GlobalPosition = Dune::FieldVector<Scalar, dimWorld>;
126 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
127 using PhaseVector = Dune::FieldVector<Scalar, getPropValue<TypeTag, Properties::NumPhases>()>;
148 void getSource(
EntryType& sourceEntry,
const Element& elementI,
const CellData& cellDataI,
const bool first);
150 void getStorage(
EntryType& storageEntry,
const Element& elementI,
const CellData& cellDataI,
const bool first);
152 void getFlux(
EntryType& entries,
const Intersection& intersection,
const CellData& cellDataI,
const bool first);
175 Dune::dinfo <<
" Warning: Regulating Boundary Permeability requires correct subface indices on reference elements!"
189 static constexpr int pressureType = getPropValue<TypeTag, Properties::PressureFormulation>();
192 Implementation &asImp_()
193 {
return *
static_cast<Implementation *
>(
this);}
196 const Implementation &asImp_()
const
197 {
return *
static_cast<const Implementation *
>(
this);}
211template<
class TypeTag>
213 const Element& elementI,
214 const CellData& cellDataI,
219 Scalar
volume = elementI.geometry().volume();
222 PrimaryVariables source(NAN);
223 problem().source(source, elementI);
227 source[contiWEqIdx] /= cellDataI.density(wPhaseIdx);
228 source[contiNEqIdx] /= cellDataI.density(nPhaseIdx);
233 const GlobalPosition& globalPos = elementI.geometry().center();
235 if (!cellDataI.hasVolumeDerivatives())
236 asImp_().volumeDerivatives(globalPos, elementI);
238 source[contiWEqIdx] *= cellDataI.dv(wCompIdx);
239 source[contiNEqIdx] *= cellDataI.dv(nCompIdx);
241 sourceEntry[rhs] =
volume * (source[contiWEqIdx] + source[contiNEqIdx]);
261template<
class TypeTag>
263 const Element& elementI,
264 const CellData& cellDataI,
269 int eIdxGlobalI = problem().variables().index(elementI);
270 Scalar
volume = elementI.geometry().volume();
273 Scalar timestep_ = problem().timeManager().timeStepSize();
276 if (!first && Dune::FloatCmp::ne<Scalar, Dune::FloatCmp::absolute>(timestep_, 0.0, 1.0e-30))
278 Scalar compress_term = cellDataI.dv_dp() / timestep_;
280 storageEntry[matrix] -= compress_term*
volume;
283 storageEntry[rhs] -= cellDataI.pressure(pressureType) * compress_term *
volume;
287 if (isnan(compress_term) || isinf(compress_term))
288 DUNE_THROW(Dune::MathError,
"Compressibility term leads to NAN matrix entry at index " << eIdxGlobalI);
290 if(!getPropValue<TypeTag, Properties::EnableCompressibility>())
291 DUNE_THROW(Dune::NotImplemented,
"Compressibility is switched off???");
296 if( problem().timeManager().episodeWillBeFinished()
297 || problem().timeManager().willBeFinished())
299 problem().variables().cellData(eIdxGlobalI).errorCorrection() = 0.;
305 problem().variables().cellData(eIdxGlobalI).volumeError() /= timestep_;
306 Scalar erri = fabs(cellDataI.volumeError());
307 Scalar x_lo = ErrorTermLowerBound_;
308 Scalar x_mi = ErrorTermUpperBound_;
309 Scalar fac = ErrorTermFactor_;
311 Scalar hifac = 1.-x_mi;
313 if ((erri*timestep_ > 5e-5) && (erri > x_lo * this->maxError_))
315 if (erri <= x_mi * this->maxError_)
317 problem().variables().cellData(eIdxGlobalI).errorCorrection() =
318 fac* (1-x_mi*(lofac-1)/(x_lo-x_mi) + (lofac-1)/(x_lo-x_mi)*erri/this->maxError_)
319 * cellDataI.volumeError() *
volume;
322 problem().variables().cellData(eIdxGlobalI).errorCorrection() =
323 fac * (1 + x_mi - hifac*x_mi/(1-x_mi) + (hifac/(1-x_mi)-1)*erri/this->maxError_)
324 * cellDataI.volumeError() *
volume;
327 problem().variables().cellData(eIdxGlobalI).errorCorrection() = 0.;
353template<
class TypeTag>
355 const Intersection& intersection,
356 const CellData& cellDataI,
360 auto elementI = intersection.inside();
361 int eIdxGlobalI = problem().variables().index(elementI);
364 const GlobalPosition& globalPos = elementI.geometry().center();
367 Scalar
volume = elementI.geometry().volume();
368 Scalar perimeter = cellDataI.perimeter();
372 const GlobalPosition& gravity_ = problem().gravity();
375 DimMatrix permeabilityI(problem().spatialParams().intrinsicPermeability(elementI));
378 Scalar fractionalWI=0, fractionalNWI=0;
381 fractionalWI = cellDataI.mobility(wPhaseIdx)
382 / (cellDataI.mobility(wPhaseIdx)+ cellDataI.mobility(nPhaseIdx));
383 fractionalNWI = cellDataI.mobility(nPhaseIdx)
384 / (cellDataI.mobility(wPhaseIdx)+ cellDataI.mobility(nPhaseIdx));
388 const GlobalPosition& unitOuterNormal = intersection.centerUnitOuterNormal();
391 Scalar faceArea = intersection.geometry().volume();
394 auto neighbor = intersection.outside();
395 int eIdxGlobalJ = problem().variables().index(neighbor);
396 CellData& cellDataJ = problem().variables().cellData(eIdxGlobalJ);
399 const GlobalPosition& globalPosNeighbor = neighbor.geometry().center();
402 GlobalPosition distVec = globalPosNeighbor - globalPos;
405 Scalar dist = distVec.two_norm();
407 GlobalPosition unitDistVec(distVec);
410 DimMatrix permeabilityJ
411 = problem().spatialParams().intrinsicPermeability(neighbor);
414 DimMatrix meanPermeability(0);
421 Scalar rhoMeanW = 0.5 * (cellDataI.density(wPhaseIdx) + cellDataJ.density(wPhaseIdx));
422 Scalar rhoMeanNW = 0.5 * (cellDataI.density(nPhaseIdx) + cellDataJ.density(nPhaseIdx));
425 Scalar potentialW = 0;
426 Scalar potentialNW = 0;
431 Scalar fractionalWJ = cellDataJ.mobility(wPhaseIdx)
432 / (cellDataJ.mobility(wPhaseIdx)+ cellDataJ.mobility(nPhaseIdx));
433 Scalar fractionalNWJ = cellDataJ.mobility(nPhaseIdx)
434 / (cellDataJ.mobility(wPhaseIdx)+ cellDataJ.mobility(nPhaseIdx));
437 Scalar lambda = (cellDataI.mobility(wPhaseIdx) + cellDataJ.mobility(wPhaseIdx)) * 0.5
438 + (cellDataI.mobility(nPhaseIdx) + cellDataJ.mobility(nPhaseIdx)) * 0.5;
440 entries[0] = fabs(lambda*faceArea*fabs(
permeability*unitOuterNormal)/(dist));
442 Scalar factor = (fractionalWI + fractionalWJ) * (rhoMeanW) * 0.5
443 + (fractionalNWI + fractionalNWJ) * (rhoMeanNW) * 0.5;
444 entries[1] = factor * lambda * faceArea * fabs(unitOuterNormal*
permeability)
445 * (gravity_ * unitDistVec);
450 if (!cellDataJ.hasVolumeDerivatives())
451 asImp_().volumeDerivatives(globalPosNeighbor, neighbor);
453 Scalar dv_dC1 = (cellDataJ.dv(wPhaseIdx)
454 + cellDataI.dv(wPhaseIdx)) / 2;
455 Scalar dv_dC2 = (cellDataJ.dv(nPhaseIdx)
456 + cellDataI.dv(nPhaseIdx)) / 2;
458 Scalar graddv_dC1 = (cellDataJ.dv(wPhaseIdx)
459 - cellDataI.dv(wPhaseIdx)) / dist;
460 Scalar graddv_dC2 = (cellDataJ.dv(nPhaseIdx)
461 - cellDataI.dv(nPhaseIdx)) / dist;
472 Scalar densityW = rhoMeanW;
473 Scalar densityNW = rhoMeanNW;
475 potentialW = (cellDataI.pressure(wPhaseIdx) - cellDataJ.pressure(wPhaseIdx))/dist;
476 potentialNW = (cellDataI.pressure(nPhaseIdx) - cellDataJ.pressure(nPhaseIdx))/dist;
478 potentialW += densityW * (unitDistVec * gravity_);
479 potentialNW += densityNW * (unitDistVec * gravity_);
482 Scalar lambdaW(0.), lambdaN(0.);
483 Scalar dV_w(0.), dV_n(0.);
484 Scalar gV_w(0.), gV_n(0.);
488 const CellData* upwindWCellData(0);
489 const CellData* upwindNCellData(0);
491 upwindWCellData = &cellDataI;
492 else if (potentialW < 0.)
493 upwindWCellData = &cellDataJ;
496 if(cellDataI.isUpwindCell(intersection.indexInInside(), contiWEqIdx))
497 upwindWCellData = &cellDataI;
498 else if(cellDataJ.isUpwindCell(intersection.indexInOutside(), contiWEqIdx))
499 upwindWCellData = &cellDataJ;
504 if (potentialNW > 0.)
505 upwindNCellData = &cellDataI;
506 else if (potentialNW < 0.)
507 upwindNCellData = &cellDataJ;
510 if(cellDataI.isUpwindCell(intersection.indexInInside(), contiNEqIdx))
511 upwindNCellData = &cellDataI;
512 else if(cellDataJ.isUpwindCell(intersection.indexInOutside(), contiNEqIdx))
513 upwindNCellData = &cellDataJ;
519 if(!upwindWCellData or (cellDataI.wasRefined() && cellDataJ.wasRefined() && elementI.father() == neighbor.father()))
521 if (cellDataI.wasRefined() && cellDataJ.wasRefined())
523 problem().variables().cellData(eIdxGlobalI).setUpwindCell(intersection.indexInInside(), contiWEqIdx,
false);
524 cellDataJ.setUpwindCell(intersection.indexInOutside(), contiWEqIdx,
false);
527 Scalar averagedMassFraction[2];
528 averagedMassFraction[wCompIdx]
529 =
harmonicMean(cellDataI.massFraction(wPhaseIdx, wCompIdx), cellDataJ.massFraction(wPhaseIdx, wCompIdx));
530 averagedMassFraction[nCompIdx]
531 =
harmonicMean(cellDataI.massFraction(wPhaseIdx, nCompIdx), cellDataJ.massFraction(wPhaseIdx, nCompIdx));
532 Scalar averageDensity =
harmonicMean(cellDataI.density(wPhaseIdx), cellDataJ.density(wPhaseIdx));
535 dV_w = dv_dC1 * averagedMassFraction[wCompIdx] + dv_dC2 * averagedMassFraction[nCompIdx];
536 dV_w *= averageDensity;
537 gV_w = graddv_dC1 * averagedMassFraction[wCompIdx] + graddv_dC2 * averagedMassFraction[nCompIdx];
538 gV_w *= averageDensity;
539 lambdaW =
harmonicMean(cellDataI.mobility(wPhaseIdx), cellDataJ.mobility(wPhaseIdx));
543 dV_w = (dv_dC1 * upwindWCellData->massFraction(wPhaseIdx, wCompIdx)
544 + dv_dC2 * upwindWCellData->massFraction(wPhaseIdx, nCompIdx));
545 lambdaW = upwindWCellData->mobility(wPhaseIdx);
546 gV_w = (graddv_dC1 * upwindWCellData->massFraction(wPhaseIdx, wCompIdx)
547 + graddv_dC2 * upwindWCellData->massFraction(wPhaseIdx, nCompIdx));
548 dV_w *= upwindWCellData->density(wPhaseIdx);
549 gV_w *= upwindWCellData->density(wPhaseIdx);
552 if(!upwindNCellData or (cellDataI.wasRefined() && cellDataJ.wasRefined()))
554 if (cellDataI.wasRefined() && cellDataJ.wasRefined())
556 problem().variables().cellData(eIdxGlobalI).setUpwindCell(intersection.indexInInside(), contiNEqIdx,
false);
557 cellDataJ.setUpwindCell(intersection.indexInOutside(), contiNEqIdx,
false);
559 Scalar averagedMassFraction[2];
560 averagedMassFraction[wCompIdx]
561 =
harmonicMean(cellDataI.massFraction(nPhaseIdx, wCompIdx), cellDataJ.massFraction(nPhaseIdx, wCompIdx));
562 averagedMassFraction[nCompIdx]
563 =
harmonicMean(cellDataI.massFraction(nPhaseIdx, nCompIdx), cellDataJ.massFraction(nPhaseIdx, nCompIdx));
564 Scalar averageDensity =
harmonicMean(cellDataI.density(nPhaseIdx), cellDataJ.density(nPhaseIdx));
567 dV_n = dv_dC1 * averagedMassFraction[wCompIdx] + dv_dC2 * averagedMassFraction[nCompIdx];
568 dV_n *= averageDensity;
569 gV_n = graddv_dC1 * averagedMassFraction[wCompIdx] + graddv_dC2 * averagedMassFraction[nCompIdx];
570 gV_n *= averageDensity;
571 lambdaN =
harmonicMean(cellDataI.mobility(nPhaseIdx), cellDataJ.mobility(nPhaseIdx));
575 dV_n = (dv_dC1 * upwindNCellData->massFraction(nPhaseIdx, wCompIdx)
576 + dv_dC2 * upwindNCellData->massFraction(nPhaseIdx, nCompIdx));
577 lambdaN = upwindNCellData->mobility(nPhaseIdx);
578 gV_n = (graddv_dC1 * upwindNCellData->massFraction(nPhaseIdx, wCompIdx)
579 + graddv_dC2 * upwindNCellData->massFraction(nPhaseIdx, nCompIdx));
580 dV_n *= upwindNCellData->density(nPhaseIdx);
581 gV_n *= upwindNCellData->density(nPhaseIdx);
585 entries[matrix] = faceArea * (lambdaW * dV_w + lambdaN * dV_n)
587 if(enableVolumeIntegral)
588 entries[matrix] -=
volume * faceArea / perimeter * (lambdaW * gV_w + lambdaN * gV_n)
592 entries[rhs] = faceArea * (densityW * lambdaW * dV_w + densityNW * lambdaN * dV_n);
594 if(enableVolumeIntegral)
595 entries[rhs] -=
volume * faceArea / perimeter * (densityW * lambdaW * gV_w + densityNW * lambdaN * gV_n)
597 entries[rhs] *= (gravity_ * unitDistVec);
600 Scalar pCGradient = (cellDataI.capillaryPressure() - cellDataJ.capillaryPressure()) / dist;
603 switch (pressureType)
608 entries[rhs] += lambdaN * dV_n * fabs(
permeability * unitOuterNormal) * pCGradient * faceArea;
609 if(enableVolumeIntegral)
610 entries[rhs]-= lambdaN * gV_n * (
permeability * unitDistVec) * pCGradient *
volume * faceArea / perimeter;
616 entries[rhs] -= lambdaW * dV_w * fabs(
permeability * unitOuterNormal) * pCGradient * faceArea;
617 if(enableVolumeIntegral)
618 entries[rhs]+= lambdaW * gV_w * (
permeability * unitDistVec) * pCGradient *
volume * faceArea / perimeter;
644template<
class TypeTag>
646 const Intersection& intersection,
const CellData& cellDataI,
const bool first)
650 auto elementI = intersection.inside();
651 const GlobalPosition& globalPos = elementI.geometry().center();
654 const GlobalPosition& unitOuterNormal = intersection.centerUnitOuterNormal();
656 Scalar faceArea = intersection.geometry().volume();
659 Scalar dv_dC1 = cellDataI.dv(wCompIdx);
660 Scalar dv_dC2 = cellDataI.dv(nCompIdx);
663 const GlobalPosition& globalPosFace = intersection.geometry().center();
666 GlobalPosition distVec(globalPosFace - globalPos);
667 Scalar dist = distVec.two_norm();
668 GlobalPosition unitDistVec(distVec);
672 BoundaryTypes bcType;
673 problem().boundaryTypes(bcType, intersection);
676 PhaseVector pressBC(0.);
679 const auto fluidMatrixInteraction = problem().spatialParams().fluidMatrixInteractionAtPos(elementI.geometry().center());
682 if (bcType.isDirichlet(Indices::pressureEqIdx))
685 DimMatrix permeabilityI(problem().spatialParams().intrinsicPermeability(elementI));
687 if(regulateBoundaryPermeability)
689 int axis = intersection.indexInInside() / 2;
690 if(permeabilityI[axis][axis] < minimalBoundaryPermeability)
691 permeabilityI[axis][axis] = minimalBoundaryPermeability;
693 const GlobalPosition& gravity_ = problem().gravity();
700 FluidState BCfluidState;
703 PrimaryVariables primaryVariablesOnBoundary(NAN);
704 problem().dirichlet(primaryVariablesOnBoundary, intersection);
709 Scalar fractionalWI=0, fractionalNWI=0;
710 fractionalWI = cellDataI.mobility(wPhaseIdx)
711 / (cellDataI.mobility(wPhaseIdx)+ cellDataI.mobility(nPhaseIdx));
712 fractionalNWI = cellDataI.mobility(nPhaseIdx)
713 / (cellDataI.mobility(wPhaseIdx)+ cellDataI.mobility(nPhaseIdx));
715 Scalar lambda = cellDataI.mobility(wPhaseIdx)+cellDataI.mobility(nPhaseIdx);
716 entries[matrix] += lambda * faceArea * fabs(
permeability * unitOuterNormal) / (dist);
717 Scalar pressBoundary = primaryVariablesOnBoundary[Indices::pressureEqIdx];
718 entries[rhs] += lambda * faceArea * pressBoundary * fabs(
permeability * unitOuterNormal) / (dist);
719 Scalar rightentry = (fractionalWI * cellDataI.density(wPhaseIdx)
720 + fractionalNWI * cellDataI.density(nPhaseIdx))
721 * lambda * faceArea * fabs(unitOuterNormal *
permeability)
722 * ( unitDistVec * gravity_);
723 entries[rhs] -= rightentry;
728 problem().transportModel().evalBoundary(globalPosFace,
732 pcBound = pressBC[nPhaseIdx] - pressBC[wPhaseIdx];
735 Scalar lambdaWBound = 0.;
736 Scalar lambdaNWBound = 0.;
738 Scalar densityWBound =
740 Scalar densityNWBound =
742 Scalar viscosityWBound =
744 Scalar viscosityNWBound =
748 if(getPropValue<TypeTag, Properties::BoundaryMobility>() == Indices::satDependent)
750 lambdaWBound = BCfluidState.saturation(wPhaseIdx)
752 lambdaNWBound = BCfluidState.saturation(nPhaseIdx)
755 else if(getPropValue<TypeTag, Properties::BoundaryMobility>() == Indices::permDependent)
758 = fluidMatrixInteraction.krw(BCfluidState.saturation(wPhaseIdx)) / viscosityWBound;
760 = fluidMatrixInteraction.krn(BCfluidState.saturation(wPhaseIdx)) / viscosityNWBound;
763 Scalar rhoMeanW = 0.5 * (cellDataI.density(wPhaseIdx) + densityWBound);
764 Scalar rhoMeanNW = 0.5 * (cellDataI.density(nPhaseIdx) + densityNWBound);
766 Scalar potentialW = 0;
767 Scalar potentialNW = 0;
779 Scalar densityW=rhoMeanW;
780 Scalar densityNW=rhoMeanNW;
783 potentialW = (cellDataI.pressure(wPhaseIdx) - pressBC[wPhaseIdx])/dist;
784 potentialNW = (cellDataI.pressure(nPhaseIdx) - pressBC[nPhaseIdx])/dist;
786 potentialW += densityW * (unitDistVec * gravity_);
787 potentialNW += densityNW * (unitDistVec * gravity_);
791 Scalar lambdaW, lambdaNW;
792 Scalar densityW(0.), densityNW(0.);
796 if (potentialW >= 0.)
798 densityW = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialW, 0.0, 1.0e-30)) ? rhoMeanW : cellDataI.density(wPhaseIdx);
799 dV_w = (dv_dC1 * cellDataI.massFraction(wPhaseIdx, wCompIdx)
800 + dv_dC2 * cellDataI.massFraction(wPhaseIdx, nCompIdx));
802 lambdaW = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialW, 0.0, 1.0e-30)) ? 0.5 * (cellDataI.mobility(wPhaseIdx) + lambdaWBound)
803 : cellDataI.mobility(wPhaseIdx);
807 densityW = densityWBound;
808 dV_w = (dv_dC1 * BCfluidState.massFraction(wPhaseIdx, wCompIdx)
809 + dv_dC2 * BCfluidState.massFraction(wPhaseIdx, nCompIdx));
811 lambdaW = lambdaWBound;
813 if (potentialNW >= 0.)
815 densityNW = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialNW, 0.0, 1.0e-30)) ? rhoMeanNW : cellDataI.density(nPhaseIdx);
816 dV_n = (dv_dC1 * cellDataI.massFraction(nPhaseIdx, wCompIdx)
817 + dv_dC2 * cellDataI.massFraction(nPhaseIdx, nCompIdx));
819 lambdaNW = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialNW, 0.0, 1.0e-30)) ? 0.5 * (cellDataI.mobility(nPhaseIdx) + lambdaNWBound)
820 : cellDataI.mobility(nPhaseIdx);
824 densityNW = densityNWBound;
825 dV_n = (dv_dC1 * BCfluidState.massFraction(nPhaseIdx, wCompIdx)
826 + dv_dC2 * BCfluidState.massFraction(nPhaseIdx, nCompIdx));
828 lambdaNW = lambdaNWBound;
832 Scalar entry = (lambdaW * dV_w + lambdaNW * dV_n)
833 * (fabs(unitOuterNormal *
permeability) / dist) * faceArea;
836 Scalar rightEntry = (lambdaW * densityW * dV_w + lambdaNW * densityNW * dV_n)
837 * fabs(unitOuterNormal *
permeability) * (gravity_ * unitDistVec) * faceArea ;
841 Scalar pCGradient = (cellDataI.capillaryPressure() - pcBound) / dist;
842 switch (pressureType)
847 rightEntry += lambdaNW * dV_n * pCGradient * fabs(unitOuterNormal *
permeability) * faceArea;
853 rightEntry -= lambdaW * dV_w * pCGradient * fabs(unitOuterNormal *
permeability) * faceArea;
860 entries[matrix] += entry;
861 entries[rhs] += entry * primaryVariablesOnBoundary[Indices::pressureEqIdx];
862 entries[rhs] -= rightEntry;
869 else if(bcType.isNeumann(Indices::pressureEqIdx))
871 PrimaryVariables J(NAN);
872 problem().neumann(J, intersection);
875 J[contiWEqIdx] /= cellDataI.density(wPhaseIdx);
876 J[contiNEqIdx] /= cellDataI.density(nPhaseIdx);
880 J[contiWEqIdx] *= dv_dC1;
881 J[contiNEqIdx] *= dv_dC2;
884 entries[rhs] -= (J[contiWEqIdx] + J[contiNEqIdx]) * faceArea;
887 DUNE_THROW(Dune::NotImplemented,
"Boundary Condition neither Dirichlet nor Neumann!");
902template<
class TypeTag>
906 GlobalPosition globalPos = element.geometry().center();
909 int eIdxGlobal = problem().variables().index(element);
910 CellData& cellData = problem().variables().cellData(eIdxGlobal);
913 FluidState& fluidState = cellData.manipulateFluidState();
915 Scalar temperature_ = problem().temperatureAtPos(globalPos);
922 Scalar Z0 = cellData.massConcentration(wCompIdx)
923 / (cellData.massConcentration(wCompIdx)
924 + cellData.massConcentration(nCompIdx));
927 if(Z0 < 0. || Z0 > 1.)
929 Dune::dgrave <<
"Feed mass fraction unphysical: Z0 = " << Z0
930 <<
" at global Idx " << eIdxGlobal
931 <<
" , because totalConcentration(wCompIdx) = "
932 << cellData.totalConcentration(wCompIdx)
933 <<
" and totalConcentration(nCompIdx) = "
934 << cellData.totalConcentration(nCompIdx)<< std::endl;
938 cellData.setTotalConcentration(wCompIdx, 0.);
939 problem().transportModel().totalConcentration(wCompIdx, eIdxGlobal) = 0.;
940 Dune::dgrave <<
"Regularize totalConcentration(wCompIdx) = "
941 << cellData.totalConcentration(wCompIdx)<< std::endl;
946 cellData.setTotalConcentration(nCompIdx, 0.);
947 problem().transportModel().totalConcentration(nCompIdx,eIdxGlobal) = 0.;
948 Dune::dgrave <<
"Regularize totalConcentration(eIdxGlobal, nCompIdx) = "
949 << cellData.totalConcentration(nCompIdx)<< std::endl;
956 const auto fluidMatrixInteraction = problem().spatialParams().fluidMatrixInteractionAtPos(element.geometry().center());
957 if(getPropValue<TypeTag, Properties::EnableCapillarity>())
959 unsigned int maxiter = 6;
960 Scalar pc = cellData.capillaryPressure();
962 for (
unsigned int iter = 0; iter < maxiter; iter++)
964 switch (pressureType)
968 pressure[wPhaseIdx] = asImp_().pressure(eIdxGlobal);
969 pressure[nPhaseIdx] = asImp_().pressure(eIdxGlobal) + pc;
974 pressure[wPhaseIdx] = asImp_().pressure(eIdxGlobal) - pc;
975 pressure[nPhaseIdx] = asImp_().pressure(eIdxGlobal);
985 pc = fluidMatrixInteraction.pc(fluidState.saturation(wPhaseIdx));
987 if (fabs(oldPc-pc)<10 && iter != 0)
990 if (iter++ == maxiter)
991 Dune::dinfo << iter <<
"times iteration of pc was applied at Idx " << eIdxGlobal
992 <<
", pc delta still " << fabs(oldPc-pc) << std::endl;
997 pressure[wPhaseIdx] =
pressure[nPhaseIdx] = asImp_().pressure()[eIdxGlobal];
1002 cellData.setMobility(wPhaseIdx, fluidMatrixInteraction.krw(fluidState.saturation(wPhaseIdx))
1003 / cellData.viscosity(wPhaseIdx));
1004 cellData.setMobility(nPhaseIdx, fluidMatrixInteraction.krn(fluidState.saturation(wPhaseIdx))
1005 / cellData.viscosity(nPhaseIdx));
1008 Scalar sumConc = (cellData.totalConcentration(wCompIdx)
1009 + cellData.totalConcentration(nCompIdx));
1010 Scalar massw = sumConc * fluidState.phaseMassFraction(wPhaseIdx);
1011 Scalar massn = sumConc * fluidState.phaseMassFraction(nPhaseIdx);
1013 if (Dune::FloatCmp::eq<Scalar>((cellData.density(wPhaseIdx)*cellData.density(nPhaseIdx)), 0))
1014 DUNE_THROW(Dune::MathError,
"Sequential2p2c::postProcessUpdate: try to divide by 0 density");
1015 Scalar vol = massw / cellData.density(wPhaseIdx) + massn / cellData.density(nPhaseIdx);
1016 if (Dune::FloatCmp::ne<Scalar, Dune::FloatCmp::absolute>(problem().timeManager().timeStepSize(), 0.0, 1.0e-30))
1018 cellData.volumeError()=(vol - problem().spatialParams().porosity(element));
1021 if (isnan(cellData.volumeError()))
1023 DUNE_THROW(Dune::MathError,
"Sequential2p2c::postProcessUpdate:\n"
1024 <<
"volErr[" << eIdxGlobal <<
"] isnan: vol = " << vol
1025 <<
", massw = " << massw <<
", rho_l = " << cellData.density(wPhaseIdx)
1026 <<
", massn = " << massn <<
", rho_g = " << cellData.density(nPhaseIdx)
1027 <<
", poro = " << problem().spatialParams().
porosity(element)
1028 <<
", dt = " << problem().timeManager().timeStepSize());
1032 cellData.volumeError()=0.;
Define some often used mathematical functions.
Simplifies writing multi-file VTK datasets.
Determines the pressures and saturations of all fluid phases given the total mass of all components.
Base class for compositional pressure equations.
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:69
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:109
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:150
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 pressure(int phaseIdx) noexcept
I/O name of pressure for multiphase systems.
Definition: name.hh:34
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
Scalar volume(Shape shape, Scalar inscribedRadius)
Returns the volume of a given geometry based on the inscribed radius.
Definition: poreproperties.hh:73
Flash calculation routines for compositional sequential models.
Definition: compositionalflash.hh:51
static void concentrationFlash2p2c(FluidState &fluidState, const Scalar Z0, const PhaseVector &phasePressure, const Scalar temperature)
Definition: compositionalflash.hh:91
The finite volume model for the solution of the compositional pressure equation.
Definition: fvpressure.hh:73
bool enableVolumeIntegral
Enables the volume integral of the pressure equation.
Definition: fvpressure.hh:182
Scalar ErrorTermLowerBound_
Handling of error term: lower bound for error dampening.
Definition: fvpressure.hh:186
void getSource(EntryType &sourceEntry, const Element &elementI, const CellData &cellDataI, const bool first)
Assembles the source term.
Definition: fvpressure.hh:212
Scalar ErrorTermFactor_
Handling of error term: relaxation factor.
Definition: fvpressure.hh:185
void getFlux(EntryType &entries, const Intersection &intersection, const CellData &cellDataI, const bool first)
Get flux at an interface between two cells.
Definition: fvpressure.hh:354
void getFluxOnBoundary(EntryType &entries, const Intersection &intersection, const CellData &cellDataI, const bool first)
Get flux on Boundary.
Definition: fvpressure.hh:645
const Problem & problem() const
Definition: fvpressure.hh:141
Problem & problem_
Definition: fvpressure.hh:181
static constexpr int pressureType
gives kind of pressure used ( , , )
Definition: fvpressure.hh:189
bool regulateBoundaryPermeability
Enables regulation of permeability in the direction of a Dirichlet Boundary Condition.
Definition: fvpressure.hh:183
Problem & problem()
Definition: fvpressure.hh:137
Scalar minimalBoundaryPermeability
Minimal limit for the boundary permeability.
Definition: fvpressure.hh:184
Dune::FieldVector< Scalar, 2 > EntryType
Definition: fvpressure.hh:135
void getStorage(EntryType &storageEntry, const Element &elementI, const CellData &cellDataI, const bool first)
Assembles the storage term.
Definition: fvpressure.hh:262
Scalar ErrorTermUpperBound_
Definition: fvpressure.hh:187
FVPressure2P2C(Problem &problem)
Constructs a FVPressure2P2C object.
Definition: fvpressure.hh:163
void updateMaterialLawsInElement(const Element &elementI, bool postTimeStep)
Updates secondary variables of one cell.
Definition: fvpressure.hh:903
The finite volume model for the solution of the compositional pressure equation.
Definition: fvpressurecompositional.hh:67
Defines the properties required for the sequential 2p2c models.