3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
3dvelocityadaptive.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*****************************************************************************
4 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_FVMPFAL3DVELOCITY2P_ADAPTIVE_HH
25#define DUMUX_FVMPFAL3DVELOCITY2P_ADAPTIVE_HH
26
27#include "3dvelocity.hh"
28
29namespace Dumux {
30
53template<class TypeTag> class FvMpfaL3dVelocity2pAdaptive: public FvMpfaL3dVelocity2p<TypeTag>
54{
56
58
59 enum
60 {
61 dim = GridView::dimension, dimWorld = GridView::dimensionworld
62 };
63
66
68
70
73
75
77 using PrimaryVariables = typename SolutionTypes::PrimaryVariables;
78
80
81 using Element = typename GridView::Traits::template Codim<0>::Entity;
82 using Grid = typename GridView::Grid;
83 using IndexSet = typename GridView::IndexSet;
84
86
90 using TransmissibilityType = typename TransmissibilityCalculator::TransmissibilityType;
91
92
93 enum
94 {
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
103 };
104 enum
105 {
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>()
113 };
114 enum
115 {
116 globalCorner = 2,
117 globalEdge = 3,
118 neumannNeumann = 0,
119 dirichletDirichlet = 1,
120 dirichletNeumann = 2,
121 neumannDirichlet = 3
122 };
123 enum
124 {
125 innerEdgeFace = 2, innerSideFace = 1
126 };
127
128 using GravityVector = Dune::FieldVector<Scalar, dimWorld>;
129 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
130 using DimVector = Dune::FieldVector<Scalar, dim>;
131
132 InteractionVolumeContainer& interactionVolumes_()
133 {
134 return problem_.pressureModel().interactionVolumes();
135 }
136
137 TransmissibilityCalculator& transmissibilityCalculator_()
138 {
139 return problem_.pressureModel().transmissibilityCalculator();
140 }
141
142public:
148 FvMpfaL3dVelocity2pAdaptive(Problem& problem) :
149 ParentType(problem), problem_(problem), gravity_(problem.gravity())
150{
151 density_[wPhaseIdx] = 0.;
152 density_[nPhaseIdx] = 0.;
153 viscosity_[wPhaseIdx] = 0.;
154 viscosity_[nPhaseIdx] = 0.;
155}
156
158 void calculateHangingNodeInteractionVolumeVelocity(InteractionVolume& interactionVolume,
159 CellData & cellData1, CellData & cellData2, CellData & cellData3, CellData & cellData4,
160 CellData & cellData5, CellData & cellData6, CellData & cellData7, CellData & cellData8, int fIdx = -1);
161
164 {
166
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.);
174 density_[wPhaseIdx] = FluidSystem::density(fluidState, wPhaseIdx);
175 density_[nPhaseIdx] = FluidSystem::density(fluidState, nPhaseIdx);
176 viscosity_[wPhaseIdx] = FluidSystem::viscosity(fluidState, wPhaseIdx);
177 viscosity_[nPhaseIdx] = FluidSystem::viscosity(fluidState, nPhaseIdx);
178
179 return;
180 }
181
182private:
183 Problem& problem_;
184 const GravityVector& gravity_;
185
186 Scalar density_[numPhases];
187 Scalar viscosity_[numPhases];
188
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>();
196};
197// end of template
198
217template<class TypeTag>
219 CellData & cellData1, CellData & cellData2, CellData & cellData3, CellData & cellData4,
220 CellData & cellData5, CellData & cellData6, CellData & cellData7, CellData & cellData8, int fIdx)
221 {
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);
230
231 // cell index
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);
240
241 // pressures flux calculation
242 Dune::FieldVector<Scalar, 8> potW(0);
243 Dune::FieldVector<Scalar, 8> potNw(0);
244
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);
253
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);
262
263 // get mobilities of the phases
264 Dune::FieldVector<Scalar, numPhases> lambda1(cellData1.mobility(wPhaseIdx));
265 lambda1[nPhaseIdx] = cellData1.mobility(nPhaseIdx);
266
267 // compute total mobility of cell 1
268 Scalar lambdaTotal1 = lambda1[wPhaseIdx] + lambda1[nPhaseIdx];
269
270 // get mobilities of the phases
271 Dune::FieldVector<Scalar, numPhases> lambda2(cellData2.mobility(wPhaseIdx));
272 lambda2[nPhaseIdx] = cellData2.mobility(nPhaseIdx);
273
274 // compute total mobility of cell 2
275 Scalar lambdaTotal2 = lambda2[wPhaseIdx] + lambda2[nPhaseIdx];
276
277 // get mobilities of the phases
278 Dune::FieldVector<Scalar, numPhases> lambda3(cellData3.mobility(wPhaseIdx));
279 lambda3[nPhaseIdx] = cellData3.mobility(nPhaseIdx);
280
281 // compute total mobility of cell 3
282 Scalar lambdaTotal3 = lambda3[wPhaseIdx] + lambda3[nPhaseIdx];
283
284 // get mobilities of the phases
285 Dune::FieldVector<Scalar, numPhases> lambda4(cellData4.mobility(wPhaseIdx));
286 lambda4[nPhaseIdx] = cellData4.mobility(nPhaseIdx);
287
288 // compute total mobility of cell 4
289 Scalar lambdaTotal4 = lambda4[wPhaseIdx] + lambda4[nPhaseIdx];
290
291 // get mobilities of the phases
292 Dune::FieldVector<Scalar, numPhases> lambda5(cellData5.mobility(wPhaseIdx));
293 lambda5[nPhaseIdx] = cellData5.mobility(nPhaseIdx);
294
295 // compute total mobility of cell 5
296 Scalar lambdaTotal5 = lambda5[wPhaseIdx] + lambda5[nPhaseIdx];
297
298 // get mobilities of the phases
299 Dune::FieldVector<Scalar, numPhases> lambda6(cellData6.mobility(wPhaseIdx));
300 lambda6[nPhaseIdx] = cellData6.mobility(nPhaseIdx);
301
302 // compute total mobility of cell 6
303 Scalar lambdaTotal6 = lambda6[wPhaseIdx] + lambda6[nPhaseIdx];
304
305 // get mobilities of the phases
306 Dune::FieldVector<Scalar, numPhases> lambda7(cellData7.mobility(wPhaseIdx));
307 lambda7[nPhaseIdx] = cellData7.mobility(nPhaseIdx);
308
309 // compute total mobility of cell 7
310 Scalar lambdaTotal7 = lambda7[wPhaseIdx] + lambda7[nPhaseIdx];
311
312 // get mobilities of the phases
313 Dune::FieldVector<Scalar, numPhases> lambda8(cellData8.mobility(wPhaseIdx));
314 lambda8[nPhaseIdx] = cellData8.mobility(nPhaseIdx);
315
316 // compute total mobility of cell 8
317 Scalar lambdaTotal8 = lambda8[wPhaseIdx] + lambda8[nPhaseIdx];
318
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;
344
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;
357
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;
370
371 //flux vector
372 Dune::FieldVector<Scalar, 12> fluxW(0);
373 Dune::FieldVector<Scalar, 12> fluxNw(0);
374
375 DimVector Tu(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);
380
381 int hangingNodeType = interactionVolume.getHangingNodeType();
382
383
384
385 if (fIdx < 0 || fIdx == 0)
386 {
387 // calculate the flux through the subvolumeface 1 (subVolumeFaceIdx = 0)
388 int caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 0, 1, 2, 3,
389 4, 5);
390
391 if (caseL == 1)
392 {
393 u[0] = potW[0];
394 u[1] = potW[1];
395 u[2] = potW[2];
396 u[3] = potW[4];
397
398 T.mv(u, Tu);
399
400 fluxW[0] = Tu[0];
401 potentialDiffW0 = Tu[0];
402
403 u[0] = potNw[0];
404 u[1] = potNw[1];
405 u[2] = potNw[2];
406 u[3] = potNw[4];
407
408 T.mv(u, Tu);
409
410 fluxNw[0] = Tu[0];
411 potentialDiffNw0 = Tu[0];
412 }
413 else if (caseL == 2)
414 {
415 u[0] = potW[0];
416 u[1] = potW[1];
417 u[2] = potW[3];
418 u[3] = potW[5];
419
420 T.mv(u, Tu);
421
422 fluxW[0] = Tu[0];
423 potentialDiffW0 = Tu[0];
424
425 u[0] = potNw[0];
426 u[1] = potNw[1];
427 u[2] = potNw[3];
428 u[3] = potNw[5];
429
430 T.mv(u, Tu);
431
432 fluxNw[0] = Tu[0];
433 potentialDiffNw0 = Tu[0];
434 }
435 else if (caseL == 3)
436 {
437 u[0] = potW[0];
438 u[1] = potW[1];
439 u[2] = potW[3];
440 u[3] = potW[4];
441
442 T.mv(u, Tu);
443
444 fluxW[0] = Tu[0];
445 potentialDiffW0 = Tu[0];
446
447 u[0] = potNw[0];
448 u[1] = potNw[1];
449 u[2] = potNw[3];
450 u[3] = potNw[4];
451
452 T.mv(u, Tu);
453
454 fluxNw[0] = Tu[0];
455 potentialDiffNw0 = Tu[0];
456 }
457 else
458 {
459 u[0] = potW[0];
460 u[1] = potW[1];
461 u[2] = potW[2];
462 u[3] = potW[5];
463
464 T.mv(u, Tu);
465
466 fluxW[0] = Tu[0];
467 potentialDiffW0 = Tu[0];
468
469 u[0] = potNw[0];
470 u[1] = potNw[1];
471 u[2] = potNw[2];
472 u[3] = potNw[5];
473
474 T.mv(u, Tu);
475
476 fluxNw[0] = Tu[0];
477 potentialDiffNw0 = Tu[0];
478 }
479 }
480
481 if (fIdx < 0 || fIdx == 1)
482 {
483 // calculate the flux through the subvolumeface 2 (subVolumeFaceIdx = 1)
484 T = 0;
485 int caseL = 0;
486 if (hangingNodeType == InteractionVolume::twoSmallCells
487 || hangingNodeType == InteractionVolume::fourSmallCellsDiag)
488 {
489 useCases[0] = true;
490 useCases[1] = false;
491 useCases[2] = false;
492 useCases[3] = true;
493 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 1, 3, 0,
494 2, 5, 7, useCases);
495 }
496 else
497 {
498 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 1, 3, 0,
499 2, 5, 7);
500 }
501
502 if (caseL == 1)
503 {
504 u[0] = potW[1];
505 u[1] = potW[3];
506 u[2] = potW[0];
507 u[3] = potW[5];
508
509 T.mv(u, Tu);
510
511 fluxW[1] = Tu[0];
512 potentialDiffW1 = Tu[0];
513
514 u[0] = potNw[1];
515 u[1] = potNw[3];
516 u[2] = potNw[0];
517 u[3] = potNw[5];
518
519 T.mv(u, Tu);
520
521 fluxNw[1] = Tu[0];
522 potentialDiffNw1 = Tu[0];
523 }
524 else if (caseL == 2)
525 {
526 u[0] = potW[1];
527 u[1] = potW[3];
528 u[2] = potW[2];
529 u[3] = potW[7];
530
531 T.mv(u, Tu);
532
533 fluxW[1] = Tu[0];
534 potentialDiffW1 = Tu[0];
535
536 u[0] = potNw[1];
537 u[1] = potNw[3];
538 u[2] = potNw[2];
539 u[3] = potNw[7];
540
541 T.mv(u, Tu);
542
543 fluxNw[1] = Tu[0];
544 potentialDiffNw1 = Tu[0];
545 }
546 else if (caseL == 3)
547 {
548 u[0] = potW[1];
549 u[1] = potW[3];
550 u[2] = potW[2];
551 u[3] = potW[5];
552
553 T.mv(u, Tu);
554
555 fluxW[1] = Tu[0];
556 potentialDiffW1 = Tu[0];
557
558 u[0] = potNw[1];
559 u[1] = potNw[3];
560 u[2] = potNw[2];
561 u[3] = potNw[5];
562
563 T.mv(u, Tu);
564
565 fluxNw[1] = Tu[0];
566 potentialDiffNw1 = Tu[0];
567 }
568 else
569 {
570 u[0] = potW[1];
571 u[1] = potW[3];
572 u[2] = potW[0];
573 u[3] = potW[7];
574
575 T.mv(u, Tu);
576
577 fluxW[1] = Tu[0];
578 potentialDiffW1 = Tu[0];
579
580 u[0] = potNw[1];
581 u[1] = potNw[3];
582 u[2] = potNw[0];
583 u[3] = potNw[7];
584
585 T.mv(u, Tu);
586
587 fluxNw[1] = Tu[0];
588 potentialDiffNw1 = Tu[0];
589 }
590 }
591
592 if (fIdx < 0 || fIdx == 2)
593 {
594 // calculate the flux through the subvolumeface 3 (subVolumeFaceIdx = 2)
595 T = 0;
596 int caseL = 0;
597 if (hangingNodeType != InteractionVolume::twoSmallCells
598 && hangingNodeType != InteractionVolume::fourSmallCellsDiag)
599 {
600 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 3,
601 2, 1, 0, 7, 6);//, false);
602
603 if (caseL == 1)
604 {
605 u[0] = potW[3];
606 u[1] = potW[2];
607 u[2] = potW[1];
608 u[3] = potW[7];
609
610 T.mv(u, Tu);
611
612 fluxW[2] = Tu[0];
613 potentialDiffW2 = Tu[0];
614
615 u[0] = potNw[3];
616 u[1] = potNw[2];
617 u[2] = potNw[1];
618 u[3] = potNw[7];
619
620 T.mv(u, Tu);
621
622 fluxNw[2] = Tu[0];
623 potentialDiffNw2 = Tu[0];
624 }
625 else if (caseL == 2)
626 {
627 u[0] = potW[3];
628 u[1] = potW[2];
629 u[2] = potW[0];
630 u[3] = potW[6];
631
632 T.mv(u, Tu);
633
634 fluxW[2] = Tu[0];
635 potentialDiffW2 = Tu[0];
636
637 u[0] = potNw[3];
638 u[1] = potNw[2];
639 u[2] = potNw[0];
640 u[3] = potNw[6];
641
642 T.mv(u, Tu);
643
644 fluxNw[2] = Tu[0];
645 potentialDiffNw2 = Tu[0];
646 }
647 else if (caseL == 3)
648 {
649 u[0] = potW[3];
650 u[1] = potW[2];
651 u[2] = potW[0];
652 u[3] = potW[7];
653
654 T.mv(u, Tu);
655
656 fluxW[2] = Tu[0];
657 potentialDiffW2 = Tu[0];
658
659 u[0] = potNw[3];
660 u[1] = potNw[2];
661 u[2] = potNw[0];
662 u[3] = potNw[7];
663
664 T.mv(u, Tu);
665
666 fluxNw[2] = Tu[0];
667 potentialDiffNw2 = Tu[0];
668 }
669 else
670 {
671 u[0] = potW[3];
672 u[1] = potW[2];
673 u[2] = potW[1];
674 u[3] = potW[6];
675
676 T.mv(u, Tu);
677
678 fluxW[2] = Tu[0];
679 potentialDiffW2 = Tu[0];
680
681 u[0] = potNw[3];
682 u[1] = potNw[2];
683 u[2] = potNw[1];
684 u[3] = potNw[6];
685
686 T.mv(u, Tu);
687
688 fluxNw[2] = Tu[0];
689 potentialDiffNw2 = Tu[0];
690 }
691 }
692 }
693
694 if (fIdx < 0 || fIdx == 3)
695 {
696 // calculate the flux through the subvolumeface 4 (subVolumeFaceIdx = 3)
697 T = 0;
698 int caseL = 0;
699 if (hangingNodeType == InteractionVolume::twoSmallCells
700 || hangingNodeType == InteractionVolume::fourSmallCellsDiag)
701 {
702 useCases[0] = false;
703 useCases[1] = true;
704 useCases[2] = true;
705 useCases[3] = false;
706 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 2, 0, 3, 1, 6,
707 4, useCases);
708 }
709 else
710 {
711 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 2,
712 0, 3, 1, 6, 4);//, false);
713 }
714
715 if (caseL == 1)
716 {
717 u[0] = potW[2];
718 u[1] = potW[0];
719 u[2] = potW[3];
720 u[3] = potW[6];
721
722 T.mv(u, Tu);
723
724 fluxW[3] = Tu[0];
725 potentialDiffW3 = Tu[0];
726
727 u[0] = potNw[2];
728 u[1] = potNw[0];
729 u[2] = potNw[3];
730 u[3] = potNw[6];
731
732 T.mv(u, Tu);
733
734 fluxNw[3] = Tu[0];
735 potentialDiffNw3 = Tu[0];
736 }
737 else if (caseL == 2)
738 {
739 u[0] = potW[2];
740 u[1] = potW[0];
741 u[2] = potW[1];
742 u[3] = potW[4];
743
744 T.mv(u, Tu);
745
746 fluxW[3] = Tu[0];
747 potentialDiffW3 = Tu[0];
748
749 u[0] = potNw[2];
750 u[1] = potNw[0];
751 u[2] = potNw[1];
752 u[3] = potNw[4];
753
754 T.mv(u, Tu);
755
756 fluxNw[3] = Tu[0];
757 potentialDiffNw3 = Tu[0];
758 }
759 else if (caseL == 3)
760 {
761 u[0] = potW[2];
762 u[1] = potW[0];
763 u[2] = potW[1];
764 u[3] = potW[6];
765
766 T.mv(u, Tu);
767
768 fluxW[3] = Tu[0];
769 potentialDiffW3 = Tu[0];
770
771 u[0] = potNw[2];
772 u[1] = potNw[0];
773 u[2] = potNw[1];
774 u[3] = potNw[6];
775
776 T.mv(u, Tu);
777
778 fluxNw[3] = Tu[0];
779 potentialDiffNw3 = Tu[0];
780 }
781 else
782 {
783 u[0] = potW[2];
784 u[1] = potW[0];
785 u[2] = potW[3];
786 u[3] = potW[4];
787
788 T.mv(u, Tu);
789
790 fluxW[3] = Tu[0];
791 potentialDiffW3 = Tu[0];
792
793 u[0] = potNw[2];
794 u[1] = potNw[0];
795 u[2] = potNw[3];
796 u[3] = potNw[4];
797
798 T.mv(u, Tu);
799
800 fluxNw[3] = Tu[0];
801 potentialDiffNw3 = Tu[0];
802 }
803 }
804
805 if (fIdx < 0 || fIdx == 4)
806 {
807
808 // calculate the flux through the subvolumeface 5 (subVolumeFaceIdx = 4)
809 T = 0;
810 TSecond = 0;
811 int caseL = 0;
812 if (hangingNodeType == InteractionVolume::sixSmallCells)
813 {
814 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 5, 4, 7,
815 6, 1, 0);
816 }
817 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx6)
818 {
819 caseL = transmissibilityCalculator_().transmissibilityCaseThree(T, interactionVolume, lambda, 4, 0, 2,
820 5);
821
822 int caseLSecond = transmissibilityCalculator_().transmissibilityCaseFour(TSecond, interactionVolume, lambda, 1, 5, 3,
823 4);
824
825 caseL = transmissibilityCalculator_().chooseTransmissibility(T, TSecond, 3, 4);
826
827 if (caseL == caseLSecond)
828 T = TSecond;
829 }
830 if (hangingNodeType == InteractionVolume::sixSmallCells)
831 {
832 if (caseL == 1)
833 {
834 u[0] = potW[5];
835 u[1] = potW[4];
836 u[2] = potW[7];
837 u[3] = potW[1];
838
839 T.mv(u, Tu);
840
841 fluxW[4] = Tu[0];
842 potentialDiffW4 = Tu[0];
843
844 u[0] = potNw[5];
845 u[1] = potNw[4];
846 u[2] = potNw[7];
847 u[3] = potNw[1];
848
849 T.mv(u, Tu);
850
851 fluxNw[4] = Tu[0];
852 potentialDiffNw4 = Tu[0];
853 }
854 else if (caseL == 2)
855 {
856 u[0] = potW[5];
857 u[1] = potW[4];
858 u[2] = potW[6];
859 u[3] = potW[0];
860
861 T.mv(u, Tu);
862
863 fluxW[4] = Tu[0];
864 potentialDiffW4 = Tu[0];
865
866 u[0] = potNw[5];
867 u[1] = potNw[4];
868 u[2] = potNw[6];
869 u[3] = potNw[0];
870
871 T.mv(u, Tu);
872
873 fluxNw[4] = Tu[0];
874 potentialDiffNw4 = Tu[0];
875 }
876 else if (caseL == 3)
877 {
878 u[0] = potW[5];
879 u[1] = potW[4];
880 u[2] = potW[6];
881 u[3] = potW[1];
882
883 T.mv(u, Tu);
884
885 fluxW[4] = Tu[0];
886 potentialDiffW4 = Tu[0];
887
888 u[0] = potNw[5];
889 u[1] = potNw[4];
890 u[2] = potNw[6];
891 u[3] = potNw[1];
892
893 T.mv(u, Tu);
894
895 fluxNw[4] = Tu[0];
896 potentialDiffNw4 = Tu[0];
897 }
898 else
899 {
900 u[0] = potW[5];
901 u[1] = potW[4];
902 u[2] = potW[7];
903 u[3] = potW[0];
904
905 T.mv(u, Tu);
906
907 fluxW[4] = Tu[0];
908 potentialDiffW4 = Tu[0];
909
910 u[0] = potNw[5];
911 u[1] = potNw[4];
912 u[2] = potNw[7];
913 u[3] = potNw[0];
914
915 T.mv(u, Tu);
916
917 fluxNw[4] = Tu[0];
918 potentialDiffNw4 = Tu[0];
919 }
920 }
921 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx6)
922 {
923 if (caseL == 3)
924 {
925 u[0] = potW[4];
926 u[1] = potW[0];
927 u[2] = potW[2];
928 u[3] = potW[5];
929
930 T.mv(u, Tu);
931
932 fluxW[4] = -Tu[2];
933 potentialDiffW4 = -Tu[2];
934
935 u[0] = potNw[4];
936 u[1] = potNw[0];
937 u[2] = potNw[2];
938 u[3] = potNw[5];
939
940 T.mv(u, Tu);
941
942 fluxNw[4] = -Tu[2];
943 potentialDiffNw4 = -Tu[2];
944 }
945 else if (caseL == 4)
946 {
947 u[0] = potW[1];
948 u[1] = potW[5];
949 u[2] = potW[3];
950 u[3] = potW[4];
951
952 T.mv(u, Tu);
953
954 fluxW[4] = Tu[2];
955 potentialDiffW4 = Tu[2];
956
957 u[0] = potNw[1];
958 u[1] = potNw[5];
959 u[2] = potNw[3];
960 u[3] = potNw[4];
961
962 T.mv(u, Tu);
963
964 fluxNw[4] = Tu[2];
965 potentialDiffNw4 = Tu[2];
966 }
967 }
968 }
969
970 if (fIdx < 0 || fIdx == 5)
971 {
972 // calculate the flux through the subvolumeface 6 (subVolumeFaceIdx = 5)
973 T = 0;
974 TSecond = 0;
975 int caseL = 0;
976 if (hangingNodeType == InteractionVolume::sixSmallCells)
977 {
978 useCases[0] = false;
979 useCases[1] = true;
980 useCases[2] = true;
981 useCases[3] = false;
982
983 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 7, 5, 6, 4, 3,
984 1, useCases);
985 }
986 else if (hangingNodeType == InteractionVolume::fourSmallCellsDiag)
987 {
988 useCases[0] = true;
989 useCases[1] = false;
990 useCases[2] = false;
991 useCases[3] = true;
992 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 7, 5, 6, 4, 3,
993 1, useCases);
994 }
995 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx7)
996 {
997 caseL = transmissibilityCalculator_().transmissibilityCaseThree(T, interactionVolume, lambda, 1, 5, 7, 0);
998
999 int caseLSecond = transmissibilityCalculator_().transmissibilityCaseFour(TSecond, interactionVolume, lambda, 7, 3, 5,
1000 2);
1001
1002 caseL = transmissibilityCalculator_().chooseTransmissibility(T, TSecond, 3, 4);
1003
1004 if (caseL == caseLSecond)
1005 T = TSecond;
1006 }
1007
1008 if (hangingNodeType == InteractionVolume::sixSmallCells
1009 || hangingNodeType == InteractionVolume::fourSmallCellsDiag)
1010 {
1011 if (caseL == 1)
1012 {
1013 u[0] = potW[7];
1014 u[1] = potW[5];
1015 u[2] = potW[6];
1016 u[3] = potW[3];
1017
1018 T.mv(u, Tu);
1019
1020 fluxW[5] = Tu[0];
1021 potentialDiffW5 = Tu[0];
1022
1023 u[0] = potNw[7];
1024 u[1] = potNw[5];
1025 u[2] = potNw[6];
1026 u[3] = potNw[3];
1027
1028 T.mv(u, Tu);
1029
1030 fluxNw[5] = Tu[0];
1031 potentialDiffNw5 = Tu[0];
1032 }
1033 else if (caseL == 2)
1034 {
1035 u[0] = potW[7];
1036 u[1] = potW[5];
1037 u[2] = potW[4];
1038 u[3] = potW[1];
1039
1040 T.mv(u, Tu);
1041
1042 fluxW[5] = Tu[0];
1043 potentialDiffW5 = Tu[0];
1044
1045 u[0] = potNw[7];
1046 u[1] = potNw[5];
1047 u[2] = potNw[4];
1048 u[3] = potNw[1];
1049
1050 T.mv(u, Tu);
1051
1052 fluxNw[5] = Tu[0];
1053 potentialDiffNw5 = Tu[0];
1054 }
1055 else if (caseL == 3)
1056 {
1057 u[0] = potW[7];
1058 u[1] = potW[5];
1059 u[2] = potW[4];
1060 u[3] = potW[3];
1061
1062 T.mv(u, Tu);
1063
1064 fluxW[5] = Tu[0];
1065 potentialDiffW5 = Tu[0];
1066
1067 u[0] = potNw[7];
1068 u[1] = potNw[5];
1069 u[2] = potNw[4];
1070 u[3] = potNw[3];
1071
1072 T.mv(u, Tu);
1073
1074 fluxNw[5] = Tu[0];
1075 potentialDiffNw5 = Tu[0];
1076 }
1077 else
1078 {
1079 u[0] = potW[7];
1080 u[1] = potW[5];
1081 u[2] = potW[6];
1082 u[3] = potW[1];
1083
1084 T.mv(u, Tu);
1085
1086 fluxW[5] = Tu[0];
1087 potentialDiffW5 = Tu[0];
1088
1089 u[0] = potNw[7];
1090 u[1] = potNw[5];
1091 u[2] = potNw[6];
1092 u[3] = potNw[1];
1093
1094 T.mv(u, Tu);
1095
1096 fluxNw[5] = Tu[0];
1097 potentialDiffNw5 = Tu[0];
1098 }
1099 }
1100 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx7)
1101 {
1102 if (caseL == 3)
1103 {
1104 u[0] = potW[1];
1105 u[1] = potW[5];
1106 u[2] = potW[7];
1107 u[3] = potW[0];
1108
1109 T.mv(u, Tu);
1110
1111 fluxW[5] = -Tu[1];
1112 potentialDiffW5 = -Tu[1];
1113
1114 u[0] = potNw[1];
1115 u[1] = potNw[5];
1116 u[2] = potNw[7];
1117 u[3] = potNw[0];
1118
1119 T.mv(u, Tu);
1120
1121 fluxNw[5] = -Tu[1];
1122 potentialDiffNw5 = -Tu[1];
1123 }
1124 else
1125 {
1126 u[0] = potW[7];
1127 u[1] = potW[3];
1128 u[2] = potW[5];
1129 u[3] = potW[2];
1130
1131 T.mv(u, Tu);
1132
1133 fluxW[5] = Tu[1];
1134 potentialDiffW5 = Tu[1];
1135
1136 u[0] = potNw[7];
1137 u[1] = potNw[3];
1138 u[2] = potNw[5];
1139 u[3] = potNw[2];
1140
1141 T.mv(u, Tu);
1142
1143 fluxNw[5] = Tu[1];
1144 potentialDiffNw5 = Tu[1];
1145 }
1146 }
1147 }
1148
1149 if (fIdx < 0 || fIdx == 6)
1150 {
1151 // calculate the flux through the subvolumeface 7 (subVolumeFaceIdx = 6)
1152 T = 0;
1153 TSecond = 0;
1154 int caseL = 0;
1155 if (hangingNodeType == InteractionVolume::fourSmallCellsDiag)
1156 {
1157 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 6, 7, 4,
1158 5, 2, 3);
1159 }
1160 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx6)
1161 {
1162 caseL = transmissibilityCalculator_().transmissibilityCaseThree(T, interactionVolume, lambda, 7, 3, 1, 6);
1163
1164 int caseLSecond = transmissibilityCalculator_().transmissibilityCaseFour(TSecond, interactionVolume, lambda, 2, 6, 0,
1165 7);
1166
1167 caseL = transmissibilityCalculator_().chooseTransmissibility(T, TSecond, 3, 4);
1168
1169 if (caseL == caseLSecond)
1170 T = TSecond;
1171 }
1172
1173 if (hangingNodeType == InteractionVolume::fourSmallCellsDiag)
1174 {
1175 if (caseL == 1)
1176 {
1177 u[0] = potW[6];
1178 u[1] = potW[7];
1179 u[2] = potW[4];
1180 u[3] = potW[2];
1181
1182 T.mv(u, Tu);
1183
1184 fluxW[6] = Tu[0];
1185 potentialDiffW6 = Tu[0];
1186
1187 u[0] = potNw[6];
1188 u[1] = potNw[7];
1189 u[2] = potNw[4];
1190 u[3] = potNw[2];
1191
1192 T.mv(u, Tu);
1193
1194 fluxNw[6] = Tu[0];
1195 potentialDiffNw6 = Tu[0];
1196 }
1197 else if (caseL == 2)
1198 {
1199 u[0] = potW[6];
1200 u[1] = potW[7];
1201 u[2] = potW[5];
1202 u[3] = potW[3];
1203
1204 T.mv(u, Tu);
1205
1206 fluxW[6] = Tu[0];
1207 potentialDiffW6 = Tu[0];
1208
1209 u[0] = potNw[6];
1210 u[1] = potNw[7];
1211 u[2] = potNw[5];
1212 u[3] = potNw[3];
1213
1214 T.mv(u, Tu);
1215
1216 fluxNw[6] = Tu[0];
1217 potentialDiffNw6 = Tu[0];
1218 }
1219 else if (caseL == 3)
1220 {
1221 u[0] = potW[6];
1222 u[1] = potW[7];
1223 u[2] = potW[5];
1224 u[3] = potW[2];
1225
1226 T.mv(u, Tu);
1227
1228 fluxW[6] = Tu[0];
1229 potentialDiffW6 = Tu[0];
1230
1231 u[0] = potNw[6];
1232 u[1] = potNw[7];
1233 u[2] = potNw[5];
1234 u[3] = potNw[2];
1235
1236 T.mv(u, Tu);
1237
1238 fluxNw[6] = Tu[0];
1239 potentialDiffNw6 = Tu[0];
1240 }
1241 else
1242 {
1243 u[0] = potW[6];
1244 u[1] = potW[7];
1245 u[2] = potW[4];
1246 u[3] = potW[3];
1247
1248 T.mv(u, Tu);
1249
1250 fluxW[6] = Tu[0];
1251 potentialDiffW6 = Tu[0];
1252
1253 u[0] = potNw[6];
1254 u[1] = potNw[7];
1255 u[2] = potNw[4];
1256 u[3] = potNw[3];
1257
1258 T.mv(u, Tu);
1259
1260 fluxNw[6] = Tu[0];
1261 potentialDiffNw6 = Tu[0];
1262 }
1263 }
1264 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx6)
1265 {
1266 if (caseL == 3)
1267 {
1268 u[0] = potW[7];
1269 u[1] = potW[3];
1270 u[2] = potW[1];
1271 u[3] = potW[6];
1272
1273 T.mv(u, Tu);
1274
1275 fluxW[6] = -Tu[2];
1276 potentialDiffW6 = -Tu[2];
1277
1278 u[0] = potNw[7];
1279 u[1] = potNw[3];
1280 u[2] = potNw[1];
1281 u[3] = potNw[6];
1282
1283 T.mv(u, Tu);
1284
1285 fluxNw[6] = -Tu[2];
1286 potentialDiffNw6 = -Tu[2];
1287 }
1288 else if (caseL == 4)
1289 {
1290 u[0] = potW[2];
1291 u[1] = potW[6];
1292 u[2] = potW[0];
1293 u[3] = potW[7];
1294
1295 T.mv(u, Tu);
1296
1297 fluxW[6] = Tu[2];
1298 potentialDiffW6 = Tu[2];
1299
1300 u[0] = potNw[2];
1301 u[1] = potNw[6];
1302 u[2] = potNw[0];
1303 u[3] = potNw[7];
1304
1305 T.mv(u, Tu);
1306
1307 fluxNw[6] = Tu[2];
1308 potentialDiffNw6 = Tu[2];
1309 }
1310 }
1311 }
1312
1313 if (fIdx < 0 || fIdx == 7)
1314 {
1315 // calculate the flux through the subvolumeface 8 (subVolumeFaceIdx = 7)
1316 T = 0;
1317 TSecond = 0;
1318 int caseL = 0;
1319 if (hangingNodeType == InteractionVolume::sixSmallCells)
1320 {
1321 useCases[0] = true;
1322 useCases[1] = false;
1323 useCases[2] = false;
1324 useCases[3] = true;
1325
1326 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 4, 6, 5, 7, 0,
1327 2, useCases);
1328 }
1329 else if (hangingNodeType == InteractionVolume::fourSmallCellsDiag)
1330 {
1331 useCases[0] = false;
1332 useCases[1] = true;
1333 useCases[2] = true;
1334 useCases[3] = false;
1335 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 4, 6, 5, 7, 0,
1336 2, useCases);
1337 }
1338 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx7)
1339 {
1340 caseL = transmissibilityCalculator_().transmissibilityCaseThree(T, interactionVolume, lambda, 2, 6, 4,
1341 3);
1342
1343 int caseLSecond = transmissibilityCalculator_().transmissibilityCaseFour(TSecond, interactionVolume, lambda, 4, 0, 6,
1344 1);
1345
1346
1347
1348 caseL = transmissibilityCalculator_().chooseTransmissibility(T, TSecond, 3, 4);
1349
1350 if (caseL == caseLSecond)
1351 T = TSecond;
1352 }
1353
1354 if (hangingNodeType == InteractionVolume::sixSmallCells
1355 || hangingNodeType == InteractionVolume::fourSmallCellsDiag)
1356 {
1357
1358 if (caseL == 1)
1359 {
1360 u[0] = potW[4];
1361 u[1] = potW[6];
1362 u[2] = potW[5];
1363 u[3] = potW[0];
1364
1365 T.mv(u, Tu);
1366
1367 fluxW[7] = Tu[0];
1368 potentialDiffW7 = Tu[0];
1369
1370 u[0] = potNw[4];
1371 u[1] = potNw[6];
1372 u[2] = potNw[5];
1373 u[3] = potNw[0];
1374
1375 T.mv(u, Tu);
1376
1377 fluxNw[7] = Tu[0];
1378 potentialDiffNw7 = Tu[0];
1379 }
1380 else if (caseL == 2)
1381 {
1382 u[0] = potW[4];
1383 u[1] = potW[6];
1384 u[2] = potW[7];
1385 u[3] = potW[2];
1386
1387 T.mv(u, Tu);
1388
1389 fluxW[7] = Tu[0];
1390 potentialDiffW7 = Tu[0];
1391
1392 u[0] = potNw[4];
1393 u[1] = potNw[6];
1394 u[2] = potNw[7];
1395 u[3] = potNw[2];
1396
1397 T.mv(u, Tu);
1398
1399 fluxNw[7] = Tu[0];
1400 potentialDiffNw7 = Tu[0];
1401 }
1402 else if (caseL == 3)
1403 {
1404 u[0] = potW[4];
1405 u[1] = potW[6];
1406 u[2] = potW[7];
1407 u[3] = potW[0];
1408
1409 T.mv(u, Tu);
1410
1411 fluxW[7] = Tu[0];
1412 potentialDiffW7 = Tu[0];
1413
1414 u[0] = potNw[4];
1415 u[1] = potNw[6];
1416 u[2] = potNw[7];
1417 u[3] = potNw[0];
1418
1419 T.mv(u, Tu);
1420
1421 fluxNw[7] = Tu[0];
1422 potentialDiffNw7 = Tu[0];
1423 }
1424 else
1425 {
1426 u[0] = potW[4];
1427 u[1] = potW[6];
1428 u[2] = potW[5];
1429 u[3] = potW[2];
1430
1431 T.mv(u, Tu);
1432
1433 fluxW[7] = Tu[0];
1434 potentialDiffW7 = Tu[0];
1435
1436 u[0] = potNw[4];
1437 u[1] = potNw[6];
1438 u[2] = potNw[5];
1439 u[3] = potNw[2];
1440
1441 T.mv(u, Tu);
1442
1443 fluxNw[7] = Tu[0];
1444 potentialDiffNw7 = Tu[0];
1445 }
1446 }
1447 else if(hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx7)
1448 {
1449 if (caseL == 3)
1450 {
1451 u[0] = potW[2];
1452 u[1] = potW[6];
1453 u[2] = potW[4];
1454 u[3] = potW[3];
1455
1456 T.mv(u, Tu);
1457
1458 fluxW[7] = -Tu[1];
1459 potentialDiffW7 = -Tu[1];
1460
1461 u[0] = potNw[2];
1462 u[1] = potNw[6];
1463 u[2] = potNw[4];
1464 u[3] = potNw[3];
1465
1466 T.mv(u, Tu);
1467
1468 fluxNw[7] = -Tu[1];
1469 potentialDiffNw7 = -Tu[1];
1470 }
1471 else if (caseL == 4)
1472 {
1473 u[0] = potW[4];
1474 u[1] = potW[0];
1475 u[2] = potW[6];
1476 u[3] = potW[1];
1477
1478 T.mv(u, Tu);
1479
1480 fluxW[7] = Tu[1];
1481 potentialDiffW7 = Tu[1];
1482
1483 u[0] = potNw[4];
1484 u[1] = potNw[0];
1485 u[2] = potNw[6];
1486 u[3] = potNw[1];
1487
1488 T.mv(u, Tu);
1489
1490 fluxNw[7] = Tu[1];
1491 potentialDiffNw7 = Tu[1];
1492 }
1493 }
1494 }
1495
1496 if (fIdx < 0 || fIdx == 8)
1497 {
1498 // calculate the flux through the subvolumeface 9 (subVolumeFaceIdx = 8)
1499 T = 0;
1500 int caseL = 0;
1501 if (hangingNodeType == InteractionVolume::sixSmallCells)
1502 {
1503 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 4, 0, 6,
1504 2, 5, 1);
1505 }
1506 else if (hangingNodeType == InteractionVolume::twoSmallCells
1507 || hangingNodeType == InteractionVolume::fourSmallCellsFace)
1508 {
1509 caseL = transmissibilityCalculator_().transmissibilityCaseTwo(T, interactionVolume, lambda, 4, 0, 2,
1510 1);
1511 }
1512 else if (hangingNodeType == InteractionVolume::fourSmallCellsDiag
1513 || (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx7))
1514 {
1515 useCases[0] = false;
1516 useCases[1] = true;
1517 useCases[2] = false;
1518 useCases[3] = true;
1519
1520 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 4, 0, 6,
1521 2, 5, 1, useCases);
1522 }
1523 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx6)
1524 {
1525 useCases[0] = false;
1526 useCases[1] = true;
1527 useCases[2] = true;
1528 useCases[3] = false;
1529
1530 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 4, 0, 6,
1531 2, 5, 1, useCases);
1532 }
1533
1534 if (caseL == 1)
1535 {
1536 u[0] = potW[4];
1537 u[1] = potW[0];
1538 u[2] = potW[6];
1539 u[3] = potW[5];
1540
1541 T.mv(u, Tu);
1542
1543 fluxW[8] = Tu[0];
1544 potentialDiffW8 = Tu[0];
1545
1546 u[0] = potNw[4];
1547 u[1] = potNw[0];
1548 u[2] = potNw[6];
1549 u[3] = potNw[5];
1550
1551 T.mv(u, Tu);
1552
1553 fluxNw[8] = Tu[0];
1554 potentialDiffNw8 = Tu[0];
1555 }
1556 else if (caseL == 2)
1557 {
1558 u[0] = potW[4];
1559 u[1] = potW[0];
1560 u[2] = potW[2];
1561 u[3] = potW[1];
1562
1563 T.mv(u, Tu);
1564
1565 fluxW[8] = Tu[0];
1566 potentialDiffW8 = Tu[0];
1567
1568 u[0] = potNw[4];
1569 u[1] = potNw[0];
1570 u[2] = potNw[2];
1571 u[3] = potNw[1];
1572
1573 T.mv(u, Tu);
1574
1575 fluxNw[8] = Tu[0];
1576 potentialDiffNw8 = Tu[0];
1577 }
1578 else if (caseL == 3)
1579 {
1580 u[0] = potW[4];
1581 u[1] = potW[0];
1582 u[2] = potW[2];
1583 u[3] = potW[5];
1584
1585 T.mv(u, Tu);
1586
1587 fluxW[8] = Tu[0];
1588 potentialDiffW8 = Tu[0];
1589
1590 u[0] = potNw[4];
1591 u[1] = potNw[0];
1592 u[2] = potNw[2];
1593 u[3] = potNw[5];
1594
1595 T.mv(u, Tu);
1596
1597 fluxNw[8] = Tu[0];
1598 potentialDiffNw8 = Tu[0];
1599 }
1600 else
1601 {
1602 u[0] = potW[4];
1603 u[1] = potW[0];
1604 u[2] = potW[6];
1605 u[3] = potW[1];
1606
1607 T.mv(u, Tu);
1608
1609 fluxW[8] = Tu[0];
1610 potentialDiffW8 = Tu[0];
1611
1612 u[0] = potNw[4];
1613 u[1] = potNw[0];
1614 u[2] = potNw[6];
1615 u[3] = potNw[1];
1616
1617 T.mv(u, Tu);
1618
1619 fluxNw[8] = Tu[0];
1620 potentialDiffNw8 = Tu[0];
1621 }
1622 }
1623
1624 if (fIdx < 0 || fIdx == 9)
1625 {
1626 // calculate the flux through the subvolumeface 10 (subVolumeFaceIdx = 9)
1627 T = 0;
1628 int caseL = 0;
1629 if (hangingNodeType == InteractionVolume::sixSmallCells)
1630 {
1631 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 1, 5, 3,
1632 7, 0, 4);
1633 }
1634 else if (hangingNodeType == InteractionVolume::twoSmallCells
1635 || hangingNodeType == InteractionVolume::fourSmallCellsFace)
1636 {
1637 caseL = transmissibilityCalculator_().transmissibilityCaseOne(T, interactionVolume, lambda, 1, 5, 3,
1638 0);
1639 }
1640 else if (hangingNodeType == InteractionVolume::fourSmallCellsDiag
1641 || (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx7))
1642 {
1643 useCases[0] = true;
1644 useCases[1] = false;
1645 useCases[2] = true;
1646 useCases[3] = false;
1647
1648 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 1, 5, 3,
1649 7, 0, 4, useCases);
1650 }
1651 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx6)
1652 {
1653 useCases[0] = true;
1654 useCases[1] = false;
1655 useCases[2] = false;
1656 useCases[3] = true;
1657
1658 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 1, 5, 3,
1659 7, 0, 4, useCases);
1660 }
1661
1662 if (caseL == 1)
1663 {
1664 u[0] = potW[1];
1665 u[1] = potW[5];
1666 u[2] = potW[3];
1667 u[3] = potW[0];
1668
1669 T.mv(u, Tu);
1670
1671 fluxW[9] = Tu[0];
1672 potentialDiffW9 = Tu[0];
1673
1674 u[0] = potNw[1];
1675 u[1] = potNw[5];
1676 u[2] = potNw[3];
1677 u[3] = potNw[0];
1678
1679 T.mv(u, Tu);
1680
1681 fluxNw[9] = Tu[0];
1682 potentialDiffNw9 = Tu[0];
1683 }
1684 else if (caseL == 2)
1685 {
1686 u[0] = potW[1];
1687 u[1] = potW[5];
1688 u[2] = potW[7];
1689 u[3] = potW[4];
1690
1691 T.mv(u, Tu);
1692
1693 fluxW[9] = Tu[0];
1694 potentialDiffW9 = Tu[0];
1695
1696 u[0] = potNw[1];
1697 u[1] = potNw[5];
1698 u[2] = potNw[7];
1699 u[3] = potNw[4];
1700
1701 T.mv(u, Tu);
1702
1703 fluxNw[9] = Tu[0];
1704 potentialDiffNw9 = Tu[0];
1705 }
1706 else if (caseL == 3)
1707 {
1708 u[0] = potW[1];
1709 u[1] = potW[5];
1710 u[2] = potW[7];
1711 u[3] = potW[0];
1712
1713 T.mv(u, Tu);
1714
1715 fluxW[9] = Tu[0];
1716 potentialDiffW9 = Tu[0];
1717
1718 u[0] = potNw[1];
1719 u[1] = potNw[5];
1720 u[2] = potNw[7];
1721 u[3] = potNw[0];
1722
1723 T.mv(u, Tu);
1724
1725 fluxNw[9] = Tu[0];
1726 potentialDiffNw9 = Tu[0];
1727 }
1728 else
1729 {
1730 u[0] = potW[1];
1731 u[1] = potW[5];
1732 u[2] = potW[3];
1733 u[3] = potW[4];
1734
1735 T.mv(u, Tu);
1736
1737 fluxW[9] = Tu[0];
1738 potentialDiffW9 = Tu[0];
1739
1740 u[0] = potNw[1];
1741 u[1] = potNw[5];
1742 u[2] = potNw[3];
1743 u[3] = potNw[4];
1744
1745 T.mv(u, Tu);
1746
1747 fluxNw[9] = Tu[0];
1748 potentialDiffNw9 = Tu[0];
1749 }
1750 }
1751
1752 if (fIdx < 0 || fIdx == 10)
1753 {
1754 // calculate the flux through the subvolumeface 11 (subVolumeFaceIdx = 10)
1755 T = 0;
1756 int caseL = 0;
1757 if (hangingNodeType == InteractionVolume::fourSmallCellsFace)
1758 {
1759 caseL = transmissibilityCalculator_().transmissibilityCaseTwo(T, interactionVolume, lambda, 7, 3, 1, 2);
1760 }
1761 else if ((hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx7)
1762 || hangingNodeType == InteractionVolume::sixSmallCells)
1763 {
1764 useCases[0] = false;
1765 useCases[1] = true;
1766 useCases[2] = false;
1767 useCases[3] = true;
1768
1769 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 7, 3, 5, 1, 6,
1770 2, useCases);
1771 }
1772 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx6)
1773 {
1774 useCases[0] = false;
1775 useCases[1] = true;
1776 useCases[2] = true;
1777 useCases[3] = false;
1778
1779 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 7, 3, 5, 1, 6,
1780 2, useCases);
1781 }
1782 else if (hangingNodeType == InteractionVolume::fourSmallCellsDiag)
1783 {
1784 useCases[0] = true;
1785 useCases[1] = false;
1786 useCases[2] = true;
1787 useCases[3] = false;
1788
1789 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 7, 3, 5, 1, 6,
1790 2, useCases);
1791 }
1792
1793 if (hangingNodeType != InteractionVolume::twoSmallCells)
1794 {
1795
1796 if (caseL == 1)
1797 {
1798 u[0] = potW[7];
1799 u[1] = potW[3];
1800 u[2] = potW[5];
1801 u[3] = potW[6];
1802
1803 T.mv(u, Tu);
1804
1805 fluxW[10] = Tu[0];
1806 potentialDiffW10 = Tu[0];
1807
1808 u[0] = potNw[7];
1809 u[1] = potNw[3];
1810 u[2] = potNw[5];
1811 u[3] = potNw[6];
1812
1813 T.mv(u, Tu);
1814
1815 fluxNw[10] = Tu[0];
1816 potentialDiffNw10 = Tu[0];
1817 }
1818 else if (caseL == 2)
1819 {
1820 u[0] = potW[7];
1821 u[1] = potW[3];
1822 u[2] = potW[1];
1823 u[3] = potW[2];
1824
1825 T.mv(u, Tu);
1826
1827 fluxW[10] = Tu[0];
1828 potentialDiffW10 = Tu[0];
1829
1830 u[0] = potNw[7];
1831 u[1] = potNw[3];
1832 u[2] = potNw[1];
1833 u[3] = potNw[2];
1834
1835 T.mv(u, Tu);
1836
1837 fluxNw[10] = Tu[0];
1838 potentialDiffNw10 = Tu[0];
1839 }
1840 else if (caseL == 3)
1841 {
1842 u[0] = potW[7];
1843 u[1] = potW[3];
1844 u[2] = potW[1];
1845 u[3] = potW[6];
1846
1847 T.mv(u, Tu);
1848
1849 fluxW[10] = Tu[0];
1850 potentialDiffW10 = Tu[0];
1851
1852 u[0] = potNw[7];
1853 u[1] = potNw[3];
1854 u[2] = potNw[1];
1855 u[3] = potNw[6];
1856
1857 T.mv(u, Tu);
1858
1859 fluxNw[10] = Tu[0];
1860 potentialDiffNw10 = Tu[0];
1861 }
1862 else
1863 {
1864 u[0] = potW[7];
1865 u[1] = potW[3];
1866 u[2] = potW[5];
1867 u[3] = potW[2];
1868
1869 T.mv(u, Tu);
1870
1871 fluxW[10] = Tu[0];
1872 potentialDiffW10 = Tu[0];
1873
1874 u[0] = potNw[7];
1875 u[1] = potNw[3];
1876 u[2] = potNw[5];
1877 u[3] = potNw[2];
1878
1879 T.mv(u, Tu);
1880
1881 fluxNw[10] = Tu[0];
1882 potentialDiffNw10 = Tu[0];
1883 }
1884 }
1885 }
1886
1887 if (fIdx < 0 || fIdx == 11)
1888 {
1889 // calculate the flux through the subvolumeface 12 (subVolumeFaceIdx = 11)
1890 T = 0;
1891 int caseL = 0;
1892 if (hangingNodeType == InteractionVolume::fourSmallCellsFace)
1893 {
1894 caseL = transmissibilityCalculator_().transmissibilityCaseOne(T, interactionVolume, lambda, 2, 6, 0, 3);
1895 }
1896 else if ((hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx7)
1897 || hangingNodeType == InteractionVolume::sixSmallCells)
1898 {
1899 useCases[0] = true;
1900 useCases[1] = false;
1901 useCases[2] = true;
1902 useCases[3] = false;
1903
1904 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 2, 6, 0, 4, 3,
1905 7, useCases);
1906 }
1907 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx6)
1908 {
1909 useCases[0] = true;
1910 useCases[1] = false;
1911 useCases[2] = false;
1912 useCases[3] = true;
1913
1914 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 2, 6, 0, 4, 3,
1915 7, useCases);
1916 }
1917 else if (hangingNodeType == InteractionVolume::fourSmallCellsDiag)
1918 {
1919 useCases[0] = false;
1920 useCases[1] = true;
1921 useCases[2] = false;
1922 useCases[3] = true;
1923
1924 caseL = transmissibilityCalculator_().transmissibility(T, interactionVolume, lambda, 2, 6, 0, 4, 3,
1925 7, useCases);
1926 }
1927
1928 if (hangingNodeType != InteractionVolume::twoSmallCells)
1929 {
1930
1931 if (caseL == 1)
1932 {
1933 u[0] = potW[2];
1934 u[1] = potW[6];
1935 u[2] = potW[0];
1936 u[3] = potW[3];
1937
1938 T.mv(u, Tu);
1939
1940 fluxW[11] = Tu[0];
1941 potentialDiffW11 = Tu[0];
1942
1943 u[0] = potNw[2];
1944 u[1] = potNw[6];
1945 u[2] = potNw[0];
1946 u[3] = potNw[3];
1947
1948 T.mv(u, Tu);
1949
1950 fluxNw[11] = Tu[0];
1951 potentialDiffNw11 = Tu[0];
1952 }
1953 else if (caseL == 2)
1954 {
1955 u[0] = potW[2];
1956 u[1] = potW[6];
1957 u[2] = potW[4];
1958 u[3] = potW[7];
1959
1960 T.mv(u, Tu);
1961
1962 fluxW[11] = Tu[0];
1963 potentialDiffW11 = Tu[0];
1964
1965 u[0] = potNw[2];
1966 u[1] = potNw[6];
1967 u[2] = potNw[4];
1968 u[3] = potNw[7];
1969
1970 T.mv(u, Tu);
1971
1972 fluxNw[11] = Tu[0];
1973 potentialDiffNw11 = Tu[0];
1974 }
1975 else if (caseL == 3)
1976 {
1977 u[0] = potW[2];
1978 u[1] = potW[6];
1979 u[2] = potW[4];
1980 u[3] = potW[3];
1981
1982 T.mv(u, Tu);
1983
1984 fluxW[11] = Tu[0];
1985 potentialDiffW11 = Tu[0];
1986
1987 u[0] = potNw[2];
1988 u[1] = potNw[6];
1989 u[2] = potNw[4];
1990 u[3] = potNw[3];
1991
1992 T.mv(u, Tu);
1993
1994 fluxNw[11] = Tu[0];
1995 potentialDiffNw11 = Tu[0];
1996 }
1997 else
1998 {
1999 u[0] = potW[2];
2000 u[1] = potW[6];
2001 u[2] = potW[0];
2002 u[3] = potW[7];
2003
2004 T.mv(u, Tu);
2005
2006 fluxW[11] = Tu[0];
2007 potentialDiffW11 = Tu[0];
2008
2009 u[0] = potNw[2];
2010 u[1] = potNw[6];
2011 u[2] = potNw[0];
2012 u[3] = potNw[7];
2013
2014 T.mv(u, Tu);
2015
2016 fluxNw[11] = Tu[0];
2017 potentialDiffNw11 = Tu[0];
2018 }
2019 }
2020 }
2021 //store potentials for further calculations (saturation, ...)
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);
2028
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);
2035
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);
2042
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);
2049
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);
2056
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);
2063
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);
2070
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);
2077
2078 // compute mobilities of subvolumeface 1 (subVolumeFaceIdx = 0)
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];
2082
2083 // compute mobilities of subvolumeface 2 (subVolumeFaceIdx = 1)
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];
2087
2088 // compute mobilities of subvolumeface 3 (subVolumeFaceIdx = 2)
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];
2092
2093 // compute mobilities of subvolumeface 4 (subVolumeFaceIdx = 3)
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];
2097
2098 // compute mobilities of subvolumeface 5 (subVolumeFaceIdx = 4)
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];
2102
2103 // compute mobilities of subvolumeface 6 (subVolumeFaceIdx = 5)
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];
2107
2108 // compute mobilities of subvolumeface 7 (subVolumeFaceIdx = 6)
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];
2112
2113 // compute mobilities of subvolumeface 8 (subVolumeFaceIdx = 7)
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];
2117
2118 // compute mobilities of subvolumeface 9 (subVolumeFaceIdx = 8)
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];
2122
2123 // compute mobilities of subvolumeface 10 (subVolumeFaceIdx = 9)
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];
2127
2128 // compute mobilities of subvolumeface 11 (subVolumeFaceIdx = 10)
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];
2132
2133 // compute mobilities of subvolumeface 12 (subVolumeFaceIdx = 11)
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];
2137
2138 for (int i = 0; i < numPhases; i++)
2139 {
2140 // evaluate parts of velocity --> always take the normal for which the flux is calculated!
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);
2165 vel21 *= -1;
2166 vel42 *= -1;
2167 vel34 *= -1;
2168 vel13 *= -1;
2169 vel56 *= -1;
2170 vel68 *= -1;
2171 vel87 *= -1;
2172 vel75 *= -1;
2173 vel15 *= -1;
2174 vel62 *= -1;
2175 vel48 *= -1;
2176 vel73 *= -1;
2177
2178 Dune::FieldVector<Scalar, 12> flux(0);
2179 switch (i)
2180 {
2181 case wPhaseIdx:
2182 {
2183 flux = fluxW;
2184 break;
2185 }
2186 case nPhaseIdx:
2187 {
2188 flux = fluxNw;
2189 break;
2190 }
2191 }
2192
2193 if (hangingNodeType == InteractionVolume::sixSmallCells)
2194 {
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));
2207 vel78 = 0.;
2208 vel87 = 0.;
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));
2219 }
2220 else if (hangingNodeType == InteractionVolume::fourSmallCellsDiag)
2221 {
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));
2226 vel43 = 0.;
2227 vel34 = 0.;
2228 vel31 *= flux[3] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx3, 2, 0));
2229 vel13 *= flux[3] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx1, 0, 1));
2230 vel65 = 0.;
2231 vel56 = 0.;
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));
2246 }
2247 else if (hangingNodeType == InteractionVolume::fourSmallCellsFace
2248 || hangingNodeType == InteractionVolume::fourSmallCellsEdge)
2249 {
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));
2266
2267 if (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx7)
2268 {
2269 vel65 = 0.;
2270 vel56 = 0.;
2271 vel86 *= flux[5] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx8, 7, 2));
2272 vel68 *= flux[5] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx6, 5, 1));
2273 vel78 = 0.;
2274 vel87 = 0.;
2275 vel57 *= flux[7] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx5, 4, 2));
2276 vel75 *= flux[7] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx7, 6, 1));
2277 }
2278 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx6)
2279 {
2280 vel65 *= flux[4] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx6, 5, 2));
2281 vel56 *= flux[4] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx5, 4, 1));
2282 vel86 = 0.;
2283 vel68 = 0.;
2284 vel78 *= flux[6] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx7, 6, 2));
2285 vel87 *= flux[6] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx8, 7, 1));
2286 vel57 = 0.;
2287 vel75 = 0.;
2288 }
2289 else
2290 {
2291 vel65 = 0.;
2292 vel56 = 0.;
2293 vel86 = 0.;
2294 vel68 = 0.;
2295 vel78 = 0.;
2296 vel87 = 0.;
2297 vel57 = 0.;
2298 vel75 = 0.;
2299 }
2300 }
2301 else if (hangingNodeType == InteractionVolume::twoSmallCells)
2302 {
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));
2307 vel43 = 0.;
2308 vel34 = 0.;
2309 vel31 *= flux[3] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx3, 2, 0));
2310 vel13 *= flux[3] / (interactionVolumes_().getRealFluxFaceArea(interactionVolume, globalIdx1, 0, 1));
2311 vel65 = 0.;
2312 vel56 = 0.;
2313 vel86 = 0.;
2314 vel68 = 0.;
2315 vel78 = 0.;
2316 vel87 = 0.;
2317 vel57 = 0.;
2318 vel75 = 0.;
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));
2323 vel84 = 0.;
2324 vel48 = 0.;
2325 vel37 = 0.;
2326 vel73 = 0.;
2327 }
2328
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;
2353
2354 vel12 *= fracFlow0;
2355 vel21 *= fracFlow0;
2356 vel24 *= fracFlow1;
2357 vel42 *= fracFlow1;
2358 vel43 *= fracFlow2;
2359 vel34 *= fracFlow2;
2360 vel31 *= fracFlow3;
2361 vel13 *= fracFlow3;
2362 vel65 *= fracFlow4;
2363 vel56 *= fracFlow4;
2364 vel86 *= fracFlow5;
2365 vel68 *= fracFlow5;
2366 vel78 *= fracFlow6;
2367 vel87 *= fracFlow6;
2368 vel57 *= fracFlow7;
2369 vel75 *= fracFlow7;
2370 vel51 *= fracFlow8;
2371 vel15 *= fracFlow8;
2372 vel26 *= fracFlow9;
2373 vel62 *= fracFlow9;
2374 vel84 *= fracFlow10;
2375 vel48 *= fracFlow10;
2376 vel37 *= fracFlow11;
2377 vel73 *= fracFlow11;
2378
2379 //store velocities
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);
2404 }
2405 //set velocity marker
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));
2430 }
2431
2432} // end namespace Dumux
2433#endif
2-d velocity calculation using a 3-d MPFA L-method.
Definition: adapt.hh:29
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