24#ifndef DUMUX_FVMPFAL2PFABOUND3DVELOCITIES2P_HH
25#define DUMUX_FVMPFAL2PFABOUND3DVELOCITIES2P_HH
27#include <dune/common/float_cmp.hh>
54 dim = GridView::dimension, dimWorld = GridView::dimensionworld
61 pw = Indices::pressureW,
62 pn = Indices::pressureNw,
63 sw = Indices::saturationW,
64 sn = Indices::saturationNw,
65 wPhaseIdx = Indices::wPhaseIdx,
66 nPhaseIdx = Indices::nPhaseIdx,
67 pressureIdx = Indices::pressureIdx,
68 saturationIdx = Indices::saturationIdx,
69 pressEqIdx = Indices::pressureEqIdx,
70 satEqIdx = Indices::satEqIdx,
71 numPhases = getPropValue<TypeTag, Properties::NumPhases>()
80 using PrimaryVariables =
typename SolutionTypes::PrimaryVariables;
84 using Intersection =
typename GridView::Intersection;
86 using Element =
typename GridView::template Codim<0>::Entity;
88 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
90 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
99 ParentType(problem), problem_(problem), velocity_(problem)
101 density_[wPhaseIdx] = 0.;
102 density_[nPhaseIdx] = 0.;
103 viscosity_[wPhaseIdx] = 0.;
104 viscosity_[nPhaseIdx] = 0.;
106 calcVelocityInTransport_ = getParam<bool>(
"MPFA.CalcVelocityInTransport");
135 const auto element = *problem_.gridView().template begin<0>();
136 FluidState fluidState;
137 fluidState.setPressure(wPhaseIdx, problem_.referencePressure(element));
138 fluidState.setPressure(nPhaseIdx, problem_.referencePressure(element));
139 fluidState.setTemperature(problem_.temperature(element));
140 fluidState.setSaturation(wPhaseIdx, 1.);
141 fluidState.setSaturation(nPhaseIdx, 0.);
148 velocity_.initialize();
176 return calcVelocityInTransport_;
189 template<
class MultiWriter>
193 velocity_.addOutputVtkFields(writer);
200 Scalar density_[numPhases];
201 Scalar viscosity_[numPhases];
202 bool calcVelocityInTransport_;
205 static const int pressureType_ = getPropValue<TypeTag, Properties::PressureFormulation>();
207 static const int saturationType_ = getPropValue<TypeTag, Properties::SaturationFormulation>();
216template<
class TypeTag>
220 for (
const auto& vertex : vertices(problem_.gridView()))
222 int vIdxGlobal = problem_.variables().index(vertex);
224 InteractionVolume& interactionVolume = this->interactionVolumes_.interactionVolume(vIdxGlobal);
227 if (interactionVolume.isInnerVolume())
231 for (
int i = 0; i < 8; i++)
233 eIdxGlobal[i] = problem_.variables().index(interactionVolume.getSubVolumeElement(i));
237 CellData & cellData1 = problem_.variables().cellData(eIdxGlobal[0]);
238 CellData & cellData2 = problem_.variables().cellData(eIdxGlobal[1]);
239 CellData & cellData3 = problem_.variables().cellData(eIdxGlobal[2]);
240 CellData & cellData4 = problem_.variables().cellData(eIdxGlobal[3]);
241 CellData & cellData5 = problem_.variables().cellData(eIdxGlobal[4]);
242 CellData & cellData6 = problem_.variables().cellData(eIdxGlobal[5]);
243 CellData & cellData7 = problem_.variables().cellData(eIdxGlobal[6]);
244 CellData & cellData8 = problem_.variables().cellData(eIdxGlobal[7]);
246 velocity_.calculateInnerInteractionVolumeVelocity(interactionVolume,
247 cellData1, cellData2, cellData3, cellData4,
248 cellData5, cellData6, cellData7, cellData8,
249 this->interactionVolumes_, this->transmissibilityCalculator_);
254 for (
int elemIdx = 0; elemIdx < 8; elemIdx++)
256 if (!interactionVolume.hasSubVolumeElement(elemIdx))
260 bool isOutside =
false;
261 for (
int fIdx = 0; fIdx < dim; fIdx++)
263 int intVolFaceIdx = interactionVolume.getFaceIndexFromSubVolume(elemIdx, fIdx);
264 if (interactionVolume.isOutsideFace(intVolFaceIdx))
275 int eIdxGlobal = problem_.variables().index(interactionVolume.getSubVolumeElement(elemIdx));
277 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
279 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 DUNE_UNUSED = problem_.variables().index(elementI);
305 int eIdxGlobalJ = problem_.variables().index(elementJ);
307 CellData& cellDataJ = problem_.variables().cellData(eIdxGlobalJ);
309 const auto refElement = referenceElement(elementI);
311 int indexInInside = intersection.indexInInside();
312 int indexInOutside = intersection.indexInOutside();
314 Dune::FieldVector<CellData, 8> cellDataTemp;
316 for (
int vIdx = 0; vIdx < numVertices; vIdx++)
318 int localVertIdx = refElement.subEntity(indexInInside, 1, vIdx, dim);
320 int vIdxGlobal = problem_.variables().index(elementI.template subEntity<dim>(localVertIdx));
322 InteractionVolume& interactionVolume = this->interactionVolumes_.interactionVolume(vIdxGlobal);
324 if (interactionVolume.isInnerVolume())
327 int localMpfaElemIdxI = 0;
328 int localMpfaElemIdxJ = 0;
331 for (
int i = 0; i < 8; i++)
333 auto elem = interactionVolume.getSubVolumeElement(i);
335 if (elem == elementI)
336 localMpfaElemIdxI = i;
337 else if (elem == elementJ)
338 localMpfaElemIdxJ = i;
340 eIdxGlobal[i] = problem_.variables().index(elem);
341 cellDataTemp[i] = problem_.variables().cellData(eIdxGlobal[i]);
344 int mpfaFaceIdx = IndexTranslator::getFaceIndexFromElements(localMpfaElemIdxI, localMpfaElemIdxJ);
348 velocity_.calculateInnerInteractionVolumeVelocity(interactionVolume,
349 cellDataTemp[0], cellDataTemp[1], cellDataTemp[2], cellDataTemp[3],
350 cellDataTemp[4], cellDataTemp[5], cellDataTemp[6], cellDataTemp[7],
351 this->interactionVolumes_, this->transmissibilityCalculator_, mpfaFaceIdx);
354 cellData.fluxData().setVelocity(wPhaseIdx, indexInInside,
355 cellDataTemp[localMpfaElemIdxI].fluxData().velocity(wPhaseIdx, indexInInside));
356 cellData.fluxData().setVelocity(nPhaseIdx, indexInInside,
357 cellDataTemp[localMpfaElemIdxI].fluxData().velocity(nPhaseIdx, indexInInside));
358 cellData.fluxData().setUpwindPotential(wPhaseIdx, indexInInside,
359 cellDataTemp[localMpfaElemIdxI].fluxData().upwindPotential(wPhaseIdx, indexInInside));
360 cellData.fluxData().setUpwindPotential(nPhaseIdx, indexInInside,
361 cellDataTemp[localMpfaElemIdxI].fluxData().upwindPotential(nPhaseIdx, indexInInside));
363 cellDataJ.fluxData().setVelocity(wPhaseIdx, indexInOutside,
364 cellDataTemp[localMpfaElemIdxJ].fluxData().velocity(wPhaseIdx, indexInOutside));
365 cellDataJ.fluxData().setVelocity(nPhaseIdx, indexInOutside,
366 cellDataTemp[localMpfaElemIdxJ].fluxData().velocity(nPhaseIdx, indexInOutside));
367 cellDataJ.fluxData().setUpwindPotential(wPhaseIdx, indexInOutside,
368 cellDataTemp[localMpfaElemIdxJ].fluxData().upwindPotential(wPhaseIdx, indexInOutside));
369 cellDataJ.fluxData().setUpwindPotential(nPhaseIdx, indexInOutside,
370 cellDataTemp[localMpfaElemIdxJ].fluxData().upwindPotential(nPhaseIdx, indexInOutside));
373 cellData.fluxData().setVelocityMarker(indexInInside);
374 cellDataJ.fluxData().setVelocityMarker(indexInOutside);
385template<
class TypeTag>
388 auto element = intersection.inside();
391 int isIndex = intersection.indexInInside();
394 const Dune::FieldVector<Scalar, dim>& unitOuterNormal = intersection.centerUnitOuterNormal();
396 BoundaryTypes bcType;
398 problem_.boundaryTypes(bcType, intersection);
399 PrimaryVariables boundValues(0.0);
401 if (bcType.isDirichlet(pressEqIdx))
403 problem_.dirichlet(boundValues, intersection);
406 const GlobalPosition& globalPosI = element.geometry().center();
409 const GlobalPosition& globalPosJ = intersection.geometry().center();
412 Scalar lambdaWI = cellData.mobility(wPhaseIdx);
413 Scalar lambdaNwI = cellData.mobility(nPhaseIdx);
416 Scalar pcI = cellData.capillaryPressure();
419 GlobalPosition distVec = globalPosJ - globalPosI;
422 Scalar dist = distVec.two_norm();
426 DimMatrix meanPermeability(0);
428 problem_.spatialParams().meanK(meanPermeability, problem_.spatialParams().intrinsicPermeability(element));
435 if (bcType.isDirichlet(satEqIdx))
437 switch (saturationType_)
441 satW = boundValues[saturationIdx];
446 satW = 1 - boundValues[saturationIdx];
453 satW = cellData.saturation(wPhaseIdx);
456 const Scalar pressBound = boundValues[pressureIdx];
461 const auto fluidMatrixInteraction = Deprecated::makePcKrSw(Scalar{}, problem_.spatialParams(), element);
463 const Scalar pcBound = fluidMatrixInteraction.pc(satW);
466 Scalar pressWBound = 0;
467 Scalar pressNwBound = 0;
468 if (pressureType_ == pw)
470 pressWBound = pressBound;
471 pressNwBound = pressBound + pcBound;
473 else if (pressureType_ == pn)
475 pressWBound = pressBound - pcBound;
476 pressNwBound = pressBound;
479 const Scalar lambdaWBound = fluidMatrixInteraction.krw(satW) / viscosity_[wPhaseIdx];
480 const Scalar lambdaNwBound = fluidMatrixInteraction.krn(satW) / viscosity_[nPhaseIdx];
482 Scalar potentialDiffW = cellData.fluxData().upwindPotential(wPhaseIdx, isIndex);
483 Scalar potentialDiffNw = cellData.fluxData().upwindPotential(nPhaseIdx, isIndex);
486 potentialDiffW = (cellData.pressure(wPhaseIdx) - pressWBound);
487 potentialDiffNw = (cellData.pressure(nPhaseIdx) - pressNwBound);
489 potentialDiffW += density_[wPhaseIdx] * (distVec * problem_.gravity());
490 potentialDiffNw += density_[nPhaseIdx] * (distVec * problem_.gravity());
493 cellData.fluxData().setUpwindPotential(wPhaseIdx, isIndex, potentialDiffW);
494 cellData.fluxData().setUpwindPotential(nPhaseIdx, isIndex, potentialDiffNw);
497 Scalar lambdaW = (potentialDiffW > 0.) ? lambdaWI : lambdaWBound;
498 lambdaW = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffW, 0.0, 1.0e-30)) ? 0.5 * (lambdaWI + lambdaWBound) : lambdaW;
499 Scalar lambdaNw = (potentialDiffNw > 0.) ? lambdaNwI : lambdaNwBound;
500 lambdaNw = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffNw, 0.0, 1.0e-30)) ? 0.5 * (lambdaNwI + lambdaNwBound) : lambdaNw;
506 Dune::FieldVector<Scalar, dimWorld> velocityW(unitOuterNormal);
507 Dune::FieldVector<Scalar, dimWorld> velocityNw(unitOuterNormal);
511 Scalar areaScaling = (unitOuterNormal * distVec);
513 Scalar gravityTermW = (problem_.gravity() * distVec) * density_[wPhaseIdx] * areaScaling;
514 Scalar gravityTermNw = (problem_.gravity() * distVec) * density_[nPhaseIdx] * areaScaling;
517 switch (pressureType_)
521 velocityW *= lambdaW * scalarPerm * ((cellData.pressure(wPhaseIdx) - pressBound) / dist + gravityTermW);
522 velocityNw *= lambdaNw * scalarPerm * ((cellData.pressure(wPhaseIdx) - pressBound) / dist + gravityTermNw)
523 + 0.5 * (lambdaNwI + lambdaNwBound) * scalarPerm * (pcI - pcBound) / dist;
528 velocityW *= lambdaW * scalarPerm * ((cellData.pressure(nPhaseIdx) - pressBound) / dist + gravityTermW)
529 - 0.5 * (lambdaWI + lambdaWBound) * scalarPerm * (pcI - pcBound) / dist;
530 velocityNw *= lambdaNw * scalarPerm * ((cellData.pressure(nPhaseIdx) - pressBound) / dist + gravityTermNw);
536 cellData.fluxData().setVelocity(wPhaseIdx, isIndex, velocityW);
537 cellData.fluxData().setVelocity(nPhaseIdx, isIndex, velocityNw);
538 cellData.fluxData().setVelocityMarker(isIndex);
542 else if (bcType.isNeumann(pressEqIdx))
544 problem_.neumann(boundValues, intersection);
546 Dune::FieldVector<Scalar, dimWorld> velocityW(unitOuterNormal);
547 Dune::FieldVector<Scalar, dimWorld> velocityNw(unitOuterNormal);
549 velocityW *= boundValues[wPhaseIdx];
550 velocityNw *= boundValues[nPhaseIdx];
552 velocityW /= density_[wPhaseIdx];
553 velocityNw /= density_[nPhaseIdx];
556 cellData.fluxData().setUpwindPotential(wPhaseIdx, isIndex, boundValues[wPhaseIdx]);
557 cellData.fluxData().setUpwindPotential(nPhaseIdx, isIndex, boundValues[nPhaseIdx]);
559 cellData.fluxData().setVelocity(wPhaseIdx, isIndex, velocityW);
560 cellData.fluxData().setVelocity(nPhaseIdx, isIndex, velocityNw);
561 cellData.fluxData().setVelocityMarker(isIndex);
565 DUNE_THROW(Dune::NotImplemented,
"No valid boundary condition type defined for pressure equation!");
3-d finite Volume-MPFAL implementation of a two-phase pressure equation.
2-d velocity calculation using a 3-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
3-d finite volume MPFA L-method discretization of a two-phase flow pressure equation of the sequentia...
Definition: 3dpressure.hh:76
void update()
Pressure update.
Definition: 3dpressure.hh:273
GetPropType< TypeTag, Properties::MPFAInteractionVolume > InteractionVolume
Type for storing interaction volume information.
Definition: 3dpressure.hh:156
void storePressureSolution()
Globally stores the pressure solution.
Definition: 3dpressure.hh:207
void addOutputVtkFields(MultiWriter &writer)
Adds pressure output to the output file.
Definition: 3dpressure.hh:386
void updateMaterialLaws()
constitutive functions are initialized and stored in the variables object
Definition: 3dpressure.hh:2484
Class for the calculation of 3d velocities from the pressure solution of an IMPES scheme using a MPFA...
Definition: 3dpressurevelocity.hh:46
void update()
Pressure and velocity update.
Definition: 3dpressurevelocity.hh:161
FvMpfaL3dPressureVelocity2p(Problem &problem)
Constructs a FvMpfaL3dPressureVelocity2p object.
Definition: 3dpressurevelocity.hh:98
bool calculateVelocityInTransport()
Indicates if velocity is reconstructed in the pressure step or in the transport step.
Definition: 3dpressurevelocity.hh:174
void updateVelocity()
Function for updating the velocity field if iterations are necessary in the transport solution.
Definition: 3dpressurevelocity.hh:118
void calculateVelocityOnBoundary(const Intersection &intersection, CellData &cellData)
Calculates the velocity at a boundary.
Definition: 3dpressurevelocity.hh:386
void initialize(bool solveTwice=true)
Initializes pressure and velocity.
Definition: 3dpressurevelocity.hh:133
void addOutputVtkFields(MultiWriter &writer)
Adds velocity output to the output file.
Definition: 3dpressurevelocity.hh:190
void calculateVelocity()
Calculates the velocities at a cell-cell interfaces for the entire simulation grid.
Definition: 3dpressurevelocity.hh:217
Class for calculating 3-d velocities from cell-wise constant pressure values.
Definition: 3dvelocity.hh:57
The finite volume base class for the solution of a pressure equation.
Definition: sequential/cellcentered/pressure.hh:49
void initialize()
Initialize pressure model.
Definition: sequential/cellcentered/pressure.hh:213