22#ifndef DUMUX_FVMPFAL2DPRESSURE2P_HH
23#define DUMUX_FVMPFAL2DPRESSURE2P_HH
72template<
class TypeTag>
80 dim = GridView::dimension, dimWorld = GridView::dimensionworld
98 using PrimaryVariables =
typename SolutionTypes::PrimaryVariables;
99 using ScalarSolutionType =
typename SolutionTypes::ScalarSolution;
105 pw = Indices::pressureW,
106 pn = Indices::pressureNw,
107 sw = Indices::saturationW,
108 sn = Indices::saturationNw
112 wPhaseIdx = Indices::wPhaseIdx,
113 nPhaseIdx = Indices::nPhaseIdx,
114 pressureIdx = Indices::pressureIdx,
115 saturationIdx = Indices::saturationIdx,
116 pressEqIdx = Indices::pressureEqIdx,
117 satEqIdx = Indices::satEqIdx,
118 numPhases = getPropValue<TypeTag, Properties::NumPhases>()
125 dirichletDirichlet = 1,
126 dirichletNeumann = 2,
130 using Element =
typename GridView::Traits::template Codim<0>::Entity;
131 using IntersectionIterator =
typename GridView::IntersectionIterator;
132 using Intersection =
typename GridView::Intersection;
134 using LocalPosition = Dune::FieldVector<Scalar, dim>;
135 using GravityVector = Dune::FieldVector<Scalar, dimWorld>;
137 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
139 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
141 using DimVector = Dune::FieldVector<Scalar, dim>;
154 using GlobalInteractionVolumeVector = std::vector<InteractionVolume>;
155 using InnerBoundaryVolumeFaces = std::vector<Dune::FieldVector<bool, 2*dim> >;
158 Intersection getNextIntersection_(
const Element&,
const IntersectionIterator&);
162 void initializeMatrix();
164 void storeInteractionVolumeInfo();
187 storeInteractionVolumeInfo();
197 const auto element = *problem_.gridView().template begin<0>();
198 FluidState fluidState;
199 fluidState.setPressure(wPhaseIdx, problem_.referencePressure(element));
200 fluidState.setPressure(nPhaseIdx, problem_.referencePressure(element));
201 fluidState.setTemperature(problem_.temperature(element));
202 fluidState.setSaturation(wPhaseIdx, 1.);
203 fluidState.setSaturation(nPhaseIdx, 0.);
214 storeInteractionVolumeInfo();
229 for (
const auto& element : elements(problem_.gridView()))
242 int eIdxGlobal = problem_.variables().index(element);
243 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
245 switch (pressureType_)
249 Scalar potW = this->
pressure()[eIdxGlobal];
251 Scalar gravityDiff = (problem_.bBoxMax() - element.geometry().center()) * gravity_;
252 Scalar potPc = cellData.capillaryPressure()
253 + gravityDiff * (density_[nPhaseIdx] - density_[wPhaseIdx]);
255 cellData.setPotential(wPhaseIdx, potW);
256 cellData.setPotential(nPhaseIdx, potW + potPc);
258 Scalar pressW = potW - gravityDiff * density_[wPhaseIdx];
260 cellData.setPressure(wPhaseIdx, pressW);
261 cellData.setPressure(nPhaseIdx, pressW + cellData.capillaryPressure());
267 Scalar potNw = this->
pressure()[eIdxGlobal];
269 Scalar gravityDiff = (problem_.bBoxMax() - element.geometry().center()) * gravity_;
270 Scalar potPc = cellData.capillaryPressure()
271 + gravityDiff * (density_[nPhaseIdx] - density_[wPhaseIdx]);
273 cellData.setPotential(nPhaseIdx, potNw);
274 cellData.setPotential(wPhaseIdx, potNw - potPc);
276 Scalar pressNw = potNw - gravityDiff * density_[nPhaseIdx];
278 cellData.setPressure(wPhaseIdx, pressNw - cellData.capillaryPressure());
279 cellData.setPressure(nPhaseIdx, pressNw);
284 cellData.fluxData().resetVelocity();
294 timeStep_ = problem_.timeManager().timeStepSize();
296 int size = problem_.gridView().size(0);
297 for (
int i = 0; i < size; i++)
301 switch (saturationType_)
304 sat = problem_.variables().cellData(i).saturation(wPhaseIdx);
307 sat = problem_.variables().cellData(i).saturation(nPhaseIdx);
312 maxError_ = max(maxError_, (sat - 1.0) / timeStep_);
316 maxError_ = max(maxError_, (-sat) / timeStep_);
338 template<
class MultiWriter>
341 int size = problem_.gridView().size(0);
342 ScalarSolutionType *potential = writer.allocateManagedBuffer(size);
346 if (pressureType_ == pw)
348 writer.attachCellData(*potential,
"wetting potential");
351 if (pressureType_ == pn)
353 writer.attachCellData(*potential,
"nonwetting potential");
356 if (vtkOutputLevel_ > 0)
358 ScalarSolutionType *
pressure = writer.allocateManagedBuffer(size);
359 ScalarSolutionType *pressureSecond = writer.allocateManagedBuffer(size);
360 ScalarSolutionType *potentialSecond = writer.allocateManagedBuffer(size);
361 ScalarSolutionType *pc = writer.allocateManagedBuffer(size);
363 for (
const auto& element : elements(problem_.gridView()))
365 int idx = problem_.variables().index(element);
366 CellData& cellData = problem_.variables().cellData(idx);
368 (*pc)[idx] = cellData.capillaryPressure();
370 if (pressureType_ == pw)
372 (*pressure)[idx] = cellData.pressure(wPhaseIdx);
373 (*potentialSecond)[idx] = cellData.potential(nPhaseIdx);
374 (*pressureSecond)[idx] = cellData.pressure(nPhaseIdx);
377 if (pressureType_ == pn)
379 (*pressure)[idx] = cellData.pressure(nPhaseIdx);
380 (*potentialSecond)[idx] = cellData.potential(wPhaseIdx);
381 (*pressureSecond)[idx] = cellData.pressure(wPhaseIdx);
385 if (pressureType_ == pw)
387 writer.attachCellData(*
pressure,
"wetting pressure");
388 writer.attachCellData(*pressureSecond,
"nonwetting pressure");
389 writer.attachCellData(*potentialSecond,
"nonwetting potential");
392 if (pressureType_ == pn)
394 writer.attachCellData(*
pressure,
"nonwetting pressure");
395 writer.attachCellData(*pressureSecond,
"wetting pressure");
396 writer.attachCellData(*potentialSecond,
"wetting potential");
399 writer.attachCellData(*pc,
"capillary pressure");
411 ParentType(problem), problem_(problem), transmissibilityCalculator_(problem),
412 gravity_(problem.gravity()),
413 maxError_(0.), timeStep_(1.)
415 if (pressureType_ != pw && pressureType_ != pn)
417 DUNE_THROW(Dune::NotImplemented,
"Pressure type not supported!");
419 if (saturationType_ != sw && saturationType_ != sn)
421 DUNE_THROW(Dune::NotImplemented,
"Saturation type not supported!");
423 if (getPropValue<TypeTag, Properties::EnableCompressibility>())
425 DUNE_THROW(Dune::NotImplemented,
"Compressibility not supported!");
429 DUNE_THROW(Dune::NotImplemented,
"Dimension not supported!");
432 ErrorTermFactor_ = getParam<Scalar>(
"Impet.ErrorTermFactor");
433 ErrorTermLowerBound_ = getParam<Scalar>(
"Impet.ErrorTermLowerBound");
434 ErrorTermUpperBound_ = getParam<Scalar>(
"Impet.ErrorTermUpperBound");
436 density_[wPhaseIdx] = 0.;
437 density_[nPhaseIdx] = 0.;
438 viscosity_[wPhaseIdx] = 0.;
439 viscosity_[nPhaseIdx] = 0.;
441 vtkOutputLevel_ = getParam<int>(
"Vtk.OutputLevel");
454 const GravityVector& gravity_;
458 Scalar ErrorTermFactor_;
459 Scalar ErrorTermLowerBound_;
460 Scalar ErrorTermUpperBound_;
462 Scalar density_[numPhases];
463 Scalar viscosity_[numPhases];
467 static constexpr Scalar threshold_ = 1e-15;
469 static const int pressureType_ = getPropValue<TypeTag, Properties::PressureFormulation>();
471 static const int saturationType_ = getPropValue<TypeTag, Properties::SaturationFormulation>();
473 static const int velocityType_ = getPropValue<TypeTag, Properties::VelocityFormulation>();
494 Scalar evaluateErrorTerm_(CellData& cellData)
499 switch (saturationType_)
502 sat = cellData.saturation(wPhaseIdx);
505 sat = cellData.saturation(nPhaseIdx);
509 Scalar error = (sat > 1.0) ? sat - 1.0 : 0.0;
517 Scalar errorAbs = abs(error);
519 if ((errorAbs * timeStep_ > 1e-6) && (errorAbs > ErrorTermLowerBound_ * maxError_)
520 && (!problem_.timeManager().willBeFinished()))
522 return ErrorTermFactor_ * error;
530template<
class TypeTag>
531typename FvMpfaL2dPressure2p<TypeTag>::Intersection
532 FvMpfaL2dPressure2p<TypeTag>::getNextIntersection_(
const Element& element,
533 const IntersectionIterator& isIt)
535 auto isItBegin = problem_.gridView().ibegin(element);
536 const auto isEndIt = problem_.gridView().iend(element);
538 auto tempIsIt = isIt;
539 auto nextIsIt = ++tempIsIt;
542 switch (getPropValue<TypeTag, Properties::GridImplementation>())
545 case GridTypeIndices::yaspGrid:
547 if (nextIsIt == isEndIt)
549 nextIsIt = isItBegin;
553 nextIsIt = ++tempIsIt;
555 if (nextIsIt == isEndIt)
557 auto tempIsItBegin = isItBegin;
558 nextIsIt = ++tempIsItBegin;
565 case GridTypeIndices::aluGrid:
566 case GridTypeIndices::ugGrid:
568 if (nextIsIt == isEndIt)
569 nextIsIt = isItBegin;
575 DUNE_THROW(Dune::NotImplemented,
"GridType can not be used with MPFAO implementation!");
584template<
class TypeTag>
585void FvMpfaL2dPressure2p<TypeTag>::initializeMatrix()
588 for (
const auto& element : elements(problem_.gridView()))
591 int eIdxGlobalI = problem_.variables().index(element);
597 const auto isEndIt = problem_.gridView().iend(element);
598 for (
auto isIt = problem_.gridView().ibegin(element); isIt != isEndIt; ++isIt)
600 const auto& intersection = *isIt;
601 auto nextIntersection = getNextIntersection_(element, isIt);
603 if (intersection.neighbor())
607 if (intersection.neighbor() && nextIntersection.neighbor())
609 for (
const auto& innerIntersection
610 : intersections(problem_.gridView(), intersection.outside()))
611 for (
const auto& innerNextIntersection
612 : intersections(problem_.gridView(), nextIntersection.outside()))
614 if (innerIntersection.neighbor() && innerNextIntersection.neighbor())
616 if (innerIntersection.outside() == innerNextIntersection.outside()
617 && innerIntersection.outside() != intersection.inside())
627 this->A_.setrowsize(eIdxGlobalI, rowSize);
632 this->A_.endrowsizes();
635 for (
const auto& element : elements(problem_.gridView()))
638 int eIdxGlobalI = problem_.variables().index(element);
641 this->A_.addindex(eIdxGlobalI, eIdxGlobalI);
644 const auto isEndIt = problem_.gridView().iend(element);
645 for (
auto isIt = problem_.gridView().ibegin(element); isIt != isEndIt; ++isIt)
647 const auto& intersection = *isIt;
648 auto nextIntersection = getNextIntersection_(element, isIt);
650 if (intersection.neighbor())
653 int eIdxGlobalJ = problem_.variables().index(intersection.outside());
657 this->A_.addindex(eIdxGlobalI, eIdxGlobalJ);
660 if (intersection.neighbor() && nextIntersection.neighbor())
662 for (
const auto& innerIntersection
663 : intersections(problem_.gridView(), intersection.outside()))
664 for (
const auto& innerNextIntersection
665 : intersections(problem_.gridView(), nextIntersection.outside()))
667 if (innerIntersection.neighbor() && innerNextIntersection.neighbor())
669 auto innerOutside = innerIntersection.outside();
671 if (innerOutside == innerNextIntersection.outside()
672 && innerOutside != intersection.inside())
674 int eIdxGlobalJ = problem_.variables().index(innerOutside);
676 this->A_.addindex(eIdxGlobalI, eIdxGlobalJ);
685 this->A_.endindices();
712template<
class TypeTag>
713void FvMpfaL2dPressure2p<TypeTag>::storeInteractionVolumeInfo()
715 BoundaryTypes bcType;
718 for (
const auto& element : elements(problem_.gridView()))
722 int eIdxGlobal1 = problem_.variables().index(element);
724 const auto isIt12End = problem_.gridView().iend(element);
725 for (
auto isIt12 = problem_.gridView().ibegin(element); isIt12 != isIt12End; ++isIt12)
727 const auto& intersection12 = *isIt12;
728 auto intersection14 = getNextIntersection_(element, isIt12);
730 int indexInInside12 = intersection12.indexInInside();
731 int indexInInside14 = intersection14.indexInInside();
735 const auto refElement = referenceElement(element);
737 GlobalPosition corner1234(0);
739 int globalVertIdx1234 = 0;
742 for (
int i = 0; i < intersection12.geometry().corners(); ++i)
744 bool finished =
false;
746 int localVertIdx12corner = refElement.subEntity(indexInInside12, 1, i, dim);
748 int globalVertIdx12corner = problem_.variables().vertexMapper().subIndex(element, localVertIdx12corner, dim);
750 for (
int j = 0; j < intersection14.geometry().corners(); ++j)
752 int localVertIdx14corner = refElement.subEntity(indexInInside14, 1, j, dim);
754 int globalVertIdx14corner = problem_.variables().vertexMapper().subIndex(element, localVertIdx14corner, dim);
756 if (globalVertIdx12corner == globalVertIdx14corner)
758 corner1234 = element.geometry().corner(localVertIdx12corner);
760 globalVertIdx1234 = globalVertIdx12corner;
773 if (interactionVolumes_[globalVertIdx1234].isStored())
779 interactionVolumes_[globalVertIdx1234].setStored();
782 interactionVolumes_[globalVertIdx1234].setCenterPosition(corner1234);
785 interactionVolumes_[globalVertIdx1234].setSubVolumeElement(element, 0);
786 interactionVolumes_[globalVertIdx1234].setIndexOnElement(intersection12.indexInInside(), 0, 0);
787 interactionVolumes_[globalVertIdx1234].setIndexOnElement(intersection14.indexInInside(), 0, 1);
790 const GlobalPosition& globalPosFace12 = intersection12.geometry().center();
793 Scalar faceVol12 = intersection12.geometry().volume() / 2.0;
796 DimVector unitOuterNormal12 = intersection12.centerUnitOuterNormal();
799 const GlobalPosition& globalPosFace41 = intersection14.geometry().center();
802 Scalar faceVol41 = intersection14.geometry().volume() / 2.0;
805 DimVector unitOuterNormal14 = intersection14.centerUnitOuterNormal();
807 interactionVolumes_[globalVertIdx1234].setNormal(unitOuterNormal12, 0, 0);
808 interactionVolumes_[globalVertIdx1234].setNormal(unitOuterNormal14, 0, 1);
810 unitOuterNormal14 *= -1;
811 unitOuterNormal12 *= -1;
812 interactionVolumes_[globalVertIdx1234].setFaceArea(faceVol12, 0, 0);
813 interactionVolumes_[globalVertIdx1234].setFaceArea(faceVol41, 0, 1);
814 interactionVolumes_[globalVertIdx1234].setFacePosition(globalPosFace12, 0, 0);
815 interactionVolumes_[globalVertIdx1234].setFacePosition(globalPosFace41, 0, 1);
818 if (intersection12.neighbor())
821 auto element2 = intersection12.outside();
823 int eIdxGlobal2 = problem_.variables().index(element2);
826 interactionVolumes_[globalVertIdx1234].setSubVolumeElement(element2, 1);
827 interactionVolumes_[globalVertIdx1234].setIndexOnElement(intersection12.indexInOutside(), 1, 1);
828 interactionVolumes_[globalVertIdx1234].setNormal(unitOuterNormal12, 1, 1);
829 interactionVolumes_[globalVertIdx1234].setFaceArea(faceVol12, 1, 1);
830 interactionVolumes_[globalVertIdx1234].setFacePosition(globalPosFace12, 1, 1);
833 if (intersection14.neighbor())
837 auto element4 = intersection14.outside();
840 interactionVolumes_[globalVertIdx1234].setSubVolumeElement(element4, 3);
841 interactionVolumes_[globalVertIdx1234].setIndexOnElement(intersection14.indexInOutside(), 3, 0);
843 interactionVolumes_[globalVertIdx1234].setNormal(unitOuterNormal14, 3, 0);
844 interactionVolumes_[globalVertIdx1234].setFaceArea(faceVol41, 3, 0);
845 interactionVolumes_[globalVertIdx1234].setFacePosition(globalPosFace41, 3, 0);
848 GlobalPosition globalPos3(0);
850 GlobalPosition globalPosFace23(0);
851 GlobalPosition globalPosFace34(0);
853 for (
const auto& intersection23
854 : intersections(problem_.gridView(), element2))
856 bool finished =
false;
858 for (
const auto& intersection43
859 : intersections(problem_.gridView(), element4))
861 if (intersection23.neighbor() && intersection43.neighbor())
863 auto element32 = intersection23.outside();
864 auto element34 = intersection43.outside();
867 if (element32 == element34 && element32 != element)
870 interactionVolumes_[globalVertIdx1234].setSubVolumeElement(element32, 2);
872 interactionVolumes_[globalVertIdx1234].setIndexOnElement(intersection23.indexInInside(), 1,
874 interactionVolumes_[globalVertIdx1234].setIndexOnElement(intersection23.indexInOutside(), 2,
876 interactionVolumes_[globalVertIdx1234].setIndexOnElement(intersection43.indexInInside(), 3,
878 interactionVolumes_[globalVertIdx1234].setIndexOnElement(intersection43.indexInOutside(), 2,
882 globalPos3 = element32.geometry().center();
884 globalPosFace23 = intersection23.geometry().center();
885 globalPosFace34 = intersection43.geometry().center();
887 Scalar faceVol23 = intersection23.geometry().volume() / 2.0;
888 Scalar faceVol34 = intersection43.geometry().volume() / 2.0;
891 DimVector unitOuterNormal23 = intersection23.centerUnitOuterNormal();
893 DimVector unitOuterNormal43 = intersection43.centerUnitOuterNormal();
895 interactionVolumes_[globalVertIdx1234].setNormal(unitOuterNormal23, 1, 0);
896 unitOuterNormal23 *= -1;
897 interactionVolumes_[globalVertIdx1234].setNormal(unitOuterNormal23, 2, 1);
898 interactionVolumes_[globalVertIdx1234].setNormal(unitOuterNormal43, 3, 1);
899 unitOuterNormal43 *= -1;
900 interactionVolumes_[globalVertIdx1234].setNormal(unitOuterNormal43, 2, 0);
901 interactionVolumes_[globalVertIdx1234].setFaceArea(faceVol23, 1, 0);
902 interactionVolumes_[globalVertIdx1234].setFaceArea(faceVol23, 2, 1);
903 interactionVolumes_[globalVertIdx1234].setFaceArea(faceVol34, 2, 0);
904 interactionVolumes_[globalVertIdx1234].setFaceArea(faceVol34, 3, 1);
905 interactionVolumes_[globalVertIdx1234].setFacePosition(globalPosFace23, 1, 0);
906 interactionVolumes_[globalVertIdx1234].setFacePosition(globalPosFace23, 2, 1);
907 interactionVolumes_[globalVertIdx1234].setFacePosition(globalPosFace34, 2, 0);
908 interactionVolumes_[globalVertIdx1234].setFacePosition(globalPosFace34, 3, 1);
926 problem_.boundaryTypes(bcType, intersection14);
927 PrimaryVariables boundValues(0.0);
929 interactionVolumes_[globalVertIdx1234].setBoundary(bcType, 3);
930 if (bcType.isNeumann(pressEqIdx))
932 problem_.neumann(boundValues, intersection14);
933 boundValues *= faceVol41;
934 interactionVolumes_[globalVertIdx1234].setNeumannCondition(boundValues, 3);
936 if (bcType.hasDirichlet())
938 problem_.dirichlet(boundValues, intersection14);
939 interactionVolumes_[globalVertIdx1234].setDirichletCondition(boundValues, 3);
946 GlobalPosition globalPosFace23(0);
949 Scalar faceVol23 = 0;
952 DimVector unitOuterNormal23(0);
954 bool finished =
false;
956 for (
const auto& intersection2
957 : intersections(problem_.gridView(), element2))
959 if (intersection2.boundary())
961 for (
int i = 0; i < intersection2.geometry().corners(); ++i)
963 int localVertIdx2corner = refElement.subEntity(intersection2.indexInInside(), dim - 1, i,
966 int globalVertIdx2corner = problem_.variables().index(
967 element2.template subEntity < dim > (localVertIdx2corner));
969 if (globalVertIdx2corner == globalVertIdx1234)
971 interactionVolumes_[globalVertIdx1234].setIndexOnElement(intersection2.indexInInside(), 1,
974 globalPosFace23 = intersection2.geometry().center();
976 faceVol23 = intersection2.geometry().volume() / 2.0;
978 unitOuterNormal23 = intersection2.centerUnitOuterNormal();
980 interactionVolumes_[globalVertIdx1234].setNormal(unitOuterNormal23, 1, 0);
981 interactionVolumes_[globalVertIdx1234].setFaceArea(faceVol23, 1, 0);
982 interactionVolumes_[globalVertIdx1234].setFacePosition(globalPosFace23, 1, 0);
984 problem_.boundaryTypes(bcType, intersection2);
987 interactionVolumes_[globalVertIdx1234].setBoundary(bcType, 1);
988 if (bcType.isNeumann(pressEqIdx))
990 problem_.neumann(boundValues, intersection2);
991 boundValues *= faceVol23;
992 interactionVolumes_[globalVertIdx1234].setNeumannCondition(boundValues, 1);
994 if (bcType.hasDirichlet())
996 problem_.dirichlet(boundValues, intersection2);
997 interactionVolumes_[globalVertIdx1234].setDirichletCondition(boundValues, 1);
1000 interactionVolumes_[globalVertIdx1234].setOutsideFace(2);
1002 innerBoundaryVolumeFaces_[eIdxGlobal1][intersection12.indexInInside()] =
true;
1003 innerBoundaryVolumeFaces_[eIdxGlobal2][intersection12.indexInOutside()] =
true;
1019 Dune::NotImplemented,
1020 "fvmpfao2pfaboundpressure2p.hh, l. 997: boundary shape not available as interaction volume shape");
1028 problem_.boundaryTypes(bcType, intersection12);
1029 PrimaryVariables boundValues(0.0);
1031 interactionVolumes_[globalVertIdx1234].setBoundary(bcType, 0);
1032 if (bcType.isNeumann(pressEqIdx))
1034 problem_.neumann(boundValues, intersection12);
1035 boundValues *= faceVol12;
1036 interactionVolumes_[globalVertIdx1234].setNeumannCondition(boundValues, 0);
1038 if (bcType.hasDirichlet())
1040 problem_.dirichlet(boundValues, intersection12);
1041 interactionVolumes_[globalVertIdx1234].setDirichletCondition(boundValues, 0);
1045 if (intersection14.boundary())
1047 problem_.boundaryTypes(bcType, intersection14);
1050 interactionVolumes_[globalVertIdx1234].setBoundary(bcType, 3);
1051 if (bcType.isNeumann(pressEqIdx))
1053 problem_.neumann(boundValues, intersection14);
1054 boundValues *= faceVol41;
1055 interactionVolumes_[globalVertIdx1234].setNeumannCondition(boundValues, 3);
1057 if (bcType.hasDirichlet())
1059 problem_.dirichlet(boundValues, intersection14);
1060 interactionVolumes_[globalVertIdx1234].setDirichletCondition(boundValues, 3);
1063 interactionVolumes_[globalVertIdx1234].setOutsideFace(1);
1064 interactionVolumes_[globalVertIdx1234].setOutsideFace(2);
1072 auto element4 = intersection14.outside();
1073 interactionVolumes_[globalVertIdx1234].setIndexOnElement(intersection14.indexInOutside(), 3, 0);
1076 interactionVolumes_[globalVertIdx1234].setSubVolumeElement(element4, 3);
1078 interactionVolumes_[globalVertIdx1234].setNormal(unitOuterNormal14, 3, 0);
1079 interactionVolumes_[globalVertIdx1234].setFaceArea(faceVol41, 3, 0);
1080 interactionVolumes_[globalVertIdx1234].setFacePosition(globalPosFace41, 3, 0);
1082 int eIdxGlobal4 = problem_.variables().index(element4);
1084 bool finished =
false;
1087 for (
const auto& intersection4
1088 : intersections(problem_.gridView(), element4))
1090 if (intersection4.boundary())
1092 for (
int i = 0; i < intersection4.geometry().corners(); ++i)
1094 int localVertIdx4corner = refElement.subEntity(intersection4.indexInInside(), dim - 1, i,
1097 int globalVertIdx4corner = problem_.variables().index(
1098 (element4).
template subEntity < dim > (localVertIdx4corner));
1100 if (globalVertIdx4corner == globalVertIdx1234)
1102 interactionVolumes_[globalVertIdx1234].setIndexOnElement(intersection4.indexInInside(), 3,
1105 const GlobalPosition& globalPosFace34 = intersection4.geometry().center();
1107 Scalar faceVol34 = intersection4.geometry().volume() / 2.0;
1109 DimVector unitOuterNormal43 = intersection4.centerUnitOuterNormal();
1111 interactionVolumes_[globalVertIdx1234].setNormal(unitOuterNormal43, 3, 1);
1112 interactionVolumes_[globalVertIdx1234].setFaceArea(faceVol34, 3, 1);
1113 interactionVolumes_[globalVertIdx1234].setFacePosition(globalPosFace34, 3, 1);
1115 problem_.boundaryTypes(bcType, intersection4);
1118 interactionVolumes_[globalVertIdx1234].setBoundary(bcType, 2);
1119 if (bcType.isNeumann(pressEqIdx))
1121 problem_.neumann(boundValues, intersection4);
1122 boundValues *= faceVol34;
1123 interactionVolumes_[globalVertIdx1234].setNeumannCondition(boundValues, 2);
1125 if (bcType.hasDirichlet())
1127 problem_.dirichlet(boundValues, intersection4);
1128 interactionVolumes_[globalVertIdx1234].setDirichletCondition(boundValues, 2);
1131 interactionVolumes_[globalVertIdx1234].setOutsideFace(1);
1133 innerBoundaryVolumeFaces_[eIdxGlobal1][intersection14.indexInInside()] =
true;
1134 innerBoundaryVolumeFaces_[eIdxGlobal4][intersection14.indexInOutside()] =
true;
1150 Dune::NotImplemented,
1151 "fvmpfao2pfaboundpressure2p.hh, l. 1164: boundary shape not available as interaction volume shape");
1164template<
class TypeTag>
1165void FvMpfaL2dPressure2p<TypeTag>::assemble()
1172 for (
const auto& vertex : vertices(problem_.gridView()))
1174 int vIdxGlobal = problem_.variables().index(vertex);
1176 InteractionVolume& interactionVolume = interactionVolumes_[vIdxGlobal];
1178 if (interactionVolume.isInnerVolume())
1181 auto element1 = interactionVolume.getSubVolumeElement(0);
1182 auto element2 = interactionVolume.getSubVolumeElement(1);
1183 auto element3 = interactionVolume.getSubVolumeElement(2);
1184 auto element4 = interactionVolume.getSubVolumeElement(3);
1187 const GlobalPosition& globalPos1 = element1.geometry().center();
1188 const GlobalPosition& globalPos2 = element2.geometry().center();
1189 const GlobalPosition& globalPos3 = element3.geometry().center();
1190 const GlobalPosition& globalPos4 = element4.geometry().center();
1193 Scalar volume1 = element1.geometry().volume();
1194 Scalar volume2 = element2.geometry().volume();
1195 Scalar volume3 = element3.geometry().volume();
1196 Scalar volume4 = element4.geometry().volume();
1199 int eIdxGlobal1 = problem_.variables().index(element1);
1200 int eIdxGlobal2 = problem_.variables().index(element2);
1201 int eIdxGlobal3 = problem_.variables().index(element3);
1202 int eIdxGlobal4 = problem_.variables().index(element4);
1205 CellData& cellData1 = problem_.variables().cellData(eIdxGlobal1);
1206 CellData& cellData2 = problem_.variables().cellData(eIdxGlobal2);
1207 CellData& cellData3 = problem_.variables().cellData(eIdxGlobal3);
1208 CellData& cellData4 = problem_.variables().cellData(eIdxGlobal4);
1211 PrimaryVariables source(0.0);
1212 problem_.source(source, element1);
1213 this->f_[eIdxGlobal1] += volume1 / (4.0) * (source[wPhaseIdx] / density_[wPhaseIdx] + source[nPhaseIdx] / density_[nPhaseIdx]);
1214 problem_.source(source, element2);
1215 this->f_[eIdxGlobal2] += volume2 / (4.0) * (source[wPhaseIdx] / density_[wPhaseIdx] + source[nPhaseIdx] / density_[nPhaseIdx]);
1216 problem_.source(source, element3);
1217 this->f_[eIdxGlobal3] += volume3 / (4.0) * (source[wPhaseIdx] / density_[wPhaseIdx] + source[nPhaseIdx] / density_[nPhaseIdx]);
1218 problem_.source(source, element4);
1219 this->f_[eIdxGlobal4] += volume4 / (4.0) * (source[wPhaseIdx] / density_[wPhaseIdx] + source[nPhaseIdx] / density_[nPhaseIdx]);
1221 this->f_[eIdxGlobal1] += evaluateErrorTerm_(cellData1) * volume1 / (4.0);
1222 this->f_[eIdxGlobal2] += evaluateErrorTerm_(cellData2) * volume2 / (4.0);
1223 this->f_[eIdxGlobal3] += evaluateErrorTerm_(cellData3) * volume3 / (4.0);
1224 this->f_[eIdxGlobal4] += evaluateErrorTerm_(cellData4) * volume4 / (4.0);
1227 Dune::FieldVector<Scalar, numPhases> lambda1(cellData1.mobility(wPhaseIdx));
1228 lambda1[nPhaseIdx] = cellData1.mobility(nPhaseIdx);
1231 Scalar lambdaTotal1 = lambda1[wPhaseIdx] + lambda1[nPhaseIdx];
1234 Dune::FieldVector<Scalar, numPhases> lambda2(cellData2.mobility(wPhaseIdx));
1235 lambda2[nPhaseIdx] = cellData2.mobility(nPhaseIdx);
1238 Scalar lambdaTotal2 = lambda2[wPhaseIdx] + lambda2[nPhaseIdx];
1241 Dune::FieldVector<Scalar, numPhases> lambda3(cellData3.mobility(wPhaseIdx));
1242 lambda3[nPhaseIdx] = cellData3.mobility(nPhaseIdx);
1245 Scalar lambdaTotal3 = lambda3[wPhaseIdx] + lambda3[nPhaseIdx];
1248 Dune::FieldVector<Scalar, numPhases> lambda4(cellData4.mobility(wPhaseIdx));
1249 lambda4[nPhaseIdx] = cellData4.mobility(nPhaseIdx);
1252 Scalar lambdaTotal4 = lambda4[wPhaseIdx] + lambda4[nPhaseIdx];
1254 std::vector<DimVector > lambda(2*dim);
1256 lambda[0][0] = lambdaTotal1;
1257 lambda[0][1] = lambdaTotal1;
1258 lambda[1][0] = lambdaTotal2;
1259 lambda[1][1] = lambdaTotal2;
1260 lambda[2][0] = lambdaTotal3;
1261 lambda[2][1] = lambdaTotal3;
1262 lambda[3][0] = lambdaTotal4;
1263 lambda[3][1] = lambdaTotal4;
1267 Dune::FieldVector<Scalar, 2 * dim> pc(0);
1268 pc[0] = cellData1.capillaryPressure();
1269 pc[1] = cellData2.capillaryPressure();
1270 pc[2] = cellData3.capillaryPressure();
1271 pc[3] = cellData4.capillaryPressure();
1273 Dune::FieldVector<Scalar, 2 * dim> gravityDiff(0);
1277 gravityDiff[0] = (problem_.bBoxMax() - globalPos1) * gravity_ * (density_[nPhaseIdx] - density_[wPhaseIdx]);
1278 gravityDiff[1] = (problem_.bBoxMax() - globalPos2) * gravity_ * (density_[nPhaseIdx] - density_[wPhaseIdx]);
1279 gravityDiff[2] = (problem_.bBoxMax() - globalPos3) * gravity_ * (density_[nPhaseIdx] - density_[wPhaseIdx]);
1280 gravityDiff[3] = (problem_.bBoxMax() - globalPos4) * gravity_ * (density_[nPhaseIdx] - density_[wPhaseIdx]);
1284 Dune::FieldVector<Scalar, 2 * dim> pcFlux(0);
1286 Scalar pcPotential12 = 0;
1287 Scalar pcPotential14 = 0;
1288 Scalar pcPotential32 = 0;
1289 Scalar pcPotential34 = 0;
1292 Dune::FieldVector<Scalar, 2 * dim - dim + 1> u(0);
1293 Dune::FieldMatrix<Scalar,dim,2*dim-dim+1> T(0);
1295 int lType = transmissibilityCalculator_.calculateTransmissibility(T, interactionVolume, lambda, 0, 1, 2, 3);
1297 if (lType == TransmissibilityCalculator::rightTriangle)
1299 if (innerBoundaryVolumeFaces_[eIdxGlobal1][interactionVolume.getIndexOnElement(0, 0)])
1303 this->A_[eIdxGlobal1][eIdxGlobal2] += T[1][0];
1304 this->A_[eIdxGlobal1][eIdxGlobal3] += T[1][1];
1305 this->A_[eIdxGlobal1][eIdxGlobal1] += T[1][2];
1307 this->A_[eIdxGlobal2][eIdxGlobal2] -= T[1][0];
1308 this->A_[eIdxGlobal2][eIdxGlobal3] -= T[1][1];
1309 this->A_[eIdxGlobal2][eIdxGlobal1] -= T[1][2];
1318 pcPotential12 = Tu[1];
1322 if (innerBoundaryVolumeFaces_[eIdxGlobal1][interactionVolume.getIndexOnElement(0, 0)])
1326 this->A_[eIdxGlobal1][eIdxGlobal1] += T[1][0];
1327 this->A_[eIdxGlobal1][eIdxGlobal4] += T[1][1];
1328 this->A_[eIdxGlobal1][eIdxGlobal2] += T[1][2];
1330 this->A_[eIdxGlobal2][eIdxGlobal1] -= T[1][0];
1331 this->A_[eIdxGlobal2][eIdxGlobal4] -= T[1][1];
1332 this->A_[eIdxGlobal2][eIdxGlobal2] -= T[1][2];
1341 pcPotential12 = Tu[1];
1344 lType = transmissibilityCalculator_.calculateTransmissibility(T, interactionVolume, lambda, 1, 2, 3, 0);
1346 if (lType == TransmissibilityCalculator::rightTriangle)
1348 if (innerBoundaryVolumeFaces_[eIdxGlobal2][interactionVolume.getIndexOnElement(1, 0)])
1352 this->A_[eIdxGlobal2][eIdxGlobal3] += T[1][0];
1353 this->A_[eIdxGlobal2][eIdxGlobal4] += T[1][1];
1354 this->A_[eIdxGlobal2][eIdxGlobal2] += T[1][2];
1356 this->A_[eIdxGlobal3][eIdxGlobal3] -= T[1][0];
1357 this->A_[eIdxGlobal3][eIdxGlobal4] -= T[1][1];
1358 this->A_[eIdxGlobal3][eIdxGlobal2] -= T[1][2];
1367 pcPotential32 = -Tu[1];
1371 if (innerBoundaryVolumeFaces_[eIdxGlobal2][interactionVolume.getIndexOnElement(1, 0)])
1375 this->A_[eIdxGlobal2][eIdxGlobal2] += T[1][0];
1376 this->A_[eIdxGlobal2][eIdxGlobal1] += T[1][1];
1377 this->A_[eIdxGlobal2][eIdxGlobal3] += T[1][2];
1379 this->A_[eIdxGlobal3][eIdxGlobal2] -= T[1][0];
1380 this->A_[eIdxGlobal3][eIdxGlobal1] -= T[1][1];
1381 this->A_[eIdxGlobal3][eIdxGlobal3] -= T[1][2];
1390 pcPotential32 = -Tu[1];
1393 lType = transmissibilityCalculator_.calculateTransmissibility(T, interactionVolume, lambda, 2, 3, 0, 1);
1395 if (lType == TransmissibilityCalculator::rightTriangle)
1397 if (innerBoundaryVolumeFaces_[eIdxGlobal3][interactionVolume.getIndexOnElement(2, 0)])
1401 this->A_[eIdxGlobal3][eIdxGlobal4] += T[1][0];
1402 this->A_[eIdxGlobal3][eIdxGlobal1] += T[1][1];
1403 this->A_[eIdxGlobal3][eIdxGlobal3] += T[1][2];
1405 this->A_[eIdxGlobal4][eIdxGlobal4] -= T[1][0];
1406 this->A_[eIdxGlobal4][eIdxGlobal1] -= T[1][1];
1407 this->A_[eIdxGlobal4][eIdxGlobal3] -= T[1][2];
1416 pcPotential34 = Tu[1];
1420 if (innerBoundaryVolumeFaces_[eIdxGlobal3][interactionVolume.getIndexOnElement(2, 0)])
1424 this->A_[eIdxGlobal3][eIdxGlobal3] += T[1][0];
1425 this->A_[eIdxGlobal3][eIdxGlobal2] += T[1][1];
1426 this->A_[eIdxGlobal3][eIdxGlobal4] += T[1][2];
1428 this->A_[eIdxGlobal4][eIdxGlobal3] -= T[1][0];
1429 this->A_[eIdxGlobal4][eIdxGlobal2] -= T[1][1];
1430 this->A_[eIdxGlobal4][eIdxGlobal4] -= T[1][2];
1439 pcPotential34 = Tu[1];
1442 lType = transmissibilityCalculator_.calculateTransmissibility(T, interactionVolume, lambda, 3, 0, 1, 2);
1444 if (lType == TransmissibilityCalculator::rightTriangle)
1446 if (innerBoundaryVolumeFaces_[eIdxGlobal4][interactionVolume.getIndexOnElement(3, 0)])
1450 this->A_[eIdxGlobal4][eIdxGlobal1] += T[1][0];
1451 this->A_[eIdxGlobal4][eIdxGlobal2] += T[1][1];
1452 this->A_[eIdxGlobal4][eIdxGlobal4] += T[1][2];
1454 this->A_[eIdxGlobal1][eIdxGlobal1] -= T[1][0];
1455 this->A_[eIdxGlobal1][eIdxGlobal2] -= T[1][1];
1456 this->A_[eIdxGlobal1][eIdxGlobal4] -= T[1][2];
1465 pcPotential14 = -Tu[1];
1469 if (innerBoundaryVolumeFaces_[eIdxGlobal4][interactionVolume.getIndexOnElement(3, 0)])
1473 this->A_[eIdxGlobal4][eIdxGlobal4] += T[1][0];
1474 this->A_[eIdxGlobal4][eIdxGlobal3] += T[1][1];
1475 this->A_[eIdxGlobal4][eIdxGlobal1] += T[1][2];
1477 this->A_[eIdxGlobal1][eIdxGlobal4] -= T[1][0];
1478 this->A_[eIdxGlobal1][eIdxGlobal3] -= T[1][1];
1479 this->A_[eIdxGlobal1][eIdxGlobal1] -= T[1][2];
1488 pcPotential14 = -Tu[1];
1491 if (pc[0] == 0 && pc[1] == 0 && pc[2] == 0 && pc[3] == 0)
1497 Dune::FieldVector<Scalar, numPhases> lambda12Upw(0.0);
1498 lambda12Upw[wPhaseIdx] = (pcPotential12 >= 0) ? lambda1[wPhaseIdx] : lambda2[wPhaseIdx];
1499 lambda12Upw[nPhaseIdx] = (pcPotential12 >= 0) ? lambda1[nPhaseIdx] : lambda2[nPhaseIdx];
1502 Dune::FieldVector<Scalar, numPhases> lambda14Upw(0.0);
1503 lambda14Upw[wPhaseIdx] = (pcPotential14 >= 0) ? lambda1[wPhaseIdx] : lambda4[wPhaseIdx];
1504 lambda14Upw[nPhaseIdx] = (pcPotential14 >= 0) ? lambda1[nPhaseIdx] : lambda4[nPhaseIdx];
1507 Dune::FieldVector<Scalar, numPhases> lambda32Upw(0.0);
1508 lambda32Upw[wPhaseIdx] = (pcPotential32 >= 0) ? lambda3[wPhaseIdx] : lambda2[wPhaseIdx];
1509 lambda32Upw[nPhaseIdx] = (pcPotential32 >= 0) ? lambda3[nPhaseIdx] : lambda2[nPhaseIdx];
1512 Dune::FieldVector<Scalar, numPhases> lambda34Upw(0.0);
1513 lambda34Upw[wPhaseIdx] = (pcPotential34 >= 0) ? lambda3[wPhaseIdx] : lambda4[wPhaseIdx];
1514 lambda34Upw[nPhaseIdx] = (pcPotential34 >= 0) ? lambda3[nPhaseIdx] : lambda4[nPhaseIdx];
1516 for (
int i = 0; i < numPhases; i++)
1518 Scalar lambdaT12 = lambda12Upw[wPhaseIdx] + lambda12Upw[nPhaseIdx];
1519 Scalar lambdaT14 = lambda14Upw[wPhaseIdx] + lambda14Upw[nPhaseIdx];
1520 Scalar lambdaT32 = lambda32Upw[wPhaseIdx] + lambda32Upw[nPhaseIdx];
1521 Scalar lambdaT34 = lambda34Upw[wPhaseIdx] + lambda34Upw[nPhaseIdx];
1522 Scalar fracFlow12 = (lambdaT12 > threshold_) ? lambda12Upw[i] / (lambdaT12) : 0.0;
1523 Scalar fracFlow14 = (lambdaT14 > threshold_) ? lambda14Upw[i] / (lambdaT14) : 0.0;
1524 Scalar fracFlow32 = (lambdaT32 > threshold_) ? lambda32Upw[i] / (lambdaT32) : 0.0;
1525 Scalar fracFlow34 = (lambdaT34 > threshold_) ? lambda34Upw[i] / (lambdaT34) : 0.0;
1527 Dune::FieldVector<Scalar, 2 * dim> pcFluxReal(pcFlux);
1529 pcFluxReal[0] *= fracFlow12;
1530 pcFluxReal[1] *= fracFlow32;
1531 pcFluxReal[2] *= fracFlow34;
1532 pcFluxReal[3] *= fracFlow14;
1537 switch (pressureType_)
1544 this->f_[eIdxGlobal1] -= (pcFluxReal[0] - pcFluxReal[3]);
1545 this->f_[eIdxGlobal2] -= (pcFluxReal[1] - pcFluxReal[0]);
1546 this->f_[eIdxGlobal3] -= (pcFluxReal[2] - pcFluxReal[1]);
1547 this->f_[eIdxGlobal4] -= (pcFluxReal[3] - pcFluxReal[2]);
1557 this->f_[eIdxGlobal1] += (pcFluxReal[0] - pcFluxReal[3]);
1558 this->f_[eIdxGlobal2] += (pcFluxReal[1] - pcFluxReal[0]);
1559 this->f_[eIdxGlobal3] += (pcFluxReal[2] - pcFluxReal[1]);
1560 this->f_[eIdxGlobal4] += (pcFluxReal[3] - pcFluxReal[2]);
1571 for (
int elemIdx = 0; elemIdx < 2 * dim; elemIdx++)
1573 bool isOutside =
false;
1574 for (
int fIdx = 0; fIdx < dim; fIdx++)
1576 int intVolFaceIdx = interactionVolume.getFaceIndexFromSubVolume(elemIdx, fIdx);
1577 if (interactionVolume.isOutsideFace(intVolFaceIdx))
1588 auto element = interactionVolume.getSubVolumeElement(elemIdx);
1591 const GlobalPosition& globalPos = element.geometry().center();
1594 Scalar volume = element.geometry().volume();
1597 int eIdxGlobal = problem_.variables().index(element);
1600 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
1603 DimMatrix
permeability(problem_.spatialParams().intrinsicPermeability(element));
1606 PrimaryVariables source(0);
1607 problem_.source(source, element);
1608 this->f_[eIdxGlobal] += volume / (4.0) * (source[wPhaseIdx] / density_[wPhaseIdx] + source[nPhaseIdx] / density_[nPhaseIdx]);
1610 this->f_[eIdxGlobal] += evaluateErrorTerm_(cellData) * volume / (4.0);
1613 Dune::FieldVector<Scalar, numPhases> lambda(cellData.mobility(wPhaseIdx));
1614 lambda[nPhaseIdx] = cellData.mobility(nPhaseIdx);
1616 Scalar pc = cellData.capillaryPressure();
1618 Scalar gravityDiff = (problem_.bBoxMax() - globalPos) * gravity_ * (density_[nPhaseIdx] - density_[wPhaseIdx]);
1622 for (
int fIdx = 0; fIdx < dim; fIdx++)
1624 int intVolFaceIdx = interactionVolume.getFaceIndexFromSubVolume(elemIdx, fIdx);
1626 if (interactionVolume.isBoundaryFace(intVolFaceIdx))
1629 if (interactionVolume.getBoundaryType(intVolFaceIdx).isDirichlet(pressEqIdx))
1631 int boundaryFaceIdx = interactionVolume.getIndexOnElement(elemIdx, fIdx);
1633 const auto refElement = referenceElement(element);
1635 const LocalPosition& localPos = refElement.position(boundaryFaceIdx, 1);
1637 const GlobalPosition& globalPosFace = element.geometry().global(localPos);
1639 DimVector distVec(globalPosFace - globalPos);
1640 Scalar dist = distVec.two_norm();
1641 DimVector unitDistVec(distVec);
1642 unitDistVec /= dist;
1644 Scalar faceArea = interactionVolume.getFaceArea(elemIdx, fIdx);
1647 Scalar satWBound = cellData.saturation(wPhaseIdx);
1649 if (interactionVolume.getBoundaryType(intVolFaceIdx).isDirichlet(satEqIdx))
1651 Scalar satBound = interactionVolume.getDirichletValues(intVolFaceIdx)[saturationIdx];
1652 switch (saturationType_)
1656 satWBound = satBound;
1661 satWBound = 1 - satBound;
1671 const auto fluidMatrixInteraction = Deprecated::makePcKrSw(Scalar{}, problem_.spatialParams(), element);
1673 Scalar pcBound = fluidMatrixInteraction.pc(satWBound);
1675 Scalar gravityDiffBound = (problem_.bBoxMax() - globalPosFace) * gravity_
1676 * (density_[nPhaseIdx] - density_[wPhaseIdx]);
1678 pcBound += gravityDiffBound;
1680 Dune::FieldVector<Scalar, numPhases> lambdaBound(fluidMatrixInteraction.krw(satWBound));
1681 lambdaBound[nPhaseIdx] = fluidMatrixInteraction.krn(satWBound);
1682 lambdaBound[wPhaseIdx] /= viscosity_[wPhaseIdx];
1683 lambdaBound[nPhaseIdx] /= viscosity_[nPhaseIdx];
1685 Scalar potentialBound = interactionVolume.getDirichletValues(intVolFaceIdx)[pressureIdx];
1686 Scalar gdeltaZ = (problem_.bBoxMax()-globalPosFace) * gravity_;
1689 Scalar potentialDiffW = 0;
1690 Scalar potentialDiffNw = 0;
1691 switch (pressureType_)
1695 potentialBound += density_[wPhaseIdx]*gdeltaZ;
1696 potentialDiffW = (cellData.potential(wPhaseIdx) - potentialBound) / dist;
1697 potentialDiffNw = (cellData.potential(nPhaseIdx) - potentialBound - pcBound) / dist;
1702 potentialBound += density_[nPhaseIdx]*gdeltaZ;
1703 potentialDiffW = (cellData.potential(wPhaseIdx) - potentialBound + pcBound) / dist;
1704 potentialDiffNw = (cellData.potential(nPhaseIdx) - potentialBound) / dist;
1709 Scalar lambdaTotal = (potentialDiffW >= 0.) ? lambda[wPhaseIdx] : lambdaBound[wPhaseIdx];
1710 lambdaTotal += (potentialDiffNw >= 0.) ? lambda[nPhaseIdx] : lambdaBound[nPhaseIdx];
1712 DimVector permTimesNormal(0);
1716 Scalar entry = lambdaTotal * (unitDistVec * permTimesNormal) / dist * faceArea;
1721 switch (pressureType_)
1726 DimVector pcGradient = unitDistVec;
1727 pcGradient *= (pc - pcBound) / dist;
1730 pcFlux = 0.5 * (lambda[nPhaseIdx] + lambdaBound[nPhaseIdx])
1731 * (permTimesNormal * pcGradient) * faceArea;
1738 DimVector pcGradient = unitDistVec;
1739 pcGradient *= (pc - pcBound) / dist;
1742 pcFlux = 0.5 * (lambda[wPhaseIdx] + lambdaBound[wPhaseIdx])
1743 * (permTimesNormal * pcGradient) * faceArea;
1751 this->A_[eIdxGlobal][eIdxGlobal] += entry;
1752 this->f_[eIdxGlobal] += entry * potentialBound;
1754 if (pc == 0 && pcBound == 0)
1759 for (
int i = 0; i < numPhases; i++)
1761 switch (pressureType_)
1768 this->f_[eIdxGlobal] -= pcFlux;
1777 this->f_[eIdxGlobal] += pcFlux;
1785 else if (interactionVolume.getBoundaryType(intVolFaceIdx).isNeumann(pressEqIdx))
1787 Scalar J = interactionVolume.getNeumannValues(intVolFaceIdx)[wPhaseIdx] / density_[wPhaseIdx];
1788 J += interactionVolume.getNeumannValues(intVolFaceIdx)[nPhaseIdx] / density_[nPhaseIdx];
1790 this->f_[eIdxGlobal] -= J;
1794 std::cout <<
"interactionVolume.getBoundaryType(intVolFaceIdx).isNeumann(pressEqIdx)"
1795 << interactionVolume.getBoundaryType(intVolFaceIdx).isNeumann(pressEqIdx) <<
"\n";
1796 DUNE_THROW(Dune::NotImplemented,
1797 "No valid boundary condition type defined for pressure equation!");
1808 if (problem_.gridView().comm().size() > 1)
1811 for (
const auto& element : elements(problem_.gridView()))
1813 if (element.partitionType() == Dune::InteriorEntity)
1817 int eIdxGlobalI = problem_.variables().index(element);
1819 this->A_[eIdxGlobalI] = 0.0;
1820 this->A_[eIdxGlobalI][eIdxGlobalI] = 1.0;
1821 this->f_[eIdxGlobalI] = this->
pressure()[eIdxGlobalI];
1833template<
class TypeTag>
1837 for (
const auto& element : elements(problem_.gridView()))
1839 int eIdxGlobal = problem_.variables().index(element);
1841 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
1843 Scalar satW = cellData.saturation(wPhaseIdx);
1848 const auto fluidMatrixInteraction = Deprecated::makePcKrSw(Scalar{}, problem_.spatialParams(), element);
1850 const Scalar pc = fluidMatrixInteraction.pc(satW);
1852 cellData.setCapillaryPressure(pc);
1855 const Scalar mobilityW = fluidMatrixInteraction.krw(satW) / viscosity_[wPhaseIdx];
1856 const Scalar mobilityNw = fluidMatrixInteraction.krn(satW) / viscosity_[nPhaseIdx];
1859 cellData.setMobility(wPhaseIdx, mobilityW);
1860 cellData.setMobility(nPhaseIdx, mobilityNw);
1863 cellData.setFracFlowFunc(wPhaseIdx, mobilityW / (mobilityW + mobilityNw));
1864 cellData.setFracFlowFunc(nPhaseIdx, mobilityNw / (mobilityW + mobilityNw));
Provides methods for transmissibility calculation 2-d.
Class including the information of an interaction volume of a MPFA L-method that does not change with...
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 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
Finite volume MPFA L-method discretization of a two-phase flow pressure equation of the sequential IM...
Definition: lmethod/2dpressure.hh:74
InnerBoundaryVolumeFaces innerBoundaryVolumeFaces_
Vector marking faces which intersect the boundary.
Definition: lmethod/2dpressure.hh:451
void storePressureSolution(const Element &element)
Stores the pressure solution of a cell.
Definition: lmethod/2dpressure.hh:240
void update()
Pressure update.
Definition: lmethod/2dpressure.hh:290
void initialize()
Initializes the pressure model.
Definition: lmethod/2dpressure.hh:193
void updateInteractionVolumeInfo()
Updates interaction volumes.
Definition: lmethod/2dpressure.hh:179
FvMpfaL2dPressure2p(Problem &problem)
Constructs a FvMpfaL2dPressure2p object.
Definition: lmethod/2dpressure.hh:410
GlobalInteractionVolumeVector interactionVolumes_
Global Vector of interaction volumes.
Definition: lmethod/2dpressure.hh:450
void storePressureSolution()
Globally stores the pressure solution.
Definition: lmethod/2dpressure.hh:227
void addOutputVtkFields(MultiWriter &writer)
Adds pressure output to the output file.
Definition: lmethod/2dpressure.hh:339
FvMpfaL2dTransmissibilityCalculator< TypeTag > TransmissibilityCalculator
Definition: lmethod/2dpressure.hh:151
void updateMaterialLaws()
constitutive functions are initialized and stored in the variables object
Definition: lmethod/2dpressure.hh:1834
Provides methods for transmissibility calculation in 2-d.
Definition: 2dtransmissibilitycalculator.hh:44
Class including the information of an interaction volume of a MPFA L-method that does not change with...
Definition: linteractionvolume.hh:41
The finite volume base class for the solution of a pressure equation.
Definition: sequential/cellcentered/pressure.hh:49
void initialize()
Initialize pressure model.
Definition: sequential/cellcentered/pressure.hh:213
PressureSolution & pressure()
Returns the vector containing the pressure solution.
Definition: sequential/cellcentered/pressure.hh:120
void solve()
Solves the global system of equations to get the spatial distribution of the pressure.
Definition: sequential/cellcentered/pressure.hh:527
Specifies the properties for immiscible 2p diffusion/pressure models.
Properties for a MPFA method.
Finite Volume Diffusion Model.