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