3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
saturation.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_FVSATURATION2P_HH
25#define DUMUX_FVSATURATION2P_HH
26
29
30namespace Dumux {
31
74template<class TypeTag>
75class FVSaturation2P: public FVTransport<TypeTag>
76{
79
81
82 enum
83 {
84 dim = GridView::dimension, dimWorld = GridView::dimensionworld
85 };
86
89
93
95
97
100
102
104 using PrimaryVariables = typename SolutionTypes::PrimaryVariables;
105
107
108 enum
109 {
110 pw = Indices::pressureW,
111 pn = Indices::pressureNw,
112 pGlobal = Indices::pressureGlobal,
113 vw = Indices::velocityW,
114 vn = Indices::velocityNw,
115 vt = Indices::velocityTotal,
116 sw = Indices::saturationW,
117 sn = Indices::saturationNw
118 };
119 enum
120 {
121 wPhaseIdx = Indices::wPhaseIdx,
122 nPhaseIdx = Indices::nPhaseIdx,
123 saturationIdx = Indices::saturationIdx,
124 satEqIdx = Indices::satEqIdx,
125 numPhases = getPropValue<TypeTag, Properties::NumPhases>()
126 };
127
129
130 using Element = typename GridView::Traits::template Codim<0>::Entity;
131 using Intersection = typename GridView::Intersection;
132
133 using GlobalPosition = Dune::FieldVector<Scalar, dimWorld>;
134 using DimVector = Dune::FieldVector<Scalar, dim>;
135
136protected:
137 CapillaryFlux& capillaryFlux()
138 {
139 return *capillaryFlux_;
140 }
141
142 const CapillaryFlux& capillaryFlux() const
143 {
144 return *capillaryFlux_;
145 }
146
147 GravityFlux& gravityFlux()
148 {
149 return *gravityFlux_;
150 }
151
152 const GravityFlux& gravityFlux() const
153 {
154 return *gravityFlux_;
155 }
156
157public:
158 Velocity& velocity()
159 {
160 return *velocity_;
161 }
162
163 Velocity& velocity() const
164 {
165 return *velocity_;
166 }
167
169 void getFlux(Scalar& update, const Intersection& intersection, CellData& cellDataI);
170
172 void getFluxOnBoundary(Scalar& update, const Intersection& intersection, CellData& cellDataI);
173
175 void getSource(Scalar& update, const Element& element, CellData& cellDataI);
176
178 void initialize();
179
181 void updateMaterialLaws();
182
183
190 void getTransportedQuantity(TransportSolutionType& transportedQuantity)
191 {
192 int size = problem_.gridView().size(0);
193 transportedQuantity.resize(size);
194 for (int i = 0; i < size; i++)
195 {
196 switch (saturationType_)
197 {
198 case sw:
199 {
200 transportedQuantity[i] = problem_.variables().cellData(i).saturation(wPhaseIdx);
201 break;
202 }
203 case sn:
204 {
205 transportedQuantity[i] = problem_.variables().cellData(i).saturation(nPhaseIdx);
206 break;
207 }
208 }
209 }
210 }
211
218 void setTransportedQuantity(TransportSolutionType& transportedQuantity)
219 {
220 int size = problem_.gridView().size(0);
221 for (int i = 0; i < size; i++)
222 {
223 CellData& cellData = problem_.variables().cellData(i);
224 switch (saturationType_)
225 {
226 case sw:
227 {
228 Scalar sat = transportedQuantity[i];
229
230 cellData.setSaturation(wPhaseIdx, sat);
231 cellData.setSaturation(nPhaseIdx, 1 - sat);
232 break;
233 }
234 case sn:
235 {
236 Scalar sat = transportedQuantity[i];
237
238 cellData.setSaturation(wPhaseIdx,1 -sat);
239 cellData.setSaturation(nPhaseIdx, sat);
240 break;
241 }
242 }
243 }
244 }
245
252 template<class DataEntry>
253 bool inPhysicalRange(DataEntry& entry)
254 {
255 return (entry > -1e-6 && entry < 1.0 + 1e-6);
256 }
257
263 void updateTransportedQuantity(TransportSolutionType& updateVec)
264 {
265 if (this->enableLocalTimeStepping())
266 this->innerUpdate(updateVec);
267 else
268 asImp_().updateSaturationSolution(updateVec);
269 }
270
277 void updateTransportedQuantity(TransportSolutionType& updateVec, Scalar dt)
278 {
279 asImp_().updateSaturationSolution(updateVec, dt);
280 }
281
287 void updateSaturationSolution(TransportSolutionType& updateVec)
288 {
289 Scalar dt = problem_.timeManager().timeStepSize();
290 int size = problem_.gridView().size(0);
291 for (int i = 0; i < size; i++)
292 {
293 asImp_().updateSaturationSolution(i, updateVec[i][0], dt);
294 }
295 }
296
303 void updateSaturationSolution(TransportSolutionType& updateVec, Scalar dt)
304 {
305 int size = problem_.gridView().size(0);
306 for (int i = 0; i < size; i++)
307 {
308 asImp_().updateSaturationSolution(i, updateVec[i][0], dt);
309 }
310 }
311
321 void updateSaturationSolution(int eIdxGlobal, Scalar update, Scalar dt)
322 {
323 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
324
325 switch (saturationType_)
326 {
327 case sw:
328 {
329 Scalar sat = cellData.saturation(wPhaseIdx) + dt*update;
330
331 cellData.setSaturation(wPhaseIdx, sat);
332 cellData.setSaturation(nPhaseIdx, 1 - sat);
333 break;
334 }
335 case sn:
336 {
337 Scalar sat = cellData.saturation(nPhaseIdx) + dt*update;
338
339 cellData.setSaturation(wPhaseIdx,1 -sat);
340 cellData.setSaturation(nPhaseIdx, sat);
341 break;
342 }
343 }
344 }
345
358 template<class MultiWriter>
359 void addOutputVtkFields(MultiWriter &writer)
360 {
361 int size = problem_.gridView().size(0);
362 TransportSolutionType *saturation = writer.allocateManagedBuffer(size);
363 TransportSolutionType *saturationSecond = 0;
364
365 if (vtkOutputLevel_ > 0)
366 {
367 saturationSecond = writer.allocateManagedBuffer(size);
368 }
369
370 for (int i = 0; i < size; i++)
371 {
372 CellData& cellData = problem_.variables().cellData(i);
373
374 if (saturationType_ == sw)
375 {
376 (*saturation)[i] = cellData.saturation(wPhaseIdx);
377 if (vtkOutputLevel_ > 0)
378 {
379 (*saturationSecond)[i] = cellData.saturation(nPhaseIdx);
380 }
381 }
382 else if (saturationType_ == sn)
383 {
384 (*saturation)[i] = cellData.saturation(nPhaseIdx);
385 if (vtkOutputLevel_ > 0)
386 {
387 (*saturationSecond)[i] = cellData.saturation(wPhaseIdx);
388 }
389 }
390 }
391
392 if (saturationType_ == sw)
393 {
394 writer.attachCellData(*saturation, "wetting saturation");
395 if (vtkOutputLevel_ > 0)
396 {
397 writer.attachCellData(*saturationSecond, "nonwetting saturation");
398 }
399 }
400 else if (saturationType_ == sn)
401 {
402 writer.attachCellData(*saturation, "nonwetting saturation");
403 if (vtkOutputLevel_ > 0)
404 {
405 writer.attachCellData(*saturationSecond, "wetting saturation");
406 }
407 }
408
409 if (velocity().calculateVelocityInTransport())
410 velocity().addOutputVtkFields(writer);
411
412 return;
413 }
414
420 void serializeEntity(std::ostream &outstream, const Element &element)
421 {
422 int eIdxGlobal = problem_.variables().index(element);
423
424 Scalar sat = 0.0;
425 switch (saturationType_)
426 {
427 case sw:
428 sat = problem_.variables().cellData(eIdxGlobal).saturation(wPhaseIdx);
429 break;
430 case sn:
431 sat = problem_.variables().cellData(eIdxGlobal).saturation(nPhaseIdx);
432 break;
433 }
434
435 outstream << sat;
436 }
437
443 void deserializeEntity(std::istream &instream, const Element &element)
444 {
445 int eIdxGlobal = problem_.variables().index(element);
446
447 Scalar sat = 0.;
448 instream >> sat;
449
450 switch (saturationType_)
451 {
452 case sw:
453 problem_.variables().cellData(eIdxGlobal).setSaturation(wPhaseIdx, sat);
454 problem_.variables().cellData(eIdxGlobal).setSaturation(nPhaseIdx, 1-sat);
455 break;
456 case sn:
457 problem_.variables().cellData(eIdxGlobal).setSaturation(nPhaseIdx, sat);
458 problem_.variables().cellData(eIdxGlobal).setSaturation(wPhaseIdx, 1-sat);
459 break;
460 }
461 }
462
463
469 FVSaturation2P(Problem& problem) :
470 ParentType(problem), problem_(problem), threshold_(1e-6),
471 switchNormals_(getParam<bool>("Impet.SwitchNormals"))
472 {
473 if (compressibility_ && velocityType_ == vt)
474 {
475 DUNE_THROW(Dune::NotImplemented,
476 "Total velocity - global pressure - model cannot be used with compressible fluids!");
477 }
478 if (saturationType_ != sw && saturationType_ != sn)
479 {
480 DUNE_THROW(Dune::NotImplemented, "Saturation type not supported!");
481 }
482 if (pressureType_ != pw && pressureType_ != pn && pressureType_ != pGlobal)
483 {
484 DUNE_THROW(Dune::NotImplemented, "Pressure type not supported!");
485 }
486 if (velocityType_ != vw && velocityType_ != vn && velocityType_ != vt)
487 {
488 DUNE_THROW(Dune::NotImplemented, "Velocity type not supported!");
489 }
490
491 capillaryFlux_ = std::make_shared<CapillaryFlux>(problem);
492 gravityFlux_ = std::make_shared<GravityFlux>(problem);
493 velocity_ = std::make_shared<Velocity>(problem);
494
495 vtkOutputLevel_ = getParam<int>("Vtk.OutputLevel");
496 porosityThreshold_ = getParam<Scalar>("Impet.PorosityThreshold");
497 }
498
499private:
501 Implementation &asImp_()
502 { return *static_cast<Implementation *>(this); }
503
505 const Implementation &asImp_() const
506 { return *static_cast<const Implementation *>(this); }
507
508 Problem& problem_;
509 std::shared_ptr<Velocity> velocity_;
510 std::shared_ptr<CapillaryFlux> capillaryFlux_;
511 std::shared_ptr<GravityFlux> gravityFlux_;
512
513 int vtkOutputLevel_;
514 Scalar porosityThreshold_;
515
516 static const bool compressibility_ = getPropValue<TypeTag, Properties::EnableCompressibility>();
517 static const int saturationType_ = getPropValue<TypeTag, Properties::SaturationFormulation>();
518 static const int velocityType_ = getPropValue<TypeTag, Properties::VelocityFormulation>();
519 static const int pressureType_ = getPropValue<TypeTag, Properties::PressureFormulation>();
520
521 Scalar density_[numPhases];
522 Scalar viscosity_[numPhases];
523
524 const Scalar threshold_;
525 const bool switchNormals_;
526};
527
539template<class TypeTag>
540void FVSaturation2P<TypeTag>::getFlux(Scalar& update, const Intersection& intersection, CellData& cellDataI)
541{
542 auto elementI = intersection.inside();
543 auto elementJ = intersection.outside();
544
545 const CellData& cellDataJ = problem_.variables().cellData(problem_.variables().index(elementJ));
546
547 // get global coordinates of cell centers
548 const GlobalPosition& globalPosI = elementI.geometry().center();
549 const GlobalPosition& globalPosJ = elementJ.geometry().center();
550
551 // cell volume, assume linear map here
552 Scalar volume = elementI.geometry().volume();
553 using std::max;
554 Scalar porosity = max(problem_.spatialParams().porosity(elementI), porosityThreshold_);
555
556 if (compressibility_)
557 {
558 viscosity_[wPhaseIdx] = cellDataI.viscosity(wPhaseIdx);
559 viscosity_[nPhaseIdx] = cellDataI.viscosity(nPhaseIdx);
560 }
561
562 // local number of faces
563 int isIndex = intersection.indexInInside();
564
565 GlobalPosition unitOuterNormal = intersection.centerUnitOuterNormal();
566 if (switchNormals_)
567 unitOuterNormal *= -1.0;
568
569 Scalar faceArea = intersection.geometry().volume();
570
571 if (velocity().calculateVelocityInTransport() && !cellDataI.fluxData().haveVelocity(isIndex))
572 velocity().calculateVelocity(intersection, cellDataI);
573
574 //get velocity*normalvector*facearea/(volume*porosity)
575 Scalar factorW = (cellDataI.fluxData().velocity(wPhaseIdx, isIndex) * unitOuterNormal) * faceArea;
576 Scalar factorNw = (cellDataI.fluxData().velocity(nPhaseIdx, isIndex) * unitOuterNormal) * faceArea;
577 Scalar factorTotal = factorW + factorNw;
578
579 // distance vector between barycenters
580 GlobalPosition distVec = globalPosJ - globalPosI;
581 // compute distance between cell centers
582 Scalar dist = distVec.two_norm();
583
584 bool takeNeighbor = (elementI.level() < elementJ.level());
585 //get phase potentials
586 bool upwindWI =
587 (takeNeighbor) ? !cellDataJ.fluxData().isUpwindCell(wPhaseIdx, intersection.indexInOutside()) :
588 cellDataI.fluxData().isUpwindCell(wPhaseIdx, isIndex);
589 bool upwindNwI =
590 (takeNeighbor) ? !cellDataJ.fluxData().isUpwindCell(nPhaseIdx, intersection.indexInOutside()) :
591 cellDataI.fluxData().isUpwindCell(nPhaseIdx, isIndex);
592
593 Scalar lambdaW = 0;
594 Scalar lambdaNw = 0;
595
596 //upwinding of lambda dependend on the phase potential gradients
597 if (upwindWI)
598 {
599 lambdaW = cellDataI.mobility(wPhaseIdx);
600 if (compressibility_)
601 {
602 lambdaW /= cellDataI.density(wPhaseIdx);
603 } //divide by density because lambda is saved as lambda*density
604 }
605 else
606 {
607 lambdaW = cellDataJ.mobility(wPhaseIdx);
608 if (compressibility_)
609 {
610 lambdaW /= cellDataJ.density(wPhaseIdx);
611 } //divide by density because lambda is saved as lambda*density
612 }
613
614 if (upwindNwI)
615 {
616 lambdaNw = cellDataI.mobility(nPhaseIdx);
617 if (compressibility_)
618 {
619 lambdaNw /= cellDataI.density(nPhaseIdx);
620 } //divide by density because lambda is saved as lambda*density
621 }
622 else
623 {
624 lambdaNw = cellDataJ.mobility(nPhaseIdx);
625 if (compressibility_)
626 {
627 lambdaNw /= cellDataJ.density(nPhaseIdx);
628 } //divide by density because lambda is saved as lambda*density
629 }
630
631 switch (velocityType_)
632 {
633 case vt:
634 {
635 //add cflFlux for time-stepping
636 this->evalCflFluxFunction().addFlux(lambdaW, lambdaNw, viscosity_[wPhaseIdx],
637 viscosity_[nPhaseIdx], factorTotal, intersection);
638
639 //determine phase saturations from primary saturation variable
640 Scalar satWI = cellDataI.saturation(wPhaseIdx);
641 Scalar satWJ = cellDataJ.saturation(wPhaseIdx);
642
643 Scalar pcI = cellDataI.capillaryPressure();
644 Scalar pcJ = cellDataJ.capillaryPressure();
645
646 // calculate the saturation gradient
647 GlobalPosition pcGradient = unitOuterNormal;
648 pcGradient *= (pcI - pcJ) / dist;
649
650 // get the diffusive part
651 DimVector flux(0.);
652 capillaryFlux().getFlux(flux, intersection, satWI, satWJ, pcGradient);
653 Scalar capillaryFlux = (flux * unitOuterNormal * faceArea);
654
655 flux = 0.0;
656 gravityFlux().getFlux(flux, intersection, satWI, satWJ);
657 Scalar gravityFlux = (flux * unitOuterNormal * faceArea);
658
659 switch (saturationType_)
660 {
661 case sw:
662 {
663 //vt*fw
664 factorTotal *= lambdaW / (lambdaW + lambdaNw);
665 break;
666 }
667 case sn:
668 {
669 //vt*fn
670 factorTotal *= lambdaNw / (lambdaW + lambdaNw);
671 capillaryFlux *= -1;
672 gravityFlux *= -1;
673 break;
674 }
675 }
676 factorTotal -= capillaryFlux;
677 factorTotal += gravityFlux;
678
679 //add cflFlux for time-stepping
680 this->evalCflFluxFunction().addFlux(lambdaW, lambdaNw, viscosity_[wPhaseIdx],
681 viscosity_[nPhaseIdx], 10 * capillaryFlux, intersection);
682 this->evalCflFluxFunction().addFlux(lambdaW, lambdaNw, viscosity_[wPhaseIdx],
683 viscosity_[nPhaseIdx], 10 * gravityFlux, intersection);
684
685 break;
686 }
687 default:
688 {
689 if (compressibility_)
690 {
691 factorW /= cellDataI.density(wPhaseIdx);
692 factorNw /= cellDataI.density(nPhaseIdx);
693 }
694
695 //add cflFlux for time-stepping
696 this->evalCflFluxFunction().addFlux(lambdaW, lambdaNw, viscosity_[wPhaseIdx],
697 viscosity_[nPhaseIdx], factorW, intersection, wPhaseIdx);
698 this->evalCflFluxFunction().addFlux(lambdaW, lambdaNw, viscosity_[wPhaseIdx],
699 viscosity_[nPhaseIdx], factorNw, intersection, nPhaseIdx);
700
701 break;
702 }
703 }
704
705 switch (velocityType_)
706 {
707 case vt:
708 update -= factorTotal / (volume * porosity); //-:v>0, if flow leaves the cell
709 break;
710 default:
711 switch (saturationType_)
712 {
713 case sw:
714 update -= factorW / (volume * porosity);//-:v>0, if flow leaves the cell
715 break;
716 case sn:
717 update -= factorNw / (volume * porosity); //-:v>0, if flow leaves the cell
718 break;
719 }
720 break;
721 }
722}
723
732template<class TypeTag>
733void FVSaturation2P<TypeTag>::getFluxOnBoundary(Scalar& update, const Intersection& intersection, CellData& cellDataI)
734{
735 auto elementI = intersection.inside();
736
737 // get global coordinates of cell centers
738 const GlobalPosition& globalPosI = elementI.geometry().center();
739
740 // center of face in global coordinates
741 const GlobalPosition& globalPosJ = intersection.geometry().center();
742
743 // cell volume, assume linear map here
744 Scalar volume = elementI.geometry().volume();
745 using std::max;
746 Scalar porosity = max(problem_.spatialParams().porosity(elementI), porosityThreshold_);
747
748 if (compressibility_)
749 {
750 viscosity_[wPhaseIdx] = cellDataI.viscosity(wPhaseIdx);
751 viscosity_[nPhaseIdx] = cellDataI.viscosity(nPhaseIdx);
752 }
753
754 // local number of faces
755 int isIndex = intersection.indexInInside();
756
757 GlobalPosition unitOuterNormal = intersection.centerUnitOuterNormal();
758 if (switchNormals_)
759 unitOuterNormal *= -1.0;
760
761 Scalar faceArea = intersection.geometry().volume();
762
763 if (velocity().calculateVelocityInTransport())
764 velocity().calculateVelocityOnBoundary(intersection, cellDataI);
765
766 //get velocity*normalvector*facearea/(volume*porosity)
767 Scalar factorW = (cellDataI.fluxData().velocity(wPhaseIdx, isIndex) * unitOuterNormal) * faceArea;
768 Scalar factorNw = (cellDataI.fluxData().velocity(nPhaseIdx, isIndex) * unitOuterNormal) * faceArea;
769 Scalar factorTotal = factorW + factorNw;
770
771 // distance vector between barycenters
772 GlobalPosition distVec = globalPosJ - globalPosI;
773
774 // compute distance between cell centers
775 Scalar dist = distVec.two_norm();
776
777 //get boundary type
778 BoundaryTypes bcType;
779 problem_.boundaryTypes(bcType, intersection);
780 PrimaryVariables boundValues(0.0);
781
782 if (bcType.isDirichlet(satEqIdx))
783 {
784 problem_.dirichlet(boundValues, intersection);
785
786 Scalar satBound = boundValues[saturationIdx];
787
788 //determine phase saturations from primary saturation variable
789 Scalar satWI = cellDataI.saturation(wPhaseIdx);
790 Scalar satWBound = 0;
791 switch (saturationType_)
792 {
793 case sw:
794 {
795 satWBound = satBound;
796 break;
797 }
798 case sn:
799 {
800 satWBound = 1 - satBound;
801 break;
802 }
803 }
804
805 const auto fluidMatrixInteraction = problem_.spatialParams().fluidMatrixInteractionAtPos(elementI.geometry().center());
806 const Scalar pcBound = fluidMatrixInteraction.pc(satWBound);
807
808 Scalar lambdaW = 0;
809 Scalar lambdaNw = 0;
810
811 //upwinding of lambda dependend on the phase potential gradients
812 if (cellDataI.fluxData().isUpwindCell(wPhaseIdx, isIndex))
813 {
814 lambdaW = cellDataI.mobility(wPhaseIdx);
815 if (compressibility_)
816 {
817 lambdaW /= cellDataI.density(wPhaseIdx);
818 } //divide by density because lambda is saved as lambda*density
819 }
820 else
821 {
822 if (compressibility_)
823 {
824 lambdaW = fluidMatrixInteraction.krw(satWBound)
825 / FluidSystem::viscosity(cellDataI.fluidState(), wPhaseIdx);
826 }
827 else
828 {
829 lambdaW = fluidMatrixInteraction.krw(satWBound)
830 / viscosity_[wPhaseIdx];
831 }
832 }
833
834 if (cellDataI.fluxData().isUpwindCell(nPhaseIdx, isIndex))
835 {
836 lambdaNw = cellDataI.mobility(nPhaseIdx);
837 if (compressibility_)
838 {
839 lambdaNw /= cellDataI.density(nPhaseIdx);
840 } //divide by density because lambda is saved as lambda*density
841 }
842 else
843 {
844 if (compressibility_)
845 {
846 lambdaNw = fluidMatrixInteraction.krn(satWBound)
847 / FluidSystem::viscosity(cellDataI.fluidState(), nPhaseIdx);
848 }
849 else
850 {
851 lambdaNw = fluidMatrixInteraction.krn(satWBound)
852 / viscosity_[nPhaseIdx];
853 }
854 }
855
856 switch (velocityType_)
857 {
858 case vt:
859 {
860 this->evalCflFluxFunction().addFlux(lambdaW, lambdaNw, viscosity_[wPhaseIdx],
861 viscosity_[nPhaseIdx], factorTotal, intersection);
862
863 Scalar pcI = cellDataI.capillaryPressure();
864
865 // calculate the saturation gradient
866 GlobalPosition pcGradient = unitOuterNormal;
867 pcGradient *= (pcI - pcBound) / dist;
868
869 // get the diffusive part -> give 1-sat because sat = S_n and lambda = lambda(S_w) and pc = pc(S_w)
870 DimVector flux(0.);
871 capillaryFlux().getFlux(flux, intersection, satWI, satWBound, pcGradient);
872 Scalar capillaryFlux = flux * unitOuterNormal * faceArea;
873
874 flux = 0.0;
875 gravityFlux().getFlux(flux, intersection, satWI, satWBound);
876 Scalar gravityFlux = flux * unitOuterNormal * faceArea;
877
878 switch (saturationType_)
879 {
880 case sw:
881 {
882 //vt*fw
883 factorTotal *= lambdaW / (lambdaW + lambdaNw);
884 break;
885 }
886 case sn:
887 {
888 //vt*fn
889 factorTotal *= lambdaNw / (lambdaW + lambdaNw);
890 capillaryFlux *= -1; //add cflFlux for time-stepping
891 this->evalCflFluxFunction().addFlux(lambdaW, lambdaNw, viscosity_[wPhaseIdx],
892 viscosity_[nPhaseIdx], factorW, intersection, wPhaseIdx);
893 this->evalCflFluxFunction().addFlux(lambdaW, lambdaNw, viscosity_[wPhaseIdx],
894 viscosity_[nPhaseIdx], factorNw, intersection, nPhaseIdx);
895 gravityFlux *= -1;
896 break;
897 }
898 }
899 //vt*fw
900 factorTotal -= capillaryFlux;
901 factorTotal += gravityFlux;
902
903 //add cflFlux for time-stepping
904 this->evalCflFluxFunction().addFlux(lambdaW, lambdaNw, viscosity_[wPhaseIdx],
905 viscosity_[nPhaseIdx], 10 * capillaryFlux, intersection);
906 this->evalCflFluxFunction().addFlux(lambdaW, lambdaNw, viscosity_[wPhaseIdx],
907 viscosity_[nPhaseIdx], 10 * gravityFlux, intersection);
908
909 break;
910 }
911 default:
912 {
913 if (compressibility_)
914 {
915 factorW /= cellDataI.density(wPhaseIdx);
916 factorNw /= cellDataI.density(nPhaseIdx);
917 }
918
919 //add cflFlux for time-stepping
920 this->evalCflFluxFunction().addFlux(lambdaW, lambdaNw, viscosity_[wPhaseIdx],
921 viscosity_[nPhaseIdx], factorW, intersection, wPhaseIdx);
922 this->evalCflFluxFunction().addFlux(lambdaW, lambdaNw, viscosity_[wPhaseIdx],
923 viscosity_[nPhaseIdx], factorNw, intersection, nPhaseIdx);
924
925 break;
926 }
927 }
928 } //end dirichlet boundary
929
930 if (bcType.isNeumann(satEqIdx))
931 {
932 problem_.neumann(boundValues, intersection);
933 factorW = boundValues[wPhaseIdx];
934 factorNw = boundValues[nPhaseIdx];
935 factorW *= faceArea;
936 factorNw *= faceArea;
937
938 //get mobilities
939 Scalar lambdaW, lambdaNw;
940
941 lambdaW = cellDataI.mobility(wPhaseIdx);
942 lambdaNw = cellDataI.mobility(nPhaseIdx);
943 if (compressibility_)
944 {
945 lambdaW /= cellDataI.density(wPhaseIdx);
946 lambdaNw /= cellDataI.density(nPhaseIdx);
947 factorW /= cellDataI.density(wPhaseIdx);
948 factorNw /= cellDataI.density(nPhaseIdx);
949 }
950 else
951 {
952 factorW /= density_[wPhaseIdx];
953 factorNw /= density_[nPhaseIdx];
954 }
955
956 switch (velocityType_)
957 {
958 case vt:
959 {
960 //add cflFlux for time-stepping
961 this->evalCflFluxFunction().addFlux(lambdaW, lambdaNw, viscosity_[wPhaseIdx],
962 viscosity_[nPhaseIdx], factorW + factorNw, intersection);
963 break;
964 }
965 default:
966 {
967 this->evalCflFluxFunction().addFlux(lambdaW, lambdaNw, viscosity_[wPhaseIdx],
968 viscosity_[nPhaseIdx], factorW, intersection, wPhaseIdx);
969 this->evalCflFluxFunction().addFlux(lambdaW, lambdaNw, viscosity_[wPhaseIdx],
970 viscosity_[nPhaseIdx], factorNw, intersection, nPhaseIdx);
971
972 break;
973 }
974 }
975
976 } //end neumann boundary
977 if (bcType.isOutflow(satEqIdx))
978 {
979 //get mobilities
980 Scalar lambdaW = cellDataI.mobility(wPhaseIdx);
981 Scalar lambdaNw = cellDataI.mobility(nPhaseIdx);
982 if (compressibility_)
983 {
984 lambdaW /= cellDataI.density(wPhaseIdx);
985 lambdaNw /= cellDataI.density(nPhaseIdx);
986 }
987
988 if (velocityType_ == vt)
989 {
990 switch (saturationType_)
991 {
992 case sw:
993 {
994 //vt*fw
995 factorTotal *= lambdaW / (lambdaW + lambdaNw);
996 this->evalCflFluxFunction().addFlux(lambdaW, lambdaNw, viscosity_[wPhaseIdx],
997 viscosity_[nPhaseIdx], factorTotal, intersection);
998 break;
999 }
1000 case sn:
1001 {
1002 //vt*fn
1003 factorTotal *= lambdaNw / (lambdaW + lambdaNw);
1004 this->evalCflFluxFunction().addFlux(lambdaW, lambdaNw, viscosity_[wPhaseIdx],
1005 viscosity_[nPhaseIdx], factorTotal, intersection);
1006 break;
1007 }
1008 }
1009 }
1010 else
1011 {
1012 this->evalCflFluxFunction().addFlux(lambdaW, lambdaNw, viscosity_[wPhaseIdx],
1013 viscosity_[nPhaseIdx], factorW, intersection, wPhaseIdx);
1014 this->evalCflFluxFunction().addFlux(lambdaW, lambdaNw, viscosity_[wPhaseIdx],
1015 viscosity_[nPhaseIdx], factorNw, intersection, nPhaseIdx);
1016 }
1017 }
1018 switch (velocityType_)
1019 {
1020 case vt:
1021 update -= factorTotal / (volume * porosity); //-:v>0, if flow leaves the cell
1022 break;
1023 default:
1024 switch (saturationType_)
1025 {
1026 case sw:
1027 update -= factorW / (volume * porosity); //-:v>0, if flow leaves the cell
1028 break;
1029 case sn:
1030 update -= factorNw / (volume * porosity); //-:v>0, if flow leaves the cell
1031 break;
1032 }
1033 break;
1034 }
1035}
1036
1044template<class TypeTag>
1045void FVSaturation2P<TypeTag>::getSource(Scalar& update, const Element& element, CellData& cellDataI)
1046{
1047 // cell volume, assume linear map here
1048 Scalar volume = element.geometry().volume();
1049
1050 using std::max;
1051 Scalar porosity = max(problem_.spatialParams().porosity(element), porosityThreshold_);
1052
1053 if (compressibility_)
1054 {
1055 viscosity_[wPhaseIdx] = cellDataI.viscosity(wPhaseIdx);
1056 viscosity_[nPhaseIdx] = cellDataI.viscosity(nPhaseIdx);
1057 }
1058
1059 PrimaryVariables sourceVec(0.0);
1060 problem_.source(sourceVec, element);
1061
1062 if (compressibility_)
1063 {
1064 sourceVec[wPhaseIdx] /= cellDataI.density(wPhaseIdx);
1065 sourceVec[nPhaseIdx] /= cellDataI.density(nPhaseIdx);
1066 }
1067 else
1068 {
1069 sourceVec[wPhaseIdx] /= density_[wPhaseIdx];
1070 sourceVec[nPhaseIdx] /= density_[nPhaseIdx];
1071 }
1072
1073 //get mobilities
1074 Scalar lambdaW = cellDataI.mobility(wPhaseIdx);
1075 Scalar lambdaNw = cellDataI.mobility(nPhaseIdx);
1076 if (compressibility_)
1077 {
1078 lambdaW /= cellDataI.density(wPhaseIdx);
1079 lambdaNw /= cellDataI.density(nPhaseIdx);
1080 }
1081
1082 switch (saturationType_)
1083 {
1084 case sw:
1085 {
1086 if (sourceVec[wPhaseIdx] < 0 && cellDataI.saturation(wPhaseIdx) < threshold_)
1087 sourceVec[wPhaseIdx] = 0.0;
1088
1089 update += sourceVec[wPhaseIdx] / porosity;
1090 break;
1091 }
1092 case sn:
1093 {
1094 if (sourceVec[nPhaseIdx] < 0 && cellDataI.saturation(nPhaseIdx) < threshold_)
1095 sourceVec[nPhaseIdx] = 0.0;
1096
1097 update += sourceVec[nPhaseIdx] / porosity;
1098 break;
1099 }
1100 }
1101
1102 switch (velocityType_)
1103 {
1104 case vt:
1105 {
1106 //add cflFlux for time-stepping
1107 this->evalCflFluxFunction().addFlux(lambdaW, lambdaNw, viscosity_[wPhaseIdx], viscosity_[nPhaseIdx],
1108 (sourceVec[wPhaseIdx] + sourceVec[nPhaseIdx]) * -1 * volume, element);
1109 break;
1110 }
1111 default:
1112 {
1113 //add cflFlux for time-stepping
1114 this->evalCflFluxFunction().addFlux(lambdaW, lambdaNw, viscosity_[wPhaseIdx],
1115 viscosity_[nPhaseIdx], sourceVec[wPhaseIdx] * -1 * volume, element,
1116 wPhaseIdx);
1117 this->evalCflFluxFunction().addFlux(lambdaW, lambdaNw, viscosity_[wPhaseIdx],
1118 viscosity_[nPhaseIdx], sourceVec[nPhaseIdx] * -1 * volume, element,
1119 nPhaseIdx);
1120 break;
1121 }
1122 }
1123}
1124
1126template<class TypeTag>
1128{
1129 ParentType::initialize();
1130
1131 if (!compressibility_)
1132 {
1133 const auto element = *problem_.gridView().template begin<0>();
1134 FluidState fluidState;
1135 fluidState.setPressure(wPhaseIdx, problem_.referencePressure(element));
1136 fluidState.setPressure(nPhaseIdx, problem_.referencePressure(element));
1137 fluidState.setTemperature(problem_.temperature(element));
1138 fluidState.setSaturation(wPhaseIdx, 1.);
1139 fluidState.setSaturation(nPhaseIdx, 0.);
1140 density_[wPhaseIdx] = FluidSystem::density(fluidState, wPhaseIdx);
1141 density_[nPhaseIdx] = FluidSystem::density(fluidState, nPhaseIdx);
1142 viscosity_[wPhaseIdx] = FluidSystem::viscosity(fluidState, wPhaseIdx);
1143 viscosity_[nPhaseIdx] = FluidSystem::viscosity(fluidState, nPhaseIdx);
1144 }
1145
1146 // iterate through leaf grid an evaluate c0 at cell center
1147 for (const auto& element : elements(problem_.gridView()))
1148 {
1149 PrimaryVariables initSol(0.0);
1150 problem_.initial(initSol, element);
1151
1152 int eIdxGlobal = problem_.variables().index(element);
1153
1154 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
1155
1156 switch (saturationType_)
1157 {
1158 case sw:
1159 {
1160 cellData.setSaturation(wPhaseIdx, initSol[saturationIdx]);
1161 cellData.setSaturation(nPhaseIdx, 1 - initSol[saturationIdx]);
1162 break;
1163 }
1164 case sn:
1165 {
1166 cellData.setSaturation(wPhaseIdx,1 -initSol[saturationIdx]);
1167 cellData.setSaturation(nPhaseIdx, initSol[saturationIdx]);
1168
1169 break;
1170 }
1171 }
1172 }
1173
1174 velocity_->initialize();
1175 capillaryFlux_->initialize();
1176 gravityFlux_->initialize();
1177
1178 return;
1179}
1180
1186template<class TypeTag>
1188{
1189 // iterate through leaf grid an evaluate c0 at cell center
1190 for (const auto& element : elements(problem_.gridView()))
1191 {
1192 int eIdxGlobal = problem_.variables().index(element);
1193
1194 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
1195
1196 //determine phase saturations from primary saturation variable
1197 Scalar satW = cellData.saturation(wPhaseIdx);
1198
1199 const auto fluidMatrixInteraction = problem_.spatialParams().fluidMatrixInteractionAtPos(element.geometry().center());
1200 const Scalar pc = fluidMatrixInteraction.pc(satW);
1201
1202 cellData.setCapillaryPressure(pc);
1203
1204 // initialize mobilities
1205 const Scalar mobilityW = fluidMatrixInteraction.krw(satW) / viscosity_[wPhaseIdx];
1206 const Scalar mobilityNw = fluidMatrixInteraction.krn(satW) / viscosity_[nPhaseIdx];
1207
1208 // initialize mobilities
1209 cellData.setMobility(wPhaseIdx, mobilityW);
1210 cellData.setMobility(nPhaseIdx, mobilityNw);
1211
1212 //initialize fractional flow functions
1213 cellData.setFracFlowFunc(wPhaseIdx, mobilityW / (mobilityW + mobilityNw));
1214 cellData.setFracFlowFunc(nPhaseIdx, mobilityNw / (mobilityW + mobilityNw));
1215 }
1216 return;
1217}
1218
1219} // end namespace Dumux
1220#endif
Finite Volume discretization of a transport equation.
T getParam(Args &&... args)
A free function to get a parameter from the parameter tree singleton.
Definition: parameters.hh:348
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 saturation(int phaseIdx) noexcept
I/O name of saturation for multiphase systems.
Definition: name.hh:43
std::string viscosity(int phaseIdx) noexcept
I/O name of viscosity for multiphase systems.
Definition: name.hh:74
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
The finite volume discretization of a saturation transport equation.
Definition: saturation.hh:76
const CapillaryFlux & capillaryFlux() const
Definition: saturation.hh:142
void getFlux(Scalar &update, const Intersection &intersection, CellData &cellDataI)
Function which calculates the flux update.
Definition: saturation.hh:540
Velocity & velocity()
Definition: saturation.hh:158
GravityFlux & gravityFlux()
Definition: saturation.hh:147
void updateTransportedQuantity(TransportSolutionType &updateVec)
Updates the primary transport variable.
Definition: saturation.hh:263
void updateSaturationSolution(int eIdxGlobal, Scalar update, Scalar dt)
Updates the saturation solution of a cell.
Definition: saturation.hh:321
void getFluxOnBoundary(Scalar &update, const Intersection &intersection, CellData &cellDataI)
Function which calculates the boundary flux update.
Definition: saturation.hh:733
void updateSaturationSolution(TransportSolutionType &updateVec)
Globally updates the saturation solution.
Definition: saturation.hh:287
CapillaryFlux & capillaryFlux()
Definition: saturation.hh:137
Velocity & velocity() const
Definition: saturation.hh:163
void updateSaturationSolution(TransportSolutionType &updateVec, Scalar dt)
Globally updates the saturation solution.
Definition: saturation.hh:303
void initialize()
Sets the initial solution.
Definition: saturation.hh:1127
void getTransportedQuantity(TransportSolutionType &transportedQuantity)
Writes the current values of the primary transport variable into the transportedQuantity-vector (come...
Definition: saturation.hh:190
const GravityFlux & gravityFlux() const
Definition: saturation.hh:152
void serializeEntity(std::ostream &outstream, const Element &element)
Function for serialization of the primary transport variable.
Definition: saturation.hh:420
void deserializeEntity(std::istream &instream, const Element &element)
Function for deserialization of the primary transport variable.
Definition: saturation.hh:443
void updateMaterialLaws()
Update the values of the material laws and constitutive relations.
Definition: saturation.hh:1187
void addOutputVtkFields(MultiWriter &writer)
Adds saturation output to the output file.
Definition: saturation.hh:359
FVSaturation2P(Problem &problem)
Constructs a FVSaturation2P object.
Definition: saturation.hh:469
void getSource(Scalar &update, const Element &element, CellData &cellDataI)
Function which calculates the source update.
Definition: saturation.hh:1045
void setTransportedQuantity(TransportSolutionType &transportedQuantity)
Writes the current values of the primary transport variable into the variable container.
Definition: saturation.hh:218
void updateTransportedQuantity(TransportSolutionType &updateVec, Scalar dt)
Updates the primary transport variable.
Definition: saturation.hh:277
bool inPhysicalRange(DataEntry &entry)
Check if saturation is in physical range.
Definition: saturation.hh:253
The finite volume discretization of a transport equation.
Definition: transport.hh:52
bool enableLocalTimeStepping()
Definition: transport.hh:202
void update(const Scalar t, Scalar &dt, TransportSolutionType &updateVec, bool impet=false)
Calculate the update vector.
Definition: transport.hh:270
void innerUpdate(TransportSolutionType &updateVec)
Definition: transport.hh:564
Specifies the properties for immiscible 2p transport.