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