24#ifndef DUMUX_FV2DTRANSPORT2P2C_ADAPTIVE_HH
25#define DUMUX_FV2DTRANSPORT2P2C_ADAPTIVE_HH
27#include <dune/grid/common/gridenums.hh>
28#include <dune/common/float_cmp.hh>
57template<
class TypeTag>
72 dim = GridView::dimension, dimWorld = GridView::dimensionworld,
73 NumPhases = getPropValue<TypeTag, Properties::NumPhases>()
77 pw = Indices::pressureW,
78 pn = Indices::pressureN
82 wPhaseIdx = Indices::wPhaseIdx, nPhaseIdx = Indices::nPhaseIdx,
83 wCompIdx = Indices::wPhaseIdx, nCompIdx = Indices::nPhaseIdx,
84 contiWEqIdx=Indices::contiWEqIdx, contiNEqIdx=Indices::contiNEqIdx
87 using IntersectionIterator =
typename GridView::IntersectionIterator;
89 using GlobalPosition = Dune::FieldVector<Scalar, dimWorld>;
90 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
91 using TransmissivityMatrix = Dune::FieldVector<Scalar,dim+1>;
92 using PhaseVector = Dune::FieldVector<Scalar, NumPhases>;
102 virtual void update(
const Scalar t, Scalar& dt, TransportSolutionType& updateVec,
bool impes =
false);
104 void getMpfaFlux(Dune::FieldVector<Scalar, 2>&, Dune::FieldVector<Scalar, 2>&,
105 const IntersectionIterator&, CellData&);
119 enableMPFA = getParam<bool>(
"GridAdapt.EnableMultiPointFluxApproximation");
131 static const int pressureType = getPropValue<TypeTag, Properties::PressureFormulation>();
153template<
class TypeTag>
156 this->impet_ = impet;
159 this->averagedFaces_ = 0;
162 unsigned int size_ = problem_.gridView().size(0);
163 updateVec.resize(getPropValue<TypeTag, Properties::NumComponents>());
164 updateVec[wCompIdx].resize(size_);
165 updateVec[nCompIdx].resize(size_);
166 updateVec[wCompIdx] = 0;
167 updateVec[nCompIdx] = 0;
169 if (this->enableLocalTimeStepping())
171 if (this->timeStepData_.size() != size_)
172 this->timeStepData_.resize(size_);
176 if(this->totalConcentration_.size() != size_)
178 this->totalConcentration_[wCompIdx].resize(size_);
179 this->totalConcentration_[nCompIdx].resize(size_);
182 for (
int i = 0; i< problem().gridView().size(0); i++)
184 CellData& cellDataI = problem().variables().cellData(i);
185 for(
int compIdx = 0; compIdx < getPropValue<TypeTag, Properties::NumComponents>(); compIdx++)
187 this->totalConcentration_[compIdx][i]
188 = cellDataI.totalConcentration(compIdx);
194 int restrictingCell = -1;
196 Dune::FieldVector<Scalar, 2> entries(0.), timestepFlux(0.);
198 for (
const auto& element : elements(problem().gridView()))
201 int globalIdxI = problem().variables().index(element);
202 CellData& cellDataI = problem().variables().cellData(globalIdxI);
205 double sumfactorin = 0;
206 double sumfactorout = 0;
208 if (this->enableLocalTimeStepping())
210 LocalTimesteppingData& localData = this->timeStepData_[globalIdxI];
211 for (
int i = 0; i < 2*dim; i++)
213 if (localData.faceTargetDt[i] < this->accumulatedDt_ + this->dtThreshold_)
215 localData.faceFluxes[i] = 0.0;
221 const auto isEndIt = problem().gridView().iend(element);
222 for (
auto isIt = problem().gridView().ibegin(element); isIt != isEndIt; ++isIt)
224 const auto& intersection = *isIt;
226 int indexInInside = intersection.indexInInside();
229 if (intersection.neighbor())
231 if (enableMPFA && (intersection.outside().level() != element.level()))
232 getMpfaFlux(entries, timestepFlux, isIt, cellDataI);
234 this->getFlux(entries, timestepFlux, intersection, cellDataI);
238 if (intersection.boundary())
240 this->getFluxOnBoundary(entries, timestepFlux, intersection, cellDataI);
243 if (this->enableLocalTimeStepping())
245 LocalTimesteppingData& localData = this->timeStepData_[globalIdxI];
247 if (localData.faceTargetDt[indexInInside] < this->accumulatedDt_ + this->dtThreshold_)
249 localData.faceFluxes[indexInInside] += entries;
255 updateVec[wCompIdx][globalIdxI] += entries[wCompIdx];
256 updateVec[nCompIdx][globalIdxI] += entries[nCompIdx];
260 sumfactorin += timestepFlux[0];
261 sumfactorout += timestepFlux[1];
265 if (this->enableLocalTimeStepping())
267 LocalTimesteppingData& localData = this->timeStepData_[globalIdxI];
268 for (
int i=0; i < 2*dim; i++)
270 updateVec[wCompIdx][globalIdxI] += localData.faceFluxes[i][wCompIdx];
271 updateVec[nCompIdx][globalIdxI] += localData.faceFluxes[i][nCompIdx];
276 PrimaryVariables q(NAN);
277 problem().source(q, element);
278 updateVec[wCompIdx][globalIdxI] += q[contiWEqIdx];
279 updateVec[nCompIdx][globalIdxI] += q[contiNEqIdx];
283 sumfactorin = max(sumfactorin,sumfactorout)
284 / problem().spatialParams().porosity(element);
287 if (this->enableLocalTimeStepping())
289 this->timeStepData_[globalIdxI].dt = 1./sumfactorin;
290 if ( 1./sumfactorin < dt)
293 restrictingCell= globalIdxI;
298 if ( 1./sumfactorin < dt)
301 restrictingCell= globalIdxI;
309 using ElementMapper =
typename SolutionTypes::ElementMapper;
311 for (
int i = 0; i < updateVec.size(); i++)
313 DataHandle dataHandle(problem_.variables().elementMapper(), updateVec[i]);
314 problem_.gridView().template communicate<DataHandle>(dataHandle,
315 Dune::InteriorBorder_All_Interface,
316 Dune::ForwardCommunication);
318 dt = problem_.gridView().comm().min(dt);
323 Dune::dinfo <<
"Timestep restricted by CellIdx " << restrictingCell <<
" leads to dt = "
324 <<dt * getParam<Scalar>(
"Impet.CFLFactor")<< std::endl;
325 if(this->averagedFaces_ != 0)
326 Dune::dinfo <<
" Averageing done for " << this->averagedFaces_ <<
" faces. "<< std::endl;
351template<
class TypeTag>
353 const IntersectionIterator& intersectionIterator, CellData& cellDataI)
358 auto elementI = intersectionIterator->inside();
359 int globalIdxI = problem().variables().index(elementI);
362 const GlobalPosition globalPos = elementI.geometry().center();
363 const GlobalPosition& gravity_ = problem().gravity();
365 Scalar volume = elementI.geometry().volume();
368 Scalar pressI = problem().pressureModel().pressure(globalIdxI);
369 Scalar pcI = cellDataI.capillaryPressure();
374 const auto fluidMatrixInteraction = Deprecated::makePcKrSw(Scalar{}, problem().spatialParams(), elementI);
376 PhaseVector SmobI(0.);
378 SmobI[wPhaseIdx] = max((cellDataI.saturation(wPhaseIdx)
379 - fluidMatrixInteraction.pcSwCurve().effToAbsParams().swr())
381 SmobI[nPhaseIdx] = max((cellDataI.saturation(nPhaseIdx)
382 - fluidMatrixInteraction.pcSwCurve().effToAbsParams().snr())
385 Scalar densityWI (0.), densityNWI(0.);
386 densityWI= cellDataI.density(wPhaseIdx);
387 densityNWI = cellDataI.density(nPhaseIdx);
389 PhaseVector potential(0.);
392 auto neighbor = intersectionIterator->outside();
393 int globalIdxJ = problem().variables().index(neighbor);
394 CellData& cellDataJ = problem().variables().cellData(globalIdxJ);
397 const GlobalPosition& globalPosNeighbor = neighbor.geometry().center();
400 GlobalPosition distVec = globalPosNeighbor - globalPos;
402 Scalar dist = distVec.two_norm();
404 GlobalPosition unitDistVec(distVec);
408 Scalar densityWJ (0.), densityNWJ(0.);
409 densityWJ = cellDataJ.density(wPhaseIdx);
410 densityNWJ = cellDataJ.density(nPhaseIdx);
413 double densityW_mean = (densityWI + densityWJ) * 0.5;
414 double densityNW_mean = (densityNWI + densityNWJ) * 0.5;
416 double pressJ = problem().pressureModel().pressure(globalIdxJ);
417 Scalar pcJ = cellDataJ.capillaryPressure();
421 GlobalPosition globalPos3(0.);
423 TransmissivityMatrix T(0.);
424 IntersectionIterator additionalIsIt = intersectionIterator;
425 TransmissivityMatrix additionalT(0.);
428 = problem().variables().getMpfaData(*intersectionIterator, additionalIsIt, T, additionalT,
429 globalPos3, globalIdx3);
430 if (halfedgesStored == 0)
431 halfedgesStored = problem().pressureModel().computeTransmissibilities(intersectionIterator,additionalIsIt,
432 T,additionalT, globalPos3, globalIdx3 );
435 Scalar press3 = problem().pressureModel().pressure(globalIdx3);
436 CellData& cellData3 = problem().variables().cellData(globalIdx3);
437 Scalar pc3 = cellData3.capillaryPressure();
438 Scalar temp1 = globalPos * gravity_;
439 Scalar temp2 = globalPosNeighbor * gravity_;
440 Scalar temp3 = globalPos3 * gravity_;
443 potential[wPhaseIdx] += (pressI-temp1*densityW_mean) * T[2]
444 +(pressJ-temp2*densityW_mean) * T[0]
445 +(press3- temp3*densityW_mean) * T[1];
446 potential[nPhaseIdx] += (pressI+pcI-temp1*densityNW_mean) * T[2]
447 +(pressJ+pcJ-temp2*densityNW_mean) * T[0]
448 +(press3+pc3- temp3*densityNW_mean) * T[1];
450 if(halfedgesStored == 2)
452 int AdditionalIdx = problem().variables().index(additionalIsIt->outside());
453 CellData& cellDataAdditional = problem().variables().cellData(AdditionalIdx);
454 potential[wPhaseIdx] += (pressI-temp1*densityW_mean) * additionalT[2]
455 +(pressJ-temp2*densityW_mean) * additionalT[0]
456 +(problem().pressureModel().pressure(AdditionalIdx)
457 -(additionalIsIt->outside().geometry().center()*gravity_*densityW_mean)
459 potential[nPhaseIdx] += (pressI+pcI-temp1*densityNW_mean) * additionalT[2]
460 +(pressJ+pcJ-temp2*densityNW_mean) * additionalT[0]
461 +(problem().pressureModel().pressure(AdditionalIdx)
462 + cellDataAdditional.capillaryPressure()
463 -(additionalIsIt->outside().geometry().center()*gravity_*densityNW_mean)
467 else if(pressureType==pn)
469 potential[wPhaseIdx] += (pressI-pcI-temp1*densityW_mean) * T[2]
470 + (pressJ-pcJ-temp2*densityW_mean) * T[0]
471 + (press3-pc3- temp3*densityW_mean) * T[1];
472 potential[nPhaseIdx] += (pressI-temp1*densityNW_mean) * T[2]
473 + (pressJ-temp2*densityNW_mean) * T[0]
474 + (press3-temp3*densityNW_mean) * T[1];
476 if(halfedgesStored == 2)
478 int AdditionalIdx = problem().variables().index(additionalIsIt->outside());
479 CellData& cellDataAdditional = problem().variables().cellData(AdditionalIdx);
481 potential[wPhaseIdx] += (pressI-pcI-temp1*densityW_mean) * additionalT[2]
482 +(pressJ-pcJ-temp2*densityW_mean) * additionalT[0]
483 +(problem().pressureModel().pressure(AdditionalIdx)
484 - cellDataAdditional.capillaryPressure()
485 -(additionalIsIt->outside().geometry().center()*gravity_*densityW_mean)
487 potential[nPhaseIdx] += (pressI-temp1*densityNW_mean) * additionalT[2]
488 +(pressJ-temp2*densityNW_mean) * additionalT[0]
489 +(problem().pressureModel().pressure(AdditionalIdx)
490 -(additionalIsIt->outside().geometry().center()*gravity_*densityNW_mean))
497 Dune::FieldVector<bool, NumPhases> doUpwinding(
true);
498 PhaseVector lambda(0.);
499 for(
int phaseIdx = 0; phaseIdx < NumPhases; phaseIdx++)
502 if(phaseIdx == wPhaseIdx)
503 contiEqIdx = contiWEqIdx;
505 contiEqIdx = contiNEqIdx;
507 if(!this->impet_ or !this->restrictFluxInTransport_)
509 if(potential[phaseIdx] > 0.)
511 lambda[phaseIdx] = cellDataI.mobility(phaseIdx);
512 if(elementI.level()>neighbor.level())
513 cellDataI.setUpwindCell(intersectionIterator->indexInInside(), contiEqIdx,
true);
515 else if(potential[phaseIdx] < 0.)
517 lambda[phaseIdx] = cellDataJ.mobility(phaseIdx);
518 if(elementI.level()>neighbor.level())
519 cellDataI.setUpwindCell(intersectionIterator->indexInInside(), contiEqIdx,
false);
523 doUpwinding[phaseIdx] =
false;
524 if(elementI.level()>neighbor.level())
525 cellDataI.setUpwindCell(intersectionIterator->indexInInside(), contiEqIdx,
false);
527 cellDataJ.setUpwindCell(intersectionIterator->indexInOutside(), contiEqIdx,
false);
532 bool cellIwasUpwindCell;
534 if(elementI.level()>neighbor.level())
535 cellIwasUpwindCell = cellDataI.isUpwindCell(intersectionIterator->indexInInside(), contiEqIdx);
537 cellIwasUpwindCell = !cellDataJ.isUpwindCell(intersectionIterator->indexInOutside(), contiEqIdx);
539 if (potential[phaseIdx] > 0. && cellIwasUpwindCell)
540 lambda[phaseIdx] = cellDataI.mobility(phaseIdx);
541 else if (potential[phaseIdx] < 0. && !cellIwasUpwindCell)
542 lambda[phaseIdx] = cellDataJ.mobility(phaseIdx);
544 else if(this->restrictFluxInTransport_ == 2)
545 doUpwinding[phaseIdx] =
false;
550 if (potential[phaseIdx] > 0. && Dune::FloatCmp::ne<Scalar, Dune::FloatCmp::absolute>(cellDataJ.mobility(phaseIdx), 0.0, 1.0e-30))
551 lambda[phaseIdx] = cellDataI.mobility(phaseIdx);
553 else if (potential[phaseIdx] < 0. && Dune::FloatCmp::ne<Scalar, Dune::FloatCmp::absolute>(cellDataI.mobility(phaseIdx), 0.0, 1.0e-30))
554 lambda[phaseIdx] = cellDataJ.mobility(phaseIdx);
556 doUpwinding[phaseIdx] =
false;
560 if(!doUpwinding[phaseIdx])
563 if(cellDataI.mobility(phaseIdx)+cellDataJ.mobility(phaseIdx)==0.)
565 potential[phaseIdx] = 0;
570 fluxEntries[wCompIdx] -= potential[phaseIdx] / volume
571 *
harmonicMean(cellDataI.massFraction(phaseIdx, wCompIdx) * cellDataI.mobility(phaseIdx) * cellDataI.density(phaseIdx),
572 cellDataJ.massFraction(phaseIdx, wCompIdx) * cellDataJ.mobility(phaseIdx) * cellDataJ.density(phaseIdx));
573 fluxEntries[nCompIdx] -= potential[phaseIdx] / volume
574 *
harmonicMean(cellDataI.massFraction(phaseIdx, nCompIdx) * cellDataI.mobility(phaseIdx) * cellDataI.density(phaseIdx),
575 cellDataJ.massFraction(phaseIdx, nCompIdx) * cellDataJ.mobility(phaseIdx) * cellDataJ.density(phaseIdx));
579 timestepFlux[0] += max(0.,
580 - potential[phaseIdx] / volume
581 *
harmonicMean(cellDataI.mobility(phaseIdx),cellDataJ.mobility(phaseIdx)));
583 timestepFlux[1] += max(0.,
584 potential[phaseIdx] / volume
585 *
harmonicMean(cellDataI.mobility(phaseIdx),cellDataJ.mobility(phaseIdx))/SmobI[phaseIdx]);
588 this->averagedFaces_++;
589 #if DUNE_MINIMAL_DEBUG_LEVEL < 3
591 if(globalIdxI > globalIdxJ)
592 Dune::dinfo <<
"harmonicMean flux of phase" << phaseIdx <<
" used from cell" << globalIdxI<<
" into " << globalIdxJ
593 <<
" ; TE upwind I = "<< cellDataI.isUpwindCell(intersectionIterator->indexInInside(), contiEqIdx)
594 <<
" but pot = "<< potential[phaseIdx] <<
" \n";
598 potential[phaseIdx] = 0;
605 double velocityJIw = max((-lambda[wPhaseIdx] * potential[wPhaseIdx]) / volume, 0.0);
606 double velocityIJw = max(( lambda[wPhaseIdx] * potential[wPhaseIdx]) / volume, 0.0);
607 double velocityJIn = max((-lambda[nPhaseIdx] * potential[nPhaseIdx]) / volume, 0.0);
608 double velocityIJn = max(( lambda[nPhaseIdx] * potential[nPhaseIdx]) / volume, 0.0);
611 timestepFlux[0] += velocityJIw + velocityJIn;
613 double foutw = velocityIJw/SmobI[wPhaseIdx];
614 double foutn = velocityIJn/SmobI[nPhaseIdx];
617 if (isnan(foutw) || isinf(foutw) || foutw < 0) foutw = 0;
618 if (isnan(foutn) || isinf(foutn) || foutn < 0) foutn = 0;
619 timestepFlux[1] += foutw + foutn;
621 fluxEntries[wCompIdx] +=
622 velocityJIw * cellDataJ.massFraction(wPhaseIdx, wCompIdx) * densityWJ
623 - velocityIJw * cellDataI.massFraction(wPhaseIdx, wCompIdx) * densityWI
624 + velocityJIn * cellDataJ.massFraction(nPhaseIdx, wCompIdx) * densityNWJ
625 - velocityIJn * cellDataI.massFraction(nPhaseIdx, wCompIdx) * densityNWI;
626 fluxEntries[nCompIdx] +=
627 velocityJIw * cellDataJ.massFraction(wPhaseIdx, nCompIdx) * densityWJ
628 - velocityIJw * cellDataI.massFraction(wPhaseIdx, nCompIdx) * densityWI
629 + velocityJIn * cellDataJ.massFraction(nPhaseIdx, nCompIdx) * densityNWJ
630 - velocityIJn * cellDataI.massFraction(nPhaseIdx, nCompIdx) * densityNWI;
Contains a class to exchange entries of a vector.
Defines the properties required for the adaptive sequential 2p2c models.
Finite volume discretization of the component transport equation.
Define some often used mathematical functions.
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
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
A data handle class to exchange entries of a vector.
Definition: vectorcommdatahandle.hh:78
Compositional Transport step in a Finite Volume discretization for a adaptive 2D-grid.
Definition: fv2dtransportadaptive.hh:59
virtual ~FV2dTransport2P2CAdaptive()
Definition: fv2dtransportadaptive.hh:122
FV2dTransport2P2CAdaptive(Problem &problem)
Constructs a FV2dTransport2P2CAdaptive object.
Definition: fv2dtransportadaptive.hh:116
static const int pressureType
gives kind of pressure used ( , , )
Definition: fv2dtransportadaptive.hh:131
Problem & problem_
Definition: fv2dtransportadaptive.hh:127
virtual void update(const Scalar t, Scalar &dt, TransportSolutionType &updateVec, bool impes=false)
Calculate the update vector and determine timestep size.
Definition: fv2dtransportadaptive.hh:154
bool enableMPFA
Definition: fv2dtransportadaptive.hh:129
void getMpfaFlux(Dune::FieldVector< Scalar, 2 > &, Dune::FieldVector< Scalar, 2 > &, const IntersectionIterator &, CellData &)
Compute flux over an irregular interface using a mpfa method.
Definition: fv2dtransportadaptive.hh:352
Compositional transport step in a Finite Volume discretization.
Definition: fvtransport.hh:62
Definition: fvtransport.hh:114