22#ifndef DUMUX_MPFAL2DPRESSUREVELOCITIES2P_ADAPTIVE_HH
23#define DUMUX_MPFAL2DPRESSUREVELOCITIES2P_ADAPTIVE_HH
25#include <dune/common/float_cmp.hh>
57 dim = GridView::dimension, dimWorld = GridView::dimensionworld
65 using PrimaryVariables =
typename SolutionTypes::PrimaryVariables;
74 using Grid =
typename GridView::Grid;
75 using IndexSet =
typename GridView::IndexSet;
76 using Intersection =
typename GridView::Intersection;
78 using Geometry =
typename Grid::template Codim<0>::Entity::Geometry;
79 using JacobianTransposed =
typename Geometry::JacobianTransposed ;
85 pw = Indices::pressureW,
86 pn = Indices::pressureNw,
87 sw = Indices::saturationW,
88 sn = Indices::saturationNw
92 wPhaseIdx = Indices::wPhaseIdx,
93 nPhaseIdx = Indices::nPhaseIdx,
94 pressureIdx = Indices::pressureIdx,
95 saturationIdx = Indices::saturationIdx,
96 pressEqIdx = Indices::pressureEqIdx,
97 satEqIdx = Indices::satEqIdx,
98 numPhases = getPropValue<TypeTag, Properties::NumPhases>()
101 using Element =
typename GridView::template Codim<0>::Entity;
103 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
104 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
105 using DimVector = Dune::FieldVector<Scalar, dim>;
114 ParentType(problem), problem_(problem), velocity_(problem)
116 density_[wPhaseIdx] = 0.;
117 density_[nPhaseIdx] = 0.;
118 viscosity_[wPhaseIdx] = 0.;
119 viscosity_[nPhaseIdx] = 0.;
121 calcVelocityInTransport_ = getParam<bool>(
"MPFA.CalcVelocityInTransport");
149 const auto element = *problem_.gridView().template begin<0>();
150 FluidState fluidState;
151 fluidState.setPressure(wPhaseIdx, problem_.referencePressure(element));
152 fluidState.setPressure(nPhaseIdx, problem_.referencePressure(element));
153 fluidState.setTemperature(problem_.temperature(element));
154 fluidState.setSaturation(wPhaseIdx, 1.);
155 fluidState.setSaturation(nPhaseIdx, 0.);
162 velocity_.initialize();
188 return calcVelocityInTransport_;
201 template<
class MultiWriter>
205 velocity_.addOutputVtkFields(writer);
212 Scalar density_[numPhases];
213 Scalar viscosity_[numPhases];
214 bool calcVelocityInTransport_;
217 static const int pressureType_ = getPropValue<TypeTag, Properties::PressureFormulation>();
219 static const int saturationType_ = getPropValue<TypeTag, Properties::SaturationFormulation>();
228template<
class TypeTag>
232 for (
const auto& vertex : vertices(problem_.gridView()))
234 int vIdxGlobal = problem_.variables().index(vertex);
248 int eIdxGlobal1 = problem_.variables().index(element1);
249 int eIdxGlobal2 = problem_.variables().index(element2);
250 int eIdxGlobal3 = problem_.variables().index(element3);
251 int eIdxGlobal4 = problem_.variables().index(element4);
254 CellData& cellData1 = problem_.variables().cellData(eIdxGlobal1);
255 CellData& cellData2 = problem_.variables().cellData(eIdxGlobal2);
256 CellData& cellData3 = problem_.variables().cellData(eIdxGlobal3);
257 CellData& cellData4 = problem_.variables().cellData(eIdxGlobal4);
259 velocity_.calculateInnerInteractionVolumeVelocity(interactionVolume, cellData1, cellData2, cellData3,
260 cellData4, this->innerBoundaryVolumeFaces_);
270 int eIdxGlobal1 = problem_.variables().index(element1);
271 int eIdxGlobal2 = problem_.variables().index(element2);
272 int eIdxGlobal4 = problem_.variables().index(element4);
275 CellData& cellData1 = problem_.variables().cellData(eIdxGlobal1);
276 CellData& cellData2 = problem_.variables().cellData(eIdxGlobal2);
277 CellData& cellData4 = problem_.variables().cellData(eIdxGlobal4);
279 velocity_.calculateHangingNodeInteractionVolumeVelocity(interactionVolume, cellData1, cellData2,
280 cellData4, this->innerBoundaryVolumeFaces_);
285 DUNE_THROW(Dune::NotImplemented,
"Unknown interactionvolume type!");
292 for (
int elemIdx = 0; elemIdx < 2 * dim; elemIdx++)
294 bool isOutside =
false;
295 for (
int fIdx = 0; fIdx < dim; fIdx++)
309 int eIdxGlobal = problem_.variables().index(interactionVolume.
getSubVolumeElement(elemIdx));
311 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
313 velocity_.calculateBoundaryInteractionVolumeVelocity(interactionVolume, cellData, elemIdx);
330template<
class TypeTag>
333 int numVertices = intersection.geometry().corners();
335 auto elementI = intersection.inside();
336 auto elementJ = intersection.outside();
338 int levelI = elementI.level();
339 int levelJ = elementJ.level();
341 int eIdxGlobalI = problem_.variables().index(elementI);
342 int eIdxGlobalJ = problem_.variables().index(elementJ);
344 CellData& cellDataJ = problem_.variables().cellData(eIdxGlobalJ);
346 const auto refElement = referenceElement(elementI);
348 int indexInInside = intersection.indexInInside();
349 int indexInOutside = intersection.indexInOutside();
351 int fIdx = indexInInside;
354 fIdx = indexInOutside;
356 std::vector<CellData> cellDataTemp(0);
358 if (levelI != levelJ)
361 cellData.fluxData().setVelocity(wPhaseIdx, indexInInside, vel);
362 cellData.fluxData().setVelocity(nPhaseIdx, indexInInside, vel);
363 cellData.fluxData().setUpwindPotential(wPhaseIdx, indexInInside, 0);
364 cellData.fluxData().setUpwindPotential(nPhaseIdx, indexInInside, 0);
366 cellDataJ.fluxData().setVelocity(wPhaseIdx, indexInOutside, vel);
367 cellDataJ.fluxData().setVelocity(nPhaseIdx, indexInOutside, vel);
368 cellDataJ.fluxData().setUpwindPotential(wPhaseIdx, indexInOutside, 0);
369 cellDataJ.fluxData().setUpwindPotential(nPhaseIdx, indexInOutside, 0);
372 for (
int vIdx = 0; vIdx < numVertices; vIdx++)
374 int localVertIdx = refElement.subEntity(fIdx, dim - 1, vIdx, dim);
377 if (levelI >= levelJ)
379 vIdxGlobal = problem_.variables().index(elementI.template subEntity<dim>(localVertIdx));
383 vIdxGlobal = problem_.variables().index(elementJ.template subEntity<dim>(localVertIdx));
391 std::vector<int> eIdxGlobal(0);
400 eIdxGlobal.resize(4);
402 eIdxGlobal[0] = problem_.variables().index(element1);
403 eIdxGlobal[1] = problem_.variables().index(element2);
404 eIdxGlobal[2] = problem_.variables().index(element3);
405 eIdxGlobal[3] = problem_.variables().index(element4);
408 cellDataTemp.resize(4);
410 cellDataTemp[0] = problem_.variables().cellData(eIdxGlobal[0]);
411 cellDataTemp[1] = problem_.variables().cellData(eIdxGlobal[1]);
412 cellDataTemp[2] = problem_.variables().cellData(eIdxGlobal[2]);
413 cellDataTemp[3] = problem_.variables().cellData(eIdxGlobal[3]);
415 velocity_.calculateInnerInteractionVolumeVelocity(interactionVolume, cellDataTemp[0], cellDataTemp[1],
416 cellDataTemp[2], cellDataTemp[3], this->innerBoundaryVolumeFaces_);
424 eIdxGlobal.resize(3);
426 eIdxGlobal[0] = problem_.variables().index(element1);
427 eIdxGlobal[1] = problem_.variables().index(element2);
428 eIdxGlobal[2] = problem_.variables().index(element4);
431 cellDataTemp.resize(3);
433 cellDataTemp[0] = problem_.variables().cellData(eIdxGlobal[0]);
434 cellDataTemp[1] = problem_.variables().cellData(eIdxGlobal[1]);
435 cellDataTemp[2] = problem_.variables().cellData(eIdxGlobal[2]);
438 velocity_.calculateHangingNodeInteractionVolumeVelocity(interactionVolume, cellDataTemp[0], cellDataTemp[1],
439 cellDataTemp[2], this->innerBoundaryVolumeFaces_);
443 DUNE_THROW(Dune::NotImplemented,
"Unknown interactionvolume type!");
446 int size = cellDataTemp.size();
447 for (
int i = 0; i < size; i++)
449 if (eIdxGlobal[i] == eIdxGlobalI)
451 if (levelI >= levelJ)
453 cellData.fluxData().setVelocity(wPhaseIdx, indexInInside,
454 cellDataTemp[i].fluxData().velocity(wPhaseIdx, indexInInside));
455 cellData.fluxData().setVelocity(nPhaseIdx, indexInInside,
456 cellDataTemp[i].fluxData().velocity(nPhaseIdx, indexInInside));
457 cellData.fluxData().setUpwindPotential(wPhaseIdx, indexInInside,
458 cellDataTemp[i].fluxData().upwindPotential(wPhaseIdx, indexInInside));
459 cellData.fluxData().setUpwindPotential(nPhaseIdx, indexInInside,
460 cellDataTemp[i].fluxData().upwindPotential(nPhaseIdx, indexInInside));
464 cellDataJ.fluxData().setVelocity(wPhaseIdx, indexInOutside,
465 cellDataTemp[i].fluxData().velocity(wPhaseIdx, indexInInside));
466 cellDataJ.fluxData().setVelocity(nPhaseIdx, indexInOutside,
467 cellDataTemp[i].fluxData().velocity(nPhaseIdx, indexInInside));
468 cellDataJ.fluxData().setUpwindPotential(wPhaseIdx, indexInOutside,
469 cellDataTemp[i].fluxData().upwindPotential(wPhaseIdx, indexInInside));
470 cellDataJ.fluxData().setUpwindPotential(nPhaseIdx, indexInOutside,
471 cellDataTemp[i].fluxData().upwindPotential(nPhaseIdx, indexInInside));
477 else if (eIdxGlobal[i] == eIdxGlobalJ)
479 if (levelJ >= levelI)
481 cellDataJ.fluxData().setVelocity(wPhaseIdx, indexInOutside,
482 cellDataTemp[i].fluxData().velocity(wPhaseIdx, indexInOutside));
483 cellDataJ.fluxData().setVelocity(nPhaseIdx, indexInOutside,
484 cellDataTemp[i].fluxData().velocity(nPhaseIdx, indexInOutside));
485 cellDataJ.fluxData().setUpwindPotential(wPhaseIdx, indexInOutside,
486 cellDataTemp[i].fluxData().upwindPotential(wPhaseIdx, indexInOutside));
487 cellDataJ.fluxData().setUpwindPotential(nPhaseIdx, indexInOutside,
488 cellDataTemp[i].fluxData().upwindPotential(nPhaseIdx, indexInOutside));
492 cellData.fluxData().setVelocity(wPhaseIdx, indexInInside,
493 cellDataTemp[i].fluxData().velocity(wPhaseIdx, indexInOutside));
494 cellData.fluxData().setVelocity(nPhaseIdx, indexInInside,
495 cellDataTemp[i].fluxData().velocity(nPhaseIdx, indexInOutside));
496 cellData.fluxData().setUpwindPotential(wPhaseIdx, indexInInside,
497 cellDataTemp[i].fluxData().upwindPotential(wPhaseIdx, indexInOutside));
498 cellData.fluxData().setUpwindPotential(nPhaseIdx, indexInInside,
499 cellDataTemp[i].fluxData().upwindPotential(nPhaseIdx, indexInOutside));
507 if (levelI == levelJ)
509 cellData.fluxData().setVelocityMarker(indexInInside);
510 cellDataJ.fluxData().setVelocityMarker(indexInOutside);
522template<
class TypeTag>
525 auto element = intersection.inside();
528 int isIndex = intersection.indexInInside();
531 const Dune::FieldVector<Scalar, dim>& unitOuterNormal = intersection.centerUnitOuterNormal();
533 BoundaryTypes bcType;
535 problem_.boundaryTypes(bcType, intersection);
536 PrimaryVariables boundValues(0.0);
538 if (bcType.isDirichlet(pressEqIdx))
540 problem_.dirichlet(boundValues, intersection);
543 const GlobalPosition& globalPosI = element.geometry().center();
546 const GlobalPosition& globalPosJ = intersection.geometry().center();
549 Scalar lambdaWI = cellData.mobility(wPhaseIdx);
550 Scalar lambdaNwI = cellData.mobility(nPhaseIdx);
553 Scalar pcI = cellData.capillaryPressure();
556 GlobalPosition distVec = globalPosJ - globalPosI;
559 Scalar dist = distVec.two_norm();
563 DimMatrix meanPermeability(0);
565 problem_.spatialParams().meanK(meanPermeability, problem_.spatialParams().intrinsicPermeability(element));
572 if (bcType.isDirichlet(satEqIdx))
574 switch (saturationType_)
578 satW = boundValues[saturationIdx];
583 satW = 1 - boundValues[saturationIdx];
590 satW = cellData.saturation(wPhaseIdx);
593 const Scalar pressBound = boundValues[pressureIdx];
598 const auto fluidMatrixInteraction = Deprecated::makePcKrSw(Scalar{}, problem_.spatialParams(), element);
600 const Scalar pcBound = fluidMatrixInteraction.pc(satW);
603 Scalar pressWBound = 0;
604 Scalar pressNwBound = 0;
605 if (pressureType_ == pw)
607 pressWBound = pressBound;
608 pressNwBound = pressBound + pcBound;
610 else if (pressureType_ == pn)
612 pressWBound = pressBound - pcBound;
613 pressNwBound = pressBound;
616 const Scalar lambdaWBound = fluidMatrixInteraction.krw(satW) / viscosity_[wPhaseIdx];
617 const Scalar lambdaNwBound = fluidMatrixInteraction.krn(satW) / viscosity_[nPhaseIdx];
619 Scalar potentialDiffW = cellData.fluxData().upwindPotential(wPhaseIdx, isIndex);
620 Scalar potentialDiffNw = cellData.fluxData().upwindPotential(nPhaseIdx, isIndex);
623 potentialDiffW = (cellData.pressure(wPhaseIdx) - pressWBound);
624 potentialDiffNw = (cellData.pressure(nPhaseIdx) - pressNwBound);
626 potentialDiffW += density_[wPhaseIdx] * (distVec * problem_.gravity());
627 potentialDiffNw += density_[nPhaseIdx] * (distVec * problem_.gravity());
630 cellData.fluxData().setUpwindPotential(wPhaseIdx, isIndex, potentialDiffW);
631 cellData.fluxData().setUpwindPotential(nPhaseIdx, isIndex, potentialDiffNw);
634 Scalar lambdaW = (potentialDiffW > 0.) ? lambdaWI : lambdaWBound;
635 lambdaW = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffW, 0.0, 1.0e-30)) ? 0.5 * (lambdaWI + lambdaWBound) : lambdaW;
636 Scalar lambdaNw = (potentialDiffNw > 0.) ? lambdaNwI : lambdaNwBound;
637 lambdaNw = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffNw, 0.0, 1.0e-30)) ? 0.5 * (lambdaNwI + lambdaNwBound) : lambdaNw;
643 Dune::FieldVector<Scalar, dimWorld> velocityW(unitOuterNormal);
644 Dune::FieldVector<Scalar, dimWorld> velocityNw(unitOuterNormal);
648 Scalar areaScaling = (unitOuterNormal * distVec);
651 Scalar gravityTermW = (problem_.gravity() * distVec) * density_[wPhaseIdx] * areaScaling;
652 Scalar gravityTermNw = (problem_.gravity() * distVec) * density_[nPhaseIdx] * areaScaling;
655 switch (pressureType_)
659 velocityW *= lambdaW * scalarPerm * ((cellData.pressure(wPhaseIdx) - pressBound) / dist + gravityTermW);
660 velocityNw *= lambdaNw * scalarPerm * ((cellData.pressure(wPhaseIdx) - pressBound) / dist + gravityTermNw)
661 + 0.5 * (lambdaNwI + lambdaNwBound) * scalarPerm * (pcI - pcBound) / dist;
666 velocityW *= lambdaW * scalarPerm * ((cellData.pressure(nPhaseIdx) - pressBound) / dist + gravityTermW)
667 - 0.5 * (lambdaWI + lambdaWBound) * scalarPerm * (pcI - pcBound) / dist;
668 velocityNw *= lambdaNw * scalarPerm * ((cellData.pressure(nPhaseIdx) - pressBound) / dist + gravityTermNw);
674 cellData.fluxData().setVelocity(wPhaseIdx, isIndex, velocityW);
675 cellData.fluxData().setVelocity(nPhaseIdx, isIndex, velocityNw);
676 cellData.fluxData().setVelocityMarker(isIndex);
680 else if (bcType.isNeumann(pressEqIdx))
682 problem_.neumann(boundValues, intersection);
684 Dune::FieldVector<Scalar, dimWorld> velocityW(unitOuterNormal);
685 Dune::FieldVector<Scalar, dimWorld> velocityNw(unitOuterNormal);
687 velocityW *= boundValues[wPhaseIdx];
688 velocityNw *= boundValues[nPhaseIdx];
690 velocityW /= density_[wPhaseIdx];
691 velocityNw /= density_[nPhaseIdx];
694 cellData.fluxData().setUpwindPotential(wPhaseIdx, isIndex, boundValues[wPhaseIdx]);
695 cellData.fluxData().setUpwindPotential(nPhaseIdx, isIndex, boundValues[nPhaseIdx]);
697 cellData.fluxData().setVelocity(wPhaseIdx, isIndex, velocityW);
698 cellData.fluxData().setVelocity(nPhaseIdx, isIndex, velocityNw);
699 cellData.fluxData().setVelocityMarker(isIndex);
703 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:77
FVMPFALInteractionVolume< TypeTag > InteractionVolume
Type of the interaction volume objects.
Definition: 2dpressureadaptive.hh:154
void addOutputVtkFields(MultiWriter &writer)
Adds pressure output to the output file.
Definition: 2dpressureadaptive.hh:376
void storePressureSolution()
Globally stores the pressure solution.
Definition: 2dpressureadaptive.hh:233
void initialize()
Initializes the pressure model.
Definition: 2dpressureadaptive.hh:202
void updateMaterialLaws()
constitutive functions are initialized and stored in the variables object
Definition: 2dpressureadaptive.hh:2552
void update()
Pressure update.
Definition: 2dpressureadaptive.hh:300
Class for the calculation of velocities from the pressure solution of an IMPES scheme using a grid ad...
Definition: 2dpressurevelocityadaptive.hh:50
void addOutputVtkFields(MultiWriter &writer)
Adds velocity output to the output file.
Definition: 2dpressurevelocityadaptive.hh:202
void updateVelocity()
Function for updating the velocity field if iterations are necessary in the transport solution.
Definition: 2dpressurevelocityadaptive.hh:132
bool calculateVelocityInTransport()
Indicates if velocity is reconstructed in the pressure step or in the transport step.
Definition: 2dpressurevelocityadaptive.hh:186
FvMpfaL2dPressureVelocity2pAdaptive(Problem &problem)
Constructs a FvMpfaL2dPressureVelocity2pAdaptive object.
Definition: 2dpressurevelocityadaptive.hh:113
void initialize()
Initializes pressure and velocity.
Definition: 2dpressurevelocityadaptive.hh:147
void calculateVelocityOnBoundary(const Intersection &intersection, CellData &cellData)
Calculates the velocity at a boundary.
Definition: 2dpressurevelocityadaptive.hh:523
void calculateVelocity()
Calculates the velocities at all cell-cell interfaces.
Definition: 2dpressurevelocityadaptive.hh:229
void update()
Pressure and velocity update.
Definition: 2dpressurevelocityadaptive.hh:174
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