3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
2p/sequential/diffusion/cellcentered/velocity.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*****************************************************************************
4 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_FVVELOCITY2P_HH
25#define DUMUX_FVVELOCITY2P_HH
26
27#include <dune/common/float_cmp.hh>
28#include <dune/grid/common/gridenums.hh>
30
32
33namespace Dumux {
34
61template<class TypeTag>
63{
67
69
71
74
77 using PrimaryVariables = typename SolutionTypes::PrimaryVariables;
79
80 using Element = typename GridView::Traits::template Codim<0>::Entity;
81 using Intersection = typename GridView::Intersection;
82
83 using Geometry = typename Element::Geometry;
84 using JacobianTransposed = typename Geometry::JacobianTransposed;
85
86 enum
87 {
88 dim = GridView::dimension, dimWorld = GridView::dimensionworld
89 };
90
91 enum
92 {
93 pw = Indices::pressureW,
94 pn = Indices::pressureNw,
95 pGlobal = Indices::pressureGlobal,
96 vw = Indices::velocityW,
97 vn = Indices::velocityNw,
98 vt = Indices::velocityTotal,
99 sw = Indices::saturationW,
100 sn = Indices::saturationNw,
101 pressureIdx = Indices::pressureIdx,
102 saturationIdx = Indices::saturationIdx,
103 eqIdxPress = Indices::pressureEqIdx,
104 eqIdxSat = Indices::satEqIdx
105 };
106 enum
107 {
108 wPhaseIdx = Indices::wPhaseIdx, nPhaseIdx = Indices::nPhaseIdx, numPhases = getPropValue<TypeTag, Properties::NumPhases>()
109 };
110
111 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
112 using GravityVector = Dune::FieldVector<Scalar, dimWorld>;
113 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
114
115public:
120 FVVelocity2P(Problem& problem) :
121 problem_(problem), gravity_(problem.gravity())
122 {
123 if (getPropValue<TypeTag, Properties::EnableCompressibility>() && velocityType_ == vt)
124 {
125 DUNE_THROW(Dune::NotImplemented,
126 "Total velocity - global pressure - model cannot be used with compressible fluids!");
127 }
128 if (velocityType_ != vw && velocityType_ != vn && velocityType_ != vt)
129 {
130 DUNE_THROW(Dune::NotImplemented, "Velocity type not supported!");
131 }
132
133 density_[wPhaseIdx] = 0.;
134 density_[nPhaseIdx] = 0.;
135 viscosity_[wPhaseIdx] = 0.;
136 viscosity_[nPhaseIdx] = 0.;
137
138 vtkOutputLevel_ = getParam<int>("Vtk.OutputLevel");
139 }
140
143 {
144 if (!compressibility_)
145 {
146 const auto element = *problem_.gridView().template begin<0> ();
147 FluidState fluidState;
148 fluidState.setPressure(wPhaseIdx, problem_.referencePressure(element));
149 fluidState.setPressure(nPhaseIdx, problem_.referencePressure(element));
150 fluidState.setTemperature(problem_.temperature(element));
151 fluidState.setSaturation(wPhaseIdx, 1.);
152 fluidState.setSaturation(nPhaseIdx, 0.);
153 density_[wPhaseIdx] = FluidSystem::density(fluidState, wPhaseIdx);
154 density_[nPhaseIdx] = FluidSystem::density(fluidState, nPhaseIdx);
155 viscosity_[wPhaseIdx] = FluidSystem::viscosity(fluidState, wPhaseIdx);
156 viscosity_[nPhaseIdx] = FluidSystem::viscosity(fluidState, nPhaseIdx);
157 }
158 }
159
161 void calculateVelocity(const Intersection&, CellData&);
162
164 void calculateVelocityOnBoundary(const Intersection&, CellData&);
165
172 {
173 return true;
174 }
175
184 template<class MultiWriter>
185 void addOutputVtkFields(MultiWriter &writer)
186 {
187 if (vtkOutputLevel_ > 0)
188 {
189 Dune::BlockVector < Dune::FieldVector<Scalar, dim> > &velocity = *(writer.template allocateManagedBuffer<Scalar,
190 dim>(problem_.gridView().size(0)));
191 Dune::BlockVector < Dune::FieldVector<Scalar, dim> > &velocitySecondPhase =
192 *(writer.template allocateManagedBuffer<Scalar, dim>(problem_.gridView().size(0)));
193
194 // compute update vector
195 for (const auto& element : elements(problem_.gridView()))
196 {
197 // cell index
198 int eIdxGlobal = problem_.variables().index(element);
199
200 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
201
202 const typename Element::Geometry& geometry = element.geometry();
203
204 // get corresponding reference element
205 const auto refElement = referenceElement(geometry);
206
207 const int numberOfFaces=refElement.size(1);
208 std::vector<Scalar> fluxW(numberOfFaces,0);
209 std::vector<Scalar> fluxNw(numberOfFaces,0);
210
211 // run through all intersections with neighbors and boundary
212 for (const auto& intersection : intersections(problem_.gridView(), element))
213 {
214 int isIndex = intersection.indexInInside();
215
216 fluxW[isIndex] += intersection.geometry().volume()
217 * (intersection.centerUnitOuterNormal() * cellData.fluxData().velocity(wPhaseIdx, isIndex));
218 fluxNw[isIndex] += intersection.geometry().volume()
219 * (intersection.centerUnitOuterNormal() * cellData.fluxData().velocity(nPhaseIdx, isIndex));
220 }
221
222 // calculate velocity on reference element as the Raviart-Thomas-0
223 // interpolant of the fluxes
224 Dune::FieldVector<Scalar, dim> refVelocity;
225 // simplices
226 if (refElement.type().isSimplex()) {
227 for (int dimIdx = 0; dimIdx < dim; dimIdx++)
228 {
229 refVelocity[dimIdx] = -fluxW[dim - 1 - dimIdx];
230 for (int fIdx = 0; fIdx < dim + 1; fIdx++)
231 {
232 refVelocity[dimIdx] += fluxW[fIdx]/(dim + 1);
233 }
234 }
235 }
236 // cubes
237 else if (refElement.type().isCube()){
238 for (int i = 0; i < dim; i++)
239 refVelocity[i] = 0.5 * (fluxW[2*i + 1] - fluxW[2*i]);
240 }
241 // 3D prism and pyramids
242 else {
243 DUNE_THROW(Dune::NotImplemented, "velocity output for prism/pyramid not implemented");
244 }
245
246 const Dune::FieldVector<Scalar, dim> localPos =
247 refElement.position(0, 0);
248
249 // get the transposed Jacobian of the element mapping
250 const JacobianTransposed jacobianT =
251 geometry.jacobianTransposed(localPos);
252
253 // calculate the element velocity by the Piola transformation
254 Dune::FieldVector < Scalar, dim > elementVelocity(0);
255 jacobianT.umtv(refVelocity, elementVelocity);
256 elementVelocity /= geometry.integrationElement(localPos);
257
258 velocity[eIdxGlobal] = elementVelocity;
259
260 // calculate velocity on reference element as the Raviart-Thomas-0
261 // interpolant of the fluxes
262 // simplices
263 if (refElement.type().isSimplex()) {
264 for (int dimIdx = 0; dimIdx < dim; dimIdx++)
265 {
266 refVelocity[dimIdx] = -fluxNw[dim - 1 - dimIdx];
267 for (int fIdx = 0; fIdx < dim + 1; fIdx++)
268 {
269 refVelocity[dimIdx] += fluxNw[fIdx]/(dim + 1);
270 }
271 }
272 }
273 // cubes
274 else if (refElement.type().isCube()){
275 for (int i = 0; i < dim; i++)
276 refVelocity[i] = 0.5 * (fluxNw[2*i + 1] - fluxNw[2*i]);
277 }
278 // 3D prism and pyramids
279 else {
280 DUNE_THROW(Dune::NotImplemented, "velocity output for prism/pyramid not implemented");
281 }
282
283 // calculate the element velocity by the Piola transformation
284 elementVelocity = 0;
285 jacobianT.umtv(refVelocity, elementVelocity);
286 elementVelocity /= geometry.integrationElement(localPos);
287
288 velocitySecondPhase[eIdxGlobal] = elementVelocity;
289 }
290
291 // switch velocities
292 if (velocityType_ == vt)
293 {
294 writer.attachCellData(velocity, "total velocity", dim);
295 }
296 else
297 {
298 writer.attachCellData(velocity, "wetting-velocity", dim);
299 writer.attachCellData(velocitySecondPhase, "nonwetting-velocity", dim);
300 }
301 }
302
303 return;
304 }
305
306private:
307 Problem& problem_;
308 const GravityVector& gravity_;
309 Scalar density_[numPhases];
310 Scalar viscosity_[numPhases];
311
312 int vtkOutputLevel_;
313
315 static const int velocityType_ = getPropValue<TypeTag, Properties::VelocityFormulation>();
316 static const bool compressibility_ = getPropValue<TypeTag, Properties::EnableCompressibility>();
318 static const int pressureType_ = getPropValue<TypeTag, Properties::PressureFormulation>();
320 static const int saturationType_ = getPropValue<TypeTag, Properties::SaturationFormulation>();
321};
322
331template<class TypeTag>
332void FVVelocity2P<TypeTag>::calculateVelocity(const Intersection& intersection, CellData& cellData)
333{
334 auto elementI = intersection.inside();
335 auto elementJ = intersection.outside();
336
337 int eIdxGlobalJ = problem_.variables().index(elementJ);
338
339 CellData& cellDataJ = problem_.variables().cellData(eIdxGlobalJ);
340
341 // get global coordinates of cell centers
342 const GlobalPosition& globalPosI = (elementI).geometry().center();
343 const GlobalPosition& globalPosJ = (elementJ).geometry().center();
344
345 // get mobilities and fractional flow factors
346 Scalar lambdaWI = cellData.mobility(wPhaseIdx);
347 Scalar lambdaNwI = cellData.mobility(nPhaseIdx);
348 Scalar lambdaWJ = cellDataJ.mobility(wPhaseIdx);
349 Scalar lambdaNwJ = cellDataJ.mobility(nPhaseIdx);
350
351 // get capillary pressure
352 Scalar pcI = cellData.capillaryPressure();
353 Scalar pcJ = cellDataJ.capillaryPressure();
354
355 // get face index
356 int isIndexI = intersection.indexInInside();
357 int isIndexJ = intersection.indexInOutside();
358
359 // get face normal
360 const Dune::FieldVector<Scalar, dim>& unitOuterNormal = intersection.centerUnitOuterNormal();
361
362 // distance vector between barycenters
363 GlobalPosition distVec = globalPosJ - globalPosI;
364
365 // compute distance between cell centers
366 Scalar dist = distVec.two_norm();
367
368 // compute vectorized permeabilities
369 DimMatrix meanPermeability(0);
370
371 problem_.spatialParams().meanK(meanPermeability, problem_.spatialParams().intrinsicPermeability(elementI),
372 problem_.spatialParams().intrinsicPermeability(elementJ));
373
374 Dune::FieldVector<Scalar, dim> permeability(0);
375 meanPermeability.mv(unitOuterNormal, permeability);
376
377 // calculate potential gradients
378 Scalar potentialDiffW = cellData.potential(wPhaseIdx) - cellDataJ.potential(wPhaseIdx);
379 Scalar potentialDiffNw = cellData.potential(nPhaseIdx) - cellDataJ.potential(nPhaseIdx);
380
381 if (compressibility_)
382 {
383 density_[wPhaseIdx] = (potentialDiffW > 0.) ? cellData.density(wPhaseIdx) : cellDataJ.density(wPhaseIdx);
384 density_[nPhaseIdx] = (potentialDiffNw > 0.) ? cellData.density(nPhaseIdx) : cellDataJ.density(nPhaseIdx);
385
386 density_[wPhaseIdx] =
387 (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffW, 0.0, 1.0e-30)) ? 0.5 * (cellData.density(wPhaseIdx) + cellDataJ.density(wPhaseIdx)) :
388 density_[wPhaseIdx];
389 density_[nPhaseIdx] =
390 (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffNw, 0.0, 1.0e-30)) ? 0.5 * (cellData.density(nPhaseIdx) + cellDataJ.density(nPhaseIdx)) :
391 density_[nPhaseIdx];
392
393 potentialDiffW = (cellData.pressure(wPhaseIdx) - cellDataJ.pressure(wPhaseIdx));
394 potentialDiffNw = (cellData.pressure(nPhaseIdx) - cellDataJ.pressure(nPhaseIdx));
395
396 potentialDiffW += density_[wPhaseIdx] * (distVec * gravity_); //delta z/delta x in unitOuterNormal[z]
397 potentialDiffNw += density_[nPhaseIdx] * (distVec * gravity_);
398 }
399
400 // store potentials for further calculations (velocity, saturation, ...)
401 cellData.fluxData().setUpwindPotential(wPhaseIdx, isIndexI, potentialDiffW);
402 cellData.fluxData().setUpwindPotential(nPhaseIdx, isIndexI, potentialDiffNw);
403
404 cellDataJ.fluxData().setUpwindPotential(wPhaseIdx, isIndexJ, -potentialDiffW);
405 cellDataJ.fluxData().setUpwindPotential(nPhaseIdx, isIndexJ, -potentialDiffNw);
406
407 // do the upwinding of the mobility depending on the phase potentials
408 Scalar lambdaW = (potentialDiffW > 0.) ? lambdaWI : lambdaWJ;
409 lambdaW = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffW, 0.0, 1.0e-30)) ? 0.5 * (lambdaWI + lambdaWJ) : lambdaW;
410 Scalar lambdaNw = (potentialDiffNw > 0.) ? lambdaNwI : lambdaNwJ;
411 lambdaNw = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffNw, 0.0, 1.0e-30)) ? 0.5 * (lambdaNwI + lambdaNwJ) : lambdaNw;
412
413 if (compressibility_)
414 {
415 density_[wPhaseIdx] = (potentialDiffW > 0.) ? cellData.density(wPhaseIdx) : cellDataJ.density(wPhaseIdx);
416 density_[nPhaseIdx] = (potentialDiffNw > 0.) ? cellData.density(nPhaseIdx) : cellDataJ.density(nPhaseIdx);
417
418 density_[wPhaseIdx] =
419 (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffW, 0.0, 1.0e-30)) ? 0.5 * (cellData.density(wPhaseIdx) + cellDataJ.density(wPhaseIdx)) :
420 density_[wPhaseIdx];
421 density_[nPhaseIdx] =
422 (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffNw, 0.0, 1.0e-30)) ? 0.5 * (cellData.density(nPhaseIdx) + cellDataJ.density(nPhaseIdx)) :
423 density_[nPhaseIdx];
424 }
425
426 Scalar scalarPerm = permeability.two_norm();
427
428 // calculate the gravity term
429 Dune::FieldVector<Scalar, dimWorld> velocityW(unitOuterNormal);
430 Dune::FieldVector<Scalar, dimWorld> velocityNw(unitOuterNormal);
431
432 // calculate unit distVec
433 distVec /= dist;
434 Scalar areaScaling = (unitOuterNormal * distVec);
435 // this treatment of g allows to account for gravity flux through faces where the face normal
436 // has no z component (e.g. parallelepiped grids)
437 Scalar gravityTermW = (gravity_ * distVec) * density_[wPhaseIdx] * areaScaling;
438 Scalar gravityTermNw = (gravity_ * distVec) * density_[nPhaseIdx] * areaScaling;
439
440 // calculate velocity depending on the pressure used -> use pc = pn - pw
441 switch (pressureType_)
442 {
443 case pw:
444 {
445 velocityW *= lambdaW * scalarPerm
446 * ((cellData.pressure(wPhaseIdx) - cellDataJ.pressure(wPhaseIdx)) / dist + gravityTermW);
447 velocityNw *= lambdaNw * scalarPerm
448 * ((cellData.pressure(wPhaseIdx) - cellDataJ.pressure(wPhaseIdx)) / dist + gravityTermNw)
449 + 0.5 * (lambdaNwI + lambdaNwJ) * scalarPerm * (pcI - pcJ) / dist;
450 break;
451 }
452 case pn:
453 {
454 velocityW *= lambdaW * scalarPerm
455 * ((cellData.pressure(nPhaseIdx) - cellDataJ.pressure(nPhaseIdx)) / dist + gravityTermW)
456 - 0.5 * (lambdaWI + lambdaWJ) * scalarPerm * (pcI - pcJ) / dist;
457 velocityNw *= lambdaNw * scalarPerm
458 * ((cellData.pressure(nPhaseIdx) - cellDataJ.pressure(nPhaseIdx)) / dist + gravityTermNw);
459 break;
460 }
461 case pGlobal:
462 {
463 velocityW *= (lambdaW + lambdaNw) * scalarPerm * (cellData.globalPressure() - cellDataJ.globalPressure()) / dist
464 + scalarPerm * (lambdaW * gravityTermW + lambdaNw * gravityTermNw);
465 velocityNw = 0;
466 break;
467 }
468 }
469
470 // store velocities
471 cellData.fluxData().setVelocity(wPhaseIdx, isIndexI, velocityW);
472 cellData.fluxData().setVelocity(nPhaseIdx, isIndexI, velocityNw);
473 cellData.fluxData().setVelocityMarker(isIndexI);
474
475 cellDataJ.fluxData().setVelocity(wPhaseIdx, isIndexJ, velocityW);
476 cellDataJ.fluxData().setVelocity(nPhaseIdx, isIndexJ, velocityNw);
477 cellDataJ.fluxData().setVelocityMarker(isIndexJ);
478
479// printvector(std::cout, cellData.fluxData().velocity(), "velocity", "row", 4, 1, 3);
480 return;
481}
482
491template<class TypeTag>
492void FVVelocity2P<TypeTag>::calculateVelocityOnBoundary(const Intersection& intersection, CellData& cellData)
493{
494 auto element = intersection.inside();
495
496 // get face index
497 int isIndex = intersection.indexInInside();
498
499 // get face normal
500 const Dune::FieldVector<Scalar, dim>& unitOuterNormal = intersection.centerUnitOuterNormal();
501
502 BoundaryTypes bcType;
503 //get boundary type
504 problem_.boundaryTypes(bcType, intersection);
505 PrimaryVariables boundValues(0.0);
506
507 if (bcType.isDirichlet(eqIdxPress))
508 {
509 problem_.dirichlet(boundValues, intersection);
510
511 // get global coordinates of cell centers
512 const GlobalPosition& globalPosI = element.geometry().center();
513
514 // center of face in global coordinates
515 const GlobalPosition& globalPosJ = intersection.geometry().center();
516
517 // get mobilities and fractional flow factors
518 Scalar lambdaWI = cellData.mobility(wPhaseIdx);
519 Scalar lambdaNwI = cellData.mobility(nPhaseIdx);
520 Scalar fractionalWI = cellData.fracFlowFunc(wPhaseIdx);
521 Scalar fractionalNwI = cellData.fracFlowFunc(nPhaseIdx);
522
523 // get capillary pressure
524 Scalar pcI = cellData.capillaryPressure();
525
526 // distance vector between barycenters
527 GlobalPosition distVec = globalPosJ - globalPosI;
528
529 // compute distance between cell centers
530 Scalar dist = distVec.two_norm();
531
532 // permeability vector at boundary
533 // compute vectorized permeabilities
534 DimMatrix meanPermeability(0);
535
536 problem_.spatialParams().meanK(meanPermeability, problem_.spatialParams().intrinsicPermeability(element));
537
538 Dune::FieldVector<Scalar, dim> permeability(0);
539 meanPermeability.mv(unitOuterNormal, permeability);
540
541 // determine saturation at the boundary -> if no saturation is known directly at the boundary use the cell saturation
542 Scalar satW = 0;
543 Scalar satNw = 0;
544 if (bcType.isDirichlet(eqIdxSat))
545 {
546 switch (saturationType_)
547 {
548 case sw:
549 {
550 satW = boundValues[saturationIdx];
551 satNw = 1 - boundValues[saturationIdx];
552 break;
553 }
554 case sn:
555 {
556 satW = 1 - boundValues[saturationIdx];
557 satNw = boundValues[saturationIdx];
558 break;
559 }
560 }
561 }
562 else
563 {
564 satW = cellData.saturation(wPhaseIdx);
565 satNw = cellData.saturation(nPhaseIdx);
566 }
567
568 const Scalar pressBound = boundValues[pressureIdx];
569
570 // old material law interface is deprecated: Replace this by
571 // const auto& fluidMatrixInteraction = spatialParams.fluidMatrixInteractionAtPos(element.geometry().center());
572 // after the release of 3.3, when the deprecated interface is no longer supported
573 const auto fluidMatrixInteraction = Deprecated::makePcKrSw(Scalar{}, problem_.spatialParams(), element);
574
575 const Scalar pcBound = fluidMatrixInteraction.pc(satW);
576
577 // determine phase pressures from primary pressure variable
578 Scalar pressWBound = 0;
579 Scalar pressNwBound = 0;
580 if (pressureType_ == pw)
581 {
582 pressWBound = pressBound;
583 pressNwBound = pressBound + pcBound;
584 }
585 else if (pressureType_ == pn)
586 {
587 pressWBound = pressBound - pcBound;
588 pressNwBound = pressBound;
589 }
590
591 // get temperature at current position
592 const Scalar temperature = problem_.temperature(element);
593
594 Scalar densityWBound = density_[wPhaseIdx];
595 Scalar densityNwBound = density_[nPhaseIdx];
596 Scalar viscosityWBound = viscosity_[wPhaseIdx];
597 Scalar viscosityNwBound = viscosity_[nPhaseIdx];
598
599 if (compressibility_)
600 {
601 FluidState fluidState;
602 fluidState.setSaturation(wPhaseIdx, satW);
603 fluidState.setSaturation(nPhaseIdx, satNw);
604 fluidState.setTemperature(temperature);
605 fluidState.setPressure(wPhaseIdx, pressWBound);
606 fluidState.setPressure(nPhaseIdx, pressNwBound);
607
608 densityWBound = FluidSystem::density(fluidState, wPhaseIdx);
609 densityNwBound = FluidSystem::density(fluidState, nPhaseIdx);
610 viscosityWBound = FluidSystem::viscosity(fluidState, wPhaseIdx) / densityWBound;
611 viscosityNwBound = FluidSystem::viscosity(fluidState, nPhaseIdx) / densityNwBound;
612 }
613
614 Scalar lambdaWBound = fluidMatrixInteraction.krw(satW)
615 / viscosityWBound;
616 Scalar lambdaNwBound = fluidMatrixInteraction.krn(satW)
617 / viscosityNwBound;
618
619 Scalar potentialDiffW = cellData.fluxData().upwindPotential(wPhaseIdx, isIndex);
620 Scalar potentialDiffNw = cellData.fluxData().upwindPotential(nPhaseIdx, isIndex);
621
622 if (compressibility_)
623 {
624 density_[wPhaseIdx] = (potentialDiffW > 0.) ? cellData.density(wPhaseIdx) : densityWBound;
625 density_[wPhaseIdx] =
626 (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffW, 0.0, 1.0e-30)) ? 0.5 * (cellData.density(wPhaseIdx) + densityWBound) : density_[wPhaseIdx];
627 density_[nPhaseIdx] = (potentialDiffNw > 0.) ? cellData.density(nPhaseIdx) : densityNwBound;
628 density_[nPhaseIdx] =
629 (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffNw, 0.0, 1.0e-30)) ? 0.5 * (cellData.density(nPhaseIdx) + densityNwBound) : density_[nPhaseIdx];
630 }
631
632 // calculate potential gradient
633 if (pressureType_ == pGlobal)
634 {
635 potentialDiffW = (cellData.globalPressure() - pressBound - fractionalNwI * (pcI - pcBound));
636 potentialDiffNw = (cellData.globalPressure() - pressBound + fractionalWI * (pcI - pcBound));
637 }
638 else
639 {
640 potentialDiffW = (cellData.pressure(wPhaseIdx) - pressWBound);
641 potentialDiffNw = (cellData.pressure(nPhaseIdx) - pressNwBound);
642 }
643
644 potentialDiffW += density_[wPhaseIdx] * (distVec * gravity_);
645 potentialDiffNw += density_[nPhaseIdx] * (distVec * gravity_);
646
647 // store potential gradients for further calculations
648 cellData.fluxData().setUpwindPotential(wPhaseIdx, isIndex, potentialDiffW);
649 cellData.fluxData().setUpwindPotential(nPhaseIdx, isIndex, potentialDiffNw);
650
651 // do the upwinding of the mobility depending on the phase potentials
652 Scalar lambdaW = (potentialDiffW > 0.) ? lambdaWI : lambdaWBound;
653 lambdaW = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffW, 0.0, 1.0e-30)) ? 0.5 * (lambdaWI + lambdaWBound) : lambdaW;
654 Scalar lambdaNw = (potentialDiffNw > 0.) ? lambdaNwI : lambdaNwBound;
655 lambdaNw = (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffNw, 0.0, 1.0e-30)) ? 0.5 * (lambdaNwI + lambdaNwBound) : lambdaNw;
656
657 if (compressibility_)
658 {
659 density_[wPhaseIdx] = (potentialDiffW > 0.) ? cellData.density(wPhaseIdx) : densityWBound;
660 density_[wPhaseIdx] =
661 (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffW, 0.0, 1.0e-30)) ? 0.5 * (cellData.density(wPhaseIdx) + densityWBound) : density_[wPhaseIdx];
662 density_[nPhaseIdx] = (potentialDiffNw > 0.) ? cellData.density(nPhaseIdx) : densityNwBound;
663 density_[nPhaseIdx] =
664 (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potentialDiffNw, 0.0, 1.0e-30)) ? 0.5 * (cellData.density(nPhaseIdx) + densityNwBound) : density_[nPhaseIdx];
665 }
666
667 const Scalar scalarPerm = permeability.two_norm();
668
669 // calculate the gravity term
670 Dune::FieldVector<Scalar, dimWorld> velocityW(unitOuterNormal);
671 Dune::FieldVector<Scalar, dimWorld> velocityNw(unitOuterNormal);
672
673 // calculate unit distVec
674 distVec /= dist;
675 Scalar areaScaling = (unitOuterNormal * distVec);
676 // this treatment of g allows to account for gravity flux through faces where the face normal
677 // has no z component (e.g. parallelepiped grids)
678 Scalar gravityTermW = (gravity_ * distVec) * density_[wPhaseIdx] * areaScaling;
679 Scalar gravityTermNw = (gravity_ * distVec) * density_[nPhaseIdx] * areaScaling;
680
681 // calculate velocity depending on the pressure used -> use pc = pn - pw
682 switch (pressureType_)
683 {
684 case pw:
685 {
686 velocityW *= lambdaW * scalarPerm * ((cellData.pressure(wPhaseIdx) - pressBound) / dist + gravityTermW);
687 velocityNw *= lambdaNw * scalarPerm * ((cellData.pressure(wPhaseIdx) - pressBound) / dist + gravityTermNw)
688 + 0.5 * (lambdaNwI + lambdaNwBound) * scalarPerm * (pcI - pcBound) / dist;
689 break;
690 }
691 case pn:
692 {
693 velocityW *= lambdaW * scalarPerm * ((cellData.pressure(nPhaseIdx) - pressBound) / dist + gravityTermW)
694 - 0.5 * (lambdaWI + lambdaWBound) * scalarPerm * (pcI - pcBound) / dist;
695 velocityNw *= lambdaNw * scalarPerm * ((cellData.pressure(nPhaseIdx) - pressBound) / dist + gravityTermNw);
696 break;
697 }
698 case pGlobal:
699 {
700 velocityW *= (lambdaW + lambdaNw) * scalarPerm * (cellData.globalPressure() - pressBound) / dist
701 + scalarPerm * (lambdaW * gravityTermW + lambdaNw * gravityTermNw);
702 velocityNw = 0;
703 break;
704 }
705 }
706
707 // store velocities
708 cellData.fluxData().setVelocity(wPhaseIdx, isIndex, velocityW);
709 cellData.fluxData().setVelocity(nPhaseIdx, isIndex, velocityNw);
710 cellData.fluxData().setVelocityMarker(isIndex);
711
712 } // end Dirichlet boundary
713
714 else if (bcType.isNeumann(eqIdxPress))
715 {
716 problem_.neumann(boundValues, intersection);
717
718 Dune::FieldVector<Scalar, dimWorld> velocityW(unitOuterNormal);
719 Dune::FieldVector<Scalar, dimWorld> velocityNw(unitOuterNormal);
720
721 velocityW *= boundValues[wPhaseIdx];
722 velocityNw *= boundValues[nPhaseIdx];
723
724 if (!compressibility_)
725 {
726 velocityW /= density_[wPhaseIdx];
727 velocityNw /= density_[nPhaseIdx];
728 }
729
730 // store potential gradients for further calculations
731 cellData.fluxData().setUpwindPotential(wPhaseIdx, isIndex, boundValues[wPhaseIdx]);
732 cellData.fluxData().setUpwindPotential(nPhaseIdx, isIndex, boundValues[nPhaseIdx]);
733
734 cellData.fluxData().setVelocity(wPhaseIdx, isIndex, velocityW);
735 cellData.fluxData().setVelocity(nPhaseIdx, isIndex, velocityNw);
736 cellData.fluxData().setVelocityMarker(isIndex);
737 } // end Neumann boundary
738 else
739 {
740 DUNE_THROW(Dune::NotImplemented, "No valid boundary condition type defined for pressure equation!");
741 }
742
743// printvector(std::cout, cellData.fluxData().velocity(), "velocity", "row", 4, 1, 3);
744 return;
745}
746} // end namespace Dumux
747#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 temperature() noexcept
I/O name of temperature for equilibrium models.
Definition: name.hh:51
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
Determines the velocity from a finite volume solution of the pressure equation of a sequential model ...
Definition: 2p/sequential/diffusion/cellcentered/velocity.hh:63
void addOutputVtkFields(MultiWriter &writer)
Adds velocity output to the output file.
Definition: 2p/sequential/diffusion/cellcentered/velocity.hh:185
void calculateVelocity(const Intersection &, CellData &)
Calculates the velocity at a cell-cell interface.
Definition: 2p/sequential/diffusion/cellcentered/velocity.hh:332
bool calculateVelocityInTransport()
Indicates if velocity is reconstructed in the pressure step or in the transport step.
Definition: 2p/sequential/diffusion/cellcentered/velocity.hh:171
FVVelocity2P(Problem &problem)
Constructs a FVVelocity2P object.
Definition: 2p/sequential/diffusion/cellcentered/velocity.hh:120
void initialize()
For initialization.
Definition: 2p/sequential/diffusion/cellcentered/velocity.hh:142
void calculateVelocityOnBoundary(const Intersection &, CellData &)
Calculates the velocity at a boundary.
Definition: 2p/sequential/diffusion/cellcentered/velocity.hh:492
Specifies the properties for immiscible 2p diffusion/pressure models.