3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
omethod/2dpressure.hh
Go to the documentation of this file.
1/*****************************************************************************
2 * See the file COPYING for full copying permissions. *
3 * *
4 * This program is free software: you can redistribute it and/or modify *
5 * it under the terms of the GNU General Public License as published by *
6 * the Free Software Foundation, either version 3 of the License, or *
7 * (at your option) any later version. *
8 * *
9 * This program is distributed in the hope that it will be useful, *
10 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
11 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
12 * GNU General Public License for more details. *
13 * *
14 * You should have received a copy of the GNU General Public License *
15 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
16 *****************************************************************************/
22#ifndef DUMUX_FVMPFAO2DPRESSURE2P_HH
23#define DUMUX_FVMPFAO2DPRESSURE2P_HH
24
25// dumux environment
30
32
33namespace Dumux {
67template<class TypeTag>
68class FvMpfaO2dPressure2p: public FVPressure<TypeTag>
69{
72
73 enum
74 {
75 dim = GridView::dimension, dimWorld = GridView::dimensionworld
76 };
77
80
82
84
87
91 using PrimaryVariables = typename SolutionTypes::PrimaryVariables;
92 using ScalarSolutionType = typename SolutionTypes::ScalarSolution;
93
95
96 enum
97 {
98 pw = Indices::pressureW,
99 pn = Indices::pressureNw,
100 sw = Indices::saturationW,
101 sn = Indices::saturationNw
102 };
103 enum
104 {
105 wPhaseIdx = Indices::wPhaseIdx,
106 nPhaseIdx = Indices::nPhaseIdx,
107 pressureIdx = Indices::pressureIdx,
108 saturationIdx = Indices::saturationIdx,
109 pressEqIdx = Indices::pressureEqIdx,
110 satEqIdx = Indices::satEqIdx,
111 numPhases = getPropValue<TypeTag, Properties::NumPhases>()
112 };
113 enum
114 {
115 globalCorner = 2,
116 globalEdge = 3,
117 neumannNeumann = 0,
118 dirichletDirichlet = 1,
119 dirichletNeumann = 2,
120 neumannDirichlet = 3
121 };
122
123 using Element = typename GridView::Traits::template Codim<0>::Entity;
124 using IntersectionIterator = typename GridView::IntersectionIterator;
125 using Intersection = typename GridView::Intersection;
126
127 using LocalPosition = Dune::FieldVector<Scalar, dim>;
128 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
129 using GravityVector = Dune::FieldVector<Scalar, dimWorld>;
130 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
131
132 using DimVector = Dune::FieldVector<Scalar, dim>;
133
135
136 using GlobalInteractionVolumeVector = std::vector<InteractionVolume>;
137 using InnerBoundaryVolumeFaces = std::vector<Dune::FieldVector<bool, 2*dim> >;
138
140 Intersection getNextIntersection_(const Element&, const IntersectionIterator&);
141
143 friend class FVPressure<TypeTag>;
144 void initializeMatrix();
145
146 void storeInteractionVolumeInfo();
147
149 void assemble();
150public:
151
153 void updateMaterialLaws();
154
161 {
162 interactionVolumes_.clear();
164
165 interactionVolumes_.resize(problem_.gridView().size(dim), InteractionVolume(problem_.gridView().grid()));
166 innerBoundaryVolumeFaces_.resize(problem_.gridView().size(0), Dune::FieldVector<bool, 2 * dim>(false));
167
168 storeInteractionVolumeInfo();
169 }
170
177 {
179
180 const auto element = *problem_.gridView().template begin<0>();
181 FluidState fluidState;
182 fluidState.setPressure(wPhaseIdx, problem_.referencePressure(element));
183 fluidState.setPressure(nPhaseIdx, problem_.referencePressure(element));
184 fluidState.setTemperature(problem_.temperature(element));
185 fluidState.setSaturation(wPhaseIdx, 1.);
186 fluidState.setSaturation(nPhaseIdx, 0.);
187 density_[wPhaseIdx] = FluidSystem::density(fluidState, wPhaseIdx);
188 density_[nPhaseIdx] = FluidSystem::density(fluidState, nPhaseIdx);
189 viscosity_[wPhaseIdx] = FluidSystem::viscosity(fluidState, wPhaseIdx);
190 viscosity_[nPhaseIdx] = FluidSystem::viscosity(fluidState, nPhaseIdx);
191
193
194 interactionVolumes_.resize(problem_.gridView().size(dim), InteractionVolume(problem_.gridView().grid()));
195 innerBoundaryVolumeFaces_.resize(problem_.gridView().size(0), Dune::FieldVector<bool, 2 * dim>(false));
196
197 storeInteractionVolumeInfo();
198
199 assemble();
200 this->solve();
201
203
204 return;
205 }
206
212 void update()
213 {
214 //error bounds for error term for incompressible models to correct unphysical saturation
215 //over/undershoots due to saturation transport
216 timeStep_ = problem_.timeManager().timeStepSize();
217 maxError_ = 0.0;
218 int size = problem_.gridView().size(0);
219 for (int i = 0; i < size; i++)
220 {
221 Scalar sat = 0;
222 using std::max;
223 switch (saturationType_)
224 {
225 case sw:
226 sat = problem_.variables().cellData(i).saturation(wPhaseIdx);
227 break;
228 case sn:
229 sat = problem_.variables().cellData(i).saturation(nPhaseIdx);
230 break;
231 }
232 if (sat > 1.0)
233 {
234 maxError_ = max(maxError_, (sat - 1.0) / timeStep_);
235 }
236 if (sat < 0.0)
237 {
238 maxError_ = max(maxError_, (-sat) / timeStep_);
239 }
240 }
241
242 assemble();
243 this->solve();
244
246
247 return;
248 }
249
254 {
255 // iterate through leaf grid an evaluate c0 at cell center
256 for (const auto& element : elements(problem_.gridView()))
257 {
258 storePressureSolution(element);
259 }
260 }
261
267 void storePressureSolution(const Element& element)
268 {
269 int eIdxGlobal = problem_.variables().index(element);
270 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
271
272 switch (pressureType_)
273 {
274 case pw:
275 {
276 Scalar potW = this->pressure()[eIdxGlobal];
277
278 Scalar gravityDiff = (problem_.bBoxMax() - element.geometry().center()) * gravity_;
279 Scalar potPc = cellData.capillaryPressure()
280 + gravityDiff * (density_[nPhaseIdx] - density_[wPhaseIdx]);
281
282 cellData.setPotential(wPhaseIdx, potW);
283 cellData.setPotential(nPhaseIdx, potW + potPc);
284
285 Scalar pressW = potW - gravityDiff * density_[wPhaseIdx];
286
287 cellData.setPressure(wPhaseIdx, pressW);
288 cellData.setPressure(nPhaseIdx, pressW + cellData.capillaryPressure());
289
290 break;
291 }
292 case pn:
293 {
294 Scalar potNw = this->pressure()[eIdxGlobal];
295
296 Scalar gravityDiff = (problem_.bBoxMax() - element.geometry().center()) * gravity_;
297 Scalar potPc = cellData.capillaryPressure()
298 + gravityDiff * (density_[nPhaseIdx] - density_[wPhaseIdx]);
299
300 cellData.setPotential(nPhaseIdx, potNw);
301 cellData.setPotential(wPhaseIdx, potNw - potPc);
302
303 Scalar pressNw = potNw - gravityDiff * density_[nPhaseIdx];
304
305 cellData.setPressure(wPhaseIdx, pressNw - cellData.capillaryPressure());
306 cellData.setPressure(nPhaseIdx, pressNw);
307
308 break;
309 }
310 }
311
312 cellData.fluxData().resetVelocity();
313 }
314
325 template<class MultiWriter>
326 void addOutputVtkFields(MultiWriter &writer)
327 {
328 int size = problem_.gridView().size(0);
329 ScalarSolutionType *potential = writer.allocateManagedBuffer(size);
330
331 (*potential) = this->pressure();
332
333 if (pressureType_ == pw)
334 {
335 writer.attachCellData(*potential, "wetting potential");
336 }
337
338 if (pressureType_ == pn)
339 {
340 writer.attachCellData(*potential, "nonwetting potential");
341 }
342
343 if (vtkOutputLevel_ > 0)
344 {
345 ScalarSolutionType *pressure = writer.allocateManagedBuffer(size);
346 ScalarSolutionType *pressureSecond = writer.allocateManagedBuffer(size);
347 ScalarSolutionType *potentialSecond = writer.allocateManagedBuffer(size);
348 ScalarSolutionType *pc = writer.allocateManagedBuffer(size);
349
350 for (const auto& element : elements(problem_.gridView()))
351 {
352 int idx = problem_.variables().index(element);
353 CellData& cellData = problem_.variables().cellData(idx);
354
355 (*pc)[idx] = cellData.capillaryPressure();
356
357 if (pressureType_ == pw)
358 {
359 (*pressure)[idx] = cellData.pressure(wPhaseIdx);
360 (*potentialSecond)[idx] = cellData.potential(nPhaseIdx);
361 (*pressureSecond)[idx] = cellData.pressure(nPhaseIdx);
362 }
363
364 if (pressureType_ == pn)
365 {
366 (*pressure)[idx] = cellData.pressure(nPhaseIdx);
367 (*potentialSecond)[idx] = cellData.potential(wPhaseIdx);
368 (*pressureSecond)[idx] = cellData.pressure(wPhaseIdx);
369 }
370 }
371
372 if (pressureType_ == pw)
373 {
374 writer.attachCellData(*pressure, "wetting pressure");
375 writer.attachCellData(*pressureSecond, "nonwetting pressure");
376 writer.attachCellData(*potentialSecond, "nonwetting potential");
377 }
378
379 if (pressureType_ == pn)
380 {
381 writer.attachCellData(*pressure, "nonwetting pressure");
382 writer.attachCellData(*pressureSecond, "wetting pressure");
383 writer.attachCellData(*potentialSecond, "wetting potential");
384 }
385
386 writer.attachCellData(*pc, "capillary pressure");
387 }
388
389 return;
390 }
391
396 FvMpfaO2dPressure2p(Problem& problem) :
397 ParentType(problem), problem_(problem), gravity_(problem.gravity()), maxError_(
398 0.), timeStep_(1.)
399 {
400 if (pressureType_ != pw && pressureType_ != pn)
401 {
402 DUNE_THROW(Dune::NotImplemented, "Pressure type not supported!");
403 }
404 if (saturationType_ != sw && saturationType_ != sn)
405 {
406 DUNE_THROW(Dune::NotImplemented, "Saturation type not supported!");
407 }
408 if (getPropValue<TypeTag, Properties::EnableCompressibility>())
409 {
410 DUNE_THROW(Dune::NotImplemented, "Compressibility not supported!");
411 }
412 if (dim != 2)
413 {
414 DUNE_THROW(Dune::NotImplemented, "Dimension not supported!");
415 }
416
417 ErrorTermFactor_ = getParam<Scalar>("Impet.ErrorTermFactor");
418 ErrorTermLowerBound_ = getParam<Scalar>("Impet.ErrorTermLowerBound");
419 ErrorTermUpperBound_ = getParam<Scalar>("Impet.ErrorTermUpperBound");
420
421 density_[wPhaseIdx] = 0.;
422 density_[nPhaseIdx] = 0.;
423 viscosity_[wPhaseIdx] = 0.;
424 viscosity_[nPhaseIdx] = 0.;
425
426 vtkOutputLevel_ = getParam<int>("Vtk.OutputLevel");
427 }
428
429private:
430 Problem& problem_;
431
448 Scalar evaluateErrorTerm_(CellData& cellData)
449 {
450 //error term for incompressible models to correct unphysical saturation over/undershoots due to saturation transport
451 // error reduction routine: volumetric error is damped and inserted to right hand side
452 Scalar sat = 0;
453 switch (saturationType_)
454 {
455 case sw:
456 sat = cellData.saturation(wPhaseIdx);
457 break;
458 case sn:
459 sat = cellData.saturation(nPhaseIdx);
460 break;
461 }
462
463 Scalar error = (sat > 1.0) ? sat - 1.0 : 0.0;
464 if (sat < 0.0)
465 {
466 error = sat;
467 }
468 error /= timeStep_;
469
470 using std::abs;
471 Scalar errorAbs = abs(error);
472
473 if ((errorAbs * timeStep_ > 1e-6) && (errorAbs > ErrorTermLowerBound_ * maxError_)
474 && (!problem_.timeManager().willBeFinished()))
475 {
476 return ErrorTermFactor_ * error;
477 }
478 return 0.0;
479 }
480
481protected:
482 GlobalInteractionVolumeVector interactionVolumes_;
483 InnerBoundaryVolumeFaces innerBoundaryVolumeFaces_;
484
485private:
486 const GravityVector& gravity_;
487
488 Scalar maxError_;
489 Scalar timeStep_;
490 Scalar ErrorTermFactor_;
491 Scalar ErrorTermLowerBound_;
492 Scalar ErrorTermUpperBound_;
493
494 Scalar density_[numPhases];
495 Scalar viscosity_[numPhases];
496
497 int vtkOutputLevel_;
498
499 static constexpr Scalar threshold_ = 1e-15;
501 static const int pressureType_ = getPropValue<TypeTag, Properties::PressureFormulation>();
503 static const int saturationType_ = getPropValue<TypeTag, Properties::SaturationFormulation>();
505 static const int velocityType_ = getPropValue<TypeTag, Properties::VelocityFormulation>();
506};
507
508// TODO doc me!
509template<class TypeTag>
510typename FvMpfaO2dPressure2p<TypeTag>::Intersection
512 const IntersectionIterator& isIt)
513{
514 auto isItBegin = problem_.gridView().ibegin(element);
515 const auto isEndIt = problem_.gridView().iend(element);
516
517 auto tempIsIt = isIt;
518 auto nextIsIt = ++tempIsIt;
519
520 // get 'nextIsIt'
521 switch (getPropValue<TypeTag, Properties::GridImplementation>())
522 {
523 // for YaspGrid
524 case GridTypeIndices::yaspGrid:
525 {
526 if (nextIsIt == isEndIt)
527 {
528 nextIsIt = isItBegin;
529 }
530 else
531 {
532 nextIsIt = ++tempIsIt;
533
534 if (nextIsIt == isEndIt)
535 {
536 auto tempIsItBegin = isItBegin;
537 nextIsIt = ++tempIsItBegin;
538 }
539 }
540
541 break;
542 }
543 // for ALUGrid and UGGrid
544 case GridTypeIndices::aluGrid:
545 case GridTypeIndices::ugGrid:
546 {
547 if (nextIsIt == isEndIt)
548 nextIsIt = isItBegin;
549
550 break;
551 }
552 default:
553 {
554 DUNE_THROW(Dune::NotImplemented, "GridType can not be used with MPFAO implementation!");
555 break;
556 }
557 }
558
559 return *nextIsIt;
560}
561
562// TODO doc me!
563template<class TypeTag>
564void FvMpfaO2dPressure2p<TypeTag>::initializeMatrix()
565{
566 // determine matrix row sizes
567 for (const auto& element : elements(problem_.gridView()))
568 {
569 // cell index
570 int eIdxGlobalI = problem_.variables().index(element);
571
572 // initialize row size
573 int rowSize = 1;
574
575 // run through all intersections with neighbors
576 const auto isEndIt = problem_.gridView().iend(element);
577 for (auto isIt = problem_.gridView().ibegin(element); isIt != isEndIt; ++isIt)
578 {
579 const auto& intersection = *isIt;
580 auto nextIntersection = getNextIntersection_(element, isIt);
581
582 if (intersection.neighbor())
583 rowSize++;
584
585 // try to find the common neighbor of intersection's and nextIntersection's outside
586 if (intersection.neighbor() && nextIntersection.neighbor())
587 {
588 for (const auto& innerIntersection
589 : intersections(problem_.gridView(), intersection.outside()))
590 for (const auto& innerNextIntersection
591 : intersections(problem_.gridView(), nextIntersection.outside()))
592 {
593 if (innerIntersection.neighbor() && innerNextIntersection.neighbor())
594 {
595 if (innerIntersection.outside() == innerNextIntersection.outside()
596 && innerIntersection.outside() != intersection.inside())
597 {
598 rowSize++;
599 }
600 }
601 }
602 }
603 } // end of intersection loop
604
605 // set number of indices in row eIdxGlobalI to rowSize
606 this->A_.setrowsize(eIdxGlobalI, rowSize);
607
608 } // end of element loop
609
610 // indicate that size of all rows is defined
611 this->A_.endrowsizes();
612
613 // determine position of matrix entries
614 for (const auto& element : elements(problem_.gridView()))
615 {
616 // cell index
617 int eIdxGlobalI = problem_.variables().index(element);
618
619 // add diagonal index
620 this->A_.addindex(eIdxGlobalI, eIdxGlobalI);
621
622 // run through all intersections with neighbors
623 const auto isEndIt = problem_.gridView().iend(element);
624 for (auto isIt = problem_.gridView().ibegin(element); isIt != isEndIt; ++isIt)
625 {
626 const auto& intersection = *isIt;
627 auto nextIntersection = getNextIntersection_(element, isIt);
628
629 if (intersection.neighbor())
630 {
631 // access neighbor
632 int eIdxGlobalJ = problem_.variables().index(intersection.outside());
633
634 // add off diagonal index
635 // add index (row,col) to the matrix
636 this->A_.addindex(eIdxGlobalI, eIdxGlobalJ);
637 }
638
639 if (intersection.neighbor() && nextIntersection.neighbor())
640 {
641 for (const auto& innerIntersection
642 : intersections(problem_.gridView(), intersection.outside()))
643 for (const auto& innerNextIntersection
644 : intersections(problem_.gridView(), nextIntersection.outside()))
645 {
646 if (innerIntersection.neighbor() && innerNextIntersection.neighbor())
647 {
648 auto innerOutside = innerIntersection.outside();
649
650 if (innerOutside == innerNextIntersection.outside()
651 && innerOutside != intersection.inside())
652 {
653 int eIdxGlobalJ = problem_.variables().index(innerOutside);
654
655 this->A_.addindex(eIdxGlobalI, eIdxGlobalJ);
656 }
657 }
658 }
659 }
660 } // end of intersection loop
661 } // end of element loop
662
663 // indicate that all indices are defined, check consistency
664 this->A_.endindices();
665
666 return;
667}
668
669// Indices used in a interaction volume of the MPFA-o method
670// ___________________________________________________
671// | | |
672// | nuxy: cell geometry | nxy: face normal |
673// | vectors (see MPFA) | |
674// | | |
675// | 4-----------3-----------3 |
676// | | --> nu43 | nu34 <-- | |
677// | | |nu41 1|--> n43 ||nu32 |
678// | | v ^ |0 ^ v| |
679// |____________4__0__|n14__|__n23_|_1__2____________|
680// | | 1 0 | 0 |
681// | | ^ |1 nu23 ^ | |
682// | | |nu14 0|--> n12 | | |
683// | | -->nu12 | nu21<-- | |
684// | 1-----------1-----------2 |
685// | elementnumber |inter- |
686// | |face- |
687// | |number |
688// |________________________|________________________|
689
690// only for 2-D general quadrilateral
691// TODO doc me!
692template<class TypeTag>
693void FvMpfaO2dPressure2p<TypeTag>::storeInteractionVolumeInfo()
694{
695 // introduce matrix R for vector rotation and R is initialized as zero matrix
696 DimMatrix R(0);
697
698 // evaluate matrix R
699 if (dim == 2)
700 {
701 R[0][1] = 1;
702 R[1][0] = -1;
703 }
704
705 BoundaryTypes bcType;
706
707 // run through all elements
708 for (const auto& element : elements(problem_.gridView()))
709 {
710 // get common geometry information for the following computation
711
712 int eIdxGlobal1 = problem_.variables().index(element);
713 // get global coordinate of cell 1 center
714 const GlobalPosition& globalPos1 = element.geometry().center();
715
716 // get absolute permeability of neighbor cell 1
717 DimMatrix K1(problem_.spatialParams().intrinsicPermeability(element));
718
719 const auto isIt12End = problem_.gridView().iend(element);
720 for (auto isIt12 = problem_.gridView().ibegin(element); isIt12 != isIt12End; ++isIt12)
721 {
722 const auto& intersection12 = *isIt12;
723 auto intersection14 = getNextIntersection_(element, isIt12);
724
725 int indexInInside12 = intersection12.indexInInside();
726 int indexInInside14 = intersection14.indexInInside();
727
728 // get the intersection node /bar^{x_3} between 'isIt12' and 'isIt14', denoted as 'corner1234'
729 // initialization of corner1234
730
731 const auto refElement = referenceElement(element);
732
733 GlobalPosition corner1234(0);
734
735 int globalVertIdx1234 = 0;
736
737 // get the global coordinate and global vertex index of corner1234
738 for (int i = 0; i < intersection12.geometry().corners(); ++i)
739 {
740 bool finished = false;
741
742 const GlobalPosition& isIt12corner = intersection12.geometry().corner(i);
743
744 int localVertIdx12corner = refElement.subEntity(indexInInside12, dim - 1, i, dim);
745
746 int globalVertIdx12corner = problem_.variables().index(element.template subEntity<dim>(localVertIdx12corner));
747
748 for (int j = 0; j < intersection14.geometry().corners(); ++j)
749 {
750 int localVertIdx14corner = refElement.subEntity(indexInInside14, dim - 1, j, dim);
751
752 int globalVertIdx14corner = problem_.variables().index(element.template subEntity<dim>(localVertIdx14corner));
753
754 if (globalVertIdx12corner == globalVertIdx14corner)
755 {
756 corner1234 = isIt12corner;
757
758 globalVertIdx1234 = globalVertIdx12corner;
759
760 finished = true;
761 break;
762 }
763 }
764
765 if (finished)
766 {
767 break;
768 }
769 }
770
771 if (interactionVolumes_[globalVertIdx1234].isStored())
772 {
773 continue;
774 }
775 else
776 {
777 interactionVolumes_[globalVertIdx1234].setStored();
778 // std::cout << "vertIdx = " << globalVertIdx1234 << "\n";
779 }
780
781 //store pointer 1
782 interactionVolumes_[globalVertIdx1234].setSubVolumeElement(element, 0);
783 interactionVolumes_[globalVertIdx1234].setIndexOnElement(intersection12.indexInInside(), 0, 0);
784 interactionVolumes_[globalVertIdx1234].setIndexOnElement(intersection14.indexInInside(), 0, 1);
785
786 // center of face in global coordinates, i.e., the midpoint of edge 'isIt12'
787 const GlobalPosition& globalPosFace12 = intersection12.geometry().center();
788
789 // get face volume
790 Scalar faceVol12 = intersection12.geometry().volume() / 2.0;
791
792 // get outer normal vector scaled with half volume of face 'isIt12'
793 DimVector unitOuterNormal12 = intersection12.centerUnitOuterNormal();
794
795 // center of face in global coordinates, i.e., the midpoint of edge 'isIt14'
796 const GlobalPosition& globalPosFace41 = intersection14.geometry().center();
797
798 // get face volume
799 Scalar faceVol41 = intersection14.geometry().volume() / 2.0;
800
801 // get outer normal vector scaled with half volume of face 'isIt14': for numbering of n see Aavatsmark, Eigestad
802 DimVector unitOuterNormal14 = intersection14.centerUnitOuterNormal();
803
804 // compute normal vectors nu14,nu12
805 DimVector nu14(0);
806 R.mv(globalPos1 - globalPosFace12, nu14);
807
808 DimVector nu12(0);
809 R.mv(globalPosFace41 - globalPos1, nu12);
810
811 interactionVolumes_[globalVertIdx1234].setPermTimesNu(nu12, K1, 0, 0);
812 interactionVolumes_[globalVertIdx1234].setPermTimesNu(nu14, K1, 0, 1);
813 interactionVolumes_[globalVertIdx1234].setNormal(unitOuterNormal12, 0, 0);
814 interactionVolumes_[globalVertIdx1234].setNormal(unitOuterNormal14, 0, 1);
815 interactionVolumes_[globalVertIdx1234].setFaceArea(faceVol12, 0, 0);
816 interactionVolumes_[globalVertIdx1234].setFaceArea(faceVol41, 0, 1);
817
818 // compute dF1, the area of quadrilateral made by normal vectors 'nu'
819 DimVector Rnu12(0);
820 R.umv(nu12, Rnu12);
821 using std::abs;
822 interactionVolumes_[globalVertIdx1234].setDF(abs(nu14 * Rnu12), 0);
823
824 // handle interior face
825 if (intersection12.neighbor())
826 {
827 // access neighbor cell 2 of 'isIt12'
828 auto element2 = intersection12.outside();
829
830 int eIdxGlobal2 = problem_.variables().index(element2);
831
832 //store pointer 2
833 interactionVolumes_[globalVertIdx1234].setSubVolumeElement(element2, 1);
834 interactionVolumes_[globalVertIdx1234].setIndexOnElement(intersection12.indexInOutside(), 1, 1);
835
836 interactionVolumes_[globalVertIdx1234].setNormal(unitOuterNormal12, 1, 1);
837 interactionVolumes_[globalVertIdx1234].setFaceArea(faceVol12, 1, 1);
838
839 // get global coordinate of neighbor cell 2 center
840 const GlobalPosition& globalPos2 = element2.geometry().center();
841
842 // get absolute permeability of neighbor cell 2
843 DimMatrix K2(problem_.spatialParams().intrinsicPermeability(element2));
844
845 // 'isIt14' is an interior face
846 if (intersection14.neighbor())
847 {
848 // neighbor cell 3
849 // access neighbor cell 3
850 auto element4 = intersection14.outside();
851
852 //store pointer 4
853 interactionVolumes_[globalVertIdx1234].setSubVolumeElement(element4, 3);
854 interactionVolumes_[globalVertIdx1234].setIndexOnElement(intersection14.indexInOutside(), 3, 0);
855
856 interactionVolumes_[globalVertIdx1234].setNormal(unitOuterNormal14, 3, 0);
857 interactionVolumes_[globalVertIdx1234].setFaceArea(faceVol41, 3, 0);
858
859 // get basic information of cell 1,2's neighbor cell 3,4
860
861 // get global coordinate of neighbor cell 4 center
862 const GlobalPosition& globalPos4 = element4.geometry().center();
863
864 // get absolute permeability of neighbor cell 2
865 DimMatrix K4(problem_.spatialParams().intrinsicPermeability(element4));
866
867 // cell 3
868 GlobalPosition globalPos3(0);
869
870 GlobalPosition globalPosFace23(0);
871 GlobalPosition globalPosFace34(0);
872
873 for (const auto& intersection2
874 : intersections(problem_.gridView(), element2))
875 {
876 bool finished = false;
877
878 for (const auto& intersection4
879 : intersections(problem_.gridView(), element4))
880 {
881 if (intersection2.neighbor() && intersection4.neighbor())
882 {
883 auto element32 = intersection2.outside();
884 auto element34 = intersection4.outside();
885
886 // find the common neighbor cell between cell 2 and cell 3, except cell 1
887 if (element32 == element34 && element32 != element)
888 {
889 //store pointer 3
890 interactionVolumes_[globalVertIdx1234].setSubVolumeElement(element32, 2);
891
892 interactionVolumes_[globalVertIdx1234].setIndexOnElement(intersection2.indexInInside(), 1,
893 0);
894 interactionVolumes_[globalVertIdx1234].setIndexOnElement(intersection2.indexInOutside(), 2,
895 1);
896 interactionVolumes_[globalVertIdx1234].setIndexOnElement(intersection4.indexInInside(), 3,
897 1);
898 interactionVolumes_[globalVertIdx1234].setIndexOnElement(intersection4.indexInOutside(), 2,
899 0);
900
901 // get global coordinate of neighbor cell 4 center
902 globalPos3 = element32.geometry().center();
903
904 globalPosFace23 = intersection2.geometry().center();
905 globalPosFace34 = intersection4.geometry().center();
906
907 Scalar faceVol23 = intersection2.geometry().volume() / 2.0;
908 Scalar faceVol34 = intersection4.geometry().volume() / 2.0;
909
910 // get outer normal vector scaled with half volume of face : for numbering of n see Aavatsmark, Eigestad
911 DimVector unitOuterNormal23 = intersection2.centerUnitOuterNormal();
912
913 DimVector unitOuterNormal43 = intersection4.centerUnitOuterNormal();
914
915 interactionVolumes_[globalVertIdx1234].setNormal(unitOuterNormal23, 1, 0);
916 interactionVolumes_[globalVertIdx1234].setNormal(unitOuterNormal23, 2, 1);
917 interactionVolumes_[globalVertIdx1234].setNormal(unitOuterNormal43, 2, 0);
918 interactionVolumes_[globalVertIdx1234].setNormal(unitOuterNormal43, 3, 1);
919 interactionVolumes_[globalVertIdx1234].setFaceArea(faceVol23, 1, 0);
920 interactionVolumes_[globalVertIdx1234].setFaceArea(faceVol23, 2, 1);
921 interactionVolumes_[globalVertIdx1234].setFaceArea(faceVol34, 2, 0);
922 interactionVolumes_[globalVertIdx1234].setFaceArea(faceVol34, 3, 1);
923
924 // get absolute permeability of neighbor cell 2
925 DimMatrix K3(
926 problem_.spatialParams().intrinsicPermeability(element32));
927
928 // compute normal vectors nu23, nu21; nu32, nu34; nu41, nu43;
929 DimVector nu23(0);
930 R.umv(globalPosFace12 - globalPos2, nu23);
931
932 DimVector nu21(0);
933 R.umv(globalPosFace23 - globalPos2, nu21);
934
935 DimVector nu32(0);
936 R.umv(globalPosFace34 - globalPos3, nu32);
937
938 DimVector nu34(0);
939 R.umv(globalPos3 - globalPosFace23, nu34);
940
941 DimVector nu41(0);
942 R.umv(globalPos4 - globalPosFace34, nu41);
943
944 DimVector nu43(0);
945 R.umv(globalPos4 - globalPosFace41, nu43);
946
947 interactionVolumes_[globalVertIdx1234].setPermTimesNu(nu23, K2, 1, 0);
948 interactionVolumes_[globalVertIdx1234].setPermTimesNu(nu21, K2, 1, 1);
949 interactionVolumes_[globalVertIdx1234].setPermTimesNu(nu34, K3, 2, 0);
950 interactionVolumes_[globalVertIdx1234].setPermTimesNu(nu32, K3, 2, 1);
951 interactionVolumes_[globalVertIdx1234].setPermTimesNu(nu41, K4, 3, 0);
952 interactionVolumes_[globalVertIdx1234].setPermTimesNu(nu43, K4, 3, 1);
953
954 // compute dF2, dF3, dF4 i.e., the area of quadrilateral made by normal vectors 'nu'
955 DimVector Rnu21(0);
956 R.umv(nu21, Rnu21);
957 interactionVolumes_[globalVertIdx1234].setDF(abs(nu23 * Rnu21), 1);
958
959 DimVector Rnu34(0);
960 R.umv(nu34, Rnu34);
961 interactionVolumes_[globalVertIdx1234].setDF(abs(nu32 * Rnu34), 2);
962
963 DimVector Rnu43(0);
964 R.umv(nu43, Rnu43);
965 interactionVolumes_[globalVertIdx1234].setDF(abs(nu41 * Rnu43), 3);
966
967 finished = true;
968
969 break;
970 }
971 }
972 }
973 if (finished)
974 {
975 break;
976 }
977 }
978 }
979
980 // 'isIt14' is on the boundary
981 else
982 {
983 problem_.boundaryTypes(bcType, intersection14);
984 PrimaryVariables boundValues(0.0);
985
986 interactionVolumes_[globalVertIdx1234].setBoundary(bcType, 3);
987 if (bcType.isNeumann(pressEqIdx))
988 {
989 problem_.neumann(boundValues, intersection14);
990 boundValues *= faceVol41;
991 interactionVolumes_[globalVertIdx1234].setNeumannCondition(boundValues, 3);
992 }
993 if (bcType.hasDirichlet())
994 {
995 problem_.dirichlet(boundValues, intersection14);
996 interactionVolumes_[globalVertIdx1234].setDirichletCondition(boundValues, 3);
997 }
998
999 // get common geometry information for the following computation
1000 // get the information of the face 'isIt23' between cell2 and cell4 (locally numbered)
1001
1002 // center of face in global coordinates, i.e., the midpoint of edge 'isIt23'
1003 GlobalPosition globalPosFace23(0);
1004
1005 // get face volume
1006 Scalar faceVol23 = 0;
1007
1008 // get outer normal vector scaled with half volume of face 'isIt23'
1009 DimVector unitOuterNormal23(0);
1010
1011 bool finished = false;
1012
1013 for (const auto& intersection2
1014 : intersections(problem_.gridView(), element2))
1015 {
1016 if (intersection2.boundary())
1017 {
1018 for (int i = 0; i < intersection2.geometry().corners(); ++i)
1019 {
1020 int localVertIdx2corner = refElement.subEntity(intersection2.indexInInside(), dim - 1, i,
1021 dim);
1022
1023 int globalVertIdx2corner = problem_.variables().index(element2.template subEntity<dim>(localVertIdx2corner));
1024
1025 if (globalVertIdx2corner == globalVertIdx1234)
1026 {
1027 interactionVolumes_[globalVertIdx1234].setIndexOnElement(intersection2.indexInInside(), 1,
1028 0);
1029
1030 globalPosFace23 = intersection2.geometry().center();
1031
1032 faceVol23 = intersection2.geometry().volume() / 2.0;
1033
1034 unitOuterNormal23 = intersection2.centerUnitOuterNormal();
1035
1036 interactionVolumes_[globalVertIdx1234].setNormal(unitOuterNormal23, 1, 0);
1037 interactionVolumes_[globalVertIdx1234].setFaceArea(faceVol23, 1, 0);
1038
1039 problem_.boundaryTypes(bcType, intersection2);
1040 boundValues = 0.0;
1041
1042 interactionVolumes_[globalVertIdx1234].setBoundary(bcType, 1);
1043 if (bcType.isNeumann(pressEqIdx))
1044 {
1045 problem_.neumann(boundValues, intersection2);
1046 boundValues *= faceVol23;
1047 interactionVolumes_[globalVertIdx1234].setNeumannCondition(boundValues, 1);
1048 }
1049 if (bcType.hasDirichlet())
1050 {
1051 problem_.dirichlet(boundValues, intersection2);
1052 interactionVolumes_[globalVertIdx1234].setDirichletCondition(boundValues, 1);
1053 }
1054
1055 interactionVolumes_[globalVertIdx1234].setOutsideFace(2);
1056
1057 innerBoundaryVolumeFaces_[eIdxGlobal1][intersection12.indexInInside()] = true;
1058 innerBoundaryVolumeFaces_[eIdxGlobal2][intersection12.indexInOutside()] = true;
1059
1060 // compute normal vectors nu23, nu21;
1061 DimVector nu23(0);
1062 R.umv(globalPosFace12 - globalPos2, nu23);
1063
1064 DimVector nu21(0);
1065 R.umv(globalPosFace23 - globalPos2, nu21);
1066
1067 interactionVolumes_[globalVertIdx1234].setPermTimesNu(nu23, K2, 1, 0);
1068 interactionVolumes_[globalVertIdx1234].setPermTimesNu(nu21, K2, 1, 1);
1069
1070 // compute dF2 i.e., the area of quadrilateral made by normal vectors 'nu'
1071 DimVector Rnu21(0);
1072 R.umv(nu21, Rnu21);
1073 interactionVolumes_[globalVertIdx1234].setDF(abs(nu23 * Rnu21), 1);
1074
1075 finished = true;
1076
1077 break;
1078 }
1079 }
1080 }
1081 if (finished)
1082 {
1083 break;
1084 }
1085 }
1086 if (!finished)
1087 {
1088 DUNE_THROW(
1089 Dune::NotImplemented,
1090 "fvmpfao2pfaboundpressure2p.hh, l. 997: boundary shape not available as interaction volume shape");
1091 }
1092 }
1093 }
1094
1095 // handle boundary face 'isIt12'
1096 else
1097 {
1098 problem_.boundaryTypes(bcType, *isIt12);
1099 PrimaryVariables boundValues(0.0);
1100
1101 interactionVolumes_[globalVertIdx1234].setBoundary(bcType, 0);
1102 if (bcType.isNeumann(pressEqIdx))
1103 {
1104 problem_.neumann(boundValues, *isIt12);
1105 boundValues *= faceVol12;
1106 interactionVolumes_[globalVertIdx1234].setNeumannCondition(boundValues, 0);
1107 }
1108 if (bcType.hasDirichlet())
1109 {
1110 problem_.dirichlet(boundValues, *isIt12);
1111 interactionVolumes_[globalVertIdx1234].setDirichletCondition(boundValues, 0);
1112 }
1113
1114 // 'isIt14' is on boundary
1115 if (intersection14.boundary())
1116 {
1117 problem_.boundaryTypes(bcType, intersection14);
1118 boundValues = 0.0;
1119
1120 interactionVolumes_[globalVertIdx1234].setBoundary(bcType, 3);
1121 if (bcType.isNeumann(pressEqIdx))
1122 {
1123 problem_.neumann(boundValues, intersection14);
1124 boundValues *= faceVol41;
1125 interactionVolumes_[globalVertIdx1234].setNeumannCondition(boundValues, 3);
1126 }
1127 if (bcType.hasDirichlet())
1128 {
1129 problem_.dirichlet(boundValues, intersection14);
1130 interactionVolumes_[globalVertIdx1234].setDirichletCondition(boundValues, 3);
1131 }
1132
1133 interactionVolumes_[globalVertIdx1234].setOutsideFace(1);
1134 interactionVolumes_[globalVertIdx1234].setOutsideFace(2);
1135 }
1136
1137 // 'isIt14' is inside
1138 else
1139 {
1140 // neighbor cell 3
1141 // access neighbor cell 3
1142 auto element4 = intersection14.outside();
1143 interactionVolumes_[globalVertIdx1234].setIndexOnElement(intersection14.indexInOutside(), 3, 0);
1144
1145 //store pointer 4
1146 interactionVolumes_[globalVertIdx1234].setSubVolumeElement(element4, 3);
1147
1148 interactionVolumes_[globalVertIdx1234].setNormal(unitOuterNormal14, 3, 0);
1149 interactionVolumes_[globalVertIdx1234].setFaceArea(faceVol41, 3, 0);
1150
1151 // get global coordinate of neighbor cell 3 center
1152 const GlobalPosition& globalPos4 = element4.geometry().center();
1153
1154 int eIdxGlobal4 = problem_.variables().index(element4);
1155
1156 bool finished = false;
1157
1158 // get the information of the face 'isIt34' between cell3 and cell4 (locally numbered)
1159 for (const auto& intersection4
1160 : intersections(problem_.gridView(), element4))
1161 {
1162 if (intersection4.boundary())
1163 {
1164 for (int i = 0; i < intersection4.geometry().corners(); ++i)
1165 {
1166 int localVertIdx4corner = refElement.subEntity(intersection4.indexInInside(), dim - 1, i,
1167 dim);
1168
1169 int globalVertIdx4corner = problem_.variables().index(element4.template subEntity<dim>(localVertIdx4corner));
1170
1171 if (globalVertIdx4corner == globalVertIdx1234)
1172 {
1173 interactionVolumes_[globalVertIdx1234].setIndexOnElement(intersection4.indexInInside(), 3,
1174 1);
1175
1176 const GlobalPosition& globalPosFace34 = intersection4.geometry().center();
1177
1178 Scalar faceVol34 = intersection4.geometry().volume() / 2.0;
1179
1180 DimVector unitOuterNormal43 = intersection4.centerUnitOuterNormal();
1181
1182 interactionVolumes_[globalVertIdx1234].setNormal(unitOuterNormal43, 3, 1);
1183 interactionVolumes_[globalVertIdx1234].setFaceArea(faceVol34, 3, 1);
1184
1185 problem_.boundaryTypes(bcType, intersection4);
1186 boundValues = 0.0;
1187
1188 interactionVolumes_[globalVertIdx1234].setBoundary(bcType, 2);
1189 if (bcType.isNeumann(pressEqIdx))
1190 {
1191 problem_.neumann(boundValues, intersection4);
1192 boundValues *= faceVol34;
1193 interactionVolumes_[globalVertIdx1234].setNeumannCondition(boundValues, 2);
1194 }
1195 if (bcType.hasDirichlet())
1196 {
1197 problem_.dirichlet(boundValues, intersection4);
1198 interactionVolumes_[globalVertIdx1234].setDirichletCondition(boundValues, 2);
1199 }
1200
1201 interactionVolumes_[globalVertIdx1234].setOutsideFace(1);
1202
1203 innerBoundaryVolumeFaces_[eIdxGlobal1][intersection14.indexInInside()] = true;
1204 innerBoundaryVolumeFaces_[eIdxGlobal4][intersection14.indexInOutside()] = true;
1205
1206 // get absolute permeability of neighbor cell 2
1207 DimMatrix K4(
1208 problem_.spatialParams().intrinsicPermeability(element4));
1209
1210 // compute normal vectors nu41, nu43;
1211 DimVector nu41(0);
1212 R.umv(globalPos4 - globalPosFace34, nu41);
1213
1214 DimVector nu43(0);
1215 R.umv(globalPos4 - globalPosFace41, nu43);
1216
1217 interactionVolumes_[globalVertIdx1234].setPermTimesNu(nu41, K4, 3, 0);
1218 interactionVolumes_[globalVertIdx1234].setPermTimesNu(nu43, K4, 3, 1);
1219
1220 // compute dF1, dF3 i.e., the area of quadrilateral made by normal vectors 'nu'
1221 DimVector Rnu43(0);
1222 R.umv(nu43, Rnu43);
1223 interactionVolumes_[globalVertIdx1234].setDF(abs(nu41 * Rnu43), 3);
1224
1225 finished = true;
1226
1227 break;
1228 }
1229 }
1230 }
1231 if (finished)
1232 {
1233 break;
1234 }
1235 }
1236 if (!finished)
1237 {
1238 DUNE_THROW(
1239 Dune::NotImplemented,
1240 "fvmpfao2pfaboundpressure2p.hh, l. 1164: boundary shape not available as interaction volume shape");
1241 }
1242 }
1243 }
1244
1245 } // end all intersections
1246 } // end grid traversal
1247
1248 return;
1249}
1250
1251// only for 2-D general quadrilateral
1252//TODO doc me!
1253template<class TypeTag>
1254void FvMpfaO2dPressure2p<TypeTag>::assemble()
1255{
1256 // initialization: set global matrix this->A_ to zero
1257 this->A_ = 0;
1258 this->f_ = 0;
1259
1260 // run through all vertices
1261 for (const auto& vertex : vertices(problem_.gridView()))
1262 {
1263 int vIdxGlobal = problem_.variables().index(vertex);
1264
1265 InteractionVolume& interactionVolume = interactionVolumes_[vIdxGlobal];
1266
1267 if (interactionVolume.isInnerVolume())
1268 {
1269
1270 auto element1 = interactionVolume.getSubVolumeElement(0);
1271 auto element2 = interactionVolume.getSubVolumeElement(1);
1272 auto element3 = interactionVolume.getSubVolumeElement(2);
1273 auto element4 = interactionVolume.getSubVolumeElement(3);
1274
1275 // get global coordinate of cell centers
1276 const GlobalPosition& globalPos1 = element1.geometry().center();
1277 const GlobalPosition& globalPos2 = element2.geometry().center();
1278 const GlobalPosition& globalPos3 = element3.geometry().center();
1279 const GlobalPosition& globalPos4 = element4.geometry().center();
1280
1281 // cell volumes
1282 Scalar volume1 = element1.geometry().volume();
1283 Scalar volume2 = element2.geometry().volume();
1284 Scalar volume3 = element3.geometry().volume();
1285 Scalar volume4 = element4.geometry().volume();
1286
1287 // cell index
1288 int eIdxGlobal1 = problem_.variables().index(element1);
1289 int eIdxGlobal2 = problem_.variables().index(element2);
1290 int eIdxGlobal3 = problem_.variables().index(element3);
1291 int eIdxGlobal4 = problem_.variables().index(element4);
1292
1293 //get the cell Data
1294 CellData& cellData1 = problem_.variables().cellData(eIdxGlobal1);
1295 CellData& cellData2 = problem_.variables().cellData(eIdxGlobal2);
1296 CellData& cellData3 = problem_.variables().cellData(eIdxGlobal3);
1297 CellData& cellData4 = problem_.variables().cellData(eIdxGlobal4);
1298
1299 // evaluate right hand side
1300 PrimaryVariables source(0.0);
1301 problem_.source(source, element1);
1302 this->f_[eIdxGlobal1] += volume1 / (4.0)
1303 * (source[wPhaseIdx] / density_[wPhaseIdx] + source[nPhaseIdx] / density_[nPhaseIdx]);
1304 problem_.source(source, element2);
1305 this->f_[eIdxGlobal2] += volume2 / (4.0)
1306 * (source[wPhaseIdx] / density_[wPhaseIdx] + source[nPhaseIdx] / density_[nPhaseIdx]);
1307 problem_.source(source, element3);
1308 this->f_[eIdxGlobal3] += volume3 / (4.0)
1309 * (source[wPhaseIdx] / density_[wPhaseIdx] + source[nPhaseIdx] / density_[nPhaseIdx]);
1310 problem_.source(source, element4);
1311 this->f_[eIdxGlobal4] += volume4 / (4.0)
1312 * (source[wPhaseIdx] / density_[wPhaseIdx] + source[nPhaseIdx] / density_[nPhaseIdx]);
1313
1314 this->f_[eIdxGlobal1] += evaluateErrorTerm_(cellData1) * volume1 / (4.0);
1315 this->f_[eIdxGlobal2] += evaluateErrorTerm_(cellData2) * volume2 / (4.0);
1316 this->f_[eIdxGlobal3] += evaluateErrorTerm_(cellData3) * volume3 / (4.0);
1317 this->f_[eIdxGlobal4] += evaluateErrorTerm_(cellData4) * volume4 / (4.0);
1318
1319 //get mobilities of the phases
1320 Dune::FieldVector<Scalar, numPhases> lambda1(cellData1.mobility(wPhaseIdx));
1321 lambda1[nPhaseIdx] = cellData1.mobility(nPhaseIdx);
1322
1323 //compute total mobility of cell 1
1324 Scalar lambdaTotal1 = lambda1[wPhaseIdx] + lambda1[nPhaseIdx];
1325
1326 //get mobilities of the phases
1327 Dune::FieldVector<Scalar, numPhases> lambda2(cellData2.mobility(wPhaseIdx));
1328 lambda2[nPhaseIdx] = cellData2.mobility(nPhaseIdx);
1329
1330 //compute total mobility of cell 1
1331 Scalar lambdaTotal2 = lambda2[wPhaseIdx] + lambda2[nPhaseIdx];
1332
1333 //get mobilities of the phases
1334 Dune::FieldVector<Scalar, numPhases> lambda3(cellData3.mobility(wPhaseIdx));
1335 lambda3[nPhaseIdx] = cellData3.mobility(nPhaseIdx);
1336
1337 //compute total mobility of cell 1
1338 Scalar lambdaTotal3 = lambda3[wPhaseIdx] + lambda3[nPhaseIdx];
1339
1340 //get mobilities of the phases
1341 Dune::FieldVector<Scalar, numPhases> lambda4(cellData4.mobility(wPhaseIdx));
1342 lambda4[nPhaseIdx] = cellData4.mobility(nPhaseIdx);
1343
1344 //compute total mobility of cell 1
1345 Scalar lambdaTotal4 = lambda4[wPhaseIdx] + lambda4[nPhaseIdx];
1346
1347 Scalar gn12nu14 = interactionVolume.getNtkrkNu_df(lambdaTotal1, 0, 0, 1);
1348 Scalar gn12nu12 = interactionVolume.getNtkrkNu_df(lambdaTotal1, 0, 0, 0);
1349 Scalar gn14nu14 = interactionVolume.getNtkrkNu_df(lambdaTotal1, 0, 1, 1);
1350 Scalar gn14nu12 = interactionVolume.getNtkrkNu_df(lambdaTotal1, 0, 1, 0);
1351 Scalar gn12nu23 = interactionVolume.getNtkrkNu_df(lambdaTotal2, 1, 1, 0);
1352 Scalar gn12nu21 = interactionVolume.getNtkrkNu_df(lambdaTotal2, 1, 1, 1);
1353 Scalar gn23nu23 = interactionVolume.getNtkrkNu_df(lambdaTotal2, 1, 0, 0);
1354 Scalar gn23nu21 = interactionVolume.getNtkrkNu_df(lambdaTotal2, 1, 0, 1);
1355 Scalar gn43nu32 = interactionVolume.getNtkrkNu_df(lambdaTotal3, 2, 0, 1);
1356 Scalar gn43nu34 = interactionVolume.getNtkrkNu_df(lambdaTotal3, 2, 0, 0);
1357 Scalar gn23nu32 = interactionVolume.getNtkrkNu_df(lambdaTotal3, 2, 1, 1);
1358 Scalar gn23nu34 = interactionVolume.getNtkrkNu_df(lambdaTotal3, 2, 1, 0);
1359 Scalar gn43nu41 = interactionVolume.getNtkrkNu_df(lambdaTotal4, 3, 1, 0);
1360 Scalar gn43nu43 = interactionVolume.getNtkrkNu_df(lambdaTotal4, 3, 1, 1);
1361 Scalar gn14nu41 = interactionVolume.getNtkrkNu_df(lambdaTotal4, 3, 0, 0);
1362 Scalar gn14nu43 = interactionVolume.getNtkrkNu_df(lambdaTotal4, 3, 0, 1);
1363
1364 // compute transmissibility matrix T = CA^{-1}B+F
1365 Dune::FieldMatrix<Scalar, 2 * dim, 2 * dim> C(0), F(0), A(0), B(0);
1366
1367 // evaluate matrix C, F, A, B
1368 C[0][0] = -gn12nu12;
1369 C[0][3] = -gn12nu14;
1370 C[1][0] = gn23nu21;
1371 C[1][1] = -gn23nu23;
1372 C[2][1] = gn43nu32;
1373 C[2][2] = gn43nu34;
1374 C[3][2] = -gn14nu43;
1375 C[3][3] = gn14nu41;
1376
1377 F[0][0] = gn12nu12 + gn12nu14;
1378 F[1][1] = -gn23nu21 + gn23nu23;
1379 F[2][2] = -gn43nu34 - gn43nu32;
1380 F[3][3] = gn14nu43 - gn14nu41;
1381
1382 A[0][0] = gn12nu12 + gn12nu21;
1383 A[0][1] = -gn12nu23;
1384 A[0][3] = gn12nu14;
1385 A[1][0] = -gn23nu21;
1386 A[1][1] = gn23nu23 + gn23nu32;
1387 A[1][2] = gn23nu34;
1388 A[2][1] = -gn43nu32;
1389 A[2][2] = -gn43nu34 - gn43nu43;
1390 A[2][3] = gn43nu41;
1391 A[3][0] = -gn14nu12;
1392 A[3][2] = gn14nu43;
1393 A[3][3] = -gn14nu41 - gn14nu14;
1394
1395 B[0][0] = gn12nu12 + gn12nu14;
1396 B[0][1] = gn12nu21 - gn12nu23;
1397 B[1][1] = -gn23nu21 + gn23nu23;
1398 B[1][2] = gn23nu34 + gn23nu32;
1399 B[2][2] = -gn43nu34 - gn43nu32;
1400 B[2][3] = -gn43nu43 + gn43nu41;
1401 B[3][0] = -gn14nu12 - gn14nu14;
1402 B[3][3] = gn14nu43 - gn14nu41;
1403
1404 // compute T
1405 A.invert();
1406 F += C.rightmultiply(B.leftmultiply(A));
1407 Dune::FieldMatrix<Scalar, 2 * dim, 2 * dim> T(F);
1408
1409 // std::cout<<"Tpress = "<<T<<"\n";
1410
1411 // assemble the global matrix this->A_ and right hand side f
1412 this->A_[eIdxGlobal1][eIdxGlobal1] += T[0][0] + T[3][0];
1413 this->A_[eIdxGlobal1][eIdxGlobal2] += T[0][1] + T[3][1];
1414 this->A_[eIdxGlobal1][eIdxGlobal3] += T[0][2] + T[3][2];
1415 this->A_[eIdxGlobal1][eIdxGlobal4] += T[0][3] + T[3][3];
1416
1417 this->A_[eIdxGlobal2][eIdxGlobal1] += -T[0][0] + T[1][0];
1418 this->A_[eIdxGlobal2][eIdxGlobal2] += -T[0][1] + T[1][1];
1419 this->A_[eIdxGlobal2][eIdxGlobal3] += -T[0][2] + T[1][2];
1420 this->A_[eIdxGlobal2][eIdxGlobal4] += -T[0][3] + T[1][3];
1421
1422 this->A_[eIdxGlobal3][eIdxGlobal1] -= T[1][0] + T[2][0];
1423 this->A_[eIdxGlobal3][eIdxGlobal2] -= T[1][1] + T[2][1];
1424 this->A_[eIdxGlobal3][eIdxGlobal3] -= T[1][2] + T[2][2];
1425 this->A_[eIdxGlobal3][eIdxGlobal4] -= T[1][3] + T[2][3];
1426
1427 this->A_[eIdxGlobal4][eIdxGlobal1] += T[2][0] - T[3][0];
1428 this->A_[eIdxGlobal4][eIdxGlobal2] += T[2][1] - T[3][1];
1429 this->A_[eIdxGlobal4][eIdxGlobal3] += T[2][2] - T[3][2];
1430 this->A_[eIdxGlobal4][eIdxGlobal4] += T[2][3] - T[3][3];
1431
1432 if (innerBoundaryVolumeFaces_[eIdxGlobal1][interactionVolume.getIndexOnElement(0, 0)])
1433 {
1434 this->A_[eIdxGlobal1][eIdxGlobal1] += T[0][0];
1435 this->A_[eIdxGlobal1][eIdxGlobal2] += T[0][1];
1436 this->A_[eIdxGlobal1][eIdxGlobal3] += T[0][2];
1437 this->A_[eIdxGlobal1][eIdxGlobal4] += T[0][3];
1438 }
1439 if (innerBoundaryVolumeFaces_[eIdxGlobal1][interactionVolume.getIndexOnElement(0, 1)])
1440 {
1441 this->A_[eIdxGlobal1][eIdxGlobal1] += T[3][0];
1442 this->A_[eIdxGlobal1][eIdxGlobal2] += T[3][1];
1443 this->A_[eIdxGlobal1][eIdxGlobal3] += T[3][2];
1444 this->A_[eIdxGlobal1][eIdxGlobal4] += T[3][3];
1445
1446 }
1447 if (innerBoundaryVolumeFaces_[eIdxGlobal2][interactionVolume.getIndexOnElement(1, 0)])
1448 {
1449 this->A_[eIdxGlobal2][eIdxGlobal1] += T[1][0];
1450 this->A_[eIdxGlobal2][eIdxGlobal2] += T[1][1];
1451 this->A_[eIdxGlobal2][eIdxGlobal3] += T[1][2];
1452 this->A_[eIdxGlobal2][eIdxGlobal4] += T[1][3];
1453 }
1454 if (innerBoundaryVolumeFaces_[eIdxGlobal2][interactionVolume.getIndexOnElement(1, 1)])
1455 {
1456 this->A_[eIdxGlobal2][eIdxGlobal1] += -T[0][0];
1457 this->A_[eIdxGlobal2][eIdxGlobal2] += -T[0][1];
1458 this->A_[eIdxGlobal2][eIdxGlobal3] += -T[0][2];
1459 this->A_[eIdxGlobal2][eIdxGlobal4] += -T[0][3];
1460 }
1461 if (innerBoundaryVolumeFaces_[eIdxGlobal3][interactionVolume.getIndexOnElement(2, 0)])
1462 {
1463 this->A_[eIdxGlobal3][eIdxGlobal1] -= T[2][0];
1464 this->A_[eIdxGlobal3][eIdxGlobal2] -= T[2][1];
1465 this->A_[eIdxGlobal3][eIdxGlobal3] -= T[2][2];
1466 this->A_[eIdxGlobal3][eIdxGlobal4] -= T[2][3];
1467 }
1468 if (innerBoundaryVolumeFaces_[eIdxGlobal3][interactionVolume.getIndexOnElement(2, 1)])
1469 {
1470 this->A_[eIdxGlobal3][eIdxGlobal1] -= T[1][0];
1471 this->A_[eIdxGlobal3][eIdxGlobal2] -= T[1][1];
1472 this->A_[eIdxGlobal3][eIdxGlobal3] -= T[1][2];
1473 this->A_[eIdxGlobal3][eIdxGlobal4] -= T[1][3];
1474 }
1475 if (innerBoundaryVolumeFaces_[eIdxGlobal4][interactionVolume.getIndexOnElement(3, 0)])
1476 {
1477 this->A_[eIdxGlobal4][eIdxGlobal1] += -T[3][0];
1478 this->A_[eIdxGlobal4][eIdxGlobal2] += -T[3][1];
1479 this->A_[eIdxGlobal4][eIdxGlobal3] += -T[3][2];
1480 this->A_[eIdxGlobal4][eIdxGlobal4] += -T[3][3];
1481 }
1482 if (innerBoundaryVolumeFaces_[eIdxGlobal4][interactionVolume.getIndexOnElement(3, 1)])
1483 {
1484 this->A_[eIdxGlobal4][eIdxGlobal1] += T[2][0];
1485 this->A_[eIdxGlobal4][eIdxGlobal2] += T[2][1];
1486 this->A_[eIdxGlobal4][eIdxGlobal3] += T[2][2];
1487 this->A_[eIdxGlobal4][eIdxGlobal4] += T[2][3];
1488 }
1489
1490 //add capillary pressure and gravity terms to right-hand-side
1491 //calculate capillary pressure velocity
1492 Dune::FieldVector<Scalar, 2 * dim> pc(0);
1493 pc[0] = cellData1.capillaryPressure();
1494 pc[1] = cellData2.capillaryPressure();
1495 pc[2] = cellData3.capillaryPressure();
1496 pc[3] = cellData4.capillaryPressure();
1497
1498 Dune::FieldVector<Scalar, 2 * dim> gravityDiff(0);
1499
1500// std::cout<<"maxPos = "<<problem_.bBoxMax()<<"\n";
1501
1502 gravityDiff[0] = (problem_.bBoxMax() - globalPos1) * gravity_ * (density_[nPhaseIdx] - density_[wPhaseIdx]);
1503 gravityDiff[1] = (problem_.bBoxMax() - globalPos2) * gravity_ * (density_[nPhaseIdx] - density_[wPhaseIdx]);
1504 gravityDiff[2] = (problem_.bBoxMax() - globalPos3) * gravity_ * (density_[nPhaseIdx] - density_[wPhaseIdx]);
1505 gravityDiff[3] = (problem_.bBoxMax() - globalPos4) * gravity_ * (density_[nPhaseIdx] - density_[wPhaseIdx]);
1506
1507 pc += gravityDiff;
1508
1509 if (pc[0] == 0 && pc[1] == 0 && pc[2] == 0 && pc[3] == 0)
1510 {
1511 continue;
1512 }
1513
1514 Dune::FieldVector<Scalar, 2 * dim> pcFlux(0);
1515
1516 T.mv(pc, pcFlux);
1517
1518 // std::cout<<"pcFlux = "<<pcFlux<<"\n";
1519
1520 Scalar pcPotential12 = pcFlux[0];
1521 Scalar pcPotential14 = pcFlux[3];
1522 Scalar pcPotential32 = -pcFlux[1];
1523 Scalar pcPotential34 = -pcFlux[2];
1524
1525 // std::cout<<"pcPotential12 = "<<pcPotential12<<"\n";
1526
1527 //compute mobilities of face 1
1528 Dune::FieldVector<Scalar, numPhases> lambda12Upw(0.0);
1529 lambda12Upw[wPhaseIdx] = (pcPotential12 >= 0) ? lambda1[wPhaseIdx] : lambda2[wPhaseIdx];
1530 lambda12Upw[nPhaseIdx] = (pcPotential12 >= 0) ? lambda1[nPhaseIdx] : lambda2[nPhaseIdx];
1531
1532 //compute mobilities of face 4
1533 Dune::FieldVector<Scalar, numPhases> lambda14Upw(0.0);
1534 lambda14Upw[wPhaseIdx] = (pcPotential14 >= 0) ? lambda1[wPhaseIdx] : lambda4[wPhaseIdx];
1535 lambda14Upw[nPhaseIdx] = (pcPotential14 >= 0) ? lambda1[nPhaseIdx] : lambda4[nPhaseIdx];
1536
1537 //compute mobilities of face 2
1538 Dune::FieldVector<Scalar, numPhases> lambda32Upw(0.0);
1539 lambda32Upw[wPhaseIdx] = (pcPotential32 >= 0) ? lambda3[wPhaseIdx] : lambda2[wPhaseIdx];
1540 lambda32Upw[nPhaseIdx] = (pcPotential32 >= 0) ? lambda3[nPhaseIdx] : lambda2[nPhaseIdx];
1541
1542 //compute mobilities of face 3
1543 Dune::FieldVector<Scalar, numPhases> lambda34Upw(0.0);
1544 lambda34Upw[wPhaseIdx] = (pcPotential34 >= 0) ? lambda3[wPhaseIdx] : lambda4[wPhaseIdx];
1545 lambda34Upw[nPhaseIdx] = (pcPotential34 >= 0) ? lambda3[nPhaseIdx] : lambda4[nPhaseIdx];
1546
1547 for (int i = 0; i < numPhases; i++)
1548 {
1549 Scalar lambdaT12 = lambda12Upw[wPhaseIdx] + lambda12Upw[nPhaseIdx];
1550 Scalar lambdaT14 = lambda14Upw[wPhaseIdx] + lambda14Upw[nPhaseIdx];
1551 Scalar lambdaT32 = lambda32Upw[wPhaseIdx] + lambda32Upw[nPhaseIdx];
1552 Scalar lambdaT34 = lambda34Upw[wPhaseIdx] + lambda34Upw[nPhaseIdx];
1553 Scalar fracFlow12 = (lambdaT12 > threshold_) ? lambda12Upw[i] / (lambdaT12) : 0.0;
1554 Scalar fracFlow14 = (lambdaT14 > threshold_) ? lambda14Upw[i] / (lambdaT14) : 0.0;
1555 Scalar fracFlow32 = (lambdaT32 > threshold_) ? lambda32Upw[i] / (lambdaT32) : 0.0;
1556 Scalar fracFlow34 = (lambdaT34 > threshold_) ? lambda34Upw[i] / (lambdaT34) : 0.0;
1557
1558 Dune::FieldVector<Scalar, 2 * dim> pcFluxReal(pcFlux);
1559
1560 pcFluxReal[0] *= fracFlow12;
1561 pcFluxReal[1] *= fracFlow32;
1562 pcFluxReal[2] *= fracFlow34;
1563 pcFluxReal[3] *= fracFlow14;
1564
1565// using std::isnan;
1566// if (isnan(pcFluxReal.two_norm()))
1567// std::cout<<"pcFlux = "<<pcFlux<<"\n";
1568
1569 switch (pressureType_)
1570 {
1571 case pw:
1572 {
1573 if (i == nPhaseIdx)
1574 {
1575 //add capillary pressure term to right hand side
1576 this->f_[eIdxGlobal1] -= (pcFluxReal[0] + pcFluxReal[3]);
1577 this->f_[eIdxGlobal2] -= (pcFluxReal[1] - pcFluxReal[0]);
1578 this->f_[eIdxGlobal3] -= (-pcFluxReal[2] - pcFluxReal[1]);
1579 this->f_[eIdxGlobal4] -= (-pcFluxReal[3] + pcFluxReal[2]);
1580
1581 if (innerBoundaryVolumeFaces_[eIdxGlobal1][interactionVolume.getIndexOnElement(0, 0)])
1582 {
1583 this->f_[eIdxGlobal1] -= pcFluxReal[0];
1584 }
1585 if (innerBoundaryVolumeFaces_[eIdxGlobal1][interactionVolume.getIndexOnElement(0, 1)])
1586 {
1587 this->f_[eIdxGlobal1] -= pcFluxReal[3];
1588 }
1589 if (innerBoundaryVolumeFaces_[eIdxGlobal2][interactionVolume.getIndexOnElement(1, 0)])
1590 {
1591 this->f_[eIdxGlobal2] -= pcFluxReal[1];
1592 }
1593 if (innerBoundaryVolumeFaces_[eIdxGlobal2][interactionVolume.getIndexOnElement(1, 1)])
1594 {
1595 this->f_[eIdxGlobal2] += pcFluxReal[0];
1596 }
1597 if (innerBoundaryVolumeFaces_[eIdxGlobal3][interactionVolume.getIndexOnElement(2, 0)])
1598 {
1599 this->f_[eIdxGlobal3] += pcFluxReal[2];
1600 }
1601 if (innerBoundaryVolumeFaces_[eIdxGlobal3][interactionVolume.getIndexOnElement(2, 1)])
1602 {
1603 this->f_[eIdxGlobal3] += pcFluxReal[1];
1604 }
1605 if (innerBoundaryVolumeFaces_[eIdxGlobal4][interactionVolume.getIndexOnElement(3, 0)])
1606 {
1607 this->f_[eIdxGlobal4] += pcFluxReal[3];
1608 }
1609 if (innerBoundaryVolumeFaces_[eIdxGlobal4][interactionVolume.getIndexOnElement(3, 1)])
1610 {
1611 this->f_[eIdxGlobal4] -= pcFluxReal[2];
1612 }
1613 }
1614 break;
1615 }
1616 case pn:
1617 {
1618 if (i == wPhaseIdx)
1619 {
1620 //add capillary pressure term to right hand side
1621 this->f_[eIdxGlobal1] += (pcFluxReal[0] + pcFluxReal[1]);
1622 this->f_[eIdxGlobal2] += (pcFluxReal[1] - pcFluxReal[0]);
1623 this->f_[eIdxGlobal3] += (-pcFluxReal[2] - pcFluxReal[1]);
1624 this->f_[eIdxGlobal4] += (-pcFluxReal[3] + pcFluxReal[2]);
1625
1626 if (innerBoundaryVolumeFaces_[eIdxGlobal1][interactionVolume.getIndexOnElement(0, 0)])
1627 {
1628 this->f_[eIdxGlobal1] += pcFluxReal[0];
1629 }
1630 if (innerBoundaryVolumeFaces_[eIdxGlobal1][interactionVolume.getIndexOnElement(0, 1)])
1631 {
1632 this->f_[eIdxGlobal1] += pcFluxReal[3];
1633 }
1634 if (innerBoundaryVolumeFaces_[eIdxGlobal2][interactionVolume.getIndexOnElement(1, 0)])
1635 {
1636 this->f_[eIdxGlobal2] += pcFluxReal[1];
1637 }
1638 if (innerBoundaryVolumeFaces_[eIdxGlobal2][interactionVolume.getIndexOnElement(1, 1)])
1639 {
1640 this->f_[eIdxGlobal2] -= pcFluxReal[0];
1641 }
1642 if (innerBoundaryVolumeFaces_[eIdxGlobal3][interactionVolume.getIndexOnElement(2, 0)])
1643 {
1644 this->f_[eIdxGlobal3] -= pcFluxReal[2];
1645 }
1646 if (innerBoundaryVolumeFaces_[eIdxGlobal3][interactionVolume.getIndexOnElement(2, 1)])
1647 {
1648 this->f_[eIdxGlobal3] -= pcFluxReal[1];
1649 }
1650 if (innerBoundaryVolumeFaces_[eIdxGlobal4][interactionVolume.getIndexOnElement(3, 0)])
1651 {
1652 this->f_[eIdxGlobal4] -= pcFluxReal[3];
1653 }
1654 if (innerBoundaryVolumeFaces_[eIdxGlobal4][interactionVolume.getIndexOnElement(3, 1)])
1655 {
1656 this->f_[eIdxGlobal4] += pcFluxReal[2];
1657 }
1658 }
1659 break;
1660 }
1661 }
1662 }
1663 }
1664
1665 // at least one face on boundary!
1666 else
1667 {
1668 for (int elemIdx = 0; elemIdx < 2 * dim; elemIdx++)
1669 {
1670 bool isOutside = false;
1671 for (int fIdx = 0; fIdx < dim; fIdx++)
1672 {
1673 int intVolFaceIdx = interactionVolume.getFaceIndexFromSubVolume(elemIdx, fIdx);
1674 if (interactionVolume.isOutsideFace(intVolFaceIdx))
1675 {
1676 isOutside = true;
1677 break;
1678 }
1679 }
1680 if (isOutside)
1681 {
1682 continue;
1683 }
1684
1685 auto element = interactionVolume.getSubVolumeElement(elemIdx);
1686
1687 // get global coordinate of cell centers
1688 const GlobalPosition& globalPos = element.geometry().center();
1689
1690 // cell volumes
1691 Scalar volume = element.geometry().volume();
1692
1693 // cell index
1694 int eIdxGlobal = problem_.variables().index(element);
1695
1696 //get the cell Data
1697 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
1698
1699 //permeability vector at boundary
1700 DimMatrix permeability(problem_.spatialParams().intrinsicPermeability(element));
1701
1702 // evaluate right hand side
1703 PrimaryVariables source(0);
1704 problem_.source(source, element);
1705 this->f_[eIdxGlobal] += volume / (4.0)
1706 * (source[wPhaseIdx] / density_[wPhaseIdx] + source[nPhaseIdx] / density_[nPhaseIdx]);
1707
1708 this->f_[eIdxGlobal] += evaluateErrorTerm_(cellData) * volume / (4.0);
1709
1710 //get mobilities of the phases
1711 Dune::FieldVector<Scalar, numPhases> lambda(cellData.mobility(wPhaseIdx));
1712 lambda[nPhaseIdx] = cellData.mobility(nPhaseIdx);
1713
1714 Scalar pc = cellData.capillaryPressure();
1715
1716 Scalar gravityDiff = (problem_.bBoxMax() - globalPos) * gravity_
1717 * (density_[nPhaseIdx] - density_[wPhaseIdx]);
1718
1719 pc += gravityDiff; //minus because of gravity definition!
1720
1721 for (int fIdx = 0; fIdx < dim; fIdx++)
1722 {
1723 int intVolFaceIdx = interactionVolume.getFaceIndexFromSubVolume(elemIdx, fIdx);
1724
1725 if (interactionVolume.isBoundaryFace(intVolFaceIdx))
1726 {
1727
1728 if (interactionVolume.getBoundaryType(intVolFaceIdx).isDirichlet(pressEqIdx))
1729 {
1730 int boundaryFaceIdx = interactionVolume.getIndexOnElement(elemIdx, fIdx);
1731
1732 const auto refElement = referenceElement(element);
1733
1734 const LocalPosition& localPos = refElement.position(boundaryFaceIdx, 1);
1735
1736 const GlobalPosition& globalPosFace = element.geometry().global(localPos);
1737
1738 DimVector distVec(globalPosFace - globalPos);
1739 Scalar dist = distVec.two_norm();
1740 DimVector unitDistVec(distVec);
1741 unitDistVec /= dist;
1742
1743 Scalar faceArea = interactionVolume.getFaceArea(elemIdx, fIdx);
1744
1745 // get pc and lambda at the boundary
1746 Scalar satWBound = cellData.saturation(wPhaseIdx);
1747 //check boundary sat at face 1
1748 if (interactionVolume.getBoundaryType(intVolFaceIdx).isDirichlet(satEqIdx))
1749 {
1750 Scalar satBound = interactionVolume.getDirichletValues(intVolFaceIdx)[saturationIdx];
1751 switch (saturationType_)
1752 {
1753 case sw:
1754 {
1755 satWBound = satBound;
1756 break;
1757 }
1758 case sn:
1759 {
1760 satWBound = 1 - satBound;
1761 break;
1762 }
1763 }
1764
1765 }
1766
1767 // old material law interface is deprecated: Replace this by
1768 // const auto& fluidMatrixInteraction = spatialParams.fluidMatrixInteractionAtPos(element.geometry().center());
1769 // after the release of 3.3, when the deprecated interface is no longer supported
1770 const auto fluidMatrixInteraction = Deprecated::makePcKrSw(Scalar{}, problem_.spatialParams(), element);
1771
1772 Scalar pcBound = fluidMatrixInteraction.pc(satWBound);
1773
1774 Scalar gravityDiffBound = (problem_.bBoxMax() - globalPosFace) * gravity_
1775 * (density_[nPhaseIdx] - density_[wPhaseIdx]);
1776
1777 pcBound += gravityDiffBound;
1778
1779 Dune::FieldVector<Scalar, numPhases> lambdaBound(fluidMatrixInteraction.krw(satWBound));
1780 lambdaBound[nPhaseIdx] = fluidMatrixInteraction.krn(satWBound);
1781 lambdaBound[wPhaseIdx] /= viscosity_[wPhaseIdx];
1782 lambdaBound[nPhaseIdx] /= viscosity_[nPhaseIdx];
1783
1784 Scalar potentialBound = interactionVolume.getDirichletValues(intVolFaceIdx)[pressureIdx];
1785 Scalar gdeltaZ = (problem_.bBoxMax()-globalPosFace) * gravity_;
1786
1787 //calculate potential gradients
1788 Scalar potentialDiffW = 0;
1789 Scalar potentialDiffNw = 0;
1790 switch (pressureType_)
1791 {
1792 case pw:
1793 {
1794 potentialBound += density_[wPhaseIdx]*gdeltaZ;
1795 potentialDiffW = (cellData.potential(wPhaseIdx) - potentialBound) / dist;
1796 potentialDiffNw = (cellData.potential(nPhaseIdx) - potentialBound - pcBound) / dist;
1797 break;
1798 }
1799 case pn:
1800 {
1801 potentialBound += density_[nPhaseIdx]*gdeltaZ;
1802 potentialDiffW = (cellData.potential(wPhaseIdx) - potentialBound + pcBound) / dist;
1803 potentialDiffNw = (cellData.potential(nPhaseIdx) - potentialBound) / dist;
1804 break;
1805 }
1806 }
1807
1808 Scalar lambdaTotal = (potentialDiffW >= 0.) ? lambda[wPhaseIdx] : lambdaBound[wPhaseIdx];
1809 lambdaTotal += (potentialDiffNw >= 0.) ? lambda[nPhaseIdx] : lambdaBound[nPhaseIdx];
1810
1811 DimVector permTimesNormal(0);
1812 permeability.mv(unitDistVec, permTimesNormal);
1813
1814 //calculate current matrix entry
1815 Scalar entry = lambdaTotal * (unitDistVec * permTimesNormal) / dist * faceArea;
1816
1817 //calculate right hand side
1818 Scalar pcFlux = 0;
1819
1820 switch (pressureType_)
1821 {
1822 case pw:
1823 {
1824 // calculate capillary pressure gradient
1825 DimVector pcGradient = unitDistVec;
1826 pcGradient *= (pc - pcBound) / dist;
1827
1828 //add capillary pressure term to right hand side
1829 pcFlux = 0.5 * (lambda[nPhaseIdx] + lambdaBound[nPhaseIdx])
1830 * (permTimesNormal * pcGradient) * faceArea;
1831
1832 break;
1833 }
1834 case pn:
1835 {
1836 // calculate capillary pressure gradient
1837 DimVector pcGradient = unitDistVec;
1838 pcGradient *= (pc - pcBound) / dist;
1839
1840 //add capillary pressure term to right hand side
1841 pcFlux = 0.5 * (lambda[wPhaseIdx] + lambdaBound[wPhaseIdx])
1842 * (permTimesNormal * pcGradient) * faceArea;
1843
1844 break;
1845
1846 }
1847 }
1848
1849 // set diagonal entry and right hand side entry
1850 this->A_[eIdxGlobal][eIdxGlobal] += entry;
1851 this->f_[eIdxGlobal] += entry * potentialBound;
1852
1853 if (pc == 0 && pcBound == 0)
1854 {
1855 continue;
1856 }
1857
1858 for (int i = 0; i < numPhases; i++)
1859 {
1860 switch (pressureType_)
1861 {
1862 case pw:
1863 {
1864 if (i == nPhaseIdx)
1865 {
1866 //add capillary pressure term to right hand side
1867 this->f_[eIdxGlobal] -= pcFlux;
1868 }
1869 break;
1870 }
1871 case pn:
1872 {
1873 if (i == wPhaseIdx)
1874 {
1875 //add capillary pressure term to right hand side
1876 this->f_[eIdxGlobal] += pcFlux;
1877 }
1878 break;
1879 }
1880 }
1881 }
1882
1883 }
1884 else if (interactionVolume.getBoundaryType(intVolFaceIdx).isNeumann(pressEqIdx))
1885 {
1886 Scalar J = interactionVolume.getNeumannValues(intVolFaceIdx)[wPhaseIdx]
1887 / density_[wPhaseIdx];
1888 J += interactionVolume.getNeumannValues(intVolFaceIdx)[nPhaseIdx] / density_[nPhaseIdx];
1889 this->f_[eIdxGlobal] -= J;
1890 }
1891 else
1892 {
1893 std::cout << "interactionVolume.getBoundaryType(intVolFaceIdx).isNeumann(pressEqIdx)"
1894 << interactionVolume.getBoundaryType(intVolFaceIdx).isNeumann(pressEqIdx) << "\n";
1895 DUNE_THROW(Dune::NotImplemented,
1896 "No valid boundary condition type defined for pressure equation!");
1897 }
1898 }
1899 }
1900 }
1901
1902 } // end boundaries
1903
1904 } // end vertex iterator
1905
1906 // only do more if we have more than one process
1907 if (problem_.gridView().comm().size() > 1)
1908 {
1909 // set ghost and overlap element entries
1910 for (const auto& element : elements(problem_.gridView()))
1911 {
1912 if (element.partitionType() == Dune::InteriorEntity)
1913 continue;
1914
1915 // get the global index of the cell
1916 int eIdxGlobalI = problem_.variables().index(element);
1917
1918 this->A_[eIdxGlobalI] = 0.0;
1919 this->A_[eIdxGlobalI][eIdxGlobalI] = 1.0;
1920 this->f_[eIdxGlobalI] = this->pressure()[eIdxGlobalI];
1921 }
1922 }
1923
1924 return;
1925}
1926
1927
1933template<class TypeTag>
1935{
1936 // iterate through leaf grid an evaluate c0 at cell center
1937 for (const auto& element : elements(problem_.gridView()))
1938 {
1939 int eIdxGlobal = problem_.variables().index(element);
1940
1941 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
1942
1943 const Scalar satW = cellData.saturation(wPhaseIdx);
1944
1945 // old material law interface is deprecated: Replace this by
1946 // const auto& fluidMatrixInteraction = spatialParams.fluidMatrixInteractionAtPos(element.geometry().center());
1947 // after the release of 3.3, when the deprecated interface is no longer supported
1948 const auto fluidMatrixInteraction = Deprecated::makePcKrSw(Scalar{}, problem_.spatialParams(), element);
1949
1950 const Scalar pc = fluidMatrixInteraction.pc(satW);
1951
1952 cellData.setCapillaryPressure(pc);
1953
1954 // initialize mobilities
1955 const Scalar mobilityW = fluidMatrixInteraction.krw(satW) / viscosity_[wPhaseIdx];
1956 const Scalar mobilityNw = fluidMatrixInteraction.krn(satW) / viscosity_[nPhaseIdx];
1957
1958 // initialize mobilities
1959 cellData.setMobility(wPhaseIdx, mobilityW);
1960 cellData.setMobility(nPhaseIdx, mobilityNw);
1961
1962 //initialize fractional flow functions
1963 cellData.setFracFlowFunc(wPhaseIdx, mobilityW / (mobilityW + mobilityNw));
1964 cellData.setFracFlowFunc(nPhaseIdx, mobilityNw / (mobilityW + mobilityNw));
1965 }
1966 return;
1967}
1968
1969} // end namespace Dumux
1970#endif
Helpers for deprecation.
Class including the information of an interaction volume of a MPFA O-method that does not change with...
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 permeability() noexcept
I/O name of permeability.
Definition: name.hh:143
std::string pressure(int phaseIdx) noexcept
I/O name of pressure for multiphase systems.
Definition: name.hh:34
std::string density(int phaseIdx) noexcept
I/O name of density for multiphase systems.
Definition: name.hh:65
Finite volume MPFA O-method discretization of a two-phase flow pressure equation of the sequential IM...
Definition: omethod/2dpressure.hh:69
void updateMaterialLaws()
Constitutive functions are initialized and stored in the variables object.
Definition: omethod/2dpressure.hh:1934
GlobalInteractionVolumeVector interactionVolumes_
Definition: omethod/2dpressure.hh:482
void storePressureSolution(const Element &element)
Stores the pressure solution of a cell.
Definition: omethod/2dpressure.hh:267
InnerBoundaryVolumeFaces innerBoundaryVolumeFaces_
Definition: omethod/2dpressure.hh:483
void updateInteractionVolumeInfo()
Updates interaction volumes.
Definition: omethod/2dpressure.hh:160
void addOutputVtkFields(MultiWriter &writer)
Adds pressure output to the output file.
Definition: omethod/2dpressure.hh:326
void update()
Pressure update.
Definition: omethod/2dpressure.hh:212
void initialize()
Initializes the pressure model.
Definition: omethod/2dpressure.hh:176
void storePressureSolution()
Globally stores the pressure solution.
Definition: omethod/2dpressure.hh:253
FvMpfaO2dPressure2p(Problem &problem)
Constructs a FvMpfaO2dPressure2p object.
Definition: omethod/2dpressure.hh:396
Class including the information of an interaction volume of a MPFA O-method that does not change with...
Definition: ointeractionvolume.hh:38
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
void solve()
Solves the global system of equations to get the spatial distribution of the pressure.
Definition: sequential/cellcentered/pressure.hh:527
Specifies the properties for immiscible 2p diffusion/pressure models.
Properties for a MPFA method.
Finite Volume Diffusion Model.