3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
3dpressurevelocity.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*****************************************************************************
4 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_FVMPFAL2PFABOUND3DVELOCITIES2P_HH
25#define DUMUX_FVMPFAL2PFABOUND3DVELOCITIES2P_HH
26
27#include <dune/common/float_cmp.hh>
28#include "3dpressure.hh"
29#include "3dvelocity.hh"
30
31namespace Dumux {
43template<class TypeTag> class FvMpfaL3dPressureVelocity2p: public FvMpfaL3dPressure2p<TypeTag>
44{
46
49
50 enum
51 {
52 dim = GridView::dimension, dimWorld = GridView::dimensionworld
53 };
54
56
57 enum
58 {
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>()
70 };
71
78 using PrimaryVariables = typename SolutionTypes::PrimaryVariables;
80
82 using Intersection = typename GridView::Intersection;
83
84 using Element = typename GridView::template Codim<0>::Entity;
85
86 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
87
88 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
89
90public:
96 FvMpfaL3dPressureVelocity2p(Problem& problem) :
97 ParentType(problem), problem_(problem), velocity_(problem)
98 {
99 density_[wPhaseIdx] = 0.;
100 density_[nPhaseIdx] = 0.;
101 viscosity_[wPhaseIdx] = 0.;
102 viscosity_[nPhaseIdx] = 0.;
103
104 calcVelocityInTransport_ = getParam<bool>("MPFA.CalcVelocityInTransport");
105 }
106
107 void calculateVelocity();
108
109public:
110
111 // Calculates the velocity at a cell-cell interface.
112 void calculateVelocity(const Intersection&, CellData&);
113 void calculateVelocityOnBoundary(const Intersection& intersection, CellData& cellData);
114
117 {
118 this->updateMaterialLaws();
119
120 this->storePressureSolution();
121
124 }
125
131 void initialize(bool solveTwice = true)
132 {
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.);
140 density_[wPhaseIdx] = FluidSystem::density(fluidState, wPhaseIdx);
141 density_[nPhaseIdx] = FluidSystem::density(fluidState, nPhaseIdx);
142 viscosity_[wPhaseIdx] = FluidSystem::viscosity(fluidState, wPhaseIdx);
143 viscosity_[nPhaseIdx] = FluidSystem::viscosity(fluidState, nPhaseIdx);
144
146 velocity_.initialize();
147
150
151 return;
152 }
153
159 void update()
160 {
162
165 }
166
173 {
174 return calcVelocityInTransport_;
175 }
176
187 template<class MultiWriter>
188 void addOutputVtkFields(MultiWriter &writer)
189 {
191 velocity_.addOutputVtkFields(writer);
192 }
193
194private:
195 Problem& problem_;
197
198 Scalar density_[numPhases];
199 Scalar viscosity_[numPhases];
200 bool calcVelocityInTransport_;
201
203 static const int pressureType_ = getPropValue<TypeTag, Properties::PressureFormulation>();
205 static const int saturationType_ = getPropValue<TypeTag, Properties::SaturationFormulation>();
206};
207// end of template
208
214template<class TypeTag>
216{
217 // run through all vertices
218 for (const auto& vertex : vertices(problem_.gridView()))
219 {
220 int vIdxGlobal = problem_.variables().index(vertex);
221
222 InteractionVolume& interactionVolume = this->interactionVolumes_.interactionVolume(vIdxGlobal);
223
224 // inner interactionvolume
225 if (interactionVolume.isInnerVolume())
226 {
227 // cell index
228 int eIdxGlobal[8];
229 for (int i = 0; i < 8; i++)
230 {
231 eIdxGlobal[i] = problem_.variables().index(interactionVolume.getSubVolumeElement(i));
232 }
233
234 //get the cell Data
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]);
243
244 velocity_.calculateInnerInteractionVolumeVelocity(interactionVolume,
245 cellData1, cellData2, cellData3, cellData4,
246 cellData5, cellData6, cellData7, cellData8,
247 this->interactionVolumes_, this->transmissibilityCalculator_);
248 }
249 // at least one face on boundary! (boundary interactionvolume)
250 else
251 {
252 for (int elemIdx = 0; elemIdx < 8; elemIdx++)
253 {
254 if (!interactionVolume.hasSubVolumeElement(elemIdx))
255 {
256 continue;
257 }
258 bool isOutside = false;
259 for (int fIdx = 0; fIdx < dim; fIdx++)
260 {
261 int intVolFaceIdx = interactionVolume.getFaceIndexFromSubVolume(elemIdx, fIdx);
262 if (interactionVolume.isOutsideFace(intVolFaceIdx))
263 {
264 isOutside = true;
265 break;
266 }
267 }
268 if (isOutside)
269 {
270 continue;
271 }
272 // cell index
273 int eIdxGlobal = problem_.variables().index(interactionVolume.getSubVolumeElement(elemIdx));
274 //get the cell Data
275 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
276
277 velocity_.calculateBoundaryInteractionVolumeVelocity(interactionVolume, cellData, elemIdx);
278 }
279 } // end boundaries
280
281 } // end vertex iterator
282
283 return;
284}
285
294template<class TypeTag>
295void FvMpfaL3dPressureVelocity2p<TypeTag>::calculateVelocity(const Intersection& intersection, CellData& cellData)
296{
297 int numVertices = intersection.geometry().corners();
298
299 auto elementI = intersection.inside();
300 auto elementJ = intersection.outside();
301
302 [[maybe_unused]] int eIdxGlobalI = problem_.variables().index(elementI);
303 int eIdxGlobalJ = problem_.variables().index(elementJ);
304
305 CellData& cellDataJ = problem_.variables().cellData(eIdxGlobalJ);
306
307 const auto refElement = referenceElement(elementI);
308
309 int indexInInside = intersection.indexInInside();
310 int indexInOutside = intersection.indexInOutside();
311
312 Dune::FieldVector<CellData, 8> cellDataTemp;
313
314 for (int vIdx = 0; vIdx < numVertices; vIdx++)
315 {
316 int localVertIdx = refElement.subEntity(indexInInside, 1, vIdx, dim);
317
318 int vIdxGlobal = problem_.variables().index(elementI.template subEntity<dim>(localVertIdx));
319
320 InteractionVolume& interactionVolume = this->interactionVolumes_.interactionVolume(vIdxGlobal);
321
322 if (interactionVolume.isInnerVolume())
323 {
324 // cell index
325 int localMpfaElemIdxI = 0;
326 int localMpfaElemIdxJ = 0;
327
328 int eIdxGlobal[8];
329 for (int i = 0; i < 8; i++)
330 {
331 auto elem = interactionVolume.getSubVolumeElement(i);
332
333 if (elem == elementI)
334 localMpfaElemIdxI = i;
335 else if (elem == elementJ)
336 localMpfaElemIdxJ = i;
337
338 eIdxGlobal[i] = problem_.variables().index(elem);
339 cellDataTemp[i] = problem_.variables().cellData(eIdxGlobal[i]);
340 }
341
342 int mpfaFaceIdx = IndexTranslator::getFaceIndexFromElements(localMpfaElemIdxI, localMpfaElemIdxJ);
343
344 // std::cout<<"mpfaFaceIdx = "<<mpfaFaceIdx<<"\n";
345
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);
350
351
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));
360
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));
369 }
370 }
371 cellData.fluxData().setVelocityMarker(indexInInside);
372 cellDataJ.fluxData().setVelocityMarker(indexInOutside);
373}
374
383template<class TypeTag>
384void FvMpfaL3dPressureVelocity2p<TypeTag>::calculateVelocityOnBoundary(const Intersection& intersection, CellData& cellData)
385{
386 auto element = intersection.inside();
387
388 //get face index
389 int isIndex = intersection.indexInInside();
390
391 //get face normal
392 const Dune::FieldVector<Scalar, dim>& unitOuterNormal = intersection.centerUnitOuterNormal();
393
394 BoundaryTypes bcType;
395 //get boundary type
396 problem_.boundaryTypes(bcType, intersection);
397 PrimaryVariables boundValues(0.0);
398
399 if (bcType.isDirichlet(pressEqIdx))
400 {
401 problem_.dirichlet(boundValues, intersection);
402
403 // get global coordinates of cell centers
404 const GlobalPosition& globalPosI = element.geometry().center();
405
406 // center of face in global coordinates
407 const GlobalPosition& globalPosJ = intersection.geometry().center();
408
409 // get mobilities and fractional flow factors
410 Scalar lambdaWI = cellData.mobility(wPhaseIdx);
411 Scalar lambdaNwI = cellData.mobility(nPhaseIdx);
412
413 // get capillary pressure
414 Scalar pcI = cellData.capillaryPressure();
415
416 // distance vector between barycenters
417 GlobalPosition distVec = globalPosJ - globalPosI;
418
419 // compute distance between cell centers
420 Scalar dist = distVec.two_norm();
421
422 //permeability vector at boundary
423 // compute vectorized permeabilities
424 DimMatrix meanPermeability(0);
425
426 problem_.spatialParams().meanK(meanPermeability, problem_.spatialParams().intrinsicPermeability(element));
427
428 Dune::FieldVector<Scalar, dim> permeability(0);
429 meanPermeability.mv(unitOuterNormal, permeability);
430
431 //determine saturation at the boundary -> if no saturation is known directly at the boundary use the cell saturation
432 Scalar satW = 0;
433 if (bcType.isDirichlet(satEqIdx))
434 {
435 switch (saturationType_)
436 {
437 case sw:
438 {
439 satW = boundValues[saturationIdx];
440 break;
441 }
442 case sn:
443 {
444 satW = 1 - boundValues[saturationIdx];
445 break;
446 }
447 }
448 }
449 else
450 {
451 satW = cellData.saturation(wPhaseIdx);
452 }
453
454 const Scalar pressBound = boundValues[pressureIdx];
455
456 const auto fluidMatrixInteraction = problem_.spatialParams().fluidMatrixInteractionAtPos(element.geometry().center());
457 const Scalar pcBound = fluidMatrixInteraction.pc(satW);
458
459 //determine phase pressures from primary pressure variable
460 Scalar pressWBound = 0;
461 Scalar pressNwBound = 0;
462 if (pressureType_ == pw)
463 {
464 pressWBound = pressBound;
465 pressNwBound = pressBound + pcBound;
466 }
467 else if (pressureType_ == pn)
468 {
469 pressWBound = pressBound - pcBound;
470 pressNwBound = pressBound;
471 }
472
473 const Scalar lambdaWBound = fluidMatrixInteraction.krw(satW) / viscosity_[wPhaseIdx];
474 const Scalar lambdaNwBound = fluidMatrixInteraction.krn(satW) / viscosity_[nPhaseIdx];
475
476 Scalar potentialDiffW = cellData.fluxData().upwindPotential(wPhaseIdx, isIndex);
477 Scalar potentialDiffNw = cellData.fluxData().upwindPotential(nPhaseIdx, isIndex);
478
479 //calculate potential gradient
480 potentialDiffW = (cellData.pressure(wPhaseIdx) - pressWBound);
481 potentialDiffNw = (cellData.pressure(nPhaseIdx) - pressNwBound);
482
483 potentialDiffW += density_[wPhaseIdx] * (distVec * problem_.gravity());
484 potentialDiffNw += density_[nPhaseIdx] * (distVec * problem_.gravity());
485
486 //store potential gradients for further calculations
487 cellData.fluxData().setUpwindPotential(wPhaseIdx, isIndex, potentialDiffW);
488 cellData.fluxData().setUpwindPotential(nPhaseIdx, isIndex, potentialDiffNw);
489
490 //do the upwinding of the mobility depending on the phase potentials
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;
495
496
497 Scalar scalarPerm = permeability.two_norm();
498
499 //calculate the gravity term
500 Dune::FieldVector<Scalar, dimWorld> velocityW(unitOuterNormal);
501 Dune::FieldVector<Scalar, dimWorld> velocityNw(unitOuterNormal);
502
503 //calculate unit distVec
504 distVec /= dist;
505 Scalar areaScaling = (unitOuterNormal * distVec);
506 //this treatment of g allows to account for gravity flux through faces where the face normal has no z component (e.g. parallelepiped grids)
507 Scalar gravityTermW = (problem_.gravity() * distVec) * density_[wPhaseIdx] * areaScaling;
508 Scalar gravityTermNw = (problem_.gravity() * distVec) * density_[nPhaseIdx] * areaScaling;
509
510 //calculate velocity depending on the pressure used -> use pc = pn - pw
511 switch (pressureType_)
512 {
513 case pw:
514 {
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;
518 break;
519 }
520 case pn:
521 {
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);
525 break;
526 }
527 }
528
529 //store velocities
530 cellData.fluxData().setVelocity(wPhaseIdx, isIndex, velocityW);
531 cellData.fluxData().setVelocity(nPhaseIdx, isIndex, velocityNw);
532 cellData.fluxData().setVelocityMarker(isIndex);
533
534 } //end dirichlet boundary
535
536 else if (bcType.isNeumann(pressEqIdx))
537 {
538 problem_.neumann(boundValues, intersection);
539
540 Dune::FieldVector<Scalar, dimWorld> velocityW(unitOuterNormal);
541 Dune::FieldVector<Scalar, dimWorld> velocityNw(unitOuterNormal);
542
543 velocityW *= boundValues[wPhaseIdx];
544 velocityNw *= boundValues[nPhaseIdx];
545
546 velocityW /= density_[wPhaseIdx];
547 velocityNw /= density_[nPhaseIdx];
548
549 //store potential gradients for further calculations
550 cellData.fluxData().setUpwindPotential(wPhaseIdx, isIndex, boundValues[wPhaseIdx]);
551 cellData.fluxData().setUpwindPotential(nPhaseIdx, isIndex, boundValues[nPhaseIdx]);
552
553 cellData.fluxData().setVelocity(wPhaseIdx, isIndex, velocityW);
554 cellData.fluxData().setVelocity(nPhaseIdx, isIndex, velocityNw);
555 cellData.fluxData().setVelocityMarker(isIndex);
556 } //end neumann boundary
557 else
558 {
559 DUNE_THROW(Dune::NotImplemented, "No valid boundary condition type defined for pressure equation!");
560 }
561}
562
563} // end namespace Dumux
564#endif
3-d finite Volume-MPFAL implementation of a two-phase pressure equation.
2-d velocity calculation using a 3-d MPFA L-method.
Definition: adapt.hh:29
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