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