24#ifndef DUMUX_FVMPFAL3DVELOCITY2P_ADAPTIVE_HH
25#define DUMUX_FVMPFAL3DVELOCITY2P_ADAPTIVE_HH
61 dim = GridView::dimension, dimWorld = GridView::dimensionworld
77 using PrimaryVariables =
typename SolutionTypes::PrimaryVariables;
81 using Element =
typename GridView::Traits::template Codim<0>::Entity;
82 using Grid =
typename GridView::Grid;
83 using IndexSet =
typename GridView::IndexSet;
95 pw = Indices::pressureW,
96 pn = Indices::pressureNw,
97 pGlobal = Indices::pressureGlobal,
98 sw = Indices::saturationW,
99 sn = Indices::saturationNw,
100 vw = Indices::velocityW,
101 vn = Indices::velocityNw,
102 vt = Indices::velocityTotal
106 wPhaseIdx = Indices::wPhaseIdx,
107 nPhaseIdx = Indices::nPhaseIdx,
108 pressureIdx = Indices::pressureIdx,
109 saturationIdx = Indices::saturationIdx,
110 pressureEqIdx = Indices::pressureEqIdx,
111 satEqIdx = Indices::satEqIdx,
112 numPhases = getPropValue<TypeTag, Properties::NumPhases>()
119 dirichletDirichlet = 1,
120 dirichletNeumann = 2,
125 innerEdgeFace = 2, innerSideFace = 1
128 using GravityVector = Dune::FieldVector<Scalar, dimWorld>;
129 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
130 using DimVector = Dune::FieldVector<Scalar, dim>;
132 InteractionVolumeContainer& interactionVolumes_()
134 return problem_.pressureModel().interactionVolumes();
139 return problem_.pressureModel().transmissibilityCalculator();
149 ParentType(problem), problem_(problem), gravity_(problem.gravity())
151 density_[wPhaseIdx] = 0.;
152 density_[nPhaseIdx] = 0.;
153 viscosity_[wPhaseIdx] = 0.;
154 viscosity_[nPhaseIdx] = 0.;
159 CellData & cellData1, CellData & cellData2, CellData & cellData3, CellData & cellData4,
160 CellData & cellData5, CellData & cellData6, CellData & cellData7, CellData & cellData8,
int fIdx = -1);
167 const auto element = *problem_.gridView().template begin<0>();
168 FluidState fluidState;
169 fluidState.setPressure(wPhaseIdx, problem_.referencePressure(element));
170 fluidState.setPressure(nPhaseIdx, problem_.referencePressure(element));
171 fluidState.setTemperature(problem_.temperature(element));
172 fluidState.setSaturation(wPhaseIdx, 1.);
173 fluidState.setSaturation(nPhaseIdx, 0.);
184 const GravityVector& gravity_;
186 Scalar density_[numPhases];
187 Scalar viscosity_[numPhases];
189 static constexpr Scalar threshold_ = 1e-15;
191 static const int velocityType_ = getPropValue<TypeTag, Properties::VelocityFormulation>();
193 static const int pressureType_ = getPropValue<TypeTag, Properties::PressureFormulation>();
195 static const int saturationType_ = getPropValue<TypeTag, Properties::SaturationFormulation>();
217template<
class TypeTag>
219 CellData & cellData1, CellData & cellData2, CellData & cellData3, CellData & cellData4,
220 CellData & cellData5, CellData & cellData6, CellData & cellData7, CellData & cellData8,
int fIdx)
222 auto element1 = interactionVolume.getSubVolumeElement(0);
223 auto element2 = interactionVolume.getSubVolumeElement(1);
224 auto element3 = interactionVolume.getSubVolumeElement(2);
225 auto element4 = interactionVolume.getSubVolumeElement(3);
226 auto element5 = interactionVolume.getSubVolumeElement(4);
227 auto element6 = interactionVolume.getSubVolumeElement(5);
228 auto element7 = interactionVolume.getSubVolumeElement(6);
229 auto element8 = interactionVolume.getSubVolumeElement(7);
232 int globalIdx1 = problem_.variables().index(element1);
233 int globalIdx2 = problem_.variables().index(element2);
234 int globalIdx3 = problem_.variables().index(element3);
235 int globalIdx4 = problem_.variables().index(element4);
236 int globalIdx5 = problem_.variables().index(element5);
237 int globalIdx6 = problem_.variables().index(element6);
238 int globalIdx7 = problem_.variables().index(element7);
239 int globalIdx8 = problem_.variables().index(element8);
242 Dune::FieldVector<Scalar, 8> potW(0);
243 Dune::FieldVector<Scalar, 8> potNw(0);
245 potW[0] = cellData1.potential(wPhaseIdx);
246 potW[1] = cellData2.potential(wPhaseIdx);
247 potW[2] = cellData3.potential(wPhaseIdx);
248 potW[3] = cellData4.potential(wPhaseIdx);
249 potW[4] = cellData5.potential(wPhaseIdx);
250 potW[5] = cellData6.potential(wPhaseIdx);
251 potW[6] = cellData7.potential(wPhaseIdx);
252 potW[7] = cellData8.potential(wPhaseIdx);
254 potNw[0] = cellData1.potential(nPhaseIdx);
255 potNw[1] = cellData2.potential(nPhaseIdx);
256 potNw[2] = cellData3.potential(nPhaseIdx);
257 potNw[3] = cellData4.potential(nPhaseIdx);
258 potNw[4] = cellData5.potential(nPhaseIdx);
259 potNw[5] = cellData6.potential(nPhaseIdx);
260 potNw[6] = cellData7.potential(nPhaseIdx);
261 potNw[7] = cellData8.potential(nPhaseIdx);
264 Dune::FieldVector<Scalar, numPhases> lambda1(cellData1.mobility(wPhaseIdx));
265 lambda1[nPhaseIdx] = cellData1.mobility(nPhaseIdx);
268 Scalar lambdaTotal1 = lambda1[wPhaseIdx] + lambda1[nPhaseIdx];
271 Dune::FieldVector<Scalar, numPhases> lambda2(cellData2.mobility(wPhaseIdx));
272 lambda2[nPhaseIdx] = cellData2.mobility(nPhaseIdx);
275 Scalar lambdaTotal2 = lambda2[wPhaseIdx] + lambda2[nPhaseIdx];
278 Dune::FieldVector<Scalar, numPhases> lambda3(cellData3.mobility(wPhaseIdx));
279 lambda3[nPhaseIdx] = cellData3.mobility(nPhaseIdx);
282 Scalar lambdaTotal3 = lambda3[wPhaseIdx] + lambda3[nPhaseIdx];
285 Dune::FieldVector<Scalar, numPhases> lambda4(cellData4.mobility(wPhaseIdx));
286 lambda4[nPhaseIdx] = cellData4.mobility(nPhaseIdx);
289 Scalar lambdaTotal4 = lambda4[wPhaseIdx] + lambda4[nPhaseIdx];
292 Dune::FieldVector<Scalar, numPhases> lambda5(cellData5.mobility(wPhaseIdx));
293 lambda5[nPhaseIdx] = cellData5.mobility(nPhaseIdx);
296 Scalar lambdaTotal5 = lambda5[wPhaseIdx] + lambda5[nPhaseIdx];
299 Dune::FieldVector<Scalar, numPhases> lambda6(cellData6.mobility(wPhaseIdx));
300 lambda6[nPhaseIdx] = cellData6.mobility(nPhaseIdx);
303 Scalar lambdaTotal6 = lambda6[wPhaseIdx] + lambda6[nPhaseIdx];
306 Dune::FieldVector<Scalar, numPhases> lambda7(cellData7.mobility(wPhaseIdx));
307 lambda7[nPhaseIdx] = cellData7.mobility(nPhaseIdx);
310 Scalar lambdaTotal7 = lambda7[wPhaseIdx] + lambda7[nPhaseIdx];
313 Dune::FieldVector<Scalar, numPhases> lambda8(cellData8.mobility(wPhaseIdx));
314 lambda8[nPhaseIdx] = cellData8.mobility(nPhaseIdx);
317 Scalar lambdaTotal8 = lambda8[wPhaseIdx] + lambda8[nPhaseIdx];
319 std::vector<DimVector> lambda(8);
320 lambda[0][0] = lambdaTotal1;
321 lambda[0][1] = lambdaTotal1;
322 lambda[0][2] = lambdaTotal1;
323 lambda[1][0] = lambdaTotal2;
324 lambda[1][1] = lambdaTotal2;
325 lambda[1][2] = lambdaTotal2;
326 lambda[2][0] = lambdaTotal3;
327 lambda[2][1] = lambdaTotal3;
328 lambda[2][2] = lambdaTotal3;
329 lambda[3][0] = lambdaTotal4;
330 lambda[3][1] = lambdaTotal4;
331 lambda[3][2] = lambdaTotal4;
332 lambda[4][0] = lambdaTotal5;
333 lambda[4][1] = lambdaTotal5;
334 lambda[4][2] = lambdaTotal5;
335 lambda[5][0] = lambdaTotal6;
336 lambda[5][1] = lambdaTotal6;
337 lambda[5][2] = lambdaTotal6;
338 lambda[6][0] = lambdaTotal7;
339 lambda[6][1] = lambdaTotal7;
340 lambda[6][2] = lambdaTotal7;
341 lambda[7][0] = lambdaTotal8;
342 lambda[7][1] = lambdaTotal8;
343 lambda[7][2] = lambdaTotal8;
345 Scalar potentialDiffW0 = 0;
346 Scalar potentialDiffW1 = 0;
347 Scalar potentialDiffW2 = 0;
348 Scalar potentialDiffW3 = 0;
349 Scalar potentialDiffW4 = 0;
350 Scalar potentialDiffW5 = 0;
351 Scalar potentialDiffW6 = 0;
352 Scalar potentialDiffW7 = 0;
353 Scalar potentialDiffW8 = 0;
354 Scalar potentialDiffW9 = 0;
355 Scalar potentialDiffW10 = 0;
356 Scalar potentialDiffW11 = 0;
358 Scalar potentialDiffNw0 = 0;
359 Scalar potentialDiffNw1 = 0;
360 Scalar potentialDiffNw2 = 0;
361 Scalar potentialDiffNw3 = 0;
362 Scalar potentialDiffNw4 = 0;
363 Scalar potentialDiffNw5 = 0;
364 Scalar potentialDiffNw6 = 0;
365 Scalar potentialDiffNw7 = 0;
366 Scalar potentialDiffNw8 = 0;
367 Scalar potentialDiffNw9 = 0;
368 Scalar potentialDiffNw10 = 0;
369 Scalar potentialDiffNw11 = 0;
372 Dune::FieldVector<Scalar, 12> fluxW(0);
373 Dune::FieldVector<Scalar, 12> fluxNw(0);
376 Dune::FieldVector<Scalar, 2 * dim - dim + 1> u(0);
377 TransmissibilityType T(0);
378 TransmissibilityType TSecond(0);
379 Dune::FieldVector<bool, 4> useCases(
false);
381 int hangingNodeType = interactionVolume.getHangingNodeType();
385 if (fIdx < 0 || fIdx == 0)
388 int caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 0, 1, 2, 3,
401 potentialDiffW0 = Tu[0];
411 potentialDiffNw0 = Tu[0];
423 potentialDiffW0 = Tu[0];
433 potentialDiffNw0 = Tu[0];
445 potentialDiffW0 = Tu[0];
455 potentialDiffNw0 = Tu[0];
467 potentialDiffW0 = Tu[0];
477 potentialDiffNw0 = Tu[0];
481 if (fIdx < 0 || fIdx == 1)
486 if (hangingNodeType == InteractionVolume::twoSmallCells
487 || hangingNodeType == InteractionVolume::fourSmallCellsDiag)
493 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 1, 3, 0,
498 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 1, 3, 0,
512 potentialDiffW1 = Tu[0];
522 potentialDiffNw1 = Tu[0];
534 potentialDiffW1 = Tu[0];
544 potentialDiffNw1 = Tu[0];
556 potentialDiffW1 = Tu[0];
566 potentialDiffNw1 = Tu[0];
578 potentialDiffW1 = Tu[0];
588 potentialDiffNw1 = Tu[0];
592 if (fIdx < 0 || fIdx == 2)
597 if (hangingNodeType != InteractionVolume::twoSmallCells
598 && hangingNodeType != InteractionVolume::fourSmallCellsDiag)
600 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 3,
613 potentialDiffW2 = Tu[0];
623 potentialDiffNw2 = Tu[0];
635 potentialDiffW2 = Tu[0];
645 potentialDiffNw2 = Tu[0];
657 potentialDiffW2 = Tu[0];
667 potentialDiffNw2 = Tu[0];
679 potentialDiffW2 = Tu[0];
689 potentialDiffNw2 = Tu[0];
694 if (fIdx < 0 || fIdx == 3)
699 if (hangingNodeType == InteractionVolume::twoSmallCells
700 || hangingNodeType == InteractionVolume::fourSmallCellsDiag)
706 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 2, 0, 3, 1, 6,
711 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 2,
725 potentialDiffW3 = Tu[0];
735 potentialDiffNw3 = Tu[0];
747 potentialDiffW3 = Tu[0];
757 potentialDiffNw3 = Tu[0];
769 potentialDiffW3 = Tu[0];
779 potentialDiffNw3 = Tu[0];
791 potentialDiffW3 = Tu[0];
801 potentialDiffNw3 = Tu[0];
805 if (fIdx < 0 || fIdx == 4)
812 if (hangingNodeType == InteractionVolume::sixSmallCells)
814 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 5, 4, 7,
817 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx6)
819 caseL = transmissibilityCalculator_().transmissibilityCaseThree(T, interactionVolume, lambda, 4, 0, 2,
822 int caseLSecond = transmissibilityCalculator_().transmissibilityCaseFour(TSecond, interactionVolume, lambda, 1, 5, 3,
825 caseL = transmissibilityCalculator_().chooseTransmissibility(T, TSecond, 3, 4);
827 if (caseL == caseLSecond)
830 if (hangingNodeType == InteractionVolume::sixSmallCells)
842 potentialDiffW4 = Tu[0];
852 potentialDiffNw4 = Tu[0];
864 potentialDiffW4 = Tu[0];
874 potentialDiffNw4 = Tu[0];
886 potentialDiffW4 = Tu[0];
896 potentialDiffNw4 = Tu[0];
908 potentialDiffW4 = Tu[0];
918 potentialDiffNw4 = Tu[0];
921 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx6)
933 potentialDiffW4 = -Tu[2];
943 potentialDiffNw4 = -Tu[2];
955 potentialDiffW4 = Tu[2];
965 potentialDiffNw4 = Tu[2];
970 if (fIdx < 0 || fIdx == 5)
976 if (hangingNodeType == InteractionVolume::sixSmallCells)
983 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 7, 5, 6, 4, 3,
986 else if (hangingNodeType == InteractionVolume::fourSmallCellsDiag)
992 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 7, 5, 6, 4, 3,
995 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx7)
997 caseL = transmissibilityCalculator_().transmissibilityCaseThree(T, interactionVolume, lambda, 1, 5, 7, 0);
999 int caseLSecond = transmissibilityCalculator_().transmissibilityCaseFour(TSecond, interactionVolume, lambda, 7, 3, 5,
1002 caseL = transmissibilityCalculator_().chooseTransmissibility(T, TSecond, 3, 4);
1004 if (caseL == caseLSecond)
1008 if (hangingNodeType == InteractionVolume::sixSmallCells
1009 || hangingNodeType == InteractionVolume::fourSmallCellsDiag)
1021 potentialDiffW5 = Tu[0];
1031 potentialDiffNw5 = Tu[0];
1033 else if (caseL == 2)
1043 potentialDiffW5 = Tu[0];
1053 potentialDiffNw5 = Tu[0];
1055 else if (caseL == 3)
1065 potentialDiffW5 = Tu[0];
1075 potentialDiffNw5 = Tu[0];
1087 potentialDiffW5 = Tu[0];
1097 potentialDiffNw5 = Tu[0];
1100 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx7)
1112 potentialDiffW5 = -Tu[1];
1122 potentialDiffNw5 = -Tu[1];
1134 potentialDiffW5 = Tu[1];
1144 potentialDiffNw5 = Tu[1];
1149 if (fIdx < 0 || fIdx == 6)
1155 if (hangingNodeType == InteractionVolume::fourSmallCellsDiag)
1157 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 6, 7, 4,
1160 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx6)
1162 caseL = transmissibilityCalculator_().transmissibilityCaseThree(T, interactionVolume, lambda, 7, 3, 1, 6);
1164 int caseLSecond = transmissibilityCalculator_().transmissibilityCaseFour(TSecond, interactionVolume, lambda, 2, 6, 0,
1167 caseL = transmissibilityCalculator_().chooseTransmissibility(T, TSecond, 3, 4);
1169 if (caseL == caseLSecond)
1173 if (hangingNodeType == InteractionVolume::fourSmallCellsDiag)
1185 potentialDiffW6 = Tu[0];
1195 potentialDiffNw6 = Tu[0];
1197 else if (caseL == 2)
1207 potentialDiffW6 = Tu[0];
1217 potentialDiffNw6 = Tu[0];
1219 else if (caseL == 3)
1229 potentialDiffW6 = Tu[0];
1239 potentialDiffNw6 = Tu[0];
1251 potentialDiffW6 = Tu[0];
1261 potentialDiffNw6 = Tu[0];
1264 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx6)
1276 potentialDiffW6 = -Tu[2];
1286 potentialDiffNw6 = -Tu[2];
1288 else if (caseL == 4)
1298 potentialDiffW6 = Tu[2];
1308 potentialDiffNw6 = Tu[2];
1313 if (fIdx < 0 || fIdx == 7)
1319 if (hangingNodeType == InteractionVolume::sixSmallCells)
1322 useCases[1] =
false;
1323 useCases[2] =
false;
1326 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 4, 6, 5, 7, 0,
1329 else if (hangingNodeType == InteractionVolume::fourSmallCellsDiag)
1331 useCases[0] =
false;
1334 useCases[3] =
false;
1335 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 4, 6, 5, 7, 0,
1338 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx7)
1340 caseL = transmissibilityCalculator_().transmissibilityCaseThree(T, interactionVolume, lambda, 2, 6, 4,
1343 int caseLSecond = transmissibilityCalculator_().transmissibilityCaseFour(TSecond, interactionVolume, lambda, 4, 0, 6,
1348 caseL = transmissibilityCalculator_().chooseTransmissibility(T, TSecond, 3, 4);
1350 if (caseL == caseLSecond)
1354 if (hangingNodeType == InteractionVolume::sixSmallCells
1355 || hangingNodeType == InteractionVolume::fourSmallCellsDiag)
1368 potentialDiffW7 = Tu[0];
1378 potentialDiffNw7 = Tu[0];
1380 else if (caseL == 2)
1390 potentialDiffW7 = Tu[0];
1400 potentialDiffNw7 = Tu[0];
1402 else if (caseL == 3)
1412 potentialDiffW7 = Tu[0];
1422 potentialDiffNw7 = Tu[0];
1434 potentialDiffW7 = Tu[0];
1444 potentialDiffNw7 = Tu[0];
1447 else if(hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx7)
1459 potentialDiffW7 = -Tu[1];
1469 potentialDiffNw7 = -Tu[1];
1471 else if (caseL == 4)
1481 potentialDiffW7 = Tu[1];
1491 potentialDiffNw7 = Tu[1];
1496 if (fIdx < 0 || fIdx == 8)
1501 if (hangingNodeType == InteractionVolume::sixSmallCells)
1503 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 4, 0, 6,
1506 else if (hangingNodeType == InteractionVolume::twoSmallCells
1507 || hangingNodeType == InteractionVolume::fourSmallCellsFace)
1509 caseL = transmissibilityCalculator_().transmissibilityCaseTwo(T, interactionVolume, lambda, 4, 0, 2,
1512 else if (hangingNodeType == InteractionVolume::fourSmallCellsDiag
1513 || (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx7))
1515 useCases[0] =
false;
1517 useCases[2] =
false;
1520 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 4, 0, 6,
1523 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx6)
1525 useCases[0] =
false;
1528 useCases[3] =
false;
1530 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 4, 0, 6,
1544 potentialDiffW8 = Tu[0];
1554 potentialDiffNw8 = Tu[0];
1556 else if (caseL == 2)
1566 potentialDiffW8 = Tu[0];
1576 potentialDiffNw8 = Tu[0];
1578 else if (caseL == 3)
1588 potentialDiffW8 = Tu[0];
1598 potentialDiffNw8 = Tu[0];
1610 potentialDiffW8 = Tu[0];
1620 potentialDiffNw8 = Tu[0];
1624 if (fIdx < 0 || fIdx == 9)
1629 if (hangingNodeType == InteractionVolume::sixSmallCells)
1631 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 1, 5, 3,
1634 else if (hangingNodeType == InteractionVolume::twoSmallCells
1635 || hangingNodeType == InteractionVolume::fourSmallCellsFace)
1637 caseL = transmissibilityCalculator_().transmissibilityCaseOne(T, interactionVolume, lambda, 1, 5, 3,
1640 else if (hangingNodeType == InteractionVolume::fourSmallCellsDiag
1641 || (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx7))
1644 useCases[1] =
false;
1646 useCases[3] =
false;
1648 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 1, 5, 3,
1651 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx6)
1654 useCases[1] =
false;
1655 useCases[2] =
false;
1658 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 1, 5, 3,
1672 potentialDiffW9 = Tu[0];
1682 potentialDiffNw9 = Tu[0];
1684 else if (caseL == 2)
1694 potentialDiffW9 = Tu[0];
1704 potentialDiffNw9 = Tu[0];
1706 else if (caseL == 3)
1716 potentialDiffW9 = Tu[0];
1726 potentialDiffNw9 = Tu[0];
1738 potentialDiffW9 = Tu[0];
1748 potentialDiffNw9 = Tu[0];
1752 if (fIdx < 0 || fIdx == 10)
1757 if (hangingNodeType == InteractionVolume::fourSmallCellsFace)
1759 caseL = transmissibilityCalculator_().transmissibilityCaseTwo(T, interactionVolume, lambda, 7, 3, 1, 2);
1761 else if ((hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx7)
1762 || hangingNodeType == InteractionVolume::sixSmallCells)
1764 useCases[0] =
false;
1766 useCases[2] =
false;
1769 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 7, 3, 5, 1, 6,
1772 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx6)
1774 useCases[0] =
false;
1777 useCases[3] =
false;
1779 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 7, 3, 5, 1, 6,
1782 else if (hangingNodeType == InteractionVolume::fourSmallCellsDiag)
1785 useCases[1] =
false;
1787 useCases[3] =
false;
1789 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 7, 3, 5, 1, 6,
1793 if (hangingNodeType != InteractionVolume::twoSmallCells)
1806 potentialDiffW10 = Tu[0];
1816 potentialDiffNw10 = Tu[0];
1818 else if (caseL == 2)
1828 potentialDiffW10 = Tu[0];
1838 potentialDiffNw10 = Tu[0];
1840 else if (caseL == 3)
1850 potentialDiffW10 = Tu[0];
1860 potentialDiffNw10 = Tu[0];
1872 potentialDiffW10 = Tu[0];
1882 potentialDiffNw10 = Tu[0];
1887 if (fIdx < 0 || fIdx == 11)
1892 if (hangingNodeType == InteractionVolume::fourSmallCellsFace)
1894 caseL = transmissibilityCalculator_().transmissibilityCaseOne(T, interactionVolume, lambda, 2, 6, 0, 3);
1896 else if ((hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx7)
1897 || hangingNodeType == InteractionVolume::sixSmallCells)
1900 useCases[1] =
false;
1902 useCases[3] =
false;
1904 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 2, 6, 0, 4, 3,
1907 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx6)
1910 useCases[1] =
false;
1911 useCases[2] =
false;
1914 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 2, 6, 0, 4, 3,
1917 else if (hangingNodeType == InteractionVolume::fourSmallCellsDiag)
1919 useCases[0] =
false;
1921 useCases[2] =
false;
1924 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 2, 6, 0, 4, 3,
1928 if (hangingNodeType != InteractionVolume::twoSmallCells)
1941 potentialDiffW11 = Tu[0];
1951 potentialDiffNw11 = Tu[0];
1953 else if (caseL == 2)
1963 potentialDiffW11 = Tu[0];
1973 potentialDiffNw11 = Tu[0];
1975 else if (caseL == 3)
1985 potentialDiffW11 = Tu[0];
1995 potentialDiffNw11 = Tu[0];
2007 potentialDiffW11 = Tu[0];
2017 potentialDiffNw11 = Tu[0];
2022 cellData1.fluxData().addUpwindPotential(wPhaseIdx, interactionVolume.getIndexOnElement(0, 0), potentialDiffW0);
2023 cellData1.fluxData().addUpwindPotential(nPhaseIdx, interactionVolume.getIndexOnElement(0, 0), potentialDiffNw0);
2024 cellData1.fluxData().addUpwindPotential(wPhaseIdx, interactionVolume.getIndexOnElement(0, 1), -potentialDiffW3);
2025 cellData1.fluxData().addUpwindPotential(nPhaseIdx, interactionVolume.getIndexOnElement(0, 1), -potentialDiffNw3);
2026 cellData1.fluxData().addUpwindPotential(wPhaseIdx, interactionVolume.getIndexOnElement(0, 2), -potentialDiffW8);
2027 cellData1.fluxData().addUpwindPotential(nPhaseIdx, interactionVolume.getIndexOnElement(0, 2), -potentialDiffNw8);
2029 cellData2.fluxData().addUpwindPotential(wPhaseIdx, interactionVolume.getIndexOnElement(1, 0), potentialDiffW1);
2030 cellData2.fluxData().addUpwindPotential(nPhaseIdx, interactionVolume.getIndexOnElement(1, 0), potentialDiffNw1);
2031 cellData2.fluxData().addUpwindPotential(wPhaseIdx, interactionVolume.getIndexOnElement(1, 1), -potentialDiffW0);
2032 cellData2.fluxData().addUpwindPotential(nPhaseIdx, interactionVolume.getIndexOnElement(1, 1), -potentialDiffNw0);
2033 cellData2.fluxData().addUpwindPotential(wPhaseIdx, interactionVolume.getIndexOnElement(1, 2), potentialDiffW9);
2034 cellData2.fluxData().addUpwindPotential(nPhaseIdx, interactionVolume.getIndexOnElement(1, 2), potentialDiffNw9);
2036 cellData3.fluxData().addUpwindPotential(wPhaseIdx, interactionVolume.getIndexOnElement(2, 0), potentialDiffW3);
2037 cellData3.fluxData().addUpwindPotential(nPhaseIdx, interactionVolume.getIndexOnElement(2, 0), potentialDiffNw3);
2038 cellData3.fluxData().addUpwindPotential(wPhaseIdx, interactionVolume.getIndexOnElement(2, 1), -potentialDiffW2);
2039 cellData3.fluxData().addUpwindPotential(nPhaseIdx, interactionVolume.getIndexOnElement(2, 1), -potentialDiffNw2);
2040 cellData3.fluxData().addUpwindPotential(wPhaseIdx, interactionVolume.getIndexOnElement(2, 2), potentialDiffW11);
2041 cellData3.fluxData().addUpwindPotential(nPhaseIdx, interactionVolume.getIndexOnElement(2, 2), potentialDiffNw11);
2043 cellData4.fluxData().addUpwindPotential(wPhaseIdx, interactionVolume.getIndexOnElement(3, 0), potentialDiffW2);
2044 cellData4.fluxData().addUpwindPotential(nPhaseIdx, interactionVolume.getIndexOnElement(3, 0), potentialDiffNw2);
2045 cellData4.fluxData().addUpwindPotential(wPhaseIdx, interactionVolume.getIndexOnElement(3, 1), -potentialDiffW1);
2046 cellData4.fluxData().addUpwindPotential(nPhaseIdx, interactionVolume.getIndexOnElement(3, 1), -potentialDiffNw1);
2047 cellData4.fluxData().addUpwindPotential(wPhaseIdx, interactionVolume.getIndexOnElement(3, 2), -potentialDiffW10);
2048 cellData4.fluxData().addUpwindPotential(nPhaseIdx, interactionVolume.getIndexOnElement(3, 2), -potentialDiffNw10);
2050 cellData5.fluxData().addUpwindPotential(wPhaseIdx, interactionVolume.getIndexOnElement(4, 0), potentialDiffW8);
2051 cellData5.fluxData().addUpwindPotential(nPhaseIdx, interactionVolume.getIndexOnElement(4, 0), potentialDiffNw8);
2052 cellData5.fluxData().addUpwindPotential(wPhaseIdx, interactionVolume.getIndexOnElement(4, 1), -potentialDiffW4);
2053 cellData5.fluxData().addUpwindPotential(nPhaseIdx, interactionVolume.getIndexOnElement(4, 1), -potentialDiffNw4);
2054 cellData5.fluxData().addUpwindPotential(wPhaseIdx, interactionVolume.getIndexOnElement(4, 2), potentialDiffW7);
2055 cellData5.fluxData().addUpwindPotential(nPhaseIdx, interactionVolume.getIndexOnElement(4, 2), potentialDiffNw7);
2057 cellData6.fluxData().addUpwindPotential(wPhaseIdx, interactionVolume.getIndexOnElement(5, 0), -potentialDiffW9);
2058 cellData6.fluxData().addUpwindPotential(nPhaseIdx, interactionVolume.getIndexOnElement(5, 0), -potentialDiffNw9);
2059 cellData6.fluxData().addUpwindPotential(wPhaseIdx, interactionVolume.getIndexOnElement(5, 1), -potentialDiffW5);
2060 cellData6.fluxData().addUpwindPotential(nPhaseIdx, interactionVolume.getIndexOnElement(5, 1), -potentialDiffNw5);
2061 cellData6.fluxData().addUpwindPotential(wPhaseIdx, interactionVolume.getIndexOnElement(5, 2), potentialDiffW4);
2062 cellData6.fluxData().addUpwindPotential(nPhaseIdx, interactionVolume.getIndexOnElement(5, 2), potentialDiffNw4);
2064 cellData7.fluxData().addUpwindPotential(wPhaseIdx, interactionVolume.getIndexOnElement(6, 0), -potentialDiffW11);
2065 cellData7.fluxData().addUpwindPotential(nPhaseIdx, interactionVolume.getIndexOnElement(6, 0), -potentialDiffNw11);
2066 cellData7.fluxData().addUpwindPotential(wPhaseIdx, interactionVolume.getIndexOnElement(6, 1), -potentialDiffW7);
2067 cellData7.fluxData().addUpwindPotential(nPhaseIdx, interactionVolume.getIndexOnElement(6, 1), -potentialDiffNw7);
2068 cellData7.fluxData().addUpwindPotential(wPhaseIdx, interactionVolume.getIndexOnElement(6, 2), potentialDiffW6);
2069 cellData7.fluxData().addUpwindPotential(nPhaseIdx, interactionVolume.getIndexOnElement(6, 2), potentialDiffNw6);
2071 cellData8.fluxData().addUpwindPotential(wPhaseIdx, interactionVolume.getIndexOnElement(7, 0), potentialDiffW10);
2072 cellData8.fluxData().addUpwindPotential(nPhaseIdx, interactionVolume.getIndexOnElement(7, 0), potentialDiffNw10);
2073 cellData8.fluxData().addUpwindPotential(wPhaseIdx, interactionVolume.getIndexOnElement(7, 1), -potentialDiffW6);
2074 cellData8.fluxData().addUpwindPotential(nPhaseIdx, interactionVolume.getIndexOnElement(7, 1), -potentialDiffNw6);
2075 cellData8.fluxData().addUpwindPotential(wPhaseIdx, interactionVolume.getIndexOnElement(7, 2), potentialDiffW5);
2076 cellData8.fluxData().addUpwindPotential(nPhaseIdx, interactionVolume.getIndexOnElement(7, 2), potentialDiffNw5);
2079 Dune::FieldVector<Scalar, numPhases> lambda0Upw(0.0);
2080 lambda0Upw[wPhaseIdx] = (potentialDiffW0 >= 0) ? lambda1[wPhaseIdx] : lambda2[wPhaseIdx];
2081 lambda0Upw[nPhaseIdx] = (potentialDiffNw0 >= 0) ? lambda1[nPhaseIdx] : lambda2[nPhaseIdx];
2084 Dune::FieldVector<Scalar, numPhases> lambda1Upw(0.0);
2085 lambda1Upw[wPhaseIdx] = (potentialDiffW1 >= 0) ? lambda2[wPhaseIdx] : lambda4[wPhaseIdx];
2086 lambda1Upw[nPhaseIdx] = (potentialDiffNw1 >= 0) ? lambda2[nPhaseIdx] : lambda4[nPhaseIdx];
2089 Dune::FieldVector<Scalar, numPhases> lambda2Upw(0.0);
2090 lambda2Upw[wPhaseIdx] = (potentialDiffW2 >= 0) ? lambda4[wPhaseIdx] : lambda3[wPhaseIdx];
2091 lambda2Upw[nPhaseIdx] = (potentialDiffNw2 >= 0) ? lambda4[nPhaseIdx] : lambda3[nPhaseIdx];
2094 Dune::FieldVector<Scalar, numPhases> lambda3Upw(0.0);
2095 lambda3Upw[wPhaseIdx] = (potentialDiffW3 >= 0) ? lambda3[wPhaseIdx] : lambda1[wPhaseIdx];
2096 lambda3Upw[nPhaseIdx] = (potentialDiffNw3 >= 0) ? lambda3[nPhaseIdx] : lambda1[nPhaseIdx];
2099 Dune::FieldVector<Scalar, numPhases> lambda4Upw(0.0);
2100 lambda4Upw[wPhaseIdx] = (potentialDiffW4 >= 0) ? lambda6[wPhaseIdx] : lambda5[wPhaseIdx];
2101 lambda4Upw[nPhaseIdx] = (potentialDiffNw4 >= 0) ? lambda6[nPhaseIdx] : lambda5[nPhaseIdx];
2104 Dune::FieldVector<Scalar, numPhases> lambda5Upw(0.0);
2105 lambda5Upw[wPhaseIdx] = (potentialDiffW5 >= 0) ? lambda8[wPhaseIdx] : lambda6[wPhaseIdx];
2106 lambda5Upw[nPhaseIdx] = (potentialDiffNw5 >= 0) ? lambda8[nPhaseIdx] : lambda6[nPhaseIdx];
2109 Dune::FieldVector<Scalar, numPhases> lambda6Upw(0.0);
2110 lambda6Upw[wPhaseIdx] = (potentialDiffW6 >= 0) ? lambda7[wPhaseIdx] : lambda8[wPhaseIdx];
2111 lambda6Upw[nPhaseIdx] = (potentialDiffNw6 >= 0) ? lambda7[nPhaseIdx] : lambda8[nPhaseIdx];
2114 Dune::FieldVector<Scalar, numPhases> lambda7Upw(0.0);
2115 lambda7Upw[wPhaseIdx] = (potentialDiffW7 >= 0) ? lambda5[wPhaseIdx] : lambda7[wPhaseIdx];
2116 lambda7Upw[nPhaseIdx] = (potentialDiffNw7 >= 0) ? lambda5[nPhaseIdx] : lambda7[nPhaseIdx];
2119 Dune::FieldVector<Scalar, numPhases> lambda8Upw(0.0);
2120 lambda8Upw[wPhaseIdx] = (potentialDiffW8 >= 0) ? lambda5[wPhaseIdx] : lambda1[wPhaseIdx];
2121 lambda8Upw[nPhaseIdx] = (potentialDiffNw8 >= 0) ? lambda5[nPhaseIdx] : lambda1[nPhaseIdx];
2124 Dune::FieldVector<Scalar, numPhases> lambda9Upw(0.0);
2125 lambda9Upw[wPhaseIdx] = (potentialDiffW9 >= 0) ? lambda2[wPhaseIdx] : lambda6[wPhaseIdx];
2126 lambda9Upw[nPhaseIdx] = (potentialDiffNw9 >= 0) ? lambda2[nPhaseIdx] : lambda6[nPhaseIdx];
2129 Dune::FieldVector<Scalar, numPhases> lambda10Upw(0.0);
2130 lambda10Upw[wPhaseIdx] = (potentialDiffW10 >= 0) ? lambda8[wPhaseIdx] : lambda4[wPhaseIdx];
2131 lambda10Upw[nPhaseIdx] = (potentialDiffNw10 >= 0) ? lambda8[nPhaseIdx] : lambda4[nPhaseIdx];
2134 Dune::FieldVector<Scalar, numPhases> lambda11Upw(0.0);
2135 lambda11Upw[wPhaseIdx] = (potentialDiffW11 >= 0) ? lambda3[wPhaseIdx] : lambda7[wPhaseIdx];
2136 lambda11Upw[nPhaseIdx] = (potentialDiffNw11 >= 0) ? lambda3[nPhaseIdx] : lambda7[nPhaseIdx];
2138 for (
int i = 0; i < numPhases; i++)
2141 DimVector vel12 = interactionVolume.getNormal(0, 0);
2142 DimVector vel21 = interactionVolume.getNormal(1, 1);
2143 DimVector vel24 = interactionVolume.getNormal(1, 0);
2144 DimVector vel42 = interactionVolume.getNormal(3, 1);
2145 DimVector vel43 = interactionVolume.getNormal(3, 0);
2146 DimVector vel34 = interactionVolume.getNormal(2, 1);
2147 DimVector vel31 = interactionVolume.getNormal(2, 0);
2148 DimVector vel13 = interactionVolume.getNormal(0, 1);
2149 DimVector vel65 = interactionVolume.getNormal(5, 2);
2150 DimVector vel56 = interactionVolume.getNormal(4, 1);
2151 DimVector vel86 = interactionVolume.getNormal(7, 2);
2152 DimVector vel68 = interactionVolume.getNormal(5, 1);
2153 DimVector vel78 = interactionVolume.getNormal(6, 2);
2154 DimVector vel87 = interactionVolume.getNormal(7, 1);
2155 DimVector vel57 = interactionVolume.getNormal(4, 2);
2156 DimVector vel75 = interactionVolume.getNormal(6, 1);
2157 DimVector vel51 = interactionVolume.getNormal(4, 0);
2158 DimVector vel15 = interactionVolume.getNormal(0, 2);
2159 DimVector vel26 = interactionVolume.getNormal(1, 2);
2160 DimVector vel62 = interactionVolume.getNormal(5, 0);
2161 DimVector vel84 = interactionVolume.getNormal(7, 0);
2162 DimVector vel48 = interactionVolume.getNormal(3, 2);
2163 DimVector vel37 = interactionVolume.getNormal(2, 2);
2164 DimVector vel73 = interactionVolume.getNormal(6, 0);
2178 Dune::FieldVector<Scalar, 12> flux(0);
2193 if (hangingNodeType == InteractionVolume::sixSmallCells)
2195 vel12 *= flux[0] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx1, 0, 0));
2196 vel21 *= flux[0] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx2, 1, 1));
2197 vel24 *= flux[1] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx2, 1, 0));
2198 vel42 *= flux[1] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx4, 3, 1));
2199 vel43 *= flux[2] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx4, 3, 0));
2200 vel34 *= flux[2] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx3, 2, 1));
2201 vel31 *= flux[3] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx3, 2, 0));
2202 vel13 *= flux[3] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx1, 0, 1));
2203 vel65 *= flux[4] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx6, 5, 2));
2204 vel56 *= flux[4] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx5, 4, 1));
2205 vel86 *= flux[5] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx8, 7, 2));
2206 vel68 *= flux[5] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx6, 5, 1));
2209 vel57 *= flux[7] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx5, 4, 2));
2210 vel75 *= flux[7] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx7, 6, 1));
2211 vel51 *= flux[8] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx5, 4, 0));
2212 vel15 *= flux[8] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx1, 0, 2));
2213 vel26 *= flux[9] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx2, 1, 2));
2214 vel62 *= flux[9] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx6, 5, 0));
2215 vel84 *= flux[10] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx8, 7, 0));
2216 vel48 *= flux[10] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx4, 3, 2));
2217 vel37 *= flux[11] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx3, 2, 2));
2218 vel73 *= flux[11] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx7, 6, 0));
2220 else if (hangingNodeType == InteractionVolume::fourSmallCellsDiag)
2222 vel12 *= flux[0] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx1, 0, 0));
2223 vel21 *= flux[0] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx2, 1, 1));
2224 vel24 *= flux[1] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx2, 1, 0));
2225 vel42 *= flux[1] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx4, 3, 1));
2228 vel31 *= flux[3] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx3, 2, 0));
2229 vel13 *= flux[3] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx1, 0, 1));
2232 vel86 *= flux[5] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx8, 7, 2));
2233 vel68 *= flux[5] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx6, 5, 1));
2234 vel78 *= flux[6] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx7, 6, 2));
2235 vel87 *= flux[6] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx8, 7, 1));
2236 vel57 *= flux[7] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx5, 4, 2));
2237 vel75 *= flux[7] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx7, 6, 1));
2238 vel51 *= flux[8] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx5, 4, 0));
2239 vel15 *= flux[8] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx1, 0, 2));
2240 vel26 *= flux[9] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx2, 1, 2));
2241 vel62 *= flux[9] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx6, 5, 0));
2242 vel84 *= flux[10] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx8, 7, 0));
2243 vel48 *= flux[10] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx4, 3, 2));
2244 vel37 *= flux[11] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx3, 2, 2));
2245 vel73 *= flux[11] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx7, 6, 0));
2247 else if (hangingNodeType == InteractionVolume::fourSmallCellsFace
2248 || hangingNodeType == InteractionVolume::fourSmallCellsEdge)
2250 vel12 *= flux[0] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx1, 0, 0));
2251 vel21 *= flux[0] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx2, 1, 1));
2252 vel24 *= flux[1] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx2, 1, 0));
2253 vel42 *= flux[1] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx4, 3, 1));
2254 vel43 *= flux[2] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx4, 3, 0));
2255 vel34 *= flux[2] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx3, 2, 1));
2256 vel31 *= flux[3] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx3, 2, 0));
2257 vel13 *= flux[3] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx1, 0, 1));
2258 vel51 *= flux[8] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx5, 4, 0));
2259 vel15 *= flux[8] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx1, 0, 2));
2260 vel26 *= flux[9] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx2, 1, 2));
2261 vel62 *= flux[9] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx6, 5, 0));
2262 vel84 *= flux[10] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx8, 7, 0));
2263 vel48 *= flux[10] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx4, 3, 2));
2264 vel37 *= flux[11] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx3, 2, 2));
2265 vel73 *= flux[11] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx7, 6, 0));
2267 if (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx7)
2271 vel86 *= flux[5] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx8, 7, 2));
2272 vel68 *= flux[5] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx6, 5, 1));
2275 vel57 *= flux[7] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx5, 4, 2));
2276 vel75 *= flux[7] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx7, 6, 1));
2278 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx6)
2280 vel65 *= flux[4] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx6, 5, 2));
2281 vel56 *= flux[4] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx5, 4, 1));
2284 vel78 *= flux[6] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx7, 6, 2));
2285 vel87 *= flux[6] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx8, 7, 1));
2301 else if (hangingNodeType == InteractionVolume::twoSmallCells)
2303 vel12 *= flux[0] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx1, 0, 0));
2304 vel21 *= flux[0] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx2, 1, 1));
2305 vel24 *= flux[1] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx2, 1, 0));
2306 vel42 *= flux[1] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx4, 3, 1));
2309 vel31 *= flux[3] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx3, 2, 0));
2310 vel13 *= flux[3] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx1, 0, 1));
2319 vel51 *= flux[8] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx5, 4, 0));
2320 vel15 *= flux[8] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx1, 0, 2));
2321 vel26 *= flux[9] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx2, 1, 2));
2322 vel62 *= flux[9] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx6, 5, 0));
2329 Scalar lambdaT0 = lambda0Upw[wPhaseIdx] + lambda0Upw[nPhaseIdx];
2330 Scalar lambdaT1 = lambda1Upw[wPhaseIdx] + lambda1Upw[nPhaseIdx];
2331 Scalar lambdaT2 = lambda2Upw[wPhaseIdx] + lambda2Upw[nPhaseIdx];
2332 Scalar lambdaT3 = lambda3Upw[wPhaseIdx] + lambda3Upw[nPhaseIdx];
2333 Scalar lambdaT4 = lambda4Upw[wPhaseIdx] + lambda4Upw[nPhaseIdx];
2334 Scalar lambdaT5 = lambda5Upw[wPhaseIdx] + lambda5Upw[nPhaseIdx];
2335 Scalar lambdaT6 = lambda6Upw[wPhaseIdx] + lambda6Upw[nPhaseIdx];
2336 Scalar lambdaT7 = lambda7Upw[wPhaseIdx] + lambda7Upw[nPhaseIdx];
2337 Scalar lambdaT8 = lambda8Upw[wPhaseIdx] + lambda8Upw[nPhaseIdx];
2338 Scalar lambdaT9 = lambda9Upw[wPhaseIdx] + lambda9Upw[nPhaseIdx];
2339 Scalar lambdaT10 = lambda10Upw[wPhaseIdx] + lambda10Upw[nPhaseIdx];
2340 Scalar lambdaT11 = lambda11Upw[wPhaseIdx] + lambda11Upw[nPhaseIdx];
2341 Scalar fracFlow0 = (lambdaT0 > threshold_) ? lambda0Upw[i] / (lambdaT0) : 0.0;
2342 Scalar fracFlow1 = (lambdaT1 > threshold_) ? lambda1Upw[i] / (lambdaT1) : 0.0;
2343 Scalar fracFlow2 = (lambdaT2 > threshold_) ? lambda2Upw[i] / (lambdaT2) : 0.0;
2344 Scalar fracFlow3 = (lambdaT3 > threshold_) ? lambda3Upw[i] / (lambdaT3) : 0.0;
2345 Scalar fracFlow4 = (lambdaT4 > threshold_) ? lambda4Upw[i] / (lambdaT4) : 0.0;
2346 Scalar fracFlow5 = (lambdaT5 > threshold_) ? lambda5Upw[i] / (lambdaT5) : 0.0;
2347 Scalar fracFlow6 = (lambdaT6 > threshold_) ? lambda6Upw[i] / (lambdaT6) : 0.0;
2348 Scalar fracFlow7 = (lambdaT7 > threshold_) ? lambda7Upw[i] / (lambdaT7) : 0.0;
2349 Scalar fracFlow8 = (lambdaT8 > threshold_) ? lambda8Upw[i] / (lambdaT8) : 0.0;
2350 Scalar fracFlow9 = (lambdaT9 > threshold_) ? lambda9Upw[i] / (lambdaT9) : 0.0;
2351 Scalar fracFlow10 = (lambdaT10 > threshold_) ? lambda10Upw[i] / (lambdaT10) : 0.0;
2352 Scalar fracFlow11 = (lambdaT11 > threshold_) ? lambda11Upw[i] / (lambdaT11) : 0.0;
2374 vel84 *= fracFlow10;
2375 vel48 *= fracFlow10;
2376 vel37 *= fracFlow11;
2377 vel73 *= fracFlow11;
2380 cellData1.fluxData().addVelocity(i, interactionVolume.getIndexOnElement(0, 0), vel12);
2381 cellData1.fluxData().addVelocity(i, interactionVolume.getIndexOnElement(0, 1), vel13);
2382 cellData1.fluxData().addVelocity(i, interactionVolume.getIndexOnElement(0, 2), vel15);
2383 cellData2.fluxData().addVelocity(i, interactionVolume.getIndexOnElement(1, 0), vel24);
2384 cellData2.fluxData().addVelocity(i, interactionVolume.getIndexOnElement(1, 1), vel21);
2385 cellData2.fluxData().addVelocity(i, interactionVolume.getIndexOnElement(1, 2), vel26);
2386 cellData3.fluxData().addVelocity(i, interactionVolume.getIndexOnElement(2, 0), vel31);
2387 cellData3.fluxData().addVelocity(i, interactionVolume.getIndexOnElement(2, 1), vel34);
2388 cellData3.fluxData().addVelocity(i, interactionVolume.getIndexOnElement(2, 2), vel37);
2389 cellData4.fluxData().addVelocity(i, interactionVolume.getIndexOnElement(3, 0), vel43);
2390 cellData4.fluxData().addVelocity(i, interactionVolume.getIndexOnElement(3, 1), vel42);
2391 cellData4.fluxData().addVelocity(i, interactionVolume.getIndexOnElement(3, 2), vel48);
2392 cellData5.fluxData().addVelocity(i, interactionVolume.getIndexOnElement(4, 0), vel51);
2393 cellData5.fluxData().addVelocity(i, interactionVolume.getIndexOnElement(4, 1), vel56);
2394 cellData5.fluxData().addVelocity(i, interactionVolume.getIndexOnElement(4, 2), vel57);
2395 cellData6.fluxData().addVelocity(i, interactionVolume.getIndexOnElement(5, 0), vel62);
2396 cellData6.fluxData().addVelocity(i, interactionVolume.getIndexOnElement(5, 1), vel68);
2397 cellData6.fluxData().addVelocity(i, interactionVolume.getIndexOnElement(5, 2), vel65);
2398 cellData7.fluxData().addVelocity(i, interactionVolume.getIndexOnElement(6, 0), vel73);
2399 cellData7.fluxData().addVelocity(i, interactionVolume.getIndexOnElement(6, 1), vel75);
2400 cellData7.fluxData().addVelocity(i, interactionVolume.getIndexOnElement(6, 2), vel78);
2401 cellData8.fluxData().addVelocity(i, interactionVolume.getIndexOnElement(7, 0), vel84);
2402 cellData8.fluxData().addVelocity(i, interactionVolume.getIndexOnElement(7, 1), vel87);
2403 cellData8.fluxData().addVelocity(i, interactionVolume.getIndexOnElement(7, 2), vel86);
2406 cellData1.fluxData().setVelocityMarker(interactionVolume.getIndexOnElement(0, 0));
2407 cellData1.fluxData().setVelocityMarker(interactionVolume.getIndexOnElement(0, 1));
2408 cellData1.fluxData().setVelocityMarker(interactionVolume.getIndexOnElement(0, 2));
2409 cellData2.fluxData().setVelocityMarker(interactionVolume.getIndexOnElement(1, 0));
2410 cellData2.fluxData().setVelocityMarker(interactionVolume.getIndexOnElement(1, 1));
2411 cellData2.fluxData().setVelocityMarker(interactionVolume.getIndexOnElement(1, 2));
2412 cellData3.fluxData().setVelocityMarker(interactionVolume.getIndexOnElement(2, 0));
2413 cellData3.fluxData().setVelocityMarker(interactionVolume.getIndexOnElement(2, 1));
2414 cellData3.fluxData().setVelocityMarker(interactionVolume.getIndexOnElement(2, 2));
2415 cellData4.fluxData().setVelocityMarker(interactionVolume.getIndexOnElement(3, 0));
2416 cellData4.fluxData().setVelocityMarker(interactionVolume.getIndexOnElement(3, 1));
2417 cellData4.fluxData().setVelocityMarker(interactionVolume.getIndexOnElement(3, 2));
2418 cellData5.fluxData().setVelocityMarker(interactionVolume.getIndexOnElement(4, 0));
2419 cellData5.fluxData().setVelocityMarker(interactionVolume.getIndexOnElement(4, 1));
2420 cellData5.fluxData().setVelocityMarker(interactionVolume.getIndexOnElement(4, 2));
2421 cellData6.fluxData().setVelocityMarker(interactionVolume.getIndexOnElement(5, 0));
2422 cellData6.fluxData().setVelocityMarker(interactionVolume.getIndexOnElement(5, 1));
2423 cellData6.fluxData().setVelocityMarker(interactionVolume.getIndexOnElement(5, 2));
2424 cellData7.fluxData().setVelocityMarker(interactionVolume.getIndexOnElement(6, 0));
2425 cellData7.fluxData().setVelocityMarker(interactionVolume.getIndexOnElement(6, 1));
2426 cellData7.fluxData().setVelocityMarker(interactionVolume.getIndexOnElement(6, 2));
2427 cellData8.fluxData().setVelocityMarker(interactionVolume.getIndexOnElement(7, 0));
2428 cellData8.fluxData().setVelocityMarker(interactionVolume.getIndexOnElement(7, 1));
2429 cellData8.fluxData().setVelocityMarker(interactionVolume.getIndexOnElement(7, 2));
2-d velocity calculation using a 3-d MPFA L-method.
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type GetProp
get the type of a property (equivalent to old macro GET_PROP(...))
Definition: propertysystem.hh:140
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
std::string viscosity(int phaseIdx) noexcept
I/O name of viscosity for multiphase systems.
Definition: name.hh:74
std::string density(int phaseIdx) noexcept
I/O name of density for multiphase systems.
Definition: name.hh:65
Provides methods for transmissibility calculation in 3-d.
Definition: 3dtransmissibilitycalculator.hh:51
Dune::FieldMatrix< Scalar, dim, 2 *dim - dim+1 > TransmissibilityType
Type of the transmissibility matrix.
Definition: 3dtransmissibilitycalculator.hh:78
Class for calculating 3-d velocities from cell-wise constant pressure values.
Definition: 3dvelocity.hh:57
void initialize(bool solveTwice=true)
Initializes the velocity model.
Definition: 3dvelocity.hh:160
Class for calculating 3-d velocities from cell-wise constant pressure values.
Definition: 3dvelocityadaptive.hh:54
FvMpfaL3dVelocity2pAdaptive(Problem &problem)
Constructs a FvMpfaL3dVelocity2pAdaptive object.
Definition: 3dvelocityadaptive.hh:148
void calculateHangingNodeInteractionVolumeVelocity(InteractionVolume &interactionVolume, CellData &cellData1, CellData &cellData2, CellData &cellData3, CellData &cellData4, CellData &cellData5, CellData &cellData6, CellData &cellData7, CellData &cellData8, int fIdx=-1)
Calculates velocities for flux faces of a hanging node interaction volume.
Definition: 3dvelocityadaptive.hh:218
void initialize()
Initializes the velocity model.
Definition: 3dvelocityadaptive.hh:163