3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
lmethod/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_MPFAL2DPRESSUREVELOCITY2P_HH
23#define DUMUX_MPFAL2DPRESSUREVELOCITY2P_HH
24
25#include <dune/common/float_cmp.hh>
26
27#include "2dpressure.hh"
28#include "2dvelocity.hh"
29
31
32namespace Dumux {
33
47template<class TypeTag> class FvMpfaL2dPressureVelocity2p: public FvMpfaL2dPressure2p<TypeTag>
48{
50
52
53 enum
54 {
55 dim = GridView::dimension, dimWorld = GridView::dimensionworld
56 };
57
63 using PrimaryVariables = typename SolutionTypes::PrimaryVariables;
64
66
68
71
72 using IndexSet = typename GridView::IndexSet;
73 using Intersection = typename GridView::Intersection;
74
76
77 enum
78 {
79 pw = Indices::pressureW,
80 pn = Indices::pressureNw,
81 sw = Indices::saturationW,
82 sn = Indices::saturationNw,
83 wPhaseIdx = Indices::wPhaseIdx,
84 nPhaseIdx = Indices::nPhaseIdx,
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 FvMpfaL2dPressureVelocity2p(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
198private:
199 Problem& problem_;
201
202 Scalar density_[numPhases];
203 Scalar viscosity_[numPhases];
204 bool calcVelocityInTransport_;
205
207 static const int pressureType_ = getPropValue<TypeTag, Properties::PressureFormulation>();
209 static const int saturationType_ = getPropValue<TypeTag, Properties::SaturationFormulation>();
210};
211// end of template
212
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,
249 cellData3, cellData4, this->innerBoundaryVolumeFaces_);
250
251 }
252 // at least one face on boundary!
253 else
254 {
255 for (int elemIdx = 0; elemIdx < 2 * dim; elemIdx++)
256 {
257 bool isOutside = false;
258 for (int fIdx = 0; fIdx < dim; fIdx++)
259 {
260 int intVolFaceIdx = interactionVolume.getFaceIndexFromSubVolume(elemIdx, fIdx);
261 if (interactionVolume.isOutsideFace(intVolFaceIdx))
262 {
263 isOutside = true;
264 break;
265 }
266 }
267 if (isOutside)
268 {
269 continue;
270 }
271 // cell index
272 int eIdxGlobal = problem_.variables().index(interactionVolume.getSubVolumeElement(elemIdx));
273 //get the cell Data
274 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
275
276 velocity_.calculateBoundaryInteractionVolumeVelocity(interactionVolume, cellData, elemIdx);
277 }
278
279 } // end boundaries
280
281 } // end vertex iterator
282
283// }
284 return;
285} // end method calcTotalVelocity
286
295template<class TypeTag>
296void FvMpfaL2dPressureVelocity2p<TypeTag>::calculateVelocity(const Intersection& intersection, CellData& cellData)
297{
298 int numVertices = intersection.geometry().corners();
299
300 auto elementI = intersection.inside();
301 auto elementJ = intersection.outside();
302
303 int eIdxGlobalI = problem_.variables().index(elementI);
304 int eIdxGlobalJ = problem_.variables().index(elementJ);
305
306 CellData& cellDataJ = problem_.variables().cellData(eIdxGlobalJ);
307
308 const auto refElement = referenceElement(elementI);
309
310 int indexInInside = intersection.indexInInside();
311 int indexInOutside = intersection.indexInOutside();
312
313 Dune::FieldVector<CellData, 4> cellDataTemp;
314
315 for (int vIdx = 0; vIdx < numVertices; vIdx++)
316 {
317 int localVertIdx = refElement.subEntity(indexInInside, dim - 1, vIdx, dim);
318
319 int vIdxGlobal = problem_.variables().index(elementI.template subEntity<dim>(localVertIdx));
320
321 InteractionVolume& interactionVolume = this->interactionVolumes_[vIdxGlobal];
322
323 if (interactionVolume.isInnerVolume())
324 {
325
326 auto element1 = interactionVolume.getSubVolumeElement(0);
327 auto element2 = interactionVolume.getSubVolumeElement(1);
328 auto element3 = interactionVolume.getSubVolumeElement(2);
329 auto element4 = interactionVolume.getSubVolumeElement(3);
330
331 // cell index
332 int eIdxGlobal[4];
333 eIdxGlobal[0] = problem_.variables().index(element1);
334 eIdxGlobal[1] = problem_.variables().index(element2);
335 eIdxGlobal[2] = problem_.variables().index(element3);
336 eIdxGlobal[3] = problem_.variables().index(element4);
337
338 //get the cell Data
339 cellDataTemp[0] = problem_.variables().cellData(eIdxGlobal[0]);
340 cellDataTemp[1] = problem_.variables().cellData(eIdxGlobal[1]);
341 cellDataTemp[2] = problem_.variables().cellData(eIdxGlobal[2]);
342 cellDataTemp[3] = problem_.variables().cellData(eIdxGlobal[3]);
343
344 velocity_.calculateInnerInteractionVolumeVelocity(interactionVolume, cellDataTemp[0], cellDataTemp[1],
345 cellDataTemp[2], cellDataTemp[3], this->innerBoundaryVolumeFaces_);
346
347 for (int i = 0; i < 4; i++)
348 {
349 if (eIdxGlobal[i] == eIdxGlobalI)
350 {
351 cellData.fluxData().setVelocity(wPhaseIdx, indexInInside,
352 cellDataTemp[i].fluxData().velocity(wPhaseIdx, indexInInside));
353 cellData.fluxData().setVelocity(nPhaseIdx, indexInInside,
354 cellDataTemp[i].fluxData().velocity(nPhaseIdx, indexInInside));
355 cellData.fluxData().setUpwindPotential(wPhaseIdx, indexInInside,
356 cellDataTemp[i].fluxData().upwindPotential(wPhaseIdx, indexInInside));
357 cellData.fluxData().setUpwindPotential(nPhaseIdx, indexInInside,
358 cellDataTemp[i].fluxData().upwindPotential(nPhaseIdx, indexInInside));
359 }
360 else if (eIdxGlobal[i] == eIdxGlobalJ)
361 {
362 cellDataJ.fluxData().setVelocity(wPhaseIdx, indexInOutside,
363 cellDataTemp[i].fluxData().velocity(wPhaseIdx, indexInOutside));
364 cellDataJ.fluxData().setVelocity(nPhaseIdx, indexInOutside,
365 cellDataTemp[i].fluxData().velocity(nPhaseIdx, indexInOutside));
366 cellDataJ.fluxData().setUpwindPotential(wPhaseIdx, indexInOutside,
367 cellDataTemp[i].fluxData().upwindPotential(wPhaseIdx, indexInOutside));
368 cellDataJ.fluxData().setUpwindPotential(nPhaseIdx, indexInOutside,
369 cellDataTemp[i].fluxData().upwindPotential(nPhaseIdx, indexInOutside));
370 }
371 }
372 }
373 }
374 cellData.fluxData().setVelocityMarker(indexInInside);
375 cellDataJ.fluxData().setVelocityMarker(indexInOutside);
376}
377
386template<class TypeTag>
387void FvMpfaL2dPressureVelocity2p<TypeTag>::calculateVelocityOnBoundary(const Intersection& intersection, CellData& cellData)
388{
389 auto element = intersection.inside();
390
391 //get face index
392 int isIndex = intersection.indexInInside();
393
394 //get face normal
395 const Dune::FieldVector<Scalar, dim>& unitOuterNormal = intersection.centerUnitOuterNormal();
396
397 BoundaryTypes bcType;
398 //get boundary type
399 problem_.boundaryTypes(bcType, intersection);
400 PrimaryVariables boundValues(0.0);
401
402 if (bcType.isDirichlet(pressEqIdx))
403 {
404 problem_.dirichlet(boundValues, intersection);
405
406 // get global coordinates of cell centers
407 const GlobalPosition& globalPosI = element.geometry().center();
408
409 // center of face in global coordinates
410 const GlobalPosition& globalPosJ = intersection.geometry().center();
411
412 // get mobilities and fractional flow factors
413 Scalar lambdaWI = cellData.mobility(wPhaseIdx);
414 Scalar lambdaNwI = cellData.mobility(nPhaseIdx);
415
416 // get capillary pressure
417 Scalar pcI = cellData.capillaryPressure();
418
419 // distance vector between barycenters
420 GlobalPosition distVec = globalPosJ - globalPosI;
421
422 // compute distance between cell centers
423 Scalar dist = distVec.two_norm();
424
425 //permeability vector at boundary
426 // compute vectorized permeabilities
427 DimMatrix meanPermeability(0);
428
429 problem_.spatialParams().meanK(meanPermeability, problem_.spatialParams().intrinsicPermeability(element));
430
431 Dune::FieldVector<Scalar, dim> permeability(0);
432 meanPermeability.mv(unitOuterNormal, permeability);
433
434 //determine saturation at the boundary -> if no saturation is known directly at the boundary use the cell saturation
435 Scalar satW = 0;
436 if (bcType.isDirichlet(satEqIdx))
437 {
438 switch (saturationType_)
439 {
440 case sw:
441 {
442 satW = boundValues[saturationIdx];
443 break;
444 }
445 case sn:
446 {
447 satW = 1 - boundValues[saturationIdx];
448 break;
449 }
450 }
451 }
452 else
453 {
454 satW = cellData.saturation(wPhaseIdx);
455 }
456
457 const Scalar pressBound = boundValues[pressureIdx];
458
459 // old material law interface is deprecated: Replace this by
460 // const auto& fluidMatrixInteraction = spatialParams.fluidMatrixInteractionAtPos(element.geometry().center());
461 // after the release of 3.3, when the deprecated interface is no longer supported
462 const auto fluidMatrixInteraction = Deprecated::makePcKrSw(Scalar{}, problem_.spatialParams(), element);
463
464 const Scalar pcBound = fluidMatrixInteraction.pc(satW);
465
466 //determine phase pressures from primary pressure variable
467 Scalar pressWBound = 0;
468 Scalar pressNwBound = 0;
469 if (pressureType_ == pw)
470 {
471 pressWBound = pressBound;
472 pressNwBound = pressBound + pcBound;
473 }
474 else if (pressureType_ == pn)
475 {
476 pressWBound = pressBound - pcBound;
477 pressNwBound = pressBound;
478 }
479
480 const Scalar lambdaWBound = fluidMatrixInteraction.krw(satW) / viscosity_[wPhaseIdx];
481 const Scalar lambdaNwBound = fluidMatrixInteraction.krn(satW) / viscosity_[nPhaseIdx];
482
483 Scalar potentialDiffW = cellData.fluxData().upwindPotential(wPhaseIdx, isIndex);
484 Scalar potentialDiffNw = cellData.fluxData().upwindPotential(nPhaseIdx, isIndex);
485
486 //calculate potential gradient
487 potentialDiffW = (cellData.pressure(wPhaseIdx) - pressWBound);
488 potentialDiffNw = (cellData.pressure(nPhaseIdx) - pressNwBound);
489
490 potentialDiffW += density_[wPhaseIdx] * (distVec * problem_.gravity());
491 potentialDiffNw += density_[nPhaseIdx] * (distVec * problem_.gravity());
492
493 //store potential gradients for further calculations
494 cellData.fluxData().setUpwindPotential(wPhaseIdx, isIndex, potentialDiffW);
495 cellData.fluxData().setUpwindPotential(nPhaseIdx, isIndex, potentialDiffNw);
496
497 //do the upwinding of the mobility depending on the phase potentials
498 Scalar lambdaW = (potentialDiffW > 0.) ? lambdaWI : lambdaWBound;
499 lambdaW = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffW, 0.0, 1.0e-30)) ? 0.5 * (lambdaWI + lambdaWBound) : lambdaW;
500 Scalar lambdaNw = (potentialDiffNw > 0.) ? lambdaNwI : lambdaNwBound;
501 lambdaNw = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffNw, 0.0, 1.0e-30)) ? 0.5 * (lambdaNwI + lambdaNwBound) : lambdaNw;
502
503
504 Scalar scalarPerm = permeability.two_norm();
505
506 //calculate the gravity term
507 Dune::FieldVector<Scalar, dimWorld> velocityW(unitOuterNormal);
508 Dune::FieldVector<Scalar, dimWorld> velocityNw(unitOuterNormal);
509
510 //calculate unit distVec
511 distVec /= dist;
512 Scalar areaScaling = (unitOuterNormal * distVec);
513 //this treatment of g allows to account for gravity flux through faces where the face normal has no z component (e.g. parallelepiped grids)
514 Scalar gravityTermW = (problem_.gravity() * distVec) * density_[wPhaseIdx] * areaScaling;
515 Scalar gravityTermNw = (problem_.gravity() * distVec) * density_[nPhaseIdx] * areaScaling;
516
517 //calculate velocity depending on the pressure used -> use pc = pn - pw
518 switch (pressureType_)
519 {
520 case pw:
521 {
522 velocityW *= lambdaW * scalarPerm * ((cellData.pressure(wPhaseIdx) - pressBound) / dist + gravityTermW);
523 velocityNw *= lambdaNw * scalarPerm * ((cellData.pressure(wPhaseIdx) - pressBound) / dist + gravityTermNw)
524 + 0.5 * (lambdaNwI + lambdaNwBound) * scalarPerm * (pcI - pcBound) / dist;
525 break;
526 }
527 case pn:
528 {
529 velocityW *= lambdaW * scalarPerm * ((cellData.pressure(nPhaseIdx) - pressBound) / dist + gravityTermW)
530 - 0.5 * (lambdaWI + lambdaWBound) * scalarPerm * (pcI - pcBound) / dist;
531 velocityNw *= lambdaNw * scalarPerm * ((cellData.pressure(nPhaseIdx) - pressBound) / dist + gravityTermNw);
532 break;
533 }
534 }
535
536 //store velocities
537 cellData.fluxData().setVelocity(wPhaseIdx, isIndex, velocityW);
538 cellData.fluxData().setVelocity(nPhaseIdx, isIndex, velocityNw);
539 cellData.fluxData().setVelocityMarker(isIndex);
540
541 } //end dirichlet boundary
542
543 else if (bcType.isNeumann(pressEqIdx))
544 {
545 problem_.neumann(boundValues, intersection);
546
547 Dune::FieldVector<Scalar, dimWorld> velocityW(unitOuterNormal);
548 Dune::FieldVector<Scalar, dimWorld> velocityNw(unitOuterNormal);
549
550 velocityW *= boundValues[wPhaseIdx];
551 velocityNw *= boundValues[nPhaseIdx];
552
553 velocityW /= density_[wPhaseIdx];
554 velocityNw /= density_[nPhaseIdx];
555
556 //store potential gradients for further calculations
557 cellData.fluxData().setUpwindPotential(wPhaseIdx, isIndex, boundValues[wPhaseIdx]);
558 cellData.fluxData().setUpwindPotential(nPhaseIdx, isIndex, boundValues[nPhaseIdx]);
559
560 cellData.fluxData().setVelocity(wPhaseIdx, isIndex, velocityW);
561 cellData.fluxData().setVelocity(nPhaseIdx, isIndex, velocityNw);
562 cellData.fluxData().setVelocityMarker(isIndex);
563 } //end neumann boundary
564 else
565 {
566 DUNE_THROW(Dune::NotImplemented, "No valid boundary condition type defined for pressure equation!");
567 }
568}
569
570} // end namespace Dumux
571#endif
Helpers for deprecation.
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 L-method discretization of a two-phase flow pressure equation of the sequential IM...
Definition: lmethod/2dpressure.hh:74
void update()
Pressure update.
Definition: lmethod/2dpressure.hh:290
void initialize()
Initializes the pressure model.
Definition: lmethod/2dpressure.hh:193
FVMPFALInteractionVolume< TypeTag > InteractionVolume
Type of the interaction volume objects.
Definition: lmethod/2dpressure.hh:150
void storePressureSolution()
Globally stores the pressure solution.
Definition: lmethod/2dpressure.hh:227
void addOutputVtkFields(MultiWriter &writer)
Adds pressure output to the output file.
Definition: lmethod/2dpressure.hh:339
void updateMaterialLaws()
constitutive functions are initialized and stored in the variables object
Definition: lmethod/2dpressure.hh:1834
Class for the calculation of velocities from the pressure solution of an IMPES scheme using a MPFA L-...
Definition: lmethod/2dpressurevelocity.hh:48
void updateVelocity()
Function for updating the velocity field if iterations are necessary in the transport solution.
Definition: lmethod/2dpressurevelocity.hh:122
void calculateVelocityOnBoundary(const Intersection &intersection, CellData &cellData)
Calculates the velocity at a boundary.
Definition: lmethod/2dpressurevelocity.hh:387
FvMpfaL2dPressureVelocity2p(Problem &problem)
Constructs a FvMpfaL2dPressureVelocity2p object.
Definition: lmethod/2dpressurevelocity.hh:103
void initialize()
Initializes pressure and velocity.
Definition: lmethod/2dpressurevelocity.hh:137
void addOutputVtkFields(MultiWriter &writer)
Adds velocity output to the output file.
Definition: lmethod/2dpressurevelocity.hh:192
bool calculateVelocityInTransport()
Indicates if velocity is reconstructed in the pressure step or in the transport step.
Definition: lmethod/2dpressurevelocity.hh:176
void calculateVelocity()
Calculates the velocities at all cell-cell interfaces.
Definition: lmethod/2dpressurevelocity.hh:220
void update()
Pressure and velocity update.
Definition: lmethod/2dpressurevelocity.hh:164
Class for calculating 2-d velocities from cell-wise constant pressure values.
Definition: lmethod/2dvelocity.hh:60
Class including the information of an interaction volume of a MPFA L-method that does not change with...
Definition: linteractionvolume.hh:41
bool isInnerVolume()
Returns true if the interaction volume is completely inside the model domain.
Definition: linteractionvolume.hh:316
Element getSubVolumeElement(int subVolumeIdx)
Get an element of the interaction volume.
Definition: linteractionvolume.hh:281
int getFaceIndexFromSubVolume(int subVolumeIdx, int subVolumeFaceIdx)
Map from local interaction volume numbering on element to numbering on interaction volume.
Definition: linteractionvolume.hh:270
bool isOutsideFace(int subVolumeFaceIdx)
Returns true if an interaction volume flux face is outside the model domain.
Definition: linteractionvolume.hh:301
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.