24#ifndef DUMUX_FV3DPRESSURE2P2C_ADAPTIVE_HH
25#define DUMUX_FV3DPRESSURE2P2C_ADAPTIVE_HH
28#include <dune/istl/bvector.hh>
29#include <dune/istl/operators.hh>
30#include <dune/istl/solvers.hh>
31#include <dune/istl/preconditioners.hh>
47template<
class TypeTag>
49template<
class TypeTag>
103 dim = GridView::dimension, dimWorld = GridView::dimensionworld,
104 NumPhases = getPropValue<TypeTag, Properties::NumPhases>(), NumComponents = getPropValue<TypeTag, Properties::NumComponents>()
108 pw = Indices::pressureW,
109 pn = Indices::pressureN,
110 pGlobal = Indices::pressureGlobal,
111 Sw = Indices::saturationW,
112 Sn = Indices::saturationN
116 wPhaseIdx = Indices::wPhaseIdx, nPhaseIdx = Indices::nPhaseIdx,
117 wCompIdx = Indices::wPhaseIdx, nCompIdx = Indices::nPhaseIdx,
118 contiWEqIdx = Indices::contiWEqIdx, contiNEqIdx = Indices::contiNEqIdx
126 using Vertex =
typename GridView::Traits::template Codim<dim>::Entity;
127 using Element =
typename GridView::Traits::template Codim<0>::Entity;
129 using Grid =
typename GridView::Grid;
130 using Intersection =
typename GridView::Intersection;
131 using IntersectionIterator =
typename GridView::IntersectionIterator;
134 using GlobalPosition = Dune::FieldVector<Scalar, dimWorld>;
135 using TransmissivityMatrix = Dune::FieldVector<Scalar,dim+1>;
136 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
137 using PhaseVector = Dune::FieldVector<Scalar, getPropValue<TypeTag, Properties::NumPhases>()>;
138 using ComponentVector = Dune::FieldVector<Scalar, getPropValue<TypeTag, Properties::NumComponents>()>;
147 using InteractionVolume =
typename InteractionVolumeContainer::InteractionVolume;
155 const Problem& problem()
const
170 new InteractionVolumeContainer(problem());
174 asImp_().initializeMatrix();
187 new InteractionVolumeContainer(problem());
197 void getMpfaFlux(
const IntersectionIterator&,
const CellData&);
199 void get1pMpfaFlux(
const IntersectionIterator&,
const CellData&);
206 TransmissivityMatrix&,
215 int gridSize = problem().gridView().size(0);
218 for (
const auto& element : elements(problem().gridView()))
221 int eIdxGlobalI = problem().variables().index(element);
224 if (element.partitionType() == Dune::InteriorEntity)
227 = problem().variables().cellData(eIdxGlobalI).pressure(this->
pressureType);
238 problem().gridView().template communicate<DataHandle>(dataHandle,
239 Dune::InteriorBorder_All_Interface,
240 Dune::ForwardCommunication);
253 enableMPFA = getParam<bool>(
"GridAdapt.EnableMultiPointFluxApproximation");
260 Implementation &asImp_()
261 {
return *
static_cast<Implementation *
>(
this);}
264 const Implementation &asImp_()
const
265 {
return *
static_cast<const Implementation *
>(
this);}
267 int searchCommonVertex_(
const Intersection& is, Vertex& vertex)
271 int localIdxLarge = 0;
272 for(localIdxLarge = 0; localIdxLarge < is.inside().subEntities(dim); ++localIdxLarge)
274 auto vLarge = is.inside().template subEntity<dim>(localIdxLarge);
277 for(
int verticeSmall = 0; verticeSmall<is.outside().subEntities(dim); ++verticeSmall)
279 auto vSmall = is.outside().template subEntity<dim>(verticeSmall);
281 if(problem().variables().index(vSmall) == problem().variables().index(vLarge) )
284 return localIdxLarge;
293 InteractionVolume& interactionVolume,
294 const int& subVolumeFaceIdx,
295 bool properFluxDirection,
296 Element& additional2,
297 Element& additional3,
298 TransmissivityMatrix& additionalT);
312template<
class TypeTag>
315 int gridSize_ = problem().gridView().size(0);
317 this->A_.setSize (gridSize_,gridSize_);
318 this->f_.resize(gridSize_);
319 irregularCellMap_.clear();
321 this->initializeMatrixRowSize();
322 this->A_.endrowsizes();
323 this->initializeMatrixIndices();
324 this->A_.endindices();
328template<
class TypeTag>
332 for (
const auto& element : elements(problem().gridView()))
335 int eIdxGlobalI = problem().variables().index(element);
336 CellData& cellDataI = problem().variables().cellData(eIdxGlobalI);
342 cellDataI.perimeter() = 0;
345 std::vector<int> foundAdditionals;
347 int numberOfIntersections = 0;
349 for (
const auto& intersection : intersections(problem().gridView(), element))
351 cellDataI.perimeter() += intersection.geometry().volume();
352 numberOfIntersections++;
353 if (intersection.neighbor())
358 if (enableMPFA && (element.level() < intersection.outside().level()))
362 int intersectionID = problem().grid().localIdSet().subId(element,
363 intersection.indexInInside(), 1);
365 int eIdxGlobalJ = problem().variables().index(intersection.outside());
368 irregularCellMap_[intersectionID].push_back(eIdxGlobalJ);
372 cellDataI.fluxData().resize(numberOfIntersections);
373 this->A_.incrementrowsize(eIdxGlobalI, rowSize);
377 if(enableMPFA && maxInteractionVolumes>1)
380 std::multimap<int, int> addionalRelations;
381 using IntPair = std::pair<int,int>;
382 std::pair<std::multimap<int,int>::iterator,std::multimap<int,int>::iterator> range;
383 std::multimap<int,int>::iterator rangeIt;
386 for (
const auto& element : elements(problem().gridView()))
389 int eIdxGlobalI = problem().variables().index(element);
391 const auto isEndIt = problem().gridView().iend(element);
392 for (
auto isIt = problem().gridView().ibegin(element); isIt != isEndIt; ++isIt)
394 const auto& intersection = *isIt;
396 if (intersection.neighbor())
399 int eIdxGlobalJ = problem().variables().index(intersection.outside());
402 if (intersection.outside().level() > element.level())
406 GlobalPosition globalPos3(0.);
408 GlobalPosition globalPos4(0.);
410 TransmissivityMatrix T(0.);
412 int interactionRegions
413 = problem().variables().getMpfaData3D(intersection, T,
414 globalPos3, eIdxGlobal3, globalPos4, eIdxGlobal4 );
415 if (interactionRegions == 0)
416 interactionRegions = problem().pressureModel().computeTransmissibilities(isIt,T,
417 globalPos3, eIdxGlobal3, globalPos4, eIdxGlobal4 );
418 if(!interactionRegions)
419 Dune::dgrave <<
"something went wrong getting mpfa data on cell " << eIdxGlobalI << std::endl;
420 if (interactionRegions == 1)
424 for (
int cocumber=1; cocumber<interactionRegions; cocumber++ )
426 problem().variables().getMpfaData3D(intersection, T,
427 globalPos3, eIdxGlobal3, globalPos4, eIdxGlobal4, cocumber);
429 int additionalIdx2 = eIdxGlobal3;
430 int additionalIdx3 = eIdxGlobal4;
432 bool addIndex =
true;
435 bool additional2isNeighbor(
false), additional3isNeighbor(
false);
437 for (
const auto& checkIntersection
438 : intersections(problem().gridView(), element))
440 if (checkIntersection.neighbor())
442 if(additionalIdx2==problem().variables().index(checkIntersection.outside()))
443 additional2isNeighbor =
true;
444 if(additionalIdx3 == problem().variables().index(checkIntersection.outside()))
445 additional3isNeighbor =
true;
450 if(!additional2isNeighbor)
453 IntPair intPair(eIdxGlobalI,additionalIdx2);
454 if(eIdxGlobalI > additionalIdx2)
457 swap(intPair.first, intPair.second);
459 range = addionalRelations.equal_range(intPair.first);
460 for (rangeIt=range.first; range.first!=range.second
461 && rangeIt!=range.second; ++rangeIt)
462 if((*rangeIt).second == intPair.second)
466 this->A_.incrementrowsize(eIdxGlobalI);
468 this->A_.incrementrowsize(additionalIdx2);
470 addionalRelations.insert(intPair);
475 if(!additional3isNeighbor)
479 IntPair intPair(eIdxGlobalI,additionalIdx3);
480 if(eIdxGlobalI > additionalIdx3)
483 swap(intPair.first, intPair.second);
485 range = addionalRelations.equal_range(intPair.first);
486 for (rangeIt=range.first; range.first!=range.second
487 && rangeIt!=range.second; ++rangeIt)
488 if((*rangeIt).second == intPair.second)
492 this->A_.incrementrowsize(eIdxGlobalI);
494 this->A_.incrementrowsize(additionalIdx3);
496 addionalRelations.insert(intPair);
501 additional2isNeighbor = additional3isNeighbor =
false;
503 for (
const auto& checkIntersection
504 : intersections(problem().gridView(), intersection.outside()))
506 if (checkIntersection.neighbor())
508 if(additionalIdx2 == problem().variables().index(checkIntersection.outside()))
509 additional2isNeighbor =
true;
510 if(additionalIdx3 == problem().variables().index(checkIntersection.outside()))
511 additional3isNeighbor =
true;
516 if(!additional2isNeighbor)
520 IntPair intPair(eIdxGlobalJ,additionalIdx2);
521 if(eIdxGlobalJ > additionalIdx2)
524 swap(intPair.first, intPair.second);
526 range = addionalRelations.equal_range(intPair.first);
527 for (rangeIt=range.first; range.first!=range.second
528 && rangeIt!=range.second; ++rangeIt)
529 if((*rangeIt).second == intPair.second)
533 this->A_.incrementrowsize(eIdxGlobalJ);
535 this->A_.incrementrowsize(additionalIdx2);
537 addionalRelations.insert(intPair);
542 if(!additional3isNeighbor)
546 IntPair intPair(eIdxGlobalJ,additionalIdx3);
547 if(eIdxGlobalJ > additionalIdx3)
550 swap(intPair.first, intPair.second);
552 range = addionalRelations.equal_range(intPair.first);
553 for (rangeIt=range.first; range.first!=range.second
554 && rangeIt!=range.second; ++rangeIt)
555 if((*rangeIt).second == intPair.second)
559 this->A_.incrementrowsize(eIdxGlobalJ);
561 this->A_.incrementrowsize(additionalIdx3);
563 addionalRelations.insert(intPair);
577template<
class TypeTag>
581 for (
const auto& element : elements(problem().gridView()))
584 int eIdxGlobalI = problem().variables().index(element);
587 this->A_.addindex(eIdxGlobalI, eIdxGlobalI);
590 const auto isEndIt = problem().gridView().iend(element);
591 for (
auto isIt = problem().gridView().ibegin(element); isIt != isEndIt; ++isIt)
593 const auto& intersection = *isIt;
595 if (intersection.neighbor())
598 int eIdxGlobalJ = problem().variables().index(intersection.outside());
601 this->A_.addindex(eIdxGlobalI, eIdxGlobalJ);
604 if (enableMPFA && (element.level() < intersection.outside().level()))
607 GlobalPosition globalPos3(0.);
609 GlobalPosition globalPos4(0.);
611 TransmissivityMatrix T(0.);
612 TransmissivityMatrix additionalT(0.);
614 int interactionRegions
615 = problem().variables().getMpfaData3D(intersection, T, globalPos3, eIdxGlobal3, globalPos4, eIdxGlobal4 );
616 if (interactionRegions == 0)
617 interactionRegions = problem().pressureModel().computeTransmissibilities(isIt,T,
618 globalPos3, eIdxGlobal3, globalPos4, eIdxGlobal4 );
620 for (
int cocumber=1; cocumber<interactionRegions; cocumber++ )
622 problem().variables().getMpfaData3D(intersection, T,
623 globalPos3, eIdxGlobal3, globalPos4, eIdxGlobal4, cocumber);
626 this->A_.addindex(eIdxGlobalI, eIdxGlobal3);
627 this->A_.addindex(eIdxGlobal3, eIdxGlobalI);
628 this->A_.addindex(eIdxGlobalI, eIdxGlobal4);
629 this->A_.addindex(eIdxGlobal4, eIdxGlobalI);
630 this->A_.addindex(eIdxGlobalJ, eIdxGlobal3);
631 this->A_.addindex(eIdxGlobal3, eIdxGlobalJ);
632 this->A_.addindex(eIdxGlobalJ, eIdxGlobal4);
633 this->A_.addindex(eIdxGlobal4, eIdxGlobalJ);
642template<
class TypeTag>
647 BaseType::assemble(
true);
655 for (
const auto& element : elements(problem().gridView()))
658 int eIdxGlobalI = problem().variables().index(element);
661 if (element.partitionType() == Dune::InteriorEntity)
664 CellData& cellDataI = problem().variables().cellData(eIdxGlobalI);
666 Dune::FieldVector<Scalar, 2> entries(0.);
669#ifndef noMultiphysics
670 if(cellDataI.subdomain() != 2)
671 problem().pressureModel().get1pSource(entries,element, cellDataI);
674 problem().pressureModel().getSource(entries,element, cellDataI, first);
676 this->f_[eIdxGlobalI] += entries[rhs];
680 const auto isEndIt = problem().gridView().iend(element);
681 for (
auto isIt = problem().gridView().ibegin(element); isIt != isEndIt; ++isIt)
683 const auto& intersection = *isIt;
686 if (intersection.neighbor())
688 auto neighbor = intersection.outside();
689 int eIdxGlobalJ = problem().variables().index(neighbor);
692 bool haveSameLevel = (element.level() == neighbor.level());
696 if (getPropValue<TypeTag, Properties::VisitFacesOnlyOnce>()
697 && (eIdxGlobalI > eIdxGlobalJ) && haveSameLevel
698 && neighbor.partitionType() == Dune::InteriorEntity)
703 if(!haveSameLevel && enableMPFA)
705 if (cellDataI.subdomain() != 2
706 || problem().variables().cellData(eIdxGlobalJ).subdomain() != 2)
708 asImp_().get1pMpfaFlux(isIt, cellDataI);
712 asImp_().getMpfaFlux(isIt, cellDataI);
717 CellData cellDataJ = problem().variables().cellData(eIdxGlobalJ);
718 if (cellDataI.subdomain() != 2
719 || problem().variables().cellData(eIdxGlobalJ).subdomain() != 2)
721 asImp_().get1pFlux(entries, intersection, cellDataI);
725 asImp_().getFlux(entries, intersection, cellDataI, first);
729 this->f_[eIdxGlobalI] -= entries[rhs];
732 this->A_[eIdxGlobalI][eIdxGlobalI] += entries[matrix];
735 this->A_[eIdxGlobalI][eIdxGlobalJ] -= entries[matrix];
738 if (getPropValue<TypeTag, Properties::VisitFacesOnlyOnce>()
739 && neighbor.partitionType() == Dune::InteriorEntity)
741 this->f_[eIdxGlobalJ] += entries[rhs];
742 this->A_[eIdxGlobalJ][eIdxGlobalJ] += entries[matrix];
743 this->A_[eIdxGlobalJ][eIdxGlobalI] -= entries[matrix];
752 if (cellDataI.subdomain() != 2)
753 asImp_().get1pFluxOnBoundary(entries, intersection, cellDataI);
755 asImp_().getFluxOnBoundary(entries, intersection, cellDataI, first);
758 this->f_[eIdxGlobalI] += entries[rhs];
760 this->A_[eIdxGlobalI][eIdxGlobalI] += entries[matrix];
765 if (cellDataI.subdomain() != 2)
766 asImp_().get1pStorage(entries, element, cellDataI);
768 asImp_().getStorage(entries, element, cellDataI, first);
770 this->f_[eIdxGlobalI] += entries[rhs];
772 this->A_[eIdxGlobalI][eIdxGlobalI] += entries[matrix];
777 this->A_[eIdxGlobalI] = 0.0;
778 this->A_[eIdxGlobalI][eIdxGlobalI] = 1.0;
779 this->f_[eIdxGlobalI] = this->
pressure()[eIdxGlobalI];
803template<
class TypeTag>
805 const CellData& cellDataI)
807 const auto& intersection = *isIt;
810 auto elementI = intersection.inside();
811 int eIdxGlobalI = problem().variables().index(elementI);
814 const GlobalPosition& globalPos = elementI.geometry().center();
817 Scalar volume = elementI.geometry().volume();
818 Scalar perimeter = cellDataI.perimeter();
821 DimMatrix permeabilityI(problem().spatialParams().intrinsicPermeability(elementI));
824 auto neighbor = intersection.outside();
825 int eIdxGlobalJ = problem().variables().index(neighbor);
826 CellData& cellDataJ = problem().variables().cellData(eIdxGlobalJ);
829 const GlobalPosition& globalPosNeighbor = neighbor.geometry().center();
832 GlobalPosition distVec = globalPosNeighbor - globalPos;
835 Scalar dist = distVec.two_norm();
837 GlobalPosition unitDistVec(distVec);
840 DimMatrix permeabilityJ
841 = problem().spatialParams().intrinsicPermeability(neighbor);
844 DimMatrix meanPermeability(0);
851 PhaseVector rhoMean(0.);
852 for (
int phaseIdx=0; phaseIdx<NumPhases; phaseIdx++)
853 rhoMean[phaseIdx] =0.5 * (cellDataI.density(phaseIdx)+cellDataJ.density(phaseIdx));
856 Dune::FieldVector<Scalar,2> potential(0.);
859 if (!cellDataJ.hasVolumeDerivatives())
860 asImp_().volumeDerivatives(globalPosNeighbor, neighbor);
862 ComponentVector dv_dC(0.), graddv_dC(0.);
863 for (
int compIdx = 0; compIdx < NumComponents; ++compIdx)
865 dv_dC[compIdx]= (cellDataJ.dv(compIdx)
866 + cellDataI.dv(compIdx)) * 0.5;
867 graddv_dC[compIdx] = (cellDataJ.dv(compIdx)
868 - cellDataI.dv(compIdx)) / dist;
882 GlobalPosition globalPos3(0.);
884 GlobalPosition globalPos4(0.);
886 TransmissivityMatrix T(0.);
889 GlobalPosition globalPosAdditional3(0.);
890 int eIdxGlobalAdditional3=-1;
891 GlobalPosition globalPosAdditional4(0.);
892 int eIdxGlobalAdditional4=-1;
894 TransmissivityMatrix additionalT(0.);
896 int interactionRegions
897 = problem().variables().getMpfaData3D(intersection, T,
898 globalPos3, eIdxGlobal3, globalPos4, eIdxGlobal4 );
899 if (interactionRegions == 0)
900 interactionRegions = problem().pressureModel().computeTransmissibilities(isIt,T,
901 globalPos3, eIdxGlobal3, globalPos4, eIdxGlobal4 );
902 if(!interactionRegions)
903 Dune::dgrave <<
"something went wrong getting mpfa data on cell " << eIdxGlobalI << std::endl;
906 CellData& cellData3 = problem().variables().cellData(eIdxGlobal3);
907 CellData& cellData4 = problem().variables().cellData(eIdxGlobal4);
908 Scalar temp1 = globalPos * this->gravity_;
909 Scalar temp2 = globalPosNeighbor * this->gravity_;
910 Scalar temp3 = globalPos3 * this->gravity_;
911 Scalar temp4 = globalPos4 * this->gravity_;
913 for (
int phaseIdx = 0; phaseIdx < NumPhases; ++phaseIdx)
915 potential[phaseIdx] = (cellDataI.pressure(phaseIdx)-temp1*rhoMean[phaseIdx]) * T[0]
916 + (cellDataJ.pressure(phaseIdx)-temp2*rhoMean[phaseIdx]) * T[1]
917 + (cellData3.pressure(phaseIdx)-temp3*rhoMean[phaseIdx]) * T[2]
918 + (cellData4.pressure(phaseIdx)-temp4*rhoMean[phaseIdx]) * T[3];
922 if(interactionRegions != 1)
924 for(
int banana = 1; banana < interactionRegions; banana ++)
927 problem().variables().getMpfaData3D(intersection, additionalT,
928 globalPosAdditional3, eIdxGlobalAdditional3,
929 globalPosAdditional4, eIdxGlobalAdditional4 ,
932 Scalar gravityContributionAdditonal
933 = temp1 * additionalT[0] + temp2 * additionalT[1]
934 + globalPosAdditional3*this->gravity_ * additionalT[2]
935 + globalPosAdditional4*this->gravity_ * additionalT[3];
936 CellData& cellDataA3 = problem().variables().cellData(eIdxGlobalAdditional3);
937 CellData& cellDataA4 = problem().variables().cellData(eIdxGlobalAdditional4);
939 for (
int phaseIdx = 0; phaseIdx < NumPhases; ++phaseIdx)
941 potential[phaseIdx] += cellDataI.pressure(phaseIdx) * additionalT[0]
942 + cellDataJ.pressure(phaseIdx) * additionalT[1]
943 + cellDataA3.pressure(phaseIdx) * additionalT[2]
944 + cellDataA4.pressure(phaseIdx) * additionalT[3];
945 potential[phaseIdx] -= gravityContributionAdditonal * rhoMean[phaseIdx];
951 PhaseVector lambda(0.), dV(0.), gV(0.);
954 std::vector<const CellData*> upwindCellData(NumPhases);
956 for (
int phaseIdx = 0; phaseIdx < NumPhases; ++phaseIdx)
958 int eqIdx = phaseIdx + 1;
959 if (potential[phaseIdx] > 0.)
960 upwindCellData[phaseIdx] = &cellDataI;
961 else if (potential[phaseIdx] < 0.)
962 upwindCellData[phaseIdx] = &cellDataJ;
965 if(cellDataI.isUpwindCell(intersection.indexInInside(), eqIdx))
966 upwindCellData[phaseIdx] = &cellDataI;
967 else if(cellDataJ.isUpwindCell(intersection.indexInOutside(), eqIdx))
968 upwindCellData[phaseIdx] = &cellDataJ;
974 if(!upwindCellData[phaseIdx])
976 Scalar averagedMassFraction[NumComponents];
977 for (
int compIdx = 0; compIdx < NumComponents; ++compIdx)
978 averagedMassFraction[compIdx]
979 =
harmonicMean(cellDataI.massFraction(phaseIdx, compIdx), cellDataJ.massFraction(phaseIdx, compIdx));
980 Scalar averageDensity =
harmonicMean(cellDataI.density(phaseIdx), cellDataJ.density(phaseIdx));
982 for (
int compIdx = 0; compIdx < NumComponents; ++compIdx)
984 dV[phaseIdx] += dv_dC[compIdx] * averagedMassFraction[compIdx];
985 gV[phaseIdx] += graddv_dC[compIdx] * averagedMassFraction[compIdx];
987 dV[phaseIdx] *= averageDensity;
988 gV[phaseIdx] *= averageDensity;
989 lambda[phaseIdx] =
harmonicMean(cellDataI.mobility(phaseIdx), cellDataJ.mobility(phaseIdx));
993 for (
int compIdx = 0; compIdx < NumComponents; ++compIdx)
995 dV[phaseIdx] += dv_dC[compIdx] * upwindCellData[phaseIdx]->massFraction(phaseIdx, compIdx);
996 gV[phaseIdx] += graddv_dC[compIdx] * upwindCellData[phaseIdx]->massFraction(phaseIdx, compIdx);
998 lambda[phaseIdx] = upwindCellData[phaseIdx]->mobility(phaseIdx);
999 dV[phaseIdx] *= upwindCellData[phaseIdx]->density(phaseIdx);
1000 gV[phaseIdx] *= upwindCellData[phaseIdx]->density(phaseIdx);
1006 this->A_[eIdxGlobalI][eIdxGlobalI] += (lambda[wPhaseIdx] * dV[wPhaseIdx] + lambda[nPhaseIdx] * dV[nPhaseIdx]) * T[0];
1007 this->A_[eIdxGlobalI][eIdxGlobalJ] += (lambda[wPhaseIdx] * dV[wPhaseIdx] + lambda[nPhaseIdx] * dV[nPhaseIdx]) * T[1];
1008 this->A_[eIdxGlobalI][eIdxGlobal3] += (lambda[wPhaseIdx] * dV[wPhaseIdx] + lambda[nPhaseIdx] * dV[nPhaseIdx]) * T[2];
1009 this->A_[eIdxGlobalI][eIdxGlobal4] += (lambda[wPhaseIdx] * dV[wPhaseIdx] + lambda[nPhaseIdx] * dV[nPhaseIdx]) * T[3];
1013 Scalar gravityContribution = temp1 * T[0] + temp2 * T[1] + temp3 * T[2] + temp4 * T[3];
1014 this->f_[eIdxGlobalI] += (rhoMean[wPhaseIdx] * lambda[wPhaseIdx] * dV[wPhaseIdx]
1015 + rhoMean[nPhaseIdx] * lambda[nPhaseIdx] * dV[nPhaseIdx]) * gravityContribution;
1018 Scalar weightingFactor = volume / perimeter;
1019 if(enableVolumeIntegral_)
1022 this->A_[eIdxGlobalI][eIdxGlobalI] -=
1023 weightingFactor * (lambda[wPhaseIdx] * gV[wPhaseIdx] + lambda[nPhaseIdx] * gV[nPhaseIdx]) * T[0];
1024 this->A_[eIdxGlobalI][eIdxGlobalJ] -=
1025 weightingFactor * (lambda[wPhaseIdx] * gV[wPhaseIdx] + lambda[nPhaseIdx] * gV[nPhaseIdx]) * T[1];
1026 this->A_[eIdxGlobalI][eIdxGlobal3] -=
1027 weightingFactor * (lambda[wPhaseIdx] * gV[wPhaseIdx] + lambda[nPhaseIdx] * gV[nPhaseIdx]) * T[2];
1028 this->A_[eIdxGlobalI][eIdxGlobal4] -=
1029 weightingFactor * (lambda[wPhaseIdx] * gV[wPhaseIdx] + lambda[nPhaseIdx] * gV[nPhaseIdx]) * T[3];
1032 this->f_[eIdxGlobalI] -= weightingFactor * gravityContribution *
1033 (rhoMean[wPhaseIdx] * lambda[wPhaseIdx] * gV[wPhaseIdx] + rhoMean[nPhaseIdx] * lambda[nPhaseIdx] * gV[nPhaseIdx]);
1037 Scalar pcGradient = cellDataI.capillaryPressure() * T[0]
1038 + cellDataJ.capillaryPressure() * T[1]
1039 + cellData3.capillaryPressure() * T[2]
1040 + cellData4.capillaryPressure() * T[3];
1042 if (this->pressureType == pw)
1043 pcGradient *= + lambda[nPhaseIdx] * dV[nPhaseIdx]
1044 - enableVolumeIntegral_ * weightingFactor * lambda[nPhaseIdx] * gV[nPhaseIdx];
1045 else if (this->pressureType == pn)
1046 pcGradient *= - lambda[wPhaseIdx] * dV[wPhaseIdx]
1047 + enableVolumeIntegral_ * weightingFactor * lambda[wPhaseIdx] * gV[wPhaseIdx];
1049 this->f_[eIdxGlobalI] += pcGradient;
1052 if(interactionRegions != 1)
1054 for(
int banana = 1; banana < interactionRegions; banana ++)
1057 problem().variables().getMpfaData3D(intersection, additionalT,
1058 globalPosAdditional3, eIdxGlobalAdditional3,
1059 globalPosAdditional4, eIdxGlobalAdditional4 ,
1062 Scalar gravityContributionAdditonal
1063 = temp1 * additionalT[0] + temp2 * additionalT[1]
1064 + globalPosAdditional3*this->gravity_ * additionalT[2]
1065 + globalPosAdditional4*this->gravity_ * additionalT[3];
1066 CellData& cellDataA3 = problem().variables().cellData(eIdxGlobalAdditional3);
1067 CellData& cellDataA4 = problem().variables().cellData(eIdxGlobalAdditional4);
1071 this->A_[eIdxGlobalI][eIdxGlobalI] +=
1072 (lambda[wPhaseIdx] * dV[wPhaseIdx] + lambda[nPhaseIdx] * dV[nPhaseIdx]) * additionalT[0];
1073 this->A_[eIdxGlobalI][eIdxGlobalJ] +=
1074 (lambda[wPhaseIdx] * dV[wPhaseIdx] + lambda[nPhaseIdx] * dV[nPhaseIdx]) * additionalT[1];
1075 this->A_[eIdxGlobalI][eIdxGlobalAdditional3] +=
1076 (lambda[wPhaseIdx] * dV[wPhaseIdx] + lambda[nPhaseIdx] * dV[nPhaseIdx]) * additionalT[2];
1077 this->A_[eIdxGlobalI][eIdxGlobalAdditional4] +=
1078 (lambda[wPhaseIdx] * dV[wPhaseIdx] + lambda[nPhaseIdx] * dV[nPhaseIdx]) * additionalT[3];
1082 this->f_[eIdxGlobalI] += (rhoMean[wPhaseIdx] * lambda[wPhaseIdx] * dV[wPhaseIdx]
1083 + rhoMean[nPhaseIdx] * lambda[nPhaseIdx] * dV[nPhaseIdx]) * gravityContributionAdditonal;
1086 if(enableVolumeIntegral_)
1089 this->A_[eIdxGlobalI][eIdxGlobalI] -=
1090 weightingFactor * (lambda[wPhaseIdx] * gV[wPhaseIdx] + lambda[nPhaseIdx] * gV[nPhaseIdx]) * additionalT[0];
1091 this->A_[eIdxGlobalI][eIdxGlobalJ] -=
1092 weightingFactor * (lambda[wPhaseIdx] * gV[wPhaseIdx] + lambda[nPhaseIdx] * gV[nPhaseIdx]) * additionalT[1];
1093 this->A_[eIdxGlobalI][eIdxGlobalAdditional3] -=
1094 weightingFactor * (lambda[wPhaseIdx] * gV[wPhaseIdx] + lambda[nPhaseIdx] * gV[nPhaseIdx]) * additionalT[2];
1095 this->A_[eIdxGlobalI][eIdxGlobalAdditional4] -=
1096 weightingFactor * (lambda[wPhaseIdx] * gV[wPhaseIdx] + lambda[nPhaseIdx] * gV[nPhaseIdx]) * additionalT[3];
1099 this->f_[eIdxGlobalI] -= weightingFactor * gravityContribution *
1100 (rhoMean[wPhaseIdx] * lambda[wPhaseIdx] * gV[wPhaseIdx] + rhoMean[nPhaseIdx] * lambda[nPhaseIdx] * gV[nPhaseIdx]);
1104 Scalar addPcGradient = cellDataI.capillaryPressure() * additionalT[0]
1105 + cellDataJ.capillaryPressure() * additionalT[1]
1106 + cellDataA3.capillaryPressure() * additionalT[2]
1107 + cellDataA4.capillaryPressure() * additionalT[3];
1109 if (this->pressureType == pw)
1110 addPcGradient *= + lambda[nPhaseIdx] * dV[nPhaseIdx]
1111 - enableVolumeIntegral_ * weightingFactor * lambda[nPhaseIdx] * gV[nPhaseIdx];
1112 else if (this->pressureType == pn)
1113 addPcGradient *= - lambda[wPhaseIdx] * dV[wPhaseIdx]
1114 + enableVolumeIntegral_ * weightingFactor * lambda[wPhaseIdx] * gV[wPhaseIdx];
1116 this->f_[eIdxGlobalI] += addPcGradient;
1140template<
class TypeTag>
1142 const CellData& cellDataI)
1144 const auto& intersection = *isIt;
1147 auto elementI = intersection.inside();
1148 int eIdxGlobalI = problem().variables().index(elementI);
1151 const GlobalPosition& globalPos = elementI.geometry().center();
1154 auto neighbor = intersection.outside();
1155 int eIdxGlobalJ = problem().variables().index(neighbor);
1156 CellData& cellDataJ = problem().variables().cellData(eIdxGlobalJ);
1159 const GlobalPosition& globalPosNeighbor = neighbor.geometry().center();
1164 int phaseIdx = min(cellDataI.subdomain(), cellDataJ.subdomain());
1167 Scalar rhoMean = 0.5 * (cellDataI.density(phaseIdx) + cellDataJ.density(phaseIdx));
1170 Scalar potential = 0.;
1174 GlobalPosition globalPos3(0.);
1176 GlobalPosition globalPos4(0.);
1178 TransmissivityMatrix T(0.);
1181 GlobalPosition globalPosAdditional3(0.);
1182 int eIdxGlobalAdditional3=-1;
1183 GlobalPosition globalPosAdditional4(0.);
1184 int eIdxGlobalAdditional4=-1;
1186 TransmissivityMatrix additionalT(0.);
1188 int interactionRegions
1189 = problem().variables().getMpfaData3D(intersection, T,
1190 globalPos3, eIdxGlobal3, globalPos4, eIdxGlobal4 );
1191 if (interactionRegions == 0)
1192 interactionRegions = problem().pressureModel().computeTransmissibilities(isIt,T,
1193 globalPos3, eIdxGlobal3, globalPos4, eIdxGlobal4 );
1194 if(!interactionRegions)
1195 Dune::dgrave <<
"something went wrong getting mpfa data on cell " << eIdxGlobalI << std::endl;
1198 CellData& cellData3 = problem().variables().cellData(eIdxGlobal3);
1199 CellData& cellData4 = problem().variables().cellData(eIdxGlobal4);
1200 Scalar temp1 = globalPos * this->gravity_;
1201 Scalar temp2 = globalPosNeighbor * this->gravity_;
1202 Scalar temp3 = globalPos3 * this->gravity_;
1203 Scalar temp4 = globalPos4 * this->gravity_;
1205 potential = (cellDataI.pressure(phaseIdx)-temp1*rhoMean) * T[0]
1206 + (cellDataJ.pressure(phaseIdx)-temp2*rhoMean) * T[1]
1207 + (cellData3.pressure(phaseIdx)-temp3*rhoMean) * T[2]
1208 + (cellData4.pressure(phaseIdx)-temp4*rhoMean) * T[3];
1211 if(interactionRegions != 1)
1213 for(
int banana = 1; banana < interactionRegions; banana ++)
1216 problem().variables().getMpfaData3D(intersection, additionalT,
1217 globalPosAdditional3, eIdxGlobalAdditional3,
1218 globalPosAdditional4, eIdxGlobalAdditional4 ,
1221 Scalar gravityContributionAdditonal
1222 = temp1 * additionalT[0] + temp2 * additionalT[1]
1223 + globalPosAdditional3*this->gravity_ * additionalT[2]
1224 + globalPosAdditional4*this->gravity_ * additionalT[3];
1225 CellData& cellDataA3 = problem().variables().cellData(eIdxGlobalAdditional3);
1226 CellData& cellDataA4 = problem().variables().cellData(eIdxGlobalAdditional4);
1228 potential += cellDataI.pressure(phaseIdx) * additionalT[0]
1229 + cellDataJ.pressure(phaseIdx) * additionalT[1]
1230 + cellDataA3.pressure(phaseIdx) * additionalT[2]
1231 + cellDataA4.pressure(phaseIdx) * additionalT[3];
1232 potential -= gravityContributionAdditonal * rhoMean;
1240 if (potential >= 0.)
1241 lambda = cellDataI.mobility(phaseIdx);
1243 lambda = cellDataJ.mobility(phaseIdx);
1247 this->A_[eIdxGlobalI][eIdxGlobalI] += lambda * T[0];
1248 this->A_[eIdxGlobalI][eIdxGlobalJ] += lambda * T[1];
1249 this->A_[eIdxGlobalI][eIdxGlobal3] += lambda * T[2];
1250 this->A_[eIdxGlobalI][eIdxGlobal4] += lambda * T[3];
1253 Scalar gravityContribution = temp1 * T[0] + temp2 * T[1] + temp3 * T[2] + temp4 * T[3];
1254 this->f_[eIdxGlobalI] += lambda * rhoMean * gravityContribution;
1257 if(interactionRegions != 1)
1259 for(
int banana = 1; banana < interactionRegions; banana ++)
1262 problem().variables().getMpfaData3D(intersection, additionalT,
1263 globalPosAdditional3, eIdxGlobalAdditional3,
1264 globalPosAdditional4, eIdxGlobalAdditional4 ,
1267 Scalar gravityContributionAdditonal
1268 = temp1 * additionalT[0] + temp2 * additionalT[1]
1269 + globalPosAdditional3*this->gravity_ * additionalT[2]
1270 + globalPosAdditional4*this->gravity_ * additionalT[3];
1274 this->A_[eIdxGlobalI][eIdxGlobalI] += lambda * additionalT[0];
1275 this->A_[eIdxGlobalI][eIdxGlobalJ] += lambda * additionalT[1];
1276 this->A_[eIdxGlobalI][eIdxGlobalAdditional3] += lambda * additionalT[2];
1277 this->A_[eIdxGlobalI][eIdxGlobalAdditional4] += lambda * additionalT[3];
1280 this->f_[eIdxGlobalI] += lambda * rhoMean* gravityContributionAdditonal;
1293template<
class TypeTag>
1297#ifdef noMultiphysics
1302 if(!fromPostTimestep)
1304 Scalar maxError = 0.;
1306 for (
const auto& element : elements(problem().gridView()))
1308 int eIdxGlobal = problem().variables().index(element);
1309 CellData& cellData = problem().variables().cellData(eIdxGlobal);
1311 if(cellData.fluidStateType() == 0)
1312 problem().pressureModel().updateMaterialLawsInElement(element, fromPostTimestep);
1314 problem().pressureModel().update1pMaterialLawsInElement(element, cellData, fromPostTimestep);
1316 maxError = max(maxError, fabs(cellData.volumeError()));
1318 this->maxError_ = maxError/problem().timeManager().timeStepSize();
1347template <
class TypeTag>
1349 TransmissivityMatrix& T,
1350 GlobalPosition& globalPos4,
1352 GlobalPosition& globalPos6,
1355 const auto& intersection = *isIt;
1358 auto element = intersection.inside();
1359 auto neighbor = intersection.outside();
1360 GlobalPosition globalPos1 = element.geometry().center();
1361 GlobalPosition globalPos2 = neighbor.geometry().center();
1362 DimMatrix K1(problem().spatialParams().intrinsicPermeability(element));
1363 DimMatrix K2(problem().spatialParams().intrinsicPermeability(neighbor));
1366 int intersectionID = 0;
1367 if(intersection.inside().level() < intersection.outside().level())
1368 intersectionID = problem().grid().localIdSet().subId(element,
1369 intersection.indexInInside(), 1);
1371 DUNE_THROW(Dune::NotImplemented,
" ABORT, transmiss calculated from wrong side!!");
1373 std::vector<int> localIrregularCells = irregularCellMap_[intersectionID];
1377 GlobalPosition globalPosFace12 = intersection.geometry().center();
1378 GlobalPosition outerNormaln12 = intersection.centerUnitOuterNormal();
1384 const auto isEndIt = problem().gridView().iend(neighbor);
1385 for (
auto isIt2 = problem().gridView().ibegin(neighbor); isIt2 != isEndIt; ++isIt2)
1387 const auto& intersection2 = *isIt2;
1390 if(!(intersection2.neighbor()) || intersection2.outside() == element)
1393 int currentNeighbor = problem().variables().index(intersection2.outside());
1396 if (find(localIrregularCells.begin(), localIrregularCells.end(),
1397 currentNeighbor) != localIrregularCells.end() && face24==isIt)
1399 else if (find(localIrregularCells.begin(), localIrregularCells.end(),
1400 currentNeighbor) != localIrregularCells.end() && face26==isIt)
1404 GlobalPosition vectorProduct =
crossProduct(face24->centerUnitOuterNormal(), intersection2.centerUnitOuterNormal());
1407 if((vectorProduct * outerNormaln12) > 0.)
1419 globalPos4 = face24->outside().geometry().center();
1420 eIdxGlobal4 = problem().variables().index(face24->outside());
1421 GlobalPosition outerNormaln24 = face24->centerUnitOuterNormal();
1423 DimMatrix K4(problem().spatialParams().intrinsicPermeability(face24->outside()));
1426 globalPos6 = face26->outside().geometry().center();
1427 eIdxGlobal6 = problem().variables().index(face26->outside());
1428 GlobalPosition outerNormaln26 = face26->centerUnitOuterNormal();
1430 DimMatrix K6(problem().spatialParams().intrinsicPermeability(face26->outside()));
1433 int localFace12 = intersection.indexInOutside();
1434 int localFace24 = face24->indexInInside();
1435 int localFace26 = face26->indexInInside();
1437 const auto refElement = referenceElement(neighbor);
1442 for(
int nectarine=0; nectarine < refElement.size(localFace12, 1, dim-1); nectarine++)
1445 int localEdgeOn12 = refElement.subEntity(localFace12, 1, nectarine, dim-1);
1447 for(
int plum = 0; plum < refElement.size(localFace26, 1,dim-1); plum++)
1450 if(refElement.subEntity(localFace12, 1, nectarine, dim-1)
1451 == refElement.subEntity(localFace26, 1, plum, dim-1))
1453 edge1226 = localEdgeOn12;
1458 GlobalPosition edgeCoord1226 =
1459 neighbor.geometry().global(refElement.position(edge1226, dim-1));
1464 for(
int nectarine=0; nectarine < refElement.size(localFace12, 1, dim-1); nectarine++)
1467 int localEdgeOn12 = refElement.subEntity(localFace12, 1, nectarine, dim-1);
1469 for(
int plum = 0; plum < refElement.size(localFace24, 1, dim-1); plum++)
1471 int localEdge24 = refElement.subEntity(localFace24, 1, plum, dim-1);
1472 if(localEdgeOn12 == localEdge24)
1474 edge1224 = localEdgeOn12;
1479 GlobalPosition edgeCoord1224 =
1480 neighbor.geometry().global(refElement.position(edge1224, dim-1));
1485 for(
int nectarine=0; nectarine < refElement.size(localFace24, 1, dim-1); nectarine++)
1488 int localEdgeOn24 = refElement.subEntity(localFace24, 1, nectarine, dim-1);
1490 for(
int plum = 0; plum < refElement.size(localFace26, 1, dim-1); plum++)
1492 int localEdge26 = refElement.subEntity(localFace26, 1, plum, dim-1);
1493 if(localEdgeOn24 == localEdge26)
1495 edge2426 = localEdgeOn24;
1500 GlobalPosition edgeCoord2426 =
1501 neighbor.geometry().global(refElement.position(edge2426, dim-1));
1505 GlobalPosition globalPosFace24 = face24->geometry().center();
1506 GlobalPosition globalPosFace26 = face26->geometry().center();
1510 Scalar subFaceArea12 =
crossProduct(edgeCoord1226-globalPosFace12, edgeCoord1224-globalPosFace12).two_norm();
1513 Scalar subFaceArea24 =
crossProduct(edgeCoord1224-globalPosFace24, edgeCoord2426-globalPosFace24).two_norm();
1516 Scalar subFaceArea26 =
crossProduct(edgeCoord1226-globalPosFace26, edgeCoord2426-globalPosFace26).two_norm();
1519 GlobalPosition nu11C2 =
crossProduct(edgeCoord1226-globalPos1, globalPosFace12-globalPos1);
1520 GlobalPosition nu12C2 =
crossProduct(edgeCoord1224-globalPos1, edgeCoord1226-globalPos1);
1521 GlobalPosition nu13C2 =
crossProduct(globalPosFace12-globalPos1, edgeCoord1224-globalPos1);
1523 GlobalPosition nu21C2 =
crossProduct(globalPosFace26-globalPos2, globalPosFace24-globalPos2);
1524 GlobalPosition nu22C2 =
crossProduct(globalPosFace12-globalPos2, globalPosFace26-globalPos2);
1525 GlobalPosition nu23C2 =
crossProduct(globalPosFace24-globalPos2, globalPosFace12-globalPos2);
1527 GlobalPosition nu41C2 =
crossProduct(edgeCoord2426-globalPos4, edgeCoord1224-globalPos4);
1528 GlobalPosition nu42C2 =
crossProduct(globalPosFace24-globalPos4, edgeCoord2426-globalPos4);
1529 GlobalPosition nu43C2 =
crossProduct(edgeCoord1224-globalPos4, globalPosFace24-globalPos4);
1531 GlobalPosition nu61C2 =
crossProduct(globalPosFace26-globalPos6, edgeCoord1226-globalPos6);
1532 GlobalPosition nu62C2 =
crossProduct(edgeCoord2426-globalPos6, globalPosFace26-globalPos6);
1533 GlobalPosition nu63C2 =
crossProduct(edgeCoord1226-globalPos6, edgeCoord2426-globalPos6);
1536 Scalar T1C2 = (globalPosFace12-globalPos1) *
crossProduct(edgeCoord1224-globalPos1, edgeCoord1226-globalPos1);
1537 Scalar T2C2 = (globalPosFace12-globalPos2) *
crossProduct(globalPosFace26-globalPos2, globalPosFace24-globalPos2);
1538 Scalar T4C2 = (globalPosFace24-globalPos4) *
crossProduct(edgeCoord2426-globalPos4, edgeCoord1224-globalPos4);
1539 Scalar T6C2 = (globalPosFace26-globalPos6) *
crossProduct(edgeCoord1226-globalPos6, edgeCoord2426-globalPos6);
1542 GlobalPosition K1nu11C2(0);
1543 K1.mv(nu11C2, K1nu11C2);
1544 GlobalPosition K1nu12C2(0);
1545 K1.mv(nu12C2, K1nu12C2);
1546 GlobalPosition K1nu13C2(0);
1547 K1.mv(nu13C2, K1nu13C2);
1549 GlobalPosition K2nu21C2(0);
1550 K2.mv(nu21C2, K2nu21C2);
1551 GlobalPosition K2nu22C2(0);
1552 K2.mv(nu22C2, K2nu22C2);
1553 GlobalPosition K2nu23C2(0);
1554 K2.mv(nu23C2, K2nu23C2);
1556 GlobalPosition K4nu41C2(0);
1557 K4.mv(nu41C2, K4nu41C2);
1558 GlobalPosition K4nu42C2(0);
1559 K4.mv(nu42C2, K4nu42C2);
1560 GlobalPosition K4nu43C2(0);
1561 K4.mv(nu43C2, K4nu43C2);
1563 GlobalPosition K6nu61C2(0);
1564 K6.mv(nu61C2, K6nu61C2);
1565 GlobalPosition K6nu62C2(0);
1566 K6.mv(nu62C2, K6nu62C2);
1567 GlobalPosition K6nu63C2(0);
1568 K6.mv(nu63C2, K6nu63C2);
1572 Scalar omega111C2 = (outerNormaln12 * K1nu11C2) * subFaceArea12/T1C2;
1573 Scalar omega112C2 = (outerNormaln12 * K1nu12C2) * subFaceArea12/T1C2;
1574 Scalar omega113C2 = (outerNormaln12 * K1nu13C2) * subFaceArea12/T1C2;
1576 Scalar omega121C2 = (outerNormaln12 * K2nu21C2) * subFaceArea12/T2C2;
1577 Scalar omega122C2 = (outerNormaln12 * K2nu22C2) * subFaceArea12/T2C2;
1578 Scalar omega123C2 = (outerNormaln12 * K2nu23C2) * subFaceArea12/T2C2;
1580 Scalar omega221C2 = (outerNormaln24 * K2nu21C2) * subFaceArea24/T2C2;
1581 Scalar omega222C2 = (outerNormaln24 * K2nu22C2) * subFaceArea24/T2C2;
1582 Scalar omega223C2 = (outerNormaln24 * K2nu23C2) * subFaceArea24/T2C2;
1584 Scalar omega241C2 = (outerNormaln24 * K4nu41C2) * subFaceArea24/T4C2;
1585 Scalar omega242C2 = (outerNormaln24 * K4nu42C2) * subFaceArea24/T4C2;
1586 Scalar omega243C2 = (outerNormaln24 * K4nu43C2) * subFaceArea24/T4C2;
1588 Scalar omega321C2 = (outerNormaln26 * K2nu21C2) * subFaceArea26/T2C2;
1589 Scalar omega322C2 = (outerNormaln26 * K2nu22C2) * subFaceArea26/T2C2;
1590 Scalar omega323C2 = (outerNormaln26 * K2nu23C2) * subFaceArea26/T2C2;
1592 Scalar omega361C2 = (outerNormaln26 * K6nu61C2) * subFaceArea26/T6C2;
1593 Scalar omega362C2 = (outerNormaln26 * K6nu62C2) * subFaceArea26/T6C2;
1594 Scalar omega363C2 = (outerNormaln26 * K6nu63C2) * subFaceArea26/T6C2;
1596 Scalar r211C2 = (nu21C2 * (edgeCoord1224-globalPos2))/T2C2;
1597 Scalar r212C2 = (nu21C2 * (edgeCoord1226-globalPos2))/T2C2;
1598 Scalar r213C2 = (nu21C2 * (edgeCoord2426-globalPos2))/T2C2;
1600 Scalar r221C2 = (nu22C2 * (edgeCoord1224-globalPos2))/T2C2;
1601 Scalar r222C2 = (nu22C2 * (edgeCoord1226-globalPos2))/T2C2;
1602 Scalar r223C2 = (nu22C2 * (edgeCoord2426-globalPos2))/T2C2;
1604 Scalar r231C2 = (nu23C2 * (edgeCoord1224-globalPos2))/T2C2;
1605 Scalar r232C2 = (nu23C2 * (edgeCoord1226-globalPos2))/T2C2;
1606 Scalar r233C2 = (nu23C2 * (edgeCoord2426-globalPos2))/T2C2;
1611 DimMatrix C(0), A(0);
1612 Dune::FieldMatrix<Scalar,dim,dim+1> D(0), B(0);
1615 C[0][0] = -omega121C2;
1616 C[0][1] = -omega122C2;
1617 C[0][2] = -omega123C2;
1618 C[1][0] = -omega221C2;
1619 C[1][1] = -omega222C2;
1620 C[1][2] = -omega223C2;
1621 C[2][0] = -omega321C2;
1622 C[2][1] = -omega322C2;
1623 C[2][2] = -omega323C2;
1625 D[0][1] = omega121C2 + omega122C2 + omega123C2;
1626 D[1][1] = omega221C2 + omega222C2 + omega223C2;
1627 D[2][1] = omega321C2 + omega322C2 + omega323C2;
1629 A[0][0] = omega121C2 - omega112C2 - omega111C2*r211C2 - omega113C2*r212C2;
1630 A[0][1] = omega122C2 - omega111C2*r221C2 - omega113C2*r222C2;
1631 A[0][2] = omega123C2 - omega111C2*r231C2 - omega113C2*r232C2;
1632 A[1][0] = omega221C2 - omega242C2*r211C2 - omega243C2*r213C2;
1633 A[1][1] = omega222C2 - omega241C2 - omega242C2*r221C2 - omega243C2*r223C2;
1634 A[1][2] = omega223C2 - omega242C2*r231C2 - omega243C2*r233C2;
1635 A[2][0] = omega321C2 - omega361C2*r213C2 - omega362C2*r212C2;
1636 A[2][1] = omega322C2 - omega361C2*r223C2 - omega362C2*r222C2;
1637 A[2][2] = omega323C2 - omega363C2 - omega361C2*r233C2 - omega362C2*r232C2;
1639 B[0][0] = -omega111C2 - omega112C2 - omega113C2;
1640 B[0][1] = omega121C2 + omega122C2 + omega123C2 + omega111C2*(1.0 - r211C2 - r221C2 -r231C2)
1641 + omega113C2*(1.0 - r212C2 - r222C2 - r232C2);
1642 B[1][1] = omega221C2 + omega222C2 + omega223C2 + omega242C2*(1.0 - r211C2 - r221C2 - r231C2)
1643 + omega243C2*(1.0 - r213C2 - r223C2 - r233C2);
1644 B[1][2] = -omega241C2 - omega242C2 - omega243C2;
1645 B[2][1] = omega321C2 + omega322C2 + omega323C2 + omega361C2*(1.0 - r213C2 - r223C2 - r233C2)
1646 + omega362C2*(1.0 - r212C2 -r222C2 -r232C2);
1647 B[2][3] = -omega361C2 - omega362C2 -omega363C2;
1656 D += B.leftmultiply(C.rightmultiply(A));
1660 if(maxInteractionVolumes ==1)
1662 T *= intersection.geometry().volume()/subFaceArea12 ;
1665 problem().variables().storeMpfaData3D(intersection, T, globalPos4, eIdxGlobal4, globalPos6, eIdxGlobal6);
1670 int countInteractionRegions = 1;
1672 TransmissivityMatrix additionalT(0.);
1674 auto outerCorner = intersection.inside().template subEntity<dim>(0);
1676 auto additional2 = intersection.inside();
1677 auto additional3 = intersection.inside();
1682 int localIdxLarge = searchCommonVertex_(intersection, outerCorner);
1686 if (problem().gridView().comm().size() > 1)
1687 if(outerCorner.partitionType() != Dune::InteriorEntity)
1692 int vIdxGlobal = problem().variables().index(outerCorner);
1693 InteractionVolume& interactionVolume
1694 = interactionVolumesContainer_->interactionVolume(vIdxGlobal);
1697 if(interactionVolume.isBoundaryInteractionVolume())
1701 int subVolumeFaceIdx = -1;
1702 bool properFluxDirection =
true;
1706 int hangingNodeType = interactionVolume.getHangingNodeType();
1708 if(hangingNodeType == InteractionVolume::noHangingNode)
1709 subVolumeFaceIdx = interactionVolumesContainer_->getMpfaCase8cells(isIt, localIdxLarge, interactionVolume, properFluxDirection);
1710 else if(hangingNodeType == InteractionVolume::sixSmallCells)
1711 subVolumeFaceIdx = interactionVolumesContainer_->getMpfaCase6cells(isIt, interactionVolume, properFluxDirection);
1713 Dune::dgrave <<
"HangingType " << hangingNodeType <<
" not implemented " << std::endl;
1715 caseL = this->transmissibilityAdapter_(isIt, interactionVolume, subVolumeFaceIdx,
1716 properFluxDirection, additional2, additional3, additionalT);
1722 problem().variables().storeMpfaData3D(intersection, additionalT,
1723 additional2.geometry().center(), problem().variables().index(additional2),
1724 additional3.geometry().center(), problem().variables().index(additional3),
1726 countInteractionRegions++;
1730 if(maxInteractionVolumes>2)
1733 std::vector<int> diagonal;
1734 for(
int verticeSmall = 0; verticeSmall < intersection.outside().subEntities(dim); ++verticeSmall)
1736 auto vSmall = intersection.outside().template subEntity<dim>(verticeSmall);
1740 if (problem().gridView().comm().size() > 1)
1741 if(vSmall.partitionType() != Dune::InteriorEntity)
1746 GlobalPosition vertexOnElement
1747 = refElement.position(verticeSmall, dim);
1749 for (
int indexOnFace = 0; indexOnFace < 4; indexOnFace++)
1752 GlobalPosition vertexOnInterface
1753 = intersection.geometryInOutside().corner(indexOnFace);
1755 if(vSmall != outerCorner
1756 && ((vertexOnInterface - vertexOnElement).two_norm()<1e-5))
1758 int vIdxGlobal2 = problem().variables().index(vSmall);
1760 InteractionVolume& interactionVolume2
1761 = interactionVolumesContainer_->interactionVolume(vIdxGlobal2);
1762 if(interactionVolume2.isBoundaryInteractionVolume())
1765 int hangingNodeType = interactionVolume2.getHangingNodeType();
1767 properFluxDirection =
true;
1769 if(hangingNodeType != InteractionVolume::fourSmallCellsFace)
1771 diagonal.push_back(problem().variables().index(vSmall));
1773 if(hangingNodeType == InteractionVolume::noHangingNode)
1776 Dune::dgrave <<
" args, noHangingNode on additional interaction region";
1779 else if(hangingNodeType == InteractionVolume::sixSmallCells)
1780 subVolumeFaceIdx = interactionVolumesContainer_->getMpfaCase6cells(isIt,
1781 interactionVolume2, properFluxDirection);
1783 subVolumeFaceIdx = interactionVolumesContainer_->getMpfaCase2or4cells(isIt,
1784 interactionVolume2, properFluxDirection);
1787 caseL = this->transmissibilityAdapter_(isIt, interactionVolume2, subVolumeFaceIdx,
1788 properFluxDirection, additional2, additional3, additionalT);
1794 problem().variables().storeMpfaData3D(intersection, additionalT,
1795 additional2.geometry().center(), problem().variables().index(additional2),
1796 additional3.geometry().center(), problem().variables().index(additional3),
1797 countInteractionRegions);
1798 countInteractionRegions++;
1808 problem().variables().storeMpfaData3D(intersection, T, globalPos4, eIdxGlobal4, globalPos6, eIdxGlobal6);
1811 Scalar weight = intersection.geometry().volume()/subFaceArea12;
1812 weight /=
static_cast<Scalar
>(countInteractionRegions);
1815 problem().variables().performTransmissitivityWeighting(intersection, weight);
1817 return countInteractionRegions;
1840template<
class TypeTag>
1842 InteractionVolume& interactionVolume,
1843 const int& subVolumeFaceIdx,
1844 bool properFluxDirection,
1845 Element& additional2,
1846 Element& additional3,
1847 TransmissivityMatrix& additionalT)
1850 if(interactionVolume.isBoundaryInteractionVolume())
1854 int hangingNodeType = interactionVolume.getHangingNodeType();
1857 Dune::FieldMatrix<Scalar, dim, 2 * dim - dim + 1> T(0);
1859 Dune::FieldVector<Scalar, dim> unity(1.);
1860 std::vector<Dune::FieldVector<Scalar, dim> > lambda
1861 = {unity, unity, unity, unity, unity, unity, unity, unity};
1862 Dune::FieldVector<bool, 4> useCases(
false);
1865 switch(subVolumeFaceIdx)
1868 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
1869 lambda, 0, 1, 2, 3, 4, 5);
1877 additional2 = interactionVolume.getSubVolumeElement(2);
1878 additional3 = interactionVolume.getSubVolumeElement(4);
1881 additional2 = interactionVolume.getSubVolumeElement(3);
1882 additional3 = interactionVolume.getSubVolumeElement(5);
1885 additional2 = interactionVolume.getSubVolumeElement(3);
1886 additional3 = interactionVolume.getSubVolumeElement(4);
1889 additional2 = interactionVolume.getSubVolumeElement(2);
1890 additional3 = interactionVolume.getSubVolumeElement(5);
1897 if (hangingNodeType == InteractionVolume::twoSmallCells
1898 || hangingNodeType == InteractionVolume::fourSmallCellsDiag )
1901 useCases[1] =
false;
1902 useCases[2] =
false;
1903 if(hangingNodeType != InteractionVolume::twoSmallCells)
1905 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
1906 lambda, 1, 3, 0, 2, 5, 7, useCases);
1909 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
1910 lambda, 1, 3, 0, 2, 5, 7);
1919 additional2 = interactionVolume.getSubVolumeElement(0);
1920 additional3 = interactionVolume.getSubVolumeElement(5);
1924 additional2 = interactionVolume.getSubVolumeElement(2);
1925 additional3 = interactionVolume.getSubVolumeElement(7);
1929 additional2 = interactionVolume.getSubVolumeElement(2);
1930 additional3 = interactionVolume.getSubVolumeElement(5);
1934 additional2 = interactionVolume.getSubVolumeElement(0);
1935 additional3 = interactionVolume.getSubVolumeElement(7);
1943 assert (hangingNodeType != InteractionVolume::twoSmallCells
1944 && hangingNodeType != InteractionVolume::fourSmallCellsDiag);
1946 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
1947 lambda, 3, 2, 1, 0, 7, 6);
1956 additional2 = interactionVolume.getSubVolumeElement(1);
1957 additional3 = interactionVolume.getSubVolumeElement(7);
1961 additional2 = interactionVolume.getSubVolumeElement(0);
1962 additional3 = interactionVolume.getSubVolumeElement(6);
1966 additional2 = interactionVolume.getSubVolumeElement(0);
1967 additional3 = interactionVolume.getSubVolumeElement(7);
1971 additional2 = interactionVolume.getSubVolumeElement(1);
1972 additional3 = interactionVolume.getSubVolumeElement(6);
1979 if (hangingNodeType == InteractionVolume::twoSmallCells
1980 || hangingNodeType == InteractionVolume::fourSmallCellsDiag)
1982 useCases[0] =
false;
1984 if(hangingNodeType != InteractionVolume::twoSmallCells)
1986 useCases[3] =
false;
1987 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T,
1988 interactionVolume, lambda, 2, 0, 3, 1, 6, 4, useCases);
1991 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
1992 lambda, 2, 0, 3, 1, 6, 4);
2001 additional2 = interactionVolume.getSubVolumeElement(3);
2002 additional3 = interactionVolume.getSubVolumeElement(6);
2006 additional2 = interactionVolume.getSubVolumeElement(1);
2007 additional3 = interactionVolume.getSubVolumeElement(4);
2011 additional2 = interactionVolume.getSubVolumeElement(1);
2012 additional3 = interactionVolume.getSubVolumeElement(6);
2016 additional2 = interactionVolume.getSubVolumeElement(3);
2017 additional3 = interactionVolume.getSubVolumeElement(4);
2024 if (hangingNodeType == InteractionVolume::fourSmallCellsEdge)
2026 assert(problem().variables().index(interactionVolume.getSubVolumeElement(4))
2027 != problem().variables().index(interactionVolume.getSubVolumeElement(5)));
2028 Dune::dgrave <<
" SubVolumeFace4 in hanging node type 3 should be modelled by"
2029 <<
" Tpfa!!" <<std::endl;
2033 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
2034 lambda, 5, 4, 7, 6, 1, 0);
2042 additional2 = interactionVolume.getSubVolumeElement(7);
2043 additional3 = interactionVolume.getSubVolumeElement(1);
2046 additional2 = interactionVolume.getSubVolumeElement(6);
2047 additional3 = interactionVolume.getSubVolumeElement(0);
2050 additional2 = interactionVolume.getSubVolumeElement(6);
2051 additional3 = interactionVolume.getSubVolumeElement(1);
2054 additional2 = interactionVolume.getSubVolumeElement(7);
2055 additional3 = interactionVolume.getSubVolumeElement(0);
2062 if (hangingNodeType == InteractionVolume::fourSmallCellsEdge)
2064 assert (problem().variables().index(interactionVolume.getSubVolumeElement(4))
2065 != problem().variables().index(interactionVolume.getSubVolumeElement(6)));
2066 Dune::dgrave <<
" SubVolumeFace5 in hanging node type 3 should be modelled by"
2067 <<
" Tpfa!!" <<std::endl;
2070 else if (hangingNodeType == InteractionVolume::sixSmallCells)
2072 useCases[0] =
false;
2075 useCases[3] =
false;
2077 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
2078 lambda, 7, 5, 6, 4, 3, 1, useCases);
2080 else if (hangingNodeType == InteractionVolume::fourSmallCellsDiag)
2083 useCases[1] =
false;
2084 useCases[2] =
false;
2086 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
2087 lambda, 7, 5, 6, 4, 3, 1, useCases);
2090 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
2091 lambda, 7, 5, 6, 4, 3, 1);
2100 additional2 = interactionVolume.getSubVolumeElement(6);
2101 additional3 = interactionVolume.getSubVolumeElement(3);
2105 additional2 = interactionVolume.getSubVolumeElement(4);
2106 additional3 = interactionVolume.getSubVolumeElement(1);
2110 additional2 = interactionVolume.getSubVolumeElement(4);
2111 additional3 = interactionVolume.getSubVolumeElement(3);
2115 additional2 = interactionVolume.getSubVolumeElement(6);
2116 additional3 = interactionVolume.getSubVolumeElement(1);
2123 if (hangingNodeType == InteractionVolume::fourSmallCellsEdge)
2125 assert (problem().variables().index(interactionVolume.getSubVolumeElement(4))
2126 != problem().variables().index(interactionVolume.getSubVolumeElement(5)));
2127 Dune::dgrave <<
" SubVolumeFace6 in hanging node type 3 should be modelled by"
2128 <<
" Tpfa!!" <<std::endl;
2132 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
2133 lambda, 6, 7, 4, 5, 2, 3);
2142 additional2 = interactionVolume.getSubVolumeElement(4);
2143 additional3 = interactionVolume.getSubVolumeElement(2);
2147 additional2 = interactionVolume.getSubVolumeElement(5);
2148 additional3 = interactionVolume.getSubVolumeElement(3);
2152 additional2 = interactionVolume.getSubVolumeElement(5);
2153 additional3 = interactionVolume.getSubVolumeElement(2);
2157 additional2 = interactionVolume.getSubVolumeElement(4);
2158 additional3 = interactionVolume.getSubVolumeElement(3);
2166 if (hangingNodeType == InteractionVolume::fourSmallCellsEdge)
2168 assert (problem().variables().index(interactionVolume.getSubVolumeElement(4))
2169 != problem().variables().index(interactionVolume.getSubVolumeElement(6)));
2170 Dune::dgrave <<
" SubVolumeFace5 in hanging node type 3 should be modelled by"
2171 <<
" Tpfa!!" <<std::endl;
2174 else if (hangingNodeType == InteractionVolume::sixSmallCells)
2177 useCases[1] =
false;
2178 useCases[2] =
false;
2181 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
2182 lambda, 4, 6, 5, 7, 0, 2, useCases);
2184 else if (hangingNodeType == InteractionVolume::fourSmallCellsDiag)
2186 useCases[0] =
false;
2189 useCases[3] =
false;
2190 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
2191 lambda, 4, 6, 5, 7, 0, 2, useCases);
2194 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
2195 lambda, 4, 6, 5, 7, 0, 2);
2204 additional2 = interactionVolume.getSubVolumeElement(5);
2205 additional3 = interactionVolume.getSubVolumeElement(0);
2209 additional2 = interactionVolume.getSubVolumeElement(7);
2210 additional3 = interactionVolume.getSubVolumeElement(2);
2214 additional2 = interactionVolume.getSubVolumeElement(7);
2215 additional3 = interactionVolume.getSubVolumeElement(0);
2219 additional2 = interactionVolume.getSubVolumeElement(5);
2220 additional3 = interactionVolume.getSubVolumeElement(2);
2228 if(hangingNodeType == InteractionVolume::noHangingNode
2229 || hangingNodeType == InteractionVolume::sixSmallCells)
2231 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
2232 lambda, 4, 0, 6, 2, 5, 1);
2234 else if (hangingNodeType == InteractionVolume::twoSmallCells
2235 || hangingNodeType == InteractionVolume::fourSmallCellsFace)
2237 caseL = mpfal3DTransmissibilityCalculator_.transmissibilityCaseTwo(T, interactionVolume,
2238 lambda, 4, 0, 2, 1);
2240 else if (hangingNodeType == InteractionVolume::fourSmallCellsDiag
2241 || (hangingNodeType == InteractionVolume::fourSmallCellsEdge
2242 && (problem().variables().index(interactionVolume.getSubVolumeElement(4))
2243 != problem().variables().index(interactionVolume.getSubVolumeElement(6)))) )
2245 useCases[0] =
false;
2247 useCases[2] =
false;
2250 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
2251 lambda, 4, 0, 6, 2, 5, 1, useCases);
2253 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge)
2255 useCases[0] =
false;
2258 useCases[3] =
false;
2260 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
2261 lambda, 4, 0, 6, 2, 5, 1, useCases);
2271 additional2 = interactionVolume.getSubVolumeElement(6);
2272 additional3 = interactionVolume.getSubVolumeElement(5);
2276 additional2 = interactionVolume.getSubVolumeElement(2);
2277 additional3 = interactionVolume.getSubVolumeElement(1);
2281 additional2 = interactionVolume.getSubVolumeElement(2);
2282 additional3 = interactionVolume.getSubVolumeElement(5);
2286 additional2 = interactionVolume.getSubVolumeElement(6);
2287 additional3 = interactionVolume.getSubVolumeElement(1);
2295 if(hangingNodeType == InteractionVolume::noHangingNode
2296 || hangingNodeType == InteractionVolume::sixSmallCells)
2298 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
2299 lambda, 1, 5, 3, 7, 0, 4);
2301 else if (hangingNodeType == InteractionVolume::twoSmallCells || hangingNodeType == InteractionVolume::fourSmallCellsFace)
2303 caseL = mpfal3DTransmissibilityCalculator_.transmissibilityCaseOne(T, interactionVolume,
2304 lambda, 1, 5, 3, 0);
2306 else if (hangingNodeType == InteractionVolume::fourSmallCellsDiag
2307 || (hangingNodeType == InteractionVolume::fourSmallCellsEdge
2308 &&(problem().variables().index(interactionVolume.getSubVolumeElement(4))
2309 != problem().variables().index(interactionVolume.getSubVolumeElement(6))) ))
2312 useCases[1] =
false;
2314 useCases[3] =
false;
2316 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
2317 lambda, 1, 5, 3, 7, 0, 4, useCases);
2319 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge
2320 &&(problem().variables().index(interactionVolume.getSubVolumeElement(4))
2321 != problem().variables().index(interactionVolume.getSubVolumeElement(5))) )
2324 useCases[1] =
false;
2325 useCases[2] =
false;
2328 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
2329 lambda, 1, 5, 3, 7, 0, 4, useCases);
2332 Dune::dgrave <<
" Missing case for subVolFaceIdx 9!!" <<std::endl;
2341 additional2 = interactionVolume.getSubVolumeElement(3);
2342 additional3 = interactionVolume.getSubVolumeElement(0);
2346 additional2 = interactionVolume.getSubVolumeElement(7);
2347 additional3 = interactionVolume.getSubVolumeElement(4);
2351 additional2 = interactionVolume.getSubVolumeElement(7);
2352 additional3 = interactionVolume.getSubVolumeElement(0);
2356 additional2 = interactionVolume.getSubVolumeElement(3);
2357 additional3 = interactionVolume.getSubVolumeElement(4);
2365 if (hangingNodeType == InteractionVolume::fourSmallCellsFace)
2366 caseL = mpfal3DTransmissibilityCalculator_.transmissibilityCaseTwo(T, interactionVolume,
2367 lambda, 7, 3, 1, 2);
2368 else if ((hangingNodeType == InteractionVolume::fourSmallCellsEdge
2369 &&(problem().variables().index(interactionVolume.getSubVolumeElement(4))
2370 != problem().variables().index(interactionVolume.getSubVolumeElement(6))) )
2371 || hangingNodeType == InteractionVolume::sixSmallCells)
2373 useCases[0] =
false;
2375 useCases[2] =
false;
2378 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
2379 lambda, 7, 3, 5, 1, 6, 2, useCases);
2381 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge)
2383 useCases[0] =
false;
2386 useCases[3] =
false;
2388 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
2389 lambda, 7, 3, 5, 1, 6, 2, useCases);
2391 else if (hangingNodeType == InteractionVolume::fourSmallCellsDiag)
2394 useCases[1] =
false;
2396 useCases[3] =
false;
2398 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
2399 lambda, 7, 3, 5, 1, 6, 2, useCases);
2402 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
2403 lambda, 7, 3, 5, 1, 6, 2);
2412 additional2 = interactionVolume.getSubVolumeElement(5);
2413 additional3 = interactionVolume.getSubVolumeElement(6);
2417 additional2 = interactionVolume.getSubVolumeElement(1);
2418 additional3 = interactionVolume.getSubVolumeElement(2);
2422 additional2 = interactionVolume.getSubVolumeElement(1);
2423 additional3 = interactionVolume.getSubVolumeElement(6);
2427 additional2 = interactionVolume.getSubVolumeElement(5);
2428 additional3 = interactionVolume.getSubVolumeElement(2);
2436 if(!interactionVolume.isHangingNodeVolume())
2437 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
2438 lambda, 2, 6, 0, 4, 3, 7);
2439 else if (hangingNodeType == InteractionVolume::fourSmallCellsFace)
2440 caseL = mpfal3DTransmissibilityCalculator_.transmissibilityCaseOne(T, interactionVolume,
2441 lambda, 2, 6, 0, 3);
2442 else if ((hangingNodeType == InteractionVolume::fourSmallCellsEdge
2443 && (problem().variables().index(interactionVolume.getSubVolumeElement(4))
2444 != problem().variables().index(interactionVolume.getSubVolumeElement(6))))
2445 || hangingNodeType == InteractionVolume::sixSmallCells)
2448 useCases[1] =
false;
2450 useCases[3] =
false;
2452 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
2453 lambda, 2, 6, 0, 4, 3, 7, useCases);
2455 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge
2456 && (problem().variables().index(interactionVolume.getSubVolumeElement(4))
2457 != problem().variables().index(interactionVolume.getSubVolumeElement(5))) )
2460 useCases[1] =
false;
2461 useCases[2] =
false;
2464 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
2465 lambda, 2, 6, 0, 4, 3, 7, useCases);
2467 else if (hangingNodeType == InteractionVolume::fourSmallCellsDiag)
2469 useCases[0] =
false;
2471 useCases[2] =
false;
2474 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
2475 lambda, 2, 6, 0, 4, 3, 7, useCases);
2485 additional2 = interactionVolume.getSubVolumeElement(0);
2486 additional3 = interactionVolume.getSubVolumeElement(3);
2490 additional2 = interactionVolume.getSubVolumeElement(4);
2491 additional3 = interactionVolume.getSubVolumeElement(7);
2495 additional2 = interactionVolume.getSubVolumeElement(4);
2496 additional3 = interactionVolume.getSubVolumeElement(3);
2500 additional2 = interactionVolume.getSubVolumeElement(0);
2501 additional3 = interactionVolume.getSubVolumeElement(7);
2508 if(!properFluxDirection)
2515 swap(additionalT[0], additionalT[1]);
Define some often used mathematical functions.
Simplifies writing multi-file VTK datasets.
Contains a class to exchange entries of a vector.
Provides methods for transmissibility calculation in 3-d.
Defines the properties required for the adaptive sequential 2p2c models.
Interaction volume container for compositional adaptive 3-d (using MPFA L-method).
Finite volume 2p2c pressure model with multi-physics.
constexpr Scalar harmonicMean(Scalar x, Scalar y, Scalar wx=1.0, Scalar wy=1.0) noexcept
Calculate the (weighted) harmonic mean of two scalar values.
Definition: math.hh:68
void harmonicMeanMatrix(Dune::FieldMatrix< Scalar, m, n > &K, const Dune::FieldMatrix< Scalar, m, n > &Ki, const Dune::FieldMatrix< Scalar, m, n > &Kj)
Calculate the harmonic mean of a fixed-size matrix.
Definition: math.hh:108
Dune::FieldVector< Scalar, 3 > crossProduct(const Dune::FieldVector< Scalar, 3 > &vec1, const Dune::FieldVector< Scalar, 3 > &vec2)
Cross product of two vectors in three-dimensional Euclidean space.
Definition: math.hh:640
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 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
A data handle class to exchange entries of a vector.
Definition: vectorcommdatahandle.hh:78
Provides methods for transmissibility calculation in 3-d.
Definition: 3dtransmissibilitycalculator.hh:51
The finite volume model for the solution of the compositional pressure equation.
Definition: fv3dpressureadaptive.hh:81
void update()
Definition: fv3dpressureadaptive.hh:162
void adaptPressure()
Adapt primary variables vector after adapting the grid.
Definition: fv3dpressureadaptive.hh:213
bool enableMPFA
Enables mpfa on hanging nodes (on by default)
Definition: fv3dpressureadaptive.hh:302
int transmissibilityAdapter_(const IntersectionIterator &isIt, InteractionVolume &interactionVolume, const int &subVolumeFaceIdx, bool properFluxDirection, Element &additional2, Element &additional3, TransmissivityMatrix &additionalT)
Adapter to use the general implementation of the mpfa-l for the compositional models.
Definition: fv3dpressureadaptive.hh:1841
void initializeMatrix()
initializes the matrix to store the system of equations
Definition: fv3dpressureadaptive.hh:313
void initializeMatrixRowSize()
Initialize the row sizes of the sparse global matrix.
Definition: fv3dpressureadaptive.hh:329
FvMpfaL3dTransmissibilityCalculator< TypeTag > mpfal3DTransmissibilityCalculator_
The common implementation to calculate the Transmissibility with the mpfa-L-method.
Definition: fv3dpressureadaptive.hh:308
bool enableVolumeIntegral_
Calculates the volume integral (on by default)
Definition: fv3dpressureadaptive.hh:301
void initializeMatrixIndices()
Determine position of matrix entries.
Definition: fv3dpressureadaptive.hh:578
int computeTransmissibilities(const IntersectionIterator &, TransmissivityMatrix &, GlobalPosition &, int &, GlobalPosition &, int &)
Computes the transmissibility coefficients for the MPFA-l method in 3D.
Definition: fv3dpressureadaptive.hh:1348
void getMpfaFlux(const IntersectionIterator &, const CellData &)
Compute flux through an irregular interface using a mpfa method.
Definition: fv3dpressureadaptive.hh:804
void assemble(bool first)
function which assembles the system of equations to be solved
Definition: fv3dpressureadaptive.hh:643
FV3dPressure2P2CAdaptive(Problem &problem)
Constructs a FVPressure2P2C object.
Definition: fv3dpressureadaptive.hh:248
std::map< int, std::vector< int > > irregularCellMap_
Container to store all cell's Indice with a hanging node.
Definition: fv3dpressureadaptive.hh:300
void initialize(bool solveTwice=false)
Definition: fv3dpressureadaptive.hh:182
void get1pMpfaFlux(const IntersectionIterator &, const CellData &)
Compute single-phase flux through an irregular interface using a mpfa method.
Definition: fv3dpressureadaptive.hh:1141
InteractionVolumeContainer * interactionVolumesContainer_
A pointer to the adaptive interaction volumes container.
Definition: fv3dpressureadaptive.hh:306
void updateMaterialLaws(bool fromPostTimestep=false)
Definition: fv3dpressureadaptive.hh:1294
int maxInteractionVolumes
Maximum number of interaction volumes considered (4 by default)
Definition: fv3dpressureadaptive.hh:303
Interaction volume container for compositional adaptive 3-d (using MPFA L-method) Model.
Definition: fvmpfal3dinteractionvolumecontaineradaptive.hh:46
bool enableVolumeIntegral
Enables the volume integral of the pressure equation.
Definition: fvpressure.hh:184
Problem & problem_
Definition: fvpressure.hh:183
void updateMaterialLaws(bool postTimeStep=false)
Updates secondary variables.
Definition: fvpressurecompositional.hh:706
void update()
Compositional pressure solution routine: update estimate for secants, assemble, solve.
Definition: fvpressurecompositional.hh:130
The finite volume model for the solution of the compositional pressure equation.
Definition: fvpressuremultiphysics.hh:72
void updateMaterialLaws(bool postTimeStep=false)
constitutive functions are updated once if new concentrations are calculated and stored in the variab...
Definition: fvpressuremultiphysics.hh:793
static constexpr int pressureType
gives kind of pressure used ( , , )
Definition: fvpressuremultiphysics.hh:234
typename SolutionTypes::ElementMapper ElementMapper
Definition: fvpressuremultiphysics.hh:226
Class including the information of a 3d interaction volume of an adaptive MPFA L-method that does not...
Definition: linteractionvolume3dadaptive.hh:192
Type of the data container for one interaction volume.
Definition: porousmediumflow/sequential/cellcentered/mpfa/properties.hh:102
Type of the data container for one interaction volume.
Definition: porousmediumflow/sequential/cellcentered/mpfa/properties.hh:104
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
@ rhs
index for the right hand side entry
Definition: sequential/cellcentered/pressure.hh:88
@ matrix
index for the global matrix entry
Definition: sequential/cellcentered/pressure.hh:89
Properties for a MPFA method.