3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
fv2dtransportadaptive.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_FV2DTRANSPORT2P2C_ADAPTIVE_HH
25#define DUMUX_FV2DTRANSPORT2P2C_ADAPTIVE_HH
26
27#include <dune/grid/common/gridenums.hh>
28#include <dune/common/float_cmp.hh>
29
30#include <dumux/common/math.hh>
32
33#include "adaptiveproperties.hh"
34#include "fvtransport.hh"
35
36namespace Dumux {
55template<class TypeTag>
57{
61
63
65
67
68 enum
69 {
70 dim = GridView::dimension, dimWorld = GridView::dimensionworld,
71 NumPhases = getPropValue<TypeTag, Properties::NumPhases>()
72 };
73 enum
74 {
75 pw = Indices::pressureW,
76 pn = Indices::pressureN
77 };
78 enum
79 {
80 wPhaseIdx = Indices::wPhaseIdx, nPhaseIdx = Indices::nPhaseIdx,
81 wCompIdx = Indices::wPhaseIdx, nCompIdx = Indices::nPhaseIdx,
82 contiWEqIdx=Indices::contiWEqIdx, contiNEqIdx=Indices::contiNEqIdx
83 };
84
85 using IntersectionIterator = typename GridView::IntersectionIterator;
86
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>;
92
94 Problem& problem()
95 { return problem_; }
96
97 using LocalTimesteppingData = typename FVTransport2P2C<TypeTag>::LocalTimesteppingData;
98
99public:
100 virtual void update(const Scalar t, Scalar& dt, TransportSolutionType& updateVec, bool impes = false);
101
102 void getMpfaFlux(Dune::FieldVector<Scalar, 2>&, Dune::FieldVector<Scalar, 2>&,
103 const IntersectionIterator&, CellData&);
104
114 FV2dTransport2P2CAdaptive(Problem& problem) : FVTransport2P2C<TypeTag>(problem),
115 problem_(problem)
116 {
117 enableMPFA = getParam<bool>("GridAdapt.EnableMultiPointFluxApproximation");
118 }
119
121 {
122 }
123
124protected:
125 Problem& problem_;
126
129 static const int pressureType = getPropValue<TypeTag, Properties::PressureFormulation>();
130};
131
151template<class TypeTag>
152void FV2dTransport2P2CAdaptive<TypeTag>::update(const Scalar t, Scalar& dt, TransportSolutionType& updateVec, bool impet)
153{
154 this->impet_ = impet;
155 // initialize dt very large
156 dt = 1E100;
157 this->averagedFaces_ = 0;
158
159 // resize update vector and set to zero
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;
166
167 if (this->enableLocalTimeStepping())
168 {
169 if (this->timeStepData_.size() != size_)
170 this->timeStepData_.resize(size_);
171 }
172
173 //also resize PV vector if necessary
174 if(this->totalConcentration_.size() != size_)
175 {
176 this->totalConcentration_[wCompIdx].resize(size_);
177 this->totalConcentration_[nCompIdx].resize(size_);
178 // copy data //TODO: remove this, remove PM Transport Vector!!
179 // loop through all elements
180 for (int i = 0; i< problem().gridView().size(0); i++)
181 {
182 CellData& cellDataI = problem().variables().cellData(i);
183 for(int compIdx = 0; compIdx < getPropValue<TypeTag, Properties::NumComponents>(); compIdx++)
184 {
185 this->totalConcentration_[compIdx][i]
186 = cellDataI.totalConcentration(compIdx);
187 }
188 }
189 }
190
191 // Cell which restricts time step size
192 int restrictingCell = -1;
193
194 Dune::FieldVector<Scalar, 2> entries(0.), timestepFlux(0.);
195 // compute update vector
196 for (const auto& element : elements(problem().gridView()))
197 {
198 // get cell infos
199 int globalIdxI = problem().variables().index(element);
200 CellData& cellDataI = problem().variables().cellData(globalIdxI);
201
202 // some variables for time step calculation
203 double sumfactorin = 0;
204 double sumfactorout = 0;
205
206 if (this->enableLocalTimeStepping())
207 {
208 LocalTimesteppingData& localData = this->timeStepData_[globalIdxI];
209 for (int i = 0; i < 2*dim; i++)
210 {
211 if (localData.faceTargetDt[i] < this->accumulatedDt_ + this->dtThreshold_)
212 {
213 localData.faceFluxes[i] = 0.0;
214 }
215 }
216 }
217
218 // run through all intersections with neighbors and boundary
219 const auto isEndIt = problem().gridView().iend(element);
220 for (auto isIt = problem().gridView().ibegin(element); isIt != isEndIt; ++isIt)
221 {
222 const auto& intersection = *isIt;
223
224 int indexInInside = intersection.indexInInside();
225
226 // handle interior face
227 if (intersection.neighbor())
228 {
229 if (enableMPFA && (intersection.outside().level() != element.level()))
230 getMpfaFlux(entries, timestepFlux, isIt, cellDataI);
231 else
232 this->getFlux(entries, timestepFlux, intersection, cellDataI);
233 }
234
235 // Boundary Face
236 if (intersection.boundary())
237 {
238 this->getFluxOnBoundary(entries, timestepFlux, intersection, cellDataI);
239 }
240
241 if (this->enableLocalTimeStepping())
242 {
243 LocalTimesteppingData& localData = this->timeStepData_[globalIdxI];
244
245 if (localData.faceTargetDt[indexInInside] < this->accumulatedDt_ + this->dtThreshold_)
246 {
247 localData.faceFluxes[indexInInside] += entries;
248 }
249 }
250 else
251 {
252 // add to update vector
253 updateVec[wCompIdx][globalIdxI] += entries[wCompIdx];
254 updateVec[nCompIdx][globalIdxI] += entries[nCompIdx];
255 }
256
257 // for time step calculation
258 sumfactorin += timestepFlux[0];
259 sumfactorout += timestepFlux[1];
260
261 }// end all intersections
262
263 if (this->enableLocalTimeStepping())
264 {
265 LocalTimesteppingData& localData = this->timeStepData_[globalIdxI];
266 for (int i=0; i < 2*dim; i++)
267 {
268 updateVec[wCompIdx][globalIdxI] += localData.faceFluxes[i][wCompIdx];
269 updateVec[nCompIdx][globalIdxI] += localData.faceFluxes[i][nCompIdx];
270 }
271 }
272
273 /*********** Handle source term ***************/
274 PrimaryVariables q(NAN);
275 problem().source(q, element);
276 updateVec[wCompIdx][globalIdxI] += q[contiWEqIdx];
277 updateVec[nCompIdx][globalIdxI] += q[contiNEqIdx];
278
279 // account for porosity
280 using std::max;
281 sumfactorin = max(sumfactorin,sumfactorout)
282 / problem().spatialParams().porosity(element);
283
284 //calculate time step
285 if (this->enableLocalTimeStepping())
286 {
287 this->timeStepData_[globalIdxI].dt = 1./sumfactorin;
288 if ( 1./sumfactorin < dt)
289 {
290 dt = 1./sumfactorin;
291 restrictingCell= globalIdxI;
292 }
293 }
294 else
295 {
296 if ( 1./sumfactorin < dt)
297 {
298 dt = 1./sumfactorin;
299 restrictingCell= globalIdxI;
300 }
301 }
302 } // end grid traversal
303
304#if HAVE_MPI
305 // communicate updated values
307 using ElementMapper = typename SolutionTypes::ElementMapper;
309 for (int i = 0; i < updateVec.size(); i++)
310 {
311 DataHandle dataHandle(problem_.variables().elementMapper(), updateVec[i]);
312 problem_.gridView().template communicate<DataHandle>(dataHandle,
313 Dune::InteriorBorder_All_Interface,
314 Dune::ForwardCommunication);
315 }
316 dt = problem_.gridView().comm().min(dt);
317#endif
318
319 if(impet)
320 {
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;
325 }
326 return;
327}
328
349template<class TypeTag>
350void FV2dTransport2P2CAdaptive<TypeTag>::getMpfaFlux(Dune::FieldVector<Scalar, 2>& fluxEntries, Dune::FieldVector<Scalar, 2>& timestepFlux,
351 const IntersectionIterator& intersectionIterator, CellData& cellDataI)
352{
353 fluxEntries = 0.;
354 timestepFlux = 0.;
355 // cell information
356 auto elementI = intersectionIterator->inside();
357 int globalIdxI = problem().variables().index(elementI);
358
359 // get position
360 const GlobalPosition globalPos = elementI.geometry().center();
361 const GlobalPosition& gravity_ = problem().gravity();
362 // cell volume, assume linear map here
363 Scalar volume = elementI.geometry().volume();
364
365 // get values of cell I
366 Scalar pressI = problem().pressureModel().pressure(globalIdxI);
367 Scalar pcI = cellDataI.capillaryPressure();
368
369 const auto fluidMatrixInteraction = problem().spatialParams().fluidMatrixInteractionAtPos(elementI.geometry().center());
370
371 PhaseVector SmobI(0.);
372 using std::max;
373 SmobI[wPhaseIdx] = max((cellDataI.saturation(wPhaseIdx)
374 - fluidMatrixInteraction.pcSwCurve().effToAbsParams().swr())
375 , 1e-2);
376 SmobI[nPhaseIdx] = max((cellDataI.saturation(nPhaseIdx)
377 - fluidMatrixInteraction.pcSwCurve().effToAbsParams().snr())
378 , 1e-2);
379
380 Scalar densityWI (0.), densityNWI(0.);
381 densityWI= cellDataI.density(wPhaseIdx);
382 densityNWI = cellDataI.density(nPhaseIdx);
383
384 PhaseVector potential(0.);
385
386 // access neighbor
387 auto neighbor = intersectionIterator->outside();
388 int globalIdxJ = problem().variables().index(neighbor);
389 CellData& cellDataJ = problem().variables().cellData(globalIdxJ);
390
391 // neighbor cell center in global coordinates
392 const GlobalPosition& globalPosNeighbor = neighbor.geometry().center();
393
394 // distance vector between barycenters
395 GlobalPosition distVec = globalPosNeighbor - globalPos;
396 // compute distance between cell centers
397 Scalar dist = distVec.two_norm();
398
399 GlobalPosition unitDistVec(distVec);
400 unitDistVec /= dist;
401
402 // phase densities in neighbor
403 Scalar densityWJ (0.), densityNWJ(0.);
404 densityWJ = cellDataJ.density(wPhaseIdx);
405 densityNWJ = cellDataJ.density(nPhaseIdx);
406
407 // average phase densities with central weighting
408 double densityW_mean = (densityWI + densityWJ) * 0.5;
409 double densityNW_mean = (densityNWI + densityNWJ) * 0.5;
410
411 double pressJ = problem().pressureModel().pressure(globalIdxJ);
412 Scalar pcJ = cellDataJ.capillaryPressure();
413
414 // determine potentials for upwind
416 GlobalPosition globalPos3(0.);
417 int globalIdx3=-1;
418 TransmissivityMatrix T(0.);
419 IntersectionIterator additionalIsIt = intersectionIterator;
420 TransmissivityMatrix additionalT(0.);
421
422 int halfedgesStored
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 );
428
429 // acess cell 3 and prepare mpfa
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_;
436 if(pressureType==pw)
437 {
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];
444 // second half edge, if there is one
445 if(halfedgesStored == 2)
446 {
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)
453 ) * additionalT[1];
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)
459 ) * additionalT[1];
460 }
461 }
462 else if(pressureType==pn)
463 {
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];
470 // second half edge, if there is one
471 if(halfedgesStored == 2)
472 {
473 int AdditionalIdx = problem().variables().index(additionalIsIt->outside());
474 CellData& cellDataAdditional = problem().variables().cellData(AdditionalIdx);
475
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)
481 ) * additionalT[1];
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))
486 * additionalT[1];
487 }
488 }
489 //end of mpfa specific stuff
490
491 // determine upwinding direction, perform upwinding if possible
492 Dune::FieldVector<bool, NumPhases> doUpwinding(true);
493 PhaseVector lambda(0.);
494 for(int phaseIdx = 0; phaseIdx < NumPhases; phaseIdx++)
495 {
496 int contiEqIdx = 0;
497 if(phaseIdx == wPhaseIdx)
498 contiEqIdx = contiWEqIdx;
499 else
500 contiEqIdx = contiNEqIdx;
501
502 if(!this->impet_ or !this->restrictFluxInTransport_) // perform a strict uwpind scheme
503 {
504 if(potential[phaseIdx] > 0.)
505 {
506 lambda[phaseIdx] = cellDataI.mobility(phaseIdx);
507 if(elementI.level()>neighbor.level())
508 cellDataI.setUpwindCell(intersectionIterator->indexInInside(), contiEqIdx, true);
509 }
510 else if(potential[phaseIdx] < 0.)
511 {
512 lambda[phaseIdx] = cellDataJ.mobility(phaseIdx);
513 if(elementI.level()>neighbor.level())
514 cellDataI.setUpwindCell(intersectionIterator->indexInInside(), contiEqIdx, false);
515 }
516 else
517 {
518 doUpwinding[phaseIdx] = false;
519 if(elementI.level()>neighbor.level())
520 cellDataI.setUpwindCell(intersectionIterator->indexInInside(), contiEqIdx, false);
521 else
522 cellDataJ.setUpwindCell(intersectionIterator->indexInOutside(), contiEqIdx, false);
523 }
524 }
525 else // Transport after PE with check on flow direction
526 {
527 bool cellIwasUpwindCell;
528 //get the information from smaller (higher level) cell, as its IS is unique
529 if(elementI.level()>neighbor.level())
530 cellIwasUpwindCell = cellDataI.isUpwindCell(intersectionIterator->indexInInside(), contiEqIdx);
531 else // reverse neighbors information gathered
532 cellIwasUpwindCell = !cellDataJ.isUpwindCell(intersectionIterator->indexInOutside(), contiEqIdx);
533
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);
538 // potential direction does not coincide with that of P.E.
539 else if(this->restrictFluxInTransport_ == 2) // perform central averaging for all direction changes
540 doUpwinding[phaseIdx] = false;
541 else // i.e. restrictFluxInTransport == 1
542 {
543 //check if harmonic weighting is necessary
544 // check if outflow induce neglected (i.e. mob=0) phase flux
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);
547 // check if inflow induce neglected phase flux
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);
550 else
551 doUpwinding[phaseIdx] = false;
552 }
553
554 // do not perform upwinding if so desired
555 if(!doUpwinding[phaseIdx])
556 {
557 //a) no flux if there wont be any flux regardless how to average/upwind
558 if(cellDataI.mobility(phaseIdx)+cellDataJ.mobility(phaseIdx)==0.)
559 {
560 potential[phaseIdx] = 0;
561 continue;
562 }
563
564 //b) perform harmonic averaging
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));
571 // c) timestep control
572 // for timestep control : influx
573 using std::max;
574 timestepFlux[0] += max(0.,
575 - potential[phaseIdx] / volume
576 * harmonicMean(cellDataI.mobility(phaseIdx),cellDataJ.mobility(phaseIdx)));
577 // outflux
578 timestepFlux[1] += max(0.,
579 potential[phaseIdx] / volume
580 * harmonicMean(cellDataI.mobility(phaseIdx),cellDataJ.mobility(phaseIdx))/SmobI[phaseIdx]);
581
582 //d) output (only for one side)
583 this->averagedFaces_++;
584 #if DUNE_MINIMAL_DEBUG_LEVEL < 3
585 // verbose (only for one side)
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";
590 #endif
591
592 //e) stop further standard calculations
593 potential[phaseIdx] = 0;
594 }
595 }
596 }
597
598 // calculate and standardized velocity
599 using std::max;
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);
604
605 // for timestep control
606 timestepFlux[0] += velocityJIw + velocityJIn;
607
608 double foutw = velocityIJw/SmobI[wPhaseIdx];
609 double foutn = velocityIJn/SmobI[nPhaseIdx];
610 using std::isnan;
611 using std::isinf;
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;
615
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;
626 return;
627}
628
629} // end namespace Dumux
630#endif
Define some often used mathematical functions.
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.
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
Definition: adapt.hh:29
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