22#ifndef DUMUX_MPFAL2DPRESSUREVELOCITIES2P_ADAPTIVE_HH
23#define DUMUX_MPFAL2DPRESSUREVELOCITIES2P_ADAPTIVE_HH
25#include <dune/common/float_cmp.hh>
55 dim = GridView::dimension, dimWorld = GridView::dimensionworld
63 using PrimaryVariables =
typename SolutionTypes::PrimaryVariables;
65 using ReferenceElements = Dune::ReferenceElements<Scalar, dim>;
68 using MaterialLaw =
typename SpatialParams::MaterialLaw;
75 using Grid =
typename GridView::Grid;
76 using IndexSet =
typename GridView::IndexSet;
77 using Intersection =
typename GridView::Intersection;
79 using Geometry =
typename Grid::template Codim<0>::Entity::Geometry;
80 using JacobianTransposed =
typename Geometry::JacobianTransposed ;
86 pw = Indices::pressureW,
87 pn = Indices::pressureNw,
88 sw = Indices::saturationW,
89 sn = Indices::saturationNw
93 wPhaseIdx = Indices::wPhaseIdx,
94 nPhaseIdx = Indices::nPhaseIdx,
95 pressureIdx = Indices::pressureIdx,
96 saturationIdx = Indices::saturationIdx,
97 pressEqIdx = Indices::pressureEqIdx,
98 satEqIdx = Indices::satEqIdx,
99 numPhases = getPropValue<TypeTag, Properties::NumPhases>()
102 using Element =
typename GridView::template Codim<0>::Entity;
104 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
105 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
106 using DimVector = Dune::FieldVector<Scalar, dim>;
115 ParentType(problem), problem_(problem), velocity_(problem)
117 density_[wPhaseIdx] = 0.;
118 density_[nPhaseIdx] = 0.;
119 viscosity_[wPhaseIdx] = 0.;
120 viscosity_[nPhaseIdx] = 0.;
122 calcVelocityInTransport_ = getParam<bool>(
"MPFA.CalcVelocityInTransport");
150 const auto element = *problem_.gridView().template begin<0>();
151 FluidState fluidState;
152 fluidState.setPressure(wPhaseIdx, problem_.referencePressure(element));
153 fluidState.setPressure(nPhaseIdx, problem_.referencePressure(element));
154 fluidState.setTemperature(problem_.temperature(element));
155 fluidState.setSaturation(wPhaseIdx, 1.);
156 fluidState.setSaturation(nPhaseIdx, 0.);
163 velocity_.initialize();
189 return calcVelocityInTransport_;
202 template<
class MultiWriter>
206 velocity_.addOutputVtkFields(writer);
213 Scalar density_[numPhases];
214 Scalar viscosity_[numPhases];
215 bool calcVelocityInTransport_;
218 static const int pressureType_ = getPropValue<TypeTag, Properties::PressureFormulation>();
220 static const int saturationType_ = getPropValue<TypeTag, Properties::SaturationFormulation>();
229template<
class TypeTag>
233 for (
const auto& vertex : vertices(problem_.gridView()))
235 int vIdxGlobal = problem_.variables().index(vertex);
249 int eIdxGlobal1 = problem_.variables().index(element1);
250 int eIdxGlobal2 = problem_.variables().index(element2);
251 int eIdxGlobal3 = problem_.variables().index(element3);
252 int eIdxGlobal4 = problem_.variables().index(element4);
255 CellData& cellData1 = problem_.variables().cellData(eIdxGlobal1);
256 CellData& cellData2 = problem_.variables().cellData(eIdxGlobal2);
257 CellData& cellData3 = problem_.variables().cellData(eIdxGlobal3);
258 CellData& cellData4 = problem_.variables().cellData(eIdxGlobal4);
260 velocity_.calculateInnerInteractionVolumeVelocity(interactionVolume, cellData1, cellData2, cellData3,
261 cellData4, this->innerBoundaryVolumeFaces_);
271 int eIdxGlobal1 = problem_.variables().index(element1);
272 int eIdxGlobal2 = problem_.variables().index(element2);
273 int eIdxGlobal4 = problem_.variables().index(element4);
276 CellData& cellData1 = problem_.variables().cellData(eIdxGlobal1);
277 CellData& cellData2 = problem_.variables().cellData(eIdxGlobal2);
278 CellData& cellData4 = problem_.variables().cellData(eIdxGlobal4);
280 velocity_.calculateHangingNodeInteractionVolumeVelocity(interactionVolume, cellData1, cellData2,
281 cellData4, this->innerBoundaryVolumeFaces_);
286 DUNE_THROW(Dune::NotImplemented,
"Unknown interactionvolume type!");
293 for (
int elemIdx = 0; elemIdx < 2 * dim; elemIdx++)
295 bool isOutside =
false;
296 for (
int fIdx = 0; fIdx < dim; fIdx++)
310 int eIdxGlobal = problem_.variables().index(interactionVolume.
getSubVolumeElement(elemIdx));
312 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
314 velocity_.calculateBoundaryInteractionVolumeVelocity(interactionVolume, cellData, elemIdx);
331template<
class TypeTag>
334 int numVertices = intersection.geometry().corners();
336 auto elementI = intersection.inside();
337 auto elementJ = intersection.outside();
339 int levelI = elementI.level();
340 int levelJ = elementJ.level();
342 int eIdxGlobalI = problem_.variables().index(elementI);
343 int eIdxGlobalJ = problem_.variables().index(elementJ);
345 CellData& cellDataJ = problem_.variables().cellData(eIdxGlobalJ);
347 const auto referenceElement = ReferenceElements::general(elementI.type());
349 int indexInInside = intersection.indexInInside();
350 int indexInOutside = intersection.indexInOutside();
352 int fIdx = indexInInside;
355 fIdx = indexInOutside;
357 std::vector<CellData> cellDataTemp(0);
359 if (levelI != levelJ)
362 cellData.fluxData().setVelocity(wPhaseIdx, indexInInside, vel);
363 cellData.fluxData().setVelocity(nPhaseIdx, indexInInside, vel);
364 cellData.fluxData().setUpwindPotential(wPhaseIdx, indexInInside, 0);
365 cellData.fluxData().setUpwindPotential(nPhaseIdx, indexInInside, 0);
367 cellDataJ.fluxData().setVelocity(wPhaseIdx, indexInOutside, vel);
368 cellDataJ.fluxData().setVelocity(nPhaseIdx, indexInOutside, vel);
369 cellDataJ.fluxData().setUpwindPotential(wPhaseIdx, indexInOutside, 0);
370 cellDataJ.fluxData().setUpwindPotential(nPhaseIdx, indexInOutside, 0);
373 for (
int vIdx = 0; vIdx < numVertices; vIdx++)
375 int localVertIdx = referenceElement.subEntity(fIdx, dim - 1, vIdx, dim);
378 if (levelI >= levelJ)
380 vIdxGlobal = problem_.variables().index(elementI.template subEntity<dim>(localVertIdx));
384 vIdxGlobal = problem_.variables().index(elementJ.template subEntity<dim>(localVertIdx));
392 std::vector<int> eIdxGlobal(0);
401 eIdxGlobal.resize(4);
403 eIdxGlobal[0] = problem_.variables().index(element1);
404 eIdxGlobal[1] = problem_.variables().index(element2);
405 eIdxGlobal[2] = problem_.variables().index(element3);
406 eIdxGlobal[3] = problem_.variables().index(element4);
409 cellDataTemp.resize(4);
411 cellDataTemp[0] = problem_.variables().cellData(eIdxGlobal[0]);
412 cellDataTemp[1] = problem_.variables().cellData(eIdxGlobal[1]);
413 cellDataTemp[2] = problem_.variables().cellData(eIdxGlobal[2]);
414 cellDataTemp[3] = problem_.variables().cellData(eIdxGlobal[3]);
416 velocity_.calculateInnerInteractionVolumeVelocity(interactionVolume, cellDataTemp[0], cellDataTemp[1],
417 cellDataTemp[2], cellDataTemp[3], this->innerBoundaryVolumeFaces_);
425 eIdxGlobal.resize(3);
427 eIdxGlobal[0] = problem_.variables().index(element1);
428 eIdxGlobal[1] = problem_.variables().index(element2);
429 eIdxGlobal[2] = problem_.variables().index(element4);
432 cellDataTemp.resize(3);
434 cellDataTemp[0] = problem_.variables().cellData(eIdxGlobal[0]);
435 cellDataTemp[1] = problem_.variables().cellData(eIdxGlobal[1]);
436 cellDataTemp[2] = problem_.variables().cellData(eIdxGlobal[2]);
439 velocity_.calculateHangingNodeInteractionVolumeVelocity(interactionVolume, cellDataTemp[0], cellDataTemp[1],
440 cellDataTemp[2], this->innerBoundaryVolumeFaces_);
444 DUNE_THROW(Dune::NotImplemented,
"Unknown interactionvolume type!");
447 int size = cellDataTemp.size();
448 for (
int i = 0; i < size; i++)
450 if (eIdxGlobal[i] == eIdxGlobalI)
452 if (levelI >= levelJ)
454 cellData.fluxData().setVelocity(wPhaseIdx, indexInInside,
455 cellDataTemp[i].fluxData().velocity(wPhaseIdx, indexInInside));
456 cellData.fluxData().setVelocity(nPhaseIdx, indexInInside,
457 cellDataTemp[i].fluxData().velocity(nPhaseIdx, indexInInside));
458 cellData.fluxData().setUpwindPotential(wPhaseIdx, indexInInside,
459 cellDataTemp[i].fluxData().upwindPotential(wPhaseIdx, indexInInside));
460 cellData.fluxData().setUpwindPotential(nPhaseIdx, indexInInside,
461 cellDataTemp[i].fluxData().upwindPotential(nPhaseIdx, indexInInside));
465 cellDataJ.fluxData().setVelocity(wPhaseIdx, indexInOutside,
466 cellDataTemp[i].fluxData().velocity(wPhaseIdx, indexInInside));
467 cellDataJ.fluxData().setVelocity(nPhaseIdx, indexInOutside,
468 cellDataTemp[i].fluxData().velocity(nPhaseIdx, indexInInside));
469 cellDataJ.fluxData().setUpwindPotential(wPhaseIdx, indexInOutside,
470 cellDataTemp[i].fluxData().upwindPotential(wPhaseIdx, indexInInside));
471 cellDataJ.fluxData().setUpwindPotential(nPhaseIdx, indexInOutside,
472 cellDataTemp[i].fluxData().upwindPotential(nPhaseIdx, indexInInside));
478 else if (eIdxGlobal[i] == eIdxGlobalJ)
480 if (levelJ >= levelI)
482 cellDataJ.fluxData().setVelocity(wPhaseIdx, indexInOutside,
483 cellDataTemp[i].fluxData().velocity(wPhaseIdx, indexInOutside));
484 cellDataJ.fluxData().setVelocity(nPhaseIdx, indexInOutside,
485 cellDataTemp[i].fluxData().velocity(nPhaseIdx, indexInOutside));
486 cellDataJ.fluxData().setUpwindPotential(wPhaseIdx, indexInOutside,
487 cellDataTemp[i].fluxData().upwindPotential(wPhaseIdx, indexInOutside));
488 cellDataJ.fluxData().setUpwindPotential(nPhaseIdx, indexInOutside,
489 cellDataTemp[i].fluxData().upwindPotential(nPhaseIdx, indexInOutside));
493 cellData.fluxData().setVelocity(wPhaseIdx, indexInInside,
494 cellDataTemp[i].fluxData().velocity(wPhaseIdx, indexInOutside));
495 cellData.fluxData().setVelocity(nPhaseIdx, indexInInside,
496 cellDataTemp[i].fluxData().velocity(nPhaseIdx, indexInOutside));
497 cellData.fluxData().setUpwindPotential(wPhaseIdx, indexInInside,
498 cellDataTemp[i].fluxData().upwindPotential(wPhaseIdx, indexInOutside));
499 cellData.fluxData().setUpwindPotential(nPhaseIdx, indexInInside,
500 cellDataTemp[i].fluxData().upwindPotential(nPhaseIdx, indexInOutside));
508 if (levelI == levelJ)
510 cellData.fluxData().setVelocityMarker(indexInInside);
511 cellDataJ.fluxData().setVelocityMarker(indexInOutside);
523template<
class TypeTag>
526 auto element = intersection.inside();
529 int isIndex = intersection.indexInInside();
532 const Dune::FieldVector<Scalar, dim>& unitOuterNormal = intersection.centerUnitOuterNormal();
534 BoundaryTypes bcType;
536 problem_.boundaryTypes(bcType, intersection);
537 PrimaryVariables boundValues(0.0);
539 if (bcType.isDirichlet(pressEqIdx))
541 problem_.dirichlet(boundValues, intersection);
544 const GlobalPosition& globalPosI = element.geometry().center();
547 const GlobalPosition& globalPosJ = intersection.geometry().center();
550 Scalar lambdaWI = cellData.mobility(wPhaseIdx);
551 Scalar lambdaNwI = cellData.mobility(nPhaseIdx);
554 Scalar pcI = cellData.capillaryPressure();
557 GlobalPosition distVec = globalPosJ - globalPosI;
560 Scalar dist = distVec.two_norm();
564 DimMatrix meanPermeability(0);
566 problem_.spatialParams().meanK(meanPermeability, problem_.spatialParams().intrinsicPermeability(element));
573 if (bcType.isDirichlet(satEqIdx))
575 switch (saturationType_)
579 satW = boundValues[saturationIdx];
584 satW = 1 - boundValues[saturationIdx];
591 satW = cellData.saturation(wPhaseIdx);
594 Scalar pressBound = boundValues[pressureIdx];
595 Scalar pcBound = MaterialLaw::pc(problem_.spatialParams().materialLawParams(element), satW);
598 Scalar pressWBound = 0;
599 Scalar pressNwBound = 0;
600 if (pressureType_ == pw)
602 pressWBound = pressBound;
603 pressNwBound = pressBound + pcBound;
605 else if (pressureType_ == pn)
607 pressWBound = pressBound - pcBound;
608 pressNwBound = pressBound;
611 Scalar lambdaWBound = MaterialLaw::krw(problem_.spatialParams().materialLawParams(element), satW)
612 / viscosity_[wPhaseIdx];
613 Scalar lambdaNwBound = MaterialLaw::krn(problem_.spatialParams().materialLawParams(element), satW)
614 / viscosity_[nPhaseIdx];
616 Scalar potentialDiffW = cellData.fluxData().upwindPotential(wPhaseIdx, isIndex);
617 Scalar potentialDiffNw = cellData.fluxData().upwindPotential(nPhaseIdx, isIndex);
620 potentialDiffW = (cellData.pressure(wPhaseIdx) - pressWBound);
621 potentialDiffNw = (cellData.pressure(nPhaseIdx) - pressNwBound);
623 potentialDiffW += density_[wPhaseIdx] * (distVec * problem_.gravity());
624 potentialDiffNw += density_[nPhaseIdx] * (distVec * problem_.gravity());
627 cellData.fluxData().setUpwindPotential(wPhaseIdx, isIndex, potentialDiffW);
628 cellData.fluxData().setUpwindPotential(nPhaseIdx, isIndex, potentialDiffNw);
631 Scalar lambdaW = (potentialDiffW > 0.) ? lambdaWI : lambdaWBound;
632 lambdaW = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffW, 0.0, 1.0e-30)) ? 0.5 * (lambdaWI + lambdaWBound) : lambdaW;
633 Scalar lambdaNw = (potentialDiffNw > 0.) ? lambdaNwI : lambdaNwBound;
634 lambdaNw = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffNw, 0.0, 1.0e-30)) ? 0.5 * (lambdaNwI + lambdaNwBound) : lambdaNw;
640 Dune::FieldVector<Scalar, dimWorld> velocityW(unitOuterNormal);
641 Dune::FieldVector<Scalar, dimWorld> velocityNw(unitOuterNormal);
645 Scalar areaScaling = (unitOuterNormal * distVec);
648 Scalar gravityTermW = (problem_.gravity() * distVec) * density_[wPhaseIdx] * areaScaling;
649 Scalar gravityTermNw = (problem_.gravity() * distVec) * density_[nPhaseIdx] * areaScaling;
652 switch (pressureType_)
656 velocityW *= lambdaW * scalarPerm * ((cellData.pressure(wPhaseIdx) - pressBound) / dist + gravityTermW);
657 velocityNw *= lambdaNw * scalarPerm * ((cellData.pressure(wPhaseIdx) - pressBound) / dist + gravityTermNw)
658 + 0.5 * (lambdaNwI + lambdaNwBound) * scalarPerm * (pcI - pcBound) / dist;
663 velocityW *= lambdaW * scalarPerm * ((cellData.pressure(nPhaseIdx) - pressBound) / dist + gravityTermW)
664 - 0.5 * (lambdaWI + lambdaWBound) * scalarPerm * (pcI - pcBound) / dist;
665 velocityNw *= lambdaNw * scalarPerm * ((cellData.pressure(nPhaseIdx) - pressBound) / dist + gravityTermNw);
671 cellData.fluxData().setVelocity(wPhaseIdx, isIndex, velocityW);
672 cellData.fluxData().setVelocity(nPhaseIdx, isIndex, velocityNw);
673 cellData.fluxData().setVelocityMarker(isIndex);
677 else if (bcType.isNeumann(pressEqIdx))
679 problem_.neumann(boundValues, intersection);
681 Dune::FieldVector<Scalar, dimWorld> velocityW(unitOuterNormal);
682 Dune::FieldVector<Scalar, dimWorld> velocityNw(unitOuterNormal);
684 velocityW *= boundValues[wPhaseIdx];
685 velocityNw *= boundValues[nPhaseIdx];
687 velocityW /= density_[wPhaseIdx];
688 velocityNw /= density_[nPhaseIdx];
691 cellData.fluxData().setUpwindPotential(wPhaseIdx, isIndex, boundValues[wPhaseIdx]);
692 cellData.fluxData().setUpwindPotential(nPhaseIdx, isIndex, boundValues[nPhaseIdx]);
694 cellData.fluxData().setVelocity(wPhaseIdx, isIndex, velocityW);
695 cellData.fluxData().setVelocity(nPhaseIdx, isIndex, velocityNw);
696 cellData.fluxData().setVelocityMarker(isIndex);
700 DUNE_THROW(Dune::NotImplemented,
"No valid boundary condition type defined for pressure equation!");
Grid adaptive finite volume MPFA L-method discretization of a two-phase pressure equation of the sequ...
Velocity calculation on adaptive grids using a 2-d MPFA L-method.
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
std::string viscosity(int phaseIdx) noexcept
I/O name of viscosity for multiphase systems.
Definition: name.hh:74
std::string permeability() noexcept
I/O name of permeability.
Definition: name.hh:143
std::string density(int phaseIdx) noexcept
I/O name of density for multiphase systems.
Definition: name.hh:65
Grid adaptive finite volume MPFA L-method discretization of a two-phase flow pressure equation of the...
Definition: 2dpressureadaptive.hh:75
FVMPFALInteractionVolume< TypeTag > InteractionVolume
Type of the interaction volume objects.
Definition: 2dpressureadaptive.hh:155
void addOutputVtkFields(MultiWriter &writer)
Adds pressure output to the output file.
Definition: 2dpressureadaptive.hh:377
void storePressureSolution()
Globally stores the pressure solution.
Definition: 2dpressureadaptive.hh:234
void initialize()
Initializes the pressure model.
Definition: 2dpressureadaptive.hh:203
void updateMaterialLaws()
constitutive functions are initialized and stored in the variables object
Definition: 2dpressureadaptive.hh:2552
void update()
Pressure update.
Definition: 2dpressureadaptive.hh:301
Class for the calculation of velocities from the pressure solution of an IMPES scheme using a grid ad...
Definition: 2dpressurevelocityadaptive.hh:48
void addOutputVtkFields(MultiWriter &writer)
Adds velocity output to the output file.
Definition: 2dpressurevelocityadaptive.hh:203
void updateVelocity()
Function for updating the velocity field if iterations are necessary in the transport solution.
Definition: 2dpressurevelocityadaptive.hh:133
bool calculateVelocityInTransport()
Indicates if velocity is reconstructed in the pressure step or in the transport step.
Definition: 2dpressurevelocityadaptive.hh:187
FvMpfaL2dPressureVelocity2pAdaptive(Problem &problem)
Constructs a FvMpfaL2dPressureVelocity2pAdaptive object.
Definition: 2dpressurevelocityadaptive.hh:114
void initialize()
Initializes pressure and velocity.
Definition: 2dpressurevelocityadaptive.hh:148
void calculateVelocityOnBoundary(const Intersection &intersection, CellData &cellData)
Calculates the velocity at a boundary.
Definition: 2dpressurevelocityadaptive.hh:524
void calculateVelocity()
Calculates the velocities at all cell-cell interfaces.
Definition: 2dpressurevelocityadaptive.hh:230
void update()
Pressure and velocity update.
Definition: 2dpressurevelocityadaptive.hh:175
Class for calculating 2-d velocities from cell-wise constant pressure values.
Definition: 2dvelocityadaptive.hh:54
Class including the information of an interaction volume of a MPFA L-method that does not change with...
Definition: linteractionvolume.hh:41
bool isInnerVolume()
Returns true if the interaction volume is completely inside the model domain.
Definition: linteractionvolume.hh:316
int getElementNumber()
Get number of stored elements.
Definition: linteractionvolume.hh:246
Element getSubVolumeElement(int subVolumeIdx)
Get an element of the interaction volume.
Definition: linteractionvolume.hh:281
int getFaceIndexFromSubVolume(int subVolumeIdx, int subVolumeFaceIdx)
Map from local interaction volume numbering on element to numbering on interaction volume.
Definition: linteractionvolume.hh:270
bool isOutsideFace(int subVolumeFaceIdx)
Returns true if an interaction volume flux face is outside the model domain.
Definition: linteractionvolume.hh:301
The finite volume base class for the solution of a pressure equation.
Definition: sequential/cellcentered/pressure.hh:49