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>
55template<
class TypeTag>
70 dim = GridView::dimension, dimWorld = GridView::dimensionworld,
71 NumPhases = getPropValue<TypeTag, Properties::NumPhases>()
75 pw = Indices::pressureW,
76 pn = Indices::pressureN
80 wPhaseIdx = Indices::wPhaseIdx, nPhaseIdx = Indices::nPhaseIdx,
81 wCompIdx = Indices::wPhaseIdx, nCompIdx = Indices::nPhaseIdx,
82 contiWEqIdx=Indices::contiWEqIdx, contiNEqIdx=Indices::contiNEqIdx
85 using IntersectionIterator =
typename GridView::IntersectionIterator;
87 using GlobalPosition = Dune::FieldVector<Scalar, dimWorld>;
88 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
89 using TransmissivityMatrix = Dune::FieldVector<Scalar,dim+1>;
90 using PhaseVector = Dune::FieldVector<Scalar, NumPhases>;
100 virtual void update(
const Scalar t, Scalar& dt, TransportSolutionType& updateVec,
bool impes =
false);
102 void getMpfaFlux(Dune::FieldVector<Scalar, 2>&, Dune::FieldVector<Scalar, 2>&,
103 const IntersectionIterator&, CellData&);
117 enableMPFA = getParam<bool>(
"GridAdapt.EnableMultiPointFluxApproximation");
129 static const int pressureType = getPropValue<TypeTag, Properties::PressureFormulation>();
151template<
class TypeTag>
154 this->impet_ = impet;
157 this->averagedFaces_ = 0;
160 unsigned int size_ = problem_.gridView().size(0);
161 updateVec.resize(getPropValue<TypeTag, Properties::NumComponents>());
162 updateVec[wCompIdx].resize(size_);
163 updateVec[nCompIdx].resize(size_);
164 updateVec[wCompIdx] = 0;
165 updateVec[nCompIdx] = 0;
167 if (this->enableLocalTimeStepping())
169 if (this->timeStepData_.size() != size_)
170 this->timeStepData_.resize(size_);
174 if(this->totalConcentration_.size() != size_)
176 this->totalConcentration_[wCompIdx].resize(size_);
177 this->totalConcentration_[nCompIdx].resize(size_);
180 for (
int i = 0; i< problem().gridView().size(0); i++)
182 CellData& cellDataI = problem().variables().cellData(i);
183 for(
int compIdx = 0; compIdx < getPropValue<TypeTag, Properties::NumComponents>(); compIdx++)
185 this->totalConcentration_[compIdx][i]
186 = cellDataI.totalConcentration(compIdx);
192 int restrictingCell = -1;
194 Dune::FieldVector<Scalar, 2> entries(0.), timestepFlux(0.);
196 for (
const auto& element : elements(problem().gridView()))
199 int globalIdxI = problem().variables().index(element);
200 CellData& cellDataI = problem().variables().cellData(globalIdxI);
203 double sumfactorin = 0;
204 double sumfactorout = 0;
206 if (this->enableLocalTimeStepping())
208 LocalTimesteppingData& localData = this->timeStepData_[globalIdxI];
209 for (
int i = 0; i < 2*dim; i++)
211 if (localData.faceTargetDt[i] < this->accumulatedDt_ + this->dtThreshold_)
213 localData.faceFluxes[i] = 0.0;
219 const auto isEndIt = problem().gridView().iend(element);
220 for (
auto isIt = problem().gridView().ibegin(element); isIt != isEndIt; ++isIt)
222 const auto& intersection = *isIt;
224 int indexInInside = intersection.indexInInside();
227 if (intersection.neighbor())
229 if (enableMPFA && (intersection.outside().level() != element.level()))
230 getMpfaFlux(entries, timestepFlux, isIt, cellDataI);
232 this->getFlux(entries, timestepFlux, intersection, cellDataI);
236 if (intersection.boundary())
238 this->getFluxOnBoundary(entries, timestepFlux, intersection, cellDataI);
241 if (this->enableLocalTimeStepping())
243 LocalTimesteppingData& localData = this->timeStepData_[globalIdxI];
245 if (localData.faceTargetDt[indexInInside] < this->accumulatedDt_ + this->dtThreshold_)
247 localData.faceFluxes[indexInInside] += entries;
253 updateVec[wCompIdx][globalIdxI] += entries[wCompIdx];
254 updateVec[nCompIdx][globalIdxI] += entries[nCompIdx];
258 sumfactorin += timestepFlux[0];
259 sumfactorout += timestepFlux[1];
263 if (this->enableLocalTimeStepping())
265 LocalTimesteppingData& localData = this->timeStepData_[globalIdxI];
266 for (
int i=0; i < 2*dim; i++)
268 updateVec[wCompIdx][globalIdxI] += localData.faceFluxes[i][wCompIdx];
269 updateVec[nCompIdx][globalIdxI] += localData.faceFluxes[i][nCompIdx];
274 PrimaryVariables q(NAN);
275 problem().source(q, element);
276 updateVec[wCompIdx][globalIdxI] += q[contiWEqIdx];
277 updateVec[nCompIdx][globalIdxI] += q[contiNEqIdx];
281 sumfactorin = max(sumfactorin,sumfactorout)
282 / problem().spatialParams().porosity(element);
285 if (this->enableLocalTimeStepping())
287 this->timeStepData_[globalIdxI].dt = 1./sumfactorin;
288 if ( 1./sumfactorin < dt)
291 restrictingCell= globalIdxI;
296 if ( 1./sumfactorin < dt)
299 restrictingCell= globalIdxI;
307 using ElementMapper =
typename SolutionTypes::ElementMapper;
309 for (
int i = 0; i < updateVec.size(); i++)
311 DataHandle dataHandle(problem_.variables().elementMapper(), updateVec[i]);
312 problem_.gridView().template communicate<DataHandle>(dataHandle,
313 Dune::InteriorBorder_All_Interface,
314 Dune::ForwardCommunication);
316 dt = problem_.gridView().comm().min(dt);
321 Dune::dinfo <<
"Timestep restricted by CellIdx " << restrictingCell <<
" leads to dt = "
322 <<dt * getParam<Scalar>(
"Impet.CFLFactor")<< std::endl;
323 if(this->averagedFaces_ != 0)
324 Dune::dinfo <<
" Averageing done for " << this->averagedFaces_ <<
" faces. "<< std::endl;
349template<
class TypeTag>
351 const IntersectionIterator& intersectionIterator, CellData& cellDataI)
356 auto elementI = intersectionIterator->inside();
357 int globalIdxI = problem().variables().index(elementI);
360 const GlobalPosition globalPos = elementI.geometry().center();
361 const GlobalPosition& gravity_ = problem().gravity();
363 Scalar
volume = elementI.geometry().volume();
366 Scalar pressI = problem().pressureModel().pressure(globalIdxI);
367 Scalar pcI = cellDataI.capillaryPressure();
369 const auto fluidMatrixInteraction = problem().spatialParams().fluidMatrixInteractionAtPos(elementI.geometry().center());
371 PhaseVector SmobI(0.);
373 SmobI[wPhaseIdx] = max((cellDataI.saturation(wPhaseIdx)
374 - fluidMatrixInteraction.pcSwCurve().effToAbsParams().swr())
376 SmobI[nPhaseIdx] = max((cellDataI.saturation(nPhaseIdx)
377 - fluidMatrixInteraction.pcSwCurve().effToAbsParams().snr())
380 Scalar densityWI (0.), densityNWI(0.);
381 densityWI= cellDataI.density(wPhaseIdx);
382 densityNWI = cellDataI.density(nPhaseIdx);
384 PhaseVector potential(0.);
387 auto neighbor = intersectionIterator->outside();
388 int globalIdxJ = problem().variables().index(neighbor);
389 CellData& cellDataJ = problem().variables().cellData(globalIdxJ);
392 const GlobalPosition& globalPosNeighbor = neighbor.geometry().center();
395 GlobalPosition distVec = globalPosNeighbor - globalPos;
397 Scalar dist = distVec.two_norm();
399 GlobalPosition unitDistVec(distVec);
403 Scalar densityWJ (0.), densityNWJ(0.);
404 densityWJ = cellDataJ.density(wPhaseIdx);
405 densityNWJ = cellDataJ.density(nPhaseIdx);
408 double densityW_mean = (densityWI + densityWJ) * 0.5;
409 double densityNW_mean = (densityNWI + densityNWJ) * 0.5;
411 double pressJ = problem().pressureModel().pressure(globalIdxJ);
412 Scalar pcJ = cellDataJ.capillaryPressure();
416 GlobalPosition globalPos3(0.);
418 TransmissivityMatrix T(0.);
419 IntersectionIterator additionalIsIt = intersectionIterator;
420 TransmissivityMatrix additionalT(0.);
423 = problem().variables().getMpfaData(*intersectionIterator, additionalIsIt, T, additionalT,
424 globalPos3, globalIdx3);
425 if (halfedgesStored == 0)
426 halfedgesStored = problem().pressureModel().computeTransmissibilities(intersectionIterator,additionalIsIt,
427 T,additionalT, globalPos3, globalIdx3 );
430 Scalar press3 = problem().pressureModel().pressure(globalIdx3);
431 CellData& cellData3 = problem().variables().cellData(globalIdx3);
432 Scalar pc3 = cellData3.capillaryPressure();
433 Scalar temp1 = globalPos * gravity_;
434 Scalar temp2 = globalPosNeighbor * gravity_;
435 Scalar temp3 = globalPos3 * gravity_;
438 potential[wPhaseIdx] += (pressI-temp1*densityW_mean) * T[2]
439 +(pressJ-temp2*densityW_mean) * T[0]
440 +(press3- temp3*densityW_mean) * T[1];
441 potential[nPhaseIdx] += (pressI+pcI-temp1*densityNW_mean) * T[2]
442 +(pressJ+pcJ-temp2*densityNW_mean) * T[0]
443 +(press3+pc3- temp3*densityNW_mean) * T[1];
445 if(halfedgesStored == 2)
447 int AdditionalIdx = problem().variables().index(additionalIsIt->outside());
448 CellData& cellDataAdditional = problem().variables().cellData(AdditionalIdx);
449 potential[wPhaseIdx] += (pressI-temp1*densityW_mean) * additionalT[2]
450 +(pressJ-temp2*densityW_mean) * additionalT[0]
451 +(problem().pressureModel().pressure(AdditionalIdx)
452 -(additionalIsIt->outside().geometry().center()*gravity_*densityW_mean)
454 potential[nPhaseIdx] += (pressI+pcI-temp1*densityNW_mean) * additionalT[2]
455 +(pressJ+pcJ-temp2*densityNW_mean) * additionalT[0]
456 +(problem().pressureModel().pressure(AdditionalIdx)
457 + cellDataAdditional.capillaryPressure()
458 -(additionalIsIt->outside().geometry().center()*gravity_*densityNW_mean)
462 else if(pressureType==pn)
464 potential[wPhaseIdx] += (pressI-pcI-temp1*densityW_mean) * T[2]
465 + (pressJ-pcJ-temp2*densityW_mean) * T[0]
466 + (press3-pc3- temp3*densityW_mean) * T[1];
467 potential[nPhaseIdx] += (pressI-temp1*densityNW_mean) * T[2]
468 + (pressJ-temp2*densityNW_mean) * T[0]
469 + (press3-temp3*densityNW_mean) * T[1];
471 if(halfedgesStored == 2)
473 int AdditionalIdx = problem().variables().index(additionalIsIt->outside());
474 CellData& cellDataAdditional = problem().variables().cellData(AdditionalIdx);
476 potential[wPhaseIdx] += (pressI-pcI-temp1*densityW_mean) * additionalT[2]
477 +(pressJ-pcJ-temp2*densityW_mean) * additionalT[0]
478 +(problem().pressureModel().pressure(AdditionalIdx)
479 - cellDataAdditional.capillaryPressure()
480 -(additionalIsIt->outside().geometry().center()*gravity_*densityW_mean)
482 potential[nPhaseIdx] += (pressI-temp1*densityNW_mean) * additionalT[2]
483 +(pressJ-temp2*densityNW_mean) * additionalT[0]
484 +(problem().pressureModel().pressure(AdditionalIdx)
485 -(additionalIsIt->outside().geometry().center()*gravity_*densityNW_mean))
492 Dune::FieldVector<bool, NumPhases> doUpwinding(
true);
493 PhaseVector lambda(0.);
494 for(
int phaseIdx = 0; phaseIdx < NumPhases; phaseIdx++)
497 if(phaseIdx == wPhaseIdx)
498 contiEqIdx = contiWEqIdx;
500 contiEqIdx = contiNEqIdx;
502 if(!this->impet_ or !this->restrictFluxInTransport_)
504 if(potential[phaseIdx] > 0.)
506 lambda[phaseIdx] = cellDataI.mobility(phaseIdx);
507 if(elementI.level()>neighbor.level())
508 cellDataI.setUpwindCell(intersectionIterator->indexInInside(), contiEqIdx,
true);
510 else if(potential[phaseIdx] < 0.)
512 lambda[phaseIdx] = cellDataJ.mobility(phaseIdx);
513 if(elementI.level()>neighbor.level())
514 cellDataI.setUpwindCell(intersectionIterator->indexInInside(), contiEqIdx,
false);
518 doUpwinding[phaseIdx] =
false;
519 if(elementI.level()>neighbor.level())
520 cellDataI.setUpwindCell(intersectionIterator->indexInInside(), contiEqIdx,
false);
522 cellDataJ.setUpwindCell(intersectionIterator->indexInOutside(), contiEqIdx,
false);
527 bool cellIwasUpwindCell;
529 if(elementI.level()>neighbor.level())
530 cellIwasUpwindCell = cellDataI.isUpwindCell(intersectionIterator->indexInInside(), contiEqIdx);
532 cellIwasUpwindCell = !cellDataJ.isUpwindCell(intersectionIterator->indexInOutside(), contiEqIdx);
534 if (potential[phaseIdx] > 0. && cellIwasUpwindCell)
535 lambda[phaseIdx] = cellDataI.mobility(phaseIdx);
536 else if (potential[phaseIdx] < 0. && !cellIwasUpwindCell)
537 lambda[phaseIdx] = cellDataJ.mobility(phaseIdx);
539 else if(this->restrictFluxInTransport_ == 2)
540 doUpwinding[phaseIdx] =
false;
545 if (potential[phaseIdx] > 0. && Dune::FloatCmp::ne<Scalar, Dune::FloatCmp::absolute>(cellDataJ.mobility(phaseIdx), 0.0, 1.0e-30))
546 lambda[phaseIdx] = cellDataI.mobility(phaseIdx);
548 else if (potential[phaseIdx] < 0. && Dune::FloatCmp::ne<Scalar, Dune::FloatCmp::absolute>(cellDataI.mobility(phaseIdx), 0.0, 1.0e-30))
549 lambda[phaseIdx] = cellDataJ.mobility(phaseIdx);
551 doUpwinding[phaseIdx] =
false;
555 if(!doUpwinding[phaseIdx])
558 if(cellDataI.mobility(phaseIdx)+cellDataJ.mobility(phaseIdx)==0.)
560 potential[phaseIdx] = 0;
565 fluxEntries[wCompIdx] -= potential[phaseIdx] /
volume
566 *
harmonicMean(cellDataI.massFraction(phaseIdx, wCompIdx) * cellDataI.mobility(phaseIdx) * cellDataI.density(phaseIdx),
567 cellDataJ.massFraction(phaseIdx, wCompIdx) * cellDataJ.mobility(phaseIdx) * cellDataJ.density(phaseIdx));
568 fluxEntries[nCompIdx] -= potential[phaseIdx] /
volume
569 *
harmonicMean(cellDataI.massFraction(phaseIdx, nCompIdx) * cellDataI.mobility(phaseIdx) * cellDataI.density(phaseIdx),
570 cellDataJ.massFraction(phaseIdx, nCompIdx) * cellDataJ.mobility(phaseIdx) * cellDataJ.density(phaseIdx));
574 timestepFlux[0] += max(0.,
575 - potential[phaseIdx] /
volume
576 *
harmonicMean(cellDataI.mobility(phaseIdx),cellDataJ.mobility(phaseIdx)));
578 timestepFlux[1] += max(0.,
579 potential[phaseIdx] /
volume
580 *
harmonicMean(cellDataI.mobility(phaseIdx),cellDataJ.mobility(phaseIdx))/SmobI[phaseIdx]);
583 this->averagedFaces_++;
584 #if DUNE_MINIMAL_DEBUG_LEVEL < 3
586 if(globalIdxI > globalIdxJ)
587 Dune::dinfo <<
"harmonicMean flux of phase" << phaseIdx <<
" used from cell" << globalIdxI<<
" into " << globalIdxJ
588 <<
" ; TE upwind I = "<< cellDataI.isUpwindCell(intersectionIterator->indexInInside(), contiEqIdx)
589 <<
" but pot = "<< potential[phaseIdx] <<
" \n";
593 potential[phaseIdx] = 0;
600 double velocityJIw = max((-lambda[wPhaseIdx] * potential[wPhaseIdx]) /
volume, 0.0);
601 double velocityIJw = max(( lambda[wPhaseIdx] * potential[wPhaseIdx]) /
volume, 0.0);
602 double velocityJIn = max((-lambda[nPhaseIdx] * potential[nPhaseIdx]) /
volume, 0.0);
603 double velocityIJn = max(( lambda[nPhaseIdx] * potential[nPhaseIdx]) /
volume, 0.0);
606 timestepFlux[0] += velocityJIw + velocityJIn;
608 double foutw = velocityIJw/SmobI[wPhaseIdx];
609 double foutn = velocityIJn/SmobI[nPhaseIdx];
612 if (isnan(foutw) || isinf(foutw) || foutw < 0) foutw = 0;
613 if (isnan(foutn) || isinf(foutn) || foutn < 0) foutn = 0;
614 timestepFlux[1] += foutw + foutn;
616 fluxEntries[wCompIdx] +=
617 velocityJIw * cellDataJ.massFraction(wPhaseIdx, wCompIdx) * densityWJ
618 - velocityIJw * cellDataI.massFraction(wPhaseIdx, wCompIdx) * densityWI
619 + velocityJIn * cellDataJ.massFraction(nPhaseIdx, wCompIdx) * densityNWJ
620 - velocityIJn * cellDataI.massFraction(nPhaseIdx, wCompIdx) * densityNWI;
621 fluxEntries[nCompIdx] +=
622 velocityJIw * cellDataJ.massFraction(wPhaseIdx, nCompIdx) * densityWJ
623 - velocityIJw * cellDataI.massFraction(wPhaseIdx, nCompIdx) * densityWI
624 + velocityJIn * cellDataJ.massFraction(nPhaseIdx, nCompIdx) * densityNWJ
625 - velocityIJn * cellDataI.massFraction(nPhaseIdx, nCompIdx) * densityNWI;
Contains a class to exchange entries of a vector.
Finite volume discretization of the component transport equation.
Defines the properties required for the adaptive sequential 2p2c models.
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:69
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type GetProp
get the type of a property
Definition: propertysystem.hh:141
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:150
Scalar volume(Shape shape, Scalar inscribedRadius)
Returns the volume of a given geometry based on the inscribed radius.
Definition: poreproperties.hh:73
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:57
virtual ~FV2dTransport2P2CAdaptive()
Definition: fv2dtransportadaptive.hh:120
FV2dTransport2P2CAdaptive(Problem &problem)
Constructs a FV2dTransport2P2CAdaptive object.
Definition: fv2dtransportadaptive.hh:114
static const int pressureType
gives kind of pressure used ( , , )
Definition: fv2dtransportadaptive.hh:129
Problem & problem_
Definition: fv2dtransportadaptive.hh:125
virtual void update(const Scalar t, Scalar &dt, TransportSolutionType &updateVec, bool impes=false)
Calculate the update vector and determine timestep size.
Definition: fv2dtransportadaptive.hh:152
bool enableMPFA
Definition: fv2dtransportadaptive.hh:127
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:350
Compositional transport step in a Finite Volume discretization.
Definition: fvtransport.hh:60
Definition: fvtransport.hh:112