3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
3dpressurevelocityadaptive.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2/*****************************************************************************
3 * See the file COPYING for full copying permissions. *
4 * *
5 * This program is free software: you can redistribute it and/or modify *
6 * it under the terms of the GNU General Public License as published by *
7 * the Free Software Foundation, either version 3 of the License, or *
8 * (at your option) any later version. *
9 * *
10 * This program is distributed in the hope that it will be useful, *
11 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
12 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
13 * GNU General Public License for more details. *
14 * *
15 * You should have received a copy of the GNU General Public License *
16 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
17 *****************************************************************************/
23#ifndef DUMUX_FVMPFALPRESSUREVELOCITY2P_ADAPTIVE_HH
24#define DUMUX_FVMPFALPRESSUREVELOCITY2P_ADAPTIVE_HH
25
26#include <dune/common/float_cmp.hh>
27#include "3dpressureadaptive.hh"
28#include "3dvelocityadaptive.hh"
29
31
32namespace Dumux {
33
47template<class TypeTag> class FvMpfaL3dPressureVelocity2pAdaptive: public FvMpfaL3dPressure2pAdaptive<TypeTag>
48{
51
53
54 enum
55 {
56 dim = GridView::dimension, dimWorld = GridView::dimensionworld
57 };
58
60
61 enum
62 {
63 pw = Indices::pressureW,
64 pn = Indices::pressureNw,
65 sw = Indices::saturationW,
66 sn = Indices::saturationNw,
67 wPhaseIdx = Indices::wPhaseIdx,
68 nPhaseIdx = Indices::nPhaseIdx,
69 pressureIdx = Indices::pressureIdx,
70 saturationIdx = Indices::saturationIdx,
71 pressEqIdx = Indices::pressureEqIdx,
72 satEqIdx = Indices::satEqIdx,
73 numPhases = getPropValue<TypeTag, Properties::NumPhases>()
74 };
75
82 using PrimaryVariables = typename SolutionTypes::PrimaryVariables;
84
86 using Intersection = typename GridView::Intersection;
87
88 using Element = typename GridView::template Codim<0>::Entity;
89
90 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
91 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
92 using DimVector = Dune::FieldVector<Scalar, dim>;
93
94public:
101 ParentType(problem), problem_(problem), velocity_(problem)
102{
103 density_[wPhaseIdx] = 0.;
104 density_[nPhaseIdx] = 0.;
105 viscosity_[wPhaseIdx] = 0.;
106 viscosity_[nPhaseIdx] = 0.;
107
108 calcVelocityInTransport_ = getParam<bool>("MPFA.CalcVelocityInTransport");
109}
110
111 void calculateVelocity();
112
114 void calculateVelocity(const Intersection&, CellData&);
115 void calculateVelocityOnBoundary(const Intersection& intersection, CellData& cellData);
116
119 {
120 this->updateMaterialLaws();
121
122 this->storePressureSolution();
123
126 }
127
133 void initialize(bool solveTwice = true)
134 {
135 const auto element = *problem_.gridView().template begin<0>();
136 FluidState fluidState;
137 fluidState.setPressure(wPhaseIdx, problem_.referencePressure(element));
138 fluidState.setPressure(nPhaseIdx, problem_.referencePressure(element));
139 fluidState.setTemperature(problem_.temperature(element));
140 fluidState.setSaturation(wPhaseIdx, 1.);
141 fluidState.setSaturation(nPhaseIdx, 0.);
142 density_[wPhaseIdx] = FluidSystem::density(fluidState, wPhaseIdx);
143 density_[nPhaseIdx] = FluidSystem::density(fluidState, nPhaseIdx);
144 viscosity_[wPhaseIdx] = FluidSystem::viscosity(fluidState, wPhaseIdx);
145 viscosity_[nPhaseIdx] = FluidSystem::viscosity(fluidState, nPhaseIdx);
146
148 velocity_.initialize();
149
152
153 return;
154 }
155
161 void update()
162 {
164
167 }
168
175 {
176 return calcVelocityInTransport_;
177 }
178
189 template<class MultiWriter>
190 void addOutputVtkFields(MultiWriter &writer)
191 {
193 velocity_.addOutputVtkFields(writer);
194 }
195
196private:
197 Problem& problem_;
199
200 Scalar density_[numPhases];
201 Scalar viscosity_[numPhases];
202 bool calcVelocityInTransport_;
203
205 static const int pressureType_ = getPropValue<TypeTag, Properties::PressureFormulation>();
207 static const int saturationType_ = getPropValue<TypeTag, Properties::SaturationFormulation>();
208};
209// end of template
210
216template<class TypeTag>
218{
219 // run through all vertices
220 for (const auto& vertex : vertices(problem_.gridView()))
221 {
222 int vIdxGlobal = problem_.variables().index(vertex);
223
224 InteractionVolume& interactionVolume = this->interactionVolumes_.interactionVolume(vIdxGlobal);
225
226 // inner interactionvolume
227 if (interactionVolume.isInnerVolume())
228 {
229 // cell index
230 int eIdxGlobal[8];
231 for (int i = 0; i < 8; i++)
232 {
233 eIdxGlobal[i] = problem_.variables().index(interactionVolume.getSubVolumeElement(i));
234 }
235
236 //get the cell Data
237 CellData & cellData1 = problem_.variables().cellData(eIdxGlobal[0]);
238 CellData & cellData2 = problem_.variables().cellData(eIdxGlobal[1]);
239 CellData & cellData3 = problem_.variables().cellData(eIdxGlobal[2]);
240 CellData & cellData4 = problem_.variables().cellData(eIdxGlobal[3]);
241 CellData & cellData5 = problem_.variables().cellData(eIdxGlobal[4]);
242 CellData & cellData6 = problem_.variables().cellData(eIdxGlobal[5]);
243 CellData & cellData7 = problem_.variables().cellData(eIdxGlobal[6]);
244 CellData & cellData8 = problem_.variables().cellData(eIdxGlobal[7]);
245
246 if (!interactionVolume.isHangingNodeVolume())
247 {
248 velocity_.calculateInnerInteractionVolumeVelocity(interactionVolume,
249 cellData1, cellData2, cellData3, cellData4,
250 cellData5, cellData6, cellData7, cellData8,
251 this->interactionVolumes_, this->transmissibilityCalculator_);
252 }
253 else
254 {
255 velocity_.calculateHangingNodeInteractionVolumeVelocity(interactionVolume,
256 cellData1, cellData2, cellData3, cellData4,
257 cellData5, cellData6, cellData7, cellData8);
258 }
259 }
260 // at least one face on boundary! (boundary interactionvolume)
261 else
262 {
263 for (int elemIdx = 0; elemIdx < 8; elemIdx++)
264 {
265 if (!interactionVolume.hasSubVolumeElement(elemIdx))
266 {
267 continue;
268 }
269 bool isOutside = false;
270 for (int fIdx = 0; fIdx < dim; fIdx++)
271 {
272 int intVolFaceIdx = interactionVolume.getFaceIndexFromSubVolume(elemIdx, fIdx);
273 if (interactionVolume.isOutsideFace(intVolFaceIdx))
274 {
275 isOutside = true;
276 break;
277 }
278 }
279 if (isOutside)
280 {
281 continue;
282 }
283
284 // cell index
285 int eIdxGlobal = problem_.variables().index(interactionVolume.getSubVolumeElement(elemIdx));
286 //get the cell Data
287 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
288
289 velocity_.calculateBoundaryInteractionVolumeVelocity(interactionVolume, cellData, elemIdx);
290 }
291 } // end boundaries
292
293 } // end vertex iterator
294
295 return;
296}
297
306template<class TypeTag>
307void FvMpfaL3dPressureVelocity2pAdaptive<TypeTag>::calculateVelocity(const Intersection& intersection, CellData& cellData)
308{
309 auto elementI = intersection.inside();
310 auto elementJ = intersection.outside();
311
312 int levelI = elementI.level();
313 int levelJ = elementJ.level();
314
315 int eIdxGlobalI = problem_.variables().index(elementI);
316 int eIdxGlobalJ = problem_.variables().index(elementJ);
317
318 CellData& cellDataJ = problem_.variables().cellData(eIdxGlobalJ);
319
320 int indexInInside = intersection.indexInInside();
321 int indexInOutside = intersection.indexInOutside();
322
323 Dune::FieldVector<CellData, 8> cellDataTemp;
324
325 if (levelI != levelJ)
326 {
327 DimVector vel(0);
328 cellData.fluxData().setVelocity(wPhaseIdx, indexInInside, vel);
329 cellData.fluxData().setVelocity(nPhaseIdx, indexInInside, vel);
330 cellData.fluxData().setUpwindPotential(wPhaseIdx, indexInInside, 0);
331 cellData.fluxData().setUpwindPotential(nPhaseIdx, indexInInside, 0);
332
333 cellDataJ.fluxData().setVelocity(wPhaseIdx, indexInOutside, vel);
334 cellDataJ.fluxData().setVelocity(nPhaseIdx, indexInOutside, vel);
335 cellDataJ.fluxData().setUpwindPotential(wPhaseIdx, indexInOutside, 0);
336 cellDataJ.fluxData().setUpwindPotential(nPhaseIdx, indexInOutside, 0);
337 }
338
339 std::set<int> vIdxGlobal;
340
341 if (levelI >= levelJ)
342 {
343 vIdxGlobal = this->interactionVolumes_.faceVerticeIndices(eIdxGlobalI, indexInInside);
344 }
345 else
346 {
347 vIdxGlobal = this->interactionVolumes_.faceVerticeIndices(eIdxGlobalJ, indexInOutside);
348 }
349
350 std::set<int>::iterator itEnd = vIdxGlobal.end();
351
352 for (std::set<int>::iterator vIdxIt = vIdxGlobal.begin(); vIdxIt != itEnd; ++vIdxIt)
353 {
354 InteractionVolume& interactionVolume = this->interactionVolumes_.interactionVolume(*vIdxIt);
355
356 if (interactionVolume.isInnerVolume())
357 {
358 // cell index
359 std::vector<std::pair<int,int> > localMpfaElemIdx(0);
360
361 int eIdxGlobal[8];
362 for (int i = 0; i < 8; i++)
363 {
364 auto elem = interactionVolume.getSubVolumeElement(i);
365
366 if (interactionVolume.isHangingNodeVolume())
367 {
368 if (elem == elementI)
369 {
370 for (int j = 0; j < 3; j++)
371 {
372 for (int k = 0; k < 8; k++)
373 {
374 if (interactionVolume.getSubVolumeElement(k) == elementJ &&
375 IndexTranslator::getFaceIndexFromSubVolume(i,j)
376 == IndexTranslator::getFaceIndexFromElements(i,k))
377 {
378 std::pair<int,int> localIdx = std::make_pair(i, k);
379 localMpfaElemIdx.push_back(localIdx);
380 }
381 }
382 }
383 }
384 }
385 else
386 {
387 if (localMpfaElemIdx.size() != 1)
388 localMpfaElemIdx.resize(1);
389
390 if (elem == elementI)
391 localMpfaElemIdx[0].first = i;
392 else if (elem == elementJ)
393 localMpfaElemIdx[0].second = i;
394 }
395
396 eIdxGlobal[i] = problem_.variables().index(elem);
397 cellDataTemp[i] = problem_.variables().cellData(eIdxGlobal[i]);
398 }
399
400 int size = localMpfaElemIdx.size();
401
402 // if (size > 1)
403 // std::cout<<"size = "<<size<<"\n";
404
405 for (int i = 0; i < size; i++)
406 {
407 int mpfaFaceIdx = IndexTranslator::getFaceIndexFromElements(localMpfaElemIdx[i].first, localMpfaElemIdx[i].second);
408
409 if (!interactionVolume.isHangingNodeVolume())
410 {
411 velocity_.calculateInnerInteractionVolumeVelocity(interactionVolume,
412 cellDataTemp[0], cellDataTemp[1], cellDataTemp[2], cellDataTemp[3],
413 cellDataTemp[4], cellDataTemp[5], cellDataTemp[6], cellDataTemp[7],
414 this->interactionVolumes_, this->transmissibilityCalculator_, mpfaFaceIdx);
415 }
416 else
417 {
418 velocity_.calculateHangingNodeInteractionVolumeVelocity(interactionVolume,
419 cellDataTemp[0], cellDataTemp[1], cellDataTemp[2], cellDataTemp[3],
420 cellDataTemp[4], cellDataTemp[5], cellDataTemp[6], cellDataTemp[7], mpfaFaceIdx);
421 }
422
423 if (levelI >= levelJ)
424 {
425 cellData.fluxData().setVelocity(wPhaseIdx, indexInInside,
426 cellDataTemp[localMpfaElemIdx[i].first].fluxData().velocity(wPhaseIdx, indexInInside));
427 cellData.fluxData().setVelocity(nPhaseIdx, indexInInside,
428 cellDataTemp[localMpfaElemIdx[i].first].fluxData().velocity(nPhaseIdx, indexInInside));
429 cellData.fluxData().setUpwindPotential(wPhaseIdx, indexInInside,
430 cellDataTemp[localMpfaElemIdx[i].first].fluxData().upwindPotential(wPhaseIdx, indexInInside));
431 cellData.fluxData().setUpwindPotential(nPhaseIdx, indexInInside,
432 cellDataTemp[localMpfaElemIdx[i].first].fluxData().upwindPotential(nPhaseIdx, indexInInside));
433 }
434 if (levelJ >= levelI)
435 {
436 cellDataJ.fluxData().setVelocity(wPhaseIdx, indexInOutside,
437 cellDataTemp[localMpfaElemIdx[i].second].fluxData().velocity(wPhaseIdx, indexInOutside));
438 cellDataJ.fluxData().setVelocity(nPhaseIdx, indexInOutside,
439 cellDataTemp[localMpfaElemIdx[i].second].fluxData().velocity(nPhaseIdx, indexInOutside));
440 cellDataJ.fluxData().setUpwindPotential(wPhaseIdx, indexInOutside,
441 cellDataTemp[localMpfaElemIdx[i].second].fluxData().upwindPotential(wPhaseIdx, indexInOutside));
442 cellDataJ.fluxData().setUpwindPotential(nPhaseIdx, indexInOutside,
443 cellDataTemp[localMpfaElemIdx[i].second].fluxData().upwindPotential(nPhaseIdx, indexInOutside));
444 }
445
446 if (size > 1)
447 {
448 for (int j = 0; j < 8; j++)
449 {
450 cellDataTemp[j] = problem_.variables().cellData(eIdxGlobal[j]);
451 }
452 }
453 }
454 }
455 }
456
457 if (levelI == levelJ)
458 {
459 cellData.fluxData().setVelocityMarker(indexInInside);
460 cellDataJ.fluxData().setVelocityMarker(indexInOutside);
461 }
462 else if (levelI > levelJ)
463 {
464 cellDataJ.fluxData().setVelocity(wPhaseIdx, indexInOutside, cellData.fluxData().velocity(wPhaseIdx, indexInInside));
465 cellDataJ.fluxData().setVelocity(nPhaseIdx, indexInOutside, cellData.fluxData().velocity(nPhaseIdx, indexInInside));
466 cellDataJ.fluxData().setUpwindPotential(wPhaseIdx, indexInOutside, -1*cellData.fluxData().upwindPotential(wPhaseIdx, indexInInside));
467 cellDataJ.fluxData().setUpwindPotential(nPhaseIdx, indexInOutside, -1*cellData.fluxData().upwindPotential(nPhaseIdx, indexInInside));
468 }
469 else if (levelJ > levelI)
470 {
471 cellData.fluxData().setVelocity(wPhaseIdx, indexInInside, cellDataJ.fluxData().velocity(wPhaseIdx, indexInOutside));
472 cellData.fluxData().setVelocity(nPhaseIdx, indexInInside, cellDataJ.fluxData().velocity(nPhaseIdx, indexInOutside));
473 cellData.fluxData().setUpwindPotential(wPhaseIdx, indexInInside, -1*cellDataJ.fluxData().upwindPotential(wPhaseIdx, indexInOutside));
474 cellData.fluxData().setUpwindPotential(nPhaseIdx, indexInInside, -1*cellDataJ.fluxData().upwindPotential(nPhaseIdx, indexInOutside));
475 }
476}
477
486template<class TypeTag>
487void FvMpfaL3dPressureVelocity2pAdaptive<TypeTag>::calculateVelocityOnBoundary(const Intersection& intersection, CellData& cellData)
488{
489 auto element = intersection.inside();
490
491 //get face index
492 int isIndex = intersection.indexInInside();
493
494 //get face normal
495 const Dune::FieldVector<Scalar, dim>& unitOuterNormal = intersection.centerUnitOuterNormal();
496
497 BoundaryTypes bcType;
498 //get boundary type
499 problem_.boundaryTypes(bcType, intersection);
500 PrimaryVariables boundValues(0.0);
501
502 if (bcType.isDirichlet(pressEqIdx))
503 {
504 problem_.dirichlet(boundValues, intersection);
505
506 // get global coordinates of cell centers
507 const GlobalPosition& globalPosI = element.geometry().center();
508
509 // center of face in global coordinates
510 const GlobalPosition& globalPosJ = intersection.geometry().center();
511
512 // get mobilities and fractional flow factors
513 Scalar lambdaWI = cellData.mobility(wPhaseIdx);
514 Scalar lambdaNwI = cellData.mobility(nPhaseIdx);
515
516 // get capillary pressure
517 Scalar pcI = cellData.capillaryPressure();
518
519 // distance vector between barycenters
520 GlobalPosition distVec = globalPosJ - globalPosI;
521
522 // compute distance between cell centers
523 Scalar dist = distVec.two_norm();
524
525 //permeability vector at boundary
526 // compute vectorized permeabilities
527 DimMatrix meanPermeability(0);
528
529 problem_.spatialParams().meanK(meanPermeability, problem_.spatialParams().intrinsicPermeability(element));
530
531 Dune::FieldVector<Scalar, dim> permeability(0);
532 meanPermeability.mv(unitOuterNormal, permeability);
533
534 //determine saturation at the boundary -> if no saturation is known directly at the boundary use the cell saturation
535 Scalar satW = 0;
536 if (bcType.isDirichlet(satEqIdx))
537 {
538 switch (saturationType_)
539 {
540 case sw:
541 {
542 satW = boundValues[saturationIdx];
543 break;
544 }
545 case sn:
546 {
547 satW = 1 - boundValues[saturationIdx];
548 break;
549 }
550 }
551 }
552 else
553 {
554 satW = cellData.saturation(wPhaseIdx);
555 }
556
557 const Scalar pressBound = boundValues[pressureIdx];
558
559 // old material law interface is deprecated: Replace this by
560 // const auto& fluidMatrixInteraction = spatialParams.fluidMatrixInteractionAtPos(element.geometry().center());
561 // after the release of 3.3, when the deprecated interface is no longer supported
562 const auto fluidMatrixInteraction = Deprecated::makePcKrSw(Scalar{}, problem_.spatialParams(), element);
563
564 const Scalar pcBound = fluidMatrixInteraction.pc(satW);
565
566 //determine phase pressures from primary pressure variable
567 Scalar pressWBound = 0;
568 Scalar pressNwBound = 0;
569 if (pressureType_ == pw)
570 {
571 pressWBound = pressBound;
572 pressNwBound = pressBound + pcBound;
573 }
574 else if (pressureType_ == pn)
575 {
576 pressWBound = pressBound - pcBound;
577 pressNwBound = pressBound;
578 }
579
580 const Scalar lambdaWBound = fluidMatrixInteraction.krw(satW) / viscosity_[wPhaseIdx];
581 const Scalar lambdaNwBound = fluidMatrixInteraction.krn(satW) / viscosity_[nPhaseIdx];
582
583 Scalar potentialDiffW = cellData.fluxData().upwindPotential(wPhaseIdx, isIndex);
584 Scalar potentialDiffNw = cellData.fluxData().upwindPotential(nPhaseIdx, isIndex);
585
586 //calculate potential gradient
587 potentialDiffW = (cellData.pressure(wPhaseIdx) - pressWBound);
588 potentialDiffNw = (cellData.pressure(nPhaseIdx) - pressNwBound);
589
590 potentialDiffW += density_[wPhaseIdx] * (distVec * problem_.gravity());
591 potentialDiffNw += density_[nPhaseIdx] * (distVec * problem_.gravity());
592
593 //store potential gradients for further calculations
594 cellData.fluxData().setUpwindPotential(wPhaseIdx, isIndex, potentialDiffW);
595 cellData.fluxData().setUpwindPotential(nPhaseIdx, isIndex, potentialDiffNw);
596
597 //do the upwinding of the mobility depending on the phase potentials
598 Scalar lambdaW = (potentialDiffW > 0.) ? lambdaWI : lambdaWBound;
599 lambdaW = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffW, 0.0, 1.0e-30)) ? 0.5 * (lambdaWI + lambdaWBound) : lambdaW;
600 Scalar lambdaNw = (potentialDiffNw > 0.) ? lambdaNwI : lambdaNwBound;
601 lambdaNw = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffNw, 0.0, 1.0e-30)) ? 0.5 * (lambdaNwI + lambdaNwBound) : lambdaNw;
602
603
604 Scalar scalarPerm = permeability.two_norm();
605
606 //calculate the gravity term
607 Dune::FieldVector<Scalar, dimWorld> velocityW(unitOuterNormal);
608 Dune::FieldVector<Scalar, dimWorld> velocityNw(unitOuterNormal);
609
610 //calculate unit distVec
611 distVec /= dist;
612 Scalar areaScaling = (unitOuterNormal * distVec);
613 //this treatment of g allows to account for gravity flux through faces where the face normal has no z component (e.g. parallelepiped grids)
614 Scalar gravityTermW = (problem_.gravity() * distVec) * density_[wPhaseIdx] * areaScaling;
615 Scalar gravityTermNw = (problem_.gravity() * distVec) * density_[nPhaseIdx] * areaScaling;
616
617 //calculate velocity depending on the pressure used -> use pc = pn - pw
618 switch (pressureType_)
619 {
620 case pw:
621 {
622 velocityW *= lambdaW * scalarPerm * ((cellData.pressure(wPhaseIdx) - pressBound) / dist + gravityTermW);
623 velocityNw *= lambdaNw * scalarPerm * ((cellData.pressure(wPhaseIdx) - pressBound) / dist + gravityTermNw)
624 + 0.5 * (lambdaNwI + lambdaNwBound) * scalarPerm * (pcI - pcBound) / dist;
625 break;
626 }
627 case pn:
628 {
629 velocityW *= lambdaW * scalarPerm * ((cellData.pressure(nPhaseIdx) - pressBound) / dist + gravityTermW)
630 - 0.5 * (lambdaWI + lambdaWBound) * scalarPerm * (pcI - pcBound) / dist;
631 velocityNw *= lambdaNw * scalarPerm * ((cellData.pressure(nPhaseIdx) - pressBound) / dist + gravityTermNw);
632 break;
633 }
634 }
635
636 //store velocities
637 cellData.fluxData().setVelocity(wPhaseIdx, isIndex, velocityW);
638 cellData.fluxData().setVelocity(nPhaseIdx, isIndex, velocityNw);
639 cellData.fluxData().setVelocityMarker(isIndex);
640
641 } //end dirichlet boundary
642
643 else if (bcType.isNeumann(pressEqIdx))
644 {
645 problem_.neumann(boundValues, intersection);
646
647 Dune::FieldVector<Scalar, dimWorld> velocityW(unitOuterNormal);
648 Dune::FieldVector<Scalar, dimWorld> velocityNw(unitOuterNormal);
649
650 velocityW *= boundValues[wPhaseIdx];
651 velocityNw *= boundValues[nPhaseIdx];
652
653 velocityW /= density_[wPhaseIdx];
654 velocityNw /= density_[nPhaseIdx];
655
656 //store potential gradients for further calculations
657 cellData.fluxData().setUpwindPotential(wPhaseIdx, isIndex, boundValues[wPhaseIdx]);
658 cellData.fluxData().setUpwindPotential(nPhaseIdx, isIndex, boundValues[nPhaseIdx]);
659
660 cellData.fluxData().setVelocity(wPhaseIdx, isIndex, velocityW);
661 cellData.fluxData().setVelocity(nPhaseIdx, isIndex, velocityNw);
662 cellData.fluxData().setVelocityMarker(isIndex);
663 } //end neumann boundary
664 else
665 {
666 DUNE_THROW(Dune::NotImplemented, "No valid boundary condition type defined for pressure equation!");
667 }
668}
669
670} // end namespace Dumux
671#endif
Helpers for deprecation.
3-d finite Volume-MPFAL implementation of a two-phase pressure equation on h-adaptive grids
3-d velocity calculation on adaptive grids using a 3-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
GetPropType< TypeTag, Properties::MPFAInteractionVolume > InteractionVolume
Type for storing interaction volume information.
Definition: 3dpressure.hh:156
void storePressureSolution()
Globally stores the pressure solution.
Definition: 3dpressure.hh:207
void addOutputVtkFields(MultiWriter &writer)
Adds pressure output to the output file.
Definition: 3dpressure.hh:386
void updateMaterialLaws()
constitutive functions are initialized and stored in the variables object
Definition: 3dpressure.hh:2484
3-d finite volume MPFA L-method discretization of a two-phase flow pressure equation of the sequentia...
Definition: 3dpressureadaptive.hh:78
void update()
updates the pressure field (analog to update function in IMPET)
Definition: 3dpressureadaptive.hh:194
Class for the calculation of 3d velocities from the pressure solution of an IMPES scheme using a grid...
Definition: 3dpressurevelocityadaptive.hh:48
void calculateVelocity()
Calculates the velocities at a cell-cell interfaces for the entire simulation grid.
Definition: 3dpressurevelocityadaptive.hh:217
FvMpfaL3dPressureVelocity2pAdaptive(Problem &problem)
Constructs a FvMpfaL3dPressureVelocity2pAdaptive object.
Definition: 3dpressurevelocityadaptive.hh:100
void updateVelocity()
Function for updating the velocity field if iterations are necessary in the transport solution.
Definition: 3dpressurevelocityadaptive.hh:118
void calculateVelocityOnBoundary(const Intersection &intersection, CellData &cellData)
Calculates the velocity at a boundary.
Definition: 3dpressurevelocityadaptive.hh:487
bool calculateVelocityInTransport()
Indicates if velocity is reconstructed in the pressure step or in the transport step.
Definition: 3dpressurevelocityadaptive.hh:174
void initialize(bool solveTwice=true)
Initializes pressure and velocity.
Definition: 3dpressurevelocityadaptive.hh:133
void addOutputVtkFields(MultiWriter &writer)
Adds velocity output to the output file.
Definition: 3dpressurevelocityadaptive.hh:190
void update()
Pressure and velocity update.
Definition: 3dpressurevelocityadaptive.hh:161
Class for calculating 3-d velocities from cell-wise constant pressure values.
Definition: 3dvelocityadaptive.hh:54
The finite volume base class for the solution of a pressure equation.
Definition: sequential/cellcentered/pressure.hh:49
void initialize()
Initialize pressure model.
Definition: sequential/cellcentered/pressure.hh:213