3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
omethod/2dpressurevelocity.hh
Go to the documentation of this file.
1/*****************************************************************************
2 * See the file COPYING for full copying permissions. *
3 * *
4 * This program is free software: you can redistribute it and/or modify *
5 * it under the terms of the GNU General Public License as published by *
6 * the Free Software Foundation, either version 3 of the License, or *
7 * (at your option) any later version. *
8 * *
9 * This program is distributed in the hope that it will be useful, *
10 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
11 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
12 * GNU General Public License for more details. *
13 * *
14 * You should have received a copy of the GNU General Public License *
15 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
16 *****************************************************************************/
22#ifndef DUMUX_MPFAO2DPRESSUREVELOCITIES2P_HH
23#define DUMUX_MPFAO2DPRESSUREVELOCITIES2P_HH
24
25#include <dune/common/float_cmp.hh>
26
27#include "2dpressure.hh"
28#include "2dvelocity.hh"
29
30namespace Dumux {
44template<class TypeTag> class FvMpfaO2dPressureVelocity2p: public FvMpfaO2dPressure2p<TypeTag>
45{
47
49
50 enum
51 {
52 dim = GridView::dimension, dimWorld = GridView::dimensionworld
53 };
54
60 using PrimaryVariables = typename SolutionTypes::PrimaryVariables;
61
63
66
67 using Vertex = typename GridView::Traits::template Codim<dim>::Entity;
68 using Intersection = typename GridView::Intersection;
69
71
73 enum
74 {
75 wPhaseIdx = Indices::wPhaseIdx,
76 nPhaseIdx = Indices::nPhaseIdx,
77 pw = Indices::pressureW,
78 pn = Indices::pressureNw,
79 vw = Indices::velocityW,
80 vn = Indices::velocityNw,
81 sw = Indices::saturationW,
82 sn = Indices::saturationNw,
83 pressureIdx = Indices::pressureIdx,
84 saturationIdx = Indices::saturationIdx,
85 pressEqIdx = Indices::pressureEqIdx,
86 satEqIdx = Indices::satEqIdx,
87 numPhases = getPropValue<TypeTag, Properties::NumPhases>()
88 };
89
90 using Element = typename GridView::template Codim<0>::Entity;
91
92 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
93 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
94
95public:
101 FvMpfaO2dPressureVelocity2p(Problem& problem) :
102 ParentType(problem), problem_(problem), velocity_(problem)
103 {
104 density_[wPhaseIdx] = 0.;
105 density_[nPhaseIdx] = 0.;
106 viscosity_[wPhaseIdx] = 0.;
107 viscosity_[nPhaseIdx] = 0.;
108
109 calcVelocityInTransport_ = getParam<bool>("MPFA.CalcVelocityInTransport");
110 }
111
113 void calculateVelocity();
114
116 void calculateVelocity(const Intersection&, CellData&);
117 void calculateVelocityOnBoundary(const Intersection& intersection, CellData& cellData);
118
121 {
122 this->updateMaterialLaws();
123
124 this->storePressureSolution();
125
128 }
129
136 {
137 const auto element = *problem_.gridView().template begin<0>();
138 FluidState fluidState;
139 fluidState.setPressure(wPhaseIdx, problem_.referencePressure(element));
140 fluidState.setPressure(nPhaseIdx, problem_.referencePressure(element));
141 fluidState.setTemperature(problem_.temperature(element));
142 fluidState.setSaturation(wPhaseIdx, 1.);
143 fluidState.setSaturation(nPhaseIdx, 0.);
144 density_[wPhaseIdx] = FluidSystem::density(fluidState, wPhaseIdx);
145 density_[nPhaseIdx] = FluidSystem::density(fluidState, nPhaseIdx);
146 viscosity_[wPhaseIdx] = FluidSystem::viscosity(fluidState, wPhaseIdx);
147 viscosity_[nPhaseIdx] = FluidSystem::viscosity(fluidState, nPhaseIdx);
148
150 velocity_.initialize();
153
154 return;
155 }
156
162 void update()
163 {
167 }
168
175 {
176 return calcVelocityInTransport_;
177 }
178
189 template<class MultiWriter>
190 void addOutputVtkFields(MultiWriter &writer)
191 {
193 velocity_.addOutputVtkFields(writer);
194 }
195
197 {
198 return velocity_;
199 }
200
201private:
202 Problem& problem_;
204 Scalar density_[numPhases];
205 Scalar viscosity_[numPhases];
206 bool calcVelocityInTransport_;
208 static const int pressureType_ = getPropValue<TypeTag, Properties::PressureFormulation>();
210 static const int saturationType_ = getPropValue<TypeTag, Properties::SaturationFormulation>();
211};
212// end of template
213
219template<class TypeTag>
221{
222 // run through all elements
223 for (const auto& vertex : vertices(problem_.gridView()))
224 {
225 int vIdxGlobal = problem_.variables().index(vertex);
226
227 InteractionVolume& interactionVolume = this->interactionVolumes_[vIdxGlobal];
228
229 if (interactionVolume.isInnerVolume())
230 {
231 auto element1 = interactionVolume.getSubVolumeElement(0);
232 auto element2 = interactionVolume.getSubVolumeElement(1);
233 auto element3 = interactionVolume.getSubVolumeElement(2);
234 auto element4 = interactionVolume.getSubVolumeElement(3);
235
236 // cell index
237 int eIdxGlobal1 = problem_.variables().index(element1);
238 int eIdxGlobal2 = problem_.variables().index(element2);
239 int eIdxGlobal3 = problem_.variables().index(element3);
240 int eIdxGlobal4 = problem_.variables().index(element4);
241
242 //get the cell Data
243 CellData& cellData1 = problem_.variables().cellData(eIdxGlobal1);
244 CellData& cellData2 = problem_.variables().cellData(eIdxGlobal2);
245 CellData& cellData3 = problem_.variables().cellData(eIdxGlobal3);
246 CellData& cellData4 = problem_.variables().cellData(eIdxGlobal4);
247
248 velocity_.calculateInnerInteractionVolumeVelocity(interactionVolume, cellData1, cellData2, cellData3,
249 cellData4, this->innerBoundaryVolumeFaces_);
250 }
251 // at least one face on boundary!
252 else
253 {
254 for (int elemIdx = 0; elemIdx < 2 * dim; elemIdx++)
255 {
256 bool isOutside = false;
257 for (int fIdx = 0; fIdx < dim; fIdx++)
258 {
259 int intVolFaceIdx = interactionVolume.getFaceIndexFromSubVolume(elemIdx, fIdx);
260 if (interactionVolume.isOutsideFace(intVolFaceIdx))
261 {
262 isOutside = true;
263 break;
264 }
265 }
266 if (isOutside)
267 {
268 continue;
269 }
270 // cell index
271 int eIdxGlobal = problem_.variables().index(interactionVolume.getSubVolumeElement(elemIdx));
272 //get the cell Data
273 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
274
275 velocity_.calculateBoundaryInteractionVolumeVelocity(interactionVolume, cellData, elemIdx);
276 }
277 } // end boundaries
278
279 } // end vertex iterator
280
281 return;
282} // end method calcTotalVelocity
283
292template<class TypeTag>
293void FvMpfaO2dPressureVelocity2p<TypeTag>::calculateVelocity(const Intersection& intersection, CellData& cellData)
294{
295 int numVertices = intersection.geometry().corners();
296
297 auto elementI = intersection.inside();
298 auto elementJ = intersection.outside();
299
300 int eIdxGlobalI = problem_.variables().index(elementI);
301 int eIdxGlobalJ = problem_.variables().index(elementJ);
302
303 CellData& cellDataJ = problem_.variables().cellData(eIdxGlobalJ);
304
305 const auto refElement = referenceElement(elementI);
306
307 int indexInInside = intersection.indexInInside();
308 int indexInOutside = intersection.indexInOutside();
309
310 Dune::FieldVector<CellData, 4> cellDataTemp;
311
312 for (int vIdx = 0; vIdx < numVertices; vIdx++)
313 {
314 int localVertIdx = refElement.subEntity(indexInInside, dim - 1, vIdx, dim);
315
316 int vIdxGlobal = problem_.variables().index(elementI.template subEntity<dim>(localVertIdx));
317
318 InteractionVolume& interactionVolume = this->interactionVolumes_[vIdxGlobal];
319
320 if (interactionVolume.isInnerVolume())
321 {
322
323 auto element1 = interactionVolume.getSubVolumeElement(0);
324 auto element2 = interactionVolume.getSubVolumeElement(1);
325 auto element3 = interactionVolume.getSubVolumeElement(2);
326 auto element4 = interactionVolume.getSubVolumeElement(3);
327
328 // cell index
329 int eIdxGlobal[4];
330 eIdxGlobal[0] = problem_.variables().index(element1);
331 eIdxGlobal[1] = problem_.variables().index(element2);
332 eIdxGlobal[2] = problem_.variables().index(element3);
333 eIdxGlobal[3] = problem_.variables().index(element4);
334
335 //get the cell Data
336 cellDataTemp[0] = problem_.variables().cellData(eIdxGlobal[0]);
337 cellDataTemp[1] = problem_.variables().cellData(eIdxGlobal[1]);
338 cellDataTemp[2] = problem_.variables().cellData(eIdxGlobal[2]);
339 cellDataTemp[3] = problem_.variables().cellData(eIdxGlobal[3]);
340
341 velocity_.calculateInnerInteractionVolumeVelocity(interactionVolume, cellDataTemp[0], cellDataTemp[1],
342 cellDataTemp[2], cellDataTemp[3], this->innerBoundaryVolumeFaces_);
343
344 for (int i = 0; i < 4; i++)
345 {
346 if (eIdxGlobal[i] == eIdxGlobalI)
347 {
348 cellData.fluxData().setVelocity(wPhaseIdx, indexInInside,
349 cellDataTemp[i].fluxData().velocity(wPhaseIdx, indexInInside));
350 cellData.fluxData().setVelocity(nPhaseIdx, indexInInside,
351 cellDataTemp[i].fluxData().velocity(nPhaseIdx, indexInInside));
352 cellData.fluxData().setUpwindPotential(wPhaseIdx, indexInInside,
353 cellDataTemp[i].fluxData().upwindPotential(wPhaseIdx, indexInInside));
354 cellData.fluxData().setUpwindPotential(nPhaseIdx, indexInInside,
355 cellDataTemp[i].fluxData().upwindPotential(nPhaseIdx, indexInInside));
356 }
357 else if (eIdxGlobal[i] == eIdxGlobalJ)
358 {
359 cellDataJ.fluxData().setVelocity(wPhaseIdx, indexInOutside,
360 cellDataTemp[i].fluxData().velocity(wPhaseIdx, indexInOutside));
361 cellDataJ.fluxData().setVelocity(nPhaseIdx, indexInOutside,
362 cellDataTemp[i].fluxData().velocity(nPhaseIdx, indexInOutside));
363 cellDataJ.fluxData().setUpwindPotential(wPhaseIdx, indexInOutside,
364 cellDataTemp[i].fluxData().upwindPotential(wPhaseIdx, indexInOutside));
365 cellDataJ.fluxData().setUpwindPotential(nPhaseIdx, indexInOutside,
366 cellDataTemp[i].fluxData().upwindPotential(nPhaseIdx, indexInOutside));
367 }
368 }
369 }
370 }
371
372
373 cellData.fluxData().setVelocityMarker(indexInInside);
374 cellDataJ.fluxData().setVelocityMarker(indexInOutside);
375}
376
385template<class TypeTag>
386void FvMpfaO2dPressureVelocity2p<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 const auto fluidMatrixInteraction = problem_.spatialParams().fluidMatrixInteractionAtPos(element.geometry().center());
459 const Scalar pcBound = fluidMatrixInteraction.pc(satW);
460
461 //determine phase pressures from primary pressure variable
462 Scalar pressWBound = 0;
463 Scalar pressNwBound = 0;
464 if (pressureType_ == pw)
465 {
466 pressWBound = pressBound;
467 pressNwBound = pressBound + pcBound;
468 }
469 else if (pressureType_ == pn)
470 {
471 pressWBound = pressBound - pcBound;
472 pressNwBound = pressBound;
473 }
474
475 const Scalar lambdaWBound = fluidMatrixInteraction.krw(satW) / viscosity_[wPhaseIdx];
476 const Scalar lambdaNwBound = fluidMatrixInteraction.krn(satW) / viscosity_[nPhaseIdx];
477
478 Scalar potentialDiffW = cellData.fluxData().upwindPotential(wPhaseIdx, isIndex);
479 Scalar potentialDiffNw = cellData.fluxData().upwindPotential(nPhaseIdx, isIndex);
480
481 //calculate potential gradient
482 potentialDiffW = (cellData.pressure(wPhaseIdx) - pressWBound);
483 potentialDiffNw = (cellData.pressure(nPhaseIdx) - pressNwBound);
484
485 potentialDiffW += density_[wPhaseIdx] * (distVec * problem_.gravity());
486 potentialDiffNw += density_[nPhaseIdx] * (distVec * problem_.gravity());
487
488 //store potential gradients for further calculations
489 cellData.fluxData().setUpwindPotential(wPhaseIdx, isIndex, potentialDiffW);
490 cellData.fluxData().setUpwindPotential(nPhaseIdx, isIndex, potentialDiffNw);
491
492 //do the upwinding of the mobility depending on the phase potentials
493 Scalar lambdaW = (potentialDiffW > 0.) ? lambdaWI : lambdaWBound;
494 lambdaW = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffW, 0.0, 1.0e-30)) ? 0.5 * (lambdaWI + lambdaWBound) : lambdaW;
495 Scalar lambdaNw = (potentialDiffNw > 0.) ? lambdaNwI : lambdaNwBound;
496 lambdaNw = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffNw, 0.0, 1.0e-30)) ? 0.5 * (lambdaNwI + lambdaNwBound) : lambdaNw;
497
498
499 Scalar scalarPerm = permeability.two_norm();
500
501 //calculate the gravity term
502 Dune::FieldVector<Scalar, dimWorld> velocityW(unitOuterNormal);
503 Dune::FieldVector<Scalar, dimWorld> velocityNw(unitOuterNormal);
504
505 //calculate unit distVec
506 distVec /= dist;
507 Scalar areaScaling = (unitOuterNormal * distVec);
508 //this treatment of g allows to account for gravity flux through faces where the face normal
509 //has no z component (e.g. parallelepiped grids)
510 Scalar gravityTermW = (problem_.gravity() * distVec) * density_[wPhaseIdx] * areaScaling;
511 Scalar gravityTermNw = (problem_.gravity() * distVec) * density_[nPhaseIdx] * areaScaling;
512
513 //calculate velocity depending on the pressure used -> use pc = pn - pw
514 switch (pressureType_)
515 {
516 case pw:
517 {
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;
521 break;
522 }
523 case pn:
524 {
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);
528 break;
529 }
530 }
531
532 //store velocities
533 cellData.fluxData().setVelocity(wPhaseIdx, isIndex, velocityW);
534 cellData.fluxData().setVelocity(nPhaseIdx, isIndex, velocityNw);
535 cellData.fluxData().setVelocityMarker(isIndex);
536
537 } //end dirichlet boundary
538
539 else if (bcType.isNeumann(pressEqIdx))
540 {
541 problem_.neumann(boundValues, intersection);
542
543 Dune::FieldVector<Scalar, dimWorld> velocityW(unitOuterNormal);
544 Dune::FieldVector<Scalar, dimWorld> velocityNw(unitOuterNormal);
545
546 velocityW *= boundValues[wPhaseIdx];
547 velocityNw *= boundValues[nPhaseIdx];
548
549 velocityW /= density_[wPhaseIdx];
550 velocityNw /= density_[nPhaseIdx];
551
552 //store potential gradients for further calculations
553 cellData.fluxData().setUpwindPotential(wPhaseIdx, isIndex, boundValues[wPhaseIdx]);
554 cellData.fluxData().setUpwindPotential(nPhaseIdx, isIndex, boundValues[nPhaseIdx]);
555
556 cellData.fluxData().setVelocity(wPhaseIdx, isIndex, velocityW);
557 cellData.fluxData().setVelocity(nPhaseIdx, isIndex, velocityNw);
558 cellData.fluxData().setVelocityMarker(isIndex);
559 } //end neumann boundary
560 else
561 {
562 DUNE_THROW(Dune::NotImplemented, "No valid boundary condition type defined for pressure equation!");
563 }
564}
565
566} // end namespace Dumux
567#endif
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
Finite volume MPFA O-method discretization of a two-phase flow pressure equation of the sequential IM...
Definition: omethod/2dpressure.hh:67
void updateMaterialLaws()
Constitutive functions are initialized and stored in the variables object.
Definition: omethod/2dpressure.hh:1928
void addOutputVtkFields(MultiWriter &writer)
Adds pressure output to the output file.
Definition: omethod/2dpressure.hh:324
void update()
Pressure update.
Definition: omethod/2dpressure.hh:210
void initialize()
Initializes the pressure model.
Definition: omethod/2dpressure.hh:174
void storePressureSolution()
Globally stores the pressure solution.
Definition: omethod/2dpressure.hh:251
Class for the calculation of velocities from the pressure solution of an IMPES scheme using a MPFA O-...
Definition: omethod/2dpressurevelocity.hh:45
bool calculateVelocityInTransport()
Indicates if velocity is reconstructed in the pressure step or in the transport step.
Definition: omethod/2dpressurevelocity.hh:174
void calculateVelocityOnBoundary(const Intersection &intersection, CellData &cellData)
Calculates the velocity at a boundary.
Definition: omethod/2dpressurevelocity.hh:386
void addOutputVtkFields(MultiWriter &writer)
Adds velocity output to the output file.
Definition: omethod/2dpressurevelocity.hh:190
void update()
Pressure and velocity update.
Definition: omethod/2dpressurevelocity.hh:162
void calculateVelocity()
Calculates the velocities at all cell-cell interfaces.
Definition: omethod/2dpressurevelocity.hh:220
FvMpfaO2dPressureVelocity2p(Problem &problem)
Constructs a FvMpfaO2dPressureVelocity2p object.
Definition: omethod/2dpressurevelocity.hh:101
void updateVelocity()
Function for updating the velocity field if iterations are necessary in the transport solution.
Definition: omethod/2dpressurevelocity.hh:120
void initialize()
Initializes pressure and velocity.
Definition: omethod/2dpressurevelocity.hh:135
FvMpfaO2dVelocity2P< TypeTag > velocity()
Definition: omethod/2dpressurevelocity.hh:196
Class for calculating velocities from cell-wise constant pressure values.
Definition: omethod/2dvelocity.hh:57
Class including the information of an interaction volume of a MPFA O-method that does not change with...
Definition: ointeractionvolume.hh:38
bool isInnerVolume()
Returns true if the interaction volume is completely inside the model domain.
Definition: ointeractionvolume.hh:291
bool isOutsideFace(int subVolumeFaceIdx)
Returns true if an interaction volume flux face is outside the model domain.
Definition: ointeractionvolume.hh:276
int getFaceIndexFromSubVolume(int subVolumeIdx, int subVolumeFaceIdx)
Map from local interaction volume numbering on element to numbering on interaction volume.
Definition: ointeractionvolume.hh:245
Element getSubVolumeElement(int subVolumeIdx)
Get an element of the interaction volume.
Definition: ointeractionvolume.hh:256
The finite volume base class for the solution of a pressure equation.
Definition: sequential/cellcentered/pressure.hh:49
Finite volume MPFA O-method discretization of a two-phase pressure equation of the sequential IMPES m...
Velocity calculation using a 2-d MPFA O-method.