22#ifndef DUMUX_MPFAL2DPRESSUREVELOCITY2P_HH
23#define DUMUX_MPFAL2DPRESSUREVELOCITY2P_HH
25#include <dune/common/float_cmp.hh>
53 dim = GridView::dimension, dimWorld = GridView::dimensionworld
59 using BoundaryTypes =
typename GET_PROP_TYPE(TypeTag, BoundaryTypes);
60 using SolutionTypes =
typename GET_PROP(TypeTag, SolutionTypes);
61 using PrimaryVariables =
typename SolutionTypes::PrimaryVariables;
63 using ReferenceElements = Dune::ReferenceElements<Scalar, dim>;
65 using SpatialParams =
typename GET_PROP_TYPE(TypeTag, SpatialParams);
66 using MaterialLaw =
typename SpatialParams::MaterialLaw;
68 using Indices =
typename GET_PROP_TYPE(TypeTag, ModelTraits)::Indices;
70 using FluidSystem =
typename GET_PROP_TYPE(TypeTag, FluidSystem);
71 using FluidState =
typename GET_PROP_TYPE(TypeTag, FluidState);
73 using IndexSet =
typename GridView::IndexSet;
74 using Intersection =
typename GridView::Intersection;
80 pw = Indices::pressureW,
81 pn = Indices::pressureNw,
82 sw = Indices::saturationW,
83 sn = Indices::saturationNw,
84 wPhaseIdx = Indices::wPhaseIdx,
85 nPhaseIdx = Indices::nPhaseIdx,
86 pressureIdx = Indices::pressureIdx,
87 saturationIdx = Indices::saturationIdx,
88 pressEqIdx = Indices::pressureEqIdx,
89 satEqIdx = Indices::satEqIdx,
93 using Element =
typename GridView::template Codim<0>::Entity;
95 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
96 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
105 ParentType(problem), problem_(problem), velocity_(problem)
107 density_[wPhaseIdx] = 0.;
108 density_[nPhaseIdx] = 0.;
109 viscosity_[wPhaseIdx] = 0.;
110 viscosity_[nPhaseIdx] = 0.;
112 calcVelocityInTransport_ = getParam<bool>(
"MPFA.CalcVelocityInTransport");
140 const auto element = *problem_.gridView().template begin<0>();
141 FluidState fluidState;
142 fluidState.setPressure(wPhaseIdx, problem_.referencePressure(element));
143 fluidState.setPressure(nPhaseIdx, problem_.referencePressure(element));
144 fluidState.setTemperature(problem_.temperature(element));
145 fluidState.setSaturation(wPhaseIdx, 1.);
146 fluidState.setSaturation(nPhaseIdx, 0.);
153 velocity_.initialize();
179 return calcVelocityInTransport_;
192 template<
class MultiWriter>
196 velocity_.addOutputVtkFields(writer);
203 Scalar density_[numPhases];
204 Scalar viscosity_[numPhases];
205 bool calcVelocityInTransport_;
208 static const int pressureType_ =
GET_PROP_VALUE(TypeTag, PressureFormulation);
220template<
class TypeTag>
224 for (
const auto& vertex : vertices(problem_.gridView()))
226 int vIdxGlobal = problem_.variables().index(vertex);
238 int eIdxGlobal1 = problem_.variables().index(element1);
239 int eIdxGlobal2 = problem_.variables().index(element2);
240 int eIdxGlobal3 = problem_.variables().index(element3);
241 int eIdxGlobal4 = problem_.variables().index(element4);
244 CellData& cellData1 = problem_.variables().cellData(eIdxGlobal1);
245 CellData& cellData2 = problem_.variables().cellData(eIdxGlobal2);
246 CellData& cellData3 = problem_.variables().cellData(eIdxGlobal3);
247 CellData& cellData4 = problem_.variables().cellData(eIdxGlobal4);
249 velocity_.calculateInnerInteractionVolumeVelocity(interactionVolume, cellData1, cellData2,
250 cellData3, cellData4, this->innerBoundaryVolumeFaces_);
256 for (
int elemIdx = 0; elemIdx < 2 * dim; elemIdx++)
258 bool isOutside =
false;
259 for (
int fIdx = 0; fIdx < dim; fIdx++)
273 int eIdxGlobal = problem_.variables().index(interactionVolume.
getSubVolumeElement(elemIdx));
275 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
277 velocity_.calculateBoundaryInteractionVolumeVelocity(interactionVolume, cellData, elemIdx);
296template<
class TypeTag>
299 int numVertices = intersection.geometry().corners();
301 auto elementI = intersection.inside();
302 auto elementJ = intersection.outside();
304 int eIdxGlobalI = problem_.variables().index(elementI);
305 int eIdxGlobalJ = problem_.variables().index(elementJ);
307 CellData& cellDataJ = problem_.variables().cellData(eIdxGlobalJ);
309 const auto referenceElement = ReferenceElements::general(elementI.type());
311 int indexInInside = intersection.indexInInside();
312 int indexInOutside = intersection.indexInOutside();
314 Dune::FieldVector<CellData, 4> cellDataTemp;
316 for (
int vIdx = 0; vIdx < numVertices; vIdx++)
318 int localVertIdx = referenceElement.subEntity(indexInInside, dim - 1, vIdx, dim);
320 int vIdxGlobal = problem_.variables().index(elementI.template subEntity<dim>(localVertIdx));
334 eIdxGlobal[0] = problem_.variables().index(element1);
335 eIdxGlobal[1] = problem_.variables().index(element2);
336 eIdxGlobal[2] = problem_.variables().index(element3);
337 eIdxGlobal[3] = problem_.variables().index(element4);
340 cellDataTemp[0] = problem_.variables().cellData(eIdxGlobal[0]);
341 cellDataTemp[1] = problem_.variables().cellData(eIdxGlobal[1]);
342 cellDataTemp[2] = problem_.variables().cellData(eIdxGlobal[2]);
343 cellDataTemp[3] = problem_.variables().cellData(eIdxGlobal[3]);
345 velocity_.calculateInnerInteractionVolumeVelocity(interactionVolume, cellDataTemp[0], cellDataTemp[1],
346 cellDataTemp[2], cellDataTemp[3], this->innerBoundaryVolumeFaces_);
348 for (
int i = 0; i < 4; i++)
350 if (eIdxGlobal[i] == eIdxGlobalI)
352 cellData.fluxData().setVelocity(wPhaseIdx, indexInInside,
353 cellDataTemp[i].fluxData().velocity(wPhaseIdx, indexInInside));
354 cellData.fluxData().setVelocity(nPhaseIdx, indexInInside,
355 cellDataTemp[i].fluxData().velocity(nPhaseIdx, indexInInside));
356 cellData.fluxData().setUpwindPotential(wPhaseIdx, indexInInside,
357 cellDataTemp[i].fluxData().upwindPotential(wPhaseIdx, indexInInside));
358 cellData.fluxData().setUpwindPotential(nPhaseIdx, indexInInside,
359 cellDataTemp[i].fluxData().upwindPotential(nPhaseIdx, indexInInside));
361 else if (eIdxGlobal[i] == eIdxGlobalJ)
363 cellDataJ.fluxData().setVelocity(wPhaseIdx, indexInOutside,
364 cellDataTemp[i].fluxData().velocity(wPhaseIdx, indexInOutside));
365 cellDataJ.fluxData().setVelocity(nPhaseIdx, indexInOutside,
366 cellDataTemp[i].fluxData().velocity(nPhaseIdx, indexInOutside));
367 cellDataJ.fluxData().setUpwindPotential(wPhaseIdx, indexInOutside,
368 cellDataTemp[i].fluxData().upwindPotential(wPhaseIdx, indexInOutside));
369 cellDataJ.fluxData().setUpwindPotential(nPhaseIdx, indexInOutside,
370 cellDataTemp[i].fluxData().upwindPotential(nPhaseIdx, indexInOutside));
375 cellData.fluxData().setVelocityMarker(indexInInside);
376 cellDataJ.fluxData().setVelocityMarker(indexInOutside);
387template<
class TypeTag>
390 auto element = intersection.inside();
393 int isIndex = intersection.indexInInside();
396 const Dune::FieldVector<Scalar, dim>& unitOuterNormal = intersection.centerUnitOuterNormal();
398 BoundaryTypes bcType;
400 problem_.boundaryTypes(bcType, intersection);
401 PrimaryVariables boundValues(0.0);
403 if (bcType.isDirichlet(pressEqIdx))
405 problem_.dirichlet(boundValues, intersection);
408 const GlobalPosition& globalPosI = element.geometry().center();
411 const GlobalPosition& globalPosJ = intersection.geometry().center();
414 Scalar lambdaWI = cellData.mobility(wPhaseIdx);
415 Scalar lambdaNwI = cellData.mobility(nPhaseIdx);
418 Scalar pcI = cellData.capillaryPressure();
421 GlobalPosition distVec = globalPosJ - globalPosI;
424 Scalar dist = distVec.two_norm();
428 DimMatrix meanPermeability(0);
430 problem_.spatialParams().meanK(meanPermeability, problem_.spatialParams().intrinsicPermeability(element));
437 if (bcType.isDirichlet(satEqIdx))
439 switch (saturationType_)
443 satW = boundValues[saturationIdx];
448 satW = 1 - boundValues[saturationIdx];
455 satW = cellData.saturation(wPhaseIdx);
458 Scalar pressBound = boundValues[pressureIdx];
459 Scalar pcBound = MaterialLaw::pc(problem_.spatialParams().materialLawParams(element), satW);
462 Scalar pressWBound = 0;
463 Scalar pressNwBound = 0;
464 if (pressureType_ == pw)
466 pressWBound = pressBound;
467 pressNwBound = pressBound + pcBound;
469 else if (pressureType_ == pn)
471 pressWBound = pressBound - pcBound;
472 pressNwBound = pressBound;
475 Scalar lambdaWBound = MaterialLaw::krw(problem_.spatialParams().materialLawParams(element), satW)
476 / viscosity_[wPhaseIdx];
477 Scalar lambdaNwBound = MaterialLaw::krn(problem_.spatialParams().materialLawParams(element), satW)
478 / viscosity_[nPhaseIdx];
480 Scalar potentialDiffW = cellData.fluxData().upwindPotential(wPhaseIdx, isIndex);
481 Scalar potentialDiffNw = cellData.fluxData().upwindPotential(nPhaseIdx, isIndex);
484 potentialDiffW = (cellData.pressure(wPhaseIdx) - pressWBound);
485 potentialDiffNw = (cellData.pressure(nPhaseIdx) - pressNwBound);
487 potentialDiffW += density_[wPhaseIdx] * (distVec * problem_.gravity());
488 potentialDiffNw += density_[nPhaseIdx] * (distVec * problem_.gravity());
491 cellData.fluxData().setUpwindPotential(wPhaseIdx, isIndex, potentialDiffW);
492 cellData.fluxData().setUpwindPotential(nPhaseIdx, isIndex, potentialDiffNw);
495 Scalar lambdaW = (potentialDiffW > 0.) ? lambdaWI : lambdaWBound;
496 lambdaW = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffW, 0.0, 1.0e-30)) ? 0.5 * (lambdaWI + lambdaWBound) : lambdaW;
497 Scalar lambdaNw = (potentialDiffNw > 0.) ? lambdaNwI : lambdaNwBound;
498 lambdaNw = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffNw, 0.0, 1.0e-30)) ? 0.5 * (lambdaNwI + lambdaNwBound) : lambdaNw;
504 Dune::FieldVector<Scalar, dimWorld> velocityW(unitOuterNormal);
505 Dune::FieldVector<Scalar, dimWorld> velocityNw(unitOuterNormal);
509 Scalar areaScaling = (unitOuterNormal * distVec);
511 Scalar gravityTermW = (problem_.gravity() * distVec) * density_[wPhaseIdx] * areaScaling;
512 Scalar gravityTermNw = (problem_.gravity() * distVec) * density_[nPhaseIdx] * areaScaling;
515 switch (pressureType_)
519 velocityW *= lambdaW * scalarPerm * ((cellData.pressure(wPhaseIdx) - pressBound) / dist + gravityTermW);
520 velocityNw *= lambdaNw * scalarPerm * ((cellData.pressure(wPhaseIdx) - pressBound) / dist + gravityTermNw)
521 + 0.5 * (lambdaNwI + lambdaNwBound) * scalarPerm * (pcI - pcBound) / dist;
526 velocityW *= lambdaW * scalarPerm * ((cellData.pressure(nPhaseIdx) - pressBound) / dist + gravityTermW)
527 - 0.5 * (lambdaWI + lambdaWBound) * scalarPerm * (pcI - pcBound) / dist;
528 velocityNw *= lambdaNw * scalarPerm * ((cellData.pressure(nPhaseIdx) - pressBound) / dist + gravityTermNw);
534 cellData.fluxData().setVelocity(wPhaseIdx, isIndex, velocityW);
535 cellData.fluxData().setVelocity(nPhaseIdx, isIndex, velocityNw);
536 cellData.fluxData().setVelocityMarker(isIndex);
540 else if (bcType.isNeumann(pressEqIdx))
542 problem_.neumann(boundValues, intersection);
544 Dune::FieldVector<Scalar, dimWorld> velocityW(unitOuterNormal);
545 Dune::FieldVector<Scalar, dimWorld> velocityNw(unitOuterNormal);
547 velocityW *= boundValues[wPhaseIdx];
548 velocityNw *= boundValues[nPhaseIdx];
550 velocityW /= density_[wPhaseIdx];
551 velocityNw /= density_[nPhaseIdx];
554 cellData.fluxData().setUpwindPotential(wPhaseIdx, isIndex, boundValues[wPhaseIdx]);
555 cellData.fluxData().setUpwindPotential(nPhaseIdx, isIndex, boundValues[nPhaseIdx]);
557 cellData.fluxData().setVelocity(wPhaseIdx, isIndex, velocityW);
558 cellData.fluxData().setVelocity(nPhaseIdx, isIndex, velocityNw);
559 cellData.fluxData().setVelocityMarker(isIndex);
563 DUNE_THROW(Dune::NotImplemented,
"No valid boundary condition type defined for pressure equation!");
#define GET_PROP_VALUE(TypeTag, PropTagName)
Definition: propertysystemmacros.hh:282
#define GET_PROP(TypeTag, PropTagName)
Definition: propertysystemmacros.hh:281
#define GET_PROP_TYPE(TypeTag, PropTagName)
Definition: propertysystemmacros.hh:283
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
Property tag SaturationFormulation
The formulation of the saturation model.
Definition: porousmediumflow/2p/sequential/properties.hh:53
Property tag NumPhases
Number of phases in the system.
Definition: porousmediumflow/sequential/properties.hh:69
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
Finite volume MPFA L-method discretization of a two-phase flow pressure equation of the sequential IM...
Definition: lmethod/2dpressure.hh:72
void update()
Pressure update.
Definition: lmethod/2dpressure.hh:291
void initialize()
Initializes the pressure model.
Definition: lmethod/2dpressure.hh:194
FVMPFALInteractionVolume< TypeTag > InteractionVolume
Type of the interaction volume objects.
Definition: lmethod/2dpressure.hh:151
void storePressureSolution()
Globally stores the pressure solution.
Definition: lmethod/2dpressure.hh:228
void addOutputVtkFields(MultiWriter &writer)
Adds pressure output to the output file.
Definition: lmethod/2dpressure.hh:340
void updateMaterialLaws()
constitutive functions are initialized and stored in the variables object
Definition: lmethod/2dpressure.hh:1838
Class for the calculation of velocities from the pressure solution of an IMPES scheme using a MPFA L-...
Definition: lmethod/2dpressurevelocity.hh:46
void updateVelocity()
Function for updating the velocity field if iterations are necessary in the transport solution.
Definition: lmethod/2dpressurevelocity.hh:123
void calculateVelocityOnBoundary(const Intersection &intersection, CellData &cellData)
Calculates the velocity at a boundary.
Definition: lmethod/2dpressurevelocity.hh:388
FvMpfaL2dPressureVelocity2p(Problem &problem)
Constructs a FvMpfaL2dPressureVelocity2p object.
Definition: lmethod/2dpressurevelocity.hh:104
void initialize()
Initializes pressure and velocity.
Definition: lmethod/2dpressurevelocity.hh:138
void addOutputVtkFields(MultiWriter &writer)
Adds velocity output to the output file.
Definition: lmethod/2dpressurevelocity.hh:193
bool calculateVelocityInTransport()
Indicates if velocity is reconstructed in the pressure step or in the transport step.
Definition: lmethod/2dpressurevelocity.hh:177
void calculateVelocity()
Calculates the velocities at all cell-cell interfaces.
Definition: lmethod/2dpressurevelocity.hh:221
void update()
Pressure and velocity update.
Definition: lmethod/2dpressurevelocity.hh:165
Class for calculating 2-d velocities from cell-wise constant pressure values.
Definition: lmethod/2dvelocity.hh:58
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
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:48
Finite volume MPFA O-method discretization of a two-phase pressure equation of the sequential IMPES m...
Velocity calculation using a 2-d MPFA O-method.