3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
2dpressurevelocityadaptive.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_MPFAL2DPRESSUREVELOCITIES2P_ADAPTIVE_HH
23#define DUMUX_MPFAL2DPRESSUREVELOCITIES2P_ADAPTIVE_HH
24
25#include <dune/common/float_cmp.hh>
26
27#include "2dpressureadaptive.hh"
28#include "2dvelocityadaptive.hh"
29
31
32namespace Dumux {
33
49template<class TypeTag> class FvMpfaL2dPressureVelocity2pAdaptive: public FvMpfaL2dPressure2pAdaptive<TypeTag>
50{
52
54
55 enum
56 {
57 dim = GridView::dimension, dimWorld = GridView::dimensionworld
58 };
59
65 using PrimaryVariables = typename SolutionTypes::PrimaryVariables;
66
68
70
73
74 using Grid = typename GridView::Grid;
75 using IndexSet = typename GridView::IndexSet;
76 using Intersection = typename GridView::Intersection;
77
78 using Geometry = typename Grid::template Codim<0>::Entity::Geometry;
79 using JacobianTransposed = typename Geometry::JacobianTransposed ;
80
82
83 enum
84 {
85 pw = Indices::pressureW,
86 pn = Indices::pressureNw,
87 sw = Indices::saturationW,
88 sn = Indices::saturationNw
89 };
90 enum
91 {
92 wPhaseIdx = Indices::wPhaseIdx,
93 nPhaseIdx = Indices::nPhaseIdx,
94 pressureIdx = Indices::pressureIdx,
95 saturationIdx = Indices::saturationIdx,
96 pressEqIdx = Indices::pressureEqIdx,
97 satEqIdx = Indices::satEqIdx,
98 numPhases = getPropValue<TypeTag, Properties::NumPhases>()
99 };
100
101 using Element = typename GridView::template Codim<0>::Entity;
102
103 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
104 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
105 using DimVector = Dune::FieldVector<Scalar, dim>;
106
107public:
114 ParentType(problem), problem_(problem), velocity_(problem)
115 {
116 density_[wPhaseIdx] = 0.;
117 density_[nPhaseIdx] = 0.;
118 viscosity_[wPhaseIdx] = 0.;
119 viscosity_[nPhaseIdx] = 0.;
120
121 calcVelocityInTransport_ = getParam<bool>("MPFA.CalcVelocityInTransport");
122 }
123
125 void calculateVelocity();
126
128 void calculateVelocity(const Intersection&, CellData&);
129 void calculateVelocityOnBoundary(const Intersection& intersection, CellData& cellData);
130
133 {
134 this->updateMaterialLaws();
135
136 this->storePressureSolution();
137
140 }
141
148 {
149 const auto element = *problem_.gridView().template begin<0>();
150 FluidState fluidState;
151 fluidState.setPressure(wPhaseIdx, problem_.referencePressure(element));
152 fluidState.setPressure(nPhaseIdx, problem_.referencePressure(element));
153 fluidState.setTemperature(problem_.temperature(element));
154 fluidState.setSaturation(wPhaseIdx, 1.);
155 fluidState.setSaturation(nPhaseIdx, 0.);
156 density_[wPhaseIdx] = FluidSystem::density(fluidState, wPhaseIdx);
157 density_[nPhaseIdx] = FluidSystem::density(fluidState, nPhaseIdx);
158 viscosity_[wPhaseIdx] = FluidSystem::viscosity(fluidState, wPhaseIdx);
159 viscosity_[nPhaseIdx] = FluidSystem::viscosity(fluidState, nPhaseIdx);
160
162 velocity_.initialize();
165
166 return;
167 }
168
174 void update()
175 {
179 }
180
187 {
188 return calcVelocityInTransport_;
189 }
190
201 template<class MultiWriter>
202 void addOutputVtkFields(MultiWriter &writer)
203 {
205 velocity_.addOutputVtkFields(writer);
206 }
207
208private:
209 Problem& problem_;
211
212 Scalar density_[numPhases];
213 Scalar viscosity_[numPhases];
214 bool calcVelocityInTransport_;
215
217 static const int pressureType_ = getPropValue<TypeTag, Properties::PressureFormulation>();
219 static const int saturationType_ = getPropValue<TypeTag, Properties::SaturationFormulation>();
220};
221// end of template
222
228template<class TypeTag>
230{
231 // run through all elements
232 for (const auto& vertex : vertices(problem_.gridView()))
233 {
234 int vIdxGlobal = problem_.variables().index(vertex);
235
236 InteractionVolume& interactionVolume = this->interactionVolumes_[vIdxGlobal];
237
238 if (interactionVolume.isInnerVolume())
239 {
240 if (interactionVolume.getElementNumber() == 4)
241 {
242 auto element1 = interactionVolume.getSubVolumeElement(0);
243 auto element2 = interactionVolume.getSubVolumeElement(1);
244 auto element3 = interactionVolume.getSubVolumeElement(2);
245 auto element4 = interactionVolume.getSubVolumeElement(3);
246
247 // cell index
248 int eIdxGlobal1 = problem_.variables().index(element1);
249 int eIdxGlobal2 = problem_.variables().index(element2);
250 int eIdxGlobal3 = problem_.variables().index(element3);
251 int eIdxGlobal4 = problem_.variables().index(element4);
252
253 //get the cell Data
254 CellData& cellData1 = problem_.variables().cellData(eIdxGlobal1);
255 CellData& cellData2 = problem_.variables().cellData(eIdxGlobal2);
256 CellData& cellData3 = problem_.variables().cellData(eIdxGlobal3);
257 CellData& cellData4 = problem_.variables().cellData(eIdxGlobal4);
258
259 velocity_.calculateInnerInteractionVolumeVelocity(interactionVolume, cellData1, cellData2, cellData3,
260 cellData4, this->innerBoundaryVolumeFaces_);
261
262 }
263 else if (interactionVolume.getElementNumber() == 3)
264 {
265 auto element1 = interactionVolume.getSubVolumeElement(0);
266 auto element2 = interactionVolume.getSubVolumeElement(1);
267 auto element4 = interactionVolume.getSubVolumeElement(3);
268
269 // cell index
270 int eIdxGlobal1 = problem_.variables().index(element1);
271 int eIdxGlobal2 = problem_.variables().index(element2);
272 int eIdxGlobal4 = problem_.variables().index(element4);
273
274 //get the cell Data
275 CellData& cellData1 = problem_.variables().cellData(eIdxGlobal1);
276 CellData& cellData2 = problem_.variables().cellData(eIdxGlobal2);
277 CellData& cellData4 = problem_.variables().cellData(eIdxGlobal4);
278
279 velocity_.calculateHangingNodeInteractionVolumeVelocity(interactionVolume, cellData1, cellData2,
280 cellData4, this->innerBoundaryVolumeFaces_);
281
282 }
283 else
284 {
285 DUNE_THROW(Dune::NotImplemented, "Unknown interactionvolume type!");
286 }
287 }
288
289 // at least one face on boundary!
290 else
291 {
292 for (int elemIdx = 0; elemIdx < 2 * dim; elemIdx++)
293 {
294 bool isOutside = false;
295 for (int fIdx = 0; fIdx < dim; fIdx++)
296 {
297 int intVolFaceIdx = interactionVolume.getFaceIndexFromSubVolume(elemIdx, fIdx);
298 if (interactionVolume.isOutsideFace(intVolFaceIdx))
299 {
300 isOutside = true;
301 break;
302 }
303 }
304 if (isOutside)
305 {
306 continue;
307 }
308 // cell index
309 int eIdxGlobal = problem_.variables().index(interactionVolume.getSubVolumeElement(elemIdx));
310 //get the cell Data
311 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
312
313 velocity_.calculateBoundaryInteractionVolumeVelocity(interactionVolume, cellData, elemIdx);
314 }
315
316 } // end boundaries
317
318 } // end vertex iterator
319 return;
320} // end method calcTotalVelocity
321
330template<class TypeTag>
331void FvMpfaL2dPressureVelocity2pAdaptive<TypeTag>::calculateVelocity(const Intersection& intersection, CellData& cellData)
332{
333 int numVertices = intersection.geometry().corners();
334
335 auto elementI = intersection.inside();
336 auto elementJ = intersection.outside();
337
338 int levelI = elementI.level();
339 int levelJ = elementJ.level();
340
341 int eIdxGlobalI = problem_.variables().index(elementI);
342 int eIdxGlobalJ = problem_.variables().index(elementJ);
343
344 CellData& cellDataJ = problem_.variables().cellData(eIdxGlobalJ);
345
346 const auto refElement = referenceElement(elementI);
347
348 int indexInInside = intersection.indexInInside();
349 int indexInOutside = intersection.indexInOutside();
350
351 int fIdx = indexInInside;
352
353 if (levelI < levelJ)
354 fIdx = indexInOutside;
355
356 std::vector<CellData> cellDataTemp(0);
357
358 if (levelI != levelJ)
359 {
360 DimVector vel(0);
361 cellData.fluxData().setVelocity(wPhaseIdx, indexInInside, vel);
362 cellData.fluxData().setVelocity(nPhaseIdx, indexInInside, vel);
363 cellData.fluxData().setUpwindPotential(wPhaseIdx, indexInInside, 0);
364 cellData.fluxData().setUpwindPotential(nPhaseIdx, indexInInside, 0);
365
366 cellDataJ.fluxData().setVelocity(wPhaseIdx, indexInOutside, vel);
367 cellDataJ.fluxData().setVelocity(nPhaseIdx, indexInOutside, vel);
368 cellDataJ.fluxData().setUpwindPotential(wPhaseIdx, indexInOutside, 0);
369 cellDataJ.fluxData().setUpwindPotential(nPhaseIdx, indexInOutside, 0);
370 }
371
372 for (int vIdx = 0; vIdx < numVertices; vIdx++)
373 {
374 int localVertIdx = refElement.subEntity(fIdx, dim - 1, vIdx, dim);
375
376 int vIdxGlobal = 0;
377 if (levelI >= levelJ)
378 {
379 vIdxGlobal = problem_.variables().index(elementI.template subEntity<dim>(localVertIdx));
380 }
381 else
382 {
383 vIdxGlobal = problem_.variables().index(elementJ.template subEntity<dim>(localVertIdx));
384 }
385
386 InteractionVolume& interactionVolume = this->interactionVolumes_[vIdxGlobal];
387
388 if (interactionVolume.isInnerVolume())
389 {
390 // cell index vector
391 std::vector<int> eIdxGlobal(0);
392
393 if (interactionVolume.getElementNumber() == 4)
394 {
395 auto element1 = interactionVolume.getSubVolumeElement(0);
396 auto element2 = interactionVolume.getSubVolumeElement(1);
397 auto element3 = interactionVolume.getSubVolumeElement(2);
398 auto element4 = interactionVolume.getSubVolumeElement(3);
399
400 eIdxGlobal.resize(4);
401
402 eIdxGlobal[0] = problem_.variables().index(element1);
403 eIdxGlobal[1] = problem_.variables().index(element2);
404 eIdxGlobal[2] = problem_.variables().index(element3);
405 eIdxGlobal[3] = problem_.variables().index(element4);
406
407 //cell Data vector
408 cellDataTemp.resize(4);
409
410 cellDataTemp[0] = problem_.variables().cellData(eIdxGlobal[0]);
411 cellDataTemp[1] = problem_.variables().cellData(eIdxGlobal[1]);
412 cellDataTemp[2] = problem_.variables().cellData(eIdxGlobal[2]);
413 cellDataTemp[3] = problem_.variables().cellData(eIdxGlobal[3]);
414
415 velocity_.calculateInnerInteractionVolumeVelocity(interactionVolume, cellDataTemp[0], cellDataTemp[1],
416 cellDataTemp[2], cellDataTemp[3], this->innerBoundaryVolumeFaces_);
417 }
418 else if (interactionVolume.getElementNumber() == 3)
419 {
420 auto element1 = interactionVolume.getSubVolumeElement(0);
421 auto element2 = interactionVolume.getSubVolumeElement(1);
422 auto element4 = interactionVolume.getSubVolumeElement(3);
423
424 eIdxGlobal.resize(3);
425
426 eIdxGlobal[0] = problem_.variables().index(element1);
427 eIdxGlobal[1] = problem_.variables().index(element2);
428 eIdxGlobal[2] = problem_.variables().index(element4);
429
430 //cell Data vector
431 cellDataTemp.resize(3);
432
433 cellDataTemp[0] = problem_.variables().cellData(eIdxGlobal[0]);
434 cellDataTemp[1] = problem_.variables().cellData(eIdxGlobal[1]);
435 cellDataTemp[2] = problem_.variables().cellData(eIdxGlobal[2]);
436
437
438 velocity_.calculateHangingNodeInteractionVolumeVelocity(interactionVolume, cellDataTemp[0], cellDataTemp[1],
439 cellDataTemp[2], this->innerBoundaryVolumeFaces_);
440 }
441 else
442 {
443 DUNE_THROW(Dune::NotImplemented, "Unknown interactionvolume type!");
444 }
445
446 int size = cellDataTemp.size();
447 for (int i = 0; i < size; i++)
448 {
449 if (eIdxGlobal[i] == eIdxGlobalI)
450 {
451 if (levelI >= levelJ)
452 {
453 cellData.fluxData().setVelocity(wPhaseIdx, indexInInside,
454 cellDataTemp[i].fluxData().velocity(wPhaseIdx, indexInInside));
455 cellData.fluxData().setVelocity(nPhaseIdx, indexInInside,
456 cellDataTemp[i].fluxData().velocity(nPhaseIdx, indexInInside));
457 cellData.fluxData().setUpwindPotential(wPhaseIdx, indexInInside,
458 cellDataTemp[i].fluxData().upwindPotential(wPhaseIdx, indexInInside));
459 cellData.fluxData().setUpwindPotential(nPhaseIdx, indexInInside,
460 cellDataTemp[i].fluxData().upwindPotential(nPhaseIdx, indexInInside));
461
462 if (levelI > levelJ)
463 {
464 cellDataJ.fluxData().setVelocity(wPhaseIdx, indexInOutside,
465 cellDataTemp[i].fluxData().velocity(wPhaseIdx, indexInInside));
466 cellDataJ.fluxData().setVelocity(nPhaseIdx, indexInOutside,
467 cellDataTemp[i].fluxData().velocity(nPhaseIdx, indexInInside));
468 cellDataJ.fluxData().setUpwindPotential(wPhaseIdx, indexInOutside,
469 cellDataTemp[i].fluxData().upwindPotential(wPhaseIdx, indexInInside));
470 cellDataJ.fluxData().setUpwindPotential(nPhaseIdx, indexInOutside,
471 cellDataTemp[i].fluxData().upwindPotential(nPhaseIdx, indexInInside));
472
473 }
474 }
475
476 }
477 else if (eIdxGlobal[i] == eIdxGlobalJ)
478 {
479 if (levelJ >= levelI)
480 {
481 cellDataJ.fluxData().setVelocity(wPhaseIdx, indexInOutside,
482 cellDataTemp[i].fluxData().velocity(wPhaseIdx, indexInOutside));
483 cellDataJ.fluxData().setVelocity(nPhaseIdx, indexInOutside,
484 cellDataTemp[i].fluxData().velocity(nPhaseIdx, indexInOutside));
485 cellDataJ.fluxData().setUpwindPotential(wPhaseIdx, indexInOutside,
486 cellDataTemp[i].fluxData().upwindPotential(wPhaseIdx, indexInOutside));
487 cellDataJ.fluxData().setUpwindPotential(nPhaseIdx, indexInOutside,
488 cellDataTemp[i].fluxData().upwindPotential(nPhaseIdx, indexInOutside));
489
490 if (levelJ > levelI)
491 {
492 cellData.fluxData().setVelocity(wPhaseIdx, indexInInside,
493 cellDataTemp[i].fluxData().velocity(wPhaseIdx, indexInOutside));
494 cellData.fluxData().setVelocity(nPhaseIdx, indexInInside,
495 cellDataTemp[i].fluxData().velocity(nPhaseIdx, indexInOutside));
496 cellData.fluxData().setUpwindPotential(wPhaseIdx, indexInInside,
497 cellDataTemp[i].fluxData().upwindPotential(wPhaseIdx, indexInOutside));
498 cellData.fluxData().setUpwindPotential(nPhaseIdx, indexInInside,
499 cellDataTemp[i].fluxData().upwindPotential(nPhaseIdx, indexInOutside));
500
501 }
502 }
503 }
504 }
505 }
506 }
507 if (levelI == levelJ)
508 {
509 cellData.fluxData().setVelocityMarker(indexInInside);
510 cellDataJ.fluxData().setVelocityMarker(indexInOutside);
511 }
512}
513
522template<class TypeTag>
523void FvMpfaL2dPressureVelocity2pAdaptive<TypeTag>::calculateVelocityOnBoundary(const Intersection& intersection, CellData& cellData)
524{
525 auto element = intersection.inside();
526
527 //get face index
528 int isIndex = intersection.indexInInside();
529
530 //get face normal
531 const Dune::FieldVector<Scalar, dim>& unitOuterNormal = intersection.centerUnitOuterNormal();
532
533 BoundaryTypes bcType;
534 //get boundary type
535 problem_.boundaryTypes(bcType, intersection);
536 PrimaryVariables boundValues(0.0);
537
538 if (bcType.isDirichlet(pressEqIdx))
539 {
540 problem_.dirichlet(boundValues, intersection);
541
542 // get global coordinates of cell centers
543 const GlobalPosition& globalPosI = element.geometry().center();
544
545 // center of face in global coordinates
546 const GlobalPosition& globalPosJ = intersection.geometry().center();
547
548 // get mobilities and fractional flow factors
549 Scalar lambdaWI = cellData.mobility(wPhaseIdx);
550 Scalar lambdaNwI = cellData.mobility(nPhaseIdx);
551
552 // get capillary pressure
553 Scalar pcI = cellData.capillaryPressure();
554
555 // distance vector between barycenters
556 GlobalPosition distVec = globalPosJ - globalPosI;
557
558 // compute distance between cell centers
559 Scalar dist = distVec.two_norm();
560
561 //permeability vector at boundary
562 // compute vectorized permeabilities
563 DimMatrix meanPermeability(0);
564
565 problem_.spatialParams().meanK(meanPermeability, problem_.spatialParams().intrinsicPermeability(element));
566
567 Dune::FieldVector<Scalar, dim> permeability(0);
568 meanPermeability.mv(unitOuterNormal, permeability);
569
570 //determine saturation at the boundary -> if no saturation is known directly at the boundary use the cell saturation
571 Scalar satW = 0;
572 if (bcType.isDirichlet(satEqIdx))
573 {
574 switch (saturationType_)
575 {
576 case sw:
577 {
578 satW = boundValues[saturationIdx];
579 break;
580 }
581 case sn:
582 {
583 satW = 1 - boundValues[saturationIdx];
584 break;
585 }
586 }
587 }
588 else
589 {
590 satW = cellData.saturation(wPhaseIdx);
591 }
592
593 const Scalar pressBound = boundValues[pressureIdx];
594
595 // old material law interface is deprecated: Replace this by
596 // const auto& fluidMatrixInteraction = spatialParams.fluidMatrixInteractionAtPos(element.geometry().center());
597 // after the release of 3.3, when the deprecated interface is no longer supported
598 const auto fluidMatrixInteraction = Deprecated::makePcKrSw(Scalar{}, problem_.spatialParams(), element);
599
600 const Scalar pcBound = fluidMatrixInteraction.pc(satW);
601
602 //determine phase pressures from primary pressure variable
603 Scalar pressWBound = 0;
604 Scalar pressNwBound = 0;
605 if (pressureType_ == pw)
606 {
607 pressWBound = pressBound;
608 pressNwBound = pressBound + pcBound;
609 }
610 else if (pressureType_ == pn)
611 {
612 pressWBound = pressBound - pcBound;
613 pressNwBound = pressBound;
614 }
615
616 const Scalar lambdaWBound = fluidMatrixInteraction.krw(satW) / viscosity_[wPhaseIdx];
617 const Scalar lambdaNwBound = fluidMatrixInteraction.krn(satW) / viscosity_[nPhaseIdx];
618
619 Scalar potentialDiffW = cellData.fluxData().upwindPotential(wPhaseIdx, isIndex);
620 Scalar potentialDiffNw = cellData.fluxData().upwindPotential(nPhaseIdx, isIndex);
621
622 //calculate potential gradient
623 potentialDiffW = (cellData.pressure(wPhaseIdx) - pressWBound);
624 potentialDiffNw = (cellData.pressure(nPhaseIdx) - pressNwBound);
625
626 potentialDiffW += density_[wPhaseIdx] * (distVec * problem_.gravity());
627 potentialDiffNw += density_[nPhaseIdx] * (distVec * problem_.gravity());
628
629 //store potential gradients for further calculations
630 cellData.fluxData().setUpwindPotential(wPhaseIdx, isIndex, potentialDiffW);
631 cellData.fluxData().setUpwindPotential(nPhaseIdx, isIndex, potentialDiffNw);
632
633 //do the upwinding of the mobility depending on the phase potentials
634 Scalar lambdaW = (potentialDiffW > 0.) ? lambdaWI : lambdaWBound;
635 lambdaW = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffW, 0.0, 1.0e-30)) ? 0.5 * (lambdaWI + lambdaWBound) : lambdaW;
636 Scalar lambdaNw = (potentialDiffNw > 0.) ? lambdaNwI : lambdaNwBound;
637 lambdaNw = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffNw, 0.0, 1.0e-30)) ? 0.5 * (lambdaNwI + lambdaNwBound) : lambdaNw;
638
639
640 Scalar scalarPerm = permeability.two_norm();
641
642 //calculate the gravity term
643 Dune::FieldVector<Scalar, dimWorld> velocityW(unitOuterNormal);
644 Dune::FieldVector<Scalar, dimWorld> velocityNw(unitOuterNormal);
645
646 //calculate unit distVec
647 distVec /= dist;
648 Scalar areaScaling = (unitOuterNormal * distVec);
649 //this treatment of g allows to account for gravity flux through faces where the face normal
650 //has no z component (e.g. parallelepiped grids)
651 Scalar gravityTermW = (problem_.gravity() * distVec) * density_[wPhaseIdx] * areaScaling;
652 Scalar gravityTermNw = (problem_.gravity() * distVec) * density_[nPhaseIdx] * areaScaling;
653
654 //calculate velocity depending on the pressure used -> use pc = pn - pw
655 switch (pressureType_)
656 {
657 case pw:
658 {
659 velocityW *= lambdaW * scalarPerm * ((cellData.pressure(wPhaseIdx) - pressBound) / dist + gravityTermW);
660 velocityNw *= lambdaNw * scalarPerm * ((cellData.pressure(wPhaseIdx) - pressBound) / dist + gravityTermNw)
661 + 0.5 * (lambdaNwI + lambdaNwBound) * scalarPerm * (pcI - pcBound) / dist;
662 break;
663 }
664 case pn:
665 {
666 velocityW *= lambdaW * scalarPerm * ((cellData.pressure(nPhaseIdx) - pressBound) / dist + gravityTermW)
667 - 0.5 * (lambdaWI + lambdaWBound) * scalarPerm * (pcI - pcBound) / dist;
668 velocityNw *= lambdaNw * scalarPerm * ((cellData.pressure(nPhaseIdx) - pressBound) / dist + gravityTermNw);
669 break;
670 }
671 }
672
673 //store velocities
674 cellData.fluxData().setVelocity(wPhaseIdx, isIndex, velocityW);
675 cellData.fluxData().setVelocity(nPhaseIdx, isIndex, velocityNw);
676 cellData.fluxData().setVelocityMarker(isIndex);
677
678 } //end dirichlet boundary
679
680 else if (bcType.isNeumann(pressEqIdx))
681 {
682 problem_.neumann(boundValues, intersection);
683
684 Dune::FieldVector<Scalar, dimWorld> velocityW(unitOuterNormal);
685 Dune::FieldVector<Scalar, dimWorld> velocityNw(unitOuterNormal);
686
687 velocityW *= boundValues[wPhaseIdx];
688 velocityNw *= boundValues[nPhaseIdx];
689
690 velocityW /= density_[wPhaseIdx];
691 velocityNw /= density_[nPhaseIdx];
692
693 //store potential gradients for further calculations
694 cellData.fluxData().setUpwindPotential(wPhaseIdx, isIndex, boundValues[wPhaseIdx]);
695 cellData.fluxData().setUpwindPotential(nPhaseIdx, isIndex, boundValues[nPhaseIdx]);
696
697 cellData.fluxData().setVelocity(wPhaseIdx, isIndex, velocityW);
698 cellData.fluxData().setVelocity(nPhaseIdx, isIndex, velocityNw);
699 cellData.fluxData().setVelocityMarker(isIndex);
700 } //end neumann boundary
701 else
702 {
703 DUNE_THROW(Dune::NotImplemented, "No valid boundary condition type defined for pressure equation!");
704 }
705}
706
707}
708// end of Dune namespace
709#endif
Helpers for deprecation.
Grid adaptive finite volume MPFA L-method discretization of a two-phase pressure equation of the sequ...
Velocity calculation on adaptive grids using a 2-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
Grid adaptive finite volume MPFA L-method discretization of a two-phase flow pressure equation of the...
Definition: 2dpressureadaptive.hh:77
FVMPFALInteractionVolume< TypeTag > InteractionVolume
Type of the interaction volume objects.
Definition: 2dpressureadaptive.hh:154
void addOutputVtkFields(MultiWriter &writer)
Adds pressure output to the output file.
Definition: 2dpressureadaptive.hh:376
void storePressureSolution()
Globally stores the pressure solution.
Definition: 2dpressureadaptive.hh:233
void initialize()
Initializes the pressure model.
Definition: 2dpressureadaptive.hh:202
void updateMaterialLaws()
constitutive functions are initialized and stored in the variables object
Definition: 2dpressureadaptive.hh:2552
void update()
Pressure update.
Definition: 2dpressureadaptive.hh:300
Class for the calculation of velocities from the pressure solution of an IMPES scheme using a grid ad...
Definition: 2dpressurevelocityadaptive.hh:50
void addOutputVtkFields(MultiWriter &writer)
Adds velocity output to the output file.
Definition: 2dpressurevelocityadaptive.hh:202
void updateVelocity()
Function for updating the velocity field if iterations are necessary in the transport solution.
Definition: 2dpressurevelocityadaptive.hh:132
bool calculateVelocityInTransport()
Indicates if velocity is reconstructed in the pressure step or in the transport step.
Definition: 2dpressurevelocityadaptive.hh:186
FvMpfaL2dPressureVelocity2pAdaptive(Problem &problem)
Constructs a FvMpfaL2dPressureVelocity2pAdaptive object.
Definition: 2dpressurevelocityadaptive.hh:113
void initialize()
Initializes pressure and velocity.
Definition: 2dpressurevelocityadaptive.hh:147
void calculateVelocityOnBoundary(const Intersection &intersection, CellData &cellData)
Calculates the velocity at a boundary.
Definition: 2dpressurevelocityadaptive.hh:523
void calculateVelocity()
Calculates the velocities at all cell-cell interfaces.
Definition: 2dpressurevelocityadaptive.hh:229
void update()
Pressure and velocity update.
Definition: 2dpressurevelocityadaptive.hh:174
Class for calculating 2-d velocities from cell-wise constant pressure values.
Definition: 2dvelocityadaptive.hh:54
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
int getElementNumber()
Get number of stored elements.
Definition: linteractionvolume.hh:246
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