22#ifndef DUMUX_FVMPFAO2DPRESSURE2P_HH
23#define DUMUX_FVMPFAO2DPRESSURE2P_HH
65template<
class TypeTag>
73 dim = GridView::dimension, dimWorld = GridView::dimensionworld
79 using ReferenceElements = Dune::ReferenceElements<Scalar, dim>;
81 using SpatialParams =
typename GET_PROP_TYPE(TypeTag, SpatialParams);
82 using MaterialLaw =
typename SpatialParams::MaterialLaw;
84 using Indices =
typename GET_PROP_TYPE(TypeTag, ModelTraits)::Indices;
86 using FluidSystem =
typename GET_PROP_TYPE(TypeTag, FluidSystem);
87 using FluidState =
typename GET_PROP_TYPE(TypeTag, FluidState);
89 using BoundaryTypes =
typename GET_PROP_TYPE(TypeTag, BoundaryTypes);
91 using SolutionTypes =
typename GET_PROP(TypeTag, SolutionTypes);
92 using PrimaryVariables =
typename SolutionTypes::PrimaryVariables;
93 using ScalarSolutionType =
typename SolutionTypes::ScalarSolution;
95 using GridTypeIndices =
typename GET_PROP_TYPE(TypeTag, GridTypeIndices);
99 pw = Indices::pressureW,
100 pn = Indices::pressureNw,
101 sw = Indices::saturationW,
102 sn = Indices::saturationNw
106 wPhaseIdx = Indices::wPhaseIdx,
107 nPhaseIdx = Indices::nPhaseIdx,
108 pressureIdx = Indices::pressureIdx,
109 saturationIdx = Indices::saturationIdx,
110 pressEqIdx = Indices::pressureEqIdx,
111 satEqIdx = Indices::satEqIdx,
119 dirichletDirichlet = 1,
120 dirichletNeumann = 2,
124 using Element =
typename GridView::Traits::template Codim<0>::Entity;
125 using IntersectionIterator =
typename GridView::IntersectionIterator;
126 using Intersection =
typename GridView::Intersection;
128 using LocalPosition = Dune::FieldVector<Scalar, dim>;
129 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
130 using GravityVector = Dune::FieldVector<Scalar, dimWorld>;
131 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
133 using DimVector = Dune::FieldVector<Scalar, dim>;
137 using GlobalInteractionVolumeVector = std::vector<InteractionVolume>;
138 using InnerBoundaryVolumeFaces = std::vector<Dune::FieldVector<bool, 2*dim> >;
141 Intersection getNextIntersection_(
const Element&,
const IntersectionIterator&);
145 void initializeMatrix();
147 void storeInteractionVolumeInfo();
169 storeInteractionVolumeInfo();
181 const auto element = *problem_.gridView().template begin<0>();
182 FluidState fluidState;
183 fluidState.setPressure(wPhaseIdx, problem_.referencePressure(element));
184 fluidState.setPressure(nPhaseIdx, problem_.referencePressure(element));
185 fluidState.setTemperature(problem_.temperature(element));
186 fluidState.setSaturation(wPhaseIdx, 1.);
187 fluidState.setSaturation(nPhaseIdx, 0.);
198 storeInteractionVolumeInfo();
217 timeStep_ = problem_.timeManager().timeStepSize();
219 int size = problem_.gridView().size(0);
220 for (
int i = 0; i < size; i++)
224 switch (saturationType_)
227 sat = problem_.variables().cellData(i).saturation(wPhaseIdx);
230 sat = problem_.variables().cellData(i).saturation(nPhaseIdx);
235 maxError_ = max(maxError_, (sat - 1.0) / timeStep_);
239 maxError_ = max(maxError_, (-sat) / timeStep_);
257 for (
const auto& element : elements(problem_.gridView()))
270 int eIdxGlobal = problem_.variables().index(element);
271 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
273 switch (pressureType_)
277 Scalar potW = this->
pressure()[eIdxGlobal];
279 Scalar gravityDiff = (problem_.bBoxMax() - element.geometry().center()) * gravity_;
280 Scalar potPc = cellData.capillaryPressure()
281 + gravityDiff * (density_[nPhaseIdx] - density_[wPhaseIdx]);
283 cellData.setPotential(wPhaseIdx, potW);
284 cellData.setPotential(nPhaseIdx, potW + potPc);
286 Scalar pressW = potW - gravityDiff * density_[wPhaseIdx];
288 cellData.setPressure(wPhaseIdx, pressW);
289 cellData.setPressure(nPhaseIdx, pressW + cellData.capillaryPressure());
295 Scalar potNw = this->
pressure()[eIdxGlobal];
297 Scalar gravityDiff = (problem_.bBoxMax() - element.geometry().center()) * gravity_;
298 Scalar potPc = cellData.capillaryPressure()
299 + gravityDiff * (density_[nPhaseIdx] - density_[wPhaseIdx]);
301 cellData.setPotential(nPhaseIdx, potNw);
302 cellData.setPotential(wPhaseIdx, potNw - potPc);
304 Scalar pressNw = potNw - gravityDiff * density_[nPhaseIdx];
306 cellData.setPressure(wPhaseIdx, pressNw - cellData.capillaryPressure());
307 cellData.setPressure(nPhaseIdx, pressNw);
313 cellData.fluxData().resetVelocity();
326 template<
class MultiWriter>
329 int size = problem_.gridView().size(0);
330 ScalarSolutionType *potential = writer.allocateManagedBuffer(size);
334 if (pressureType_ == pw)
336 writer.attachCellData(*potential,
"wetting potential");
339 if (pressureType_ == pn)
341 writer.attachCellData(*potential,
"nonwetting potential");
344 if (vtkOutputLevel_ > 0)
346 ScalarSolutionType *
pressure = writer.allocateManagedBuffer(size);
347 ScalarSolutionType *pressureSecond = writer.allocateManagedBuffer(size);
348 ScalarSolutionType *potentialSecond = writer.allocateManagedBuffer(size);
349 ScalarSolutionType *pc = writer.allocateManagedBuffer(size);
351 for (
const auto& element : elements(problem_.gridView()))
353 int idx = problem_.variables().index(element);
354 CellData& cellData = problem_.variables().cellData(idx);
356 (*pc)[idx] = cellData.capillaryPressure();
358 if (pressureType_ == pw)
360 (*pressure)[idx] = cellData.pressure(wPhaseIdx);
361 (*potentialSecond)[idx] = cellData.potential(nPhaseIdx);
362 (*pressureSecond)[idx] = cellData.pressure(nPhaseIdx);
365 if (pressureType_ == pn)
367 (*pressure)[idx] = cellData.pressure(nPhaseIdx);
368 (*potentialSecond)[idx] = cellData.potential(wPhaseIdx);
369 (*pressureSecond)[idx] = cellData.pressure(wPhaseIdx);
373 if (pressureType_ == pw)
375 writer.attachCellData(*
pressure,
"wetting pressure");
376 writer.attachCellData(*pressureSecond,
"nonwetting pressure");
377 writer.attachCellData(*potentialSecond,
"nonwetting potential");
380 if (pressureType_ == pn)
382 writer.attachCellData(*
pressure,
"nonwetting pressure");
383 writer.attachCellData(*pressureSecond,
"wetting pressure");
384 writer.attachCellData(*potentialSecond,
"wetting potential");
387 writer.attachCellData(*pc,
"capillary pressure");
398 ParentType(problem), problem_(problem), gravity_(problem.gravity()), maxError_(
401 if (pressureType_ != pw && pressureType_ != pn)
403 DUNE_THROW(Dune::NotImplemented,
"Pressure type not supported!");
405 if (saturationType_ != sw && saturationType_ != sn)
407 DUNE_THROW(Dune::NotImplemented,
"Saturation type not supported!");
411 DUNE_THROW(Dune::NotImplemented,
"Compressibility not supported!");
415 DUNE_THROW(Dune::NotImplemented,
"Dimension not supported!");
418 ErrorTermFactor_ = getParam<Scalar>(
"Impet.ErrorTermFactor");
419 ErrorTermLowerBound_ = getParam<Scalar>(
"Impet.ErrorTermLowerBound");
420 ErrorTermUpperBound_ = getParam<Scalar>(
"Impet.ErrorTermUpperBound");
422 density_[wPhaseIdx] = 0.;
423 density_[nPhaseIdx] = 0.;
424 viscosity_[wPhaseIdx] = 0.;
425 viscosity_[nPhaseIdx] = 0.;
427 vtkOutputLevel_ = getParam<int>(
"Vtk.OutputLevel");
449 Scalar evaluateErrorTerm_(
CellData& cellData)
454 switch (saturationType_)
457 sat = cellData.saturation(wPhaseIdx);
460 sat = cellData.saturation(nPhaseIdx);
464 Scalar error = (sat > 1.0) ? sat - 1.0 : 0.0;
472 Scalar errorAbs = abs(error);
474 if ((errorAbs * timeStep_ > 1e-6) && (errorAbs > ErrorTermLowerBound_ * maxError_)
475 && (!problem_.timeManager().willBeFinished()))
477 return ErrorTermFactor_ * error;
487 const GravityVector& gravity_;
491 Scalar ErrorTermFactor_;
492 Scalar ErrorTermLowerBound_;
493 Scalar ErrorTermUpperBound_;
495 Scalar density_[numPhases];
496 Scalar viscosity_[numPhases];
500 static constexpr Scalar threshold_ = 1e-15;
502 static const int pressureType_ =
GET_PROP_VALUE(TypeTag, PressureFormulation);
510template<
class TypeTag>
511typename FvMpfaO2dPressure2p<TypeTag>::Intersection
513 const IntersectionIterator& isIt)
515 auto isItBegin = problem_.gridView().ibegin(element);
516 const auto isEndIt = problem_.gridView().iend(element);
518 auto tempIsIt = isIt;
519 auto nextIsIt = ++tempIsIt;
525 case GridTypeIndices::yaspGrid:
527 if (nextIsIt == isEndIt)
529 nextIsIt = isItBegin;
533 nextIsIt = ++tempIsIt;
535 if (nextIsIt == isEndIt)
537 auto tempIsItBegin = isItBegin;
538 nextIsIt = ++tempIsItBegin;
545 case GridTypeIndices::aluGrid:
546 case GridTypeIndices::ugGrid:
548 if (nextIsIt == isEndIt)
549 nextIsIt = isItBegin;
555 DUNE_THROW(Dune::NotImplemented,
"GridType can not be used with MPFAO implementation!");
564template<
class TypeTag>
565void FvMpfaO2dPressure2p<TypeTag>::initializeMatrix()
568 for (
const auto& element : elements(problem_.gridView()))
571 int eIdxGlobalI = problem_.variables().index(element);
577 const auto isEndIt = problem_.gridView().iend(element);
578 for (
auto isIt = problem_.gridView().ibegin(element); isIt != isEndIt; ++isIt)
580 const auto& intersection = *isIt;
581 auto nextIntersection = getNextIntersection_(element, isIt);
583 if (intersection.neighbor())
587 if (intersection.neighbor() && nextIntersection.neighbor())
589 for (
const auto& innerIntersection
590 :
intersections(problem_.gridView(), intersection.outside()))
591 for (
const auto& innerNextIntersection
592 :
intersections(problem_.gridView(), nextIntersection.outside()))
594 if (innerIntersection.neighbor() && innerNextIntersection.neighbor())
596 if (innerIntersection.outside() == innerNextIntersection.outside()
597 && innerIntersection.outside() != intersection.inside())
607 this->A_.setrowsize(eIdxGlobalI, rowSize);
612 this->A_.endrowsizes();
615 for (
const auto& element : elements(problem_.gridView()))
618 int eIdxGlobalI = problem_.variables().index(element);
621 this->A_.addindex(eIdxGlobalI, eIdxGlobalI);
624 const auto isEndIt = problem_.gridView().iend(element);
625 for (
auto isIt = problem_.gridView().ibegin(element); isIt != isEndIt; ++isIt)
627 const auto& intersection = *isIt;
628 auto nextIntersection = getNextIntersection_(element, isIt);
630 if (intersection.neighbor())
633 int eIdxGlobalJ = problem_.variables().index(intersection.outside());
637 this->A_.addindex(eIdxGlobalI, eIdxGlobalJ);
640 if (intersection.neighbor() && nextIntersection.neighbor())
642 for (
const auto& innerIntersection
643 :
intersections(problem_.gridView(), intersection.outside()))
644 for (
const auto& innerNextIntersection
645 :
intersections(problem_.gridView(), nextIntersection.outside()))
647 if (innerIntersection.neighbor() && innerNextIntersection.neighbor())
649 auto innerOutside = innerIntersection.outside();
651 if (innerOutside == innerNextIntersection.outside()
652 && innerOutside != intersection.inside())
654 int eIdxGlobalJ = problem_.variables().index(innerOutside);
656 this->A_.addindex(eIdxGlobalI, eIdxGlobalJ);
665 this->A_.endindices();
693template<
class TypeTag>
694void FvMpfaO2dPressure2p<TypeTag>::storeInteractionVolumeInfo()
706 BoundaryTypes bcType;
709 for (
const auto& element : elements(problem_.gridView()))
713 int eIdxGlobal1 = problem_.variables().index(element);
715 const GlobalPosition& globalPos1 = element.geometry().center();
718 DimMatrix K1(problem_.spatialParams().intrinsicPermeability(element));
720 const auto isIt12End = problem_.gridView().iend(element);
721 for (
auto isIt12 = problem_.gridView().ibegin(element); isIt12 != isIt12End; ++isIt12)
723 const auto& intersection12 = *isIt12;
724 auto intersection14 = getNextIntersection_(element, isIt12);
726 int indexInInside12 = intersection12.indexInInside();
727 int indexInInside14 = intersection14.indexInInside();
732 const auto referenceElement = ReferenceElements::general(element.type());
734 GlobalPosition corner1234(0);
736 int globalVertIdx1234 = 0;
739 for (
int i = 0; i < intersection12.geometry().corners(); ++i)
741 bool finished =
false;
743 const GlobalPosition& isIt12corner = intersection12.geometry().corner(i);
745 int localVertIdx12corner = referenceElement.subEntity(indexInInside12, dim - 1, i, dim);
747 int globalVertIdx12corner = problem_.variables().index(element.template subEntity<dim>(localVertIdx12corner));
749 for (
int j = 0; j < intersection14.geometry().corners(); ++j)
751 int localVertIdx14corner = referenceElement.subEntity(indexInInside14, dim - 1, j, dim);
753 int globalVertIdx14corner = problem_.variables().index(element.template subEntity<dim>(localVertIdx14corner));
755 if (globalVertIdx12corner == globalVertIdx14corner)
757 corner1234 = isIt12corner;
759 globalVertIdx1234 = globalVertIdx12corner;
772 if (interactionVolumes_[globalVertIdx1234].isStored())
778 interactionVolumes_[globalVertIdx1234].setStored();
783 interactionVolumes_[globalVertIdx1234].setSubVolumeElement(element, 0);
784 interactionVolumes_[globalVertIdx1234].setIndexOnElement(intersection12.indexInInside(), 0, 0);
785 interactionVolumes_[globalVertIdx1234].setIndexOnElement(intersection14.indexInInside(), 0, 1);
788 const GlobalPosition& globalPosFace12 = intersection12.geometry().center();
791 Scalar faceVol12 = intersection12.geometry().volume() / 2.0;
794 DimVector unitOuterNormal12 = intersection12.centerUnitOuterNormal();
797 const GlobalPosition& globalPosFace41 = intersection14.geometry().center();
800 Scalar faceVol41 = intersection14.geometry().volume() / 2.0;
803 DimVector unitOuterNormal14 = intersection14.centerUnitOuterNormal();
807 R.mv(globalPos1 - globalPosFace12, nu14);
810 R.mv(globalPosFace41 - globalPos1, nu12);
812 interactionVolumes_[globalVertIdx1234].setPermTimesNu(nu12, K1, 0, 0);
813 interactionVolumes_[globalVertIdx1234].setPermTimesNu(nu14, K1, 0, 1);
814 interactionVolumes_[globalVertIdx1234].setNormal(unitOuterNormal12, 0, 0);
815 interactionVolumes_[globalVertIdx1234].setNormal(unitOuterNormal14, 0, 1);
816 interactionVolumes_[globalVertIdx1234].setFaceArea(faceVol12, 0, 0);
817 interactionVolumes_[globalVertIdx1234].setFaceArea(faceVol41, 0, 1);
823 interactionVolumes_[globalVertIdx1234].setDF(abs(nu14 * Rnu12), 0);
826 if (intersection12.neighbor())
829 auto element2 = intersection12.outside();
831 int eIdxGlobal2 = problem_.variables().index(element2);
834 interactionVolumes_[globalVertIdx1234].setSubVolumeElement(element2, 1);
835 interactionVolumes_[globalVertIdx1234].setIndexOnElement(intersection12.indexInOutside(), 1, 1);
837 interactionVolumes_[globalVertIdx1234].setNormal(unitOuterNormal12, 1, 1);
838 interactionVolumes_[globalVertIdx1234].setFaceArea(faceVol12, 1, 1);
841 const GlobalPosition& globalPos2 = element2.geometry().center();
844 DimMatrix K2(problem_.spatialParams().intrinsicPermeability(element2));
847 if (intersection14.neighbor())
851 auto element4 = intersection14.outside();
854 interactionVolumes_[globalVertIdx1234].setSubVolumeElement(element4, 3);
855 interactionVolumes_[globalVertIdx1234].setIndexOnElement(intersection14.indexInOutside(), 3, 0);
857 interactionVolumes_[globalVertIdx1234].setNormal(unitOuterNormal14, 3, 0);
858 interactionVolumes_[globalVertIdx1234].setFaceArea(faceVol41, 3, 0);
863 const GlobalPosition& globalPos4 = element4.geometry().center();
866 DimMatrix K4(problem_.spatialParams().intrinsicPermeability(element4));
869 GlobalPosition globalPos3(0);
871 GlobalPosition globalPosFace23(0);
872 GlobalPosition globalPosFace34(0);
874 for (
const auto& intersection2
877 bool finished =
false;
879 for (
const auto& intersection4
882 if (intersection2.neighbor() && intersection4.neighbor())
884 auto element32 = intersection2.outside();
885 auto element34 = intersection4.outside();
888 if (element32 == element34 && element32 != element)
891 interactionVolumes_[globalVertIdx1234].setSubVolumeElement(element32, 2);
893 interactionVolumes_[globalVertIdx1234].setIndexOnElement(intersection2.indexInInside(), 1,
895 interactionVolumes_[globalVertIdx1234].setIndexOnElement(intersection2.indexInOutside(), 2,
897 interactionVolumes_[globalVertIdx1234].setIndexOnElement(intersection4.indexInInside(), 3,
899 interactionVolumes_[globalVertIdx1234].setIndexOnElement(intersection4.indexInOutside(), 2,
903 globalPos3 = element32.geometry().center();
905 globalPosFace23 = intersection2.geometry().center();
906 globalPosFace34 = intersection4.geometry().center();
908 Scalar faceVol23 = intersection2.geometry().volume() / 2.0;
909 Scalar faceVol34 = intersection4.geometry().volume() / 2.0;
912 DimVector unitOuterNormal23 = intersection2.centerUnitOuterNormal();
914 DimVector unitOuterNormal43 = intersection4.centerUnitOuterNormal();
916 interactionVolumes_[globalVertIdx1234].setNormal(unitOuterNormal23, 1, 0);
917 interactionVolumes_[globalVertIdx1234].setNormal(unitOuterNormal23, 2, 1);
918 interactionVolumes_[globalVertIdx1234].setNormal(unitOuterNormal43, 2, 0);
919 interactionVolumes_[globalVertIdx1234].setNormal(unitOuterNormal43, 3, 1);
920 interactionVolumes_[globalVertIdx1234].setFaceArea(faceVol23, 1, 0);
921 interactionVolumes_[globalVertIdx1234].setFaceArea(faceVol23, 2, 1);
922 interactionVolumes_[globalVertIdx1234].setFaceArea(faceVol34, 2, 0);
923 interactionVolumes_[globalVertIdx1234].setFaceArea(faceVol34, 3, 1);
927 problem_.spatialParams().intrinsicPermeability(element32));
931 R.umv(globalPosFace12 - globalPos2, nu23);
934 R.umv(globalPosFace23 - globalPos2, nu21);
937 R.umv(globalPosFace34 - globalPos3, nu32);
940 R.umv(globalPos3 - globalPosFace23, nu34);
943 R.umv(globalPos4 - globalPosFace34, nu41);
946 R.umv(globalPos4 - globalPosFace41, nu43);
948 interactionVolumes_[globalVertIdx1234].setPermTimesNu(nu23, K2, 1, 0);
949 interactionVolumes_[globalVertIdx1234].setPermTimesNu(nu21, K2, 1, 1);
950 interactionVolumes_[globalVertIdx1234].setPermTimesNu(nu34, K3, 2, 0);
951 interactionVolumes_[globalVertIdx1234].setPermTimesNu(nu32, K3, 2, 1);
952 interactionVolumes_[globalVertIdx1234].setPermTimesNu(nu41, K4, 3, 0);
953 interactionVolumes_[globalVertIdx1234].setPermTimesNu(nu43, K4, 3, 1);
958 interactionVolumes_[globalVertIdx1234].setDF(abs(nu23 * Rnu21), 1);
962 interactionVolumes_[globalVertIdx1234].setDF(abs(nu32 * Rnu34), 2);
966 interactionVolumes_[globalVertIdx1234].setDF(abs(nu41 * Rnu43), 3);
984 problem_.boundaryTypes(bcType, intersection14);
985 PrimaryVariables boundValues(0.0);
987 interactionVolumes_[globalVertIdx1234].setBoundary(bcType, 3);
988 if (bcType.isNeumann(pressEqIdx))
990 problem_.neumann(boundValues, intersection14);
991 boundValues *= faceVol41;
992 interactionVolumes_[globalVertIdx1234].setNeumannCondition(boundValues, 3);
994 if (bcType.hasDirichlet())
996 problem_.dirichlet(boundValues, intersection14);
997 interactionVolumes_[globalVertIdx1234].setDirichletCondition(boundValues, 3);
1004 GlobalPosition globalPosFace23(0);
1007 Scalar faceVol23 = 0;
1010 DimVector unitOuterNormal23(0);
1012 bool finished =
false;
1014 for (
const auto& intersection2
1017 if (intersection2.boundary())
1019 for (
int i = 0; i < intersection2.geometry().corners(); ++i)
1021 int localVertIdx2corner = referenceElement.subEntity(intersection2.indexInInside(), dim - 1, i,
1024 int globalVertIdx2corner = problem_.variables().index(element2.template subEntity<dim>(localVertIdx2corner));
1026 if (globalVertIdx2corner == globalVertIdx1234)
1028 interactionVolumes_[globalVertIdx1234].setIndexOnElement(intersection2.indexInInside(), 1,
1031 globalPosFace23 = intersection2.geometry().center();
1033 faceVol23 = intersection2.geometry().volume() / 2.0;
1035 unitOuterNormal23 = intersection2.centerUnitOuterNormal();
1037 interactionVolumes_[globalVertIdx1234].setNormal(unitOuterNormal23, 1, 0);
1038 interactionVolumes_[globalVertIdx1234].setFaceArea(faceVol23, 1, 0);
1040 problem_.boundaryTypes(bcType, intersection2);
1043 interactionVolumes_[globalVertIdx1234].setBoundary(bcType, 1);
1044 if (bcType.isNeumann(pressEqIdx))
1046 problem_.neumann(boundValues, intersection2);
1047 boundValues *= faceVol23;
1048 interactionVolumes_[globalVertIdx1234].setNeumannCondition(boundValues, 1);
1050 if (bcType.hasDirichlet())
1052 problem_.dirichlet(boundValues, intersection2);
1053 interactionVolumes_[globalVertIdx1234].setDirichletCondition(boundValues, 1);
1056 interactionVolumes_[globalVertIdx1234].setOutsideFace(2);
1058 innerBoundaryVolumeFaces_[eIdxGlobal1][intersection12.indexInInside()] =
true;
1059 innerBoundaryVolumeFaces_[eIdxGlobal2][intersection12.indexInOutside()] =
true;
1063 R.umv(globalPosFace12 - globalPos2, nu23);
1066 R.umv(globalPosFace23 - globalPos2, nu21);
1068 interactionVolumes_[globalVertIdx1234].setPermTimesNu(nu23, K2, 1, 0);
1069 interactionVolumes_[globalVertIdx1234].setPermTimesNu(nu21, K2, 1, 1);
1074 interactionVolumes_[globalVertIdx1234].setDF(abs(nu23 * Rnu21), 1);
1090 Dune::NotImplemented,
1091 "fvmpfao2pfaboundpressure2p.hh, l. 997: boundary shape not available as interaction volume shape");
1099 problem_.boundaryTypes(bcType, *isIt12);
1100 PrimaryVariables boundValues(0.0);
1102 interactionVolumes_[globalVertIdx1234].setBoundary(bcType, 0);
1103 if (bcType.isNeumann(pressEqIdx))
1105 problem_.neumann(boundValues, *isIt12);
1106 boundValues *= faceVol12;
1107 interactionVolumes_[globalVertIdx1234].setNeumannCondition(boundValues, 0);
1109 if (bcType.hasDirichlet())
1111 problem_.dirichlet(boundValues, *isIt12);
1112 interactionVolumes_[globalVertIdx1234].setDirichletCondition(boundValues, 0);
1116 if (intersection14.boundary())
1118 problem_.boundaryTypes(bcType, intersection14);
1121 interactionVolumes_[globalVertIdx1234].setBoundary(bcType, 3);
1122 if (bcType.isNeumann(pressEqIdx))
1124 problem_.neumann(boundValues, intersection14);
1125 boundValues *= faceVol41;
1126 interactionVolumes_[globalVertIdx1234].setNeumannCondition(boundValues, 3);
1128 if (bcType.hasDirichlet())
1130 problem_.dirichlet(boundValues, intersection14);
1131 interactionVolumes_[globalVertIdx1234].setDirichletCondition(boundValues, 3);
1134 interactionVolumes_[globalVertIdx1234].setOutsideFace(1);
1135 interactionVolumes_[globalVertIdx1234].setOutsideFace(2);
1143 auto element4 = intersection14.outside();
1144 interactionVolumes_[globalVertIdx1234].setIndexOnElement(intersection14.indexInOutside(), 3, 0);
1147 interactionVolumes_[globalVertIdx1234].setSubVolumeElement(element4, 3);
1149 interactionVolumes_[globalVertIdx1234].setNormal(unitOuterNormal14, 3, 0);
1150 interactionVolumes_[globalVertIdx1234].setFaceArea(faceVol41, 3, 0);
1153 const GlobalPosition& globalPos4 = element4.geometry().center();
1155 int eIdxGlobal4 = problem_.variables().index(element4);
1157 bool finished =
false;
1160 for (
const auto& intersection4
1163 if (intersection4.boundary())
1165 for (
int i = 0; i < intersection4.geometry().corners(); ++i)
1167 int localVertIdx4corner = referenceElement.subEntity(intersection4.indexInInside(), dim - 1, i,
1170 int globalVertIdx4corner = problem_.variables().index(element4.template subEntity<dim>(localVertIdx4corner));
1172 if (globalVertIdx4corner == globalVertIdx1234)
1174 interactionVolumes_[globalVertIdx1234].setIndexOnElement(intersection4.indexInInside(), 3,
1177 const GlobalPosition& globalPosFace34 = intersection4.geometry().center();
1179 Scalar faceVol34 = intersection4.geometry().volume() / 2.0;
1181 DimVector unitOuterNormal43 = intersection4.centerUnitOuterNormal();
1183 interactionVolumes_[globalVertIdx1234].setNormal(unitOuterNormal43, 3, 1);
1184 interactionVolumes_[globalVertIdx1234].setFaceArea(faceVol34, 3, 1);
1186 problem_.boundaryTypes(bcType, intersection4);
1189 interactionVolumes_[globalVertIdx1234].setBoundary(bcType, 2);
1190 if (bcType.isNeumann(pressEqIdx))
1192 problem_.neumann(boundValues, intersection4);
1193 boundValues *= faceVol34;
1194 interactionVolumes_[globalVertIdx1234].setNeumannCondition(boundValues, 2);
1196 if (bcType.hasDirichlet())
1198 problem_.dirichlet(boundValues, intersection4);
1199 interactionVolumes_[globalVertIdx1234].setDirichletCondition(boundValues, 2);
1202 interactionVolumes_[globalVertIdx1234].setOutsideFace(1);
1204 innerBoundaryVolumeFaces_[eIdxGlobal1][intersection14.indexInInside()] =
true;
1205 innerBoundaryVolumeFaces_[eIdxGlobal4][intersection14.indexInOutside()] =
true;
1209 problem_.spatialParams().intrinsicPermeability(element4));
1213 R.umv(globalPos4 - globalPosFace34, nu41);
1216 R.umv(globalPos4 - globalPosFace41, nu43);
1218 interactionVolumes_[globalVertIdx1234].setPermTimesNu(nu41, K4, 3, 0);
1219 interactionVolumes_[globalVertIdx1234].setPermTimesNu(nu43, K4, 3, 1);
1224 interactionVolumes_[globalVertIdx1234].setDF(abs(nu41 * Rnu43), 3);
1240 Dune::NotImplemented,
1241 "fvmpfao2pfaboundpressure2p.hh, l. 1164: boundary shape not available as interaction volume shape");
1254template<
class TypeTag>
1255void FvMpfaO2dPressure2p<TypeTag>::assemble()
1262 for (
const auto& vertex : vertices(problem_.gridView()))
1264 int vIdxGlobal = problem_.variables().index(vertex);
1266 InteractionVolume& interactionVolume = interactionVolumes_[vIdxGlobal];
1268 if (interactionVolume.isInnerVolume())
1271 auto element1 = interactionVolume.getSubVolumeElement(0);
1272 auto element2 = interactionVolume.getSubVolumeElement(1);
1273 auto element3 = interactionVolume.getSubVolumeElement(2);
1274 auto element4 = interactionVolume.getSubVolumeElement(3);
1277 const GlobalPosition& globalPos1 = element1.geometry().center();
1278 const GlobalPosition& globalPos2 = element2.geometry().center();
1279 const GlobalPosition& globalPos3 = element3.geometry().center();
1280 const GlobalPosition& globalPos4 = element4.geometry().center();
1283 Scalar volume1 = element1.geometry().volume();
1284 Scalar volume2 = element2.geometry().volume();
1285 Scalar volume3 = element3.geometry().volume();
1286 Scalar volume4 = element4.geometry().volume();
1289 int eIdxGlobal1 = problem_.variables().index(element1);
1290 int eIdxGlobal2 = problem_.variables().index(element2);
1291 int eIdxGlobal3 = problem_.variables().index(element3);
1292 int eIdxGlobal4 = problem_.variables().index(element4);
1295 CellData& cellData1 = problem_.variables().cellData(eIdxGlobal1);
1296 CellData& cellData2 = problem_.variables().cellData(eIdxGlobal2);
1297 CellData& cellData3 = problem_.variables().cellData(eIdxGlobal3);
1298 CellData& cellData4 = problem_.variables().cellData(eIdxGlobal4);
1301 PrimaryVariables source(0.0);
1302 problem_.source(source, element1);
1303 this->f_[eIdxGlobal1] += volume1 / (4.0)
1304 * (source[wPhaseIdx] / density_[wPhaseIdx] + source[nPhaseIdx] / density_[nPhaseIdx]);
1305 problem_.source(source, element2);
1306 this->f_[eIdxGlobal2] += volume2 / (4.0)
1307 * (source[wPhaseIdx] / density_[wPhaseIdx] + source[nPhaseIdx] / density_[nPhaseIdx]);
1308 problem_.source(source, element3);
1309 this->f_[eIdxGlobal3] += volume3 / (4.0)
1310 * (source[wPhaseIdx] / density_[wPhaseIdx] + source[nPhaseIdx] / density_[nPhaseIdx]);
1311 problem_.source(source, element4);
1312 this->f_[eIdxGlobal4] += volume4 / (4.0)
1313 * (source[wPhaseIdx] / density_[wPhaseIdx] + source[nPhaseIdx] / density_[nPhaseIdx]);
1315 this->f_[eIdxGlobal1] += evaluateErrorTerm_(cellData1) * volume1 / (4.0);
1316 this->f_[eIdxGlobal2] += evaluateErrorTerm_(cellData2) * volume2 / (4.0);
1317 this->f_[eIdxGlobal3] += evaluateErrorTerm_(cellData3) * volume3 / (4.0);
1318 this->f_[eIdxGlobal4] += evaluateErrorTerm_(cellData4) * volume4 / (4.0);
1321 Dune::FieldVector<Scalar, numPhases> lambda1(cellData1.mobility(wPhaseIdx));
1322 lambda1[nPhaseIdx] = cellData1.mobility(nPhaseIdx);
1325 Scalar lambdaTotal1 = lambda1[wPhaseIdx] + lambda1[nPhaseIdx];
1328 Dune::FieldVector<Scalar, numPhases> lambda2(cellData2.mobility(wPhaseIdx));
1329 lambda2[nPhaseIdx] = cellData2.mobility(nPhaseIdx);
1332 Scalar lambdaTotal2 = lambda2[wPhaseIdx] + lambda2[nPhaseIdx];
1335 Dune::FieldVector<Scalar, numPhases> lambda3(cellData3.mobility(wPhaseIdx));
1336 lambda3[nPhaseIdx] = cellData3.mobility(nPhaseIdx);
1339 Scalar lambdaTotal3 = lambda3[wPhaseIdx] + lambda3[nPhaseIdx];
1342 Dune::FieldVector<Scalar, numPhases> lambda4(cellData4.mobility(wPhaseIdx));
1343 lambda4[nPhaseIdx] = cellData4.mobility(nPhaseIdx);
1346 Scalar lambdaTotal4 = lambda4[wPhaseIdx] + lambda4[nPhaseIdx];
1348 Scalar gn12nu14 = interactionVolume.getNtkrkNu_df(lambdaTotal1, 0, 0, 1);
1349 Scalar gn12nu12 = interactionVolume.getNtkrkNu_df(lambdaTotal1, 0, 0, 0);
1350 Scalar gn14nu14 = interactionVolume.getNtkrkNu_df(lambdaTotal1, 0, 1, 1);
1351 Scalar gn14nu12 = interactionVolume.getNtkrkNu_df(lambdaTotal1, 0, 1, 0);
1352 Scalar gn12nu23 = interactionVolume.getNtkrkNu_df(lambdaTotal2, 1, 1, 0);
1353 Scalar gn12nu21 = interactionVolume.getNtkrkNu_df(lambdaTotal2, 1, 1, 1);
1354 Scalar gn23nu23 = interactionVolume.getNtkrkNu_df(lambdaTotal2, 1, 0, 0);
1355 Scalar gn23nu21 = interactionVolume.getNtkrkNu_df(lambdaTotal2, 1, 0, 1);
1356 Scalar gn43nu32 = interactionVolume.getNtkrkNu_df(lambdaTotal3, 2, 0, 1);
1357 Scalar gn43nu34 = interactionVolume.getNtkrkNu_df(lambdaTotal3, 2, 0, 0);
1358 Scalar gn23nu32 = interactionVolume.getNtkrkNu_df(lambdaTotal3, 2, 1, 1);
1359 Scalar gn23nu34 = interactionVolume.getNtkrkNu_df(lambdaTotal3, 2, 1, 0);
1360 Scalar gn43nu41 = interactionVolume.getNtkrkNu_df(lambdaTotal4, 3, 1, 0);
1361 Scalar gn43nu43 = interactionVolume.getNtkrkNu_df(lambdaTotal4, 3, 1, 1);
1362 Scalar gn14nu41 = interactionVolume.getNtkrkNu_df(lambdaTotal4, 3, 0, 0);
1363 Scalar gn14nu43 = interactionVolume.getNtkrkNu_df(lambdaTotal4, 3, 0, 1);
1366 Dune::FieldMatrix<Scalar, 2 * dim, 2 * dim> C(0), F(0), A(0), B(0);
1369 C[0][0] = -gn12nu12;
1370 C[0][3] = -gn12nu14;
1372 C[1][1] = -gn23nu23;
1375 C[3][2] = -gn14nu43;
1378 F[0][0] = gn12nu12 + gn12nu14;
1379 F[1][1] = -gn23nu21 + gn23nu23;
1380 F[2][2] = -gn43nu34 - gn43nu32;
1381 F[3][3] = gn14nu43 - gn14nu41;
1383 A[0][0] = gn12nu12 + gn12nu21;
1384 A[0][1] = -gn12nu23;
1386 A[1][0] = -gn23nu21;
1387 A[1][1] = gn23nu23 + gn23nu32;
1389 A[2][1] = -gn43nu32;
1390 A[2][2] = -gn43nu34 - gn43nu43;
1392 A[3][0] = -gn14nu12;
1394 A[3][3] = -gn14nu41 - gn14nu14;
1396 B[0][0] = gn12nu12 + gn12nu14;
1397 B[0][1] = gn12nu21 - gn12nu23;
1398 B[1][1] = -gn23nu21 + gn23nu23;
1399 B[1][2] = gn23nu34 + gn23nu32;
1400 B[2][2] = -gn43nu34 - gn43nu32;
1401 B[2][3] = -gn43nu43 + gn43nu41;
1402 B[3][0] = -gn14nu12 - gn14nu14;
1403 B[3][3] = gn14nu43 - gn14nu41;
1407 F += C.rightmultiply(B.leftmultiply(A));
1408 Dune::FieldMatrix<Scalar, 2 * dim, 2 * dim> T(F);
1413 this->A_[eIdxGlobal1][eIdxGlobal1] += T[0][0] + T[3][0];
1414 this->A_[eIdxGlobal1][eIdxGlobal2] += T[0][1] + T[3][1];
1415 this->A_[eIdxGlobal1][eIdxGlobal3] += T[0][2] + T[3][2];
1416 this->A_[eIdxGlobal1][eIdxGlobal4] += T[0][3] + T[3][3];
1418 this->A_[eIdxGlobal2][eIdxGlobal1] += -T[0][0] + T[1][0];
1419 this->A_[eIdxGlobal2][eIdxGlobal2] += -T[0][1] + T[1][1];
1420 this->A_[eIdxGlobal2][eIdxGlobal3] += -T[0][2] + T[1][2];
1421 this->A_[eIdxGlobal2][eIdxGlobal4] += -T[0][3] + T[1][3];
1423 this->A_[eIdxGlobal3][eIdxGlobal1] -= T[1][0] + T[2][0];
1424 this->A_[eIdxGlobal3][eIdxGlobal2] -= T[1][1] + T[2][1];
1425 this->A_[eIdxGlobal3][eIdxGlobal3] -= T[1][2] + T[2][2];
1426 this->A_[eIdxGlobal3][eIdxGlobal4] -= T[1][3] + T[2][3];
1428 this->A_[eIdxGlobal4][eIdxGlobal1] += T[2][0] - T[3][0];
1429 this->A_[eIdxGlobal4][eIdxGlobal2] += T[2][1] - T[3][1];
1430 this->A_[eIdxGlobal4][eIdxGlobal3] += T[2][2] - T[3][2];
1431 this->A_[eIdxGlobal4][eIdxGlobal4] += T[2][3] - T[3][3];
1433 if (innerBoundaryVolumeFaces_[eIdxGlobal1][interactionVolume.getIndexOnElement(0, 0)])
1435 this->A_[eIdxGlobal1][eIdxGlobal1] += T[0][0];
1436 this->A_[eIdxGlobal1][eIdxGlobal2] += T[0][1];
1437 this->A_[eIdxGlobal1][eIdxGlobal3] += T[0][2];
1438 this->A_[eIdxGlobal1][eIdxGlobal4] += T[0][3];
1440 if (innerBoundaryVolumeFaces_[eIdxGlobal1][interactionVolume.getIndexOnElement(0, 1)])
1442 this->A_[eIdxGlobal1][eIdxGlobal1] += T[3][0];
1443 this->A_[eIdxGlobal1][eIdxGlobal2] += T[3][1];
1444 this->A_[eIdxGlobal1][eIdxGlobal3] += T[3][2];
1445 this->A_[eIdxGlobal1][eIdxGlobal4] += T[3][3];
1448 if (innerBoundaryVolumeFaces_[eIdxGlobal2][interactionVolume.getIndexOnElement(1, 0)])
1450 this->A_[eIdxGlobal2][eIdxGlobal1] += T[1][0];
1451 this->A_[eIdxGlobal2][eIdxGlobal2] += T[1][1];
1452 this->A_[eIdxGlobal2][eIdxGlobal3] += T[1][2];
1453 this->A_[eIdxGlobal2][eIdxGlobal4] += T[1][3];
1455 if (innerBoundaryVolumeFaces_[eIdxGlobal2][interactionVolume.getIndexOnElement(1, 1)])
1457 this->A_[eIdxGlobal2][eIdxGlobal1] += -T[0][0];
1458 this->A_[eIdxGlobal2][eIdxGlobal2] += -T[0][1];
1459 this->A_[eIdxGlobal2][eIdxGlobal3] += -T[0][2];
1460 this->A_[eIdxGlobal2][eIdxGlobal4] += -T[0][3];
1462 if (innerBoundaryVolumeFaces_[eIdxGlobal3][interactionVolume.getIndexOnElement(2, 0)])
1464 this->A_[eIdxGlobal3][eIdxGlobal1] -= T[2][0];
1465 this->A_[eIdxGlobal3][eIdxGlobal2] -= T[2][1];
1466 this->A_[eIdxGlobal3][eIdxGlobal3] -= T[2][2];
1467 this->A_[eIdxGlobal3][eIdxGlobal4] -= T[2][3];
1469 if (innerBoundaryVolumeFaces_[eIdxGlobal3][interactionVolume.getIndexOnElement(2, 1)])
1471 this->A_[eIdxGlobal3][eIdxGlobal1] -= T[1][0];
1472 this->A_[eIdxGlobal3][eIdxGlobal2] -= T[1][1];
1473 this->A_[eIdxGlobal3][eIdxGlobal3] -= T[1][2];
1474 this->A_[eIdxGlobal3][eIdxGlobal4] -= T[1][3];
1476 if (innerBoundaryVolumeFaces_[eIdxGlobal4][interactionVolume.getIndexOnElement(3, 0)])
1478 this->A_[eIdxGlobal4][eIdxGlobal1] += -T[3][0];
1479 this->A_[eIdxGlobal4][eIdxGlobal2] += -T[3][1];
1480 this->A_[eIdxGlobal4][eIdxGlobal3] += -T[3][2];
1481 this->A_[eIdxGlobal4][eIdxGlobal4] += -T[3][3];
1483 if (innerBoundaryVolumeFaces_[eIdxGlobal4][interactionVolume.getIndexOnElement(3, 1)])
1485 this->A_[eIdxGlobal4][eIdxGlobal1] += T[2][0];
1486 this->A_[eIdxGlobal4][eIdxGlobal2] += T[2][1];
1487 this->A_[eIdxGlobal4][eIdxGlobal3] += T[2][2];
1488 this->A_[eIdxGlobal4][eIdxGlobal4] += T[2][3];
1493 Dune::FieldVector<Scalar, 2 * dim> pc(0);
1494 pc[0] = cellData1.capillaryPressure();
1495 pc[1] = cellData2.capillaryPressure();
1496 pc[2] = cellData3.capillaryPressure();
1497 pc[3] = cellData4.capillaryPressure();
1499 Dune::FieldVector<Scalar, 2 * dim> gravityDiff(0);
1503 gravityDiff[0] = (problem_.bBoxMax() - globalPos1) * gravity_ * (density_[nPhaseIdx] - density_[wPhaseIdx]);
1504 gravityDiff[1] = (problem_.bBoxMax() - globalPos2) * gravity_ * (density_[nPhaseIdx] - density_[wPhaseIdx]);
1505 gravityDiff[2] = (problem_.bBoxMax() - globalPos3) * gravity_ * (density_[nPhaseIdx] - density_[wPhaseIdx]);
1506 gravityDiff[3] = (problem_.bBoxMax() - globalPos4) * gravity_ * (density_[nPhaseIdx] - density_[wPhaseIdx]);
1510 if (pc[0] == 0 && pc[1] == 0 && pc[2] == 0 && pc[3] == 0)
1515 Dune::FieldVector<Scalar, 2 * dim> pcFlux(0);
1521 Scalar pcPotential12 = pcFlux[0];
1522 Scalar pcPotential14 = pcFlux[3];
1523 Scalar pcPotential32 = -pcFlux[1];
1524 Scalar pcPotential34 = -pcFlux[2];
1529 Dune::FieldVector<Scalar, numPhases> lambda12Upw(0.0);
1530 lambda12Upw[wPhaseIdx] = (pcPotential12 >= 0) ? lambda1[wPhaseIdx] : lambda2[wPhaseIdx];
1531 lambda12Upw[nPhaseIdx] = (pcPotential12 >= 0) ? lambda1[nPhaseIdx] : lambda2[nPhaseIdx];
1534 Dune::FieldVector<Scalar, numPhases> lambda14Upw(0.0);
1535 lambda14Upw[wPhaseIdx] = (pcPotential14 >= 0) ? lambda1[wPhaseIdx] : lambda4[wPhaseIdx];
1536 lambda14Upw[nPhaseIdx] = (pcPotential14 >= 0) ? lambda1[nPhaseIdx] : lambda4[nPhaseIdx];
1539 Dune::FieldVector<Scalar, numPhases> lambda32Upw(0.0);
1540 lambda32Upw[wPhaseIdx] = (pcPotential32 >= 0) ? lambda3[wPhaseIdx] : lambda2[wPhaseIdx];
1541 lambda32Upw[nPhaseIdx] = (pcPotential32 >= 0) ? lambda3[nPhaseIdx] : lambda2[nPhaseIdx];
1544 Dune::FieldVector<Scalar, numPhases> lambda34Upw(0.0);
1545 lambda34Upw[wPhaseIdx] = (pcPotential34 >= 0) ? lambda3[wPhaseIdx] : lambda4[wPhaseIdx];
1546 lambda34Upw[nPhaseIdx] = (pcPotential34 >= 0) ? lambda3[nPhaseIdx] : lambda4[nPhaseIdx];
1548 for (
int i = 0; i < numPhases; i++)
1550 Scalar lambdaT12 = lambda12Upw[wPhaseIdx] + lambda12Upw[nPhaseIdx];
1551 Scalar lambdaT14 = lambda14Upw[wPhaseIdx] + lambda14Upw[nPhaseIdx];
1552 Scalar lambdaT32 = lambda32Upw[wPhaseIdx] + lambda32Upw[nPhaseIdx];
1553 Scalar lambdaT34 = lambda34Upw[wPhaseIdx] + lambda34Upw[nPhaseIdx];
1554 Scalar fracFlow12 = (lambdaT12 > threshold_) ? lambda12Upw[i] / (lambdaT12) : 0.0;
1555 Scalar fracFlow14 = (lambdaT14 > threshold_) ? lambda14Upw[i] / (lambdaT14) : 0.0;
1556 Scalar fracFlow32 = (lambdaT32 > threshold_) ? lambda32Upw[i] / (lambdaT32) : 0.0;
1557 Scalar fracFlow34 = (lambdaT34 > threshold_) ? lambda34Upw[i] / (lambdaT34) : 0.0;
1559 Dune::FieldVector<Scalar, 2 * dim> pcFluxReal(pcFlux);
1561 pcFluxReal[0] *= fracFlow12;
1562 pcFluxReal[1] *= fracFlow32;
1563 pcFluxReal[2] *= fracFlow34;
1564 pcFluxReal[3] *= fracFlow14;
1570 switch (pressureType_)
1577 this->f_[eIdxGlobal1] -= (pcFluxReal[0] + pcFluxReal[3]);
1578 this->f_[eIdxGlobal2] -= (pcFluxReal[1] - pcFluxReal[0]);
1579 this->f_[eIdxGlobal3] -= (-pcFluxReal[2] - pcFluxReal[1]);
1580 this->f_[eIdxGlobal4] -= (-pcFluxReal[3] + pcFluxReal[2]);
1582 if (innerBoundaryVolumeFaces_[eIdxGlobal1][interactionVolume.getIndexOnElement(0, 0)])
1584 this->f_[eIdxGlobal1] -= pcFluxReal[0];
1586 if (innerBoundaryVolumeFaces_[eIdxGlobal1][interactionVolume.getIndexOnElement(0, 1)])
1588 this->f_[eIdxGlobal1] -= pcFluxReal[3];
1590 if (innerBoundaryVolumeFaces_[eIdxGlobal2][interactionVolume.getIndexOnElement(1, 0)])
1592 this->f_[eIdxGlobal2] -= pcFluxReal[1];
1594 if (innerBoundaryVolumeFaces_[eIdxGlobal2][interactionVolume.getIndexOnElement(1, 1)])
1596 this->f_[eIdxGlobal2] += pcFluxReal[0];
1598 if (innerBoundaryVolumeFaces_[eIdxGlobal3][interactionVolume.getIndexOnElement(2, 0)])
1600 this->f_[eIdxGlobal3] += pcFluxReal[2];
1602 if (innerBoundaryVolumeFaces_[eIdxGlobal3][interactionVolume.getIndexOnElement(2, 1)])
1604 this->f_[eIdxGlobal3] += pcFluxReal[1];
1606 if (innerBoundaryVolumeFaces_[eIdxGlobal4][interactionVolume.getIndexOnElement(3, 0)])
1608 this->f_[eIdxGlobal4] += pcFluxReal[3];
1610 if (innerBoundaryVolumeFaces_[eIdxGlobal4][interactionVolume.getIndexOnElement(3, 1)])
1612 this->f_[eIdxGlobal4] -= pcFluxReal[2];
1622 this->f_[eIdxGlobal1] += (pcFluxReal[0] + pcFluxReal[1]);
1623 this->f_[eIdxGlobal2] += (pcFluxReal[1] - pcFluxReal[0]);
1624 this->f_[eIdxGlobal3] += (-pcFluxReal[2] - pcFluxReal[1]);
1625 this->f_[eIdxGlobal4] += (-pcFluxReal[3] + pcFluxReal[2]);
1627 if (innerBoundaryVolumeFaces_[eIdxGlobal1][interactionVolume.getIndexOnElement(0, 0)])
1629 this->f_[eIdxGlobal1] += pcFluxReal[0];
1631 if (innerBoundaryVolumeFaces_[eIdxGlobal1][interactionVolume.getIndexOnElement(0, 1)])
1633 this->f_[eIdxGlobal1] += pcFluxReal[3];
1635 if (innerBoundaryVolumeFaces_[eIdxGlobal2][interactionVolume.getIndexOnElement(1, 0)])
1637 this->f_[eIdxGlobal2] += pcFluxReal[1];
1639 if (innerBoundaryVolumeFaces_[eIdxGlobal2][interactionVolume.getIndexOnElement(1, 1)])
1641 this->f_[eIdxGlobal2] -= pcFluxReal[0];
1643 if (innerBoundaryVolumeFaces_[eIdxGlobal3][interactionVolume.getIndexOnElement(2, 0)])
1645 this->f_[eIdxGlobal3] -= pcFluxReal[2];
1647 if (innerBoundaryVolumeFaces_[eIdxGlobal3][interactionVolume.getIndexOnElement(2, 1)])
1649 this->f_[eIdxGlobal3] -= pcFluxReal[1];
1651 if (innerBoundaryVolumeFaces_[eIdxGlobal4][interactionVolume.getIndexOnElement(3, 0)])
1653 this->f_[eIdxGlobal4] -= pcFluxReal[3];
1655 if (innerBoundaryVolumeFaces_[eIdxGlobal4][interactionVolume.getIndexOnElement(3, 1)])
1657 this->f_[eIdxGlobal4] += pcFluxReal[2];
1669 for (
int elemIdx = 0; elemIdx < 2 * dim; elemIdx++)
1671 bool isOutside =
false;
1672 for (
int fIdx = 0; fIdx < dim; fIdx++)
1674 int intVolFaceIdx = interactionVolume.getFaceIndexFromSubVolume(elemIdx, fIdx);
1675 if (interactionVolume.isOutsideFace(intVolFaceIdx))
1686 auto element = interactionVolume.getSubVolumeElement(elemIdx);
1689 const GlobalPosition& globalPos = element.geometry().center();
1692 Scalar volume = element.geometry().volume();
1695 int eIdxGlobal = problem_.variables().index(element);
1698 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
1701 DimMatrix
permeability(problem_.spatialParams().intrinsicPermeability(element));
1704 PrimaryVariables source(0);
1705 problem_.source(source, element);
1706 this->f_[eIdxGlobal] += volume / (4.0)
1707 * (source[wPhaseIdx] / density_[wPhaseIdx] + source[nPhaseIdx] / density_[nPhaseIdx]);
1709 this->f_[eIdxGlobal] += evaluateErrorTerm_(cellData) * volume / (4.0);
1712 Dune::FieldVector<Scalar, numPhases> lambda(cellData.mobility(wPhaseIdx));
1713 lambda[nPhaseIdx] = cellData.mobility(nPhaseIdx);
1715 Scalar pc = cellData.capillaryPressure();
1717 Scalar gravityDiff = (problem_.bBoxMax() - globalPos) * gravity_
1718 * (density_[nPhaseIdx] - density_[wPhaseIdx]);
1722 for (
int fIdx = 0; fIdx < dim; fIdx++)
1724 int intVolFaceIdx = interactionVolume.getFaceIndexFromSubVolume(elemIdx, fIdx);
1726 if (interactionVolume.isBoundaryFace(intVolFaceIdx))
1729 if (interactionVolume.getBoundaryType(intVolFaceIdx).isDirichlet(pressEqIdx))
1731 int boundaryFaceIdx = interactionVolume.getIndexOnElement(elemIdx, fIdx);
1733 const auto referenceElement = ReferenceElements::general(element.type());
1735 const LocalPosition& localPos = referenceElement.position(boundaryFaceIdx, 1);
1737 const GlobalPosition& globalPosFace = element.geometry().global(localPos);
1739 DimVector distVec(globalPosFace - globalPos);
1740 Scalar dist = distVec.two_norm();
1741 DimVector unitDistVec(distVec);
1742 unitDistVec /= dist;
1744 Scalar faceArea = interactionVolume.getFaceArea(elemIdx, fIdx);
1747 Scalar satWBound = cellData.saturation(wPhaseIdx);
1749 if (interactionVolume.getBoundaryType(intVolFaceIdx).isDirichlet(satEqIdx))
1751 Scalar satBound = interactionVolume.getDirichletValues(intVolFaceIdx)[saturationIdx];
1752 switch (saturationType_)
1756 satWBound = satBound;
1761 satWBound = 1 - satBound;
1768 Scalar pcBound = MaterialLaw::pc(
1769 problem_.spatialParams().materialLawParams(element), satWBound);
1771 Scalar gravityDiffBound = (problem_.bBoxMax() - globalPosFace) * gravity_
1772 * (density_[nPhaseIdx] - density_[wPhaseIdx]);
1774 pcBound += gravityDiffBound;
1776 Dune::FieldVector<Scalar, numPhases> lambdaBound(
1777 MaterialLaw::krw(problem_.spatialParams().materialLawParams(element),
1779 lambdaBound[nPhaseIdx] = MaterialLaw::krn(
1780 problem_.spatialParams().materialLawParams(element), satWBound);
1781 lambdaBound[wPhaseIdx] /= viscosity_[wPhaseIdx];
1782 lambdaBound[nPhaseIdx] /= viscosity_[nPhaseIdx];
1784 Scalar potentialBound = interactionVolume.getDirichletValues(intVolFaceIdx)[pressureIdx];
1785 Scalar gdeltaZ = (problem_.bBoxMax()-globalPosFace) * gravity_;
1788 Scalar potentialDiffW = 0;
1789 Scalar potentialDiffNw = 0;
1790 switch (pressureType_)
1794 potentialBound += density_[wPhaseIdx]*gdeltaZ;
1795 potentialDiffW = (cellData.potential(wPhaseIdx) - potentialBound) / dist;
1796 potentialDiffNw = (cellData.potential(nPhaseIdx) - potentialBound - pcBound) / dist;
1801 potentialBound += density_[nPhaseIdx]*gdeltaZ;
1802 potentialDiffW = (cellData.potential(wPhaseIdx) - potentialBound + pcBound) / dist;
1803 potentialDiffNw = (cellData.potential(nPhaseIdx) - potentialBound) / dist;
1808 Scalar lambdaTotal = (potentialDiffW >= 0.) ? lambda[wPhaseIdx] : lambdaBound[wPhaseIdx];
1809 lambdaTotal += (potentialDiffNw >= 0.) ? lambda[nPhaseIdx] : lambdaBound[nPhaseIdx];
1811 DimVector permTimesNormal(0);
1815 Scalar entry = lambdaTotal * (unitDistVec * permTimesNormal) / dist * faceArea;
1820 switch (pressureType_)
1825 DimVector pcGradient = unitDistVec;
1826 pcGradient *= (pc - pcBound) / dist;
1829 pcFlux = 0.5 * (lambda[nPhaseIdx] + lambdaBound[nPhaseIdx])
1830 * (permTimesNormal * pcGradient) * faceArea;
1837 DimVector pcGradient = unitDistVec;
1838 pcGradient *= (pc - pcBound) / dist;
1841 pcFlux = 0.5 * (lambda[wPhaseIdx] + lambdaBound[wPhaseIdx])
1842 * (permTimesNormal * pcGradient) * faceArea;
1850 this->A_[eIdxGlobal][eIdxGlobal] += entry;
1851 this->f_[eIdxGlobal] += entry * potentialBound;
1853 if (pc == 0 && pcBound == 0)
1858 for (
int i = 0; i < numPhases; i++)
1860 switch (pressureType_)
1867 this->f_[eIdxGlobal] -= pcFlux;
1876 this->f_[eIdxGlobal] += pcFlux;
1884 else if (interactionVolume.getBoundaryType(intVolFaceIdx).isNeumann(pressEqIdx))
1886 Scalar J = interactionVolume.getNeumannValues(intVolFaceIdx)[wPhaseIdx]
1887 / density_[wPhaseIdx];
1888 J += interactionVolume.getNeumannValues(intVolFaceIdx)[nPhaseIdx] / density_[nPhaseIdx];
1889 this->f_[eIdxGlobal] -= J;
1893 std::cout <<
"interactionVolume.getBoundaryType(intVolFaceIdx).isNeumann(pressEqIdx)"
1894 << interactionVolume.getBoundaryType(intVolFaceIdx).isNeumann(pressEqIdx) <<
"\n";
1895 DUNE_THROW(Dune::NotImplemented,
1896 "No valid boundary condition type defined for pressure equation!");
1907 if (problem_.gridView().comm().size() > 1)
1910 for (
const auto& element : elements(problem_.gridView()))
1912 if (element.partitionType() == Dune::InteriorEntity)
1916 int eIdxGlobalI = problem_.variables().index(element);
1918 this->A_[eIdxGlobalI] = 0.0;
1919 this->A_[eIdxGlobalI][eIdxGlobalI] = 1.0;
1920 this->f_[eIdxGlobalI] = this->
pressure()[eIdxGlobalI];
1933template<
class TypeTag>
1937 for (
const auto& element : elements(problem_.gridView()))
1939 int eIdxGlobal = problem_.variables().index(element);
1941 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
1943 Scalar satW = cellData.saturation(wPhaseIdx);
1945 Scalar pc = MaterialLaw::pc(problem_.spatialParams().materialLawParams(element), satW);
1947 cellData.setCapillaryPressure(pc);
1950 Scalar mobilityW = MaterialLaw::krw(problem_.spatialParams().materialLawParams(element), satW)
1951 / viscosity_[wPhaseIdx];
1952 Scalar mobilityNw = MaterialLaw::krn(problem_.spatialParams().materialLawParams(element), satW)
1953 / viscosity_[nPhaseIdx];
1956 cellData.setMobility(wPhaseIdx, mobilityW);
1957 cellData.setMobility(nPhaseIdx, mobilityNw);
1960 cellData.setFracFlowFunc(wPhaseIdx, mobilityW / (mobilityW + mobilityNw));
1961 cellData.setFracFlowFunc(nPhaseIdx, mobilityNw / (mobilityW + mobilityNw));
#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
Class including the information of an interaction volume of a MPFA O-method that does not change with...
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 CellData
Defines data object to be stored.
Definition: porousmediumflow/sequential/properties.hh:72
Property tag GridImplementation
Gives kind of grid implementation in form of a GridType.
Definition: porousmediumflow/sequential/cellcentered/mpfa/properties.hh:95
Property tag EnableCompressibility
Returns whether compressibility is allowed.
Definition: porousmediumflow/2p/sequential/properties.hh:55
Property tag SaturationFormulation
The formulation of the saturation model.
Definition: porousmediumflow/2p/sequential/properties.hh:53
Property tag NumPhases
Number of phases in the system.
Definition: porousmediumflow/sequential/properties.hh:69
Property tag VelocityFormulation
The type of velocity reconstructed for the transport model.
Definition: porousmediumflow/2p/sequential/properties.hh:54
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 O-method discretization of a two-phase flow pressure equation of the sequential IM...
Definition: omethod/2dpressure.hh:67
void updateMaterialLaws()
Constitutive functions are initialized and stored in the variables object.
Definition: omethod/2dpressure.hh:1934
GlobalInteractionVolumeVector interactionVolumes_
Definition: omethod/2dpressure.hh:483
void storePressureSolution(const Element &element)
Stores the pressure solution of a cell.
Definition: omethod/2dpressure.hh:268
InnerBoundaryVolumeFaces innerBoundaryVolumeFaces_
Definition: omethod/2dpressure.hh:484
void updateInteractionVolumeInfo()
Updates interaction volumes.
Definition: omethod/2dpressure.hh:161
void addOutputVtkFields(MultiWriter &writer)
Adds pressure output to the output file.
Definition: omethod/2dpressure.hh:327
void update()
Pressure update.
Definition: omethod/2dpressure.hh:213
void initialize()
Initializes the pressure model.
Definition: omethod/2dpressure.hh:177
void storePressureSolution()
Globally stores the pressure solution.
Definition: omethod/2dpressure.hh:254
FvMpfaO2dPressure2p(Problem &problem)
Constructs a FvMpfaO2dPressure2p object.
Definition: omethod/2dpressure.hh:397
Class including the information of an interaction volume of a MPFA O-method that does not change with...
Definition: ointeractionvolume.hh:38
The finite volume base class for the solution of a pressure equation.
Definition: sequential/cellcentered/pressure.hh:48
void initialize()
Initialize pressure model.
Definition: sequential/cellcentered/pressure.hh:212
PressureSolution & pressure()
Returns the vector containing the pressure solution.
Definition: sequential/cellcentered/pressure.hh:119
void solve()
Solves the global system of equations to get the spatial distribution of the pressure.
Definition: sequential/cellcentered/pressure.hh:526
Specifies the properties for immiscible 2p diffusion/pressure models.
Properties for a MPFA method.
Finite Volume Diffusion Model.