24#ifndef DUMUX_FVMPFAL2PFABOUND3DVELOCITIES2P_HH
25#define DUMUX_FVMPFAL2PFABOUND3DVELOCITIES2P_HH
27#include <dune/common/float_cmp.hh>
52 dim = GridView::dimension, dimWorld = GridView::dimensionworld
59 pw = Indices::pressureW,
60 pn = Indices::pressureNw,
61 sw = Indices::saturationW,
62 sn = Indices::saturationNw,
63 wPhaseIdx = Indices::wPhaseIdx,
64 nPhaseIdx = Indices::nPhaseIdx,
65 pressureIdx = Indices::pressureIdx,
66 saturationIdx = Indices::saturationIdx,
67 pressEqIdx = Indices::pressureEqIdx,
68 satEqIdx = Indices::satEqIdx,
69 numPhases = getPropValue<TypeTag, Properties::NumPhases>()
78 using PrimaryVariables =
typename SolutionTypes::PrimaryVariables;
82 using Intersection =
typename GridView::Intersection;
84 using Element =
typename GridView::template Codim<0>::Entity;
86 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
88 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
97 ParentType(problem), problem_(problem), velocity_(problem)
99 density_[wPhaseIdx] = 0.;
100 density_[nPhaseIdx] = 0.;
101 viscosity_[wPhaseIdx] = 0.;
102 viscosity_[nPhaseIdx] = 0.;
104 calcVelocityInTransport_ = getParam<bool>(
"MPFA.CalcVelocityInTransport");
133 const auto element = *problem_.gridView().template begin<0>();
134 FluidState fluidState;
135 fluidState.setPressure(wPhaseIdx, problem_.referencePressure(element));
136 fluidState.setPressure(nPhaseIdx, problem_.referencePressure(element));
137 fluidState.setTemperature(problem_.temperature(element));
138 fluidState.setSaturation(wPhaseIdx, 1.);
139 fluidState.setSaturation(nPhaseIdx, 0.);
146 velocity_.initialize();
174 return calcVelocityInTransport_;
187 template<
class MultiWriter>
191 velocity_.addOutputVtkFields(writer);
198 Scalar density_[numPhases];
199 Scalar viscosity_[numPhases];
200 bool calcVelocityInTransport_;
203 static const int pressureType_ = getPropValue<TypeTag, Properties::PressureFormulation>();
205 static const int saturationType_ = getPropValue<TypeTag, Properties::SaturationFormulation>();
214template<
class TypeTag>
218 for (
const auto& vertex : vertices(problem_.gridView()))
220 int vIdxGlobal = problem_.variables().index(vertex);
222 InteractionVolume& interactionVolume = this->interactionVolumes_.interactionVolume(vIdxGlobal);
225 if (interactionVolume.isInnerVolume())
229 for (
int i = 0; i < 8; i++)
231 eIdxGlobal[i] = problem_.variables().index(interactionVolume.getSubVolumeElement(i));
235 CellData & cellData1 = problem_.variables().cellData(eIdxGlobal[0]);
236 CellData & cellData2 = problem_.variables().cellData(eIdxGlobal[1]);
237 CellData & cellData3 = problem_.variables().cellData(eIdxGlobal[2]);
238 CellData & cellData4 = problem_.variables().cellData(eIdxGlobal[3]);
239 CellData & cellData5 = problem_.variables().cellData(eIdxGlobal[4]);
240 CellData & cellData6 = problem_.variables().cellData(eIdxGlobal[5]);
241 CellData & cellData7 = problem_.variables().cellData(eIdxGlobal[6]);
242 CellData & cellData8 = problem_.variables().cellData(eIdxGlobal[7]);
244 velocity_.calculateInnerInteractionVolumeVelocity(interactionVolume,
245 cellData1, cellData2, cellData3, cellData4,
246 cellData5, cellData6, cellData7, cellData8,
247 this->interactionVolumes_, this->transmissibilityCalculator_);
252 for (
int elemIdx = 0; elemIdx < 8; elemIdx++)
254 if (!interactionVolume.hasSubVolumeElement(elemIdx))
258 bool isOutside =
false;
259 for (
int fIdx = 0; fIdx < dim; fIdx++)
261 int intVolFaceIdx = interactionVolume.getFaceIndexFromSubVolume(elemIdx, fIdx);
262 if (interactionVolume.isOutsideFace(intVolFaceIdx))
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 [[maybe_unused]]
int eIdxGlobalI = problem_.variables().index(elementI);
303 int eIdxGlobalJ = problem_.variables().index(elementJ);
305 CellData& cellDataJ = problem_.variables().cellData(eIdxGlobalJ);
307 const auto refElement = referenceElement(elementI);
309 int indexInInside = intersection.indexInInside();
310 int indexInOutside = intersection.indexInOutside();
312 Dune::FieldVector<CellData, 8> cellDataTemp;
314 for (
int vIdx = 0; vIdx < numVertices; vIdx++)
316 int localVertIdx = refElement.subEntity(indexInInside, 1, vIdx, dim);
318 int vIdxGlobal = problem_.variables().index(elementI.template subEntity<dim>(localVertIdx));
320 InteractionVolume& interactionVolume = this->interactionVolumes_.interactionVolume(vIdxGlobal);
322 if (interactionVolume.isInnerVolume())
325 int localMpfaElemIdxI = 0;
326 int localMpfaElemIdxJ = 0;
329 for (
int i = 0; i < 8; i++)
331 auto elem = interactionVolume.getSubVolumeElement(i);
333 if (elem == elementI)
334 localMpfaElemIdxI = i;
335 else if (elem == elementJ)
336 localMpfaElemIdxJ = i;
338 eIdxGlobal[i] = problem_.variables().index(elem);
339 cellDataTemp[i] = problem_.variables().cellData(eIdxGlobal[i]);
342 int mpfaFaceIdx = IndexTranslator::getFaceIndexFromElements(localMpfaElemIdxI, localMpfaElemIdxJ);
346 velocity_.calculateInnerInteractionVolumeVelocity(interactionVolume,
347 cellDataTemp[0], cellDataTemp[1], cellDataTemp[2], cellDataTemp[3],
348 cellDataTemp[4], cellDataTemp[5], cellDataTemp[6], cellDataTemp[7],
349 this->interactionVolumes_, this->transmissibilityCalculator_, mpfaFaceIdx);
352 cellData.fluxData().setVelocity(wPhaseIdx, indexInInside,
353 cellDataTemp[localMpfaElemIdxI].fluxData().velocity(wPhaseIdx, indexInInside));
354 cellData.fluxData().setVelocity(nPhaseIdx, indexInInside,
355 cellDataTemp[localMpfaElemIdxI].fluxData().velocity(nPhaseIdx, indexInInside));
356 cellData.fluxData().setUpwindPotential(wPhaseIdx, indexInInside,
357 cellDataTemp[localMpfaElemIdxI].fluxData().upwindPotential(wPhaseIdx, indexInInside));
358 cellData.fluxData().setUpwindPotential(nPhaseIdx, indexInInside,
359 cellDataTemp[localMpfaElemIdxI].fluxData().upwindPotential(nPhaseIdx, indexInInside));
361 cellDataJ.fluxData().setVelocity(wPhaseIdx, indexInOutside,
362 cellDataTemp[localMpfaElemIdxJ].fluxData().velocity(wPhaseIdx, indexInOutside));
363 cellDataJ.fluxData().setVelocity(nPhaseIdx, indexInOutside,
364 cellDataTemp[localMpfaElemIdxJ].fluxData().velocity(nPhaseIdx, indexInOutside));
365 cellDataJ.fluxData().setUpwindPotential(wPhaseIdx, indexInOutside,
366 cellDataTemp[localMpfaElemIdxJ].fluxData().upwindPotential(wPhaseIdx, indexInOutside));
367 cellDataJ.fluxData().setUpwindPotential(nPhaseIdx, indexInOutside,
368 cellDataTemp[localMpfaElemIdxJ].fluxData().upwindPotential(nPhaseIdx, indexInOutside));
371 cellData.fluxData().setVelocityMarker(indexInInside);
372 cellDataJ.fluxData().setVelocityMarker(indexInOutside);
383template<
class TypeTag>
386 auto element = intersection.inside();
389 int isIndex = intersection.indexInInside();
392 const Dune::FieldVector<Scalar, dim>& unitOuterNormal = intersection.centerUnitOuterNormal();
394 BoundaryTypes bcType;
396 problem_.boundaryTypes(bcType, intersection);
397 PrimaryVariables boundValues(0.0);
399 if (bcType.isDirichlet(pressEqIdx))
401 problem_.dirichlet(boundValues, intersection);
404 const GlobalPosition& globalPosI = element.geometry().center();
407 const GlobalPosition& globalPosJ = intersection.geometry().center();
410 Scalar lambdaWI = cellData.mobility(wPhaseIdx);
411 Scalar lambdaNwI = cellData.mobility(nPhaseIdx);
414 Scalar pcI = cellData.capillaryPressure();
417 GlobalPosition distVec = globalPosJ - globalPosI;
420 Scalar dist = distVec.two_norm();
424 DimMatrix meanPermeability(0);
426 problem_.spatialParams().meanK(meanPermeability, problem_.spatialParams().intrinsicPermeability(element));
433 if (bcType.isDirichlet(satEqIdx))
435 switch (saturationType_)
439 satW = boundValues[saturationIdx];
444 satW = 1 - boundValues[saturationIdx];
451 satW = cellData.saturation(wPhaseIdx);
454 const Scalar pressBound = boundValues[pressureIdx];
456 const auto fluidMatrixInteraction = problem_.spatialParams().fluidMatrixInteractionAtPos(element.geometry().center());
457 const Scalar pcBound = fluidMatrixInteraction.pc(satW);
460 Scalar pressWBound = 0;
461 Scalar pressNwBound = 0;
462 if (pressureType_ == pw)
464 pressWBound = pressBound;
465 pressNwBound = pressBound + pcBound;
467 else if (pressureType_ == pn)
469 pressWBound = pressBound - pcBound;
470 pressNwBound = pressBound;
473 const Scalar lambdaWBound = fluidMatrixInteraction.krw(satW) / viscosity_[wPhaseIdx];
474 const Scalar lambdaNwBound = fluidMatrixInteraction.krn(satW) / viscosity_[nPhaseIdx];
476 Scalar potentialDiffW = cellData.fluxData().upwindPotential(wPhaseIdx, isIndex);
477 Scalar potentialDiffNw = cellData.fluxData().upwindPotential(nPhaseIdx, isIndex);
480 potentialDiffW = (cellData.pressure(wPhaseIdx) - pressWBound);
481 potentialDiffNw = (cellData.pressure(nPhaseIdx) - pressNwBound);
483 potentialDiffW += density_[wPhaseIdx] * (distVec * problem_.gravity());
484 potentialDiffNw += density_[nPhaseIdx] * (distVec * problem_.gravity());
487 cellData.fluxData().setUpwindPotential(wPhaseIdx, isIndex, potentialDiffW);
488 cellData.fluxData().setUpwindPotential(nPhaseIdx, isIndex, potentialDiffNw);
491 Scalar lambdaW = (potentialDiffW > 0.) ? lambdaWI : lambdaWBound;
492 lambdaW = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffW, 0.0, 1.0e-30)) ? 0.5 * (lambdaWI + lambdaWBound) : lambdaW;
493 Scalar lambdaNw = (potentialDiffNw > 0.) ? lambdaNwI : lambdaNwBound;
494 lambdaNw = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffNw, 0.0, 1.0e-30)) ? 0.5 * (lambdaNwI + lambdaNwBound) : lambdaNw;
500 Dune::FieldVector<Scalar, dimWorld> velocityW(unitOuterNormal);
501 Dune::FieldVector<Scalar, dimWorld> velocityNw(unitOuterNormal);
505 Scalar areaScaling = (unitOuterNormal * distVec);
507 Scalar gravityTermW = (problem_.gravity() * distVec) * density_[wPhaseIdx] * areaScaling;
508 Scalar gravityTermNw = (problem_.gravity() * distVec) * density_[nPhaseIdx] * areaScaling;
511 switch (pressureType_)
515 velocityW *= lambdaW * scalarPerm * ((cellData.pressure(wPhaseIdx) - pressBound) / dist + gravityTermW);
516 velocityNw *= lambdaNw * scalarPerm * ((cellData.pressure(wPhaseIdx) - pressBound) / dist + gravityTermNw)
517 + 0.5 * (lambdaNwI + lambdaNwBound) * scalarPerm * (pcI - pcBound) / dist;
522 velocityW *= lambdaW * scalarPerm * ((cellData.pressure(nPhaseIdx) - pressBound) / dist + gravityTermW)
523 - 0.5 * (lambdaWI + lambdaWBound) * scalarPerm * (pcI - pcBound) / dist;
524 velocityNw *= lambdaNw * scalarPerm * ((cellData.pressure(nPhaseIdx) - pressBound) / dist + gravityTermNw);
530 cellData.fluxData().setVelocity(wPhaseIdx, isIndex, velocityW);
531 cellData.fluxData().setVelocity(nPhaseIdx, isIndex, velocityNw);
532 cellData.fluxData().setVelocityMarker(isIndex);
536 else if (bcType.isNeumann(pressEqIdx))
538 problem_.neumann(boundValues, intersection);
540 Dune::FieldVector<Scalar, dimWorld> velocityW(unitOuterNormal);
541 Dune::FieldVector<Scalar, dimWorld> velocityNw(unitOuterNormal);
543 velocityW *= boundValues[wPhaseIdx];
544 velocityNw *= boundValues[nPhaseIdx];
546 velocityW /= density_[wPhaseIdx];
547 velocityNw /= density_[nPhaseIdx];
550 cellData.fluxData().setUpwindPotential(wPhaseIdx, isIndex, boundValues[wPhaseIdx]);
551 cellData.fluxData().setUpwindPotential(nPhaseIdx, isIndex, boundValues[nPhaseIdx]);
553 cellData.fluxData().setVelocity(wPhaseIdx, isIndex, velocityW);
554 cellData.fluxData().setVelocity(nPhaseIdx, isIndex, velocityNw);
555 cellData.fluxData().setVelocityMarker(isIndex);
559 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
Definition: propertysystem.hh:141
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:150
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:74
void update()
Pressure update.
Definition: 3dpressure.hh:271
GetPropType< TypeTag, Properties::MPFAInteractionVolume > InteractionVolume
Type for storing interaction volume information.
Definition: 3dpressure.hh:154
void storePressureSolution()
Globally stores the pressure solution.
Definition: 3dpressure.hh:205
void addOutputVtkFields(MultiWriter &writer)
Adds pressure output to the output file.
Definition: 3dpressure.hh:384
void updateMaterialLaws()
constitutive functions are initialized and stored in the variables object
Definition: 3dpressure.hh:2478
Class for the calculation of 3d velocities from the pressure solution of an IMPES scheme using a MPFA...
Definition: 3dpressurevelocity.hh:44
void update()
Pressure and velocity update.
Definition: 3dpressurevelocity.hh:159
FvMpfaL3dPressureVelocity2p(Problem &problem)
Constructs a FvMpfaL3dPressureVelocity2p object.
Definition: 3dpressurevelocity.hh:96
bool calculateVelocityInTransport()
Indicates if velocity is reconstructed in the pressure step or in the transport step.
Definition: 3dpressurevelocity.hh:172
void updateVelocity()
Function for updating the velocity field if iterations are necessary in the transport solution.
Definition: 3dpressurevelocity.hh:116
void calculateVelocityOnBoundary(const Intersection &intersection, CellData &cellData)
Calculates the velocity at a boundary.
Definition: 3dpressurevelocity.hh:384
void initialize(bool solveTwice=true)
Initializes pressure and velocity.
Definition: 3dpressurevelocity.hh:131
void addOutputVtkFields(MultiWriter &writer)
Adds velocity output to the output file.
Definition: 3dpressurevelocity.hh:188
void calculateVelocity()
Calculates the velocities at a cell-cell interfaces for the entire simulation grid.
Definition: 3dpressurevelocity.hh:215
Class for calculating 3-d velocities from cell-wise constant pressure values.
Definition: 3dvelocity.hh:55
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