3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
fvtransport.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_FVTRANSPORT2P2C_HH
25#define DUMUX_FVTRANSPORT2P2C_HH
26
27#include <cmath>
28#include <unordered_map>
29
30#include <dune/grid/common/gridenums.hh>
31#include <dune/common/float_cmp.hh>
32
35#include <dumux/common/math.hh>
37
39
40namespace Dumux {
60template<class TypeTag>
62{
67
69
72
75
77
79
80 enum
81 {
82 dim = GridView::dimension, dimWorld = GridView::dimensionworld
83 };
84 enum
85 {
86 pw = Indices::pressureW,
87 pn = Indices::pressureN
88 };
89 enum
90 {
91 wPhaseIdx = Indices::wPhaseIdx, nPhaseIdx = Indices::nPhaseIdx,
92 wCompIdx = Indices::wPhaseIdx, nCompIdx = Indices::nPhaseIdx,
93 contiWEqIdx=Indices::contiWEqIdx, contiNEqIdx=Indices::contiNEqIdx,
94 NumPhases = getPropValue<TypeTag, Properties::NumPhases>(),
95 NumComponents = getPropValue<TypeTag, Properties::NumComponents>()
96 };
97
98 using Element = typename GridView::Traits::template Codim<0>::Entity;
99 using Intersection = typename GridView::Intersection;
100
101 using GlobalPosition = Dune::FieldVector<Scalar, dimWorld>;
102 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
103 using PhaseVector = Dune::FieldVector<Scalar, NumPhases>;
104 using ComponentVector = Dune::FieldVector<Scalar, NumComponents>;
106
107public:
109 using EntryType = Dune::FieldVector<Scalar, 2>;
110 using TimeStepFluxType = Dune::FieldVector<Scalar, 2>;
111
112protected:
114 {
115 Dune::FieldVector<EntryType, 2*dim> faceFluxes;
116 Dune::FieldVector<Scalar, 2*dim> faceTargetDt;
117 Scalar dt;
119 {}
120 };
121
123 Problem& problem()
124 { return problem_; }
125
126 void innerUpdate(TransportSolutionType& updateVec);
127
128public:
129 virtual void update(const Scalar t, Scalar& dt, TransportSolutionType& updateVec, bool impet = false);
130
131 void updateTransportedQuantity(TransportSolutionType& updateVector);
132
133 void updateTransportedQuantity(TransportSolutionType& updateVector, Scalar dt);
134
135 void updateConcentrations(TransportSolutionType& updateVector, Scalar dt);
136
137 // Function which calculates the flux update
138 void getFlux(ComponentVector& fluxEntries, EntryType& timestepFlux,
139 const Intersection& intersection, CellData& cellDataI);
140
141 // Function which calculates the boundary flux update
142 void getFluxOnBoundary(ComponentVector& fluxEntries, EntryType& timestepFlux,
143 const Intersection& intersection, const CellData& cellDataI);
144
145 void evalBoundary(GlobalPosition globalPosFace,
146 const Intersection& intersection,
147 FluidState& BCfluidState,
148 PhaseVector& pressBound);
149
150
156 {
157 // resize update vector and set to zero
158 int transportedQuantities = getPropValue<TypeTag, Properties::NumEq>() - 1; // NumEq - 1 pressure Eq
159 totalConcentration_.resize(transportedQuantities);
160 for (int eqNumber = 0; eqNumber < transportedQuantities; eqNumber++)
161 {
162 totalConcentration_[eqNumber].resize(problem().gridView().size(0));
163 totalConcentration_[eqNumber] = 0;
164 }
165 }
166
171 template<class MultiWriter>
172 void addOutputVtkFields(MultiWriter &writer)
173 {
174 if(problem().vtkOutputLevel()>3)
175 {
176 using ScalarSolutionType = typename GetProp<TypeTag, Properties::SolutionTypes>::ScalarSolution;
177 int size = problem_.gridView().size(0);
178 ScalarSolutionType *totalC1PV = writer.allocateManagedBuffer(size);
179 ScalarSolutionType *totalC2PV = writer.allocateManagedBuffer(size);
180 *totalC1PV = this->totalConcentration_[wCompIdx];
181 *totalC2PV = this->totalConcentration_[nCompIdx];
182 writer.attachCellData(*totalC1PV, "total Concentration w-Comp - vector");
183 writer.attachCellData(*totalC2PV, "total Concentration n-Comp - vector");
184 }
185 }
186
188 void serializeEntity(std::ostream &outstream, const Element &element)
189 {
190 int eIdxGlobal = problem().variables().index(element);
191 outstream << totalConcentration_[wCompIdx][eIdxGlobal]
192 << " " << totalConcentration_[nCompIdx][eIdxGlobal];
193 }
195 void deserializeEntity(std::istream &instream, const Element &element)
196 {
197 int eIdxGlobal = problem().variables().index(element);
198 CellData& cellData = problem().variables().cellData(eIdxGlobal);
199 instream >> totalConcentration_[wCompIdx][eIdxGlobal]
200 >> totalConcentration_[nCompIdx][eIdxGlobal];
201 cellData.setMassConcentration(wCompIdx, totalConcentration_[wCompIdx][eIdxGlobal]);
202 cellData.setMassConcentration(nCompIdx, totalConcentration_[nCompIdx][eIdxGlobal]);
203 }
204
213 void getTransportedQuantity(TransportSolutionType& transportedQuantity)
214 {
215 // resize update vector and set to zero
216 transportedQuantity.resize((getPropValue<TypeTag, Properties::NumEq>() - 1));
217 for(int compIdx = 0; compIdx < (getPropValue<TypeTag, Properties::NumEq>() - 1); compIdx++)
218 transportedQuantity[compIdx].resize(problem_.gridView().size(0));
219
220 transportedQuantity = totalConcentration_;
221 }
222
232 Scalar& totalConcentration(int compIdx, int eIdxGlobal)
233 {
234 return totalConcentration_[compIdx][eIdxGlobal][0];
235 }
236
237 void getSource(Scalar& update, const Element& element, CellData& cellDataI)
238 {}
239
240
246 template<class DataEntry>
247 bool inPhysicalRange(DataEntry& entry)
248 {
249 int numComp = getPropValue<TypeTag, Properties::NumEq>() - 1;
250 for(int compIdx = 0; compIdx < numComp; compIdx++)
251 {
252 if (entry[compIdx] < -1.0e-6)
253 {
254 return false;
255 }
256 }
257 return true;
258 }
259
262 {
263 return localTimeStepping_;
264 }
265
267
276 switchNormals(getParam<bool>("Impet.SwitchNormals")), accumulatedDt_(0),
277 dtThreshold_(1e-6), subCFLFactor_(1.0)
278 {
279 restrictFluxInTransport_ = getParam<int>("Impet.RestrictFluxInTransport", 0);
280 regulateBoundaryPermeability = getPropValue<TypeTag, Properties::RegulateBoundaryPermeability>();
282 minimalBoundaryPermeability = getParam<Scalar>("SpatialParams.MinBoundaryPermeability");
283
284 Scalar cFLFactor = getParam<Scalar>("Impet.CFLFactor");
285 using std::min;
286 subCFLFactor_ = min(getParam<Scalar>("Impet.SubCFLFactor"), cFLFactor);
287 verbosity_ = getParam<int>("TimeManager.SubTimestepVerbosity");
288
290
292 std::cout<<"max CFL-Number of "<<cFLFactor<<", max Sub-CFL-Number of "
293 <<subCFLFactor_<<": Enable local time-stepping!" << std::endl;
294 }
295
297 {
298 }
299
300protected:
301 TransportSolutionType totalConcentration_;
302 Problem& problem_;
303 bool impet_;
305
307 static const int pressureType = getPropValue<TypeTag, Properties::PressureFormulation>();
318 const Scalar dtThreshold_;
320 std::vector<LocalTimesteppingData> timeStepData_;
321
322 void updatedTargetDt_(Scalar &dt);
323
325 {
326 timeStepData_.clear();
327 accumulatedDt_ = 0;
328 }
329
333
334
335private:
337 Implementation &asImp_()
338 { return *static_cast<Implementation *>(this); }
339
341 const Implementation &asImp_() const
342 { return *static_cast<const Implementation *>(this); }
343};
344
345
363template<class TypeTag>
364void FVTransport2P2C<TypeTag>::update(const Scalar t, Scalar& dt,
365 TransportSolutionType& updateVec, bool impet)
366{
367 // initialize dt very large
368 dt = 1E100;
369
370 unsigned int size = problem_.gridView().size(0);
371 if (localTimeStepping_)
372 {
373 if (this->timeStepData_.size() != size)
374 this->timeStepData_.resize(size);
375 }
376 // store if we do update Estimate for flux functions
377 impet_ = impet;
378 averagedFaces_ = 0;
379
380 // resize update vector and set to zero
381 updateVec.resize(getPropValue<TypeTag, Properties::NumComponents>());
382 updateVec[wCompIdx].resize(problem_.gridView().size(0));
383 updateVec[nCompIdx].resize(problem_.gridView().size(0));
384 updateVec[wCompIdx] = 0;
385 updateVec[nCompIdx] = 0;
386
387 // Cell which restricts time step size
388 int restrictingCell = -1;
389
390 ComponentVector entries(0.);
391 EntryType timestepFlux(0.);
392 // compute update vector
393 for (const auto& element : elements(problem().gridView()))
394 {
395 // get cell infos
396 int eIdxGlobalI = problem().variables().index(element);
397 CellData& cellDataI = problem().variables().cellData(eIdxGlobalI);
398
399 // some variables for time step calculation
400 double sumfactorin = 0;
401 double sumfactorout = 0;
402
403 // run through all intersections with neighbors and boundary
404 for (const auto& intersection : intersections(problem().gridView(), element))
405 {
406 int indexInInside = intersection.indexInInside();
407
408 /****** interior face *****************/
409 if (intersection.neighbor())
410 asImp_().getFlux(entries, timestepFlux, intersection, cellDataI);
411
412 /****** Boundary Face *****************/
413 if (intersection.boundary())
414 asImp_().getFluxOnBoundary(entries, timestepFlux, intersection, cellDataI);
415
416
417 if (localTimeStepping_)
418 {
419 LocalTimesteppingData& localData = timeStepData_[eIdxGlobalI];
420 if (localData.faceTargetDt[indexInInside] < accumulatedDt_ + dtThreshold_)
421 {
422 localData.faceFluxes[indexInInside] = entries;
423 }
424 }
425 else
426 {
427 // add to update vector
428 updateVec[wCompIdx][eIdxGlobalI] += entries[wCompIdx];
429 updateVec[nCompIdx][eIdxGlobalI] += entries[nCompIdx];
430 }
431
432 // for time step calculation
433 sumfactorin += timestepFlux[0];
434 sumfactorout += timestepFlux[1];
435
436 }// end all intersections
437
438 if (localTimeStepping_)
439 {
440 LocalTimesteppingData& localData = timeStepData_[eIdxGlobalI];
441 for (int i=0; i < 2*dim; i++)
442 {
443 updateVec[wCompIdx][eIdxGlobalI] += localData.faceFluxes[i][wCompIdx];
444 updateVec[nCompIdx][eIdxGlobalI] += localData.faceFluxes[i][nCompIdx];
445 }
446 }
447
448 /*********** Handle source term ***************/
449 PrimaryVariables q(NAN);
450 problem().source(q, element);
451 updateVec[wCompIdx][eIdxGlobalI] += q[contiWEqIdx];
452 updateVec[nCompIdx][eIdxGlobalI] += q[contiNEqIdx];
453
454 // account for porosity in fluxes for time-step
455 using std::max;
456 sumfactorin = max(sumfactorin,sumfactorout)
457 / problem().spatialParams().porosity(element);
458
459 //calculate time step
460 if (localTimeStepping_)
461 {
462 timeStepData_[eIdxGlobalI].dt = 1./sumfactorin;
463 if ( 1./sumfactorin < dt)
464 {
465 dt = 1./sumfactorin;
466 restrictingCell= eIdxGlobalI;
467 }
468 }
469 else
470 {
471 if ( 1./sumfactorin < dt)
472 {
473 dt = 1./sumfactorin;
474 restrictingCell= eIdxGlobalI;
475 }
476 }
477 } // end grid traversal
478
479#if HAVE_MPI
480 // communicate updated values
482 using ElementMapper = typename SolutionTypes::ElementMapper;
484 for (int i = 0; i < updateVec.size(); i++)
485 {
486 DataHandle dataHandle(problem_.variables().elementMapper(), updateVec[i]);
487 problem_.gridView().template communicate<DataHandle>(dataHandle,
488 Dune::InteriorBorder_All_Interface,
489 Dune::ForwardCommunication);
490 }
491
492 if (localTimeStepping_)
493 {
495
496 TimeDataHandle timeDataHandle(problem_.elementMapper(), timeStepData_);
497 problem_.gridView().template communicate<TimeDataHandle>(timeDataHandle,
498 Dune::InteriorBorder_All_Interface,
499 Dune::ForwardCommunication);
500 }
501
502 dt = problem_.gridView().comm().min(dt);
503#endif
504
505 if(impet)
506 {
507 Dune::dinfo << "Timestep restricted by CellIdx " << restrictingCell << " leads to dt = "
508 <<dt * getParam<Scalar>("Impet.CFLFactor")<< std::endl;
509 if (averagedFaces_ != 0)
510 Dune::dinfo << " Averageing done for " << averagedFaces_ << " faces. "<< std::endl;
511 }
512}
519template<class TypeTag>
520void FVTransport2P2C<TypeTag>::updateTransportedQuantity(TransportSolutionType& updateVector)
521{
522 if (this->enableLocalTimeStepping())
523 this->innerUpdate(updateVector);
524 else
525 updateConcentrations(updateVector, problem().timeManager().timeStepSize());
526}
527
535template<class TypeTag>
536void FVTransport2P2C<TypeTag>::updateTransportedQuantity(TransportSolutionType& updateVector, Scalar dt)
537{
538 updateConcentrations(updateVector, dt);
539}
540
548template<class TypeTag>
549void FVTransport2P2C<TypeTag>::updateConcentrations(TransportSolutionType& updateVector, Scalar dt)
550{
551 // loop thorugh all elements
552 for (int i = 0; i< problem().gridView().size(0); i++)
553 {
554 CellData& cellDataI = problem().variables().cellData(i);
555 for(int compIdx = 0; compIdx < getPropValue<TypeTag, Properties::NumComponents>(); compIdx++)
556 {
557 totalConcentration_[compIdx][i] += (updateVector[compIdx][i]*=dt);
558 cellDataI.setMassConcentration(compIdx, totalConcentration_[compIdx][i]);
559 }
560 }
561}
562
579template<class TypeTag>
580void FVTransport2P2C<TypeTag>::getFlux(ComponentVector& fluxEntries,
581 Dune::FieldVector<Scalar, 2>& timestepFlux,
582 const Intersection& intersection,
583 CellData& cellDataI)
584{
585 fluxEntries = 0.;
586 timestepFlux = 0.;
587 // cell information
588 auto elementI = intersection.inside();
589 int eIdxGlobalI = problem().variables().index(elementI);
590
591 // get position
592 const GlobalPosition globalPos = elementI.geometry().center();
593 const GlobalPosition& gravity_ = problem().gravity();
594 // cell volume, assume linear map here
595 Scalar volume = elementI.geometry().volume();
596
597 // get values of cell I
598 Scalar pressI = problem().pressureModel().pressure(eIdxGlobalI);
599 Scalar pcI = cellDataI.capillaryPressure();
600 DimMatrix K_I(problem().spatialParams().intrinsicPermeability(elementI));
601
602 // old material law interface is deprecated: Replace this by
603 // const auto& fluidMatrixInteraction = problem().spatialParams.fluidMatrixInteractionAtPos(elementI.geometry().center());
604 // after the release of 3.3, when the deprecated interface is no longer supported
605 const auto fluidMatrixInteraction = Deprecated::makePcKrSw(Scalar{}, problem().spatialParams(), elementI);
606
607 PhaseVector SmobI(0.);
608 using std::max;
609 SmobI[wPhaseIdx] = max((cellDataI.saturation(wPhaseIdx)
610 - fluidMatrixInteraction.pcSwCurve().effToAbsParams().swr())
611 , 1e-2);
612 SmobI[nPhaseIdx] = max((cellDataI.saturation(nPhaseIdx)
613 - fluidMatrixInteraction.pcSwCurve().effToAbsParams().snr())
614 , 1e-2);
615
616 Scalar densityWI (0.), densityNWI(0.);
617 densityWI= cellDataI.density(wPhaseIdx);
618 densityNWI = cellDataI.density(nPhaseIdx);
619
620 // face properties
621 GlobalPosition unitOuterNormal = intersection.centerUnitOuterNormal();
622 if (switchNormals)
623 unitOuterNormal *= -1.0;
624 Scalar faceArea = intersection.geometry().volume();
625
626// // local interface index
627// const int indexInInside = intersection.indexInInside();
628// const int indexInOutside = intersection.indexInOutside();
629
630 // create vector for timestep and for update
631 Dune::FieldVector<Scalar, 2> factor (0.);
632 Dune::FieldVector<Scalar, 2> updFactor (0.);
633
634 PhaseVector potential(0.);
635
636 // access neighbor
637 auto neighbor = intersection.outside();
638 int eIdxGlobalJ = problem().variables().index(neighbor);
639 CellData& cellDataJ = problem().variables().cellData(eIdxGlobalJ);
640
641 // neighbor cell center in global coordinates
642 const GlobalPosition& globalPosNeighbor = neighbor.geometry().center();
643
644 // distance vector between barycenters
645 GlobalPosition distVec = globalPosNeighbor - globalPos;
646 // compute distance between cell centers
647 Scalar dist = distVec.two_norm();
648
649 GlobalPosition unitDistVec(distVec);
650 unitDistVec /= dist;
651
652 // phase densities in neighbor
653 Scalar densityWJ (0.), densityNWJ(0.);
654 densityWJ = cellDataJ.density(wPhaseIdx);
655 densityNWJ = cellDataJ.density(nPhaseIdx);
656
657 // average phase densities with central weighting
658 double densityW_mean = (densityWI + densityWJ) * 0.5;
659 double densityNW_mean = (densityNWI + densityNWJ) * 0.5;
660
661 double pressJ = problem().pressureModel().pressure(eIdxGlobalJ);
662 Scalar pcJ = cellDataJ.capillaryPressure();
663
664 // compute mean permeability
665 DimMatrix meanK_(0.);
666 harmonicMeanMatrix(meanK_,
667 K_I,
668 problem().spatialParams().intrinsicPermeability(neighbor));
669 Dune::FieldVector<Scalar,dim> K(0);
670 meanK_.umv(unitDistVec,K);
671
672 // determine potentials for upwind
673 switch (pressureType)
674 {
675 case pw:
676 {
677 potential[wPhaseIdx] = (pressI - pressJ) / (dist);
678 potential[nPhaseIdx] = (pressI - pressJ + pcI - pcJ) / (dist);
679 break;
680 }
681 case pn:
682 {
683 potential[wPhaseIdx] = (pressI - pressJ - pcI + pcJ) / (dist);
684 potential[nPhaseIdx] = (pressI - pressJ) / (dist);
685 break;
686 }
687 }
688 // add gravity term
689 potential[nPhaseIdx] += (gravity_ * unitDistVec) * densityNW_mean;
690 potential[wPhaseIdx] += (gravity_ * unitDistVec) * densityW_mean;
691
692 potential[wPhaseIdx] *= fabs(K * unitOuterNormal);
693 potential[nPhaseIdx] *= fabs(K * unitOuterNormal);
694
695 // determine upwinding direction, perform upwinding if possible
696 Dune::FieldVector<bool, NumPhases> doUpwinding(true);
697 PhaseVector lambda(0.);
698 for(int phaseIdx = 0; phaseIdx < NumPhases; phaseIdx++)
699 {
700 int contiEqIdx = 0;
701 if(phaseIdx == wPhaseIdx)
702 contiEqIdx = contiWEqIdx;
703 else
704 contiEqIdx = contiNEqIdx;
705
706 if(!impet_ || restrictFluxInTransport_==0) // perform a strict uwpind scheme
707 {
708 if (potential[phaseIdx] > 0.)
709 {
710 lambda[phaseIdx] = cellDataI.mobility(phaseIdx);
711 cellDataI.setUpwindCell(intersection.indexInInside(), contiEqIdx, true);
712
713 }
714 else if(potential[phaseIdx] < 0.)
715 {
716 lambda[phaseIdx] = cellDataJ.mobility(phaseIdx);
717 cellDataI.setUpwindCell(intersection.indexInInside(), contiEqIdx, false);
718 }
719 else
720 {
721 doUpwinding[phaseIdx] = false;
722 cellDataI.setUpwindCell(intersection.indexInInside(), contiEqIdx, false);
723 cellDataJ.setUpwindCell(intersection.indexInOutside(), contiEqIdx, false);
724 }
725 }
726 else // Transport after PE with check on flow direction
727 {
728 bool wasUpwindCell = cellDataI.isUpwindCell(intersection.indexInInside(), contiEqIdx);
729
730 if (potential[phaseIdx] > 0. && wasUpwindCell)
731 lambda[phaseIdx] = cellDataI.mobility(phaseIdx);
732 else if (potential[phaseIdx] < 0. && !wasUpwindCell)
733 lambda[phaseIdx] = cellDataJ.mobility(phaseIdx);
734 // potential direction does not coincide with that of P.E.
735 else if(restrictFluxInTransport_ == 2) // perform central averaging for all direction changes
736 doUpwinding[phaseIdx] = false;
737 else // i.e. restrictFluxInTransport == 1
738 {
739 //check if harmonic weighting is necessary
740 if (potential[phaseIdx] > 0. && (Dune::FloatCmp::ne<Scalar, Dune::FloatCmp::absolute>(cellDataJ.mobility(phaseIdx), 0.0, 1.0e-30) // check if outflow induce neglected (i.e. mob=0) phase flux
741 || (cellDataI.wasRefined() && cellDataJ.wasRefined() && elementI.father() == neighbor.father())))
742 lambda[phaseIdx] = cellDataI.mobility(phaseIdx);
743 else if (potential[phaseIdx] < 0. && (Dune::FloatCmp::ne<Scalar, Dune::FloatCmp::absolute>(cellDataI.mobility(phaseIdx), 0.0, 1.0e-30) // check if inflow induce neglected phase flux
744 || (cellDataI.wasRefined() && cellDataJ.wasRefined() && elementI.father() == neighbor.father())))
745 lambda[phaseIdx] = cellDataJ.mobility(phaseIdx);
746 else
747 doUpwinding[phaseIdx] = false;
748 }
749
750 }
751
752 // do not perform upwinding if so desired
753 if(!doUpwinding[phaseIdx])
754 {
755 //a) no flux if there wont be any flux regardless how to average/upwind
756 if(cellDataI.mobility(phaseIdx)+cellDataJ.mobility(phaseIdx)==0.)
757 {
758 potential[phaseIdx] = 0;
759 continue;
760 }
761
762 //b) perform harmonic averaging
763 fluxEntries[wCompIdx] -= potential[phaseIdx] * faceArea / volume
764 * harmonicMean(cellDataI.massFraction(phaseIdx, wCompIdx) * cellDataI.mobility(phaseIdx) * cellDataI.density(phaseIdx),
765 cellDataJ.massFraction(phaseIdx, wCompIdx) * cellDataJ.mobility(phaseIdx) * cellDataJ.density(phaseIdx));
766 fluxEntries[nCompIdx] -= potential[phaseIdx] * faceArea / volume
767 * harmonicMean(cellDataI.massFraction(phaseIdx, nCompIdx) * cellDataI.mobility(phaseIdx) * cellDataI.density(phaseIdx),
768 cellDataJ.massFraction(phaseIdx, nCompIdx) * cellDataJ.mobility(phaseIdx) * cellDataJ.density(phaseIdx));
769 // c) timestep control
770 // for timestep control : influx
771 using std::max;
772 timestepFlux[0] += max(0.,
773 - potential[phaseIdx] * faceArea / volume
774 * harmonicMean(cellDataI.mobility(phaseIdx),cellDataJ.mobility(phaseIdx)));
775 // outflux
776 timestepFlux[1] += max(0.,
777 potential[phaseIdx] * faceArea / volume
778 * harmonicMean(cellDataI.mobility(phaseIdx),cellDataJ.mobility(phaseIdx))/SmobI[phaseIdx]);
779
780 //d) output
781 if(!(cellDataI.wasRefined() && cellDataJ.wasRefined() && elementI.father() == neighbor.father())
782 && eIdxGlobalI > eIdxGlobalJ) //(only for one side)
783 {
784 averagedFaces_++;
785 #if DUNE_MINIMAL_DEBUG_LEVEL < 3
786 // verbose (only for one side)
787 if(eIdxGlobalI > eIdxGlobalJ)
788 Dune::dinfo << "harmonicMean flux of phase" << phaseIdx <<" used from cell" << eIdxGlobalI<< " into " << eIdxGlobalJ
789 << " ; TE upwind I = "<< cellDataI.isUpwindCell(intersection.indexInInside(), contiEqIdx)
790 << " but pot = "<< potential[phaseIdx] << std::endl;
791 #endif
792 }
793
794 //e) stop further standard calculations
795 potential[phaseIdx] = 0;
796 }
797 }
798
799 // calculate and standardized velocity
800 using std::max;
801 double velocityJIw = max((-lambda[wPhaseIdx] * potential[wPhaseIdx]) * faceArea / volume, 0.0);
802 double velocityIJw = max(( lambda[wPhaseIdx] * potential[wPhaseIdx]) * faceArea / volume, 0.0);
803 double velocityJIn = max((-lambda[nPhaseIdx] * potential[nPhaseIdx]) * faceArea / volume, 0.0);
804 double velocityIJn = max(( lambda[nPhaseIdx] * potential[nPhaseIdx]) * faceArea / volume, 0.0);
805
806 // for timestep control : influx
807 timestepFlux[0] += velocityJIw + velocityJIn;
808
809 double foutw = velocityIJw/SmobI[wPhaseIdx];
810 double foutn = velocityIJn/SmobI[nPhaseIdx];
811 using std::isnan;
812 using std::isinf;
813 if (isnan(foutw) || isinf(foutw) || foutw < 0) foutw = 0;
814 if (isnan(foutn) || isinf(foutn) || foutn < 0) foutn = 0;
815 timestepFlux[1] += foutw + foutn;
816
817 fluxEntries[wCompIdx] +=
818 velocityJIw * cellDataJ.massFraction(wPhaseIdx, wCompIdx) * densityWJ
819 - velocityIJw * cellDataI.massFraction(wPhaseIdx, wCompIdx) * densityWI
820 + velocityJIn * cellDataJ.massFraction(nPhaseIdx, wCompIdx) * densityNWJ
821 - velocityIJn * cellDataI.massFraction(nPhaseIdx, wCompIdx) * densityNWI;
822 fluxEntries[nCompIdx] +=
823 velocityJIw * cellDataJ.massFraction(wPhaseIdx, nCompIdx) * densityWJ
824 - velocityIJw * cellDataI.massFraction(wPhaseIdx, nCompIdx) * densityWI
825 + velocityJIn * cellDataJ.massFraction(nPhaseIdx, nCompIdx) * densityNWJ
826 - velocityIJn * cellDataI.massFraction(nPhaseIdx, nCompIdx) * densityNWI;
827
828
829// //add dispersion
830// for (int phaseIdx = 0; phaseIdx < NumPhases; ++phaseIdx)
831// {
832//
833// //calculate porous media diffusion coefficient
834// // for phases which exist in both finite volumes
835// if (cellDataI.saturation(phaseIdx) <= 0 ||
836// cellDataJ.saturation(phaseIdx) <= 0)
837// {
838// continue;
839// }
840//
841// // calculate tortuosity at the nodes i and j needed
842// // for porous media diffusion coefficient
843// Scalar poroI = problem().spatialParams().porosity(elementI);
844// Scalar poroJ = problem().spatialParams().porosity(neighbor);
845// Scalar tauI =
846// 1.0/(poroI * poroI) *
847// pow(poroI * cellDataI.saturation(phaseIdx), 7.0/3);
848// Scalar tauJ =
849// 1.0/(poroJ * poroJ) *
850// pow(poroJ * cellDataJ.saturation(phaseIdx), 7.0/3);
851// // Diffusion coefficient in the porous medium
852// // -> harmonic mean
853// Scalar porousDiffCoeff
854// = harmonicMean(poroI * cellDataI.saturation(phaseIdx) * tauI
855// * FluidSystem::binaryDiffusionCoefficient(cellDataI.fluidState(), phaseIdx, wCompIdx, nCompIdx),
856// poroJ * cellDataJ.saturation(phaseIdx) * tauJ
857// * FluidSystem::binaryDiffusionCoefficient(cellDataJ.fluidState(), phaseIdx, wCompIdx, nCompIdx));
858// Scalar averagedMolarDensity = (cellDataI.density(phaseIdx) / cellDataI.fluidState().averageMolarMass(phaseIdx)
859// +cellDataJ.density(phaseIdx) / cellDataJ.fluidState().averageMolarMass(phaseIdx))*0.5;
860//
861// for (int compIdx = 0; compIdx < NumComponents; ++compIdx)
862// {
863//
864// Scalar gradx = (cellDataJ.moleFraction(phaseIdx, compIdx)-cellDataI.moleFraction(phaseIdx, compIdx) )/ dist;
865// Scalar gradX = (cellDataJ.massFraction(phaseIdx, compIdx)-cellDataI.massFraction(phaseIdx, compIdx) )/ dist;
866//
867// fluxEntries[compIdx] += gradx * averagedMolarDensity * porousDiffCoeff * faceArea / volume;
868// }
869// }
870
871 return;
872}
873
889template<class TypeTag>
890void FVTransport2P2C<TypeTag>::getFluxOnBoundary(ComponentVector& fluxEntries,
891 Dune::FieldVector<Scalar, 2>& timestepFlux,
892 const Intersection& intersection,
893 const CellData& cellDataI)
894{
895 using std::max;
896 // cell information
897 auto elementI = intersection.inside();
898 int eIdxGlobalI = problem().variables().index(elementI);
899
900 // get position
901 const GlobalPosition globalPos = elementI.geometry().center();
902
903 // cell volume, assume linear map here
904 Scalar volume = elementI.geometry().volume();
905 const GlobalPosition& gravity_ = problem().gravity();
906 // get values of cell I
907 Scalar pressI = problem().pressureModel().pressure(eIdxGlobalI);
908 Scalar pcI = cellDataI.capillaryPressure();
909 DimMatrix K_I(problem().spatialParams().intrinsicPermeability(elementI));
910
911 if(regulateBoundaryPermeability)
912 {
913 int axis = intersection.indexInInside() / 2;
914 if(K_I[axis][axis] < minimalBoundaryPermeability)
915 K_I[axis][axis] = minimalBoundaryPermeability;
916 }
917
918 // old material law interface is deprecated: Replace this by
919 // const auto& fluidMatrixInteraction = problem().spatialParams.fluidMatrixInteractionAtPos(elementI.geometry().center());
920 // after the release of 3.3, when the deprecated interface is no longer supported
921 const auto fluidMatrixInteraction = Deprecated::makePcKrSw(Scalar{}, problem().spatialParams(), elementI);
922
923 Scalar SwmobI = max((cellDataI.saturation(wPhaseIdx)
924 - fluidMatrixInteraction.pcSwCurve().effToAbsParams().swr())
925 , 1e-2);
926 Scalar SnmobI = max((cellDataI.saturation(nPhaseIdx)
927 - fluidMatrixInteraction.pcSwCurve().effToAbsParams().snr())
928 , 1e-2);
929
930 Scalar densityWI (0.), densityNWI(0.);
931 densityWI= cellDataI.density(wPhaseIdx);
932 densityNWI = cellDataI.density(nPhaseIdx);
933
934 // face properties
935 GlobalPosition unitOuterNormal = intersection.centerUnitOuterNormal();
936 if (switchNormals)
937 unitOuterNormal *= -1.0;
938 Scalar faceArea = intersection.geometry().volume();
939
940 // create vector for timestep and for update
941 Dune::FieldVector<Scalar, 2> factor (0.);
942 Dune::FieldVector<Scalar, 2> updFactor (0.);
943
944 PhaseVector potential(0.);
945 // center of face in global coordinates
946 const GlobalPosition& globalPosFace = intersection.geometry().center();
947
948 // distance vector between barycenters
949 GlobalPosition distVec = globalPosFace - globalPos;
950 // compute distance between cell centers
951 Scalar dist = distVec.two_norm();
952
953 GlobalPosition unitDistVec(distVec);
954 unitDistVec /= dist;
955
956 //instantiate a fluid state
957 FluidState BCfluidState;
958
959 //get boundary type
960 BoundaryTypes bcTypes;
961 problem().boundaryTypes(bcTypes, intersection);
962
963 /********** Dirichlet Boundary *************/
964 if (bcTypes.isDirichlet(contiWEqIdx)) // if contiWEq is Dirichlet, so is contiNEq
965 {
966 //get dirichlet pressure boundary condition
967 PhaseVector pressBound(0.);
968 Scalar pcBound (0.);
969
970 // read boundary values
971 this->evalBoundary(globalPosFace,
972 intersection,
973 BCfluidState,
974 pressBound);
975
976 // determine fluid properties at the boundary
977 Scalar densityWBound = BCfluidState.density(wPhaseIdx);
978 Scalar densityNWBound = BCfluidState.density(nPhaseIdx);
979 Scalar viscosityWBound = FluidSystem::viscosity(BCfluidState, wPhaseIdx);
980 Scalar viscosityNWBound = FluidSystem::viscosity(BCfluidState, nPhaseIdx);
981 if(getPropValue<TypeTag, Properties::EnableCapillarity>())
982 pcBound = (BCfluidState.pressure(nPhaseIdx) - BCfluidState.pressure(wPhaseIdx));
983 // average
984 double densityW_mean = (densityWI + densityWBound) / 2;
985 double densityNW_mean = (densityNWI + densityNWBound) / 2;
986
987 // prepare K
988 Dune::FieldVector<Scalar,dim> K(0);
989 K_I.umv(unitDistVec,K);
990
991 //calculate potential gradient
992 switch (pressureType)
993 {
994 case pw:
995 {
996 potential[wPhaseIdx] = (pressI - pressBound[wPhaseIdx]) / dist;
997 potential[nPhaseIdx] = (pressI + pcI - pressBound[wPhaseIdx] - pcBound)
998 / dist;
999 break;
1000 }
1001 case pn:
1002 {
1003 potential[wPhaseIdx] = (pressI - pcI - pressBound[nPhaseIdx] + pcBound)
1004 / dist;
1005 potential[nPhaseIdx] = (pressI - pressBound[nPhaseIdx]) / dist;
1006 break;
1007 }
1008 }
1009 potential[wPhaseIdx] += (gravity_ * unitDistVec) * densityW_mean;
1010 potential[nPhaseIdx] += (gravity_ * unitDistVec) * densityNW_mean;
1011
1012 potential[wPhaseIdx] *= fabs(K * unitOuterNormal);
1013 potential[nPhaseIdx] *= fabs(K * unitOuterNormal);
1014
1015 // old material law interface is deprecated: Replace this by
1016 // const auto& fluidMatrixInteraction = problem().spatialParams.fluidMatrixInteractionAtPos(elementI.geometry().center());
1017 // after the release of 3.3, when the deprecated interface is no longer supported
1018 const auto fluidMatrixInteraction = Deprecated::makePcKrSw(Scalar{}, problem().spatialParams(), elementI);
1019
1020 // do upwinding for lambdas
1021 PhaseVector lambda(0.);
1022 if (potential[wPhaseIdx] >= 0.)
1023 lambda[wPhaseIdx] = cellDataI.mobility(wPhaseIdx);
1024 else
1025 {
1026 if(getPropValue<TypeTag, Properties::BoundaryMobility>()==Indices::satDependent)
1027 lambda[wPhaseIdx] = BCfluidState.saturation(wPhaseIdx) / viscosityWBound;
1028 else
1029 lambda[wPhaseIdx] = fluidMatrixInteraction.krw(BCfluidState.saturation(wPhaseIdx)) / viscosityWBound;
1030 }
1031 if (potential[nPhaseIdx] >= 0.)
1032 lambda[nPhaseIdx] = cellDataI.mobility(nPhaseIdx);
1033 else
1034 {
1035 if(getPropValue<TypeTag, Properties::BoundaryMobility>()==Indices::satDependent)
1036 lambda[nPhaseIdx] = BCfluidState.saturation(nPhaseIdx) / viscosityNWBound;
1037 else
1038 lambda[nPhaseIdx] = fluidMatrixInteraction.krn(BCfluidState.saturation(wPhaseIdx)) / viscosityNWBound;
1039 }
1040 // calculate and standardized velocity
1041
1042 double velocityJIw = max((-lambda[wPhaseIdx] * potential[wPhaseIdx]) * faceArea / volume, 0.0);
1043 double velocityIJw = max(( lambda[wPhaseIdx] * potential[wPhaseIdx]) * faceArea / volume, 0.0);
1044 double velocityJIn = max((-lambda[nPhaseIdx] * potential[nPhaseIdx]) * faceArea / volume, 0.0);
1045 double velocityIJn = max(( lambda[nPhaseIdx] * potential[nPhaseIdx]) * faceArea / volume, 0.0);
1046
1047 // for timestep control
1048 timestepFlux[0] = velocityJIw + velocityJIn;
1049
1050 double foutw = velocityIJw/SwmobI;
1051 double foutn = velocityIJn/SnmobI;
1052 using std::isnan;
1053 using std::isinf;
1054 if (isnan(foutw) || isinf(foutw) || foutw < 0) foutw = 0;
1055 if (isnan(foutn) || isinf(foutn) || foutn < 0) foutn = 0;
1056 timestepFlux[1] = foutw + foutn;
1057
1058 fluxEntries[wCompIdx] =
1059 + velocityJIw * BCfluidState.massFraction(wPhaseIdx, wCompIdx) * densityWBound
1060 - velocityIJw * cellDataI.massFraction(wPhaseIdx, wCompIdx) * densityWI
1061 + velocityJIn * BCfluidState.massFraction(nPhaseIdx, wCompIdx) * densityNWBound
1062 - velocityIJn * cellDataI.massFraction(nPhaseIdx, wCompIdx) * densityNWI ;
1063 fluxEntries[nCompIdx] =
1064 velocityJIw * BCfluidState.massFraction(wPhaseIdx, nCompIdx) * densityWBound
1065 - velocityIJw * cellDataI.massFraction(wPhaseIdx, nCompIdx) * densityWI
1066 + velocityJIn * BCfluidState.massFraction(nPhaseIdx, nCompIdx) * densityNWBound
1067 - velocityIJn * cellDataI.massFraction(nPhaseIdx, nCompIdx) * densityNWI ;
1068 }//end dirichlet boundary
1069 else if (bcTypes.isNeumann(contiWEqIdx))
1070 {
1071 // Convention: outflow => positive sign : has to be subtracted from update vec
1072 PrimaryVariables J(NAN);
1073 problem().neumann(J, intersection);
1074 fluxEntries[wCompIdx] = - J[contiWEqIdx] * faceArea / volume;
1075 fluxEntries[nCompIdx] = - J[contiNEqIdx] * faceArea / volume;
1076
1077 // for timestep control
1078 #define cflIgnoresNeumann
1079 #ifdef cflIgnoresNeumann
1080 timestepFlux[0] = 0;
1081 timestepFlux[1] = 0;
1082 #else
1083 double inflow = updFactor[wCompIdx] / densityW + updFactor[nCompIdx] / densityNW;
1084 if (inflow>0)
1085 {
1086 timestepFlux[0] = updFactor[wCompIdx] / densityW
1087 + updFactor[nCompIdx] / densityNW; // =factor in
1088 timestepFlux[1] = -(updFactor[wCompIdx] / densityW /SwmobI
1089 + updFactor[nCompIdx] / densityNW / SnmobI); // =factor out
1090 }
1091 else
1092 {
1093 timestepFlux[0] = -(updFactor[wCompIdx] / densityW
1094 + updFactor[nCompIdx] / densityNW); // =factor in
1095 timestepFlux[1] = updFactor[wCompIdx] / densityW /SwmobI
1096 + updFactor[nCompIdx] / densityNW / SnmobI; // =factor out
1097 }
1098 #endif
1099 }//end neumann boundary
1100 return;
1101}
1102
1103
1118template<class TypeTag>
1119void FVTransport2P2C<TypeTag>::evalBoundary(GlobalPosition globalPosFace,
1120 const Intersection& intersection,
1121 FluidState &BCfluidState,
1122 PhaseVector &pressBound)
1123{
1124 // prepare a flash solver
1126
1127 auto element = intersection.inside();
1128 // read boundary values
1129 PrimaryVariables primaryVariablesOnBoundary(0.);
1130 problem().dirichlet(primaryVariablesOnBoundary, intersection);
1131
1132 // old material law interface is deprecated: Replace this by
1133 // const auto& fluidMatrixInteraction = problem().spatialParams.fluidMatrixInteractionAtPos(element.geometry().center());
1134 // after the release of 3.3, when the deprecated interface is no longer supported
1135 const auto fluidMatrixInteraction = Deprecated::makePcKrSw(Scalar{}, problem().spatialParams(), element);
1136
1137 // read boundary type
1138 typename Indices::BoundaryFormulation bcType;
1139 problem().boundaryFormulation(bcType, intersection);
1140 if (bcType == Indices::saturation)
1141 {
1142 Scalar satBound = primaryVariablesOnBoundary[contiWEqIdx];
1143
1144 if(getPropValue<TypeTag, Properties::EnableCapillarity>())
1145 {
1146 Scalar pcBound = fluidMatrixInteraction.pc(satBound);
1147 switch (pressureType)
1148 {
1149 case pw:
1150 {
1151 pressBound[wPhaseIdx] = primaryVariablesOnBoundary[Indices::pressureEqIdx];
1152 pressBound[nPhaseIdx] = primaryVariablesOnBoundary[Indices::pressureEqIdx] + pcBound;
1153 break;
1154 }
1155 case pn:
1156 {
1157 pressBound[wPhaseIdx] = primaryVariablesOnBoundary[Indices::pressureEqIdx] - pcBound;
1158 pressBound[nPhaseIdx] = primaryVariablesOnBoundary[Indices::pressureEqIdx];
1159 break;
1160 }
1161 }
1162 }
1163 else // capillarity neglected
1164 pressBound[wPhaseIdx] = pressBound[nPhaseIdx] = primaryVariablesOnBoundary[Indices::pressureEqIdx];
1165
1166 flashSolver.saturationFlash2p2c(BCfluidState, satBound, pressBound,
1167 problem().temperatureAtPos(globalPosFace));
1168 }
1169 else if (bcType == Indices::concentration)
1170 {
1171 // saturation and hence pc and hence corresponding pressure unknown
1172 pressBound[wPhaseIdx] = pressBound[nPhaseIdx] = primaryVariablesOnBoundary[Indices::pressureEqIdx];
1173 Scalar Z0Bound = primaryVariablesOnBoundary[contiWEqIdx];
1174 flashSolver.concentrationFlash2p2c(BCfluidState, Z0Bound, pressBound,
1175 problem().temperatureAtPos(globalPosFace));
1176
1177 if(getPropValue<TypeTag, Properties::EnableCapillarity>())
1178 {
1179 Scalar pcBound = fluidMatrixInteraction.pc(BCfluidState.saturation(wPhaseIdx));
1180 int maxiter = 3;
1181 //start iteration loop
1182 for(int iter=0; iter < maxiter; iter++)
1183 {
1184 //prepare pressures to enter flash calculation
1185 switch (pressureType)
1186 {
1187 case pw:
1188 {
1189 pressBound[wPhaseIdx] = primaryVariablesOnBoundary[Indices::pressureEqIdx];
1190 pressBound[nPhaseIdx] = primaryVariablesOnBoundary[Indices::pressureEqIdx]
1191 + pcBound;
1192 break;
1193 }
1194 case pn:
1195 {
1196 pressBound[wPhaseIdx] = primaryVariablesOnBoundary[Indices::pressureEqIdx]
1197 - pcBound;
1198 pressBound[nPhaseIdx] = primaryVariablesOnBoundary[Indices::pressureEqIdx];
1199 break;
1200 }
1201 }
1202
1203 //store old pc
1204 Scalar oldPc = pcBound;
1205 //update with better pressures
1206 flashSolver.concentrationFlash2p2c(BCfluidState, Z0Bound, pressBound,
1207 problem().temperatureAtPos(globalPosFace));
1208 pcBound = fluidMatrixInteraction.pc(BCfluidState.saturation(wPhaseIdx));
1209 // TODO: get right criterion, do output for evaluation
1210 //converge criterion
1211 using std::abs;
1212 if (abs(oldPc - pcBound) < 10.0)
1213 iter = maxiter;
1214 }
1215 }
1216 }
1217 else
1218 DUNE_THROW(Dune::NotImplemented, "Boundary Formulation neither Concentration nor Saturation??");
1219}
1220
1221template<class TypeTag>
1223{
1224 dt = std::numeric_limits<Scalar>::max();
1225
1226 // update target time-step-sizes
1227 for (const auto& element : elements(problem_.gridView()))
1228 {
1229#if HAVE_MPI
1230 if (element.partitionType() != Dune::InteriorEntity)
1231 {
1232 continue;
1233 }
1234#endif
1235
1236 // cell index
1237 int eIdxGlobalI = problem_.variables().index(element);
1238
1239 LocalTimesteppingData& localDataI = timeStepData_[eIdxGlobalI];
1240
1241
1242 using FaceDt = std::unordered_map<int, Scalar>;
1243 FaceDt faceDt;
1244
1245 // run through all intersections with neighbors and boundary
1246 for (const auto& intersection : intersections(problem_.gridView(), element))
1247 {
1248 int indexInInside = intersection.indexInInside();
1249 using std::min;
1250 if (intersection.neighbor())
1251 {
1252 auto neighbor = intersection.outside();
1253 int eIdxGlobalJ = problem_.variables().index(neighbor);
1254
1255 int levelI = element.level();
1256 int levelJ = neighbor.level();
1257
1258
1259 if (eIdxGlobalI < eIdxGlobalJ && levelI <= levelJ)
1260 {
1261 LocalTimesteppingData& localDataJ = timeStepData_[eIdxGlobalJ];
1262
1263 int indexInOutside = intersection.indexInOutside();
1264
1265 if (localDataI.faceTargetDt[indexInInside] < accumulatedDt_ + dtThreshold_
1266 || localDataJ.faceTargetDt[indexInOutside] < accumulatedDt_ + dtThreshold_)
1267 {
1268 Scalar timeStep = min(localDataI.dt, localDataJ.dt);
1269
1270 if (levelI < levelJ)
1271 {
1272 typename FaceDt::iterator it = faceDt.find(indexInInside);
1273 if (it != faceDt.end())
1274 {
1275 it->second = min(it->second, timeStep);
1276 }
1277 else
1278 {
1279 faceDt.insert(std::make_pair(indexInInside, timeStep));
1280 }
1281 }
1282 else
1283 {
1284 localDataI.faceTargetDt[indexInInside] += subCFLFactor_ * timeStep;
1285 localDataJ.faceTargetDt[indexInOutside] += subCFLFactor_ * timeStep;
1286 }
1287
1288 dt = min(dt, timeStep);
1289 }
1290 }
1291 }
1292 else if (intersection.boundary())
1293 {
1294 if (localDataI.faceTargetDt[indexInInside] < accumulatedDt_ + dtThreshold_)
1295 {
1296 localDataI.faceTargetDt[indexInInside] += subCFLFactor_ * localDataI.dt;
1297 dt = min(dt, subCFLFactor_ * localDataI.dt);
1298 }
1299 }
1300 }
1301 if (faceDt.size() > 0)
1302 {
1303 typename FaceDt::iterator it = faceDt.begin();
1304 for (;it != faceDt.end();++it)
1305 {
1306 localDataI.faceTargetDt[it->first] += subCFLFactor_ * it->second;
1307 }
1308
1309 for (const auto& intersection : intersections(problem_.gridView(), element))
1310 {
1311 if (intersection.neighbor())
1312 {
1313 int indexInInside = intersection.indexInInside();
1314
1315 it = faceDt.find(indexInInside);
1316 if (it != faceDt.end())
1317 {
1318 auto neighbor = intersection.outside();
1319 int eIdxGlobalJ = problem_.variables().index(neighbor);
1320
1321 LocalTimesteppingData& localDataJ = timeStepData_[eIdxGlobalJ];
1322
1323 int indexInOutside = intersection.indexInOutside();
1324
1325 localDataJ.faceTargetDt[indexInOutside] += subCFLFactor_ * it->second;
1326 }
1327 }
1328 }
1329 }
1330 }
1331
1332#if HAVE_MPI
1333 // communicate updated values
1334 using SolutionTypes = GetProp<TypeTag, Properties::SolutionTypes>;
1335 using ElementMapper = typename SolutionTypes::ElementMapper;
1336 using TimeDataHandle = VectorCommDataHandleEqual<ElementMapper, std::vector<LocalTimesteppingData>, 0/*elementCodim*/>;
1337
1338 TimeDataHandle timeDataHandle(problem_.elementMapper(), timeStepData_);
1339 problem_.gridView().template communicate<TimeDataHandle>(timeDataHandle,
1340 Dune::InteriorBorder_All_Interface,
1341 Dune::ForwardCommunication);
1342
1343 dt = problem_.gridView().comm().min(dt);
1344#endif
1345}
1346
1347template<class TypeTag>
1348void FVTransport2P2C<TypeTag>::innerUpdate(TransportSolutionType& updateVec)
1349{
1350 if (localTimeStepping_)
1351 {
1352 Scalar realDt = problem_.timeManager().timeStepSize();
1353
1354 Scalar subDt = realDt;
1355
1356 updatedTargetDt_(subDt);
1357
1358 Scalar accumulatedDtOld = accumulatedDt_;
1359 accumulatedDt_ += subDt;
1360
1361 Scalar t = problem_.timeManager().time();
1362
1363 if (accumulatedDt_ < realDt)
1364 {
1365 using std::min;
1366 while(true)
1367 {
1368 Scalar dtCorrection = min(0.0, realDt - accumulatedDt_);
1369 subDt += dtCorrection;
1370
1371 if (verbosity_ > 0)
1372 std::cout<<" Sub-time-step size: "<<subDt<< std::endl;
1373
1374 bool stopTimeStep = false;
1375 int size = problem_.gridView().size(0);
1376 for (int i = 0; i < size; i++)
1377 {
1378 EntryType newVal(0);
1379 int transportedQuantities = getPropValue<TypeTag, Properties::NumEq>() - 1; // NumEq - 1 pressure Eq
1380 for (int eqNumber = 0; eqNumber < transportedQuantities; eqNumber++)
1381 {
1382 newVal[eqNumber] = totalConcentration_[eqNumber][i];
1383 newVal[eqNumber] += updateVec[eqNumber][i] * subDt;
1384 }
1385 if (!asImp_().inPhysicalRange(newVal))
1386 {
1387 stopTimeStep = true;
1388
1389 break;
1390 }
1391 }
1392
1393#if HAVE_MPI
1394 int rank = 0;
1395 if (stopTimeStep)
1396 rank = problem_.gridView().comm().rank();
1397
1398 rank = problem_.gridView().comm().max(rank);
1399 problem_.gridView().comm().broadcast(&stopTimeStep,1,rank);
1400#endif
1401
1402
1403 if (stopTimeStep && accumulatedDtOld > dtThreshold_)
1404 {
1405 problem_.timeManager().setTimeStepSize(accumulatedDtOld);
1406 break;
1407 }
1408 else
1409 {
1410 asImp_().updateTransportedQuantity(updateVec, subDt);
1411 }
1412
1413
1414 if (accumulatedDt_ >= realDt)
1415 {
1416 break;
1417 }
1418
1419 problem_.pressureModel().updateMaterialLaws();
1420 problem_.model().updateTransport(t, subDt, updateVec);
1421
1422 updatedTargetDt_(subDt);
1423
1424 accumulatedDtOld = accumulatedDt_;
1425 accumulatedDt_ += subDt;
1426 }
1427 }
1428 else
1429 {
1430 asImp_().updateTransportedQuantity(updateVec, realDt);
1431 }
1432
1433 resetTimeStepData_();
1434 }
1435}
1436
1437} // end namespace Dumux
1438#endif
Contains a class to exchange entries of a vector.
Determines the pressures and saturations of all fluid phases given the total mass of all components.
Define some often used mathematical functions.
Helpers for deprecation.
constexpr Scalar harmonicMean(Scalar x, Scalar y, Scalar wx=1.0, Scalar wy=1.0) noexcept
Calculate the (weighted) harmonic mean of two scalar values.
Definition: math.hh:68
void harmonicMeanMatrix(Dune::FieldMatrix< Scalar, m, n > &K, const Dune::FieldMatrix< Scalar, m, n > &Ki, const Dune::FieldMatrix< Scalar, m, n > &Kj)
Calculate the harmonic mean of a fixed-size matrix.
Definition: math.hh:108
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
Flash calculation routines for compositional sequential models.
Definition: compositionalflash.hh:51
static void saturationFlash2p2c(FluidState &fluidState, const Scalar &saturation, const PhaseVector &phasePressure, const Scalar &porosity, const Scalar &temperature)
Definition: compositionalflash.hh:219
static void concentrationFlash2p2c(FluidState &fluidState, const Scalar &Z0, const PhaseVector &phasePressure, const Scalar &porosity, const Scalar &temperature)
Definition: compositionalflash.hh:71
A data handle class to exchange entries of a vector.
Definition: vectorcommdatahandle.hh:78
Compositional transport step in a Finite Volume discretization.
Definition: fvtransport.hh:62
virtual void update(const Scalar t, Scalar &dt, TransportSolutionType &updateVec, bool impet=false)
Calculate the update vector and determine timestep size This method calculates the update vector of ...
Definition: fvtransport.hh:364
void evalBoundary(GlobalPosition globalPosFace, const Intersection &intersection, FluidState &BCfluidState, PhaseVector &pressBound)
Evaluate the boundary conditions.
Definition: fvtransport.hh:1119
void deserializeEntity(std::istream &instream, const Element &element)
Function needed for restart option of the transport model: Read in.
Definition: fvtransport.hh:195
bool regulateBoundaryPermeability
Enables regulation of permeability in the direction of a Dirichlet Boundary Condition.
Definition: fvtransport.hh:311
void updateConcentrations(TransportSolutionType &updateVector, Scalar dt)
Updates the concentrations once an update is calculated.
Definition: fvtransport.hh:549
bool localTimeStepping_
Definition: fvtransport.hh:331
Scalar subCFLFactor_
Definition: fvtransport.hh:330
int restrictFluxInTransport_
Restriction of flux on new pressure field if direction reverses from the pressure equation.
Definition: fvtransport.hh:309
FVTransport2P2C(Problem &problem)
Constructs a FVTransport2P2C object Currently, the compositional transport scheme can not be applied ...
Definition: fvtransport.hh:274
const Scalar dtThreshold_
Threshold for sub-time-stepping routine.
Definition: fvtransport.hh:318
void initialize()
Set the initial values before the first pressure equation This method is called before first pressure...
Definition: fvtransport.hh:155
void addOutputVtkFields(MultiWriter &writer)
Write transport variables into the output files.
Definition: fvtransport.hh:172
void getSource(Scalar &update, const Element &element, CellData &cellDataI)
Definition: fvtransport.hh:237
Scalar & totalConcentration(int compIdx, int eIdxGlobal)
Return the the total concentration stored in the transport vector To get real cell values,...
Definition: fvtransport.hh:232
Problem & problem_
Definition: fvtransport.hh:302
void innerUpdate(TransportSolutionType &updateVec)
Definition: fvtransport.hh:1348
int averagedFaces_
number of faces were flux was restricted
Definition: fvtransport.hh:304
void updatedTargetDt_(Scalar &dt)
Definition: fvtransport.hh:1222
void resetTimeStepData_()
Definition: fvtransport.hh:324
Dune::FieldVector< Scalar, 2 > EntryType
Definition: fvtransport.hh:109
bool impet_
indicating if we are in an estimate (false) or real impet (true) step.
Definition: fvtransport.hh:303
void serializeEntity(std::ostream &outstream, const Element &element)
Function needed for restart option of the transport model: Write out.
Definition: fvtransport.hh:188
void getTransportedQuantity(TransportSolutionType &transportedQuantity)
Return the vector of the transported quantity For an immiscible IMPES scheme, this is the saturation....
Definition: fvtransport.hh:213
void getFlux(ComponentVector &fluxEntries, EntryType &timestepFlux, const Intersection &intersection, CellData &cellDataI)
Get flux at an interface between two cells The flux through is calculated according to the underlyin...
Definition: fvtransport.hh:580
Problem & problem()
Acess function for the current problem.
Definition: fvtransport.hh:123
std::vector< LocalTimesteppingData > timeStepData_
Stores data for sub-time-stepping.
Definition: fvtransport.hh:320
int verbosity_
Definition: fvtransport.hh:332
bool inPhysicalRange(DataEntry &entry)
Function to control the abort of the transport-sub-time-stepping depending on a physical parameter ra...
Definition: fvtransport.hh:247
TransportSolutionType totalConcentration_
private vector of transported primary variables
Definition: fvtransport.hh:301
Scalar accumulatedDt_
Current time-interval in sub-time-stepping routine.
Definition: fvtransport.hh:316
virtual ~FVTransport2P2C()
Definition: fvtransport.hh:296
void updateTransportedQuantity(TransportSolutionType &updateVector)
Updates the transported quantity once an update is calculated.
Definition: fvtransport.hh:520
void getFluxOnBoundary(ComponentVector &fluxEntries, EntryType &timestepFlux, const Intersection &intersection, const CellData &cellDataI)
Get flux on Boundary.
Definition: fvtransport.hh:890
static const int pressureType
gives kind of pressure used ( , , )
Definition: fvtransport.hh:307
bool switchNormals
Definition: fvtransport.hh:314
bool enableLocalTimeStepping()
Function to check if local time stepping is activated.
Definition: fvtransport.hh:261
Dune::FieldVector< Scalar, 2 > TimeStepFluxType
Definition: fvtransport.hh:110
Scalar minimalBoundaryPermeability
Minimal limit for the boundary permeability.
Definition: fvtransport.hh:313
Definition: fvtransport.hh:114
Dune::FieldVector< EntryType, 2 *dim > faceFluxes
Definition: fvtransport.hh:115
Dune::FieldVector< Scalar, 2 *dim > faceTargetDt
Definition: fvtransport.hh:116
LocalTimesteppingData()
Definition: fvtransport.hh:118
Scalar dt
Definition: fvtransport.hh:117
Defines the properties required for the sequential 2p2c models.