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