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