3.3.0
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
37
38namespace Dumux {
57template<class TypeTag>
59{
63
65
67
69
70 enum
71 {
72 dim = GridView::dimension, dimWorld = GridView::dimensionworld,
73 NumPhases = getPropValue<TypeTag, Properties::NumPhases>()
74 };
75 enum
76 {
77 pw = Indices::pressureW,
78 pn = Indices::pressureN
79 };
80 enum
81 {
82 wPhaseIdx = Indices::wPhaseIdx, nPhaseIdx = Indices::nPhaseIdx,
83 wCompIdx = Indices::wPhaseIdx, nCompIdx = Indices::nPhaseIdx,
84 contiWEqIdx=Indices::contiWEqIdx, contiNEqIdx=Indices::contiNEqIdx
85 };
86
87 using IntersectionIterator = typename GridView::IntersectionIterator;
88
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>;
94
96 Problem& problem()
97 { return problem_; }
98
99 using LocalTimesteppingData = typename FVTransport2P2C<TypeTag>::LocalTimesteppingData;
100
101public:
102 virtual void update(const Scalar t, Scalar& dt, TransportSolutionType& updateVec, bool impes = false);
103
104 void getMpfaFlux(Dune::FieldVector<Scalar, 2>&, Dune::FieldVector<Scalar, 2>&,
105 const IntersectionIterator&, CellData&);
106
116 FV2dTransport2P2CAdaptive(Problem& problem) : FVTransport2P2C<TypeTag>(problem),
117 problem_(problem)
118 {
119 enableMPFA = getParam<bool>("GridAdapt.EnableMultiPointFluxApproximation");
120 }
121
123 {
124 }
125
126protected:
127 Problem& problem_;
128
131 static const int pressureType = getPropValue<TypeTag, Properties::PressureFormulation>();
132};
133
153template<class TypeTag>
154void FV2dTransport2P2CAdaptive<TypeTag>::update(const Scalar t, Scalar& dt, TransportSolutionType& updateVec, bool impet)
155{
156 this->impet_ = impet;
157 // initialize dt very large
158 dt = 1E100;
159 this->averagedFaces_ = 0;
160
161 // resize update vector and set to zero
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;
168
169 if (this->enableLocalTimeStepping())
170 {
171 if (this->timeStepData_.size() != size_)
172 this->timeStepData_.resize(size_);
173 }
174
175 //also resize PV vector if necessary
176 if(this->totalConcentration_.size() != size_)
177 {
178 this->totalConcentration_[wCompIdx].resize(size_);
179 this->totalConcentration_[nCompIdx].resize(size_);
180 // copy data //TODO: remove this, remove PM Transport Vector!!
181 // loop through all elements
182 for (int i = 0; i< problem().gridView().size(0); i++)
183 {
184 CellData& cellDataI = problem().variables().cellData(i);
185 for(int compIdx = 0; compIdx < getPropValue<TypeTag, Properties::NumComponents>(); compIdx++)
186 {
187 this->totalConcentration_[compIdx][i]
188 = cellDataI.totalConcentration(compIdx);
189 }
190 }
191 }
192
193 // Cell which restricts time step size
194 int restrictingCell = -1;
195
196 Dune::FieldVector<Scalar, 2> entries(0.), timestepFlux(0.);
197 // compute update vector
198 for (const auto& element : elements(problem().gridView()))
199 {
200 // get cell infos
201 int globalIdxI = problem().variables().index(element);
202 CellData& cellDataI = problem().variables().cellData(globalIdxI);
203
204 // some variables for time step calculation
205 double sumfactorin = 0;
206 double sumfactorout = 0;
207
208 if (this->enableLocalTimeStepping())
209 {
210 LocalTimesteppingData& localData = this->timeStepData_[globalIdxI];
211 for (int i = 0; i < 2*dim; i++)
212 {
213 if (localData.faceTargetDt[i] < this->accumulatedDt_ + this->dtThreshold_)
214 {
215 localData.faceFluxes[i] = 0.0;
216 }
217 }
218 }
219
220 // run through all intersections with neighbors and boundary
221 const auto isEndIt = problem().gridView().iend(element);
222 for (auto isIt = problem().gridView().ibegin(element); isIt != isEndIt; ++isIt)
223 {
224 const auto& intersection = *isIt;
225
226 int indexInInside = intersection.indexInInside();
227
228 // handle interior face
229 if (intersection.neighbor())
230 {
231 if (enableMPFA && (intersection.outside().level() != element.level()))
232 getMpfaFlux(entries, timestepFlux, isIt, cellDataI);
233 else
234 this->getFlux(entries, timestepFlux, intersection, cellDataI);
235 }
236
237 // Boundary Face
238 if (intersection.boundary())
239 {
240 this->getFluxOnBoundary(entries, timestepFlux, intersection, cellDataI);
241 }
242
243 if (this->enableLocalTimeStepping())
244 {
245 LocalTimesteppingData& localData = this->timeStepData_[globalIdxI];
246
247 if (localData.faceTargetDt[indexInInside] < this->accumulatedDt_ + this->dtThreshold_)
248 {
249 localData.faceFluxes[indexInInside] += entries;
250 }
251 }
252 else
253 {
254 // add to update vector
255 updateVec[wCompIdx][globalIdxI] += entries[wCompIdx];
256 updateVec[nCompIdx][globalIdxI] += entries[nCompIdx];
257 }
258
259 // for time step calculation
260 sumfactorin += timestepFlux[0];
261 sumfactorout += timestepFlux[1];
262
263 }// end all intersections
264
265 if (this->enableLocalTimeStepping())
266 {
267 LocalTimesteppingData& localData = this->timeStepData_[globalIdxI];
268 for (int i=0; i < 2*dim; i++)
269 {
270 updateVec[wCompIdx][globalIdxI] += localData.faceFluxes[i][wCompIdx];
271 updateVec[nCompIdx][globalIdxI] += localData.faceFluxes[i][nCompIdx];
272 }
273 }
274
275 /*********** Handle source term ***************/
276 PrimaryVariables q(NAN);
277 problem().source(q, element);
278 updateVec[wCompIdx][globalIdxI] += q[contiWEqIdx];
279 updateVec[nCompIdx][globalIdxI] += q[contiNEqIdx];
280
281 // account for porosity
282 using std::max;
283 sumfactorin = max(sumfactorin,sumfactorout)
284 / problem().spatialParams().porosity(element);
285
286 //calculate time step
287 if (this->enableLocalTimeStepping())
288 {
289 this->timeStepData_[globalIdxI].dt = 1./sumfactorin;
290 if ( 1./sumfactorin < dt)
291 {
292 dt = 1./sumfactorin;
293 restrictingCell= globalIdxI;
294 }
295 }
296 else
297 {
298 if ( 1./sumfactorin < dt)
299 {
300 dt = 1./sumfactorin;
301 restrictingCell= globalIdxI;
302 }
303 }
304 } // end grid traversal
305
306#if HAVE_MPI
307 // communicate updated values
309 using ElementMapper = typename SolutionTypes::ElementMapper;
311 for (int i = 0; i < updateVec.size(); i++)
312 {
313 DataHandle dataHandle(problem_.variables().elementMapper(), updateVec[i]);
314 problem_.gridView().template communicate<DataHandle>(dataHandle,
315 Dune::InteriorBorder_All_Interface,
316 Dune::ForwardCommunication);
317 }
318 dt = problem_.gridView().comm().min(dt);
319#endif
320
321 if(impet)
322 {
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;
327 }
328 return;
329}
330
351template<class TypeTag>
352void FV2dTransport2P2CAdaptive<TypeTag>::getMpfaFlux(Dune::FieldVector<Scalar, 2>& fluxEntries, Dune::FieldVector<Scalar, 2>& timestepFlux,
353 const IntersectionIterator& intersectionIterator, CellData& cellDataI)
354{
355 fluxEntries = 0.;
356 timestepFlux = 0.;
357 // cell information
358 auto elementI = intersectionIterator->inside();
359 int globalIdxI = problem().variables().index(elementI);
360
361 // get position
362 const GlobalPosition globalPos = elementI.geometry().center();
363 const GlobalPosition& gravity_ = problem().gravity();
364 // cell volume, assume linear map here
365 Scalar volume = elementI.geometry().volume();
366
367 // get values of cell I
368 Scalar pressI = problem().pressureModel().pressure(globalIdxI);
369 Scalar pcI = cellDataI.capillaryPressure();
370
371 // old material law interface is deprecated: Replace this by
372 // const auto& fluidMatrixInteraction = problem().spatialParams.fluidMatrixInteractionAtPos(elementI.geometry().center());
373 // after the release of 3.3, when the deprecated interface is no longer supported
374 const auto fluidMatrixInteraction = Deprecated::makePcKrSw(Scalar{}, problem().spatialParams(), elementI);
375
376 PhaseVector SmobI(0.);
377 using std::max;
378 SmobI[wPhaseIdx] = max((cellDataI.saturation(wPhaseIdx)
379 - fluidMatrixInteraction.pcSwCurve().effToAbsParams().swr())
380 , 1e-2);
381 SmobI[nPhaseIdx] = max((cellDataI.saturation(nPhaseIdx)
382 - fluidMatrixInteraction.pcSwCurve().effToAbsParams().snr())
383 , 1e-2);
384
385 Scalar densityWI (0.), densityNWI(0.);
386 densityWI= cellDataI.density(wPhaseIdx);
387 densityNWI = cellDataI.density(nPhaseIdx);
388
389 PhaseVector potential(0.);
390
391 // access neighbor
392 auto neighbor = intersectionIterator->outside();
393 int globalIdxJ = problem().variables().index(neighbor);
394 CellData& cellDataJ = problem().variables().cellData(globalIdxJ);
395
396 // neighbor cell center in global coordinates
397 const GlobalPosition& globalPosNeighbor = neighbor.geometry().center();
398
399 // distance vector between barycenters
400 GlobalPosition distVec = globalPosNeighbor - globalPos;
401 // compute distance between cell centers
402 Scalar dist = distVec.two_norm();
403
404 GlobalPosition unitDistVec(distVec);
405 unitDistVec /= dist;
406
407 // phase densities in neighbor
408 Scalar densityWJ (0.), densityNWJ(0.);
409 densityWJ = cellDataJ.density(wPhaseIdx);
410 densityNWJ = cellDataJ.density(nPhaseIdx);
411
412 // average phase densities with central weighting
413 double densityW_mean = (densityWI + densityWJ) * 0.5;
414 double densityNW_mean = (densityNWI + densityNWJ) * 0.5;
415
416 double pressJ = problem().pressureModel().pressure(globalIdxJ);
417 Scalar pcJ = cellDataJ.capillaryPressure();
418
419 // determine potentials for upwind
421 GlobalPosition globalPos3(0.);
422 int globalIdx3=-1;
423 TransmissivityMatrix T(0.);
424 IntersectionIterator additionalIsIt = intersectionIterator;
425 TransmissivityMatrix additionalT(0.);
426
427 int halfedgesStored
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 );
433
434 // acess cell 3 and prepare mpfa
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_;
441 if(pressureType==pw)
442 {
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];
449 // second half edge, if there is one
450 if(halfedgesStored == 2)
451 {
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)
458 ) * additionalT[1];
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)
464 ) * additionalT[1];
465 }
466 }
467 else if(pressureType==pn)
468 {
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];
475 // second half edge, if there is one
476 if(halfedgesStored == 2)
477 {
478 int AdditionalIdx = problem().variables().index(additionalIsIt->outside());
479 CellData& cellDataAdditional = problem().variables().cellData(AdditionalIdx);
480
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)
486 ) * additionalT[1];
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))
491 * additionalT[1];
492 }
493 }
494 //end of mpfa specific stuff
495
496 // determine upwinding direction, perform upwinding if possible
497 Dune::FieldVector<bool, NumPhases> doUpwinding(true);
498 PhaseVector lambda(0.);
499 for(int phaseIdx = 0; phaseIdx < NumPhases; phaseIdx++)
500 {
501 int contiEqIdx = 0;
502 if(phaseIdx == wPhaseIdx)
503 contiEqIdx = contiWEqIdx;
504 else
505 contiEqIdx = contiNEqIdx;
506
507 if(!this->impet_ or !this->restrictFluxInTransport_) // perform a strict uwpind scheme
508 {
509 if(potential[phaseIdx] > 0.)
510 {
511 lambda[phaseIdx] = cellDataI.mobility(phaseIdx);
512 if(elementI.level()>neighbor.level())
513 cellDataI.setUpwindCell(intersectionIterator->indexInInside(), contiEqIdx, true);
514 }
515 else if(potential[phaseIdx] < 0.)
516 {
517 lambda[phaseIdx] = cellDataJ.mobility(phaseIdx);
518 if(elementI.level()>neighbor.level())
519 cellDataI.setUpwindCell(intersectionIterator->indexInInside(), contiEqIdx, false);
520 }
521 else
522 {
523 doUpwinding[phaseIdx] = false;
524 if(elementI.level()>neighbor.level())
525 cellDataI.setUpwindCell(intersectionIterator->indexInInside(), contiEqIdx, false);
526 else
527 cellDataJ.setUpwindCell(intersectionIterator->indexInOutside(), contiEqIdx, false);
528 }
529 }
530 else // Transport after PE with check on flow direction
531 {
532 bool cellIwasUpwindCell;
533 //get the information from smaller (higher level) cell, as its IS is unique
534 if(elementI.level()>neighbor.level())
535 cellIwasUpwindCell = cellDataI.isUpwindCell(intersectionIterator->indexInInside(), contiEqIdx);
536 else // reverse neighbors information gathered
537 cellIwasUpwindCell = !cellDataJ.isUpwindCell(intersectionIterator->indexInOutside(), contiEqIdx);
538
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);
543 // potential direction does not coincide with that of P.E.
544 else if(this->restrictFluxInTransport_ == 2) // perform central averaging for all direction changes
545 doUpwinding[phaseIdx] = false;
546 else // i.e. restrictFluxInTransport == 1
547 {
548 //check if harmonic weighting is necessary
549 // check if outflow induce neglected (i.e. mob=0) phase flux
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);
552 // check if inflow induce neglected phase flux
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);
555 else
556 doUpwinding[phaseIdx] = false;
557 }
558
559 // do not perform upwinding if so desired
560 if(!doUpwinding[phaseIdx])
561 {
562 //a) no flux if there wont be any flux regardless how to average/upwind
563 if(cellDataI.mobility(phaseIdx)+cellDataJ.mobility(phaseIdx)==0.)
564 {
565 potential[phaseIdx] = 0;
566 continue;
567 }
568
569 //b) perform harmonic averaging
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));
576 // c) timestep control
577 // for timestep control : influx
578 using std::max;
579 timestepFlux[0] += max(0.,
580 - potential[phaseIdx] / volume
581 * harmonicMean(cellDataI.mobility(phaseIdx),cellDataJ.mobility(phaseIdx)));
582 // outflux
583 timestepFlux[1] += max(0.,
584 potential[phaseIdx] / volume
585 * harmonicMean(cellDataI.mobility(phaseIdx),cellDataJ.mobility(phaseIdx))/SmobI[phaseIdx]);
586
587 //d) output (only for one side)
588 this->averagedFaces_++;
589 #if DUNE_MINIMAL_DEBUG_LEVEL < 3
590 // verbose (only for one side)
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";
595 #endif
596
597 //e) stop further standard calculations
598 potential[phaseIdx] = 0;
599 }
600 }
601 }
602
603 // calculate and standardized velocity
604 using std::max;
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);
609
610 // for timestep control
611 timestepFlux[0] += velocityJIw + velocityJIn;
612
613 double foutw = velocityIJw/SmobI[wPhaseIdx];
614 double foutn = velocityIJn/SmobI[nPhaseIdx];
615 using std::isnan;
616 using std::isinf;
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;
620
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;
631 return;
632}
633
634} // end namespace Dumux
635#endif
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.
Helpers for deprecation.
constexpr Scalar harmonicMean(Scalar x, Scalar y, Scalar wx=1.0, Scalar wy=1.0) noexcept
Calculate the (weighted) harmonic mean of two scalar values.
Definition: math.hh:68
Definition: adapt.hh:29
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type GetProp
get the type of a property (equivalent to old macro GET_PROP(...))
Definition: propertysystem.hh:140
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
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