3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
2p/sequential/diffusion/cellcentered/pressure.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*****************************************************************************
4 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_FVPRESSURE2P_HH
25#define DUMUX_FVPRESSURE2P_HH
26
27#include <dune/common/float_cmp.hh>
28
29// dumux environment
32
33namespace Dumux {
112template<class TypeTag> class FVPressure2P: public FVPressure<TypeTag>
113{
115
116 //the model implementation
118
122
124
126
129
132 using PrimaryVariables = typename SolutionTypes::PrimaryVariables;
135
136 using ScalarSolutionType = typename SolutionTypes::ScalarSolution;
137
138 enum
139 {
140 dim = GridView::dimension, dimWorld = GridView::dimensionworld
141 };
142 enum
143 {
144 pw = Indices::pressureW,
145 pn = Indices::pressureNw,
146 pGlobal = Indices::pressureGlobal,
147 sw = Indices::saturationW,
148 sn = Indices::saturationNw,
149 pressureIdx = Indices::pressureIdx,
150 saturationIdx = Indices::saturationIdx,
151 eqIdxPress = Indices::pressureEqIdx,
152 eqIdxSat = Indices::satEqIdx
153 };
154 enum
155 {
156 wPhaseIdx = Indices::wPhaseIdx, nPhaseIdx = Indices::nPhaseIdx, numPhases = getPropValue<TypeTag, Properties::NumPhases>()
157 };
158
159 using Element = typename GridView::Traits::template Codim<0>::Entity;
160 using Intersection = typename GridView::Intersection;
161
162 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
163 using GravityVector = Dune::FieldVector<Scalar, dimWorld>;
164 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
165
166protected:
168 using EntryType = typename ParentType::EntryType;
169 enum
170 {
172 };
174
175public:
177 void getSource(EntryType& entry, const Element& element, const CellData& cellData, const bool first);
178
180 void getStorage(EntryType& entry, const Element& element, const CellData& cellData, const bool first);
181
183 void getFlux(EntryType& entry, const Intersection& intersection, const CellData& cellData, const bool first);
184
186 void getFluxOnBoundary(EntryType& entry,
187 const Intersection& intersection, const CellData& cellData, const bool first);
188
190 void updateMaterialLaws();
191
199 void initialize(bool solveTwice = true)
200 {
202
203 if (!compressibility_)
204 {
205 const auto element = *problem_.gridView().template begin<0>();
206 FluidState fluidState;
207 fluidState.setPressure(wPhaseIdx, problem_.referencePressure(element));
208 fluidState.setPressure(nPhaseIdx, problem_.referencePressure(element));
209 fluidState.setTemperature(problem_.temperature(element));
210 fluidState.setSaturation(wPhaseIdx, 1.);
211 fluidState.setSaturation(nPhaseIdx, 0.);
212 density_[wPhaseIdx] = FluidSystem::density(fluidState, wPhaseIdx);
213 density_[nPhaseIdx] = FluidSystem::density(fluidState, nPhaseIdx);
214 viscosity_[wPhaseIdx] = FluidSystem::viscosity(fluidState, wPhaseIdx);
215 viscosity_[nPhaseIdx] = FluidSystem::viscosity(fluidState, nPhaseIdx);
216 }
217
220
221 this->assemble(true);
222 this->solve();
223 if (solveTwice)
224 {
225 PressureSolutionVector pressureOld(this->pressure());
226
227 this->assemble(false);
228 this->solve();
229
230 PressureSolutionVector pressureDiff(pressureOld);
231 pressureDiff -= this->pressure();
232 pressureOld = this->pressure();
233 Scalar pressureNorm = pressureDiff.infinity_norm();
234 pressureNorm /= pressureOld.infinity_norm();
235 int numIter = 0;
236 while (pressureNorm > 1e-5 && numIter < 10)
237 {
239 this->assemble(false);
240 this->solve();
241
242 pressureDiff = pressureOld;
243 pressureDiff -= this->pressure();
244 pressureNorm = pressureDiff.infinity_norm();
245 pressureOld = this->pressure();
246 pressureNorm /= pressureOld.infinity_norm();
247
248 numIter++;
249 }
250 // std::cout<<"Pressure defect = "<<pressureNorm<<"; "<<
251 // numIter<<" Iterations needed for initial pressure field"<<std::endl;
252 }
253
255 }
256
262 void update()
263 {
264 timeStep_ = problem_.timeManager().timeStepSize();
265 //error bounds for error term for incompressible models
266 //to correct unphysical saturation over/undershoots due to saturation transport
267 if (!compressibility_)
268 {
269 maxError_ = 0.0;
270 int size = problem_.gridView().size(0);
271 for (int i = 0; i < size; i++)
272 {
273 Scalar sat = 0;
274 switch (saturationType_)
275 {
276 case sw:
277 sat = problem_.variables().cellData(i).saturation(wPhaseIdx);
278 break;
279 case sn:
280 sat = problem_.variables().cellData(i).saturation(nPhaseIdx);
281 break;
282 }
283 if (sat > 1.0)
284 {
285 using std::max;
286 maxError_ = max(maxError_, (sat - 1.0) / timeStep_);
287 }
288 if (sat < 0.0)
289 {
290 using std::max;
291 maxError_ = max(maxError_, (-sat) / timeStep_);
292 }
293 }
294 }
295
297
299 }
300
307 {
309
310 //reset velocities
311 int size = problem_.gridView().size(0);
312 for (int i = 0; i < size; i++)
313 {
314 CellData& cellData = problem_.variables().cellData(i);
315 cellData.fluxData().resetVelocity();
316 }
317 }
318
321 {
322 // iterate through leaf grid
323 for (const auto& element : elements(problem_.gridView()))
324 {
325 asImp_().storePressureSolution(element);
326 }
327 }
328
336 void storePressureSolution(const Element& element)
337 {
338 int eIdxGlobal = problem_.variables().index(element);
339 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
340
341 if (compressibility_)
342 {
343 density_[wPhaseIdx] = cellData.density(wPhaseIdx);
344 density_[nPhaseIdx] = cellData.density(nPhaseIdx);
345 }
346
347 switch (pressureType_)
348 {
349 case pw:
350 {
351 Scalar pressW = this->pressure()[eIdxGlobal];
352 Scalar pc = cellData.capillaryPressure();
353
354 cellData.setPressure(wPhaseIdx, pressW);
355 cellData.setPressure(nPhaseIdx, pressW + pc);
356
357 Scalar gravityDiff = (problem_.bBoxMax() - element.geometry().center()) * gravity_;
358 Scalar potW = pressW + gravityDiff * density_[wPhaseIdx];
359 Scalar potNw = pressW + pc + gravityDiff * density_[nPhaseIdx];
360
361 cellData.setPotential(wPhaseIdx, potW);
362 cellData.setPotential(nPhaseIdx, potNw);
363
364 break;
365 }
366 case pn:
367 {
368 Scalar pressNw = this->pressure()[eIdxGlobal];
369 Scalar pc = cellData.capillaryPressure();
370
371 cellData.setPressure(nPhaseIdx, pressNw);
372 cellData.setPressure(wPhaseIdx, pressNw - pc);
373
374 Scalar gravityDiff = (problem_.bBoxMax() - element.geometry().center()) * gravity_;
375 Scalar potW = pressNw - pc + gravityDiff * density_[wPhaseIdx];
376 Scalar potNw = pressNw + gravityDiff * density_[nPhaseIdx];
377
378 cellData.setPotential(wPhaseIdx, potW);
379 cellData.setPotential(nPhaseIdx, potNw);
380
381 break;
382 }
383 case pGlobal:
384 {
385 Scalar press = this->pressure()[eIdxGlobal];
386 cellData.setGlobalPressure(press);
387
388 Scalar pc = cellData.capillaryPressure();
389 Scalar gravityDiff = (problem_.bBoxMax() - element.geometry().center()) * gravity_;
390
391 //This is only an estimation!!! -> only used for uwind directions and time-step estimation!!!
392 Scalar potW = press - cellData.fracFlowFunc(nPhaseIdx) * pc + gravityDiff * density_[wPhaseIdx];
393 Scalar potNw = press - cellData.fracFlowFunc(wPhaseIdx) * pc + gravityDiff * density_[nPhaseIdx];
394
395 cellData.setPotential(wPhaseIdx, potW);
396 cellData.setPotential(nPhaseIdx, potNw);
397
398 break;
399 }
400 }
401 cellData.fluxData().resetVelocity();
402 }
403
413 template<class MultiWriter>
414 void addOutputVtkFields(MultiWriter &writer)
415 {
416 int size = problem_.gridView().size(0);
417 ScalarSolutionType *pressure = writer.allocateManagedBuffer(size);
418 ScalarSolutionType *pressureSecond = 0;
419 ScalarSolutionType *pc = 0;
420 ScalarSolutionType *potentialW = 0;
421 ScalarSolutionType *potentialNw = 0;
422 ScalarSolutionType *mobilityW = 0;
423 ScalarSolutionType *mobilityNW = 0;
424
425 if (vtkOutputLevel_ > 0)
426 {
427 pressureSecond = writer.allocateManagedBuffer(size);
428 pc = writer.allocateManagedBuffer(size);
429 mobilityW = writer.allocateManagedBuffer(size);
430 mobilityNW = writer.allocateManagedBuffer(size);
431 }
432 if (vtkOutputLevel_ > 1)
433 {
434 potentialW = writer.allocateManagedBuffer(size);
435 potentialNw = writer.allocateManagedBuffer(size);
436 }
437
438 for (int i = 0; i < size; i++)
439 {
440 CellData& cellData = problem_.variables().cellData(i);
441
442 if (pressureType_ == pw)
443 {
444 (*pressure)[i] = cellData.pressure(wPhaseIdx);
445 if (vtkOutputLevel_ > 0)
446 {
447 (*pressureSecond)[i] = cellData.pressure(nPhaseIdx);
448 }
449 }
450 else if (pressureType_ == pn)
451 {
452 (*pressure)[i] = cellData.pressure(nPhaseIdx);
453 if (vtkOutputLevel_ > 0)
454 {
455 (*pressureSecond)[i] = cellData.pressure(wPhaseIdx);
456 }
457 }
458 else if (pressureType_ == pGlobal)
459 {
460 (*pressure)[i] = cellData.globalPressure();
461 }
462 if (vtkOutputLevel_ > 0)
463 {
464 (*pc)[i] = cellData.capillaryPressure();
465 (*mobilityW)[i] = cellData.mobility(wPhaseIdx);
466 (*mobilityNW)[i] = cellData.mobility(nPhaseIdx);
467 if (compressibility_)
468 {
469 (*mobilityW)[i] = (*mobilityW)[i]/cellData.density(wPhaseIdx);
470 (*mobilityNW)[i] = (*mobilityNW)[i]/cellData.density(nPhaseIdx);
471 }
472 }
473 if (vtkOutputLevel_ > 1)
474 {
475 (*potentialW)[i] = cellData.potential(wPhaseIdx);
476 (*potentialNw)[i] = cellData.potential(nPhaseIdx);
477 }
478 }
479
480 if (pressureType_ == pw)
481 {
482 writer.attachCellData(*pressure, "wetting pressure");
483 if (vtkOutputLevel_ > 0)
484 {
485 writer.attachCellData(*pressureSecond, "nonwetting pressure");
486 }
487 }
488 else if (pressureType_ == pn)
489 {
490 writer.attachCellData(*pressure, "nonwetting pressure");
491 if (vtkOutputLevel_ > 0)
492 {
493 writer.attachCellData(*pressureSecond, "wetting pressure");
494 }
495 }
496 if (pressureType_ == pGlobal)
497 {
498
499 writer.attachCellData(*pressure, "global pressure");
500 }
501
502 if (vtkOutputLevel_ > 0)
503 {
504 writer.attachCellData(*pc, "capillary pressure");
505 writer.attachCellData(*mobilityW, "wetting mobility");
506 writer.attachCellData(*mobilityNW, "nonwetting mobility");
507 }
508 if (vtkOutputLevel_ > 1)
509 {
510 writer.attachCellData(*potentialW, "wetting potential");
511 writer.attachCellData(*potentialNw, "nonwetting potential");
512 }
513
514 if (compressibility_)
515 {
516 if (vtkOutputLevel_ > 0)
517 {
518 ScalarSolutionType *densityWetting = writer.allocateManagedBuffer(size);
519 ScalarSolutionType *densityNonwetting = writer.allocateManagedBuffer(size);
520 ScalarSolutionType *viscosityWetting = writer.allocateManagedBuffer(size);
521 ScalarSolutionType *viscosityNonwetting = writer.allocateManagedBuffer(size);
522
523 for (int i = 0; i < size; i++)
524 {
525 CellData& cellData = problem_.variables().cellData(i);
526 (*densityWetting)[i] = cellData.density(wPhaseIdx);
527 (*densityNonwetting)[i] = cellData.density(nPhaseIdx);
528 (*viscosityWetting)[i] = cellData.viscosity(wPhaseIdx);
529 (*viscosityNonwetting)[i] = cellData.viscosity(nPhaseIdx);
530 }
531
532 writer.attachCellData(*densityWetting, "wetting density");
533 writer.attachCellData(*densityNonwetting, "nonwetting density");
534 writer.attachCellData(*viscosityWetting, "wetting viscosity");
535 writer.attachCellData(*viscosityNonwetting, "nonwetting viscosity");
536 }
537 }
538 }
539
545 FVPressure2P(Problem& problem) :
546 ParentType(problem), problem_(problem), gravity_(problem.gravity()), maxError_(0.), timeStep_(1.)
547 {
548 if (pressureType_ != pw && pressureType_ != pn && pressureType_ != pGlobal)
549 {
550 DUNE_THROW(Dune::NotImplemented, "Pressure type not supported!");
551 }
552 if (pressureType_ == pGlobal && compressibility_)
553 {
554 DUNE_THROW(Dune::NotImplemented, "Compressibility not supported for global pressure!");
555 }
556 if (saturationType_ != sw && saturationType_ != sn)
557 {
558 DUNE_THROW(Dune::NotImplemented, "Saturation type not supported!");
559 }
560
561 ErrorTermFactor_ = getParam<Scalar>("Impet.ErrorTermFactor");
562 ErrorTermLowerBound_ = getParam<Scalar>("Impet.ErrorTermLowerBound");
563 ErrorTermUpperBound_ = getParam<Scalar>("Impet.ErrorTermUpperBound");
564
565 density_[wPhaseIdx] = 0.;
566 density_[nPhaseIdx] = 0.;
567 viscosity_[wPhaseIdx] = 0.;
568 viscosity_[nPhaseIdx] = 0.;
569
570 vtkOutputLevel_ = getParam<int>("Vtk.OutputLevel");
571 }
572
573private:
575 Implementation &asImp_()
576 { return *static_cast<Implementation *>(this); }
577
579 const Implementation &asImp_() const
580 { return *static_cast<const Implementation *>(this); }
581
582 Problem& problem_;
583 const GravityVector& gravity_;
584
585 Scalar maxError_;
586 Scalar timeStep_;
587 Scalar ErrorTermFactor_;
588 Scalar ErrorTermLowerBound_;
589 Scalar ErrorTermUpperBound_;
590
591 Scalar density_[numPhases];
592 Scalar viscosity_[numPhases];
593
594 int vtkOutputLevel_;
595
596 static const bool compressibility_ = getPropValue<TypeTag, Properties::EnableCompressibility>();
598 static const int pressureType_ = getPropValue<TypeTag, Properties::PressureFormulation>();
600 static const int saturationType_ = getPropValue<TypeTag, Properties::SaturationFormulation>();
601};
602
610template<class TypeTag>
611void FVPressure2P<TypeTag>::getSource(EntryType& entry, const Element& element
612 , const CellData& cellData, const bool first)
613{
614 // cell volume, assume linear map here
615 Scalar volume = element.geometry().volume();
616
617 // get sources from problem
618 PrimaryVariables sourcePhase(0.0);
619 problem_.source(sourcePhase, element);
620
621 if (!compressibility_)
622 {
623 sourcePhase[wPhaseIdx] /= density_[wPhaseIdx];
624 sourcePhase[nPhaseIdx] /= density_[nPhaseIdx];
625 }
626
627 entry[rhs] = volume * (sourcePhase[wPhaseIdx] + sourcePhase[nPhaseIdx]);
628}
629
654template<class TypeTag>
655void FVPressure2P<TypeTag>::getStorage(EntryType& entry, const Element& element
656 , const CellData& cellData, const bool first)
657{
658 //volume correction due to density differences
659 if (compressibility_ && !first)
660 {
661 // cell volume, assume linear map here
662 Scalar volume = element.geometry().volume();
663
664 Scalar porosity = problem_.spatialParams().porosity(element);
665
666 switch (saturationType_)
667 {
668 case sw:
669 {
670 entry[rhs] = -(cellData.volumeCorrection()/timeStep_ * porosity * volume
671 * (cellData.density(wPhaseIdx) - cellData.density(nPhaseIdx)));
672 break;
673 }
674 case sn:
675 {
676 entry[rhs] = -(cellData.volumeCorrection()/timeStep_ * porosity * volume
677 * (cellData.density(nPhaseIdx) - cellData.density(wPhaseIdx)));
678 break;
679 }
680 }
681 }
682 else if (!compressibility_ && !first)
683 {
684 // error term for incompressible models to correct non-physical saturation over/undershoots due to saturation transport
685 // error reduction routine: volumetric error is damped and inserted to right hand side
686 Scalar sat = 0;
687 switch (saturationType_)
688 {
689 case sw:
690 sat = cellData.saturation(wPhaseIdx);
691 break;
692 case sn:
693 sat = cellData.saturation(nPhaseIdx);
694 break;
695 }
696
697 Scalar volume = element.geometry().volume();
698
699 Scalar error = (sat > 1.0) ? sat - 1.0 : 0.0;
700 if (sat < 0.0) {error = sat;}
701 error /= timeStep_;
702
703 using std::abs;
704 Scalar errorAbs = abs(error);
705
706 if ((errorAbs*timeStep_ > 1e-6) && (errorAbs > ErrorTermLowerBound_ * maxError_) && (!problem_.timeManager().willBeFinished()))
707 {
708 entry[rhs] = ErrorTermFactor_ * error * volume;
709 }
710 }
711}
712
718template<class TypeTag>
719void FVPressure2P<TypeTag>::getFlux(EntryType& entry, const Intersection& intersection
720 , const CellData& cellData, const bool first)
721{
722 auto elementI = intersection.inside();
723 auto elementJ = intersection.outside();
724
725 const CellData& cellDataJ = problem_.variables().cellData(problem_.variables().index(elementJ));
726
727 // get global coordinates of cell centers
728 const GlobalPosition& globalPosI = elementI.geometry().center();
729 const GlobalPosition& globalPosJ = elementJ.geometry().center();
730
731 // get mobilities and fractional flow factors
732 Scalar lambdaWI = cellData.mobility(wPhaseIdx);
733 Scalar lambdaNwI = cellData.mobility(nPhaseIdx);
734 Scalar lambdaWJ = cellDataJ.mobility(wPhaseIdx);
735 Scalar lambdaNwJ = cellDataJ.mobility(nPhaseIdx);
736
737 // get capillary pressure
738 Scalar pcI = cellData.capillaryPressure();
739 Scalar pcJ = cellDataJ.capillaryPressure();
740
741 // get face normal
742 const Dune::FieldVector<Scalar, dim>& unitOuterNormal = intersection.centerUnitOuterNormal();
743
744 // get face area
745 Scalar faceArea = intersection.geometry().volume();
746
747 // distance vector between barycenters
748 GlobalPosition distVec = globalPosJ - globalPosI;
749
750 // compute distance between cell centers
751 Scalar dist = distVec.two_norm();
752
753 // compute vectorized permeabilities
754 DimMatrix meanPermeability(0);
755
756 problem_.spatialParams().meanK(meanPermeability, problem_.spatialParams().intrinsicPermeability(elementI),
757 problem_.spatialParams().intrinsicPermeability(elementJ));
758
759 Dune::FieldVector<Scalar, dim> permeability(0);
760 meanPermeability.mv(unitOuterNormal, permeability);
761
762 Scalar rhoMeanW = 0;
763 Scalar rhoMeanNw = 0;
764 if (compressibility_)
765 {
766 rhoMeanW = 0.5 * (cellData.density(wPhaseIdx) + cellDataJ.density(wPhaseIdx));
767 rhoMeanNw = 0.5 * (cellData.density(nPhaseIdx) + cellDataJ.density(nPhaseIdx));
768 }
769
770 // calculate potential gradients
771 Scalar potentialDiffW = 0;
772 Scalar potentialDiffNw = 0;
773
774 // if we are at the very first iteration we can't calculate phase potentials
775 if (!first)
776 {
777 potentialDiffW = cellData.potential(wPhaseIdx) - cellDataJ.potential(wPhaseIdx);
778 potentialDiffNw = cellData.potential(nPhaseIdx) - cellDataJ.potential(nPhaseIdx);
779
780 if (compressibility_)
781 {
782 density_[wPhaseIdx] = (potentialDiffW > 0.) ? cellData.density(wPhaseIdx) : cellDataJ.density(wPhaseIdx);
783 density_[nPhaseIdx] = (potentialDiffNw > 0.) ? cellData.density(nPhaseIdx) : cellDataJ.density(nPhaseIdx);
784
785 density_[wPhaseIdx] = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffW, 0.0, 1.0e-30)) ? rhoMeanW : density_[wPhaseIdx];
786 density_[nPhaseIdx] = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffNw, 0.0, 1.0e-30)) ? rhoMeanNw : density_[nPhaseIdx];
787
788 potentialDiffW = cellData.pressure(wPhaseIdx) - cellDataJ.pressure(wPhaseIdx);
789 potentialDiffNw = cellData.pressure(nPhaseIdx) - cellDataJ.pressure(nPhaseIdx);
790
791 potentialDiffW += density_[wPhaseIdx] * (distVec * gravity_);
792 potentialDiffNw += density_[nPhaseIdx] * (distVec * gravity_);
793 }
794 }
795
796 // do the upwinding of the mobility depending on the phase potentials
797 Scalar lambdaW = (potentialDiffW > 0.) ? lambdaWI : lambdaWJ;
798 lambdaW = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffW, 0.0, 1.0e-30)) ? 0.5 * (lambdaWI + lambdaWJ) : lambdaW;
799 Scalar lambdaNw = (potentialDiffNw > 0) ? lambdaNwI : lambdaNwJ;
800 lambdaNw = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffNw, 0.0, 1.0e-30)) ? 0.5 * (lambdaNwI + lambdaNwJ) : lambdaNw;
801
802 if (compressibility_)
803 {
804 density_[wPhaseIdx] = (potentialDiffW > 0.) ? cellData.density(wPhaseIdx) : cellDataJ.density(wPhaseIdx);
805 density_[nPhaseIdx] = (potentialDiffNw > 0.) ? cellData.density(nPhaseIdx) : cellDataJ.density(nPhaseIdx);
806
807 density_[wPhaseIdx] = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffW, 0.0, 1.0e-30)) ? rhoMeanW : density_[wPhaseIdx];
808 density_[nPhaseIdx] = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffNw, 0.0, 1.0e-30)) ? rhoMeanNw : density_[nPhaseIdx];
809 }
810
811 Scalar scalarPerm = permeability.two_norm();
812 // calculate current matrix entry
813 entry[matrix] = (lambdaW + lambdaNw) * scalarPerm / dist * faceArea;
814
815 // calculate right hand side
816 // calculate unit distVec
817 distVec /= dist;
818 Scalar areaScaling = (unitOuterNormal * distVec);
819 // this treatment of g allows to account for gravity flux through faces where the face normal has no z component (e.g. parallelepiped grids)
820 entry[rhs] = (lambdaW * density_[wPhaseIdx] + lambdaNw * density_[nPhaseIdx]) * scalarPerm * (gravity_ * distVec) * faceArea * areaScaling;
821
822 if (pressureType_ == pw)
823 {
824 // add capillary pressure term to right hand side
825 entry[rhs] += 0.5 * (lambdaNwI + lambdaNwJ) * scalarPerm * (pcI - pcJ) / dist * faceArea;
826 }
827 else if (pressureType_ == pn)
828 {
829 // add capillary pressure term to right hand side
830 entry[rhs] -= 0.5 * (lambdaWI + lambdaWJ) * scalarPerm * (pcI - pcJ) / dist * faceArea;
831 }
832}
833
842template<class TypeTag>
844const Intersection& intersection, const CellData& cellData, const bool first)
845{
846 auto element = intersection.inside();
847
848 // get global coordinates of cell centers
849 const GlobalPosition& globalPosI = element.geometry().center();
850
851 // center of face in global coordinates
852 const GlobalPosition& globalPosJ = intersection.geometry().center();
853
854 // get mobilities and fractional flow factors
855 Scalar lambdaWI = cellData.mobility(wPhaseIdx);
856 Scalar lambdaNwI = cellData.mobility(nPhaseIdx);
857 Scalar fractionalWI = cellData.fracFlowFunc(wPhaseIdx);
858 Scalar fractionalNwI = cellData.fracFlowFunc(nPhaseIdx);
859
860 // get capillary pressure
861 Scalar pcI = cellData.capillaryPressure();
862
863 // get face index
864 int isIndexI = intersection.indexInInside();
865
866 // get face normal
867 const Dune::FieldVector<Scalar, dim>& unitOuterNormal = intersection.centerUnitOuterNormal();
868
869 // get face area
870 Scalar faceArea = intersection.geometry().volume();
871
872 // distance vector between barycenters
873 GlobalPosition distVec = globalPosJ - globalPosI;
874
875 // compute distance between cell centers
876 Scalar dist = distVec.two_norm();
877
878 BoundaryTypes bcType;
879 problem_.boundaryTypes(bcType, intersection);
880 PrimaryVariables boundValues(0.0);
881
882 if (bcType.isDirichlet(eqIdxPress))
883 {
884 problem_.dirichlet(boundValues, intersection);
885
886 //permeability vector at boundary
887 // compute vectorized permeabilities
888 DimMatrix meanPermeability(0);
889
890 problem_.spatialParams().meanK(meanPermeability,
891 problem_.spatialParams().intrinsicPermeability(element));
892
893 Dune::FieldVector<Scalar, dim> permeability(0);
894 meanPermeability.mv(unitOuterNormal, permeability);
895
896 // determine saturation at the boundary -> if no saturation is known directly at the boundary use the cell saturation
897 Scalar satW = 0;
898 Scalar satNw = 0;
899 if (bcType.isDirichlet(eqIdxSat))
900 {
901 switch (saturationType_)
902 {
903 case sw:
904 {
905 satW = boundValues[saturationIdx];
906 satNw = 1 - boundValues[saturationIdx];
907 break;
908 }
909 case sn:
910 {
911 satW = 1 - boundValues[saturationIdx];
912 satNw = boundValues[saturationIdx];
913 break;
914 }
915 }
916 }
917 else
918 {
919 satW = cellData.saturation(wPhaseIdx);
920 satNw = cellData.saturation(nPhaseIdx);
921 }
922
923 const Scalar temperature = problem_.temperature(element);
924
925 // get Dirichlet pressure boundary condition
926 const Scalar pressBound = boundValues[pressureIdx];
927
928 //calculate constitutive relations depending on the kind of saturation used
929
930 const auto fluidMatrixInteraction = problem_.spatialParams().fluidMatrixInteractionAtPos(element.geometry().center());
931 const Scalar pcBound = fluidMatrixInteraction.pc(satW);
932
933 // determine phase pressures from primary pressure variable
934 Scalar pressW = 0;
935 Scalar pressNw = 0;
936 if (pressureType_ == pw)
937 {
938 pressW = pressBound;
939 pressNw = pressBound + pcBound;
940 }
941 else if (pressureType_ == pn)
942 {
943 pressW = pressBound - pcBound;
944 pressNw = pressBound;
945 }
946
947 Scalar densityWBound = density_[wPhaseIdx];
948 Scalar densityNwBound = density_[nPhaseIdx];
949 Scalar viscosityWBound = viscosity_[wPhaseIdx];
950 Scalar viscosityNwBound = viscosity_[nPhaseIdx];
951 Scalar rhoMeanW = 0;
952 Scalar rhoMeanNw = 0;
953
954 if (compressibility_)
955 {
956 FluidState fluidState;
957 fluidState.setSaturation(wPhaseIdx, satW);
958 fluidState.setSaturation(nPhaseIdx, satNw);
959 fluidState.setTemperature(temperature);
960 fluidState.setPressure(wPhaseIdx, pressW);
961 fluidState.setPressure(nPhaseIdx, pressNw);
962
963 densityWBound = FluidSystem::density(fluidState, wPhaseIdx);
964 densityNwBound = FluidSystem::density(fluidState, nPhaseIdx);
965 viscosityWBound = FluidSystem::viscosity(fluidState, wPhaseIdx) / densityWBound;
966 viscosityNwBound = FluidSystem::viscosity(fluidState, nPhaseIdx) / densityNwBound;
967
968 rhoMeanW = 0.5 * (cellData.density(wPhaseIdx) + densityWBound);
969 rhoMeanNw = 0.5 * (cellData.density(nPhaseIdx) + densityNwBound);
970 }
971
972 Scalar lambdaWBound = fluidMatrixInteraction.krw(satW)
973 / viscosityWBound;
974 Scalar lambdaNwBound = fluidMatrixInteraction.krn(satW)
975 / viscosityNwBound;
976
977 Scalar fractionalWBound = lambdaWBound / (lambdaWBound + lambdaNwBound);
978 Scalar fractionalNwBound = lambdaNwBound / (lambdaWBound + lambdaNwBound);
979
980 Scalar fMeanW = 0.5 * (fractionalWI + fractionalWBound);
981 Scalar fMeanNw = 0.5 * (fractionalNwI + fractionalNwBound);
982
983 Scalar potentialDiffW = 0;
984 Scalar potentialDiffNw = 0;
985
986 if (!first)
987 {
988 potentialDiffW = cellData.fluxData().upwindPotential(wPhaseIdx, isIndexI);
989 potentialDiffNw = cellData.fluxData().upwindPotential(nPhaseIdx, isIndexI);
990
991 if (compressibility_)
992 {
993 density_[wPhaseIdx] = (potentialDiffW > 0.) ? cellData.density(wPhaseIdx) : densityWBound;
994 density_[nPhaseIdx] = (potentialDiffNw > 0.) ? cellData.density(nPhaseIdx) : densityNwBound;
995
996 density_[wPhaseIdx] = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffW, 0.0, 1.0e-30)) ? rhoMeanW : density_[wPhaseIdx];
997 density_[nPhaseIdx] = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffNw, 0.0, 1.0e-30)) ? rhoMeanNw : density_[nPhaseIdx];
998 }
999
1000 // calculate potential gradient
1001 switch (pressureType_)
1002 {
1003 case pw:
1004 {
1005 potentialDiffW = (cellData.pressure(wPhaseIdx) - pressBound);
1006 potentialDiffNw = (cellData.pressure(nPhaseIdx) - pressBound - pcBound);
1007 break;
1008 }
1009 case pn:
1010 {
1011 potentialDiffW = (cellData.pressure(wPhaseIdx) - pressBound + pcBound);
1012 potentialDiffNw = (cellData.pressure(nPhaseIdx) - pressBound);
1013 break;
1014 }
1015 case pGlobal:
1016 {
1017 potentialDiffW = (cellData.globalPressure() - pressBound - fMeanNw * (pcI - pcBound));
1018 potentialDiffNw = (cellData.globalPressure() - pressBound + fMeanW * (pcI - pcBound));
1019 break;
1020 }
1021 }
1022
1023 potentialDiffW += density_[wPhaseIdx] * (distVec * gravity_);
1024 potentialDiffNw += density_[nPhaseIdx] * (distVec * gravity_);
1025 }
1026
1027 // do the upwinding of the mobility depending on the phase potentials
1028 Scalar lambdaW = (potentialDiffW > 0.) ? lambdaWI : lambdaWBound;
1029 lambdaW = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffW, 0.0, 1.0e-30)) ? 0.5 * (lambdaWI + lambdaWBound) : lambdaW;
1030 Scalar lambdaNw = (potentialDiffNw > 0.) ? lambdaNwI : lambdaNwBound;
1031 lambdaNw = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffNw, 0.0, 1.0e-30)) ? 0.5 * (lambdaNwI + lambdaNwBound) : lambdaNw;
1032
1033 if (compressibility_)
1034 {
1035 density_[wPhaseIdx] = (potentialDiffW > 0.) ? cellData.density(wPhaseIdx) : densityWBound;
1036 density_[wPhaseIdx] = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffW, 0.0, 1.0e-30)) ? rhoMeanW : density_[wPhaseIdx];
1037 density_[nPhaseIdx] = (potentialDiffNw > 0.) ? cellData.density(nPhaseIdx) : densityNwBound;
1038 density_[nPhaseIdx] = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffNw, 0.0, 1.0e-30)) ? rhoMeanNw : density_[nPhaseIdx];
1039 }
1040
1041 Scalar scalarPerm = permeability.two_norm();
1042 // calculate current matrix entry
1043 entry[matrix] = (lambdaW + lambdaNw) * scalarPerm / dist * faceArea;
1044 entry[rhs] = entry[matrix] * pressBound;
1045
1046 // calculate right hand side
1047 // calculate unit distVec
1048 distVec /= dist;
1049 Scalar areaScaling = (unitOuterNormal * distVec);
1050 // this treatment of g allows to account for gravity flux through faces where the face normal has no z component (e.g. parallelepiped grids)
1051 entry[rhs] -= (lambdaW * density_[wPhaseIdx] + lambdaNw * density_[nPhaseIdx]) * scalarPerm * (gravity_ * distVec)
1052 * faceArea * areaScaling;
1053
1054 if (pressureType_ == pw)
1055 {
1056 //add capillary pressure term to right hand side
1057 entry[rhs] -= 0.5 * (lambdaNwI + lambdaNwBound) * scalarPerm * (pcI - pcBound) / dist * faceArea;
1058 }
1059 else if (pressureType_ == pn)
1060 {
1061 //add capillary pressure term to right hand side
1062 entry[rhs] += 0.5 * (lambdaWI + lambdaWBound) * scalarPerm * (pcI - pcBound) / dist * faceArea;
1063 }
1064 }
1065 // set Neumann boundary condition
1066 else if (bcType.isNeumann(eqIdxPress))
1067 {
1068 problem_.neumann(boundValues, intersection);
1069
1070 if (!compressibility_)
1071 {
1072 boundValues[wPhaseIdx] /= density_[wPhaseIdx];
1073 boundValues[nPhaseIdx] /= density_[nPhaseIdx];
1074 }
1075 entry[rhs] = -(boundValues[wPhaseIdx] + boundValues[nPhaseIdx]) * faceArea;
1076 }
1077 else
1078 {
1079 DUNE_THROW(Dune::NotImplemented, "No valid boundary condition type defined for pressure equation!");
1080 }
1081}
1082
1091template<class TypeTag>
1093{
1094 for (const auto& element : elements(problem_.gridView()))
1095 {
1096 int eIdxGlobal = problem_.variables().index(element);
1097
1098 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
1099
1100 const Scalar temperature = problem_.temperature(element);
1101
1102 // determine phase saturations from primary saturation variable
1103 const Scalar satW = cellData.saturation(wPhaseIdx);
1104
1105 const auto fluidMatrixInteraction = problem_.spatialParams().fluidMatrixInteractionAtPos(element.geometry().center());
1106 const Scalar pc = fluidMatrixInteraction.pc(satW);
1107
1108 // determine phase pressures from primary pressure variable
1109 Scalar pressW = 0;
1110 Scalar pressNw = 0;
1111 if (pressureType_ == pw)
1112 {
1113 pressW = cellData.pressure(wPhaseIdx);
1114 pressNw = pressW + pc;
1115 }
1116 else if (pressureType_ == pn)
1117 {
1118 pressNw = cellData.pressure(nPhaseIdx);
1119 pressW = pressNw - pc;
1120 }
1121
1122 if (compressibility_)
1123 {
1124 FluidState& fluidState = cellData.fluidState();
1125 fluidState.setTemperature(temperature);
1126
1127 fluidState.setPressure(wPhaseIdx, pressW);
1128 fluidState.setPressure(nPhaseIdx, pressNw);
1129
1130 density_[wPhaseIdx] = FluidSystem::density(fluidState, wPhaseIdx);
1131 density_[nPhaseIdx] = FluidSystem::density(fluidState, nPhaseIdx);
1132
1133 viscosity_[wPhaseIdx]= FluidSystem::viscosity(fluidState, wPhaseIdx);
1134 viscosity_[nPhaseIdx] = FluidSystem::viscosity(fluidState, nPhaseIdx);
1135
1136 //store density
1137 fluidState.setDensity(wPhaseIdx, density_[wPhaseIdx]);
1138 fluidState.setDensity(nPhaseIdx, density_[nPhaseIdx]);
1139
1140 //store viscosity
1141 fluidState.setViscosity(wPhaseIdx, viscosity_[wPhaseIdx]);
1142 fluidState.setViscosity(nPhaseIdx, viscosity_[nPhaseIdx]);
1143 }
1144 else
1145 {
1146 cellData.setCapillaryPressure(pc);
1147
1148 if (pressureType_ != pGlobal)
1149 {
1150 cellData.setPressure(wPhaseIdx, pressW);
1151 cellData.setPressure(nPhaseIdx, pressNw);
1152 }
1153 }
1154
1155 // initialize mobilities
1156 Scalar mobilityW = fluidMatrixInteraction.krw(satW) / viscosity_[wPhaseIdx];
1157 Scalar mobilityNw = fluidMatrixInteraction.krn(satW) / viscosity_[nPhaseIdx];
1158
1159 if (compressibility_)
1160 {
1161 mobilityW *= density_[wPhaseIdx];
1162 mobilityNw *= density_[nPhaseIdx];
1163 }
1164
1165 // initialize mobilities
1166 cellData.setMobility(wPhaseIdx, mobilityW);
1167 cellData.setMobility(nPhaseIdx, mobilityNw);
1168
1169 // initialize fractional flow functions
1170 cellData.setFracFlowFunc(wPhaseIdx, mobilityW / (mobilityW + mobilityNw));
1171 cellData.setFracFlowFunc(nPhaseIdx, mobilityNw / (mobilityW + mobilityNw));
1172
1173 const Scalar gravityDiff = (problem_.bBoxMax() - element.geometry().center()) * gravity_;
1174
1175 Scalar potW = pressW + gravityDiff * density_[wPhaseIdx];
1176 Scalar potNw = pressNw + gravityDiff * density_[nPhaseIdx];
1177
1178 if (pressureType_ == pGlobal)
1179 {
1180 potW = cellData.globalPressure() - cellData.fracFlowFunc(nPhaseIdx) * pc + gravityDiff * density_[wPhaseIdx];
1181 potNw = cellData.globalPressure() - cellData.fracFlowFunc(wPhaseIdx) * pc + gravityDiff * density_[nPhaseIdx];
1182 }
1183
1184 cellData.setPotential(wPhaseIdx, potW);
1185 cellData.setPotential(nPhaseIdx, potNw);
1186 }
1187}
1188
1189} // end namespace Dumux
1190#endif
Definition: adapt.hh:29
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type GetProp
get the type of a property
Definition: propertysystem.hh:141
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:150
std::string temperature() noexcept
I/O name of temperature for equilibrium models.
Definition: name.hh:51
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 density(int phaseIdx) noexcept
I/O name of density for multiphase systems.
Definition: name.hh:65
std::string porosity() noexcept
I/O name of porosity.
Definition: name.hh:139
Scalar volume(Shape shape, Scalar inscribedRadius)
Returns the volume of a given geometry based on the inscribed radius.
Definition: poreproperties.hh:73
Finite Volume discretization of a two-phase flow pressure equation of the sequential IMPES model.
Definition: 2p/sequential/diffusion/cellcentered/pressure.hh:113
void storePressureSolution(const Element &element)
Stores the pressure solution of a cell.
Definition: 2p/sequential/diffusion/cellcentered/pressure.hh:336
void storePressureSolution()
Globally stores the pressure solution.
Definition: 2p/sequential/diffusion/cellcentered/pressure.hh:320
void getFlux(EntryType &entry, const Intersection &intersection, const CellData &cellData, const bool first)
Function which calculates the flux entry.
Definition: 2p/sequential/diffusion/cellcentered/pressure.hh:719
FVPressure2P(Problem &problem)
Constructs a FVPressure2P object.
Definition: 2p/sequential/diffusion/cellcentered/pressure.hh:545
void updateMaterialLaws()
Updates and stores constitutive relations.
Definition: 2p/sequential/diffusion/cellcentered/pressure.hh:1092
void getFluxOnBoundary(EntryType &entry, const Intersection &intersection, const CellData &cellData, const bool first)
Function which calculates the boundary flux entry.
Definition: 2p/sequential/diffusion/cellcentered/pressure.hh:843
void initialize(bool solveTwice=true)
Initializes the pressure model.
Definition: 2p/sequential/diffusion/cellcentered/pressure.hh:199
void getSource(EntryType &entry, const Element &element, const CellData &cellData, const bool first)
Function which calculates the source entry.
Definition: 2p/sequential/diffusion/cellcentered/pressure.hh:611
void getStorage(EntryType &entry, const Element &element, const CellData &cellData, const bool first)
Function which calculates the storage entry.
Definition: 2p/sequential/diffusion/cellcentered/pressure.hh:655
void addOutputVtkFields(MultiWriter &writer)
Adds pressure output to the output file.
Definition: 2p/sequential/diffusion/cellcentered/pressure.hh:414
void updateVelocity()
Velocity update.
Definition: 2p/sequential/diffusion/cellcentered/pressure.hh:306
void update()
Pressure update.
Definition: 2p/sequential/diffusion/cellcentered/pressure.hh:262
The finite volume base class for the solution of a pressure equation.
Definition: sequential/cellcentered/pressure.hh:49
void assemble(bool first)
Function which assembles the system of equations to be solved.
Definition: sequential/cellcentered/pressure.hh:401
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
@ rhs
index for the right hand side entry
Definition: sequential/cellcentered/pressure.hh:88
@ matrix
index for the global matrix entry
Definition: sequential/cellcentered/pressure.hh:89
void solve()
Solves the global system of equations to get the spatial distribution of the pressure.
Definition: sequential/cellcentered/pressure.hh:527
void update()
Pressure update.
Definition: sequential/cellcentered/pressure.hh:228
Dune::FieldVector< Scalar, 2 > EntryType
Definition: sequential/cellcentered/pressure.hh:78
Specifies the properties for immiscible 2p diffusion/pressure models.
Finite Volume Diffusion Model.