3.2-git
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#include "2dpressure.hh"
27#include "2dvelocity.hh"
28
29namespace Dumux {
43template<class TypeTag> class FvMpfaO2dPressureVelocity2p: public FvMpfaO2dPressure2p<TypeTag>
44{
46
48
49 enum
50 {
51 dim = GridView::dimension, dimWorld = GridView::dimensionworld
52 };
53
59 using PrimaryVariables = typename SolutionTypes::PrimaryVariables;
60
62 using MaterialLaw = typename SpatialParams::MaterialLaw;
63
66
67 using Vertex = typename GridView::Traits::template Codim<dim>::Entity;
68 using Intersection = typename GridView::Intersection;
69
70 using ReferenceElements = Dune::ReferenceElements<Scalar, dim>;
71
73
75 enum
76 {
77 wPhaseIdx = Indices::wPhaseIdx,
78 nPhaseIdx = Indices::nPhaseIdx,
79 pw = Indices::pressureW,
80 pn = Indices::pressureNw,
81 vw = Indices::velocityW,
82 vn = Indices::velocityNw,
83 sw = Indices::saturationW,
84 sn = Indices::saturationNw,
85 pressureIdx = Indices::pressureIdx,
86 saturationIdx = Indices::saturationIdx,
87 pressEqIdx = Indices::pressureEqIdx,
88 satEqIdx = Indices::satEqIdx,
89 numPhases = getPropValue<TypeTag, Properties::NumPhases>()
90 };
91
92 using Element = typename GridView::template Codim<0>::Entity;
93
94 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
95 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
96
97public:
103 FvMpfaO2dPressureVelocity2p(Problem& problem) :
104 ParentType(problem), problem_(problem), velocity_(problem)
105 {
106 density_[wPhaseIdx] = 0.;
107 density_[nPhaseIdx] = 0.;
108 viscosity_[wPhaseIdx] = 0.;
109 viscosity_[nPhaseIdx] = 0.;
110
111 calcVelocityInTransport_ = getParam<bool>("MPFA.CalcVelocityInTransport");
112 }
113
115 void calculateVelocity();
116
118 void calculateVelocity(const Intersection&, CellData&);
119 void calculateVelocityOnBoundary(const Intersection& intersection, CellData& cellData);
120
123 {
124 this->updateMaterialLaws();
125
126 this->storePressureSolution();
127
130 }
131
138 {
139 const auto element = *problem_.gridView().template begin<0>();
140 FluidState fluidState;
141 fluidState.setPressure(wPhaseIdx, problem_.referencePressure(element));
142 fluidState.setPressure(nPhaseIdx, problem_.referencePressure(element));
143 fluidState.setTemperature(problem_.temperature(element));
144 fluidState.setSaturation(wPhaseIdx, 1.);
145 fluidState.setSaturation(nPhaseIdx, 0.);
146 density_[wPhaseIdx] = FluidSystem::density(fluidState, wPhaseIdx);
147 density_[nPhaseIdx] = FluidSystem::density(fluidState, nPhaseIdx);
148 viscosity_[wPhaseIdx] = FluidSystem::viscosity(fluidState, wPhaseIdx);
149 viscosity_[nPhaseIdx] = FluidSystem::viscosity(fluidState, nPhaseIdx);
150
152 velocity_.initialize();
155
156 return;
157 }
158
164 void update()
165 {
169 }
170
177 {
178 return calcVelocityInTransport_;
179 }
180
191 template<class MultiWriter>
192 void addOutputVtkFields(MultiWriter &writer)
193 {
195 velocity_.addOutputVtkFields(writer);
196 }
197
199 {
200 return velocity_;
201 }
202
203private:
204 Problem& problem_;
206 Scalar density_[numPhases];
207 Scalar viscosity_[numPhases];
208 bool calcVelocityInTransport_;
210 static const int pressureType_ = getPropValue<TypeTag, Properties::PressureFormulation>();
212 static const int saturationType_ = getPropValue<TypeTag, Properties::SaturationFormulation>();
213};
214// end of template
215
221template<class TypeTag>
223{
224 // run through all elements
225 for (const auto& vertex : vertices(problem_.gridView()))
226 {
227 int vIdxGlobal = problem_.variables().index(vertex);
228
229 InteractionVolume& interactionVolume = this->interactionVolumes_[vIdxGlobal];
230
231 if (interactionVolume.isInnerVolume())
232 {
233 auto element1 = interactionVolume.getSubVolumeElement(0);
234 auto element2 = interactionVolume.getSubVolumeElement(1);
235 auto element3 = interactionVolume.getSubVolumeElement(2);
236 auto element4 = interactionVolume.getSubVolumeElement(3);
237
238 // cell index
239 int eIdxGlobal1 = problem_.variables().index(element1);
240 int eIdxGlobal2 = problem_.variables().index(element2);
241 int eIdxGlobal3 = problem_.variables().index(element3);
242 int eIdxGlobal4 = problem_.variables().index(element4);
243
244 //get the cell Data
245 CellData& cellData1 = problem_.variables().cellData(eIdxGlobal1);
246 CellData& cellData2 = problem_.variables().cellData(eIdxGlobal2);
247 CellData& cellData3 = problem_.variables().cellData(eIdxGlobal3);
248 CellData& cellData4 = problem_.variables().cellData(eIdxGlobal4);
249
250 velocity_.calculateInnerInteractionVolumeVelocity(interactionVolume, cellData1, cellData2, cellData3,
251 cellData4, this->innerBoundaryVolumeFaces_);
252 }
253 // at least one face on boundary!
254 else
255 {
256 for (int elemIdx = 0; elemIdx < 2 * dim; elemIdx++)
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} // end method calcTotalVelocity
285
294template<class TypeTag>
295void FvMpfaO2dPressureVelocity2p<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 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 referenceElement = ReferenceElements::general(elementI.type());
308
309 int indexInInside = intersection.indexInInside();
310 int indexInOutside = intersection.indexInOutside();
311
312 Dune::FieldVector<CellData, 4> cellDataTemp;
313
314 for (int vIdx = 0; vIdx < numVertices; vIdx++)
315 {
316 int localVertIdx = referenceElement.subEntity(indexInInside, dim - 1, vIdx, dim);
317
318 int vIdxGlobal = problem_.variables().index(elementI.template subEntity<dim>(localVertIdx));
319
320 InteractionVolume& interactionVolume = this->interactionVolumes_[vIdxGlobal];
321
322 if (interactionVolume.isInnerVolume())
323 {
324
325 auto element1 = interactionVolume.getSubVolumeElement(0);
326 auto element2 = interactionVolume.getSubVolumeElement(1);
327 auto element3 = interactionVolume.getSubVolumeElement(2);
328 auto element4 = interactionVolume.getSubVolumeElement(3);
329
330 // cell index
331 int eIdxGlobal[4];
332 eIdxGlobal[0] = problem_.variables().index(element1);
333 eIdxGlobal[1] = problem_.variables().index(element2);
334 eIdxGlobal[2] = problem_.variables().index(element3);
335 eIdxGlobal[3] = problem_.variables().index(element4);
336
337 //get the cell Data
338 cellDataTemp[0] = problem_.variables().cellData(eIdxGlobal[0]);
339 cellDataTemp[1] = problem_.variables().cellData(eIdxGlobal[1]);
340 cellDataTemp[2] = problem_.variables().cellData(eIdxGlobal[2]);
341 cellDataTemp[3] = problem_.variables().cellData(eIdxGlobal[3]);
342
343 velocity_.calculateInnerInteractionVolumeVelocity(interactionVolume, cellDataTemp[0], cellDataTemp[1],
344 cellDataTemp[2], cellDataTemp[3], this->innerBoundaryVolumeFaces_);
345
346 for (int i = 0; i < 4; i++)
347 {
348 if (eIdxGlobal[i] == eIdxGlobalI)
349 {
350 cellData.fluxData().setVelocity(wPhaseIdx, indexInInside,
351 cellDataTemp[i].fluxData().velocity(wPhaseIdx, indexInInside));
352 cellData.fluxData().setVelocity(nPhaseIdx, indexInInside,
353 cellDataTemp[i].fluxData().velocity(nPhaseIdx, indexInInside));
354 cellData.fluxData().setUpwindPotential(wPhaseIdx, indexInInside,
355 cellDataTemp[i].fluxData().upwindPotential(wPhaseIdx, indexInInside));
356 cellData.fluxData().setUpwindPotential(nPhaseIdx, indexInInside,
357 cellDataTemp[i].fluxData().upwindPotential(nPhaseIdx, indexInInside));
358 }
359 else if (eIdxGlobal[i] == eIdxGlobalJ)
360 {
361 cellDataJ.fluxData().setVelocity(wPhaseIdx, indexInOutside,
362 cellDataTemp[i].fluxData().velocity(wPhaseIdx, indexInOutside));
363 cellDataJ.fluxData().setVelocity(nPhaseIdx, indexInOutside,
364 cellDataTemp[i].fluxData().velocity(nPhaseIdx, indexInOutside));
365 cellDataJ.fluxData().setUpwindPotential(wPhaseIdx, indexInOutside,
366 cellDataTemp[i].fluxData().upwindPotential(wPhaseIdx, indexInOutside));
367 cellDataJ.fluxData().setUpwindPotential(nPhaseIdx, indexInOutside,
368 cellDataTemp[i].fluxData().upwindPotential(nPhaseIdx, indexInOutside));
369 }
370 }
371 }
372 }
373
374
375 cellData.fluxData().setVelocityMarker(indexInInside);
376 cellDataJ.fluxData().setVelocityMarker(indexInOutside);
377}
378
387template<class TypeTag>
388void FvMpfaO2dPressureVelocity2p<TypeTag>::calculateVelocityOnBoundary(const Intersection& intersection, CellData& cellData)
389{
390 auto element = intersection.inside();
391
392 //get face index
393 int isIndex = intersection.indexInInside();
394
395 //get face normal
396 const Dune::FieldVector<Scalar, dim>& unitOuterNormal = intersection.centerUnitOuterNormal();
397
398 BoundaryTypes bcType;
399 //get boundary type
400 problem_.boundaryTypes(bcType, intersection);
401 PrimaryVariables boundValues(0.0);
402
403 if (bcType.isDirichlet(pressEqIdx))
404 {
405 problem_.dirichlet(boundValues, intersection);
406
407 // get global coordinates of cell centers
408 const GlobalPosition& globalPosI = element.geometry().center();
409
410 // center of face in global coordinates
411 const GlobalPosition& globalPosJ = intersection.geometry().center();
412
413 // get mobilities and fractional flow factors
414 Scalar lambdaWI = cellData.mobility(wPhaseIdx);
415 Scalar lambdaNwI = cellData.mobility(nPhaseIdx);
416
417 // get capillary pressure
418 Scalar pcI = cellData.capillaryPressure();
419
420 // distance vector between barycenters
421 GlobalPosition distVec = globalPosJ - globalPosI;
422
423 // compute distance between cell centers
424 Scalar dist = distVec.two_norm();
425
426 //permeability vector at boundary
427 // compute vectorized permeabilities
428 DimMatrix meanPermeability(0);
429
430 problem_.spatialParams().meanK(meanPermeability, problem_.spatialParams().intrinsicPermeability(element));
431
432 Dune::FieldVector<Scalar, dim> permeability(0);
433 meanPermeability.mv(unitOuterNormal, permeability);
434
435 //determine saturation at the boundary -> if no saturation is known directly at the boundary use the cell saturation
436 Scalar satW = 0;
437 if (bcType.isDirichlet(satEqIdx))
438 {
439 switch (saturationType_)
440 {
441 case sw:
442 {
443 satW = boundValues[saturationIdx];
444 break;
445 }
446 case sn:
447 {
448 satW = 1 - boundValues[saturationIdx];
449 break;
450 }
451 }
452 }
453 else
454 {
455 satW = cellData.saturation(wPhaseIdx);
456 }
457
458 Scalar pressBound = boundValues[pressureIdx];
459 Scalar pcBound = MaterialLaw::pc(problem_.spatialParams().materialLawParams(element), 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 Scalar lambdaWBound = MaterialLaw::krw(problem_.spatialParams().materialLawParams(element), satW)
476 / viscosity_[wPhaseIdx];
477 Scalar lambdaNwBound = MaterialLaw::krn(problem_.spatialParams().materialLawParams(element), satW)
478 / viscosity_[nPhaseIdx];
479
480 Scalar potentialDiffW = cellData.fluxData().upwindPotential(wPhaseIdx, isIndex);
481 Scalar potentialDiffNw = cellData.fluxData().upwindPotential(nPhaseIdx, isIndex);
482
483 //calculate potential gradient
484 potentialDiffW = (cellData.pressure(wPhaseIdx) - pressWBound);
485 potentialDiffNw = (cellData.pressure(nPhaseIdx) - pressNwBound);
486
487 potentialDiffW += density_[wPhaseIdx] * (distVec * problem_.gravity());
488 potentialDiffNw += density_[nPhaseIdx] * (distVec * problem_.gravity());
489
490 //store potential gradients for further calculations
491 cellData.fluxData().setUpwindPotential(wPhaseIdx, isIndex, potentialDiffW);
492 cellData.fluxData().setUpwindPotential(nPhaseIdx, isIndex, potentialDiffNw);
493
494 //do the upwinding of the mobility depending on the phase potentials
495 Scalar lambdaW = (potentialDiffW > 0.) ? lambdaWI : lambdaWBound;
496 lambdaW = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffW, 0.0, 1.0e-30)) ? 0.5 * (lambdaWI + lambdaWBound) : lambdaW;
497 Scalar lambdaNw = (potentialDiffNw > 0.) ? lambdaNwI : lambdaNwBound;
498 lambdaNw = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffNw, 0.0, 1.0e-30)) ? 0.5 * (lambdaNwI + lambdaNwBound) : lambdaNw;
499
500
501 Scalar scalarPerm = permeability.two_norm();
502
503 //calculate the gravity term
504 Dune::FieldVector<Scalar, dimWorld> velocityW(unitOuterNormal);
505 Dune::FieldVector<Scalar, dimWorld> velocityNw(unitOuterNormal);
506
507 //calculate unit distVec
508 distVec /= dist;
509 Scalar areaScaling = (unitOuterNormal * distVec);
510 //this treatment of g allows to account for gravity flux through faces where the face normal
511 //has no z component (e.g. parallelepiped grids)
512 Scalar gravityTermW = (problem_.gravity() * distVec) * density_[wPhaseIdx] * areaScaling;
513 Scalar gravityTermNw = (problem_.gravity() * distVec) * density_[nPhaseIdx] * areaScaling;
514
515 //calculate velocity depending on the pressure used -> use pc = pn - pw
516 switch (pressureType_)
517 {
518 case pw:
519 {
520 velocityW *= lambdaW * scalarPerm * ((cellData.pressure(wPhaseIdx) - pressBound) / dist + gravityTermW);
521 velocityNw *= lambdaNw * scalarPerm * ((cellData.pressure(wPhaseIdx) - pressBound) / dist + gravityTermNw)
522 + 0.5 * (lambdaNwI + lambdaNwBound) * scalarPerm * (pcI - pcBound) / dist;
523 break;
524 }
525 case pn:
526 {
527 velocityW *= lambdaW * scalarPerm * ((cellData.pressure(nPhaseIdx) - pressBound) / dist + gravityTermW)
528 - 0.5 * (lambdaWI + lambdaWBound) * scalarPerm * (pcI - pcBound) / dist;
529 velocityNw *= lambdaNw * scalarPerm * ((cellData.pressure(nPhaseIdx) - pressBound) / dist + gravityTermNw);
530 break;
531 }
532 }
533
534 //store velocities
535 cellData.fluxData().setVelocity(wPhaseIdx, isIndex, velocityW);
536 cellData.fluxData().setVelocity(nPhaseIdx, isIndex, velocityNw);
537 cellData.fluxData().setVelocityMarker(isIndex);
538
539 } //end dirichlet boundary
540
541 else if (bcType.isNeumann(pressEqIdx))
542 {
543 problem_.neumann(boundValues, intersection);
544
545 Dune::FieldVector<Scalar, dimWorld> velocityW(unitOuterNormal);
546 Dune::FieldVector<Scalar, dimWorld> velocityNw(unitOuterNormal);
547
548 velocityW *= boundValues[wPhaseIdx];
549 velocityNw *= boundValues[nPhaseIdx];
550
551 velocityW /= density_[wPhaseIdx];
552 velocityNw /= density_[nPhaseIdx];
553
554 //store potential gradients for further calculations
555 cellData.fluxData().setUpwindPotential(wPhaseIdx, isIndex, boundValues[wPhaseIdx]);
556 cellData.fluxData().setUpwindPotential(nPhaseIdx, isIndex, boundValues[nPhaseIdx]);
557
558 cellData.fluxData().setVelocity(wPhaseIdx, isIndex, velocityW);
559 cellData.fluxData().setVelocity(nPhaseIdx, isIndex, velocityNw);
560 cellData.fluxData().setVelocityMarker(isIndex);
561 } //end neumann boundary
562 else
563 {
564 DUNE_THROW(Dune::NotImplemented, "No valid boundary condition type defined for pressure equation!");
565 }
566}
567
568} // end namespace Dumux
569#endif
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
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:1934
void addOutputVtkFields(MultiWriter &writer)
Adds pressure output to the output file.
Definition: omethod/2dpressure.hh:327
void update()
Pressure update.
Definition: omethod/2dpressure.hh:213
void initialize()
Initializes the pressure model.
Definition: omethod/2dpressure.hh:177
void storePressureSolution()
Globally stores the pressure solution.
Definition: omethod/2dpressure.hh:254
Class for the calculation of velocities from the pressure solution of an IMPES scheme using a MPFA O-...
Definition: omethod/2dpressurevelocity.hh:44
bool calculateVelocityInTransport()
Indicates if velocity is reconstructed in the pressure step or in the transport step.
Definition: omethod/2dpressurevelocity.hh:176
void calculateVelocityOnBoundary(const Intersection &intersection, CellData &cellData)
Calculates the velocity at a boundary.
Definition: omethod/2dpressurevelocity.hh:388
void addOutputVtkFields(MultiWriter &writer)
Adds velocity output to the output file.
Definition: omethod/2dpressurevelocity.hh:192
void update()
Pressure and velocity update.
Definition: omethod/2dpressurevelocity.hh:164
void calculateVelocity()
Calculates the velocities at all cell-cell interfaces.
Definition: omethod/2dpressurevelocity.hh:222
FvMpfaO2dPressureVelocity2p(Problem &problem)
Constructs a FvMpfaO2dPressureVelocity2p object.
Definition: omethod/2dpressurevelocity.hh:103
void updateVelocity()
Function for updating the velocity field if iterations are necessary in the transport solution.
Definition: omethod/2dpressurevelocity.hh:122
void initialize()
Initializes pressure and velocity.
Definition: omethod/2dpressurevelocity.hh:137
FvMpfaO2dVelocity2P< TypeTag > velocity()
Definition: omethod/2dpressurevelocity.hh:198
Class for calculating velocities from cell-wise constant pressure values.
Definition: omethod/2dvelocity.hh:56
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.