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