22#ifndef DUMUX_MPFAO2DPRESSUREVELOCITIES2P_HH
23#define DUMUX_MPFAO2DPRESSUREVELOCITIES2P_HH
25#include <dune/common/float_cmp.hh>
51 dim = GridView::dimension, dimWorld = GridView::dimensionworld
59 using PrimaryVariables =
typename SolutionTypes::PrimaryVariables;
62 using MaterialLaw =
typename SpatialParams::MaterialLaw;
67 using Vertex =
typename GridView::Traits::template Codim<dim>::Entity;
68 using Intersection =
typename GridView::Intersection;
70 using ReferenceElements = Dune::ReferenceElements<Scalar, dim>;
77 wPhaseIdx = Indices::wPhaseIdx,
78 nPhaseIdx = Indices::nPhaseIdx,
79 pw = Indices::pressureW,
80 pn = Indices::pressureNw,
81 vw = Indices::velocityW,
82 vn = Indices::velocityNw,
83 sw = Indices::saturationW,
84 sn = Indices::saturationNw,
85 pressureIdx = Indices::pressureIdx,
86 saturationIdx = Indices::saturationIdx,
87 pressEqIdx = Indices::pressureEqIdx,
88 satEqIdx = Indices::satEqIdx,
89 numPhases = getPropValue<TypeTag, Properties::NumPhases>()
92 using Element =
typename GridView::template Codim<0>::Entity;
94 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
95 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
104 ParentType(problem), problem_(problem), velocity_(problem)
106 density_[wPhaseIdx] = 0.;
107 density_[nPhaseIdx] = 0.;
108 viscosity_[wPhaseIdx] = 0.;
109 viscosity_[nPhaseIdx] = 0.;
111 calcVelocityInTransport_ = getParam<bool>(
"MPFA.CalcVelocityInTransport");
139 const auto element = *problem_.gridView().template begin<0>();
140 FluidState fluidState;
141 fluidState.setPressure(wPhaseIdx, problem_.referencePressure(element));
142 fluidState.setPressure(nPhaseIdx, problem_.referencePressure(element));
143 fluidState.setTemperature(problem_.temperature(element));
144 fluidState.setSaturation(wPhaseIdx, 1.);
145 fluidState.setSaturation(nPhaseIdx, 0.);
152 velocity_.initialize();
178 return calcVelocityInTransport_;
191 template<
class MultiWriter>
195 velocity_.addOutputVtkFields(writer);
206 Scalar density_[numPhases];
207 Scalar viscosity_[numPhases];
208 bool calcVelocityInTransport_;
210 static const int pressureType_ = getPropValue<TypeTag, Properties::PressureFormulation>();
212 static const int saturationType_ = getPropValue<TypeTag, Properties::SaturationFormulation>();
221template<
class TypeTag>
225 for (
const auto& vertex : vertices(problem_.gridView()))
227 int vIdxGlobal = problem_.variables().index(vertex);
239 int eIdxGlobal1 = problem_.variables().index(element1);
240 int eIdxGlobal2 = problem_.variables().index(element2);
241 int eIdxGlobal3 = problem_.variables().index(element3);
242 int eIdxGlobal4 = problem_.variables().index(element4);
245 CellData& cellData1 = problem_.variables().cellData(eIdxGlobal1);
246 CellData& cellData2 = problem_.variables().cellData(eIdxGlobal2);
247 CellData& cellData3 = problem_.variables().cellData(eIdxGlobal3);
248 CellData& cellData4 = problem_.variables().cellData(eIdxGlobal4);
250 velocity_.calculateInnerInteractionVolumeVelocity(interactionVolume, cellData1, cellData2, cellData3,
251 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);
294template<
class TypeTag>
297 int numVertices = intersection.geometry().corners();
299 auto elementI = intersection.inside();
300 auto elementJ = intersection.outside();
302 int eIdxGlobalI = problem_.variables().index(elementI);
303 int eIdxGlobalJ = problem_.variables().index(elementJ);
305 CellData& cellDataJ = problem_.variables().cellData(eIdxGlobalJ);
307 const auto referenceElement = ReferenceElements::general(elementI.type());
309 int indexInInside = intersection.indexInInside();
310 int indexInOutside = intersection.indexInOutside();
312 Dune::FieldVector<CellData, 4> cellDataTemp;
314 for (
int vIdx = 0; vIdx < numVertices; vIdx++)
316 int localVertIdx = referenceElement.subEntity(indexInInside, dim - 1, vIdx, dim);
318 int vIdxGlobal = problem_.variables().index(elementI.template subEntity<dim>(localVertIdx));
332 eIdxGlobal[0] = problem_.variables().index(element1);
333 eIdxGlobal[1] = problem_.variables().index(element2);
334 eIdxGlobal[2] = problem_.variables().index(element3);
335 eIdxGlobal[3] = problem_.variables().index(element4);
338 cellDataTemp[0] = problem_.variables().cellData(eIdxGlobal[0]);
339 cellDataTemp[1] = problem_.variables().cellData(eIdxGlobal[1]);
340 cellDataTemp[2] = problem_.variables().cellData(eIdxGlobal[2]);
341 cellDataTemp[3] = problem_.variables().cellData(eIdxGlobal[3]);
343 velocity_.calculateInnerInteractionVolumeVelocity(interactionVolume, cellDataTemp[0], cellDataTemp[1],
344 cellDataTemp[2], cellDataTemp[3], this->innerBoundaryVolumeFaces_);
346 for (
int i = 0; i < 4; i++)
348 if (eIdxGlobal[i] == eIdxGlobalI)
350 cellData.fluxData().setVelocity(wPhaseIdx, indexInInside,
351 cellDataTemp[i].fluxData().velocity(wPhaseIdx, indexInInside));
352 cellData.fluxData().setVelocity(nPhaseIdx, indexInInside,
353 cellDataTemp[i].fluxData().velocity(nPhaseIdx, indexInInside));
354 cellData.fluxData().setUpwindPotential(wPhaseIdx, indexInInside,
355 cellDataTemp[i].fluxData().upwindPotential(wPhaseIdx, indexInInside));
356 cellData.fluxData().setUpwindPotential(nPhaseIdx, indexInInside,
357 cellDataTemp[i].fluxData().upwindPotential(nPhaseIdx, indexInInside));
359 else if (eIdxGlobal[i] == eIdxGlobalJ)
361 cellDataJ.fluxData().setVelocity(wPhaseIdx, indexInOutside,
362 cellDataTemp[i].fluxData().velocity(wPhaseIdx, indexInOutside));
363 cellDataJ.fluxData().setVelocity(nPhaseIdx, indexInOutside,
364 cellDataTemp[i].fluxData().velocity(nPhaseIdx, indexInOutside));
365 cellDataJ.fluxData().setUpwindPotential(wPhaseIdx, indexInOutside,
366 cellDataTemp[i].fluxData().upwindPotential(wPhaseIdx, indexInOutside));
367 cellDataJ.fluxData().setUpwindPotential(nPhaseIdx, indexInOutside,
368 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);
512 Scalar gravityTermW = (problem_.gravity() * distVec) * density_[wPhaseIdx] * areaScaling;
513 Scalar gravityTermNw = (problem_.gravity() * distVec) * density_[nPhaseIdx] * areaScaling;
516 switch (pressureType_)
520 velocityW *= lambdaW * scalarPerm * ((cellData.pressure(wPhaseIdx) - pressBound) / dist + gravityTermW);
521 velocityNw *= lambdaNw * scalarPerm * ((cellData.pressure(wPhaseIdx) - pressBound) / dist + gravityTermNw)
522 + 0.5 * (lambdaNwI + lambdaNwBound) * scalarPerm * (pcI - pcBound) / dist;
527 velocityW *= lambdaW * scalarPerm * ((cellData.pressure(nPhaseIdx) - pressBound) / dist + gravityTermW)
528 - 0.5 * (lambdaWI + lambdaWBound) * scalarPerm * (pcI - pcBound) / dist;
529 velocityNw *= lambdaNw * scalarPerm * ((cellData.pressure(nPhaseIdx) - pressBound) / dist + gravityTermNw);
535 cellData.fluxData().setVelocity(wPhaseIdx, isIndex, velocityW);
536 cellData.fluxData().setVelocity(nPhaseIdx, isIndex, velocityNw);
537 cellData.fluxData().setVelocityMarker(isIndex);
541 else if (bcType.isNeumann(pressEqIdx))
543 problem_.neumann(boundValues, intersection);
545 Dune::FieldVector<Scalar, dimWorld> velocityW(unitOuterNormal);
546 Dune::FieldVector<Scalar, dimWorld> velocityNw(unitOuterNormal);
548 velocityW *= boundValues[wPhaseIdx];
549 velocityNw *= boundValues[nPhaseIdx];
551 velocityW /= density_[wPhaseIdx];
552 velocityNw /= density_[nPhaseIdx];
555 cellData.fluxData().setUpwindPotential(wPhaseIdx, isIndex, boundValues[wPhaseIdx]);
556 cellData.fluxData().setUpwindPotential(nPhaseIdx, isIndex, boundValues[nPhaseIdx]);
558 cellData.fluxData().setVelocity(wPhaseIdx, isIndex, velocityW);
559 cellData.fluxData().setVelocity(nPhaseIdx, isIndex, velocityNw);
560 cellData.fluxData().setVelocityMarker(isIndex);
564 DUNE_THROW(Dune::NotImplemented,
"No valid boundary condition type defined for pressure equation!");
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
Finite volume MPFA O-method discretization of a two-phase flow pressure equation of the sequential IM...
Definition: omethod/2dpressure.hh:67
void updateMaterialLaws()
Constitutive functions are initialized and stored in the variables object.
Definition: omethod/2dpressure.hh:1934
void addOutputVtkFields(MultiWriter &writer)
Adds pressure output to the output file.
Definition: omethod/2dpressure.hh:327
void update()
Pressure update.
Definition: omethod/2dpressure.hh:213
void initialize()
Initializes the pressure model.
Definition: omethod/2dpressure.hh:177
void storePressureSolution()
Globally stores the pressure solution.
Definition: omethod/2dpressure.hh:254
Class for the calculation of velocities from the pressure solution of an IMPES scheme using a MPFA O-...
Definition: omethod/2dpressurevelocity.hh:44
bool calculateVelocityInTransport()
Indicates if velocity is reconstructed in the pressure step or in the transport step.
Definition: omethod/2dpressurevelocity.hh:176
void calculateVelocityOnBoundary(const Intersection &intersection, CellData &cellData)
Calculates the velocity at a boundary.
Definition: omethod/2dpressurevelocity.hh:388
void addOutputVtkFields(MultiWriter &writer)
Adds velocity output to the output file.
Definition: omethod/2dpressurevelocity.hh:192
void update()
Pressure and velocity update.
Definition: omethod/2dpressurevelocity.hh:164
void calculateVelocity()
Calculates the velocities at all cell-cell interfaces.
Definition: omethod/2dpressurevelocity.hh:222
FvMpfaO2dPressureVelocity2p(Problem &problem)
Constructs a FvMpfaO2dPressureVelocity2p object.
Definition: omethod/2dpressurevelocity.hh:103
void updateVelocity()
Function for updating the velocity field if iterations are necessary in the transport solution.
Definition: omethod/2dpressurevelocity.hh:122
void initialize()
Initializes pressure and velocity.
Definition: omethod/2dpressurevelocity.hh:137
FvMpfaO2dVelocity2P< TypeTag > velocity()
Definition: omethod/2dpressurevelocity.hh:198
Class for calculating velocities from cell-wise constant pressure values.
Definition: omethod/2dvelocity.hh:56
Class including the information of an interaction volume of a MPFA O-method that does not change with...
Definition: ointeractionvolume.hh:38
bool isInnerVolume()
Returns true if the interaction volume is completely inside the model domain.
Definition: ointeractionvolume.hh:291
bool isOutsideFace(int subVolumeFaceIdx)
Returns true if an interaction volume flux face is outside the model domain.
Definition: ointeractionvolume.hh:276
int getFaceIndexFromSubVolume(int subVolumeIdx, int subVolumeFaceIdx)
Map from local interaction volume numbering on element to numbering on interaction volume.
Definition: ointeractionvolume.hh:245
Element getSubVolumeElement(int subVolumeIdx)
Get an element of the interaction volume.
Definition: ointeractionvolume.hh:256
The finite volume base class for the solution of a pressure equation.
Definition: sequential/cellcentered/pressure.hh:49
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.