24#ifndef DUMUX_FVMPFAL3DVELOCITY2P_ADAPTIVE_HH
25#define DUMUX_FVMPFAL3DVELOCITY2P_ADAPTIVE_HH
61 dim = GridView::dimension, dimWorld = GridView::dimensionworld
68 using MaterialLaw =
typename SpatialParams::MaterialLaw;
78 using PrimaryVariables =
typename SolutionTypes::PrimaryVariables;
82 using Element =
typename GridView::Traits::template Codim<0>::Entity;
83 using Grid =
typename GridView::Grid;
84 using IndexSet =
typename GridView::IndexSet;
96 pw = Indices::pressureW,
97 pn = Indices::pressureNw,
98 pGlobal = Indices::pressureGlobal,
99 sw = Indices::saturationW,
100 sn = Indices::saturationNw,
101 vw = Indices::velocityW,
102 vn = Indices::velocityNw,
103 vt = Indices::velocityTotal
107 wPhaseIdx = Indices::wPhaseIdx,
108 nPhaseIdx = Indices::nPhaseIdx,
109 pressureIdx = Indices::pressureIdx,
110 saturationIdx = Indices::saturationIdx,
111 pressureEqIdx = Indices::pressureEqIdx,
112 satEqIdx = Indices::satEqIdx,
113 numPhases = getPropValue<TypeTag, Properties::NumPhases>()
120 dirichletDirichlet = 1,
121 dirichletNeumann = 2,
126 innerEdgeFace = 2, innerSideFace = 1
129 using GravityVector = Dune::FieldVector<Scalar, dimWorld>;
130 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
131 using DimVector = Dune::FieldVector<Scalar, dim>;
133 InteractionVolumeContainer& interactionVolumes_()
135 return problem_.pressureModel().interactionVolumes();
140 return problem_.pressureModel().transmissibilityCalculator();
150 ParentType(problem), problem_(problem), gravity_(problem.gravity())
152 density_[wPhaseIdx] = 0.;
153 density_[nPhaseIdx] = 0.;
154 viscosity_[wPhaseIdx] = 0.;
155 viscosity_[nPhaseIdx] = 0.;
160 CellData & cellData1, CellData & cellData2, CellData & cellData3, CellData & cellData4,
161 CellData & cellData5, CellData & cellData6, CellData & cellData7, CellData & cellData8,
int fIdx = -1);
168 const auto element = *problem_.gridView().template begin<0>();
169 FluidState fluidState;
170 fluidState.setPressure(wPhaseIdx, problem_.referencePressure(element));
171 fluidState.setPressure(nPhaseIdx, problem_.referencePressure(element));
172 fluidState.setTemperature(problem_.temperature(element));
173 fluidState.setSaturation(wPhaseIdx, 1.);
174 fluidState.setSaturation(nPhaseIdx, 0.);
185 const GravityVector& gravity_;
187 Scalar density_[numPhases];
188 Scalar viscosity_[numPhases];
190 static constexpr Scalar threshold_ = 1e-15;
192 static const int velocityType_ = getPropValue<TypeTag, Properties::VelocityFormulation>();
194 static const int pressureType_ = getPropValue<TypeTag, Properties::PressureFormulation>();
196 static const int saturationType_ = getPropValue<TypeTag, Properties::SaturationFormulation>();
218template<
class TypeTag>
220 CellData & cellData1, CellData & cellData2, CellData & cellData3, CellData & cellData4,
221 CellData & cellData5, CellData & cellData6, CellData & cellData7, CellData & cellData8,
int fIdx)
223 auto element1 = interactionVolume.getSubVolumeElement(0);
224 auto element2 = interactionVolume.getSubVolumeElement(1);
225 auto element3 = interactionVolume.getSubVolumeElement(2);
226 auto element4 = interactionVolume.getSubVolumeElement(3);
227 auto element5 = interactionVolume.getSubVolumeElement(4);
228 auto element6 = interactionVolume.getSubVolumeElement(5);
229 auto element7 = interactionVolume.getSubVolumeElement(6);
230 auto element8 = interactionVolume.getSubVolumeElement(7);
233 int globalIdx1 = problem_.variables().index(element1);
234 int globalIdx2 = problem_.variables().index(element2);
235 int globalIdx3 = problem_.variables().index(element3);
236 int globalIdx4 = problem_.variables().index(element4);
237 int globalIdx5 = problem_.variables().index(element5);
238 int globalIdx6 = problem_.variables().index(element6);
239 int globalIdx7 = problem_.variables().index(element7);
240 int globalIdx8 = problem_.variables().index(element8);
243 Dune::FieldVector<Scalar, 8> potW(0);
244 Dune::FieldVector<Scalar, 8> potNw(0);
246 potW[0] = cellData1.potential(wPhaseIdx);
247 potW[1] = cellData2.potential(wPhaseIdx);
248 potW[2] = cellData3.potential(wPhaseIdx);
249 potW[3] = cellData4.potential(wPhaseIdx);
250 potW[4] = cellData5.potential(wPhaseIdx);
251 potW[5] = cellData6.potential(wPhaseIdx);
252 potW[6] = cellData7.potential(wPhaseIdx);
253 potW[7] = cellData8.potential(wPhaseIdx);
255 potNw[0] = cellData1.potential(nPhaseIdx);
256 potNw[1] = cellData2.potential(nPhaseIdx);
257 potNw[2] = cellData3.potential(nPhaseIdx);
258 potNw[3] = cellData4.potential(nPhaseIdx);
259 potNw[4] = cellData5.potential(nPhaseIdx);
260 potNw[5] = cellData6.potential(nPhaseIdx);
261 potNw[6] = cellData7.potential(nPhaseIdx);
262 potNw[7] = cellData8.potential(nPhaseIdx);
265 Dune::FieldVector<Scalar, numPhases> lambda1(cellData1.mobility(wPhaseIdx));
266 lambda1[nPhaseIdx] = cellData1.mobility(nPhaseIdx);
269 Scalar lambdaTotal1 = lambda1[wPhaseIdx] + lambda1[nPhaseIdx];
272 Dune::FieldVector<Scalar, numPhases> lambda2(cellData2.mobility(wPhaseIdx));
273 lambda2[nPhaseIdx] = cellData2.mobility(nPhaseIdx);
276 Scalar lambdaTotal2 = lambda2[wPhaseIdx] + lambda2[nPhaseIdx];
279 Dune::FieldVector<Scalar, numPhases> lambda3(cellData3.mobility(wPhaseIdx));
280 lambda3[nPhaseIdx] = cellData3.mobility(nPhaseIdx);
283 Scalar lambdaTotal3 = lambda3[wPhaseIdx] + lambda3[nPhaseIdx];
286 Dune::FieldVector<Scalar, numPhases> lambda4(cellData4.mobility(wPhaseIdx));
287 lambda4[nPhaseIdx] = cellData4.mobility(nPhaseIdx);
290 Scalar lambdaTotal4 = lambda4[wPhaseIdx] + lambda4[nPhaseIdx];
293 Dune::FieldVector<Scalar, numPhases> lambda5(cellData5.mobility(wPhaseIdx));
294 lambda5[nPhaseIdx] = cellData5.mobility(nPhaseIdx);
297 Scalar lambdaTotal5 = lambda5[wPhaseIdx] + lambda5[nPhaseIdx];
300 Dune::FieldVector<Scalar, numPhases> lambda6(cellData6.mobility(wPhaseIdx));
301 lambda6[nPhaseIdx] = cellData6.mobility(nPhaseIdx);
304 Scalar lambdaTotal6 = lambda6[wPhaseIdx] + lambda6[nPhaseIdx];
307 Dune::FieldVector<Scalar, numPhases> lambda7(cellData7.mobility(wPhaseIdx));
308 lambda7[nPhaseIdx] = cellData7.mobility(nPhaseIdx);
311 Scalar lambdaTotal7 = lambda7[wPhaseIdx] + lambda7[nPhaseIdx];
314 Dune::FieldVector<Scalar, numPhases> lambda8(cellData8.mobility(wPhaseIdx));
315 lambda8[nPhaseIdx] = cellData8.mobility(nPhaseIdx);
318 Scalar lambdaTotal8 = lambda8[wPhaseIdx] + lambda8[nPhaseIdx];
320 std::vector<DimVector> lambda(8);
321 lambda[0][0] = lambdaTotal1;
322 lambda[0][1] = lambdaTotal1;
323 lambda[0][2] = lambdaTotal1;
324 lambda[1][0] = lambdaTotal2;
325 lambda[1][1] = lambdaTotal2;
326 lambda[1][2] = lambdaTotal2;
327 lambda[2][0] = lambdaTotal3;
328 lambda[2][1] = lambdaTotal3;
329 lambda[2][2] = lambdaTotal3;
330 lambda[3][0] = lambdaTotal4;
331 lambda[3][1] = lambdaTotal4;
332 lambda[3][2] = lambdaTotal4;
333 lambda[4][0] = lambdaTotal5;
334 lambda[4][1] = lambdaTotal5;
335 lambda[4][2] = lambdaTotal5;
336 lambda[5][0] = lambdaTotal6;
337 lambda[5][1] = lambdaTotal6;
338 lambda[5][2] = lambdaTotal6;
339 lambda[6][0] = lambdaTotal7;
340 lambda[6][1] = lambdaTotal7;
341 lambda[6][2] = lambdaTotal7;
342 lambda[7][0] = lambdaTotal8;
343 lambda[7][1] = lambdaTotal8;
344 lambda[7][2] = lambdaTotal8;
346 Scalar potentialDiffW0 = 0;
347 Scalar potentialDiffW1 = 0;
348 Scalar potentialDiffW2 = 0;
349 Scalar potentialDiffW3 = 0;
350 Scalar potentialDiffW4 = 0;
351 Scalar potentialDiffW5 = 0;
352 Scalar potentialDiffW6 = 0;
353 Scalar potentialDiffW7 = 0;
354 Scalar potentialDiffW8 = 0;
355 Scalar potentialDiffW9 = 0;
356 Scalar potentialDiffW10 = 0;
357 Scalar potentialDiffW11 = 0;
359 Scalar potentialDiffNw0 = 0;
360 Scalar potentialDiffNw1 = 0;
361 Scalar potentialDiffNw2 = 0;
362 Scalar potentialDiffNw3 = 0;
363 Scalar potentialDiffNw4 = 0;
364 Scalar potentialDiffNw5 = 0;
365 Scalar potentialDiffNw6 = 0;
366 Scalar potentialDiffNw7 = 0;
367 Scalar potentialDiffNw8 = 0;
368 Scalar potentialDiffNw9 = 0;
369 Scalar potentialDiffNw10 = 0;
370 Scalar potentialDiffNw11 = 0;
373 Dune::FieldVector<Scalar, 12> fluxW(0);
374 Dune::FieldVector<Scalar, 12> fluxNw(0);
377 Dune::FieldVector<Scalar, 2 * dim - dim + 1> u(0);
378 TransmissibilityType T(0);
379 TransmissibilityType TSecond(0);
380 Dune::FieldVector<bool, 4> useCases(
false);
382 int hangingNodeType = interactionVolume.getHangingNodeType();
386 if (fIdx < 0 || fIdx == 0)
389 int caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 0, 1, 2, 3,
402 potentialDiffW0 = Tu[0];
412 potentialDiffNw0 = Tu[0];
424 potentialDiffW0 = Tu[0];
434 potentialDiffNw0 = Tu[0];
446 potentialDiffW0 = Tu[0];
456 potentialDiffNw0 = Tu[0];
468 potentialDiffW0 = Tu[0];
478 potentialDiffNw0 = Tu[0];
482 if (fIdx < 0 || fIdx == 1)
487 if (hangingNodeType == InteractionVolume::twoSmallCells
488 || hangingNodeType == InteractionVolume::fourSmallCellsDiag)
494 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 1, 3, 0,
499 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 1, 3, 0,
513 potentialDiffW1 = Tu[0];
523 potentialDiffNw1 = Tu[0];
535 potentialDiffW1 = Tu[0];
545 potentialDiffNw1 = Tu[0];
557 potentialDiffW1 = Tu[0];
567 potentialDiffNw1 = Tu[0];
579 potentialDiffW1 = Tu[0];
589 potentialDiffNw1 = Tu[0];
593 if (fIdx < 0 || fIdx == 2)
598 if (hangingNodeType != InteractionVolume::twoSmallCells
599 && hangingNodeType != InteractionVolume::fourSmallCellsDiag)
601 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 3,
614 potentialDiffW2 = Tu[0];
624 potentialDiffNw2 = Tu[0];
636 potentialDiffW2 = Tu[0];
646 potentialDiffNw2 = Tu[0];
658 potentialDiffW2 = Tu[0];
668 potentialDiffNw2 = Tu[0];
680 potentialDiffW2 = Tu[0];
690 potentialDiffNw2 = Tu[0];
695 if (fIdx < 0 || fIdx == 3)
700 if (hangingNodeType == InteractionVolume::twoSmallCells
701 || hangingNodeType == InteractionVolume::fourSmallCellsDiag)
707 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 2, 0, 3, 1, 6,
712 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 2,
726 potentialDiffW3 = Tu[0];
736 potentialDiffNw3 = Tu[0];
748 potentialDiffW3 = Tu[0];
758 potentialDiffNw3 = Tu[0];
770 potentialDiffW3 = Tu[0];
780 potentialDiffNw3 = Tu[0];
792 potentialDiffW3 = Tu[0];
802 potentialDiffNw3 = Tu[0];
806 if (fIdx < 0 || fIdx == 4)
813 if (hangingNodeType == InteractionVolume::sixSmallCells)
815 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 5, 4, 7,
818 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx6)
820 caseL = transmissibilityCalculator_().transmissibilityCaseThree(T, interactionVolume, lambda, 4, 0, 2,
823 int caseLSecond = transmissibilityCalculator_().transmissibilityCaseFour(TSecond, interactionVolume, lambda, 1, 5, 3,
826 caseL = transmissibilityCalculator_().chooseTransmissibility(T, TSecond, 3, 4);
828 if (caseL == caseLSecond)
831 if (hangingNodeType == InteractionVolume::sixSmallCells)
843 potentialDiffW4 = Tu[0];
853 potentialDiffNw4 = Tu[0];
865 potentialDiffW4 = Tu[0];
875 potentialDiffNw4 = Tu[0];
887 potentialDiffW4 = Tu[0];
897 potentialDiffNw4 = Tu[0];
909 potentialDiffW4 = Tu[0];
919 potentialDiffNw4 = Tu[0];
922 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx6)
934 potentialDiffW4 = -Tu[2];
944 potentialDiffNw4 = -Tu[2];
956 potentialDiffW4 = Tu[2];
966 potentialDiffNw4 = Tu[2];
971 if (fIdx < 0 || fIdx == 5)
977 if (hangingNodeType == InteractionVolume::sixSmallCells)
984 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 7, 5, 6, 4, 3,
987 else if (hangingNodeType == InteractionVolume::fourSmallCellsDiag)
993 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 7, 5, 6, 4, 3,
996 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx7)
998 caseL = transmissibilityCalculator_().transmissibilityCaseThree(T, interactionVolume, lambda, 1, 5, 7, 0);
1000 int caseLSecond = transmissibilityCalculator_().transmissibilityCaseFour(TSecond, interactionVolume, lambda, 7, 3, 5,
1003 caseL = transmissibilityCalculator_().chooseTransmissibility(T, TSecond, 3, 4);
1005 if (caseL == caseLSecond)
1009 if (hangingNodeType == InteractionVolume::sixSmallCells
1010 || hangingNodeType == InteractionVolume::fourSmallCellsDiag)
1022 potentialDiffW5 = Tu[0];
1032 potentialDiffNw5 = Tu[0];
1034 else if (caseL == 2)
1044 potentialDiffW5 = Tu[0];
1054 potentialDiffNw5 = Tu[0];
1056 else if (caseL == 3)
1066 potentialDiffW5 = Tu[0];
1076 potentialDiffNw5 = Tu[0];
1088 potentialDiffW5 = Tu[0];
1098 potentialDiffNw5 = Tu[0];
1101 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx7)
1113 potentialDiffW5 = -Tu[1];
1123 potentialDiffNw5 = -Tu[1];
1135 potentialDiffW5 = Tu[1];
1145 potentialDiffNw5 = Tu[1];
1150 if (fIdx < 0 || fIdx == 6)
1156 if (hangingNodeType == InteractionVolume::fourSmallCellsDiag)
1158 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 6, 7, 4,
1161 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx6)
1163 caseL = transmissibilityCalculator_().transmissibilityCaseThree(T, interactionVolume, lambda, 7, 3, 1, 6);
1165 int caseLSecond = transmissibilityCalculator_().transmissibilityCaseFour(TSecond, interactionVolume, lambda, 2, 6, 0,
1168 caseL = transmissibilityCalculator_().chooseTransmissibility(T, TSecond, 3, 4);
1170 if (caseL == caseLSecond)
1174 if (hangingNodeType == InteractionVolume::fourSmallCellsDiag)
1186 potentialDiffW6 = Tu[0];
1196 potentialDiffNw6 = Tu[0];
1198 else if (caseL == 2)
1208 potentialDiffW6 = Tu[0];
1218 potentialDiffNw6 = Tu[0];
1220 else if (caseL == 3)
1230 potentialDiffW6 = Tu[0];
1240 potentialDiffNw6 = Tu[0];
1252 potentialDiffW6 = Tu[0];
1262 potentialDiffNw6 = Tu[0];
1265 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx6)
1277 potentialDiffW6 = -Tu[2];
1287 potentialDiffNw6 = -Tu[2];
1289 else if (caseL == 4)
1299 potentialDiffW6 = Tu[2];
1309 potentialDiffNw6 = Tu[2];
1314 if (fIdx < 0 || fIdx == 7)
1320 if (hangingNodeType == InteractionVolume::sixSmallCells)
1323 useCases[1] =
false;
1324 useCases[2] =
false;
1327 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 4, 6, 5, 7, 0,
1330 else if (hangingNodeType == InteractionVolume::fourSmallCellsDiag)
1332 useCases[0] =
false;
1335 useCases[3] =
false;
1336 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 4, 6, 5, 7, 0,
1339 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx7)
1341 caseL = transmissibilityCalculator_().transmissibilityCaseThree(T, interactionVolume, lambda, 2, 6, 4,
1344 int caseLSecond = transmissibilityCalculator_().transmissibilityCaseFour(TSecond, interactionVolume, lambda, 4, 0, 6,
1349 caseL = transmissibilityCalculator_().chooseTransmissibility(T, TSecond, 3, 4);
1351 if (caseL == caseLSecond)
1355 if (hangingNodeType == InteractionVolume::sixSmallCells
1356 || hangingNodeType == InteractionVolume::fourSmallCellsDiag)
1369 potentialDiffW7 = Tu[0];
1379 potentialDiffNw7 = Tu[0];
1381 else if (caseL == 2)
1391 potentialDiffW7 = Tu[0];
1401 potentialDiffNw7 = Tu[0];
1403 else if (caseL == 3)
1413 potentialDiffW7 = Tu[0];
1423 potentialDiffNw7 = Tu[0];
1435 potentialDiffW7 = Tu[0];
1445 potentialDiffNw7 = Tu[0];
1448 else if(hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx7)
1460 potentialDiffW7 = -Tu[1];
1470 potentialDiffNw7 = -Tu[1];
1472 else if (caseL == 4)
1482 potentialDiffW7 = Tu[1];
1492 potentialDiffNw7 = Tu[1];
1497 if (fIdx < 0 || fIdx == 8)
1502 if (hangingNodeType == InteractionVolume::sixSmallCells)
1504 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 4, 0, 6,
1507 else if (hangingNodeType == InteractionVolume::twoSmallCells
1508 || hangingNodeType == InteractionVolume::fourSmallCellsFace)
1510 caseL = transmissibilityCalculator_().transmissibilityCaseTwo(T, interactionVolume, lambda, 4, 0, 2,
1513 else if (hangingNodeType == InteractionVolume::fourSmallCellsDiag
1514 || (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx7))
1516 useCases[0] =
false;
1518 useCases[2] =
false;
1521 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 4, 0, 6,
1524 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx6)
1526 useCases[0] =
false;
1529 useCases[3] =
false;
1531 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 4, 0, 6,
1545 potentialDiffW8 = Tu[0];
1555 potentialDiffNw8 = Tu[0];
1557 else if (caseL == 2)
1567 potentialDiffW8 = Tu[0];
1577 potentialDiffNw8 = Tu[0];
1579 else if (caseL == 3)
1589 potentialDiffW8 = Tu[0];
1599 potentialDiffNw8 = Tu[0];
1611 potentialDiffW8 = Tu[0];
1621 potentialDiffNw8 = Tu[0];
1625 if (fIdx < 0 || fIdx == 9)
1630 if (hangingNodeType == InteractionVolume::sixSmallCells)
1632 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 1, 5, 3,
1635 else if (hangingNodeType == InteractionVolume::twoSmallCells
1636 || hangingNodeType == InteractionVolume::fourSmallCellsFace)
1638 caseL = transmissibilityCalculator_().transmissibilityCaseOne(T, interactionVolume, lambda, 1, 5, 3,
1641 else if (hangingNodeType == InteractionVolume::fourSmallCellsDiag
1642 || (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx7))
1645 useCases[1] =
false;
1647 useCases[3] =
false;
1649 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 1, 5, 3,
1652 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx6)
1655 useCases[1] =
false;
1656 useCases[2] =
false;
1659 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 1, 5, 3,
1673 potentialDiffW9 = Tu[0];
1683 potentialDiffNw9 = Tu[0];
1685 else if (caseL == 2)
1695 potentialDiffW9 = Tu[0];
1705 potentialDiffNw9 = Tu[0];
1707 else if (caseL == 3)
1717 potentialDiffW9 = Tu[0];
1727 potentialDiffNw9 = Tu[0];
1739 potentialDiffW9 = Tu[0];
1749 potentialDiffNw9 = Tu[0];
1753 if (fIdx < 0 || fIdx == 10)
1758 if (hangingNodeType == InteractionVolume::fourSmallCellsFace)
1760 caseL = transmissibilityCalculator_().transmissibilityCaseTwo(T, interactionVolume, lambda, 7, 3, 1, 2);
1762 else if ((hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx7)
1763 || hangingNodeType == InteractionVolume::sixSmallCells)
1765 useCases[0] =
false;
1767 useCases[2] =
false;
1770 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 7, 3, 5, 1, 6,
1773 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx6)
1775 useCases[0] =
false;
1778 useCases[3] =
false;
1780 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 7, 3, 5, 1, 6,
1783 else if (hangingNodeType == InteractionVolume::fourSmallCellsDiag)
1786 useCases[1] =
false;
1788 useCases[3] =
false;
1790 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 7, 3, 5, 1, 6,
1794 if (hangingNodeType != InteractionVolume::twoSmallCells)
1807 potentialDiffW10 = Tu[0];
1817 potentialDiffNw10 = Tu[0];
1819 else if (caseL == 2)
1829 potentialDiffW10 = Tu[0];
1839 potentialDiffNw10 = Tu[0];
1841 else if (caseL == 3)
1851 potentialDiffW10 = Tu[0];
1861 potentialDiffNw10 = Tu[0];
1873 potentialDiffW10 = Tu[0];
1883 potentialDiffNw10 = Tu[0];
1888 if (fIdx < 0 || fIdx == 11)
1893 if (hangingNodeType == InteractionVolume::fourSmallCellsFace)
1895 caseL = transmissibilityCalculator_().transmissibilityCaseOne(T, interactionVolume, lambda, 2, 6, 0, 3);
1897 else if ((hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx7)
1898 || hangingNodeType == InteractionVolume::sixSmallCells)
1901 useCases[1] =
false;
1903 useCases[3] =
false;
1905 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 2, 6, 0, 4, 3,
1908 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx6)
1911 useCases[1] =
false;
1912 useCases[2] =
false;
1915 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 2, 6, 0, 4, 3,
1918 else if (hangingNodeType == InteractionVolume::fourSmallCellsDiag)
1920 useCases[0] =
false;
1922 useCases[2] =
false;
1925 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 2, 6, 0, 4, 3,
1929 if (hangingNodeType != InteractionVolume::twoSmallCells)
1942 potentialDiffW11 = Tu[0];
1952 potentialDiffNw11 = Tu[0];
1954 else if (caseL == 2)
1964 potentialDiffW11 = Tu[0];
1974 potentialDiffNw11 = Tu[0];
1976 else if (caseL == 3)
1986 potentialDiffW11 = Tu[0];
1996 potentialDiffNw11 = Tu[0];
2008 potentialDiffW11 = Tu[0];
2018 potentialDiffNw11 = Tu[0];
2023 cellData1.fluxData().addUpwindPotential(wPhaseIdx, interactionVolume.getIndexOnElement(0, 0), potentialDiffW0);
2024 cellData1.fluxData().addUpwindPotential(nPhaseIdx, interactionVolume.getIndexOnElement(0, 0), potentialDiffNw0);
2025 cellData1.fluxData().addUpwindPotential(wPhaseIdx, interactionVolume.getIndexOnElement(0, 1), -potentialDiffW3);
2026 cellData1.fluxData().addUpwindPotential(nPhaseIdx, interactionVolume.getIndexOnElement(0, 1), -potentialDiffNw3);
2027 cellData1.fluxData().addUpwindPotential(wPhaseIdx, interactionVolume.getIndexOnElement(0, 2), -potentialDiffW8);
2028 cellData1.fluxData().addUpwindPotential(nPhaseIdx, interactionVolume.getIndexOnElement(0, 2), -potentialDiffNw8);
2030 cellData2.fluxData().addUpwindPotential(wPhaseIdx, interactionVolume.getIndexOnElement(1, 0), potentialDiffW1);
2031 cellData2.fluxData().addUpwindPotential(nPhaseIdx, interactionVolume.getIndexOnElement(1, 0), potentialDiffNw1);
2032 cellData2.fluxData().addUpwindPotential(wPhaseIdx, interactionVolume.getIndexOnElement(1, 1), -potentialDiffW0);
2033 cellData2.fluxData().addUpwindPotential(nPhaseIdx, interactionVolume.getIndexOnElement(1, 1), -potentialDiffNw0);
2034 cellData2.fluxData().addUpwindPotential(wPhaseIdx, interactionVolume.getIndexOnElement(1, 2), potentialDiffW9);
2035 cellData2.fluxData().addUpwindPotential(nPhaseIdx, interactionVolume.getIndexOnElement(1, 2), potentialDiffNw9);
2037 cellData3.fluxData().addUpwindPotential(wPhaseIdx, interactionVolume.getIndexOnElement(2, 0), potentialDiffW3);
2038 cellData3.fluxData().addUpwindPotential(nPhaseIdx, interactionVolume.getIndexOnElement(2, 0), potentialDiffNw3);
2039 cellData3.fluxData().addUpwindPotential(wPhaseIdx, interactionVolume.getIndexOnElement(2, 1), -potentialDiffW2);
2040 cellData3.fluxData().addUpwindPotential(nPhaseIdx, interactionVolume.getIndexOnElement(2, 1), -potentialDiffNw2);
2041 cellData3.fluxData().addUpwindPotential(wPhaseIdx, interactionVolume.getIndexOnElement(2, 2), potentialDiffW11);
2042 cellData3.fluxData().addUpwindPotential(nPhaseIdx, interactionVolume.getIndexOnElement(2, 2), potentialDiffNw11);
2044 cellData4.fluxData().addUpwindPotential(wPhaseIdx, interactionVolume.getIndexOnElement(3, 0), potentialDiffW2);
2045 cellData4.fluxData().addUpwindPotential(nPhaseIdx, interactionVolume.getIndexOnElement(3, 0), potentialDiffNw2);
2046 cellData4.fluxData().addUpwindPotential(wPhaseIdx, interactionVolume.getIndexOnElement(3, 1), -potentialDiffW1);
2047 cellData4.fluxData().addUpwindPotential(nPhaseIdx, interactionVolume.getIndexOnElement(3, 1), -potentialDiffNw1);
2048 cellData4.fluxData().addUpwindPotential(wPhaseIdx, interactionVolume.getIndexOnElement(3, 2), -potentialDiffW10);
2049 cellData4.fluxData().addUpwindPotential(nPhaseIdx, interactionVolume.getIndexOnElement(3, 2), -potentialDiffNw10);
2051 cellData5.fluxData().addUpwindPotential(wPhaseIdx, interactionVolume.getIndexOnElement(4, 0), potentialDiffW8);
2052 cellData5.fluxData().addUpwindPotential(nPhaseIdx, interactionVolume.getIndexOnElement(4, 0), potentialDiffNw8);
2053 cellData5.fluxData().addUpwindPotential(wPhaseIdx, interactionVolume.getIndexOnElement(4, 1), -potentialDiffW4);
2054 cellData5.fluxData().addUpwindPotential(nPhaseIdx, interactionVolume.getIndexOnElement(4, 1), -potentialDiffNw4);
2055 cellData5.fluxData().addUpwindPotential(wPhaseIdx, interactionVolume.getIndexOnElement(4, 2), potentialDiffW7);
2056 cellData5.fluxData().addUpwindPotential(nPhaseIdx, interactionVolume.getIndexOnElement(4, 2), potentialDiffNw7);
2058 cellData6.fluxData().addUpwindPotential(wPhaseIdx, interactionVolume.getIndexOnElement(5, 0), -potentialDiffW9);
2059 cellData6.fluxData().addUpwindPotential(nPhaseIdx, interactionVolume.getIndexOnElement(5, 0), -potentialDiffNw9);
2060 cellData6.fluxData().addUpwindPotential(wPhaseIdx, interactionVolume.getIndexOnElement(5, 1), -potentialDiffW5);
2061 cellData6.fluxData().addUpwindPotential(nPhaseIdx, interactionVolume.getIndexOnElement(5, 1), -potentialDiffNw5);
2062 cellData6.fluxData().addUpwindPotential(wPhaseIdx, interactionVolume.getIndexOnElement(5, 2), potentialDiffW4);
2063 cellData6.fluxData().addUpwindPotential(nPhaseIdx, interactionVolume.getIndexOnElement(5, 2), potentialDiffNw4);
2065 cellData7.fluxData().addUpwindPotential(wPhaseIdx, interactionVolume.getIndexOnElement(6, 0), -potentialDiffW11);
2066 cellData7.fluxData().addUpwindPotential(nPhaseIdx, interactionVolume.getIndexOnElement(6, 0), -potentialDiffNw11);
2067 cellData7.fluxData().addUpwindPotential(wPhaseIdx, interactionVolume.getIndexOnElement(6, 1), -potentialDiffW7);
2068 cellData7.fluxData().addUpwindPotential(nPhaseIdx, interactionVolume.getIndexOnElement(6, 1), -potentialDiffNw7);
2069 cellData7.fluxData().addUpwindPotential(wPhaseIdx, interactionVolume.getIndexOnElement(6, 2), potentialDiffW6);
2070 cellData7.fluxData().addUpwindPotential(nPhaseIdx, interactionVolume.getIndexOnElement(6, 2), potentialDiffNw6);
2072 cellData8.fluxData().addUpwindPotential(wPhaseIdx, interactionVolume.getIndexOnElement(7, 0), potentialDiffW10);
2073 cellData8.fluxData().addUpwindPotential(nPhaseIdx, interactionVolume.getIndexOnElement(7, 0), potentialDiffNw10);
2074 cellData8.fluxData().addUpwindPotential(wPhaseIdx, interactionVolume.getIndexOnElement(7, 1), -potentialDiffW6);
2075 cellData8.fluxData().addUpwindPotential(nPhaseIdx, interactionVolume.getIndexOnElement(7, 1), -potentialDiffNw6);
2076 cellData8.fluxData().addUpwindPotential(wPhaseIdx, interactionVolume.getIndexOnElement(7, 2), potentialDiffW5);
2077 cellData8.fluxData().addUpwindPotential(nPhaseIdx, interactionVolume.getIndexOnElement(7, 2), potentialDiffNw5);
2080 Dune::FieldVector<Scalar, numPhases> lambda0Upw(0.0);
2081 lambda0Upw[wPhaseIdx] = (potentialDiffW0 >= 0) ? lambda1[wPhaseIdx] : lambda2[wPhaseIdx];
2082 lambda0Upw[nPhaseIdx] = (potentialDiffNw0 >= 0) ? lambda1[nPhaseIdx] : lambda2[nPhaseIdx];
2085 Dune::FieldVector<Scalar, numPhases> lambda1Upw(0.0);
2086 lambda1Upw[wPhaseIdx] = (potentialDiffW1 >= 0) ? lambda2[wPhaseIdx] : lambda4[wPhaseIdx];
2087 lambda1Upw[nPhaseIdx] = (potentialDiffNw1 >= 0) ? lambda2[nPhaseIdx] : lambda4[nPhaseIdx];
2090 Dune::FieldVector<Scalar, numPhases> lambda2Upw(0.0);
2091 lambda2Upw[wPhaseIdx] = (potentialDiffW2 >= 0) ? lambda4[wPhaseIdx] : lambda3[wPhaseIdx];
2092 lambda2Upw[nPhaseIdx] = (potentialDiffNw2 >= 0) ? lambda4[nPhaseIdx] : lambda3[nPhaseIdx];
2095 Dune::FieldVector<Scalar, numPhases> lambda3Upw(0.0);
2096 lambda3Upw[wPhaseIdx] = (potentialDiffW3 >= 0) ? lambda3[wPhaseIdx] : lambda1[wPhaseIdx];
2097 lambda3Upw[nPhaseIdx] = (potentialDiffNw3 >= 0) ? lambda3[nPhaseIdx] : lambda1[nPhaseIdx];
2100 Dune::FieldVector<Scalar, numPhases> lambda4Upw(0.0);
2101 lambda4Upw[wPhaseIdx] = (potentialDiffW4 >= 0) ? lambda6[wPhaseIdx] : lambda5[wPhaseIdx];
2102 lambda4Upw[nPhaseIdx] = (potentialDiffNw4 >= 0) ? lambda6[nPhaseIdx] : lambda5[nPhaseIdx];
2105 Dune::FieldVector<Scalar, numPhases> lambda5Upw(0.0);
2106 lambda5Upw[wPhaseIdx] = (potentialDiffW5 >= 0) ? lambda8[wPhaseIdx] : lambda6[wPhaseIdx];
2107 lambda5Upw[nPhaseIdx] = (potentialDiffNw5 >= 0) ? lambda8[nPhaseIdx] : lambda6[nPhaseIdx];
2110 Dune::FieldVector<Scalar, numPhases> lambda6Upw(0.0);
2111 lambda6Upw[wPhaseIdx] = (potentialDiffW6 >= 0) ? lambda7[wPhaseIdx] : lambda8[wPhaseIdx];
2112 lambda6Upw[nPhaseIdx] = (potentialDiffNw6 >= 0) ? lambda7[nPhaseIdx] : lambda8[nPhaseIdx];
2115 Dune::FieldVector<Scalar, numPhases> lambda7Upw(0.0);
2116 lambda7Upw[wPhaseIdx] = (potentialDiffW7 >= 0) ? lambda5[wPhaseIdx] : lambda7[wPhaseIdx];
2117 lambda7Upw[nPhaseIdx] = (potentialDiffNw7 >= 0) ? lambda5[nPhaseIdx] : lambda7[nPhaseIdx];
2120 Dune::FieldVector<Scalar, numPhases> lambda8Upw(0.0);
2121 lambda8Upw[wPhaseIdx] = (potentialDiffW8 >= 0) ? lambda5[wPhaseIdx] : lambda1[wPhaseIdx];
2122 lambda8Upw[nPhaseIdx] = (potentialDiffNw8 >= 0) ? lambda5[nPhaseIdx] : lambda1[nPhaseIdx];
2125 Dune::FieldVector<Scalar, numPhases> lambda9Upw(0.0);
2126 lambda9Upw[wPhaseIdx] = (potentialDiffW9 >= 0) ? lambda2[wPhaseIdx] : lambda6[wPhaseIdx];
2127 lambda9Upw[nPhaseIdx] = (potentialDiffNw9 >= 0) ? lambda2[nPhaseIdx] : lambda6[nPhaseIdx];
2130 Dune::FieldVector<Scalar, numPhases> lambda10Upw(0.0);
2131 lambda10Upw[wPhaseIdx] = (potentialDiffW10 >= 0) ? lambda8[wPhaseIdx] : lambda4[wPhaseIdx];
2132 lambda10Upw[nPhaseIdx] = (potentialDiffNw10 >= 0) ? lambda8[nPhaseIdx] : lambda4[nPhaseIdx];
2135 Dune::FieldVector<Scalar, numPhases> lambda11Upw(0.0);
2136 lambda11Upw[wPhaseIdx] = (potentialDiffW11 >= 0) ? lambda3[wPhaseIdx] : lambda7[wPhaseIdx];
2137 lambda11Upw[nPhaseIdx] = (potentialDiffNw11 >= 0) ? lambda3[nPhaseIdx] : lambda7[nPhaseIdx];
2139 for (
int i = 0; i < numPhases; i++)
2142 DimVector vel12 = interactionVolume.getNormal(0, 0);
2143 DimVector vel21 = interactionVolume.getNormal(1, 1);
2144 DimVector vel24 = interactionVolume.getNormal(1, 0);
2145 DimVector vel42 = interactionVolume.getNormal(3, 1);
2146 DimVector vel43 = interactionVolume.getNormal(3, 0);
2147 DimVector vel34 = interactionVolume.getNormal(2, 1);
2148 DimVector vel31 = interactionVolume.getNormal(2, 0);
2149 DimVector vel13 = interactionVolume.getNormal(0, 1);
2150 DimVector vel65 = interactionVolume.getNormal(5, 2);
2151 DimVector vel56 = interactionVolume.getNormal(4, 1);
2152 DimVector vel86 = interactionVolume.getNormal(7, 2);
2153 DimVector vel68 = interactionVolume.getNormal(5, 1);
2154 DimVector vel78 = interactionVolume.getNormal(6, 2);
2155 DimVector vel87 = interactionVolume.getNormal(7, 1);
2156 DimVector vel57 = interactionVolume.getNormal(4, 2);
2157 DimVector vel75 = interactionVolume.getNormal(6, 1);
2158 DimVector vel51 = interactionVolume.getNormal(4, 0);
2159 DimVector vel15 = interactionVolume.getNormal(0, 2);
2160 DimVector vel26 = interactionVolume.getNormal(1, 2);
2161 DimVector vel62 = interactionVolume.getNormal(5, 0);
2162 DimVector vel84 = interactionVolume.getNormal(7, 0);
2163 DimVector vel48 = interactionVolume.getNormal(3, 2);
2164 DimVector vel37 = interactionVolume.getNormal(2, 2);
2165 DimVector vel73 = interactionVolume.getNormal(6, 0);
2179 Dune::FieldVector<Scalar, 12> flux(0);
2194 if (hangingNodeType == InteractionVolume::sixSmallCells)
2196 vel12 *= flux[0] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx1, 0, 0));
2197 vel21 *= flux[0] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx2, 1, 1));
2198 vel24 *= flux[1] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx2, 1, 0));
2199 vel42 *= flux[1] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx4, 3, 1));
2200 vel43 *= flux[2] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx4, 3, 0));
2201 vel34 *= flux[2] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx3, 2, 1));
2202 vel31 *= flux[3] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx3, 2, 0));
2203 vel13 *= flux[3] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx1, 0, 1));
2204 vel65 *= flux[4] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx6, 5, 2));
2205 vel56 *= flux[4] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx5, 4, 1));
2206 vel86 *= flux[5] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx8, 7, 2));
2207 vel68 *= flux[5] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx6, 5, 1));
2210 vel57 *= flux[7] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx5, 4, 2));
2211 vel75 *= flux[7] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx7, 6, 1));
2212 vel51 *= flux[8] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx5, 4, 0));
2213 vel15 *= flux[8] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx1, 0, 2));
2214 vel26 *= flux[9] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx2, 1, 2));
2215 vel62 *= flux[9] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx6, 5, 0));
2216 vel84 *= flux[10] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx8, 7, 0));
2217 vel48 *= flux[10] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx4, 3, 2));
2218 vel37 *= flux[11] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx3, 2, 2));
2219 vel73 *= flux[11] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx7, 6, 0));
2221 else if (hangingNodeType == InteractionVolume::fourSmallCellsDiag)
2223 vel12 *= flux[0] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx1, 0, 0));
2224 vel21 *= flux[0] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx2, 1, 1));
2225 vel24 *= flux[1] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx2, 1, 0));
2226 vel42 *= flux[1] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx4, 3, 1));
2229 vel31 *= flux[3] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx3, 2, 0));
2230 vel13 *= flux[3] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx1, 0, 1));
2233 vel86 *= flux[5] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx8, 7, 2));
2234 vel68 *= flux[5] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx6, 5, 1));
2235 vel78 *= flux[6] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx7, 6, 2));
2236 vel87 *= flux[6] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx8, 7, 1));
2237 vel57 *= flux[7] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx5, 4, 2));
2238 vel75 *= flux[7] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx7, 6, 1));
2239 vel51 *= flux[8] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx5, 4, 0));
2240 vel15 *= flux[8] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx1, 0, 2));
2241 vel26 *= flux[9] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx2, 1, 2));
2242 vel62 *= flux[9] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx6, 5, 0));
2243 vel84 *= flux[10] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx8, 7, 0));
2244 vel48 *= flux[10] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx4, 3, 2));
2245 vel37 *= flux[11] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx3, 2, 2));
2246 vel73 *= flux[11] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx7, 6, 0));
2248 else if (hangingNodeType == InteractionVolume::fourSmallCellsFace
2249 || hangingNodeType == InteractionVolume::fourSmallCellsEdge)
2251 vel12 *= flux[0] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx1, 0, 0));
2252 vel21 *= flux[0] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx2, 1, 1));
2253 vel24 *= flux[1] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx2, 1, 0));
2254 vel42 *= flux[1] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx4, 3, 1));
2255 vel43 *= flux[2] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx4, 3, 0));
2256 vel34 *= flux[2] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx3, 2, 1));
2257 vel31 *= flux[3] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx3, 2, 0));
2258 vel13 *= flux[3] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx1, 0, 1));
2259 vel51 *= flux[8] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx5, 4, 0));
2260 vel15 *= flux[8] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx1, 0, 2));
2261 vel26 *= flux[9] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx2, 1, 2));
2262 vel62 *= flux[9] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx6, 5, 0));
2263 vel84 *= flux[10] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx8, 7, 0));
2264 vel48 *= flux[10] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx4, 3, 2));
2265 vel37 *= flux[11] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx3, 2, 2));
2266 vel73 *= flux[11] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx7, 6, 0));
2268 if (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx7)
2272 vel86 *= flux[5] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx8, 7, 2));
2273 vel68 *= flux[5] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx6, 5, 1));
2276 vel57 *= flux[7] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx5, 4, 2));
2277 vel75 *= flux[7] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx7, 6, 1));
2279 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx6)
2281 vel65 *= flux[4] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx6, 5, 2));
2282 vel56 *= flux[4] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx5, 4, 1));
2285 vel78 *= flux[6] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx7, 6, 2));
2286 vel87 *= flux[6] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx8, 7, 1));
2302 else if (hangingNodeType == InteractionVolume::twoSmallCells)
2304 vel12 *= flux[0] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx1, 0, 0));
2305 vel21 *= flux[0] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx2, 1, 1));
2306 vel24 *= flux[1] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx2, 1, 0));
2307 vel42 *= flux[1] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx4, 3, 1));
2310 vel31 *= flux[3] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx3, 2, 0));
2311 vel13 *= flux[3] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx1, 0, 1));
2320 vel51 *= flux[8] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx5, 4, 0));
2321 vel15 *= flux[8] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx1, 0, 2));
2322 vel26 *= flux[9] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx2, 1, 2));
2323 vel62 *= flux[9] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx6, 5, 0));
2330 Scalar lambdaT0 = lambda0Upw[wPhaseIdx] + lambda0Upw[nPhaseIdx];
2331 Scalar lambdaT1 = lambda1Upw[wPhaseIdx] + lambda1Upw[nPhaseIdx];
2332 Scalar lambdaT2 = lambda2Upw[wPhaseIdx] + lambda2Upw[nPhaseIdx];
2333 Scalar lambdaT3 = lambda3Upw[wPhaseIdx] + lambda3Upw[nPhaseIdx];
2334 Scalar lambdaT4 = lambda4Upw[wPhaseIdx] + lambda4Upw[nPhaseIdx];
2335 Scalar lambdaT5 = lambda5Upw[wPhaseIdx] + lambda5Upw[nPhaseIdx];
2336 Scalar lambdaT6 = lambda6Upw[wPhaseIdx] + lambda6Upw[nPhaseIdx];
2337 Scalar lambdaT7 = lambda7Upw[wPhaseIdx] + lambda7Upw[nPhaseIdx];
2338 Scalar lambdaT8 = lambda8Upw[wPhaseIdx] + lambda8Upw[nPhaseIdx];
2339 Scalar lambdaT9 = lambda9Upw[wPhaseIdx] + lambda9Upw[nPhaseIdx];
2340 Scalar lambdaT10 = lambda10Upw[wPhaseIdx] + lambda10Upw[nPhaseIdx];
2341 Scalar lambdaT11 = lambda11Upw[wPhaseIdx] + lambda11Upw[nPhaseIdx];
2342 Scalar fracFlow0 = (lambdaT0 > threshold_) ? lambda0Upw[i] / (lambdaT0) : 0.0;
2343 Scalar fracFlow1 = (lambdaT1 > threshold_) ? lambda1Upw[i] / (lambdaT1) : 0.0;
2344 Scalar fracFlow2 = (lambdaT2 > threshold_) ? lambda2Upw[i] / (lambdaT2) : 0.0;
2345 Scalar fracFlow3 = (lambdaT3 > threshold_) ? lambda3Upw[i] / (lambdaT3) : 0.0;
2346 Scalar fracFlow4 = (lambdaT4 > threshold_) ? lambda4Upw[i] / (lambdaT4) : 0.0;
2347 Scalar fracFlow5 = (lambdaT5 > threshold_) ? lambda5Upw[i] / (lambdaT5) : 0.0;
2348 Scalar fracFlow6 = (lambdaT6 > threshold_) ? lambda6Upw[i] / (lambdaT6) : 0.0;
2349 Scalar fracFlow7 = (lambdaT7 > threshold_) ? lambda7Upw[i] / (lambdaT7) : 0.0;
2350 Scalar fracFlow8 = (lambdaT8 > threshold_) ? lambda8Upw[i] / (lambdaT8) : 0.0;
2351 Scalar fracFlow9 = (lambdaT9 > threshold_) ? lambda9Upw[i] / (lambdaT9) : 0.0;
2352 Scalar fracFlow10 = (lambdaT10 > threshold_) ? lambda10Upw[i] / (lambdaT10) : 0.0;
2353 Scalar fracFlow11 = (lambdaT11 > threshold_) ? lambda11Upw[i] / (lambdaT11) : 0.0;
2375 vel84 *= fracFlow10;
2376 vel48 *= fracFlow10;
2377 vel37 *= fracFlow11;
2378 vel73 *= fracFlow11;
2381 cellData1.fluxData().addVelocity(i, interactionVolume.getIndexOnElement(0, 0), vel12);
2382 cellData1.fluxData().addVelocity(i, interactionVolume.getIndexOnElement(0, 1), vel13);
2383 cellData1.fluxData().addVelocity(i, interactionVolume.getIndexOnElement(0, 2), vel15);
2384 cellData2.fluxData().addVelocity(i, interactionVolume.getIndexOnElement(1, 0), vel24);
2385 cellData2.fluxData().addVelocity(i, interactionVolume.getIndexOnElement(1, 1), vel21);
2386 cellData2.fluxData().addVelocity(i, interactionVolume.getIndexOnElement(1, 2), vel26);
2387 cellData3.fluxData().addVelocity(i, interactionVolume.getIndexOnElement(2, 0), vel31);
2388 cellData3.fluxData().addVelocity(i, interactionVolume.getIndexOnElement(2, 1), vel34);
2389 cellData3.fluxData().addVelocity(i, interactionVolume.getIndexOnElement(2, 2), vel37);
2390 cellData4.fluxData().addVelocity(i, interactionVolume.getIndexOnElement(3, 0), vel43);
2391 cellData4.fluxData().addVelocity(i, interactionVolume.getIndexOnElement(3, 1), vel42);
2392 cellData4.fluxData().addVelocity(i, interactionVolume.getIndexOnElement(3, 2), vel48);
2393 cellData5.fluxData().addVelocity(i, interactionVolume.getIndexOnElement(4, 0), vel51);
2394 cellData5.fluxData().addVelocity(i, interactionVolume.getIndexOnElement(4, 1), vel56);
2395 cellData5.fluxData().addVelocity(i, interactionVolume.getIndexOnElement(4, 2), vel57);
2396 cellData6.fluxData().addVelocity(i, interactionVolume.getIndexOnElement(5, 0), vel62);
2397 cellData6.fluxData().addVelocity(i, interactionVolume.getIndexOnElement(5, 1), vel68);
2398 cellData6.fluxData().addVelocity(i, interactionVolume.getIndexOnElement(5, 2), vel65);
2399 cellData7.fluxData().addVelocity(i, interactionVolume.getIndexOnElement(6, 0), vel73);
2400 cellData7.fluxData().addVelocity(i, interactionVolume.getIndexOnElement(6, 1), vel75);
2401 cellData7.fluxData().addVelocity(i, interactionVolume.getIndexOnElement(6, 2), vel78);
2402 cellData8.fluxData().addVelocity(i, interactionVolume.getIndexOnElement(7, 0), vel84);
2403 cellData8.fluxData().addVelocity(i, interactionVolume.getIndexOnElement(7, 1), vel87);
2404 cellData8.fluxData().addVelocity(i, interactionVolume.getIndexOnElement(7, 2), vel86);
2407 cellData1.fluxData().setVelocityMarker(interactionVolume.getIndexOnElement(0, 0));
2408 cellData1.fluxData().setVelocityMarker(interactionVolume.getIndexOnElement(0, 1));
2409 cellData1.fluxData().setVelocityMarker(interactionVolume.getIndexOnElement(0, 2));
2410 cellData2.fluxData().setVelocityMarker(interactionVolume.getIndexOnElement(1, 0));
2411 cellData2.fluxData().setVelocityMarker(interactionVolume.getIndexOnElement(1, 1));
2412 cellData2.fluxData().setVelocityMarker(interactionVolume.getIndexOnElement(1, 2));
2413 cellData3.fluxData().setVelocityMarker(interactionVolume.getIndexOnElement(2, 0));
2414 cellData3.fluxData().setVelocityMarker(interactionVolume.getIndexOnElement(2, 1));
2415 cellData3.fluxData().setVelocityMarker(interactionVolume.getIndexOnElement(2, 2));
2416 cellData4.fluxData().setVelocityMarker(interactionVolume.getIndexOnElement(3, 0));
2417 cellData4.fluxData().setVelocityMarker(interactionVolume.getIndexOnElement(3, 1));
2418 cellData4.fluxData().setVelocityMarker(interactionVolume.getIndexOnElement(3, 2));
2419 cellData5.fluxData().setVelocityMarker(interactionVolume.getIndexOnElement(4, 0));
2420 cellData5.fluxData().setVelocityMarker(interactionVolume.getIndexOnElement(4, 1));
2421 cellData5.fluxData().setVelocityMarker(interactionVolume.getIndexOnElement(4, 2));
2422 cellData6.fluxData().setVelocityMarker(interactionVolume.getIndexOnElement(5, 0));
2423 cellData6.fluxData().setVelocityMarker(interactionVolume.getIndexOnElement(5, 1));
2424 cellData6.fluxData().setVelocityMarker(interactionVolume.getIndexOnElement(5, 2));
2425 cellData7.fluxData().setVelocityMarker(interactionVolume.getIndexOnElement(6, 0));
2426 cellData7.fluxData().setVelocityMarker(interactionVolume.getIndexOnElement(6, 1));
2427 cellData7.fluxData().setVelocityMarker(interactionVolume.getIndexOnElement(6, 2));
2428 cellData8.fluxData().setVelocityMarker(interactionVolume.getIndexOnElement(7, 0));
2429 cellData8.fluxData().setVelocityMarker(interactionVolume.getIndexOnElement(7, 1));
2430 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:55
void initialize(bool solveTwice=true)
Initializes the velocity model.
Definition: 3dvelocity.hh:161
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:149
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:219
void initialize()
Initializes the velocity model.
Definition: 3dvelocityadaptive.hh:164