version 3.8
freeflow/rans/problem.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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_RANS_PROBLEM_HH
13#define DUMUX_RANS_PROBLEM_HH
14
15#include <algorithm>
16
17#include <dune/common/fmatrix.hh>
25#include "model.hh"
26
27namespace Dumux {
28
30template<class TypeTag, TurbulenceModel turbulenceModel>
32
34template<class TypeTag>
36
45template<class TypeTag>
47{
49 using Implementation = GetPropType<TypeTag, Properties::Problem>;
50
52
54 using FVElementGeometry = typename GridGeometry::LocalView;
55 using GridView = typename GridGeometry::GridView;
56 using Element = typename GridView::template Codim<0>::Entity;
57 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
58 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
60 using PrimaryVariables = typename VolumeVariables::PrimaryVariables;
61 using CellCenterPrimaryVariables = GetPropType<TypeTag, Properties::CellCenterPrimaryVariables>;
64
65 using GlobalPosition = typename SubControlVolumeFace::GlobalPosition;
66
67 static constexpr auto dim = GridView::dimension;
68 static constexpr int numCorners = SubControlVolumeFace::numCornersPerFace;
69 using DimVector = GlobalPosition;
70 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
71
72 struct WallElementInformation
73 {
74 // store the element indices for all elements with an intersection on the wall
75 unsigned int wallElementIdx;
76 // for each wall element, store the faces normal axis
77 unsigned int wallFaceNormalAxis;
78 // for each wall element, store the location of the face center and each corner.
79 GlobalPosition wallFaceCenter;
80 std::array<GlobalPosition, numCorners> wallFaceCorners;
81 };
82
83public:
84
90 RANSProblemBase(std::shared_ptr<const GridGeometry> gridGeometry, const std::string& paramGroup = "")
92 {
93 if ( !(hasParamInGroup(this->paramGroup(), "RANS.IsFlatWallBounded")))
94 {
95 std::cout << "The parameter \"Rans.IsFlatWallBounded\" is not specified. \n"
96 << " -- Based on the grid and the boundary conditions specified by the user,"
97 << " this parameter is set to be "<< std::boolalpha << isFlatWallBounded() << "\n";
98 }
99
100 // update size and initial values of the global vectors
101 wallDistance_.resize(this->gridGeometry().elementMapper().size(), std::numeric_limits<Scalar>::max());
102 neighborIdx_.resize(this->gridGeometry().elementMapper().size());
103 velocity_.resize(this->gridGeometry().elementMapper().size(), DimVector(0.0));
104 velocityGradients_.resize(this->gridGeometry().elementMapper().size(), DimMatrix(0.0));
105 stressTensorScalarProduct_.resize(this->gridGeometry().elementMapper().size(), 0.0);
106 vorticityTensorScalarProduct_.resize(this->gridGeometry().elementMapper().size(), 0.0);
107 flowDirectionAxis_.resize(this->gridGeometry().elementMapper().size(), fixedFlowDirectionAxis_);
108 storedViscosity_.resize(this->gridGeometry().elementMapper().size(), 0.0);
109 storedDensity_.resize(this->gridGeometry().elementMapper().size(), 0.0);
110 }
111
116 {
117 std::cout << "Update static wall properties. ";
119
120 checkForWalls_();
121 findWallDistances_();
122 findNeighborIndices_();
123 }
124
130 template<class SolutionVector>
131 void updateDynamicWallProperties(const SolutionVector& curSol)
132 {
133 std::cout << "Update dynamic wall properties." << std::endl;
135 DUNE_THROW(Dune::InvalidStateException,
136 "You have to call updateStaticWallProperties() once before you call updateDynamicWallProperties().");
137
138 calculateCCVelocities_(curSol);
139 calculateCCVelocityGradients_();
140 calculateMaxMinVelocities_();
141 calculateStressTensor_();
142 calculateVorticityTensor_();
143 storeViscosities_(curSol);
144 }
145
153 bool useWallFunction(const Element& element,
154 const SubControlVolumeFace& scvf,
155 const int& eqIdx) const
156 { return false; }
157
161 template<class ElementVolumeVariables, class ElementFaceVariables>
162 FacePrimaryVariables wallFunction(const Element& element,
163 const FVElementGeometry& fvGeometry,
164 const ElementVolumeVariables& elemVolVars,
165 const ElementFaceVariables& elemFaceVars,
166 const SubControlVolumeFace& scvf,
167 const SubControlVolumeFace& lateralBoundaryFace) const
168 { return FacePrimaryVariables(0.0); }
169
173 template<class ElementVolumeVariables, class ElementFaceVariables>
174 CellCenterPrimaryVariables wallFunction(const Element& element,
175 const FVElementGeometry& fvGeometry,
176 const ElementVolumeVariables& elemVolVars,
177 const ElementFaceVariables& elemFaceVars,
178 const SubControlVolumeFace& scvf) const
179 { return CellCenterPrimaryVariables(0.0); }
180
184 bool isFlatWallBounded() const
185 {
186 static const bool hasAlignedWalls = hasAlignedWalls_();
187 return hasAlignedWalls;
188 }
189
193 const Scalar karmanConstant() const
194 { return 0.41; }
195
197 const Scalar betaOmega() const
198 { return 0.0708; }
199
205 {
206 static const Scalar turbulentPrandtlNumber
207 = getParamFromGroup<Scalar>(this->paramGroup(), "RANS.TurbulentPrandtlNumber", 1.0);
209 }
210
216 {
217 static const Scalar turbulentSchmidtNumber
218 = getParamFromGroup<Scalar>(this->paramGroup(), "RANS.TurbulentSchmidtNumber", 1.0);
220 }
221
222 int wallNormalAxis(const int elementIdx) const
223 {
224 if (!isFlatWallBounded())
225 DUNE_THROW(Dune::NotImplemented, "\n Due to grid/geometric concerns, models requiring a wallNormalAxis "
226 << "can only be used for flat wall bounded flows. "
227 << "\n If your geometry is a flat channel, "
228 << "please set the runtime parameter RANS.IsFlatWallBounded to true. \n");
229 return wallNormalAxis_[elementIdx];
230 }
231
232 int flowDirectionAxis(const int elementIdx) const
233 {
234 if (!isFlatWallBounded())
235 DUNE_THROW(Dune::NotImplemented, "\n Due to grid/geometric concerns, models requiring a flowDirectionAxis "
236 << "can only be used for flat wall bounded flows. "
237 << "\n If your geometry is a flat channel, "
238 << "please set the runtime parameter RANS.IsFlatWallBounded to true. \n");
239 return flowDirectionAxis_[elementIdx];
240 }
241
242 unsigned int wallElementIndex(const int elementIdx) const
243 {
244 if (!isFlatWallBounded())
245 DUNE_THROW(Dune::NotImplemented, "\n Due to grid/geometric concerns, models requiring a wallElementIndex "
246 << "can only be used for flat wall bounded flows. "
247 << "\n If your geometry is a flat channel, "
248 << "please set the runtime parameter RANS.IsFlatWallBounded to true. \n");
249 return wallElementIdx_[elementIdx];
250
251 }
252
253 Scalar wallDistance(const int elementIdx) const
254 { return wallDistance_[elementIdx]; }
255
256 GlobalPosition cellCenter(const int elementIdx) const
257 {
258 const auto& element = this->gridGeometry().element(elementIdx);
259 return element.geometry().center();
260 }
261
262 unsigned int neighborIndex(const int elementIdx, const int axisIdx, const int sideIdx) const
263 { return neighborIdx_[elementIdx][axisIdx][sideIdx];}
264
265 DimVector ccVelocityVector(const int elementIdx) const
266 { return velocity_[elementIdx]; }
267
268 Scalar ccVelocity(const int elementIdx, const int axisIdx) const
269 { return velocity_[elementIdx][axisIdx]; }
270
271 DimVector velocityMaximum(const int elementIdx) const
272 { return velocityMaximum_[elementIdx]; }
273
274 DimVector velocityMinimum(const int elementIdx) const
275 { return velocityMinimum_[elementIdx]; }
276
277 DimMatrix velocityGradientTensor(const int elementIdx) const
278 { return velocityGradients_[elementIdx]; }
279
280 Scalar velocityGradient(const int elementIdx, const int i, const int j) const
281 { return velocityGradients_[elementIdx][i][j]; }
282
283 Scalar stressTensorScalarProduct(const int elementIdx) const
284 { return stressTensorScalarProduct_[elementIdx]; }
285
286 Scalar vorticityTensorScalarProduct(const int elementIdx) const
287 { return vorticityTensorScalarProduct_[elementIdx]; }
288
289 Scalar storedViscosity(const int elementIdx) const
290 { return storedViscosity_[elementIdx]; }
291
292 Scalar storedDensity(const int elementIdx) const
293 { return storedDensity_[elementIdx]; }
294
295 Scalar kinematicViscosity(const int elementIdx) const
296 { return storedViscosity(elementIdx) / storedDensity(elementIdx); }
297
299
300private:
301
302 bool hasAlignedWalls_() const
303 {
304 if ( hasParamInGroup(this->paramGroup(), "RANS.IsFlatWallBounded"))
305 {
306 static const bool isFlatWallBounded = getParamFromGroup<bool>(this->paramGroup(), "RANS.IsFlatWallBounded");
307 return isFlatWallBounded;
308 }
309
310 std::vector<int> wallFaceAxis;
311 wallFaceAxis.reserve(this->gridGeometry().numBoundaryScvf());
312
313 const auto gridView = this->gridGeometry().gridView();
314 auto fvGeometry = localView(this->gridGeometry());
315 for (const auto& element : elements(gridView))
316 {
317 fvGeometry.bindElement(element);
318 for (const auto& scvf : scvfs(fvGeometry))
319 if (!scvf.boundary() && asImp_().boundaryTypes(element, scvf).hasWall()) // only search for walls at a global boundary
320 wallFaceAxis.push_back(scvf.directionIndex());
321 }
322
323 // Returns if all wall directions are the same
324 return std::all_of(wallFaceAxis.begin(), wallFaceAxis.end(), [firstDir=wallFaceAxis[0]](auto dir){ return (dir == firstDir);} ) ;
325 }
326
327 void checkForWalls_()
328 {
329 for (const auto& element : elements(this->gridGeometry().gridView()))
330 {
331 auto fvGeometry = localView(this->gridGeometry());
332 fvGeometry.bindElement(element);
333 for (auto&& scvf : scvfs(fvGeometry))
334 if (asImp_().boundaryTypes(element, scvf).hasWall())
335 return;
336 }
337 // If reached, no walls were found using the boundary types has wall function.
338 DUNE_THROW(Dune::InvalidStateException, "No walls are are specified with the setWall() function");
339 }
340
346 void findWallDistances_()
347 {
348 WallDistance wallInformation(this->gridGeometry(), WallDistance<GridGeometry>::atElementCenters,
349 [this] (const FVElementGeometry& fvGeometry, const SubControlVolumeFace& scvf)
350 { return asImp_().boundaryTypes(fvGeometry.element(), scvf).hasWall(); });
351 wallDistance_ = wallInformation.wallDistance();
352 storeWallElementAndDirectionIndex_(wallInformation.wallData());
353 }
354
355 template <class WallData>
356 void storeWallElementAndDirectionIndex_(const WallData& wallData)
357 {
358 // The wall Direction Index is used for flat quadrilateral channel problems only
359 if (!(GridGeometry::discMethod == DiscretizationMethods::staggered))
360 DUNE_THROW(Dune::NotImplemented, "The wall direction Index can only be calculated for quadrilateral structured grids");
361
362 // If isFlatWallBounded, the corresponding wall element is stored for each element
363 if (isFlatWallBounded())
364 {
365 wallNormalAxis_.resize(wallData.size());
366 wallElementIdx_.resize(wallData.size());
367
368 for (const auto& element : elements(this->gridGeometry().gridView()))
369 {
370 unsigned int elementIdx = this->gridGeometry().elementMapper().index(element);
371 wallElementIdx_[elementIdx] = wallData[elementIdx].eIdx;
372 if ( ! (hasParam("RANS.WallNormalAxis")) )
373 {
374 GlobalPosition wallOuterNormal = wallData[elementIdx].scvfOuterNormal;
375 if constexpr (dim == 2) // 2D
376 wallNormalAxis_[elementIdx] = (wallOuterNormal[0] == 1) ? 0 : 1;
377 else // 3D
378 wallNormalAxis_[elementIdx] = (wallOuterNormal[0] == 1) ? 0 : ((wallOuterNormal[1] == 1) ? 1 : 2);
379 }
380 else
381 wallNormalAxis_[elementIdx] = fixedWallNormalAxis_;
382 }
383 }
384 }
385
389 void findNeighborIndices_()
390 {
391 // search for neighbor Idxs
392 for (const auto& element : elements(this->gridGeometry().gridView()))
393 {
394 unsigned int elementIdx = this->gridGeometry().elementMapper().index(element);
395 for (unsigned int axisIdx = 0; axisIdx < dim; ++axisIdx)
396 {
397 neighborIdx_[elementIdx][axisIdx][0] = elementIdx;
398 neighborIdx_[elementIdx][axisIdx][1] = elementIdx;
399 }
400
401 for (const auto& intersection : intersections(this->gridGeometry().gridView(), element))
402 {
403 if (intersection.boundary())
404 continue;
405
406 unsigned int neighborIdx = this->gridGeometry().elementMapper().index(intersection.outside());
407 for (unsigned int axisIdx = 0; axisIdx < dim; ++axisIdx)
408 {
409 if (abs(cellCenter(elementIdx)[axisIdx] - cellCenter(neighborIdx)[axisIdx]) > 1e-8)
410 {
411 if (cellCenter(elementIdx)[axisIdx] > cellCenter(neighborIdx)[axisIdx])
412 neighborIdx_[elementIdx][axisIdx][0] = neighborIdx;
413
414 if (cellCenter(elementIdx)[axisIdx] < cellCenter(neighborIdx)[axisIdx])
415 neighborIdx_[elementIdx][axisIdx][1] = neighborIdx;
416 }
417 }
418 }
419 }
420 }
421
422 template<class SolutionVector>
423 void calculateCCVelocities_(const SolutionVector& curSol)
424 {
425 auto fvGeometry = localView(this->gridGeometry());
426 // calculate cell-center-averaged velocities
427 for (const auto& element : elements(this->gridGeometry().gridView()))
428 {
429 fvGeometry.bindElement(element);
430 unsigned int elementIdx = this->gridGeometry().elementMapper().index(element);
431
432 // calculate velocities
433 DimVector velocityTemp(0.0);
434 for (auto&& scvf : scvfs(fvGeometry))
435 {
436 const int dofIdxFace = scvf.dofIndex();
437 const auto numericalSolutionFace = curSol[GridGeometry::faceIdx()][dofIdxFace][Indices::velocity(scvf.directionIndex())];
438 velocityTemp[scvf.directionIndex()] += numericalSolutionFace;
439 }
440 for (unsigned int axisIdx = 0; axisIdx < dim; ++axisIdx)
441 velocity_[elementIdx][axisIdx] = velocityTemp[axisIdx] * 0.5; // faces are equidistant to cell center
442 }
443 }
444
445
446 void calculateCCVelocityGradients_()
447 {
448 using std::abs;
449
450 // calculate cell-center-averaged velocity gradients, maximum, and minimum values
451 auto fvGeometry = localView(this->gridGeometry());
452 for (const auto& element : elements(this->gridGeometry().gridView()))
453 {
454 const unsigned int elementIdx = this->gridGeometry().elementMapper().index(element);
455
456 for (unsigned int j = 0; j < dim; ++j)
457 {
458 for (unsigned int i = 0; i < dim; ++i)
459 {
460 const unsigned int neighborIndex0 = neighborIndex(elementIdx, j, 0);
461 const unsigned int neighborIndex1 = neighborIndex(elementIdx, j, 1);
462
463 velocityGradients_[elementIdx][i][j]
464 = (ccVelocity(neighborIndex1, i) - ccVelocity(neighborIndex0, i))
465 / (cellCenter(neighborIndex1)[j] - cellCenter(neighborIndex0)[j]);
466
467 if (abs(cellCenter(neighborIndex1)[j] - cellCenter(neighborIndex0)[j]) < 1e-8)
468 velocityGradients_[elementIdx][i][j] = 0.0;
469 }
470 }
471
472 fvGeometry.bindElement(element);
473 for (auto&& scvf : scvfs(fvGeometry))
474 {
475 // adapt calculations for Dirichlet condition
476 unsigned int axisIdx = scvf.directionIndex();
477 if (scvf.boundary())
478 {
479 for (unsigned int velIdx = 0; velIdx < dim; ++velIdx)
480 {
481 if (!asImp_().boundaryTypes(element, scvf).isDirichlet(Indices::velocity(velIdx)))
482 continue;
483
484 Scalar dirichletVelocity = asImp_().dirichlet(element, scvf)[Indices::velocity(velIdx)];
485
486 unsigned int neighborIdx = neighborIndex(elementIdx, axisIdx, 0);
487 if (scvf.center()[axisIdx] < cellCenter(elementIdx)[axisIdx])
488 neighborIdx = neighborIndex(elementIdx, axisIdx, 1);
489
490 velocityGradients_[elementIdx][velIdx][axisIdx]
491 = (ccVelocity(neighborIdx, velIdx) - dirichletVelocity)
492 / (cellCenter(neighborIdx)[axisIdx] - scvf.center()[axisIdx]);
493 }
494 }
495
496 // Calculate the BJS-velocity by accounting for all sub faces.
497 std::vector<int> bjsNumFaces(dim, 0);
498 std::vector<unsigned int> bjsNeighbor(dim, 0);
499 DimVector bjsVelocityAverage(0.0);
500 DimVector normalNormCoordinate(0.0);
501 unsigned int velCompIdx = Indices::velocity(scvf.directionIndex());
502 const int numSubFaces = scvf.pairData().size();
503 for(int localSubFaceIdx = 0; localSubFaceIdx < numSubFaces; ++localSubFaceIdx)
504 {
505 const auto& lateralFace = fvGeometry.scvf(scvf.insideScvIdx(), scvf.pairData()[localSubFaceIdx].localLateralFaceIdx);
506
507 // adapt calculations for Beavers-Joseph-Saffman condition
508 unsigned int lateralAxisIdx = lateralFace.directionIndex();
509 if (lateralFace.boundary() && (asImp_().boundaryTypes(element, lateralFace).isBeaversJoseph(Indices::velocity(velCompIdx))))
510 {
511 unsigned int neighborIdx = neighborIndex(elementIdx, lateralAxisIdx, 0);
512 if (lateralFace.center()[lateralAxisIdx] < cellCenter(elementIdx)[lateralAxisIdx])
513 neighborIdx = neighborIndex(elementIdx, lateralAxisIdx, 1);
514
515 const SubControlVolume& scv = fvGeometry.scv(scvf.insideScvIdx());
516 bjsVelocityAverage[lateralAxisIdx] += ParentType::beaversJosephVelocity(element, scv, scvf, lateralFace, ccVelocity(elementIdx, velCompIdx), 0.0);
517 if (bjsNumFaces[lateralAxisIdx] > 0 && neighborIdx != bjsNeighbor[lateralAxisIdx])
518 DUNE_THROW(Dune::InvalidStateException, "Two different neighborIdx should not occur");
519 bjsNeighbor[lateralAxisIdx] = neighborIdx;
520 normalNormCoordinate[lateralAxisIdx] = lateralFace.center()[lateralAxisIdx];
521 bjsNumFaces[lateralAxisIdx]++;
522 }
523 }
524 for (unsigned axisIdx = 0; axisIdx < dim; ++axisIdx)
525 {
526 if (bjsNumFaces[axisIdx] == 0)
527 continue;
528
529 unsigned int neighborIdx = bjsNeighbor[axisIdx];
530 bjsVelocityAverage[axisIdx] /= bjsNumFaces[axisIdx];
531
532 velocityGradients_[elementIdx][velCompIdx][axisIdx]
533 = (ccVelocity(neighborIdx, velCompIdx) - bjsVelocityAverage[axisIdx])
534 / (cellCenter(neighborIdx)[axisIdx] - normalNormCoordinate[axisIdx]);
535 }
536 }
537 }
538 }
539
540 void calculateMaxMinVelocities_()
541 {
542 using std::abs;
543 if (isFlatWallBounded())
544 {
545 // If the parameter isFlatWallBounded is set to true,
546 // the maximum/minimum velocities are calculated along a profile perpendicular to the corresponding wall face.
547
548 // re-initialize min and max values
549 velocityMaximum_.assign(this->gridGeometry().elementMapper().size(), DimVector(1e-16));
550 velocityMinimum_.assign(this->gridGeometry().elementMapper().size(), DimVector(std::numeric_limits<Scalar>::max()));
551
552 // For each profile perpendicular to the channel wall, find the max and minimum velocities
553 for (const auto& element : elements(this->gridGeometry().gridView()))
554 {
555 const unsigned int elementIdx = this->gridGeometry().elementMapper().index(element);
556 Scalar maxVelocity = 0.0;
557 const unsigned int wallElementIdx = wallElementIndex(elementIdx);
558
559 for (unsigned int axisIdx = 0; axisIdx < dim; ++axisIdx)
560 {
561 if (abs(ccVelocity(elementIdx, axisIdx)) > abs(velocityMaximum_[wallElementIdx][axisIdx]))
562 velocityMaximum_[wallElementIdx][axisIdx] = ccVelocity(elementIdx, axisIdx);
563
564 if (abs(ccVelocity(elementIdx, axisIdx)) < abs(velocityMinimum_[wallElementIdx][axisIdx]))
565 velocityMinimum_[wallElementIdx][axisIdx] = ccVelocity(elementIdx, axisIdx);
566
567 // Set the flow direction axis as the direction of the max velocity
568 if ((hasParam("RANS.FlowDirectionAxis") != 1) && (maxVelocity) < abs(ccVelocity(elementIdx, axisIdx)))
569 {
570 maxVelocity = abs(ccVelocity(elementIdx, axisIdx));
571 flowDirectionAxis_[elementIdx] = axisIdx;
572 }
573 }
574 }
575 }
576 else
577 {
578 // If the parameter isFlatWallBounded is set to false, or not set,
579 // the maximum/minimum velocities are calculated as a global max/min throughout the domain.
580
581 DimVector maxVelocity(0.0);
582 DimVector minVelocity(std::numeric_limits<Scalar>::max());
583 // Find the max and minimum velocities in the full domain
584 for (const auto& element : elements(this->gridGeometry().gridView()))
585 {
586 const unsigned int elementIdx = this->gridGeometry().elementMapper().index(element);
587
588 for (unsigned int axisIdx = 0; axisIdx < dim; ++axisIdx)
589 {
590 if (abs(ccVelocity(elementIdx, axisIdx)) > abs(maxVelocity[axisIdx]))
591 maxVelocity[axisIdx] = ccVelocity(elementIdx, axisIdx);
592
593 if (abs(ccVelocity(elementIdx, axisIdx)) < abs(minVelocity[axisIdx]))
594 minVelocity[axisIdx] = ccVelocity(elementIdx, axisIdx);
595 }
596 }
597 velocityMaximum_.assign(this->gridGeometry().elementMapper().size(), maxVelocity);
598 velocityMinimum_.assign(this->gridGeometry().elementMapper().size(), minVelocity);
599 }
600 }
601
602 void calculateStressTensor_()
603 {
604 for (const auto& element : elements(this->gridGeometry().gridView()))
605 {
606 unsigned int elementIdx = this->gridGeometry().elementMapper().index(element);
607 Dune::FieldMatrix<Scalar, GridView::dimension, GridView::dimension> stressTensor(0.0);
608 for (unsigned int j = 0; j < dim; ++j)
609 {
610 for (unsigned int i = 0; i < dim; ++i)
611 {
612 stressTensor[j][i] = 0.5 * velocityGradient(elementIdx, j, i)
613 + 0.5 * velocityGradient(elementIdx, i, j);
614 }
615 }
616 stressTensorScalarProduct_[elementIdx] = 0.0;
617 for (unsigned int j = 0; j < dim; ++j)
618 {
619 for (unsigned int i = 0; i < dim; ++i)
620 {
621 stressTensorScalarProduct_[elementIdx] += stressTensor[j][i] * stressTensor[j][i];
622 }
623 }
624 }
625 }
626
627 void calculateVorticityTensor_()
628 {
629 for (const auto& element : elements(this->gridGeometry().gridView()))
630 {
631 unsigned int elementIdx = this->gridGeometry().elementMapper().index(element);
632 Dune::FieldMatrix<Scalar, GridView::dimension, GridView::dimension> vorticityTensor(0.0);
633 for (unsigned int j = 0; j < dim; ++j)
634 {
635 for (unsigned int i = 0; i < dim; ++i)
636 {
637 vorticityTensor[j][i] = 0.5 * velocityGradient(elementIdx, j, i)
638 - 0.5 * velocityGradient(elementIdx, i, j);
639 }
640 }
641 vorticityTensorScalarProduct_[elementIdx] = 0.0;
642 for (unsigned int j = 0; j < dim; ++j)
643 {
644 for (unsigned int i = 0; i < dim; ++i)
645 {
646 vorticityTensorScalarProduct_[elementIdx] += vorticityTensor[j][i] * vorticityTensor[j][i];
647 }
648 }
649 }
650 }
651
652 template<class SolutionVector>
653 void storeViscosities_(const SolutionVector& curSol)
654 {
655 // calculate or call all secondary variables
656 auto fvGeometry = localView(this->gridGeometry());
657 for (const auto& element : elements(this->gridGeometry().gridView()))
658 {
659 unsigned int elementIdx = this->gridGeometry().elementMapper().index(element);
660 fvGeometry.bindElement(element);
661 for (auto&& scv : scvs(fvGeometry))
662 {
663 const int dofIdx = scv.dofIndex();
664 // construct a privars object from the cell center solution vector
665 const auto& cellCenterPriVars = curSol[GridGeometry::cellCenterIdx()][dofIdx];
666 PrimaryVariables priVars = makePriVarsFromCellCenterPriVars<PrimaryVariables>(cellCenterPriVars);
667 auto elemSol = elementSolution<typename GridGeometry::LocalView>(std::move(priVars));
668
669 VolumeVariables volVars;
670 volVars.update(elemSol, asImp_(), element, scv);
671 storedDensity_[elementIdx] = volVars.density();
672 storedViscosity_[elementIdx] = volVars.viscosity();
673 }
674 }
675 }
676
677 const int fixedFlowDirectionAxis_ = getParam<int>("RANS.FlowDirectionAxis", 0);
678 const int fixedWallNormalAxis_ = getParam<int>("RANS.WallNormalAxis", 1);
679
680 std::vector<unsigned int> wallNormalAxis_;
681 std::vector<unsigned int> flowDirectionAxis_;
682 std::vector<Scalar> wallDistance_;
683 std::vector<unsigned int> wallElementIdx_;
684 std::vector<std::array<std::array<unsigned int, 2>, dim>> neighborIdx_;
685
686 std::vector<DimVector> velocity_;
687 std::vector<DimVector> velocityMaximum_;
688 std::vector<DimVector> velocityMinimum_;
689 std::vector<DimMatrix> velocityGradients_;
690
691 std::vector<Scalar> stressTensorScalarProduct_;
692 std::vector<Scalar> vorticityTensorScalarProduct_;
693
694 std::vector<Scalar> storedDensity_;
695 std::vector<Scalar> storedViscosity_;
696
698 Implementation &asImp_()
699 { return *static_cast<Implementation *>(this); }
700
702 const Implementation &asImp_() const
703 { return *static_cast<const Implementation *>(this); }
704};
705
706} // end namespace Dumux
707
708#endif
Base class for all finite-volume problems.
Definition: common/fvproblem.hh:43
const std::string & paramGroup() const
The parameter group in which to retrieve runtime parameters.
Definition: common/fvproblem.hh:528
const GridGeometry & gridGeometry() const
The finite volume grid geometry.
Definition: common/fvproblem.hh:524
auto boundaryTypes(const Element &element, const SubControlVolume &scv) const
Specifies which kind of boundary condition should be used for which equation on a given boundary segm...
Definition: common/fvproblem.hh:130
Navier-Stokes staggered problem base class.
Definition: freeflow/navierstokes/staggered/problem.hh:33
const Scalar beaversJosephVelocity(const Element &element, const SubControlVolume &scv, const SubControlVolumeFace &ownScvf, const SubControlVolumeFace &faceOnPorousBoundary, const Scalar velocitySelf, const Scalar tangentialVelocityGradient) const
Returns the slip velocity at a porous boundary based on the Beavers-Joseph(-Saffman) condition.
Definition: freeflow/navierstokes/staggered/problem.hh:189
Reynolds-Averaged Navier-Stokes problem base class.
Definition: freeflow/rans/problem.hh:47
unsigned int wallElementIndex(const int elementIdx) const
Definition: freeflow/rans/problem.hh:242
int flowDirectionAxis(const int elementIdx) const
Definition: freeflow/rans/problem.hh:232
void updateDynamicWallProperties(const SolutionVector &curSol)
Update the dynamic (solution dependent) turbulence parameters.
Definition: freeflow/rans/problem.hh:131
DimMatrix velocityGradientTensor(const int elementIdx) const
Definition: freeflow/rans/problem.hh:277
Scalar vorticityTensorScalarProduct(const int elementIdx) const
Definition: freeflow/rans/problem.hh:286
Scalar ccVelocity(const int elementIdx, const int axisIdx) const
Definition: freeflow/rans/problem.hh:268
void updateStaticWallProperties()
Update the static (solution independent) relations to the walls and neighbors.
Definition: freeflow/rans/problem.hh:115
GlobalPosition cellCenter(const int elementIdx) const
Definition: freeflow/rans/problem.hh:256
FacePrimaryVariables wallFunction(const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const SubControlVolumeFace &scvf, const SubControlVolumeFace &lateralBoundaryFace) const
Returns an additional wall function momentum flux.
Definition: freeflow/rans/problem.hh:162
Scalar velocityGradient(const int elementIdx, const int i, const int j) const
Definition: freeflow/rans/problem.hh:280
DimVector velocityMaximum(const int elementIdx) const
Definition: freeflow/rans/problem.hh:271
const Scalar betaOmega() const
Returns the constant.
Definition: freeflow/rans/problem.hh:197
DimVector ccVelocityVector(const int elementIdx) const
Definition: freeflow/rans/problem.hh:265
Scalar kinematicViscosity(const int elementIdx) const
Definition: freeflow/rans/problem.hh:295
Scalar turbulentSchmidtNumber() const
Return the turbulent Schmidt number which is used to convert the eddy viscosity to an eddy diffusivi...
Definition: freeflow/rans/problem.hh:215
CellCenterPrimaryVariables wallFunction(const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const SubControlVolumeFace &scvf) const
Returns an additional wall function flux for cell-centered quantities.
Definition: freeflow/rans/problem.hh:174
DimVector velocityMinimum(const int elementIdx) const
Definition: freeflow/rans/problem.hh:274
Scalar wallDistance(const int elementIdx) const
Definition: freeflow/rans/problem.hh:253
Scalar storedDensity(const int elementIdx) const
Definition: freeflow/rans/problem.hh:292
Scalar turbulentPrandtlNumber() const
Return the turbulent Prandtl number which is used to convert the eddy viscosity to an eddy thermal c...
Definition: freeflow/rans/problem.hh:204
int wallNormalAxis(const int elementIdx) const
Definition: freeflow/rans/problem.hh:222
unsigned int neighborIndex(const int elementIdx, const int axisIdx, const int sideIdx) const
Definition: freeflow/rans/problem.hh:262
RANSProblemBase(std::shared_ptr< const GridGeometry > gridGeometry, const std::string &paramGroup="")
The constructor.
Definition: freeflow/rans/problem.hh:90
bool isFlatWallBounded() const
Returns whether a given sub control volume face is on a wall.
Definition: freeflow/rans/problem.hh:184
Scalar stressTensorScalarProduct(const int elementIdx) const
Definition: freeflow/rans/problem.hh:283
bool calledUpdateStaticWallProperties
Definition: freeflow/rans/problem.hh:298
bool useWallFunction(const Element &element, const SubControlVolumeFace &scvf, const int &eqIdx) const
Returns whether a wall function should be used at a given face.
Definition: freeflow/rans/problem.hh:153
const Scalar karmanConstant() const
Returns the Karman constant.
Definition: freeflow/rans/problem.hh:193
Scalar storedViscosity(const int elementIdx) const
Definition: freeflow/rans/problem.hh:289
forward declare
Definition: freeflow/rans/problem.hh:31
static constexpr auto atElementCenters
Definition: walldistance.hh:98
Defines all properties used in Dumux.
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:26
bool hasParamInGroup(const std::string &paramGroup, const std::string &param)
Check whether a key exists in the parameter tree with a model group prefix.
Definition: parameters.hh:165
bool hasParam(const std::string &param)
Check whether a key exists in the parameter tree.
Definition: parameters.hh:157
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
Free function to get the local view of a grid cache object.
The available discretization methods in Dumux.
constexpr Staggered staggered
Definition: method.hh:149
Definition: adapt.hh:17
Adaption of the fully implicit scheme to the tracer transport model.
The local element solution class for staggered methods.
Base class for all staggered fv problems.