27#ifndef DUMUX_FVMPFAL2PFABOUND3DPRESSURE2P_ADAPTIVE_HH
28#define DUMUX_FVMPFAL2PFABOUND3DPRESSURE2P_ADAPTIVE_HH
76template<
class TypeTag>
85 dim = GridView::dimension, dimWorld = GridView::dimensionworld
102 using PrimaryVariables =
typename SolutionTypes::PrimaryVariables;
103 using ScalarSolutionType =
typename SolutionTypes::ScalarSolution;
107 pw = Indices::pressureW,
108 pn = Indices::pressureNw,
109 pGlobal = Indices::pressureGlobal,
110 sw = Indices::saturationW,
111 sn = Indices::saturationNw,
112 vw = Indices::velocityW,
113 vn = Indices::velocityNw,
114 vt = Indices::velocityTotal
118 wPhaseIdx = Indices::wPhaseIdx,
119 nPhaseIdx = Indices::nPhaseIdx,
120 pressureIdx = Indices::pressureIdx,
121 saturationIdx = Indices::saturationIdx,
122 pressureEqIdx = Indices::pressureEqIdx,
123 satEqIdx = Indices::satEqIdx,
124 numPhases = getPropValue<TypeTag, Properties::NumPhases>()
131 dirichletDirichlet = 1,
132 dirichletNeumann = 2,
141 using Element =
typename GridView::Traits::template Codim<0>::Entity;
142 using Grid =
typename GridView::Grid;
143 using Geometry =
typename Element::Geometry;
145 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
146 using GravityVector = Dune::FieldVector<Scalar, dimWorld>;
147 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
149 using DimVector = Dune::FieldVector<Scalar, dim>;
178 const auto element = *problem_.gridView().template begin<0>();
179 FluidState fluidState;
180 fluidState.setPressure(wPhaseIdx, problem_.referencePressure(element));
181 fluidState.setPressure(nPhaseIdx, problem_.referencePressure(element));
182 fluidState.setTemperature(problem_.temperature(element));
183 fluidState.setSaturation(wPhaseIdx, 1.);
184 fluidState.setSaturation(nPhaseIdx, 0.);
196 int size = problem_.gridView().size(0);
198 if (problem_.gridAdapt().wasAdapted())
201 this->
A_.setSize(size, size);
202 this->
f_.resize(size);
207 asImp_().initializeMatrix();
219 gravity_(problem.gravity())
221 if (pressureType_ != pw && pressureType_ != pn)
223 DUNE_THROW(Dune::NotImplemented,
"Pressure type not supported!");
225 if (saturationType_ != sw && saturationType_ != sn)
227 DUNE_THROW(Dune::NotImplemented,
"Saturation type not supported!");
229 if (getPropValue<TypeTag, Properties::EnableCompressibility>())
231 DUNE_THROW(Dune::NotImplemented,
"Compressibility not supported!");
235 DUNE_THROW(Dune::NotImplemented,
"Dimension not supported!");
238 density_[wPhaseIdx] = 0.;
239 density_[nPhaseIdx] = 0.;
240 viscosity_[wPhaseIdx] = 0.;
241 viscosity_[nPhaseIdx] = 0.;
248 Implementation &asImp_()
249 {
return *
static_cast<Implementation *
>(
this);}
252 const Implementation &asImp_()
const
253 {
return *
static_cast<const Implementation *
>(
this);}
255 const GravityVector& gravity_;
257 Scalar density_[numPhases];
258 Scalar viscosity_[numPhases];
260 static constexpr Scalar threshold_ = 1e-15;
262 static const int pressureType_ = getPropValue<TypeTag, Properties::PressureFormulation>();
264 static const int saturationType_ = getPropValue<TypeTag, Properties::SaturationFormulation>();
266 static const int velocityType_ = getPropValue<TypeTag, Properties::VelocityFormulation>();
270template<
class TypeTag>
273 initializeMatrixRowSize();
274 this->A_.endrowsizes();
275 initializeMatrixIndices();
276 this->A_.endindices();
280template<
class TypeTag>
284 for (
const auto& element : elements(problem_.gridView()))
287 int globalIdxI = problem_.variables().index(element);
289 int levelI = element.level();
291 std::set<int> neighborIndices;
293 int numVertices = element.geometry().corners();
295 for (
int vIdxI = 0; vIdxI < numVertices; vIdxI++)
297 int vIdxIGlobal = problem_.variables().vertexMapper().subIndex(element, vIdxI, dim);
299 InteractionVolume& interactionVolume = this->interactionVolumes_.interactionVolume(vIdxIGlobal);
301 for (
int subVolumeIdx = 0; subVolumeIdx < InteractionVolume::subVolumeTotalNum; subVolumeIdx++)
303 if (interactionVolume.hasSubVolumeElement(subVolumeIdx))
305 auto neighbor = interactionVolume.getSubVolumeElement(subVolumeIdx);
306 int neighborIdx = problem_.variables().index(neighbor);
308 neighborIndices.insert(neighborIdx);
310 if (!interactionVolume.sameLevel())
312 if (neighbor.level() == levelI + 2)
314 for (
int vIdxJ = 0; vIdxJ < numVertices; vIdxJ++)
316 int vIdxJGlobal = problem_.variables().vertexMapper().subIndex(neighbor, vIdxJ, dim);
318 if (vIdxJGlobal != vIdxIGlobal)
321 = this->interactionVolumes_.interactionVolume(vIdxJGlobal);
323 if (interactionVolumeJ.isHangingNodeVolume())
325 std::set<int> neighborIndicesJ;
326 bool additionalEntries =
false;
328 for (
int subVolumeIdxJ = 0;
329 subVolumeIdxJ < InteractionVolume::subVolumeTotalNum; subVolumeIdxJ++)
331 if (interactionVolumeJ.hasSubVolumeElement(subVolumeIdxJ))
333 int globalIdxJJ = problem_.variables().index(interactionVolumeJ.getSubVolumeElement(subVolumeIdxJ));
335 neighborIndicesJ.insert(globalIdxJJ);
337 if (globalIdxI == globalIdxJJ)
339 additionalEntries =
true;
344 if (additionalEntries)
346 neighborIndices.insert(neighborIndicesJ.begin(), neighborIndicesJ.end());
357 this->A_.setrowsize(globalIdxI, neighborIndices.size());
364template<
class TypeTag>
368 for (
const auto& element : elements(problem_.gridView()))
371 int globalIdxI = problem_.variables().index(element);
373 int levelI = element.level();
376 this->A_.addindex(globalIdxI, globalIdxI);
378 int numVertices = element.geometry().corners();
380 for (
int vIdx = 0; vIdx < numVertices; vIdx++)
382 int vIdxGlobal = problem_.variables().vertexMapper().subIndex(element, vIdx, dim);
384 InteractionVolume& interactionVolume = this->interactionVolumes_.interactionVolume(vIdxGlobal);
385 for (
int subVolumeIdx = 0; subVolumeIdx < InteractionVolume::subVolumeTotalNum; subVolumeIdx++)
387 if (interactionVolume.hasSubVolumeElement(subVolumeIdx))
389 auto neighbor = interactionVolume.getSubVolumeElement(subVolumeIdx);
390 int globalIdxJ = problem_.variables().index(neighbor);
392 this->A_.addindex(globalIdxI, globalIdxJ);
394 if (interactionVolume.isHangingNodeVolume() && interactionVolume.getHangingNodeType() ==
395 InteractionVolume::sixSmallCells && !interactionVolume.sameLevel())
397 if (neighbor.level() == levelI-2)
399 this->A_.addindex(globalIdxJ, globalIdxI);
411template<
class TypeTag>
419 for (
const auto& vertex : vertices(problem_.gridView()))
422 if (vertex.partitionType() != Dune::InteriorEntity && vertex.partitionType() != Dune::BorderEntity)
427 int vIdxGlobal = problem_.variables().index(vertex);
429 InteractionVolume& interactionVolume = this->interactionVolumes_.interactionVolume(vIdxGlobal);
432 if (interactionVolume.isInnerVolume())
434 if (!interactionVolume.isHangingNodeVolume())
436 this->assembleInnerInteractionVolume(interactionVolume);
440 assembleHangingNodeInteractionVolume(interactionVolume);
446 this->assembleBoundaryInteractionVolume(interactionVolume);
455template<
class TypeTag>
458 auto element1 = interactionVolume.getSubVolumeElement(0);
459 auto element2 = interactionVolume.getSubVolumeElement(1);
460 auto element3 = interactionVolume.getSubVolumeElement(2);
461 auto element4 = interactionVolume.getSubVolumeElement(3);
462 auto element5 = interactionVolume.getSubVolumeElement(4);
463 auto element6 = interactionVolume.getSubVolumeElement(5);
464 auto element7 = interactionVolume.getSubVolumeElement(6);
465 auto element8 = interactionVolume.getSubVolumeElement(7);
468 const GlobalPosition& globalPos1 = element1.geometry().center();
469 const GlobalPosition& globalPos2 = element2.geometry().center();
470 const GlobalPosition& globalPos3 = element3.geometry().center();
471 const GlobalPosition& globalPos4 = element4.geometry().center();
472 const GlobalPosition& globalPos5 = element5.geometry().center();
473 const GlobalPosition& globalPos6 = element6.geometry().center();
474 const GlobalPosition& globalPos7 = element7.geometry().center();
475 const GlobalPosition& globalPos8 = element8.geometry().center();
478 Scalar volume1 = element1.geometry().volume();
479 Scalar volume2 = element2.geometry().volume();
480 Scalar volume3 = element3.geometry().volume();
481 Scalar volume4 = element4.geometry().volume();
482 Scalar volume5 = element5.geometry().volume();
483 Scalar volume6 = element6.geometry().volume();
484 Scalar volume7 DUNE_UNUSED = element7.geometry().volume();
485 Scalar volume8 DUNE_UNUSED = element8.geometry().volume();
488 int globalIdx1 = problem_.variables().index(element1);
489 int globalIdx2 = problem_.variables().index(element2);
490 int globalIdx3 = problem_.variables().index(element3);
491 int globalIdx4 = problem_.variables().index(element4);
492 int globalIdx5 = problem_.variables().index(element5);
493 int globalIdx6 = problem_.variables().index(element6);
494 int globalIdx7 = problem_.variables().index(element7);
495 int globalIdx8 = problem_.variables().index(element8);
498 CellData& cellData1 = problem_.variables().cellData(globalIdx1);
499 CellData& cellData2 = problem_.variables().cellData(globalIdx2);
500 CellData& cellData3 = problem_.variables().cellData(globalIdx3);
501 CellData& cellData4 = problem_.variables().cellData(globalIdx4);
502 CellData& cellData5 = problem_.variables().cellData(globalIdx5);
503 CellData& cellData6 = problem_.variables().cellData(globalIdx6);
504 CellData& cellData7 = problem_.variables().cellData(globalIdx7);
505 CellData& cellData8 = problem_.variables().cellData(globalIdx8);
508 Dune::FieldVector<Scalar, numPhases> lambda1(cellData1.mobility(wPhaseIdx));
509 lambda1[nPhaseIdx] = cellData1.mobility(nPhaseIdx);
512 Scalar lambdaTotal1 = lambda1[wPhaseIdx] + lambda1[nPhaseIdx];
515 Dune::FieldVector<Scalar, numPhases> lambda2(cellData2.mobility(wPhaseIdx));
516 lambda2[nPhaseIdx] = cellData2.mobility(nPhaseIdx);
519 Scalar lambdaTotal2 = lambda2[wPhaseIdx] + lambda2[nPhaseIdx];
522 Dune::FieldVector<Scalar, numPhases> lambda3(cellData3.mobility(wPhaseIdx));
523 lambda3[nPhaseIdx] = cellData3.mobility(nPhaseIdx);
526 Scalar lambdaTotal3 = lambda3[wPhaseIdx] + lambda3[nPhaseIdx];
529 Dune::FieldVector<Scalar, numPhases> lambda4(cellData4.mobility(wPhaseIdx));
530 lambda4[nPhaseIdx] = cellData4.mobility(nPhaseIdx);
533 Scalar lambdaTotal4 = lambda4[wPhaseIdx] + lambda4[nPhaseIdx];
536 Dune::FieldVector<Scalar, numPhases> lambda5(cellData5.mobility(wPhaseIdx));
537 lambda5[nPhaseIdx] = cellData5.mobility(nPhaseIdx);
540 Scalar lambdaTotal5 = lambda5[wPhaseIdx] + lambda5[nPhaseIdx];
543 Dune::FieldVector<Scalar, numPhases> lambda6(cellData6.mobility(wPhaseIdx));
544 lambda6[nPhaseIdx] = cellData6.mobility(nPhaseIdx);
547 Scalar lambdaTotal6 = lambda6[wPhaseIdx] + lambda6[nPhaseIdx];
550 Dune::FieldVector<Scalar, numPhases> lambda7(cellData7.mobility(wPhaseIdx));
551 lambda7[nPhaseIdx] = cellData7.mobility(nPhaseIdx);
554 Scalar lambdaTotal7 = lambda7[wPhaseIdx] + lambda7[nPhaseIdx];
557 Dune::FieldVector<Scalar, numPhases> lambda8(cellData8.mobility(wPhaseIdx));
558 lambda8[nPhaseIdx] = cellData8.mobility(nPhaseIdx);
561 Scalar lambdaTotal8 = lambda8[wPhaseIdx] + lambda8[nPhaseIdx];
563 std::vector<DimVector> lambda(8);
564 lambda[0][0] = lambdaTotal1;
565 lambda[0][1] = lambdaTotal1;
566 lambda[0][2] = lambdaTotal1;
567 lambda[1][0] = lambdaTotal2;
568 lambda[1][1] = lambdaTotal2;
569 lambda[1][2] = lambdaTotal2;
570 lambda[2][0] = lambdaTotal3;
571 lambda[2][1] = lambdaTotal3;
572 lambda[2][2] = lambdaTotal3;
573 lambda[3][0] = lambdaTotal4;
574 lambda[3][1] = lambdaTotal4;
575 lambda[3][2] = lambdaTotal4;
576 lambda[4][0] = lambdaTotal5;
577 lambda[4][1] = lambdaTotal5;
578 lambda[4][2] = lambdaTotal5;
579 lambda[5][0] = lambdaTotal6;
580 lambda[5][1] = lambdaTotal6;
581 lambda[5][2] = lambdaTotal6;
582 lambda[6][0] = lambdaTotal7;
583 lambda[6][1] = lambdaTotal7;
584 lambda[6][2] = lambdaTotal7;
585 lambda[7][0] = lambdaTotal8;
586 lambda[7][1] = lambdaTotal8;
587 lambda[7][2] = lambdaTotal8;
591 Dune::FieldVector<Scalar, 8> pc(0);
592 pc[0] = cellData1.capillaryPressure();
593 pc[1] = cellData2.capillaryPressure();
594 pc[2] = cellData3.capillaryPressure();
595 pc[3] = cellData4.capillaryPressure();
596 pc[4] = cellData5.capillaryPressure();
597 pc[5] = cellData6.capillaryPressure();
598 pc[6] = cellData7.capillaryPressure();
599 pc[7] = cellData8.capillaryPressure();
601 Dune::FieldVector<Scalar, 8> gravityDiff(0);
605 gravityDiff[0] = (problem_.bBoxMax() - globalPos1) * gravity_ * (density_[nPhaseIdx] - density_[wPhaseIdx]);
606 gravityDiff[1] = (problem_.bBoxMax() - globalPos2) * gravity_ * (density_[nPhaseIdx] - density_[wPhaseIdx]);
607 gravityDiff[2] = (problem_.bBoxMax() - globalPos3) * gravity_ * (density_[nPhaseIdx] - density_[wPhaseIdx]);
608 gravityDiff[3] = (problem_.bBoxMax() - globalPos4) * gravity_ * (density_[nPhaseIdx] - density_[wPhaseIdx]);
609 gravityDiff[4] = (problem_.bBoxMax() - globalPos5) * gravity_ * (density_[nPhaseIdx] - density_[wPhaseIdx]);
610 gravityDiff[5] = (problem_.bBoxMax() - globalPos6) * gravity_ * (density_[nPhaseIdx] - density_[wPhaseIdx]);
611 gravityDiff[6] = (problem_.bBoxMax() - globalPos7) * gravity_ * (density_[nPhaseIdx] - density_[wPhaseIdx]);
612 gravityDiff[7] = (problem_.bBoxMax() - globalPos8) * gravity_ * (density_[nPhaseIdx] - density_[wPhaseIdx]);
616 Dune::FieldVector<Dune::FieldVector<Scalar, 3>, 8> pcFlux(Dune::FieldVector<Scalar, 3>(0));
618 Scalar pcPotential0 = 0;
619 Scalar pcPotential1 = 0;
620 Scalar pcPotential2 = 0;
621 Scalar pcPotential3 = 0;
622 Scalar pcPotential4 = 0;
623 Scalar pcPotential5 = 0;
624 Scalar pcPotential6 = 0;
625 Scalar pcPotential7 = 0;
626 Scalar pcPotential8 = 0;
627 Scalar pcPotential9 = 0;
628 Scalar pcPotential10 = 0;
629 Scalar pcPotential11 = 0;
631 int hangingNodeType = interactionVolume.getHangingNodeType();
633 if (hangingNodeType == InteractionVolume::twoSmallCells)
636 PrimaryVariables source(0.0);
637 problem_.source(source, element1);
638 this->f_[globalIdx1] += volume1 / (8.0)
639 * (source[wPhaseIdx] / density_[wPhaseIdx] + source[nPhaseIdx] / density_[nPhaseIdx]);
640 problem_.source(source, element2);
641 this->f_[globalIdx2] += volume2 / (8.0)
642 * (source[wPhaseIdx] / density_[wPhaseIdx] + source[nPhaseIdx] / density_[nPhaseIdx]);
644 this->f_[globalIdx1] += this->evaluateErrorTerm(cellData1) * volume1 / (8.0);
645 this->f_[globalIdx2] += this->evaluateErrorTerm(cellData2) * volume2 / (8.0);
647 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge
648 || hangingNodeType == InteractionVolume::fourSmallCellsFace)
651 PrimaryVariables source(0.0);
652 problem_.source(source, element1);
653 this->f_[globalIdx1] += volume1 / (8.0)
654 * (source[wPhaseIdx] / density_[wPhaseIdx] + source[nPhaseIdx] / density_[nPhaseIdx]);
655 problem_.source(source, element2);
656 this->f_[globalIdx2] += volume2 / (8.0)
657 * (source[wPhaseIdx] / density_[wPhaseIdx] + source[nPhaseIdx] / density_[nPhaseIdx]);
658 problem_.source(source, element3);
659 this->f_[globalIdx3] += volume3 / (8.0)
660 * (source[wPhaseIdx] / density_[wPhaseIdx] + source[nPhaseIdx] / density_[nPhaseIdx]);
661 problem_.source(source, element4);
662 this->f_[globalIdx4] += volume4 / (8.0)
663 * (source[wPhaseIdx] / density_[wPhaseIdx] + source[nPhaseIdx] / density_[nPhaseIdx]);
665 this->f_[globalIdx1] += this->evaluateErrorTerm(cellData1) * volume1 / (8.0);
666 this->f_[globalIdx2] += this->evaluateErrorTerm(cellData2) * volume2 / (8.0);
667 this->f_[globalIdx3] += this->evaluateErrorTerm(cellData3) * volume3 / (8.0);
668 this->f_[globalIdx4] += this->evaluateErrorTerm(cellData4) * volume4 / (8.0);
670 else if (hangingNodeType == InteractionVolume::sixSmallCells)
673 PrimaryVariables source(0.0);
674 problem_.source(source, element1);
675 this->f_[globalIdx1] += volume1 / (8.0)
676 * (source[wPhaseIdx] / density_[wPhaseIdx] + source[nPhaseIdx] / density_[nPhaseIdx]);
677 problem_.source(source, element2);
678 this->f_[globalIdx2] += volume2 / (8.0)
679 * (source[wPhaseIdx] / density_[wPhaseIdx] + source[nPhaseIdx] / density_[nPhaseIdx]);
680 problem_.source(source, element3);
681 this->f_[globalIdx3] += volume3 / (8.0)
682 * (source[wPhaseIdx] / density_[wPhaseIdx] + source[nPhaseIdx] / density_[nPhaseIdx]);
683 problem_.source(source, element4);
684 this->f_[globalIdx4] += volume4 / (8.0)
685 * (source[wPhaseIdx] / density_[wPhaseIdx] + source[nPhaseIdx] / density_[nPhaseIdx]);
686 problem_.source(source, element4);
687 this->f_[globalIdx5] += volume5 / (8.0)
688 * (source[wPhaseIdx] / density_[wPhaseIdx] + source[nPhaseIdx] / density_[nPhaseIdx]);
689 problem_.source(source, element6);
690 this->f_[globalIdx6] += volume6 / (8.0)
691 * (source[wPhaseIdx] / density_[wPhaseIdx] + source[nPhaseIdx] / density_[nPhaseIdx]);
693 this->f_[globalIdx1] += this->evaluateErrorTerm(cellData1) * volume1 / (8.0);
694 this->f_[globalIdx2] += this->evaluateErrorTerm(cellData2) * volume2 / (8.0);
695 this->f_[globalIdx3] += this->evaluateErrorTerm(cellData3) * volume3 / (8.0);
696 this->f_[globalIdx4] += this->evaluateErrorTerm(cellData4) * volume4 / (8.0);
697 this->f_[globalIdx5] += this->evaluateErrorTerm(cellData5) * volume5 / (8.0);
698 this->f_[globalIdx6] += this->evaluateErrorTerm(cellData6) * volume6 / (8.0);
702 Dune::FieldVector<Scalar, 2 * dim - dim + 1> u(0);
705 Dune::FieldVector<bool, 4> useCases(
false);
708 int caseL = this->transmissibilityCalculator_.transmissibility(T, interactionVolume, lambda, 0, 1, 2, 3, 4,
712 T *= this->interactionVolumes_.faceAreaFactor(interactionVolume, globalIdx1, 0, 0);
713 TSecond *= this->interactionVolumes_.faceAreaFactor(interactionVolume, globalIdx2, 1, 1);
717 this->A_[globalIdx1][globalIdx1] += T[0][0];
718 this->A_[globalIdx1][globalIdx2] += T[0][1];
719 this->A_[globalIdx1][globalIdx3] += T[0][2];
720 this->A_[globalIdx1][globalIdx5] += T[0][3];
722 this->A_[globalIdx2][globalIdx1] -= TSecond[0][0];
723 this->A_[globalIdx2][globalIdx2] -= TSecond[0][1];
724 this->A_[globalIdx2][globalIdx3] -= TSecond[0][2];
725 this->A_[globalIdx2][globalIdx5] -= TSecond[0][3];
734 pcFlux[0][0] = Tu[0];
735 pcPotential0 = Tu[0];
742 this->A_[globalIdx1][globalIdx1] += T[0][0];
743 this->A_[globalIdx1][globalIdx2] += T[0][1];
744 this->A_[globalIdx1][globalIdx4] += T[0][2];
745 this->A_[globalIdx1][globalIdx6] += T[0][3];
747 this->A_[globalIdx2][globalIdx1] -= TSecond[0][0];
748 this->A_[globalIdx2][globalIdx2] -= TSecond[0][1];
749 this->A_[globalIdx2][globalIdx4] -= TSecond[0][2];
750 this->A_[globalIdx2][globalIdx6] -= TSecond[0][3];
758 pcFlux[0][0] = Tu[0];
759 pcPotential0 = Tu[0];
766 this->A_[globalIdx1][globalIdx1] += T[0][0];
767 this->A_[globalIdx1][globalIdx2] += T[0][1];
768 this->A_[globalIdx1][globalIdx4] += T[0][2];
769 this->A_[globalIdx1][globalIdx5] += T[0][3];
771 this->A_[globalIdx2][globalIdx1] -= TSecond[0][0];
772 this->A_[globalIdx2][globalIdx2] -= TSecond[0][1];
773 this->A_[globalIdx2][globalIdx4] -= TSecond[0][2];
774 this->A_[globalIdx2][globalIdx5] -= TSecond[0][3];
782 pcFlux[0][0] = Tu[0];
783 pcPotential0 = Tu[0];
790 this->A_[globalIdx1][globalIdx1] += T[0][0];
791 this->A_[globalIdx1][globalIdx2] += T[0][1];
792 this->A_[globalIdx1][globalIdx3] += T[0][2];
793 this->A_[globalIdx1][globalIdx6] += T[0][3];
795 this->A_[globalIdx2][globalIdx1] -= TSecond[0][0];
796 this->A_[globalIdx2][globalIdx2] -= TSecond[0][1];
797 this->A_[globalIdx2][globalIdx3] -= TSecond[0][2];
798 this->A_[globalIdx2][globalIdx6] -= TSecond[0][3];
806 pcFlux[0][0] = Tu[0];
807 pcPotential0 = Tu[0];
816 if (hangingNodeType == InteractionVolume::twoSmallCells
817 || hangingNodeType == InteractionVolume::fourSmallCellsDiag)
823 caseL = this->transmissibilityCalculator_.transmissibility(T, interactionVolume, lambda, 1, 3, 0,
828 caseL = this->transmissibilityCalculator_.transmissibility(T, interactionVolume, lambda, 1, 3, 0,
833 T *= this->interactionVolumes_.faceAreaFactor(interactionVolume, globalIdx2, 1, 0);
834 TSecond *= this->interactionVolumes_.faceAreaFactor(interactionVolume, globalIdx4, 3, 1);
839 this->A_[globalIdx2][globalIdx2] += T[0][0];
840 this->A_[globalIdx2][globalIdx4] += T[0][1];
841 this->A_[globalIdx2][globalIdx1] += T[0][2];
842 this->A_[globalIdx2][globalIdx6] += T[0][3];
844 this->A_[globalIdx4][globalIdx2] -= TSecond[0][0];
845 this->A_[globalIdx4][globalIdx4] -= TSecond[0][1];
846 this->A_[globalIdx4][globalIdx1] -= TSecond[0][2];
847 this->A_[globalIdx4][globalIdx6] -= TSecond[0][3];
855 pcFlux[1][0] = Tu[0];
856 pcPotential1 = Tu[0];
859 pcFlux[3][1] = Tu[0];
864 this->A_[globalIdx2][globalIdx2] += T[0][0];
865 this->A_[globalIdx2][globalIdx4] += T[0][1];
866 this->A_[globalIdx2][globalIdx3] += T[0][2];
867 this->A_[globalIdx2][globalIdx8] += T[0][3];
869 this->A_[globalIdx4][globalIdx2] -= TSecond[0][0];
870 this->A_[globalIdx4][globalIdx4] -= TSecond[0][1];
871 this->A_[globalIdx4][globalIdx3] -= TSecond[0][2];
872 this->A_[globalIdx4][globalIdx8] -= TSecond[0][3];
880 pcFlux[1][0] = Tu[0];
881 pcPotential1 = Tu[0];
884 pcFlux[3][1] = Tu[0];
888 this->A_[globalIdx2][globalIdx2] += T[0][0];
889 this->A_[globalIdx2][globalIdx4] += T[0][1];
890 this->A_[globalIdx2][globalIdx3] += T[0][2];
891 this->A_[globalIdx2][globalIdx6] += T[0][3];
893 this->A_[globalIdx4][globalIdx2] -= TSecond[0][0];
894 this->A_[globalIdx4][globalIdx4] -= TSecond[0][1];
895 this->A_[globalIdx4][globalIdx3] -= TSecond[0][2];
896 this->A_[globalIdx4][globalIdx6] -= TSecond[0][3];
904 pcFlux[1][0] = Tu[0];
905 pcPotential1 = Tu[0];
908 pcFlux[3][1] = Tu[0];
912 this->A_[globalIdx2][globalIdx2] += T[0][0];
913 this->A_[globalIdx2][globalIdx4] += T[0][1];
914 this->A_[globalIdx2][globalIdx1] += T[0][2];
915 this->A_[globalIdx2][globalIdx8] += T[0][3];
917 this->A_[globalIdx4][globalIdx2] -= TSecond[0][0];
918 this->A_[globalIdx4][globalIdx4] -= TSecond[0][1];
919 this->A_[globalIdx4][globalIdx1] -= TSecond[0][2];
920 this->A_[globalIdx4][globalIdx8] -= TSecond[0][3];
928 pcFlux[1][0] = Tu[0];
929 pcPotential1 = Tu[0];
932 pcFlux[3][1] = Tu[0];
938 if (hangingNodeType != InteractionVolume::twoSmallCells
939 && hangingNodeType != InteractionVolume::fourSmallCellsDiag)
941 caseL = this->transmissibilityCalculator_.transmissibility(T, interactionVolume, lambda, 3, 2, 1,
945 T *= this->interactionVolumes_.faceAreaFactor(interactionVolume, globalIdx4, 3, 0);
946 TSecond *= this->interactionVolumes_.faceAreaFactor(interactionVolume, globalIdx3, 2, 1);
952 this->A_[globalIdx4][globalIdx4] += T[0][0];
953 this->A_[globalIdx4][globalIdx3] += T[0][1];
954 this->A_[globalIdx4][globalIdx2] += T[0][2];
955 this->A_[globalIdx4][globalIdx8] += T[0][3];
957 this->A_[globalIdx3][globalIdx4] -= TSecond[0][0];
958 this->A_[globalIdx3][globalIdx3] -= TSecond[0][1];
959 this->A_[globalIdx3][globalIdx2] -= TSecond[0][2];
960 this->A_[globalIdx3][globalIdx8] -= TSecond[0][3];
968 pcPotential2 = Tu[0];
969 pcFlux[3][0] = Tu[0];
972 pcFlux[2][1] = Tu[0];
976 this->A_[globalIdx4][globalIdx4] += T[0][0];
977 this->A_[globalIdx4][globalIdx3] += T[0][1];
978 this->A_[globalIdx4][globalIdx1] += T[0][2];
979 this->A_[globalIdx4][globalIdx7] += T[0][3];
981 this->A_[globalIdx3][globalIdx4] -= TSecond[0][0];
982 this->A_[globalIdx3][globalIdx3] -= TSecond[0][1];
983 this->A_[globalIdx3][globalIdx1] -= TSecond[0][2];
984 this->A_[globalIdx3][globalIdx7] -= TSecond[0][3];
992 pcPotential2 = Tu[0];
993 pcFlux[3][0] = Tu[0];
996 pcFlux[2][1] = Tu[0];
1000 this->A_[globalIdx4][globalIdx4] += T[0][0];
1001 this->A_[globalIdx4][globalIdx3] += T[0][1];
1002 this->A_[globalIdx4][globalIdx1] += T[0][2];
1003 this->A_[globalIdx4][globalIdx8] += T[0][3];
1005 this->A_[globalIdx3][globalIdx4] -= TSecond[0][0];
1006 this->A_[globalIdx3][globalIdx3] -= TSecond[0][1];
1007 this->A_[globalIdx3][globalIdx1] -= TSecond[0][2];
1008 this->A_[globalIdx3][globalIdx8] -= TSecond[0][3];
1016 pcPotential2 = Tu[0];
1017 pcFlux[3][0] = Tu[0];
1020 pcFlux[2][1] = Tu[0];
1024 this->A_[globalIdx4][globalIdx4] += T[0][0];
1025 this->A_[globalIdx4][globalIdx3] += T[0][1];
1026 this->A_[globalIdx4][globalIdx2] += T[0][2];
1027 this->A_[globalIdx4][globalIdx7] += T[0][3];
1029 this->A_[globalIdx3][globalIdx4] -= TSecond[0][0];
1030 this->A_[globalIdx3][globalIdx3] -= TSecond[0][1];
1031 this->A_[globalIdx3][globalIdx2] -= TSecond[0][2];
1032 this->A_[globalIdx3][globalIdx7] -= TSecond[0][3];
1040 pcPotential2 = Tu[0];
1041 pcFlux[3][0] = Tu[0];
1044 pcFlux[2][1] = Tu[0];
1051 if (hangingNodeType == InteractionVolume::twoSmallCells
1052 || hangingNodeType == InteractionVolume::fourSmallCellsDiag)
1054 useCases[0] =
false;
1057 useCases[3] =
false;
1058 caseL = this->transmissibilityCalculator_.transmissibility(T, interactionVolume, lambda, 2, 0, 3, 1, 6,
1063 caseL = this->transmissibilityCalculator_.transmissibility(T, interactionVolume, lambda, 2, 0, 3,
1068 T *= this->interactionVolumes_.faceAreaFactor(interactionVolume, globalIdx3, 2, 0);
1069 TSecond *= this->interactionVolumes_.faceAreaFactor(interactionVolume, globalIdx1, 0, 1);
1075 this->A_[globalIdx3][globalIdx3] += T[0][0];
1076 this->A_[globalIdx3][globalIdx1] += T[0][1];
1077 this->A_[globalIdx3][globalIdx4] += T[0][2];
1078 this->A_[globalIdx3][globalIdx7] += T[0][3];
1080 this->A_[globalIdx1][globalIdx3] -= TSecond[0][0];
1081 this->A_[globalIdx1][globalIdx1] -= TSecond[0][1];
1082 this->A_[globalIdx1][globalIdx4] -= TSecond[0][2];
1083 this->A_[globalIdx1][globalIdx7] -= TSecond[0][3];
1091 pcFlux[2][0] = Tu[0];
1092 pcPotential3 = Tu[0];
1095 pcFlux[0][1] = Tu[0];
1097 else if (caseL == 2)
1099 this->A_[globalIdx3][globalIdx3] += T[0][0];
1100 this->A_[globalIdx3][globalIdx1] += T[0][1];
1101 this->A_[globalIdx3][globalIdx2] += T[0][2];
1102 this->A_[globalIdx3][globalIdx5] += T[0][3];
1104 this->A_[globalIdx1][globalIdx3] -= TSecond[0][0];
1105 this->A_[globalIdx1][globalIdx1] -= TSecond[0][1];
1106 this->A_[globalIdx1][globalIdx2] -= TSecond[0][2];
1107 this->A_[globalIdx1][globalIdx5] -= TSecond[0][3];
1115 pcFlux[2][0] = Tu[0];
1116 pcPotential3 = Tu[0];
1119 pcFlux[0][1] = Tu[0];
1121 else if (caseL == 3)
1123 this->A_[globalIdx3][globalIdx3] += T[0][0];
1124 this->A_[globalIdx3][globalIdx1] += T[0][1];
1125 this->A_[globalIdx3][globalIdx2] += T[0][2];
1126 this->A_[globalIdx3][globalIdx7] += T[0][3];
1128 this->A_[globalIdx1][globalIdx3] -= TSecond[0][0];
1129 this->A_[globalIdx1][globalIdx1] -= TSecond[0][1];
1130 this->A_[globalIdx1][globalIdx2] -= TSecond[0][2];
1131 this->A_[globalIdx1][globalIdx7] -= TSecond[0][3];
1139 pcFlux[2][0] = Tu[0];
1140 pcPotential3 = Tu[0];
1143 pcFlux[0][1] = Tu[0];
1147 this->A_[globalIdx3][globalIdx3] += T[0][0];
1148 this->A_[globalIdx3][globalIdx1] += T[0][1];
1149 this->A_[globalIdx3][globalIdx4] += T[0][2];
1150 this->A_[globalIdx3][globalIdx5] += T[0][3];
1152 this->A_[globalIdx1][globalIdx3] -= TSecond[0][0];
1153 this->A_[globalIdx1][globalIdx1] -= TSecond[0][1];
1154 this->A_[globalIdx1][globalIdx4] -= TSecond[0][2];
1155 this->A_[globalIdx1][globalIdx5] -= TSecond[0][3];
1163 pcFlux[2][0] = Tu[0];
1164 pcPotential3 = Tu[0];
1167 pcFlux[0][1] = Tu[0];
1174 if (hangingNodeType == InteractionVolume::sixSmallCells)
1176 caseL = this->transmissibilityCalculator_.transmissibility(T, interactionVolume, lambda, 5, 4, 7,
1179 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx6)
1181 caseL = this->transmissibilityCalculator_.transmissibilityCaseThree(T, interactionVolume, lambda, 4, 0, 2,5);
1183 int caseLSecond = this->transmissibilityCalculator_.transmissibilityCaseFour(TSecond, interactionVolume, lambda, 1, 5, 3,
1186 caseL = this->transmissibilityCalculator_.chooseTransmissibility(T, TSecond, 3, 4);
1188 if (caseL == caseLSecond)
1193 if (hangingNodeType == InteractionVolume::sixSmallCells)
1196 T *= this->interactionVolumes_.faceAreaFactor(interactionVolume, globalIdx6, 5, 2);
1197 TSecond *= this->interactionVolumes_.faceAreaFactor(interactionVolume, globalIdx5, 4, 1);
1203 this->A_[globalIdx6][globalIdx6] += T[0][0];
1204 this->A_[globalIdx6][globalIdx5] += T[0][1];
1205 this->A_[globalIdx6][globalIdx8] += T[0][2];
1206 this->A_[globalIdx6][globalIdx2] += T[0][3];
1208 this->A_[globalIdx5][globalIdx6] -= TSecond[0][0];
1209 this->A_[globalIdx5][globalIdx5] -= TSecond[0][1];
1210 this->A_[globalIdx5][globalIdx8] -= TSecond[0][2];
1211 this->A_[globalIdx5][globalIdx2] -= TSecond[0][3];
1219 pcFlux[5][2] = Tu[0];
1220 pcPotential4 = Tu[0];
1223 pcFlux[4][1] = Tu[0];
1225 else if (caseL == 2)
1227 this->A_[globalIdx6][globalIdx6] += T[0][0];
1228 this->A_[globalIdx6][globalIdx5] += T[0][1];
1229 this->A_[globalIdx6][globalIdx7] += T[0][2];
1230 this->A_[globalIdx6][globalIdx1] += T[0][3];
1232 this->A_[globalIdx5][globalIdx6] -= TSecond[0][0];
1233 this->A_[globalIdx5][globalIdx5] -= TSecond[0][1];
1234 this->A_[globalIdx5][globalIdx7] -= TSecond[0][2];
1235 this->A_[globalIdx5][globalIdx1] -= TSecond[0][3];
1243 pcFlux[5][2] = Tu[0];
1244 pcPotential4 = Tu[0];
1247 pcFlux[4][1] = Tu[0];
1249 else if (caseL == 3)
1251 this->A_[globalIdx6][globalIdx6] += T[0][0];
1252 this->A_[globalIdx6][globalIdx5] += T[0][1];
1253 this->A_[globalIdx6][globalIdx7] += T[0][2];
1254 this->A_[globalIdx6][globalIdx2] += T[0][3];
1256 this->A_[globalIdx5][globalIdx6] -= TSecond[0][0];
1257 this->A_[globalIdx5][globalIdx5] -= TSecond[0][1];
1258 this->A_[globalIdx5][globalIdx7] -= TSecond[0][2];
1259 this->A_[globalIdx5][globalIdx2] -= TSecond[0][3];
1267 pcFlux[5][2] = Tu[0];
1268 pcPotential4 = Tu[0];
1271 pcFlux[4][1] = Tu[0];
1275 this->A_[globalIdx6][globalIdx6] += T[0][0];
1276 this->A_[globalIdx6][globalIdx5] += T[0][1];
1277 this->A_[globalIdx6][globalIdx8] += T[0][2];
1278 this->A_[globalIdx6][globalIdx1] += T[0][3];
1280 this->A_[globalIdx5][globalIdx6] -= TSecond[0][0];
1281 this->A_[globalIdx5][globalIdx5] -= TSecond[0][1];
1282 this->A_[globalIdx5][globalIdx8] -= TSecond[0][2];
1283 this->A_[globalIdx5][globalIdx1] -= TSecond[0][3];
1291 pcFlux[5][2] = Tu[0];
1292 pcPotential4 = Tu[0];
1295 pcFlux[4][1] = Tu[0];
1298 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx6)
1301 T *= this->interactionVolumes_.faceAreaFactor(interactionVolume, globalIdx6, 5, 2);
1302 TSecond *= this->interactionVolumes_.faceAreaFactor(interactionVolume, globalIdx5, 4, 1);
1306 this->A_[globalIdx6][globalIdx5] -= T[2][0];
1307 this->A_[globalIdx6][globalIdx1] -= T[2][1];
1308 this->A_[globalIdx6][globalIdx3] -= T[2][2];
1309 this->A_[globalIdx6][globalIdx6] -= T[2][3];
1311 this->A_[globalIdx5][globalIdx5] += TSecond[2][0];
1312 this->A_[globalIdx5][globalIdx1] += TSecond[2][1];
1313 this->A_[globalIdx5][globalIdx3] += TSecond[2][2];
1314 this->A_[globalIdx5][globalIdx6] += TSecond[2][3];
1322 pcFlux[5][2] = -Tu[2];
1323 pcPotential4 = -Tu[2];
1326 pcFlux[4][1] = -Tu[2];
1330 this->A_[globalIdx6][globalIdx2] += T[2][0];
1331 this->A_[globalIdx6][globalIdx6] += T[2][1];
1332 this->A_[globalIdx6][globalIdx4] += T[2][2];
1333 this->A_[globalIdx6][globalIdx5] += T[2][3];
1335 this->A_[globalIdx5][globalIdx2] -= TSecond[2][0];
1336 this->A_[globalIdx5][globalIdx6] -= TSecond[2][1];
1337 this->A_[globalIdx5][globalIdx4] -= TSecond[2][2];
1338 this->A_[globalIdx5][globalIdx5] -= TSecond[2][3];
1346 pcFlux[5][2] = Tu[2];
1347 pcPotential4 = Tu[2];
1350 pcFlux[4][1] = Tu[2];
1358 if (hangingNodeType == InteractionVolume::sixSmallCells)
1360 useCases[0] =
false;
1363 useCases[3] =
false;
1365 caseL = this->transmissibilityCalculator_.transmissibility(T, interactionVolume, lambda, 7, 5, 6, 4, 3,
1368 else if (hangingNodeType == InteractionVolume::fourSmallCellsDiag)
1371 useCases[1] =
false;
1372 useCases[2] =
false;
1374 caseL = this->transmissibilityCalculator_.transmissibility(T, interactionVolume, lambda, 7, 5, 6, 4, 3,
1377 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx7)
1379 caseL = this->transmissibilityCalculator_.transmissibilityCaseThree(T, interactionVolume, lambda, 1, 5, 7, 0);
1381 int caseLSecond = this->transmissibilityCalculator_.transmissibilityCaseFour(TSecond, interactionVolume, lambda, 7, 3, 5, 2);
1383 caseL = this->transmissibilityCalculator_.chooseTransmissibility(T, TSecond, 3, 4);
1385 if (caseL == caseLSecond)
1389 if (hangingNodeType == InteractionVolume::sixSmallCells
1390 || hangingNodeType == InteractionVolume::fourSmallCellsDiag)
1393 T *= this->interactionVolumes_.faceAreaFactor(interactionVolume, globalIdx8, 7, 2);
1394 TSecond *= this->interactionVolumes_.faceAreaFactor(interactionVolume, globalIdx6, 5, 1);
1400 this->A_[globalIdx8][globalIdx8] += T[0][0];
1401 this->A_[globalIdx8][globalIdx6] += T[0][1];
1402 this->A_[globalIdx8][globalIdx7] += T[0][2];
1403 this->A_[globalIdx8][globalIdx4] += T[0][3];
1405 this->A_[globalIdx6][globalIdx8] -= TSecond[0][0];
1406 this->A_[globalIdx6][globalIdx6] -= TSecond[0][1];
1407 this->A_[globalIdx6][globalIdx7] -= TSecond[0][2];
1408 this->A_[globalIdx6][globalIdx4] -= TSecond[0][3];
1416 pcFlux[7][2] = Tu[0];
1417 pcPotential5 = Tu[0];
1420 pcFlux[5][1] = Tu[0];
1422 else if (caseL == 2)
1424 this->A_[globalIdx8][globalIdx8] += T[0][0];
1425 this->A_[globalIdx8][globalIdx6] += T[0][1];
1426 this->A_[globalIdx8][globalIdx5] += T[0][2];
1427 this->A_[globalIdx8][globalIdx2] += T[0][3];
1429 this->A_[globalIdx6][globalIdx8] -= TSecond[0][0];
1430 this->A_[globalIdx6][globalIdx6] -= TSecond[0][1];
1431 this->A_[globalIdx6][globalIdx5] -= TSecond[0][2];
1432 this->A_[globalIdx6][globalIdx2] -= TSecond[0][3];
1440 pcFlux[7][2] = Tu[0];
1441 pcPotential5 = Tu[0];
1444 pcFlux[5][1] = Tu[0];
1446 else if (caseL == 3)
1448 this->A_[globalIdx8][globalIdx8] += T[0][0];
1449 this->A_[globalIdx8][globalIdx6] += T[0][1];
1450 this->A_[globalIdx8][globalIdx5] += T[0][2];
1451 this->A_[globalIdx8][globalIdx4] += T[0][3];
1453 this->A_[globalIdx6][globalIdx8] -= TSecond[0][0];
1454 this->A_[globalIdx6][globalIdx6] -= TSecond[0][1];
1455 this->A_[globalIdx6][globalIdx5] -= TSecond[0][2];
1456 this->A_[globalIdx6][globalIdx4] -= TSecond[0][3];
1464 pcFlux[7][2] = Tu[0];
1465 pcPotential5 = Tu[0];
1468 pcFlux[5][1] = Tu[0];
1472 this->A_[globalIdx8][globalIdx8] += T[0][0];
1473 this->A_[globalIdx8][globalIdx6] += T[0][1];
1474 this->A_[globalIdx8][globalIdx7] += T[0][2];
1475 this->A_[globalIdx8][globalIdx2] += T[0][3];
1477 this->A_[globalIdx6][globalIdx8] -= TSecond[0][0];
1478 this->A_[globalIdx6][globalIdx6] -= TSecond[0][1];
1479 this->A_[globalIdx6][globalIdx7] -= TSecond[0][2];
1480 this->A_[globalIdx6][globalIdx2] -= TSecond[0][3];
1488 pcFlux[7][2] = Tu[0];
1489 pcPotential5 = Tu[0];
1492 pcFlux[5][1] = Tu[0];
1495 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx7)
1498 T *= this->interactionVolumes_.faceAreaFactor(interactionVolume, globalIdx8, 7, 2);
1499 TSecond *= this->interactionVolumes_.faceAreaFactor(interactionVolume, globalIdx6, 5, 1);
1503 this->A_[globalIdx8][globalIdx2] -= T[1][0];
1504 this->A_[globalIdx8][globalIdx6] -= T[1][1];
1505 this->A_[globalIdx8][globalIdx8] -= T[1][2];
1506 this->A_[globalIdx8][globalIdx1] -= T[1][3];
1508 this->A_[globalIdx6][globalIdx2] += TSecond[1][0];
1509 this->A_[globalIdx6][globalIdx6] += TSecond[1][1];
1510 this->A_[globalIdx6][globalIdx8] += TSecond[1][2];
1511 this->A_[globalIdx6][globalIdx1] += TSecond[1][3];
1519 pcFlux[7][2] = -Tu[1];
1520 pcPotential5 = -Tu[1];
1523 pcFlux[5][1] = -Tu[1];
1527 this->A_[globalIdx8][globalIdx8] += T[1][0];
1528 this->A_[globalIdx8][globalIdx4] += T[1][1];
1529 this->A_[globalIdx8][globalIdx6] += T[1][2];
1530 this->A_[globalIdx8][globalIdx3] += T[1][3];
1532 this->A_[globalIdx6][globalIdx8] -= TSecond[1][0];
1533 this->A_[globalIdx6][globalIdx4] -= TSecond[1][1];
1534 this->A_[globalIdx6][globalIdx6] -= TSecond[1][2];
1535 this->A_[globalIdx6][globalIdx3] -= TSecond[1][3];
1543 pcFlux[7][2] = Tu[1];
1544 pcPotential5 = Tu[1];
1547 pcFlux[5][1] = Tu[1];
1555 if (hangingNodeType == InteractionVolume::fourSmallCellsDiag)
1557 caseL = this->transmissibilityCalculator_.transmissibility(T, interactionVolume, lambda, 6, 7, 4,
1560 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx6)
1562 caseL = this->transmissibilityCalculator_.transmissibilityCaseThree(T, interactionVolume, lambda, 7, 3, 1, 6);
1564 int caseLSecond = this->transmissibilityCalculator_.transmissibilityCaseFour(TSecond, interactionVolume, lambda, 2, 6, 0,
1567 caseL = this->transmissibilityCalculator_.chooseTransmissibility(T, TSecond, 3, 4);
1569 if (caseL == caseLSecond)
1573 if (hangingNodeType == InteractionVolume::fourSmallCellsDiag)
1576 T *= this->interactionVolumes_.faceAreaFactor(interactionVolume, globalIdx7, 6, 2);
1577 TSecond *= this->interactionVolumes_.faceAreaFactor(interactionVolume, globalIdx8, 7, 1);
1583 this->A_[globalIdx7][globalIdx7] += T[0][0];
1584 this->A_[globalIdx7][globalIdx8] += T[0][1];
1585 this->A_[globalIdx7][globalIdx5] += T[0][2];
1586 this->A_[globalIdx7][globalIdx3] += T[0][3];
1588 this->A_[globalIdx8][globalIdx7] -= TSecond[0][0];
1589 this->A_[globalIdx8][globalIdx8] -= TSecond[0][1];
1590 this->A_[globalIdx8][globalIdx5] -= TSecond[0][2];
1591 this->A_[globalIdx8][globalIdx3] -= TSecond[0][3];
1599 pcFlux[6][2] = Tu[0];
1600 pcPotential6 = Tu[0];
1603 pcFlux[7][1] = Tu[0];
1605 else if (caseL == 2)
1607 this->A_[globalIdx7][globalIdx7] += T[0][0];
1608 this->A_[globalIdx7][globalIdx8] += T[0][1];
1609 this->A_[globalIdx7][globalIdx6] += T[0][2];
1610 this->A_[globalIdx7][globalIdx4] += T[0][3];
1612 this->A_[globalIdx8][globalIdx7] -= TSecond[0][0];
1613 this->A_[globalIdx8][globalIdx8] -= TSecond[0][1];
1614 this->A_[globalIdx8][globalIdx6] -= TSecond[0][2];
1615 this->A_[globalIdx8][globalIdx4] -= TSecond[0][3];
1623 pcFlux[6][2] = Tu[0];
1624 pcPotential6 = Tu[0];
1627 pcFlux[7][1] = Tu[0];
1629 else if (caseL == 3)
1631 this->A_[globalIdx7][globalIdx7] += T[0][0];
1632 this->A_[globalIdx7][globalIdx8] += T[0][1];
1633 this->A_[globalIdx7][globalIdx6] += T[0][2];
1634 this->A_[globalIdx7][globalIdx3] += T[0][3];
1636 this->A_[globalIdx8][globalIdx7] -= TSecond[0][0];
1637 this->A_[globalIdx8][globalIdx8] -= TSecond[0][1];
1638 this->A_[globalIdx8][globalIdx6] -= TSecond[0][2];
1639 this->A_[globalIdx8][globalIdx3] -= TSecond[0][3];
1647 pcFlux[6][2] = Tu[0];
1648 pcPotential6 = Tu[0];
1651 pcFlux[7][1] = Tu[0];
1655 this->A_[globalIdx7][globalIdx7] += T[0][0];
1656 this->A_[globalIdx7][globalIdx8] += T[0][1];
1657 this->A_[globalIdx7][globalIdx5] += T[0][2];
1658 this->A_[globalIdx7][globalIdx4] += T[0][3];
1660 this->A_[globalIdx8][globalIdx7] -= TSecond[0][0];
1661 this->A_[globalIdx8][globalIdx8] -= TSecond[0][1];
1662 this->A_[globalIdx8][globalIdx5] -= TSecond[0][2];
1663 this->A_[globalIdx8][globalIdx4] -= TSecond[0][3];
1671 pcFlux[6][2] = Tu[0];
1672 pcPotential6 = Tu[0];
1675 pcFlux[7][1] = Tu[0];
1678 else if(hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx6)
1681 T *= this->interactionVolumes_.faceAreaFactor(interactionVolume, globalIdx7, 6, 2);
1682 TSecond *= this->interactionVolumes_.faceAreaFactor(interactionVolume, globalIdx8, 7, 1);
1686 this->A_[globalIdx7][globalIdx8] -= T[2][0];
1687 this->A_[globalIdx7][globalIdx4] -= T[2][1];
1688 this->A_[globalIdx7][globalIdx2] -= T[2][2];
1689 this->A_[globalIdx7][globalIdx7] -= T[2][3];
1691 this->A_[globalIdx8][globalIdx8] += TSecond[2][0];
1692 this->A_[globalIdx8][globalIdx4] += TSecond[2][1];
1693 this->A_[globalIdx8][globalIdx2] += TSecond[2][2];
1694 this->A_[globalIdx8][globalIdx7] += TSecond[2][3];
1702 pcFlux[6][2] = -Tu[2];
1703 pcPotential6 = -Tu[2];
1706 pcFlux[7][1] = -Tu[2];
1708 else if (caseL == 4)
1710 this->A_[globalIdx7][globalIdx3] += T[2][0];
1711 this->A_[globalIdx7][globalIdx7] += T[2][1];
1712 this->A_[globalIdx7][globalIdx1] += T[2][2];
1713 this->A_[globalIdx7][globalIdx8] += T[2][3];
1715 this->A_[globalIdx8][globalIdx3] -= TSecond[2][0];
1716 this->A_[globalIdx8][globalIdx7] -= TSecond[2][1];
1717 this->A_[globalIdx8][globalIdx1] -= TSecond[2][2];
1718 this->A_[globalIdx8][globalIdx8] -= TSecond[2][3];
1726 pcFlux[6][2] = Tu[2];
1727 pcPotential6 = Tu[2];
1730 pcFlux[7][1] = Tu[2];
1738 if (hangingNodeType == InteractionVolume::sixSmallCells)
1741 useCases[1] =
false;
1742 useCases[2] =
false;
1745 caseL = this->transmissibilityCalculator_.transmissibility(T, interactionVolume, lambda, 4, 6, 5, 7, 0,
1748 else if (hangingNodeType == InteractionVolume::fourSmallCellsDiag)
1750 useCases[0] =
false;
1753 useCases[3] =
false;
1754 caseL = this->transmissibilityCalculator_.transmissibility(T, interactionVolume, lambda, 4, 6, 5, 7, 0,
1757 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx7)
1759 caseL = this->transmissibilityCalculator_.transmissibilityCaseThree(T, interactionVolume, lambda, 2, 6, 4,
1762 int caseLSecond = this->transmissibilityCalculator_.transmissibilityCaseFour(TSecond, interactionVolume, lambda, 4, 0, 6,
1765 caseL = this->transmissibilityCalculator_.chooseTransmissibility(T, TSecond, 3, 4);
1767 if (caseL == caseLSecond)
1771 if (hangingNodeType == InteractionVolume::sixSmallCells
1772 || hangingNodeType == InteractionVolume::fourSmallCellsDiag)
1775 T *= this->interactionVolumes_.faceAreaFactor(interactionVolume, globalIdx5, 4, 2);
1776 TSecond *= this->interactionVolumes_.faceAreaFactor(interactionVolume, globalIdx7, 6, 1);
1782 this->A_[globalIdx5][globalIdx5] += T[0][0];
1783 this->A_[globalIdx5][globalIdx7] += T[0][1];
1784 this->A_[globalIdx5][globalIdx6] += T[0][2];
1785 this->A_[globalIdx5][globalIdx1] += T[0][3];
1787 this->A_[globalIdx7][globalIdx5] -= TSecond[0][0];
1788 this->A_[globalIdx7][globalIdx7] -= TSecond[0][1];
1789 this->A_[globalIdx7][globalIdx6] -= TSecond[0][2];
1790 this->A_[globalIdx7][globalIdx1] -= TSecond[0][3];
1799 pcFlux[4][2] = Tu[0];
1800 pcPotential7 = Tu[0];
1803 pcFlux[6][1] = Tu[0];
1805 else if (caseL == 2)
1807 this->A_[globalIdx5][globalIdx5] += T[0][0];
1808 this->A_[globalIdx5][globalIdx7] += T[0][1];
1809 this->A_[globalIdx5][globalIdx8] += T[0][2];
1810 this->A_[globalIdx5][globalIdx3] += T[0][3];
1812 this->A_[globalIdx7][globalIdx5] -= TSecond[0][0];
1813 this->A_[globalIdx7][globalIdx7] -= TSecond[0][1];
1814 this->A_[globalIdx7][globalIdx8] -= TSecond[0][2];
1815 this->A_[globalIdx7][globalIdx3] -= TSecond[0][3];
1823 pcFlux[4][2] = Tu[0];
1824 pcPotential7 = Tu[0];
1827 pcFlux[6][1] = Tu[0];
1829 else if (caseL == 3)
1831 this->A_[globalIdx5][globalIdx5] += T[0][0];
1832 this->A_[globalIdx5][globalIdx7] += T[0][1];
1833 this->A_[globalIdx5][globalIdx8] += T[0][2];
1834 this->A_[globalIdx5][globalIdx1] += T[0][3];
1836 this->A_[globalIdx7][globalIdx5] -= TSecond[0][0];
1837 this->A_[globalIdx7][globalIdx7] -= TSecond[0][1];
1838 this->A_[globalIdx7][globalIdx8] -= TSecond[0][2];
1839 this->A_[globalIdx7][globalIdx1] -= TSecond[0][3];
1847 pcFlux[4][2] = Tu[0];
1848 pcPotential7 = Tu[0];
1851 pcFlux[6][1] = Tu[0];
1855 this->A_[globalIdx5][globalIdx5] += T[0][0];
1856 this->A_[globalIdx5][globalIdx7] += T[0][1];
1857 this->A_[globalIdx5][globalIdx6] += T[0][2];
1858 this->A_[globalIdx5][globalIdx3] += T[0][3];
1860 this->A_[globalIdx7][globalIdx5] -= TSecond[0][0];
1861 this->A_[globalIdx7][globalIdx7] -= TSecond[0][1];
1862 this->A_[globalIdx7][globalIdx6] -= TSecond[0][2];
1863 this->A_[globalIdx7][globalIdx3] -= TSecond[0][3];
1871 pcFlux[4][2] = Tu[0];
1872 pcPotential7 = Tu[0];
1875 pcFlux[6][1] = Tu[0];
1878 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx7)
1881 T *= this->interactionVolumes_.faceAreaFactor(interactionVolume, globalIdx5, 4, 2);
1882 TSecond *= this->interactionVolumes_.faceAreaFactor(interactionVolume, globalIdx7, 6, 1);
1886 this->A_[globalIdx5][globalIdx3] -= T[1][0];
1887 this->A_[globalIdx5][globalIdx7] -= T[1][1];
1888 this->A_[globalIdx5][globalIdx5] -= T[1][2];
1889 this->A_[globalIdx5][globalIdx4] -= T[1][3];
1891 this->A_[globalIdx7][globalIdx3] += TSecond[1][0];
1892 this->A_[globalIdx7][globalIdx7] += TSecond[1][1];
1893 this->A_[globalIdx7][globalIdx5] += TSecond[1][2];
1894 this->A_[globalIdx7][globalIdx4] += TSecond[1][3];
1902 pcFlux[4][2] = -Tu[1];
1903 pcPotential7 = -Tu[1];
1906 pcFlux[6][1] = -Tu[1];
1908 else if (caseL == 4)
1910 this->A_[globalIdx5][globalIdx5] += T[1][0];
1911 this->A_[globalIdx5][globalIdx1] += T[1][1];
1912 this->A_[globalIdx5][globalIdx7] += T[1][2];
1913 this->A_[globalIdx5][globalIdx2] += T[1][3];
1915 this->A_[globalIdx7][globalIdx5] -= TSecond[1][0];
1916 this->A_[globalIdx7][globalIdx1] -= TSecond[1][1];
1917 this->A_[globalIdx7][globalIdx7] -= TSecond[1][2];
1918 this->A_[globalIdx7][globalIdx2] -= TSecond[1][3];
1926 pcFlux[4][2] = Tu[1];
1927 pcPotential7 = Tu[1];
1930 pcFlux[6][1] = Tu[1];
1937 if (hangingNodeType == InteractionVolume::sixSmallCells)
1939 caseL = this->transmissibilityCalculator_.transmissibility(T, interactionVolume, lambda, 4, 0, 6,
1942 else if (hangingNodeType == InteractionVolume::twoSmallCells || hangingNodeType == InteractionVolume::fourSmallCellsFace)
1944 caseL = this->transmissibilityCalculator_.transmissibilityCaseTwo(T, interactionVolume, lambda, 4, 0, 2,
1947 else if (hangingNodeType == InteractionVolume::fourSmallCellsDiag
1948 || (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx7))
1950 useCases[0] =
false;
1952 useCases[2] =
false;
1955 caseL = this->transmissibilityCalculator_.transmissibility(T, interactionVolume, lambda, 4, 0, 6,
1958 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx6)
1960 useCases[0] =
false;
1963 useCases[3] =
false;
1965 caseL = this->transmissibilityCalculator_.transmissibility(T, interactionVolume, lambda, 4, 0, 6,
1970 T *= this->interactionVolumes_.faceAreaFactor(interactionVolume, globalIdx5, 4, 0);
1971 TSecond *= this->interactionVolumes_.faceAreaFactor(interactionVolume, globalIdx1, 0, 2);
1977 this->A_[globalIdx5][globalIdx5] += T[0][0];
1978 this->A_[globalIdx5][globalIdx1] += T[0][1];
1979 this->A_[globalIdx5][globalIdx7] += T[0][2];
1980 this->A_[globalIdx5][globalIdx6] += T[0][3];
1982 this->A_[globalIdx1][globalIdx5] -= TSecond[0][0];
1983 this->A_[globalIdx1][globalIdx1] -= TSecond[0][1];
1984 this->A_[globalIdx1][globalIdx7] -= TSecond[0][2];
1985 this->A_[globalIdx1][globalIdx6] -= TSecond[0][3];
1994 pcFlux[4][0] = Tu[0];
1995 pcPotential8 = Tu[0];
1998 pcFlux[0][2] = Tu[0];
2000 else if (caseL == 2)
2002 this->A_[globalIdx5][globalIdx5] += T[0][0];
2003 this->A_[globalIdx5][globalIdx1] += T[0][1];
2004 this->A_[globalIdx5][globalIdx3] += T[0][2];
2005 this->A_[globalIdx5][globalIdx2] += T[0][3];
2007 this->A_[globalIdx1][globalIdx5] -= TSecond[0][0];
2008 this->A_[globalIdx1][globalIdx1] -= TSecond[0][1];
2009 this->A_[globalIdx1][globalIdx3] -= TSecond[0][2];
2010 this->A_[globalIdx1][globalIdx2] -= TSecond[0][3];
2018 pcFlux[4][0] = Tu[0];
2019 pcPotential8 = Tu[0];
2022 pcFlux[0][2] = Tu[0];
2024 else if (caseL == 3)
2026 this->A_[globalIdx5][globalIdx5] += T[0][0];
2027 this->A_[globalIdx5][globalIdx1] += T[0][1];
2028 this->A_[globalIdx5][globalIdx3] += T[0][2];
2029 this->A_[globalIdx5][globalIdx6] += T[0][3];
2031 this->A_[globalIdx1][globalIdx5] -= TSecond[0][0];
2032 this->A_[globalIdx1][globalIdx1] -= TSecond[0][1];
2033 this->A_[globalIdx1][globalIdx3] -= TSecond[0][2];
2034 this->A_[globalIdx1][globalIdx6] -= TSecond[0][3];
2042 pcFlux[4][0] = Tu[0];
2043 pcPotential8 = Tu[0];
2046 pcFlux[0][2] = Tu[0];
2050 this->A_[globalIdx5][globalIdx5] += T[0][0];
2051 this->A_[globalIdx5][globalIdx1] += T[0][1];
2052 this->A_[globalIdx5][globalIdx7] += T[0][2];
2053 this->A_[globalIdx5][globalIdx2] += T[0][3];
2055 this->A_[globalIdx1][globalIdx5] -= TSecond[0][0];
2056 this->A_[globalIdx1][globalIdx1] -= TSecond[0][1];
2057 this->A_[globalIdx1][globalIdx7] -= TSecond[0][2];
2058 this->A_[globalIdx1][globalIdx2] -= TSecond[0][3];
2066 pcFlux[4][0] = Tu[0];
2067 pcPotential8 = Tu[0];
2070 pcFlux[0][2] = Tu[0];
2076 if (hangingNodeType == InteractionVolume::sixSmallCells)
2078 caseL = this->transmissibilityCalculator_.transmissibility(T, interactionVolume, lambda, 1, 5, 3,
2081 else if (hangingNodeType == InteractionVolume::twoSmallCells || hangingNodeType == InteractionVolume::fourSmallCellsFace)
2083 caseL = this->transmissibilityCalculator_.transmissibilityCaseOne(T, interactionVolume, lambda, 1, 5, 3,
2086 else if (hangingNodeType == InteractionVolume::fourSmallCellsDiag
2087 || (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx7))
2090 useCases[1] =
false;
2092 useCases[3] =
false;
2094 caseL = this->transmissibilityCalculator_.transmissibility(T, interactionVolume, lambda, 1, 5, 3,
2097 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx6)
2100 useCases[1] =
false;
2101 useCases[2] =
false;
2104 caseL = this->transmissibilityCalculator_.transmissibility(T, interactionVolume, lambda, 1, 5, 3,
2109 T *= this->interactionVolumes_.faceAreaFactor(interactionVolume, globalIdx2, 1, 2);
2110 TSecond *= this->interactionVolumes_.faceAreaFactor(interactionVolume, globalIdx6, 5, 0);
2114 this->A_[globalIdx2][globalIdx2] += T[0][0];
2115 this->A_[globalIdx2][globalIdx6] += T[0][1];
2116 this->A_[globalIdx2][globalIdx4] += T[0][2];
2117 this->A_[globalIdx2][globalIdx1] += T[0][3];
2119 this->A_[globalIdx6][globalIdx2] -= TSecond[0][0];
2120 this->A_[globalIdx6][globalIdx6] -= TSecond[0][1];
2121 this->A_[globalIdx6][globalIdx4] -= TSecond[0][2];
2122 this->A_[globalIdx6][globalIdx1] -= TSecond[0][3];
2131 pcFlux[1][2] = Tu[0];
2132 pcPotential9 = Tu[0];
2135 pcFlux[5][0] = Tu[0];
2137 else if (caseL == 2)
2139 this->A_[globalIdx2][globalIdx2] += T[0][0];
2140 this->A_[globalIdx2][globalIdx6] += T[0][1];
2141 this->A_[globalIdx2][globalIdx8] += T[0][2];
2142 this->A_[globalIdx2][globalIdx5] += T[0][3];
2144 this->A_[globalIdx6][globalIdx2] -= TSecond[0][0];
2145 this->A_[globalIdx6][globalIdx6] -= TSecond[0][1];
2146 this->A_[globalIdx6][globalIdx8] -= TSecond[0][2];
2147 this->A_[globalIdx6][globalIdx5] -= TSecond[0][3];
2156 pcFlux[1][2] = Tu[0];
2157 pcPotential9 = Tu[0];
2160 pcFlux[5][0] = Tu[0];
2162 else if (caseL == 3)
2164 this->A_[globalIdx2][globalIdx2] += T[0][0];
2165 this->A_[globalIdx2][globalIdx6] += T[0][1];
2166 this->A_[globalIdx2][globalIdx8] += T[0][2];
2167 this->A_[globalIdx2][globalIdx1] += T[0][3];
2169 this->A_[globalIdx6][globalIdx2] -= TSecond[0][0];
2170 this->A_[globalIdx6][globalIdx6] -= TSecond[0][1];
2171 this->A_[globalIdx6][globalIdx8] -= TSecond[0][2];
2172 this->A_[globalIdx6][globalIdx1] -= TSecond[0][3];
2180 pcFlux[1][2] = Tu[0];
2181 pcPotential9 = Tu[0];
2184 pcFlux[5][0] = Tu[0];
2188 this->A_[globalIdx2][globalIdx2] += T[0][0];
2189 this->A_[globalIdx2][globalIdx6] += T[0][1];
2190 this->A_[globalIdx2][globalIdx4] += T[0][2];
2191 this->A_[globalIdx2][globalIdx5] += T[0][3];
2193 this->A_[globalIdx6][globalIdx2] -= TSecond[0][0];
2194 this->A_[globalIdx6][globalIdx6] -= TSecond[0][1];
2195 this->A_[globalIdx6][globalIdx4] -= TSecond[0][2];
2196 this->A_[globalIdx6][globalIdx5] -= TSecond[0][3];
2204 pcFlux[1][2] = Tu[0];
2205 pcPotential9 = Tu[0];
2208 pcFlux[5][0] = Tu[0];
2214 if (hangingNodeType == InteractionVolume::fourSmallCellsFace)
2216 caseL = this->transmissibilityCalculator_.transmissibilityCaseTwo(T, interactionVolume, lambda, 7, 3, 1, 2);
2218 else if ((hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx7)
2219 || hangingNodeType == InteractionVolume::sixSmallCells)
2221 useCases[0] =
false;
2223 useCases[2] =
false;
2226 caseL = this->transmissibilityCalculator_.transmissibility(T, interactionVolume, lambda, 7, 3, 5, 1, 6,
2229 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx6)
2231 useCases[0] =
false;
2234 useCases[3] =
false;
2236 caseL = this->transmissibilityCalculator_.transmissibility(T, interactionVolume, lambda, 7, 3, 5, 1, 6,
2239 else if (hangingNodeType == InteractionVolume::fourSmallCellsDiag)
2242 useCases[1] =
false;
2244 useCases[3] =
false;
2246 caseL = this->transmissibilityCalculator_.transmissibility(T, interactionVolume, lambda, 7, 3, 5, 1, 6,
2251 T *= this->interactionVolumes_.faceAreaFactor(interactionVolume, globalIdx8, 7, 0);
2252 TSecond *= this->interactionVolumes_.faceAreaFactor(interactionVolume, globalIdx4, 3, 2);
2254 if (hangingNodeType != InteractionVolume::twoSmallCells)
2259 this->A_[globalIdx8][globalIdx8] += T[0][0];
2260 this->A_[globalIdx8][globalIdx4] += T[0][1];
2261 this->A_[globalIdx8][globalIdx6] += T[0][2];
2262 this->A_[globalIdx8][globalIdx7] += T[0][3];
2264 this->A_[globalIdx4][globalIdx8] -= TSecond[0][0];
2265 this->A_[globalIdx4][globalIdx4] -= TSecond[0][1];
2266 this->A_[globalIdx4][globalIdx6] -= TSecond[0][2];
2267 this->A_[globalIdx4][globalIdx7] -= TSecond[0][3];
2276 pcFlux[7][0] = Tu[0];
2277 pcPotential10 = Tu[0];
2280 pcFlux[3][2] = Tu[0];
2282 else if (caseL == 2)
2284 this->A_[globalIdx8][globalIdx8] += T[0][0];
2285 this->A_[globalIdx8][globalIdx4] += T[0][1];
2286 this->A_[globalIdx8][globalIdx2] += T[0][2];
2287 this->A_[globalIdx8][globalIdx3] += T[0][3];
2289 this->A_[globalIdx4][globalIdx8] -= TSecond[0][0];
2290 this->A_[globalIdx4][globalIdx4] -= TSecond[0][1];
2291 this->A_[globalIdx4][globalIdx2] -= TSecond[0][2];
2292 this->A_[globalIdx4][globalIdx3] -= TSecond[0][3];
2300 pcFlux[7][0] = Tu[0];
2301 pcPotential10 = Tu[0];
2304 pcFlux[3][2] = Tu[0];
2306 else if (caseL == 3)
2308 this->A_[globalIdx8][globalIdx8] += T[0][0];
2309 this->A_[globalIdx8][globalIdx4] += T[0][1];
2310 this->A_[globalIdx8][globalIdx2] += T[0][2];
2311 this->A_[globalIdx8][globalIdx7] += T[0][3];
2313 this->A_[globalIdx4][globalIdx8] -= TSecond[0][0];
2314 this->A_[globalIdx4][globalIdx4] -= TSecond[0][1];
2315 this->A_[globalIdx4][globalIdx2] -= TSecond[0][2];
2316 this->A_[globalIdx4][globalIdx7] -= TSecond[0][3];
2324 pcFlux[7][0] = Tu[0];
2325 pcPotential10 = Tu[0];
2328 pcFlux[3][2] = Tu[0];
2332 this->A_[globalIdx8][globalIdx8] += T[0][0];
2333 this->A_[globalIdx8][globalIdx4] += T[0][1];
2334 this->A_[globalIdx8][globalIdx6] += T[0][2];
2335 this->A_[globalIdx8][globalIdx3] += T[0][3];
2337 this->A_[globalIdx4][globalIdx8] -= TSecond[0][0];
2338 this->A_[globalIdx4][globalIdx4] -= TSecond[0][1];
2339 this->A_[globalIdx4][globalIdx6] -= TSecond[0][2];
2340 this->A_[globalIdx4][globalIdx3] -= TSecond[0][3];
2348 pcFlux[7][0] = Tu[0];
2349 pcPotential10 = Tu[0];
2352 pcFlux[3][2] = Tu[0];
2359 if (hangingNodeType == InteractionVolume::fourSmallCellsFace)
2361 caseL = this->transmissibilityCalculator_.transmissibilityCaseOne(T, interactionVolume, lambda, 2, 6, 0, 3);
2363 else if ((hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx7)
2364 || hangingNodeType == InteractionVolume::sixSmallCells)
2367 useCases[1] =
false;
2369 useCases[3] =
false;
2371 caseL = this->transmissibilityCalculator_.transmissibility(T, interactionVolume, lambda, 2, 6, 0, 4, 3,
2374 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge && globalIdx5 != globalIdx6)
2377 useCases[1] =
false;
2378 useCases[2] =
false;
2381 caseL = this->transmissibilityCalculator_.transmissibility(T, interactionVolume, lambda, 2, 6, 0, 4, 3,
2384 else if (hangingNodeType == InteractionVolume::fourSmallCellsDiag)
2386 useCases[0] =
false;
2388 useCases[2] =
false;
2391 caseL = this->transmissibilityCalculator_.transmissibility(T, interactionVolume, lambda, 2, 6, 0, 4, 3,
2396 T *= this->interactionVolumes_.faceAreaFactor(interactionVolume, globalIdx3, 2, 2);
2397 TSecond *= this->interactionVolumes_.faceAreaFactor(interactionVolume, globalIdx7, 6, 0);
2399 if (hangingNodeType != InteractionVolume::twoSmallCells)
2403 this->A_[globalIdx3][globalIdx3] += T[0][0];
2404 this->A_[globalIdx3][globalIdx7] += T[0][1];
2405 this->A_[globalIdx3][globalIdx1] += T[0][2];
2406 this->A_[globalIdx3][globalIdx4] += T[0][3];
2408 this->A_[globalIdx7][globalIdx3] -= TSecond[0][0];
2409 this->A_[globalIdx7][globalIdx7] -= TSecond[0][1];
2410 this->A_[globalIdx7][globalIdx1] -= TSecond[0][2];
2411 this->A_[globalIdx7][globalIdx4] -= TSecond[0][3];
2420 pcFlux[2][2] = Tu[0];
2421 pcPotential11 = Tu[0];
2424 pcFlux[6][0] = Tu[0];
2426 else if (caseL == 2)
2428 this->A_[globalIdx3][globalIdx3] += T[0][0];
2429 this->A_[globalIdx3][globalIdx7] += T[0][1];
2430 this->A_[globalIdx3][globalIdx5] += T[0][2];
2431 this->A_[globalIdx3][globalIdx8] += T[0][3];
2433 this->A_[globalIdx7][globalIdx3] -= TSecond[0][0];
2434 this->A_[globalIdx7][globalIdx7] -= TSecond[0][1];
2435 this->A_[globalIdx7][globalIdx5] -= TSecond[0][2];
2436 this->A_[globalIdx7][globalIdx8] -= TSecond[0][3];
2445 pcFlux[2][2] = Tu[0];
2446 pcPotential11 = Tu[0];
2449 pcFlux[6][0] = Tu[0];
2451 else if (caseL == 3)
2453 this->A_[globalIdx3][globalIdx3] += T[0][0];
2454 this->A_[globalIdx3][globalIdx7] += T[0][1];
2455 this->A_[globalIdx3][globalIdx5] += T[0][2];
2456 this->A_[globalIdx3][globalIdx4] += T[0][3];
2458 this->A_[globalIdx7][globalIdx3] -= TSecond[0][0];
2459 this->A_[globalIdx7][globalIdx7] -= TSecond[0][1];
2460 this->A_[globalIdx7][globalIdx5] -= TSecond[0][2];
2461 this->A_[globalIdx7][globalIdx4] -= TSecond[0][3];
2469 pcFlux[2][2] = Tu[0];
2470 pcPotential11 = Tu[0];
2473 pcFlux[6][0] = Tu[0];
2477 this->A_[globalIdx3][globalIdx3] += T[0][0];
2478 this->A_[globalIdx3][globalIdx7] += T[0][1];
2479 this->A_[globalIdx3][globalIdx1] += T[0][2];
2480 this->A_[globalIdx3][globalIdx8] += T[0][3];
2482 this->A_[globalIdx7][globalIdx3] -= TSecond[0][0];
2483 this->A_[globalIdx7][globalIdx7] -= TSecond[0][1];
2484 this->A_[globalIdx7][globalIdx1] -= TSecond[0][2];
2485 this->A_[globalIdx7][globalIdx8] -= TSecond[0][3];
2493 pcFlux[2][2] = Tu[0];
2494 pcPotential11 = Tu[0];
2497 pcFlux[6][0] = Tu[0];
2501 if (pc[0] == 0 && pc[1] == 0 && pc[2] == 0 && pc[3] == 0 && pc[4] == 0 && pc[5] == 0 && pc[6] == 0
2508 Dune::FieldVector<Scalar, numPhases> lambda0Upw(0.0);
2509 lambda0Upw[wPhaseIdx] = (pcPotential0 >= 0) ? lambda1[wPhaseIdx] : lambda2[wPhaseIdx];
2510 lambda0Upw[nPhaseIdx] = (pcPotential0 >= 0) ? lambda1[nPhaseIdx] : lambda2[nPhaseIdx];
2513 Dune::FieldVector<Scalar, numPhases> lambda1Upw(0.0);
2514 lambda1Upw[wPhaseIdx] = (pcPotential1 >= 0) ? lambda2[wPhaseIdx] : lambda4[wPhaseIdx];
2515 lambda1Upw[nPhaseIdx] = (pcPotential1 >= 0) ? lambda2[nPhaseIdx] : lambda4[nPhaseIdx];
2518 Dune::FieldVector<Scalar, numPhases> lambda2Upw(0.0);
2519 lambda2Upw[wPhaseIdx] = (pcPotential2 >= 0) ? lambda4[wPhaseIdx] : lambda3[wPhaseIdx];
2520 lambda2Upw[nPhaseIdx] = (pcPotential2 >= 0) ? lambda4[nPhaseIdx] : lambda3[nPhaseIdx];
2523 Dune::FieldVector<Scalar, numPhases> lambda3Upw(0.0);
2524 lambda3Upw[wPhaseIdx] = (pcPotential3 >= 0) ? lambda3[wPhaseIdx] : lambda1[wPhaseIdx];
2525 lambda3Upw[nPhaseIdx] = (pcPotential3 >= 0) ? lambda3[nPhaseIdx] : lambda1[nPhaseIdx];
2528 Dune::FieldVector<Scalar, numPhases> lambda4Upw(0.0);
2529 lambda4Upw[wPhaseIdx] = (pcPotential4 >= 0) ? lambda6[wPhaseIdx] : lambda5[wPhaseIdx];
2530 lambda4Upw[nPhaseIdx] = (pcPotential4 >= 0) ? lambda6[nPhaseIdx] : lambda5[nPhaseIdx];
2533 Dune::FieldVector<Scalar, numPhases> lambda5Upw(0.0);
2534 lambda5Upw[wPhaseIdx] = (pcPotential5 >= 0) ? lambda8[wPhaseIdx] : lambda6[wPhaseIdx];
2535 lambda5Upw[nPhaseIdx] = (pcPotential5 >= 0) ? lambda8[nPhaseIdx] : lambda6[nPhaseIdx];
2538 Dune::FieldVector<Scalar, numPhases> lambda6Upw(0.0);
2539 lambda6Upw[wPhaseIdx] = (pcPotential6 >= 0) ? lambda7[wPhaseIdx] : lambda8[wPhaseIdx];
2540 lambda6Upw[nPhaseIdx] = (pcPotential6 >= 0) ? lambda7[nPhaseIdx] : lambda8[nPhaseIdx];
2543 Dune::FieldVector<Scalar, numPhases> lambda7Upw(0.0);
2544 lambda7Upw[wPhaseIdx] = (pcPotential7 >= 0) ? lambda5[wPhaseIdx] : lambda7[wPhaseIdx];
2545 lambda7Upw[nPhaseIdx] = (pcPotential7 >= 0) ? lambda5[nPhaseIdx] : lambda7[nPhaseIdx];
2548 Dune::FieldVector<Scalar, numPhases> lambda8Upw(0.0);
2549 lambda8Upw[wPhaseIdx] = (pcPotential8 >= 0) ? lambda5[wPhaseIdx] : lambda1[wPhaseIdx];
2550 lambda8Upw[nPhaseIdx] = (pcPotential8 >= 0) ? lambda5[nPhaseIdx] : lambda1[nPhaseIdx];
2553 Dune::FieldVector<Scalar, numPhases> lambda9Upw(0.0);
2554 lambda9Upw[wPhaseIdx] = (pcPotential9 >= 0) ? lambda2[wPhaseIdx] : lambda6[wPhaseIdx];
2555 lambda9Upw[nPhaseIdx] = (pcPotential9 >= 0) ? lambda2[nPhaseIdx] : lambda6[nPhaseIdx];
2558 Dune::FieldVector<Scalar, numPhases> lambda10Upw(0.0);
2559 lambda10Upw[wPhaseIdx] = (pcPotential10 >= 0) ? lambda8[wPhaseIdx] : lambda4[wPhaseIdx];
2560 lambda10Upw[nPhaseIdx] = (pcPotential10 >= 0) ? lambda8[nPhaseIdx] : lambda4[nPhaseIdx];
2563 Dune::FieldVector<Scalar, numPhases> lambda11Upw(0.0);
2564 lambda11Upw[wPhaseIdx] = (pcPotential11 >= 0) ? lambda3[wPhaseIdx] : lambda7[wPhaseIdx];
2565 lambda11Upw[nPhaseIdx] = (pcPotential11 >= 0) ? lambda3[nPhaseIdx] : lambda7[nPhaseIdx];
2567 for (
int i = 0; i < numPhases; i++)
2569 Scalar lambdaT0 = lambda0Upw[wPhaseIdx] + lambda0Upw[nPhaseIdx];
2570 Scalar lambdaT1 = lambda1Upw[wPhaseIdx] + lambda1Upw[nPhaseIdx];
2571 Scalar lambdaT2 = lambda2Upw[wPhaseIdx] + lambda2Upw[nPhaseIdx];
2572 Scalar lambdaT3 = lambda3Upw[wPhaseIdx] + lambda3Upw[nPhaseIdx];
2573 Scalar lambdaT4 = lambda4Upw[wPhaseIdx] + lambda4Upw[nPhaseIdx];
2574 Scalar lambdaT5 = lambda5Upw[wPhaseIdx] + lambda5Upw[nPhaseIdx];
2575 Scalar lambdaT6 = lambda6Upw[wPhaseIdx] + lambda6Upw[nPhaseIdx];
2576 Scalar lambdaT7 = lambda7Upw[wPhaseIdx] + lambda7Upw[nPhaseIdx];
2577 Scalar lambdaT8 = lambda8Upw[wPhaseIdx] + lambda8Upw[nPhaseIdx];
2578 Scalar lambdaT9 = lambda9Upw[wPhaseIdx] + lambda9Upw[nPhaseIdx];
2579 Scalar lambdaT10 = lambda10Upw[wPhaseIdx] + lambda10Upw[nPhaseIdx];
2580 Scalar lambdaT11 = lambda11Upw[wPhaseIdx] + lambda11Upw[nPhaseIdx];
2581 Scalar fracFlow0 = (lambdaT0 > threshold_) ? lambda0Upw[i] / (lambdaT0) : 0.0;
2582 Scalar fracFlow1 = (lambdaT1 > threshold_) ? lambda1Upw[i] / (lambdaT1) : 0.0;
2583 Scalar fracFlow2 = (lambdaT2 > threshold_) ? lambda2Upw[i] / (lambdaT2) : 0.0;
2584 Scalar fracFlow3 = (lambdaT3 > threshold_) ? lambda3Upw[i] / (lambdaT3) : 0.0;
2585 Scalar fracFlow4 = (lambdaT4 > threshold_) ? lambda4Upw[i] / (lambdaT4) : 0.0;
2586 Scalar fracFlow5 = (lambdaT5 > threshold_) ? lambda5Upw[i] / (lambdaT5) : 0.0;
2587 Scalar fracFlow6 = (lambdaT6 > threshold_) ? lambda6Upw[i] / (lambdaT6) : 0.0;
2588 Scalar fracFlow7 = (lambdaT7 > threshold_) ? lambda7Upw[i] / (lambdaT7) : 0.0;
2589 Scalar fracFlow8 = (lambdaT8 > threshold_) ? lambda8Upw[i] / (lambdaT8) : 0.0;
2590 Scalar fracFlow9 = (lambdaT9 > threshold_) ? lambda9Upw[i] / (lambdaT9) : 0.0;
2591 Scalar fracFlow10 = (lambdaT10 > threshold_) ? lambda10Upw[i] / (lambdaT10) : 0.0;
2592 Scalar fracFlow11 = (lambdaT11 > threshold_) ? lambda11Upw[i] / (lambdaT11) : 0.0;
2594 switch (pressureType_)
2601 this->f_[globalIdx1] -= (fracFlow0 * pcFlux[0][0] - fracFlow3 * pcFlux[0][1] - fracFlow8 * pcFlux[0][2]);
2602 this->f_[globalIdx2] -= (fracFlow1 * pcFlux[1][0] - fracFlow0 * pcFlux[1][1] + fracFlow9 * pcFlux[1][2]);
2603 this->f_[globalIdx3] -= (fracFlow3 * pcFlux[2][0] - fracFlow2 * pcFlux[2][1] + fracFlow11 * pcFlux[2][2]);
2604 this->f_[globalIdx4] -= (fracFlow2 * pcFlux[3][0] - fracFlow1 * pcFlux[3][1] - fracFlow10 * pcFlux[3][2]);
2605 this->f_[globalIdx5] -= (fracFlow8 * pcFlux[4][0] - fracFlow4 * pcFlux[4][1] + fracFlow7 * pcFlux[4][2]);
2606 this->f_[globalIdx6] -= (-fracFlow9 * pcFlux[5][0] - fracFlow5 * pcFlux[5][1] + fracFlow4 * pcFlux[5][2]);
2607 this->f_[globalIdx7] -= (-fracFlow11 * pcFlux[6][0] - fracFlow7 * pcFlux[6][1] + fracFlow6 * pcFlux[6][2]);
2608 this->f_[globalIdx8] -= (fracFlow10 * pcFlux[7][0] - fracFlow6 * pcFlux[7][1] + fracFlow5 * pcFlux[7][2]);
2617 this->f_[globalIdx1] += (fracFlow0 * pcFlux[0][0] - fracFlow3 * pcFlux[0][1] - fracFlow8 * pcFlux[0][2]);
2618 this->f_[globalIdx2] += (fracFlow1 * pcFlux[1][0] - fracFlow0 * pcFlux[1][1] + fracFlow9 * pcFlux[1][2]);
2619 this->f_[globalIdx3] += (fracFlow3 * pcFlux[2][0] - fracFlow2 * pcFlux[2][1] + fracFlow11 * pcFlux[2][2]);
2620 this->f_[globalIdx4] += (fracFlow2 * pcFlux[3][0] - fracFlow1 * pcFlux[3][1] - fracFlow10 * pcFlux[3][2]);
2621 this->f_[globalIdx5] += (fracFlow8 * pcFlux[4][0] - fracFlow4 * pcFlux[4][1] + fracFlow7 * pcFlux[4][2]);
2622 this->f_[globalIdx6] += (-fracFlow9 * pcFlux[5][0] - fracFlow5 * pcFlux[5][1] + fracFlow4 * pcFlux[5][2]);
2623 this->f_[globalIdx7] += (-fracFlow11 * pcFlux[6][0] - fracFlow7 * pcFlux[6][1] + fracFlow6 * pcFlux[6][2]);
2624 this->f_[globalIdx8] += (fracFlow10 * pcFlux[7][0] - fracFlow6 * pcFlux[7][1] + fracFlow5 * pcFlux[7][2]);
Provides methods for transmissibility calculation in 3-d.
Interactionvolume container for 3-d MPFA L-method on an h-adaptive grid.
3-d finite Volume-MPFAL implementation of a two-phase pressure equation.
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
3-d finite volume MPFA L-method discretization of a two-phase flow pressure equation of the sequentia...
Definition: 3dpressure.hh:76
typename TransmissibilityCalculator::TransmissibilityType TransmissibilityType
Type including methods for calculation of MPFA transmissibilities.
Definition: 3dpressure.hh:154
void update()
Pressure update.
Definition: 3dpressure.hh:273
InteractionVolumeContainer interactionVolumes_
Definition: 3dpressure.hh:509
GetPropType< TypeTag, Properties::MPFAInteractionVolume > InteractionVolume
Type for storing interaction volume information.
Definition: 3dpressure.hh:156
3-d finite volume MPFA L-method discretization of a two-phase flow pressure equation of the sequentia...
Definition: 3dpressureadaptive.hh:78
void update()
updates the pressure field (analog to update function in IMPET)
Definition: 3dpressureadaptive.hh:194
void assembleHangingNodeInteractionVolume(InteractionVolume &interactionVolume)
assembles the matrix entries of one hanging node interaction volume into the global matrix
Definition: 3dpressureadaptive.hh:456
void assemble()
assembles the global matrix and rhs vector for the pressure solution
Definition: 3dpressureadaptive.hh:412
void initializeMatrix()
Initializes the sparse matrix for the pressure solution.
Definition: 3dpressureadaptive.hh:271
void initializeMatrixRowSize()
Initializes the row size of the sparse matrix for the pressure solution.
Definition: 3dpressureadaptive.hh:281
void initializeMatrixIndices()
Initializes the indices of the sparse matrix for the pressure solution.
Definition: 3dpressureadaptive.hh:365
FvMpfaL3dPressure2pAdaptive(Problem &problem)
Constructs a FvMpfaL3dPressure2pAdaptive object.
Definition: 3dpressureadaptive.hh:217
friend ParentType
Definition: 3dpressureadaptive.hh:167
void initialize(bool solveTwice=true)
Initializes the pressure model.
Definition: 3dpressureadaptive.hh:176
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
The finite volume base class for the solution of a pressure equation.
Definition: sequential/cellcentered/pressure.hh:49
void initialize()
Initialize pressure model.
Definition: sequential/cellcentered/pressure.hh:213
PressureSolution & pressure()
Returns the vector containing the pressure solution.
Definition: sequential/cellcentered/pressure.hh:120
Matrix A_
Global stiffnes matrix (sparse matrix which is build by the initializeMatrix() function)
Definition: sequential/cellcentered/pressure.hh:324
RHSVector f_
Right hand side vector.
Definition: sequential/cellcentered/pressure.hh:325