3.3.0
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
32
33namespace Dumux {
45template<class TypeTag> class FvMpfaL3dPressureVelocity2p: public FvMpfaL3dPressure2p<TypeTag>
46{
48
51
52 enum
53 {
54 dim = GridView::dimension, dimWorld = GridView::dimensionworld
55 };
56
58
59 enum
60 {
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>()
72 };
73
80 using PrimaryVariables = typename SolutionTypes::PrimaryVariables;
82
84 using Intersection = typename GridView::Intersection;
85
86 using Element = typename GridView::template Codim<0>::Entity;
87
88 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
89
90 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
91
92public:
98 FvMpfaL3dPressureVelocity2p(Problem& problem) :
99 ParentType(problem), problem_(problem), velocity_(problem)
100 {
101 density_[wPhaseIdx] = 0.;
102 density_[nPhaseIdx] = 0.;
103 viscosity_[wPhaseIdx] = 0.;
104 viscosity_[nPhaseIdx] = 0.;
105
106 calcVelocityInTransport_ = getParam<bool>("MPFA.CalcVelocityInTransport");
107 }
108
109 void calculateVelocity();
110
111public:
112
113 // Calculates the velocity at a cell-cell interface.
114 void calculateVelocity(const Intersection&, CellData&);
115 void calculateVelocityOnBoundary(const Intersection& intersection, CellData& cellData);
116
119 {
120 this->updateMaterialLaws();
121
122 this->storePressureSolution();
123
126 }
127
133 void initialize(bool solveTwice = true)
134 {
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.);
142 density_[wPhaseIdx] = FluidSystem::density(fluidState, wPhaseIdx);
143 density_[nPhaseIdx] = FluidSystem::density(fluidState, nPhaseIdx);
144 viscosity_[wPhaseIdx] = FluidSystem::viscosity(fluidState, wPhaseIdx);
145 viscosity_[nPhaseIdx] = FluidSystem::viscosity(fluidState, nPhaseIdx);
146
148 velocity_.initialize();
149
152
153 return;
154 }
155
161 void update()
162 {
164
167 }
168
175 {
176 return calcVelocityInTransport_;
177 }
178
189 template<class MultiWriter>
190 void addOutputVtkFields(MultiWriter &writer)
191 {
193 velocity_.addOutputVtkFields(writer);
194 }
195
196private:
197 Problem& problem_;
199
200 Scalar density_[numPhases];
201 Scalar viscosity_[numPhases];
202 bool calcVelocityInTransport_;
203
205 static const int pressureType_ = getPropValue<TypeTag, Properties::PressureFormulation>();
207 static const int saturationType_ = getPropValue<TypeTag, Properties::SaturationFormulation>();
208};
209// end of template
210
216template<class TypeTag>
218{
219 // run through all vertices
220 for (const auto& vertex : vertices(problem_.gridView()))
221 {
222 int vIdxGlobal = problem_.variables().index(vertex);
223
224 InteractionVolume& interactionVolume = this->interactionVolumes_.interactionVolume(vIdxGlobal);
225
226 // inner interactionvolume
227 if (interactionVolume.isInnerVolume())
228 {
229 // cell index
230 int eIdxGlobal[8];
231 for (int i = 0; i < 8; i++)
232 {
233 eIdxGlobal[i] = problem_.variables().index(interactionVolume.getSubVolumeElement(i));
234 }
235
236 //get the cell Data
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]);
245
246 velocity_.calculateInnerInteractionVolumeVelocity(interactionVolume,
247 cellData1, cellData2, cellData3, cellData4,
248 cellData5, cellData6, cellData7, cellData8,
249 this->interactionVolumes_, this->transmissibilityCalculator_);
250 }
251 // at least one face on boundary! (boundary interactionvolume)
252 else
253 {
254 for (int elemIdx = 0; elemIdx < 8; elemIdx++)
255 {
256 if (!interactionVolume.hasSubVolumeElement(elemIdx))
257 {
258 continue;
259 }
260 bool isOutside = false;
261 for (int fIdx = 0; fIdx < dim; fIdx++)
262 {
263 int intVolFaceIdx = interactionVolume.getFaceIndexFromSubVolume(elemIdx, fIdx);
264 if (interactionVolume.isOutsideFace(intVolFaceIdx))
265 {
266 isOutside = true;
267 break;
268 }
269 }
270 if (isOutside)
271 {
272 continue;
273 }
274 // cell index
275 int eIdxGlobal = problem_.variables().index(interactionVolume.getSubVolumeElement(elemIdx));
276 //get the cell Data
277 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
278
279 velocity_.calculateBoundaryInteractionVolumeVelocity(interactionVolume, cellData, elemIdx);
280 }
281 } // end boundaries
282
283 } // end vertex iterator
284
285 return;
286}
287
296template<class TypeTag>
297void FvMpfaL3dPressureVelocity2p<TypeTag>::calculateVelocity(const Intersection& intersection, CellData& cellData)
298{
299 int numVertices = intersection.geometry().corners();
300
301 auto elementI = intersection.inside();
302 auto elementJ = intersection.outside();
303
304 int eIdxGlobalI DUNE_UNUSED = problem_.variables().index(elementI);
305 int eIdxGlobalJ = problem_.variables().index(elementJ);
306
307 CellData& cellDataJ = problem_.variables().cellData(eIdxGlobalJ);
308
309 const auto refElement = referenceElement(elementI);
310
311 int indexInInside = intersection.indexInInside();
312 int indexInOutside = intersection.indexInOutside();
313
314 Dune::FieldVector<CellData, 8> cellDataTemp;
315
316 for (int vIdx = 0; vIdx < numVertices; vIdx++)
317 {
318 int localVertIdx = refElement.subEntity(indexInInside, 1, vIdx, dim);
319
320 int vIdxGlobal = problem_.variables().index(elementI.template subEntity<dim>(localVertIdx));
321
322 InteractionVolume& interactionVolume = this->interactionVolumes_.interactionVolume(vIdxGlobal);
323
324 if (interactionVolume.isInnerVolume())
325 {
326 // cell index
327 int localMpfaElemIdxI = 0;
328 int localMpfaElemIdxJ = 0;
329
330 int eIdxGlobal[8];
331 for (int i = 0; i < 8; i++)
332 {
333 auto elem = interactionVolume.getSubVolumeElement(i);
334
335 if (elem == elementI)
336 localMpfaElemIdxI = i;
337 else if (elem == elementJ)
338 localMpfaElemIdxJ = i;
339
340 eIdxGlobal[i] = problem_.variables().index(elem);
341 cellDataTemp[i] = problem_.variables().cellData(eIdxGlobal[i]);
342 }
343
344 int mpfaFaceIdx = IndexTranslator::getFaceIndexFromElements(localMpfaElemIdxI, localMpfaElemIdxJ);
345
346 // std::cout<<"mpfaFaceIdx = "<<mpfaFaceIdx<<"\n";
347
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);
352
353
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));
362
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));
371 }
372 }
373 cellData.fluxData().setVelocityMarker(indexInInside);
374 cellDataJ.fluxData().setVelocityMarker(indexInOutside);
375}
376
385template<class TypeTag>
386void FvMpfaL3dPressureVelocity2p<TypeTag>::calculateVelocityOnBoundary(const Intersection& intersection, CellData& cellData)
387{
388 auto element = intersection.inside();
389
390 //get face index
391 int isIndex = intersection.indexInInside();
392
393 //get face normal
394 const Dune::FieldVector<Scalar, dim>& unitOuterNormal = intersection.centerUnitOuterNormal();
395
396 BoundaryTypes bcType;
397 //get boundary type
398 problem_.boundaryTypes(bcType, intersection);
399 PrimaryVariables boundValues(0.0);
400
401 if (bcType.isDirichlet(pressEqIdx))
402 {
403 problem_.dirichlet(boundValues, intersection);
404
405 // get global coordinates of cell centers
406 const GlobalPosition& globalPosI = element.geometry().center();
407
408 // center of face in global coordinates
409 const GlobalPosition& globalPosJ = intersection.geometry().center();
410
411 // get mobilities and fractional flow factors
412 Scalar lambdaWI = cellData.mobility(wPhaseIdx);
413 Scalar lambdaNwI = cellData.mobility(nPhaseIdx);
414
415 // get capillary pressure
416 Scalar pcI = cellData.capillaryPressure();
417
418 // distance vector between barycenters
419 GlobalPosition distVec = globalPosJ - globalPosI;
420
421 // compute distance between cell centers
422 Scalar dist = distVec.two_norm();
423
424 //permeability vector at boundary
425 // compute vectorized permeabilities
426 DimMatrix meanPermeability(0);
427
428 problem_.spatialParams().meanK(meanPermeability, problem_.spatialParams().intrinsicPermeability(element));
429
430 Dune::FieldVector<Scalar, dim> permeability(0);
431 meanPermeability.mv(unitOuterNormal, permeability);
432
433 //determine saturation at the boundary -> if no saturation is known directly at the boundary use the cell saturation
434 Scalar satW = 0;
435 if (bcType.isDirichlet(satEqIdx))
436 {
437 switch (saturationType_)
438 {
439 case sw:
440 {
441 satW = boundValues[saturationIdx];
442 break;
443 }
444 case sn:
445 {
446 satW = 1 - boundValues[saturationIdx];
447 break;
448 }
449 }
450 }
451 else
452 {
453 satW = cellData.saturation(wPhaseIdx);
454 }
455
456 const Scalar pressBound = boundValues[pressureIdx];
457
458 // old material law interface is deprecated: Replace this by
459 // const auto& fluidMatrixInteraction = spatialParams.fluidMatrixInteractionAtPos(element.geometry().center());
460 // after the release of 3.3, when the deprecated interface is no longer supported
461 const auto fluidMatrixInteraction = Deprecated::makePcKrSw(Scalar{}, problem_.spatialParams(), element);
462
463 const Scalar pcBound = fluidMatrixInteraction.pc(satW);
464
465 //determine phase pressures from primary pressure variable
466 Scalar pressWBound = 0;
467 Scalar pressNwBound = 0;
468 if (pressureType_ == pw)
469 {
470 pressWBound = pressBound;
471 pressNwBound = pressBound + pcBound;
472 }
473 else if (pressureType_ == pn)
474 {
475 pressWBound = pressBound - pcBound;
476 pressNwBound = pressBound;
477 }
478
479 const Scalar lambdaWBound = fluidMatrixInteraction.krw(satW) / viscosity_[wPhaseIdx];
480 const Scalar lambdaNwBound = fluidMatrixInteraction.krn(satW) / viscosity_[nPhaseIdx];
481
482 Scalar potentialDiffW = cellData.fluxData().upwindPotential(wPhaseIdx, isIndex);
483 Scalar potentialDiffNw = cellData.fluxData().upwindPotential(nPhaseIdx, isIndex);
484
485 //calculate potential gradient
486 potentialDiffW = (cellData.pressure(wPhaseIdx) - pressWBound);
487 potentialDiffNw = (cellData.pressure(nPhaseIdx) - pressNwBound);
488
489 potentialDiffW += density_[wPhaseIdx] * (distVec * problem_.gravity());
490 potentialDiffNw += density_[nPhaseIdx] * (distVec * problem_.gravity());
491
492 //store potential gradients for further calculations
493 cellData.fluxData().setUpwindPotential(wPhaseIdx, isIndex, potentialDiffW);
494 cellData.fluxData().setUpwindPotential(nPhaseIdx, isIndex, potentialDiffNw);
495
496 //do the upwinding of the mobility depending on the phase potentials
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;
501
502
503 Scalar scalarPerm = permeability.two_norm();
504
505 //calculate the gravity term
506 Dune::FieldVector<Scalar, dimWorld> velocityW(unitOuterNormal);
507 Dune::FieldVector<Scalar, dimWorld> velocityNw(unitOuterNormal);
508
509 //calculate unit distVec
510 distVec /= dist;
511 Scalar areaScaling = (unitOuterNormal * distVec);
512 //this treatment of g allows to account for gravity flux through faces where the face normal has no z component (e.g. parallelepiped grids)
513 Scalar gravityTermW = (problem_.gravity() * distVec) * density_[wPhaseIdx] * areaScaling;
514 Scalar gravityTermNw = (problem_.gravity() * distVec) * density_[nPhaseIdx] * areaScaling;
515
516 //calculate velocity depending on the pressure used -> use pc = pn - pw
517 switch (pressureType_)
518 {
519 case pw:
520 {
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;
524 break;
525 }
526 case pn:
527 {
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);
531 break;
532 }
533 }
534
535 //store velocities
536 cellData.fluxData().setVelocity(wPhaseIdx, isIndex, velocityW);
537 cellData.fluxData().setVelocity(nPhaseIdx, isIndex, velocityNw);
538 cellData.fluxData().setVelocityMarker(isIndex);
539
540 } //end dirichlet boundary
541
542 else if (bcType.isNeumann(pressEqIdx))
543 {
544 problem_.neumann(boundValues, intersection);
545
546 Dune::FieldVector<Scalar, dimWorld> velocityW(unitOuterNormal);
547 Dune::FieldVector<Scalar, dimWorld> velocityNw(unitOuterNormal);
548
549 velocityW *= boundValues[wPhaseIdx];
550 velocityNw *= boundValues[nPhaseIdx];
551
552 velocityW /= density_[wPhaseIdx];
553 velocityNw /= density_[nPhaseIdx];
554
555 //store potential gradients for further calculations
556 cellData.fluxData().setUpwindPotential(wPhaseIdx, isIndex, boundValues[wPhaseIdx]);
557 cellData.fluxData().setUpwindPotential(nPhaseIdx, isIndex, boundValues[nPhaseIdx]);
558
559 cellData.fluxData().setVelocity(wPhaseIdx, isIndex, velocityW);
560 cellData.fluxData().setVelocity(nPhaseIdx, isIndex, velocityNw);
561 cellData.fluxData().setVelocityMarker(isIndex);
562 } //end neumann boundary
563 else
564 {
565 DUNE_THROW(Dune::NotImplemented, "No valid boundary condition type defined for pressure equation!");
566 }
567}
568
569} // end namespace Dumux
570#endif
Helpers for deprecation.
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 (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