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