3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
fv3dtransportadaptive.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_FV3DTRANSPORT2P2C_ADAPTIVE_HH
25#define DUMUX_FV3DTRANSPORT2P2C_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 {
52template<class TypeTag>
54{
58
60 using MaterialLaw = typename SpatialParams::MaterialLaw;
61
64
67
69
71
72 enum
73 {
74 dim = GridView::dimension, dimWorld = GridView::dimensionworld,
75 NumPhases = getPropValue<TypeTag, Properties::NumPhases>()
76 };
77 enum
78 {
79 pw = Indices::pressureW,
80 pn = Indices::pressureN,
81 Sw = Indices::saturationW,
82 Sn = Indices::saturationN
83 };
84 enum
85 {
86 wPhaseIdx = Indices::wPhaseIdx, nPhaseIdx = Indices::nPhaseIdx,
87 wCompIdx = Indices::wPhaseIdx, nCompIdx = Indices::nPhaseIdx,
88 contiWEqIdx=Indices::contiWEqIdx, contiNEqIdx=Indices::contiNEqIdx
89 };
90
91 using Element = typename GridView::Traits::template Codim<0>::Entity;
92 using Grid = typename GridView::Grid;
93 using Intersection = typename GridView::Intersection;
94 using IntersectionIterator = typename GridView::IntersectionIterator;
95
96 using TransmissivityMatrix = Dune::FieldVector<Scalar,dim+1>;
97 using GlobalPosition = Dune::FieldVector<Scalar, dimWorld>;
98 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
99 using PhaseVector = Dune::FieldVector<Scalar, NumPhases>;
101
103 Problem& problem()
104 { return problem_; }
105
106public:
107 virtual void update(const Scalar t, Scalar& dt, TransportSolutionType& updateVec,
108 bool impes = false);
109
110 void getMpfaFlux(Dune::FieldVector<Scalar, 2>&, Dune::FieldVector<Scalar, 2>&,
111 const IntersectionIterator&, CellData&);
112
122 FV3dTransport2P2CAdaptive(Problem& problem) : FVTransport2P2C<TypeTag>(problem),
123 problem_(problem), enableMPFA(false)
124 {
125 enableMPFA = getParam<bool>("GridAdapt.EnableMultiPointFluxApproximation");
126 }
127
129 {
130 }
131
132protected:
133 Problem& problem_;
134
136
138 static const int pressureType = getPropValue<TypeTag, Properties::PressureFormulation>();
139};
140
161template<class TypeTag>
162void FV3dTransport2P2CAdaptive<TypeTag>::update(const Scalar t, Scalar& dt,
163 TransportSolutionType& updateVec, bool impet)
164{
165 this->impet_ = impet;
166
167 // initialize dt very large
168 dt = 1E100;
169 this->averagedFaces_ = 0;
170
171 // resize update vector and set to zero
172 int size_ = problem_.gridView().size(0);
173 updateVec.resize(getPropValue<TypeTag, Properties::NumComponents>());
174 updateVec[wCompIdx].resize(size_);
175 updateVec[nCompIdx].resize(size_);
176 updateVec[wCompIdx] = 0;
177 updateVec[nCompIdx] = 0;
178 //also resize PV vector if necessary
179 if(this->totalConcentration_.size() != size_)
180 {
181 this->totalConcentration_[wCompIdx].resize(size_);
182 this->totalConcentration_[nCompIdx].resize(size_);
183 // copy data //TODO: remove this, remove PM Transport Vector!!
184 // loop thorugh all elements
185 for (int i = 0; i< problem().gridView().size(0); i++)
186 {
187 CellData& cellDataI = problem().variables().cellData(i);
188 for(int compIdx = 0; compIdx < getPropValue<TypeTag, Properties::NumComponents>(); compIdx++)
189 {
190 this->totalConcentration_[compIdx][i]
191 = cellDataI.totalConcentration(compIdx);
192 }
193 }
194 }
195 if (this->localTimeStepping_)
196 {
197 if (this->timeStepData_.size() != size_)
198 this->timeStepData_.resize(size_);
199 }
200
201 // Cell which restricts time step size
202 int restrictingCell = -1;
203
204 Dune::FieldVector<Scalar, 2> entries(0.), timestepFlux(0.);
205 // compute update vector
206 for (const auto& element : elements(problem().gridView()))
207 {
208 // get cell infos
209 int globalIdxI = problem().variables().index(element);
210 CellData& cellDataI = problem().variables().cellData(globalIdxI);
211
212 if(!impet && cellDataI.subdomain()!=2) // estimate only necessary in subdomain
213 continue;
214
215 // some variables for time step calculation
216 double sumfactorin = 0;
217 double sumfactorout = 0;
218
219 // run through all intersections with neighbors and boundary
220 const auto isEndIt = problem().gridView().iend(element);
221 for (auto isIt = problem().gridView().ibegin(element); isIt != isEndIt; ++isIt)
222 {
223 const auto& intersection = *isIt;
224
225 // handle interior face
226 if (intersection.neighbor())
227 {
228 if (enableMPFA && intersection.outside().level() != element.level())
229 getMpfaFlux(entries, timestepFlux, isIt, cellDataI);
230 else
231 this->getFlux(entries, timestepFlux, intersection, cellDataI);
232 }
233
234 // Boundary Face
235 if (intersection.boundary())
236 {
237 this->getFluxOnBoundary(entries, timestepFlux, intersection, cellDataI);
238 }
239
240 if (this->localTimeStepping_)
241 {
242 int indexInInside = intersection.indexInInside();
243 typename FVTransport2P2C<TypeTag>::LocalTimesteppingData& localData = this->timeStepData_[globalIdxI];
244 if (localData.faceTargetDt[indexInInside] < this->accumulatedDt_ + this->dtThreshold_)
245 {
246 localData.faceFluxes[indexInInside] = entries;
247 }
248 }
249 else
250 {
251 // add to update vector
252 updateVec[wCompIdx][globalIdxI] += entries[wCompIdx];
253 updateVec[nCompIdx][globalIdxI] += entries[nCompIdx];
254 }
255 // for time step calculation
256 sumfactorin += timestepFlux[0];
257 sumfactorout += timestepFlux[1];
258
259 }// end all intersections
260
261 if (this->localTimeStepping_)
262 {
263 typename FVTransport2P2C<TypeTag>::LocalTimesteppingData& localData = this->timeStepData_[globalIdxI];
264 for (int i=0; i < 2*dim; i++)
265 {
266 updateVec[wCompIdx][globalIdxI] += localData.faceFluxes[i][wCompIdx];
267 updateVec[nCompIdx][globalIdxI] += localData.faceFluxes[i][nCompIdx];
268 }
269 }
270
271 /*********** Handle source term ***************/
272 PrimaryVariables q(NAN);
273 problem().source(q, element);
274 updateVec[wCompIdx][globalIdxI] += q[contiWEqIdx];
275 updateVec[nCompIdx][globalIdxI] += q[contiNEqIdx];
276
277 // account for porosity
278 using std::max;
279 sumfactorin = max(sumfactorin,sumfactorout)
280 / problem().spatialParams().porosity(element);
281
282 //calculate time step
283 if (this->localTimeStepping_)
284 {
285 this->timeStepData_[globalIdxI].dt = 1./sumfactorin;
286 if ( 1./sumfactorin < dt)
287 {
288 dt = 1./sumfactorin;
289 restrictingCell= globalIdxI;
290 }
291 }
292 else
293 {
294 if ( 1./sumfactorin < dt)
295 {
296 dt = 1./sumfactorin;
297 restrictingCell= globalIdxI;
298 }
299 }
300 } // end grid traversal
301
302#if HAVE_MPI
303 // communicate updated values
305 using ElementMapper = typename SolutionTypes::ElementMapper;
307 for (int i = 0; i < updateVec.size(); i++)
308 {
309 DataHandle dataHandle(problem().variables().elementMapper(), updateVec[i]);
310 problem_.gridView().template communicate<DataHandle>(dataHandle,
311 Dune::InteriorBorder_All_Interface,
312 Dune::ForwardCommunication);
313 }
314 dt = problem().gridView().comm().min(dt);
315#endif
316
317 if(impet)
318 {
319 Dune::dinfo << "Timestep restricted by CellIdx " << restrictingCell << " leads to dt = "
320 <<dt * getParam<Scalar>("Impet.CFLFactor")<< std::endl;
321 if(this->averagedFaces_ != 0)
322 Dune::dwarn << " Averageing done for " << this->averagedFaces_ << " faces. "<< std::endl;
323 }
324 return;
325}
326
344template<class TypeTag>
345void FV3dTransport2P2CAdaptive<TypeTag>::getMpfaFlux(Dune::FieldVector<Scalar, 2>& fluxEntries,
346 Dune::FieldVector<Scalar, 2>& timestepFlux,
347 const IntersectionIterator& isIt, CellData& cellDataI)
348{
349 const auto& intersection = *isIt;
350
351 fluxEntries = 0.;
352 timestepFlux = 0.;
353 // cell information
354 auto elementI = intersection.inside();
355 int globalIdxI = problem().variables().index(elementI);
356
357 // get position
358 const GlobalPosition globalPos = elementI.geometry().center();
359 const GlobalPosition& gravity_ = problem().gravity();
360 // cell volume, assume linear map here
361 Scalar volume = elementI.geometry().volume();
362
363 // get values of cell I
364 Scalar pressI = problem().pressureModel().pressure(globalIdxI);
365 Scalar pcI = cellDataI.capillaryPressure();
366
367 PhaseVector SmobI(0.);
368 using std::max;
369 SmobI[wPhaseIdx] = max((cellDataI.saturation(wPhaseIdx)
370 - problem().spatialParams().materialLawParams(elementI).swr())
371 , 1e-2);
372 SmobI[nPhaseIdx] = max((cellDataI.saturation(nPhaseIdx)
373 - problem().spatialParams().materialLawParams(elementI).snr())
374 , 1e-2);
375
376 Scalar densityWI (0.), densityNWI(0.);
377 densityWI= cellDataI.density(wPhaseIdx);
378 densityNWI = cellDataI.density(nPhaseIdx);
379
380 PhaseVector potential(0.);
381
382 // access neighbor
383 auto neighbor = intersection.outside();
384 int globalIdxJ = problem().variables().index(neighbor);
385 CellData& cellDataJ = problem().variables().cellData(globalIdxJ);
386
387 // neighbor cell center in global coordinates
388 const GlobalPosition& globalPosNeighbor = neighbor.geometry().center();
389
390 // distance vector between barycenters
391 GlobalPosition distVec = globalPosNeighbor - globalPos;
392 // compute distance between cell centers
393 Scalar dist = distVec.two_norm();
394
395 GlobalPosition unitDistVec(distVec);
396 unitDistVec /= dist;
397
398 // phase densities in neighbor
399 Scalar densityWJ (0.), densityNWJ(0.);
400 densityWJ = cellDataJ.density(wPhaseIdx);
401 densityNWJ = cellDataJ.density(nPhaseIdx);
402
403 // average phase densities with central weighting
404 double densityW_mean = (densityWI + densityWJ) * 0.5;
405 double densityNW_mean = (densityNWI + densityNWJ) * 0.5;
406
407 double pressJ = problem().pressureModel().pressure(globalIdxJ);
408 Scalar pcJ = cellDataJ.capillaryPressure();
409
410 // determine potentials for upwind
412 GlobalPosition globalPos3(0.);
413 int globalIdx3=-1;
414 GlobalPosition globalPos4(0.);
415 int globalIdx4=-1;
416 TransmissivityMatrix T(0.);
417 TransmissivityMatrix additionalT(0.);
418 // prepare second interaction region
419 GlobalPosition globalPosAdditional3(0.);
420 int globalIdxAdditional3=-1;
421 GlobalPosition globalPosAdditional4(0.);
422 int globalIdxAdditional4=-1;
423
424 int halfedgesStored
425 = problem().variables().getMpfaData3D(intersection, T, globalPos3, globalIdx3, globalPos4, globalIdx4 );
426 if (halfedgesStored == 0)
427 halfedgesStored = problem().pressureModel().computeTransmissibilities(isIt,T,
428 globalPos3, globalIdx3, globalPos4, globalIdx4 );
429
430 // acess cell 3&4 and prepare mpfa
431 Scalar press3 = problem().pressureModel().pressure(globalIdx3);
432 CellData& cellData3 = problem().variables().cellData(globalIdx3);
433 Scalar pc3 = cellData3.capillaryPressure();
434 Scalar press4 = problem().pressureModel().pressure(globalIdx4);
435 CellData& cellData4 = problem().variables().cellData(globalIdx4);
436 Scalar pc4 = cellData4.capillaryPressure();
437 Scalar temp1 = globalPos * gravity_;
438 Scalar temp2 = globalPosNeighbor * gravity_;
439 Scalar temp3 = globalPos3 * gravity_;
440 Scalar temp4 = globalPos4 * gravity_;
441 if(pressureType==pw)
442 {
443 potential[wPhaseIdx] += (pressI-temp1*densityW_mean) * T[0]
444 +(pressJ-temp2*densityW_mean) * T[1]
445 +(press3- temp3*densityW_mean) * T[2]
446 +(press4- temp4*densityW_mean) * T[3];
447 potential[nPhaseIdx] += (pressI+pcI-temp1*densityNW_mean) * T[0]
448 +(pressJ+pcJ-temp2*densityNW_mean) * T[1]
449 +(press3+pc3- temp3*densityNW_mean) * T[2]
450 +(press4+pc4- temp4*densityNW_mean) * T[3];
451 }
452 else if(pressureType==pn)
453 {
454 potential[wPhaseIdx] += (pressI-pcI-temp1*densityW_mean) * T[0]
455 + (pressJ-pcJ-temp2*densityW_mean) * T[1]
456 + (press3-pc3- temp3*densityW_mean) * T[2]
457 + (press4-pc4- temp4*densityW_mean) * T[3];
458 potential[nPhaseIdx] += (pressI-temp1*densityNW_mean) * T[0]
459 + (pressJ-temp2*densityNW_mean) * T[1]
460 + (press3-temp3*densityNW_mean) * T[2]
461 + (press4-temp4*densityNW_mean) * T[3];
462 }
463 // regard more interaction regions, if there are more
464 if(halfedgesStored != 1)
465 {
466 for(int banana = 1; banana < halfedgesStored; banana ++)
467 {
468 // get data for second interaction region
469 problem().variables().getMpfaData3D(intersection, additionalT,
470 globalPosAdditional3, globalIdxAdditional3,
471 globalPosAdditional4, globalIdxAdditional4 ,
472 banana); // offset for second interaction region
473
474 Scalar gravityContributionAdditonal
475 = temp1 * additionalT[0] + temp2 * additionalT[1]
476 + globalPosAdditional3*gravity_ * additionalT[2]
477 + globalPosAdditional4*gravity_ * additionalT[3];
478 CellData& cellDataA3 = problem().variables().cellData(globalIdxAdditional3);
479 CellData& cellDataA4 = problem().variables().cellData(globalIdxAdditional4);
480
481 if(pressureType==pw)
482 {
483 potential[wPhaseIdx] += pressI * additionalT[0] + pressJ * additionalT[1]
484 +problem().pressureModel().pressure(globalIdxAdditional3) * additionalT[2]
485 +problem().pressureModel().pressure(globalIdxAdditional4)* additionalT[3];
486 potential[nPhaseIdx] += (pressI+pcI) * additionalT[0] + (pressJ+pcJ) * additionalT[1]
487 +(problem().pressureModel().pressure(globalIdxAdditional3)+cellDataA3.capillaryPressure()) * additionalT[2]
488 +(problem().pressureModel().pressure(globalIdxAdditional4)+cellDataA4.capillaryPressure()) * additionalT[3];
489 }
490 else if(pressureType==pn)
491 {
492 potential[wPhaseIdx] += (pressI-pcI) * additionalT[0] + (pressJ-pcJ) * additionalT[1]
493 + (problem().pressureModel().pressure(globalIdxAdditional3)-cellDataA3.capillaryPressure()) * additionalT[2]
494 + (problem().pressureModel().pressure(globalIdxAdditional4)-cellDataA4.capillaryPressure()) * additionalT[3];
495 potential[nPhaseIdx] += pressI * additionalT[0] + pressJ * additionalT[1]
496 + problem().pressureModel().pressure(globalIdxAdditional3) * additionalT[2]
497 + problem().pressureModel().pressure(globalIdxAdditional4) * additionalT[3];
498 }
499 potential[wPhaseIdx] -= gravityContributionAdditonal * densityW_mean;
500 potential[nPhaseIdx] -= gravityContributionAdditonal * densityNW_mean;
501 }
502 } //end of mpfa specific stuff
503
504 // determine upwinding direction, perform upwinding if possible
505 Dune::FieldVector<bool, NumPhases> doUpwinding(true);
506 PhaseVector lambda(0.);
507 for(int phaseIdx = 0; phaseIdx < NumPhases; phaseIdx++)
508 {
509 int contiEqIdx = 0;
510 if(phaseIdx == wPhaseIdx)
511 contiEqIdx = contiWEqIdx;
512 else
513 contiEqIdx = contiNEqIdx;
514
515 if(!this->impet_ or !this->restrictFluxInTransport_) // perform a strict uwpind scheme
516 {
517 if(potential[phaseIdx] > 0.)
518 {
519 lambda[phaseIdx] = cellDataI.mobility(phaseIdx);
520 cellDataI.setUpwindCell(intersection.indexInInside(), contiEqIdx, true);
521 }
522 else if(potential[phaseIdx] < 0.)
523 {
524 lambda[phaseIdx] = cellDataJ.mobility(phaseIdx);
525 cellDataI.setUpwindCell(intersection.indexInInside(), contiEqIdx, false);
526 }
527 else
528 {
529 doUpwinding[phaseIdx] = false;
530 cellDataI.setUpwindCell(intersection.indexInInside(), contiEqIdx, false);
531 cellDataJ.setUpwindCell(intersection.indexInOutside(), contiEqIdx, false);
532 }
533 }
534 else // Transport after PE with check on flow direction
535 {
536 bool cellIwasUpwindCell;
537 //get the information from smaller (higher level) cell, as its IS is unique
538 if(elementI.level()>neighbor.level())
539 cellIwasUpwindCell = cellDataI.isUpwindCell(intersection.indexInInside(), contiEqIdx);
540 else // reverse neighbors information gathered
541 cellIwasUpwindCell = !cellDataJ.isUpwindCell(intersection.indexInOutside(), contiEqIdx);
542
543 if (potential[phaseIdx] > 0. && cellIwasUpwindCell)
544 lambda[phaseIdx] = cellDataI.mobility(phaseIdx);
545 else if (potential[phaseIdx] < 0. && !cellIwasUpwindCell)
546 lambda[phaseIdx] = cellDataJ.mobility(phaseIdx);
547 // potential direction does not coincide with that of P.E.
548 else if(this->restrictFluxInTransport_ == 2) // perform central averaging for all direction changes
549 doUpwinding[phaseIdx] = false;
550 else // i.e. restrictFluxInTransport == 1
551 {
552 //check if harmonic weithing is necessary
553 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
554 lambda[phaseIdx] = cellDataI.mobility(phaseIdx);
555 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
556 lambda[phaseIdx] = cellDataJ.mobility(phaseIdx);
557 else
558 doUpwinding[phaseIdx] = false;
559 }
560
561 // do not perform upwinding if so desired
562 if(!doUpwinding[phaseIdx])
563 {
564 //a) no flux if there wont be any flux regardless how to average/upwind
565 if(cellDataI.mobility(phaseIdx)+cellDataJ.mobility(phaseIdx)==0.)
566 {
567 potential[phaseIdx] = 0;
568 continue;
569 }
570
571 //b) perform harmonic averaging
572 fluxEntries[wCompIdx] -= potential[phaseIdx] / volume
573 * harmonicMean(cellDataI.massFraction(phaseIdx, wCompIdx) * cellDataI.mobility(phaseIdx) * cellDataI.density(phaseIdx),
574 cellDataJ.massFraction(phaseIdx, wCompIdx) * cellDataJ.mobility(phaseIdx) * cellDataJ.density(phaseIdx));
575 fluxEntries[nCompIdx] -= potential[phaseIdx] / volume
576 * harmonicMean(cellDataI.massFraction(phaseIdx, nCompIdx) * cellDataI.mobility(phaseIdx) * cellDataI.density(phaseIdx),
577 cellDataJ.massFraction(phaseIdx, nCompIdx) * cellDataJ.mobility(phaseIdx) * cellDataJ.density(phaseIdx));
578 // c) timestep control
579 // for timestep control : influx
580 using std::max;
581 timestepFlux[0] += max(0.,
582 - potential[phaseIdx] / volume
583 * harmonicMean(cellDataI.mobility(phaseIdx),cellDataJ.mobility(phaseIdx)));
584 // outflux
585 timestepFlux[1] += max(0.,
586 potential[phaseIdx] / volume
587 * harmonicMean(cellDataI.mobility(phaseIdx),cellDataJ.mobility(phaseIdx))/SmobI[phaseIdx]);
588
589 //d) output (only for one side)
590 this->averagedFaces_++;
591 #if DUNE_MINIMAL_DEBUG_LEVEL < 3
592 // verbose (only for one side)
593 if(globalIdxI > globalIdxJ)
594 Dune::dinfo << "harmonicMean flux of phase" << phaseIdx <<" used from cell" << globalIdxI<< " into " << globalIdxJ
595 << " ; TE upwind I = "<< cellIwasUpwindCell << " but pot = "<< potential[phaseIdx] << std::endl;
596 #endif
597
598 //e) stop further standard calculations
599 potential[phaseIdx] = 0;
600 }
601 }
602 }
603
604 // calculate and standardized velocity
605 using std::max;
606 double velocityJIw = max((-lambda[wPhaseIdx] * potential[wPhaseIdx]) / volume, 0.0);
607 double velocityIJw = max(( lambda[wPhaseIdx] * potential[wPhaseIdx]) / volume, 0.0);
608 double velocityJIn = max((-lambda[nPhaseIdx] * potential[nPhaseIdx]) / volume, 0.0);
609 double velocityIJn = max(( lambda[nPhaseIdx] * potential[nPhaseIdx]) / volume, 0.0);
610
611 // for timestep control
612 timestepFlux[0] += velocityJIw + velocityJIn;
613
614 double foutw = velocityIJw/SmobI[wPhaseIdx];
615 double foutn = velocityIJn/SmobI[nPhaseIdx];
616 using std::isnan;
617 using std::isinf;
618 if (isnan(foutw) || isinf(foutw) || foutw < 0) foutw = 0;
619 if (isnan(foutn) || isinf(foutn) || foutn < 0) foutn = 0;
620 timestepFlux[1] += foutw + foutn;
621
622 fluxEntries[wCompIdx] +=
623 velocityJIw * cellDataJ.massFraction(wPhaseIdx, wCompIdx) * densityWJ
624 - velocityIJw * cellDataI.massFraction(wPhaseIdx, wCompIdx) * densityWI
625 + velocityJIn * cellDataJ.massFraction(nPhaseIdx, wCompIdx) * densityNWJ
626 - velocityIJn * cellDataI.massFraction(nPhaseIdx, wCompIdx) * densityNWI;
627 fluxEntries[nCompIdx] +=
628 velocityJIw * cellDataJ.massFraction(wPhaseIdx, nCompIdx) * densityWJ
629 - velocityIJw * cellDataI.massFraction(wPhaseIdx, nCompIdx) * densityWI
630 + velocityJIn * cellDataJ.massFraction(nPhaseIdx, nCompIdx) * densityNWJ
631 - velocityIJn * cellDataI.massFraction(nPhaseIdx, nCompIdx) * densityNWI;
632 return;
633}
634
635} // end namespace Dumux
636#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.
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:79
Compositional transport step in a finite volume discretization.
Definition: fv3dtransportadaptive.hh:54
virtual void update(const Scalar t, Scalar &dt, TransportSolutionType &updateVec, bool impes=false)
Calculate the update vector and determine timestep size.
Definition: fv3dtransportadaptive.hh:162
static const int pressureType
‍Specifies if the MPFA is used on hanging nodes
Definition: fv3dtransportadaptive.hh:138
FV3dTransport2P2CAdaptive(Problem &problem)
Constructs a FV3dTransport2P2CAdaptive object.
Definition: fv3dtransportadaptive.hh:122
Problem & problem_
Definition: fv3dtransportadaptive.hh:133
void getMpfaFlux(Dune::FieldVector< Scalar, 2 > &, Dune::FieldVector< Scalar, 2 > &, const IntersectionIterator &, CellData &)
Compute flux over an irregular interface using a mpfa method.
Definition: fv3dtransportadaptive.hh:345
virtual ~FV3dTransport2P2CAdaptive()
Definition: fv3dtransportadaptive.hh:128
bool enableMPFA
Definition: fv3dtransportadaptive.hh:135
Compositional transport step in a Finite Volume discretization.
Definition: fvtransport.hh:60
Definition: fvtransport.hh:113
Dune::FieldVector< EntryType, 2 *dim > faceFluxes
Definition: fvtransport.hh:114
Dune::FieldVector< Scalar, 2 *dim > faceTargetDt
Definition: fvtransport.hh:115