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