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;
80 using MaterialLaw =
typename SpatialParams::MaterialLaw;
82 using ReferenceElements = Dune::ReferenceElements<Scalar, dim>;
85 using Intersection =
typename GridView::Intersection;
87 using Element =
typename GridView::template Codim<0>::Entity;
89 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
91 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
100 ParentType(problem), problem_(problem), velocity_(problem)
102 density_[wPhaseIdx] = 0.;
103 density_[nPhaseIdx] = 0.;
104 viscosity_[wPhaseIdx] = 0.;
105 viscosity_[nPhaseIdx] = 0.;
107 calcVelocityInTransport_ = getParam<bool>(
"MPFA.CalcVelocityInTransport");
136 const auto element = *problem_.gridView().template begin<0>();
137 FluidState fluidState;
138 fluidState.setPressure(wPhaseIdx, problem_.referencePressure(element));
139 fluidState.setPressure(nPhaseIdx, problem_.referencePressure(element));
140 fluidState.setTemperature(problem_.temperature(element));
141 fluidState.setSaturation(wPhaseIdx, 1.);
142 fluidState.setSaturation(nPhaseIdx, 0.);
149 velocity_.initialize();
177 return calcVelocityInTransport_;
190 template<
class MultiWriter>
194 velocity_.addOutputVtkFields(writer);
201 Scalar density_[numPhases];
202 Scalar viscosity_[numPhases];
203 bool calcVelocityInTransport_;
206 static const int pressureType_ = getPropValue<TypeTag, Properties::PressureFormulation>();
208 static const int saturationType_ = getPropValue<TypeTag, Properties::SaturationFormulation>();
217template<
class TypeTag>
221 for (
const auto& vertex : vertices(problem_.gridView()))
223 int vIdxGlobal = problem_.variables().index(vertex);
225 InteractionVolume& interactionVolume = this->interactionVolumes_.interactionVolume(vIdxGlobal);
228 if (interactionVolume.isInnerVolume())
232 for (
int i = 0; i < 8; i++)
234 eIdxGlobal[i] = problem_.variables().index(interactionVolume.getSubVolumeElement(i));
238 CellData & cellData1 = problem_.variables().cellData(eIdxGlobal[0]);
239 CellData & cellData2 = problem_.variables().cellData(eIdxGlobal[1]);
240 CellData & cellData3 = problem_.variables().cellData(eIdxGlobal[2]);
241 CellData & cellData4 = problem_.variables().cellData(eIdxGlobal[3]);
242 CellData & cellData5 = problem_.variables().cellData(eIdxGlobal[4]);
243 CellData & cellData6 = problem_.variables().cellData(eIdxGlobal[5]);
244 CellData & cellData7 = problem_.variables().cellData(eIdxGlobal[6]);
245 CellData & cellData8 = problem_.variables().cellData(eIdxGlobal[7]);
247 velocity_.calculateInnerInteractionVolumeVelocity(interactionVolume,
248 cellData1, cellData2, cellData3, cellData4,
249 cellData5, cellData6, cellData7, cellData8,
250 this->interactionVolumes_, this->transmissibilityCalculator_);
255 for (
int elemIdx = 0; elemIdx < 8; elemIdx++)
257 if (!interactionVolume.hasSubVolumeElement(elemIdx))
261 bool isOutside =
false;
262 for (
int fIdx = 0; fIdx < dim; fIdx++)
264 int intVolFaceIdx = interactionVolume.getFaceIndexFromSubVolume(elemIdx, fIdx);
265 if (interactionVolume.isOutsideFace(intVolFaceIdx))
276 int eIdxGlobal = problem_.variables().index(interactionVolume.getSubVolumeElement(elemIdx));
278 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
280 velocity_.calculateBoundaryInteractionVolumeVelocity(interactionVolume, cellData, elemIdx);
297template<
class TypeTag>
300 int numVertices = intersection.geometry().corners();
302 auto elementI = intersection.inside();
303 auto elementJ = intersection.outside();
305 int eIdxGlobalI DUNE_UNUSED = problem_.variables().index(elementI);
306 int eIdxGlobalJ = problem_.variables().index(elementJ);
308 CellData& cellDataJ = problem_.variables().cellData(eIdxGlobalJ);
310 const auto referenceElement = ReferenceElements::general(elementI.type());
312 int indexInInside = intersection.indexInInside();
313 int indexInOutside = intersection.indexInOutside();
315 Dune::FieldVector<CellData, 8> cellDataTemp;
317 for (
int vIdx = 0; vIdx < numVertices; vIdx++)
319 int localVertIdx = referenceElement.subEntity(indexInInside, 1, vIdx, dim);
321 int vIdxGlobal = problem_.variables().index(elementI.template subEntity<dim>(localVertIdx));
323 InteractionVolume& interactionVolume = this->interactionVolumes_.interactionVolume(vIdxGlobal);
325 if (interactionVolume.isInnerVolume())
328 int localMpfaElemIdxI = 0;
329 int localMpfaElemIdxJ = 0;
332 for (
int i = 0; i < 8; i++)
334 auto elem = interactionVolume.getSubVolumeElement(i);
336 if (elem == elementI)
337 localMpfaElemIdxI = i;
338 else if (elem == elementJ)
339 localMpfaElemIdxJ = i;
341 eIdxGlobal[i] = problem_.variables().index(elem);
342 cellDataTemp[i] = problem_.variables().cellData(eIdxGlobal[i]);
345 int mpfaFaceIdx = IndexTranslator::getFaceIndexFromElements(localMpfaElemIdxI, localMpfaElemIdxJ);
349 velocity_.calculateInnerInteractionVolumeVelocity(interactionVolume,
350 cellDataTemp[0], cellDataTemp[1], cellDataTemp[2], cellDataTemp[3],
351 cellDataTemp[4], cellDataTemp[5], cellDataTemp[6], cellDataTemp[7],
352 this->interactionVolumes_, this->transmissibilityCalculator_, mpfaFaceIdx);
355 cellData.fluxData().setVelocity(wPhaseIdx, indexInInside,
356 cellDataTemp[localMpfaElemIdxI].fluxData().velocity(wPhaseIdx, indexInInside));
357 cellData.fluxData().setVelocity(nPhaseIdx, indexInInside,
358 cellDataTemp[localMpfaElemIdxI].fluxData().velocity(nPhaseIdx, indexInInside));
359 cellData.fluxData().setUpwindPotential(wPhaseIdx, indexInInside,
360 cellDataTemp[localMpfaElemIdxI].fluxData().upwindPotential(wPhaseIdx, indexInInside));
361 cellData.fluxData().setUpwindPotential(nPhaseIdx, indexInInside,
362 cellDataTemp[localMpfaElemIdxI].fluxData().upwindPotential(nPhaseIdx, indexInInside));
364 cellDataJ.fluxData().setVelocity(wPhaseIdx, indexInOutside,
365 cellDataTemp[localMpfaElemIdxJ].fluxData().velocity(wPhaseIdx, indexInOutside));
366 cellDataJ.fluxData().setVelocity(nPhaseIdx, indexInOutside,
367 cellDataTemp[localMpfaElemIdxJ].fluxData().velocity(nPhaseIdx, indexInOutside));
368 cellDataJ.fluxData().setUpwindPotential(wPhaseIdx, indexInOutside,
369 cellDataTemp[localMpfaElemIdxJ].fluxData().upwindPotential(wPhaseIdx, indexInOutside));
370 cellDataJ.fluxData().setUpwindPotential(nPhaseIdx, indexInOutside,
371 cellDataTemp[localMpfaElemIdxJ].fluxData().upwindPotential(nPhaseIdx, indexInOutside));
374 cellData.fluxData().setVelocityMarker(indexInInside);
375 cellDataJ.fluxData().setVelocityMarker(indexInOutside);
386template<
class TypeTag>
389 auto element = intersection.inside();
392 int isIndex = intersection.indexInInside();
395 const Dune::FieldVector<Scalar, dim>& unitOuterNormal = intersection.centerUnitOuterNormal();
397 BoundaryTypes bcType;
399 problem_.boundaryTypes(bcType, intersection);
400 PrimaryVariables boundValues(0.0);
402 if (bcType.isDirichlet(pressEqIdx))
404 problem_.dirichlet(boundValues, intersection);
407 const GlobalPosition& globalPosI = element.geometry().center();
410 const GlobalPosition& globalPosJ = intersection.geometry().center();
413 Scalar lambdaWI = cellData.mobility(wPhaseIdx);
414 Scalar lambdaNwI = cellData.mobility(nPhaseIdx);
417 Scalar pcI = cellData.capillaryPressure();
420 GlobalPosition distVec = globalPosJ - globalPosI;
423 Scalar dist = distVec.two_norm();
427 DimMatrix meanPermeability(0);
429 problem_.spatialParams().meanK(meanPermeability, problem_.spatialParams().intrinsicPermeability(element));
436 if (bcType.isDirichlet(satEqIdx))
438 switch (saturationType_)
442 satW = boundValues[saturationIdx];
447 satW = 1 - boundValues[saturationIdx];
454 satW = cellData.saturation(wPhaseIdx);
457 Scalar pressBound = boundValues[pressureIdx];
458 Scalar pcBound = MaterialLaw::pc(problem_.spatialParams().materialLawParams(element), satW);
461 Scalar pressWBound = 0;
462 Scalar pressNwBound = 0;
463 if (pressureType_ == pw)
465 pressWBound = pressBound;
466 pressNwBound = pressBound + pcBound;
468 else if (pressureType_ == pn)
470 pressWBound = pressBound - pcBound;
471 pressNwBound = pressBound;
474 Scalar lambdaWBound = MaterialLaw::krw(problem_.spatialParams().materialLawParams(element), satW)
475 / viscosity_[wPhaseIdx];
476 Scalar lambdaNwBound = MaterialLaw::krn(problem_.spatialParams().materialLawParams(element), satW)
477 / viscosity_[nPhaseIdx];
479 Scalar potentialDiffW = cellData.fluxData().upwindPotential(wPhaseIdx, isIndex);
480 Scalar potentialDiffNw = cellData.fluxData().upwindPotential(nPhaseIdx, isIndex);
483 potentialDiffW = (cellData.pressure(wPhaseIdx) - pressWBound);
484 potentialDiffNw = (cellData.pressure(nPhaseIdx) - pressNwBound);
486 potentialDiffW += density_[wPhaseIdx] * (distVec * problem_.gravity());
487 potentialDiffNw += density_[nPhaseIdx] * (distVec * problem_.gravity());
490 cellData.fluxData().setUpwindPotential(wPhaseIdx, isIndex, potentialDiffW);
491 cellData.fluxData().setUpwindPotential(nPhaseIdx, isIndex, potentialDiffNw);
494 Scalar lambdaW = (potentialDiffW > 0.) ? lambdaWI : lambdaWBound;
495 lambdaW = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffW, 0.0, 1.0e-30)) ? 0.5 * (lambdaWI + lambdaWBound) : lambdaW;
496 Scalar lambdaNw = (potentialDiffNw > 0.) ? lambdaNwI : lambdaNwBound;
497 lambdaNw = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffNw, 0.0, 1.0e-30)) ? 0.5 * (lambdaNwI + lambdaNwBound) : lambdaNw;
503 Dune::FieldVector<Scalar, dimWorld> velocityW(unitOuterNormal);
504 Dune::FieldVector<Scalar, dimWorld> velocityNw(unitOuterNormal);
508 Scalar areaScaling = (unitOuterNormal * distVec);
510 Scalar gravityTermW = (problem_.gravity() * distVec) * density_[wPhaseIdx] * areaScaling;
511 Scalar gravityTermNw = (problem_.gravity() * distVec) * density_[nPhaseIdx] * areaScaling;
514 switch (pressureType_)
518 velocityW *= lambdaW * scalarPerm * ((cellData.pressure(wPhaseIdx) - pressBound) / dist + gravityTermW);
519 velocityNw *= lambdaNw * scalarPerm * ((cellData.pressure(wPhaseIdx) - pressBound) / dist + gravityTermNw)
520 + 0.5 * (lambdaNwI + lambdaNwBound) * scalarPerm * (pcI - pcBound) / dist;
525 velocityW *= lambdaW * scalarPerm * ((cellData.pressure(nPhaseIdx) - pressBound) / dist + gravityTermW)
526 - 0.5 * (lambdaWI + lambdaWBound) * scalarPerm * (pcI - pcBound) / dist;
527 velocityNw *= lambdaNw * scalarPerm * ((cellData.pressure(nPhaseIdx) - pressBound) / dist + gravityTermNw);
533 cellData.fluxData().setVelocity(wPhaseIdx, isIndex, velocityW);
534 cellData.fluxData().setVelocity(nPhaseIdx, isIndex, velocityNw);
535 cellData.fluxData().setVelocityMarker(isIndex);
539 else if (bcType.isNeumann(pressEqIdx))
541 problem_.neumann(boundValues, intersection);
543 Dune::FieldVector<Scalar, dimWorld> velocityW(unitOuterNormal);
544 Dune::FieldVector<Scalar, dimWorld> velocityNw(unitOuterNormal);
546 velocityW *= boundValues[wPhaseIdx];
547 velocityNw *= boundValues[nPhaseIdx];
549 velocityW /= density_[wPhaseIdx];
550 velocityNw /= density_[nPhaseIdx];
553 cellData.fluxData().setUpwindPotential(wPhaseIdx, isIndex, boundValues[wPhaseIdx]);
554 cellData.fluxData().setUpwindPotential(nPhaseIdx, isIndex, boundValues[nPhaseIdx]);
556 cellData.fluxData().setVelocity(wPhaseIdx, isIndex, velocityW);
557 cellData.fluxData().setVelocity(nPhaseIdx, isIndex, velocityNw);
558 cellData.fluxData().setVelocityMarker(isIndex);
562 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:74
void update()
Pressure update.
Definition: 3dpressure.hh:272
GetPropType< TypeTag, Properties::MPFAInteractionVolume > InteractionVolume
Type for storing interaction volume information.
Definition: 3dpressure.hh:155
void storePressureSolution()
Globally stores the pressure solution.
Definition: 3dpressure.hh:206
void addOutputVtkFields(MultiWriter &writer)
Adds pressure output to the output file.
Definition: 3dpressure.hh:385
void updateMaterialLaws()
constitutive functions are initialized and stored in the variables object
Definition: 3dpressure.hh:2481
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:162
FvMpfaL3dPressureVelocity2p(Problem &problem)
Constructs a FvMpfaL3dPressureVelocity2p object.
Definition: 3dpressurevelocity.hh:99
bool calculateVelocityInTransport()
Indicates if velocity is reconstructed in the pressure step or in the transport step.
Definition: 3dpressurevelocity.hh:175
void updateVelocity()
Function for updating the velocity field if iterations are necessary in the transport solution.
Definition: 3dpressurevelocity.hh:119
void calculateVelocityOnBoundary(const Intersection &intersection, CellData &cellData)
Calculates the velocity at a boundary.
Definition: 3dpressurevelocity.hh:387
void initialize(bool solveTwice=true)
Initializes pressure and velocity.
Definition: 3dpressurevelocity.hh:134
void addOutputVtkFields(MultiWriter &writer)
Adds velocity output to the output file.
Definition: 3dpressurevelocity.hh:191
void calculateVelocity()
Calculates the velocities at a cell-cell interfaces for the entire simulation grid.
Definition: 3dpressurevelocity.hh:218
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