3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
test/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 * 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 *****************************************************************************/
27#ifndef DUMUX_PIPE_LAUFER_PROBLEM_HH
28#define DUMUX_PIPE_LAUFER_PROBLEM_HH
29
30#include <dune/grid/yaspgrid.hh>
31#include <dune/common/hybridutilities.hh>
32
38
49
50namespace Dumux {
51
52template <class TypeTag>
53class PipeLauferProblem;
54
55namespace Properties {
56
57// Create new type tags
58namespace TTag {
59// Base Typetag
60struct RANSModel { using InheritsFrom = std::tuple<StaggeredFreeFlowModel>; };
61// Isothermal Typetags
62struct PipeLauferZeroEq { using InheritsFrom = std::tuple<RANSModel, ZeroEq>; };
63struct PipeLauferOneEq { using InheritsFrom = std::tuple<RANSModel, OneEq>; };
64struct PipeLauferKOmega { using InheritsFrom = std::tuple<RANSModel, KOmega>; };
65struct PipeLauferLowReKEpsilon { using InheritsFrom = std::tuple<RANSModel, LowReKEpsilon>; };
66struct PipeLauferKEpsilon { using InheritsFrom = std::tuple<RANSModel, KEpsilon>; };
67// Non-Isothermal Typetags
68struct PipeLauferNIZeroEq { using InheritsFrom = std::tuple<RANSModel, ZeroEqNI>; };
69struct PipeLauferNIOneEq { using InheritsFrom = std::tuple<RANSModel, OneEqNI>; };
70struct PipeLauferNIKOmega { using InheritsFrom = std::tuple<RANSModel, KOmegaNI>; };
71struct PipeLauferNILowReKEpsilon { using InheritsFrom = std::tuple<RANSModel, LowReKEpsilonNI>; };
72struct PipeLauferNIKEpsilon { using InheritsFrom = std::tuple<RANSModel, KEpsilonNI>; };
73} // end namespace TTag
74
75
76// the fluid system
77template<class TypeTag>
78struct FluidSystem<TypeTag, TTag::RANSModel>
79{
82};
83
84// Set the grid type
85template<class TypeTag>
86struct Grid<TypeTag, TTag::RANSModel>
87{ using type = Dune::YaspGrid<2, Dune::TensorProductCoordinates<GetPropType<TypeTag, Properties::Scalar>, 2> >; };
88
89// Set the problem property
90template<class TypeTag>
91struct Problem<TypeTag, TTag::RANSModel>
93
94template<class TypeTag>
95struct EnableGridGeometryCache<TypeTag, TTag::RANSModel> { static constexpr bool value = true; };
96
97template<class TypeTag>
98struct EnableGridFluxVariablesCache<TypeTag, TTag::RANSModel> { static constexpr bool value = true; };
99template<class TypeTag>
100struct EnableGridVolumeVariablesCache<TypeTag, TTag::RANSModel> { static constexpr bool value = true; };
101} // end namespace Properties
102
110template <class TypeTag>
111class PipeLauferProblem : public RANSProblem<TypeTag>
112{
114
123
125 using Element = typename GridGeometry::GridView::template Codim<0>::Entity;
126 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
127 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
128 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
129 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
130
131 using TimeLoopPtr = std::shared_ptr<CheckPointTimeLoop<Scalar>>;
132
133 static constexpr auto dimWorld = GridGeometry::GridView::dimensionworld;
134
135public:
136 PipeLauferProblem(std::shared_ptr<const GridGeometry> gridGeometry)
137 : ParentType(gridGeometry), eps_(1e-6)
138 {
139 inletVelocity_ = getParam<Scalar>("Problem.InletVelocity");
140 inletTemperature_ = getParam<Scalar>("Problem.InletTemperature", 283.15);
141 wallTemperature_ = getParam<Scalar>("Problem.WallTemperature", 303.15);
142 sandGrainRoughness_ = getParam<Scalar>("Problem.SandGrainRoughness", 0.0);
143
144 FluidSystem::init();
146 FluidState fluidState;
147 fluidState.setPressure(0, 1e5);
148 fluidState.setTemperature(temperature());
149 Scalar density = FluidSystem::density(fluidState, 0);
150 Scalar kinematicViscosity = FluidSystem::viscosity(fluidState, 0) / density;
151 Scalar diameter = this->gridGeometry().bBoxMax()[1] - this->gridGeometry().bBoxMin()[1];
152
153 // ideally the viscosityTilde parameter as inflow for the Spalart-Allmaras model should be zero
154 viscosityTilde_ = 1e-3 * turbulenceProperties.viscosityTilde(inletVelocity_, diameter, kinematicViscosity);
155 turbulentKineticEnergy_ = turbulenceProperties.turbulentKineticEnergy(inletVelocity_, diameter, kinematicViscosity);
156
157 if (ModelTraits::turbulenceModel() == TurbulenceModel::komega)
158 dissipation_ = turbulenceProperties.dissipationRate(inletVelocity_, diameter, kinematicViscosity);
159 else
160 dissipation_ = turbulenceProperties.dissipation(inletVelocity_, diameter, kinematicViscosity);
161
162 if (ModelTraits::turbulenceModel() == TurbulenceModel::oneeq)
163 initializationTime_ = getParam<Scalar>("TimeLoop.Initialization", 1.0);
164 else
165 initializationTime_ = getParam<Scalar>("TimeLoop.Initialization", -1.0);
166
167 turbulenceModelName_ = turbulenceModelToString(ModelTraits::turbulenceModel());
168 std::cout << "Using the "<< turbulenceModelName_ << " Turbulence Model. \n";
169 std::cout << std::endl;
170 }
171
175 // \{
176
177 bool isOnWallAtPos(const GlobalPosition &globalPos) const
178 {
179 return globalPos[1] < this->gridGeometry().bBoxMin()[1] + eps_
180 || globalPos[1] > this->gridGeometry().bBoxMax()[1] - eps_;
181 }
182
183 Scalar sandGrainRoughnessAtPos(const GlobalPosition &globalPos) const
184 {
185 return sandGrainRoughness_;
186 }
187
189 {
190 return false;
191 }
192
196 Scalar temperature() const
197 { return inletTemperature_; }
198
204 NumEqVector sourceAtPos(const GlobalPosition &globalPos) const
205 {
206 return NumEqVector(0.0);
207 }
208 // \}
209
213 // \{
214
221 BoundaryTypes boundaryTypesAtPos(const GlobalPosition &globalPos) const
222 {
223 BoundaryTypes values;
224
225 // common boundary types for all turbulence models
226 if(isOutlet_(globalPos))
227 values.setDirichlet(Indices::pressureIdx);
228 else
229 {
230 // walls and inflow
231 values.setDirichlet(Indices::velocityXIdx);
232 values.setDirichlet(Indices::velocityYIdx);
233 }
234
235#if NONISOTHERMAL
236 if(isOutlet_(globalPos))
237 values.setOutflow(Indices::energyEqIdx);
238 else
239 values.setDirichlet(Indices::temperatureIdx);
240#endif
241
242 // turbulence model-specific boundary types
243 static constexpr auto numEq = numTurbulenceEq(ModelTraits::turbulenceModel());
244 setBcTypes_(values, globalPos, Dune::index_constant<numEq>{});
245
246 return values;
247 }
248
257 template<class Element, class FVElementGeometry, class SubControlVolume>
258 bool isDirichletCell(const Element& element,
259 const FVElementGeometry& fvGeometry,
260 const SubControlVolume& scv,
261 int pvIdx) const
262 {
263 using IsKOmegaKEpsilon = std::integral_constant<bool, (ModelTraits::turbulenceModel() == TurbulenceModel::komega
264 || ModelTraits::turbulenceModel() == TurbulenceModel::kepsilon)>;
265 return isDirichletCell_(element, fvGeometry, scv, pvIdx, IsKOmegaKEpsilon{});
266 }
267
275 PrimaryVariables dirichlet(const Element& element, const SubControlVolumeFace& scvf) const
276 {
277 const auto globalPos = scvf.ipGlobal();
278 PrimaryVariables values(initialAtPos(globalPos));
279
280#if NONISOTHERMAL
281 values[Indices::temperatureIdx] = isOnWallAtPos(globalPos) ? wallTemperature_ : inletTemperature_;
282#endif
283
284 return values;
285 }
286
294 template<bool enable = (ModelTraits::turbulenceModel() == TurbulenceModel::komega
295 || ModelTraits::turbulenceModel() == TurbulenceModel::kepsilon),
296 std::enable_if_t<!enable, int> = 0>
297 PrimaryVariables dirichlet(const Element& element, const SubControlVolume& scv) const
298 {
299 const auto globalPos = scv.center();
300 PrimaryVariables values(initialAtPos(globalPos));
301 return values;
302 }
303
311 template<bool enable = (ModelTraits::turbulenceModel() == TurbulenceModel::komega
312 || ModelTraits::turbulenceModel() == TurbulenceModel::kepsilon),
313 std::enable_if_t<enable, int> = 0>
314 PrimaryVariables dirichlet(const Element& element, const SubControlVolume& scv) const
315 {
316 using SetDirichletCellForBothTurbEq = std::integral_constant<bool, (ModelTraits::turbulenceModel() == TurbulenceModel::kepsilon)>;
317
318 return dirichletTurbulentTwoEq_(element, scv, SetDirichletCellForBothTurbEq{});
319 }
320
326 PrimaryVariables initialAtPos(const GlobalPosition& globalPos) const
327 {
328 // common initial conditions for all turbulence models
329 PrimaryVariables values(0.0);
330 values[Indices::pressureIdx] = 1.0e+5;
331 values[Indices::velocityXIdx] = time() > initializationTime_
332 ? inletVelocity_
333 : time() / initializationTime_ * inletVelocity_;
334 if (isOnWallAtPos(globalPos))
335 values[Indices::velocityXIdx] = 0.0;
336
337#if NONISOTHERMAL
338 values[Indices::temperatureIdx] = isOnWallAtPos(globalPos) ? wallTemperature_ : inletTemperature_;
339#endif
340
341 // turbulence model-specific initial conditions
342 static constexpr auto numEq = numTurbulenceEq(ModelTraits::turbulenceModel());
343 setInitialAtPos_(values, globalPos, Dune::index_constant<numEq>{});
344
345 return values;
346 }
347
348 // \}
349
350 void setTimeLoop(TimeLoopPtr timeLoop)
351 {
352 timeLoop_ = timeLoop;
353 }
354
355 Scalar time() const
356 {
357 return timeLoop_->time();
358 }
359
360private:
361 bool isInlet_(const GlobalPosition& globalPos) const
362 {
363 return globalPos[0] < this->gridGeometry().bBoxMin()[0] + eps_;
364 }
365
366 bool isOutlet_(const GlobalPosition& globalPos) const
367 {
368 return globalPos[0] > this->gridGeometry().bBoxMax()[0] - eps_;
369 }
370
372 void setInitialAtPos_(PrimaryVariables& values, const GlobalPosition &globalPos, Dune::index_constant<0>) const {}
373
375 void setInitialAtPos_(PrimaryVariables& values, const GlobalPosition &globalPos, Dune::index_constant<1>) const
376 {
377 values[Indices::viscosityTildeIdx] = viscosityTilde_;
378 if (isOnWallAtPos(globalPos))
379 values[Indices::viscosityTildeIdx] = 0.0;
380 }
381
383 void setInitialAtPos_(PrimaryVariables& values, const GlobalPosition &globalPos, Dune::index_constant<2>) const
384 {
385 values[Indices::turbulentKineticEnergyIdx] = turbulentKineticEnergy_;
386 values[Indices::dissipationIdx] = dissipation_;
387 if (isOnWallAtPos(globalPos))
388 {
389 values[Indices::turbulentKineticEnergyIdx] = 0.0;
390 values[Indices::dissipationIdx] = 0.0;
391 }
392 }
393
395 void setBcTypes_(BoundaryTypes& values, const GlobalPosition& pos, Dune::index_constant<0>) const {}
396
398 void setBcTypes_(BoundaryTypes& values, const GlobalPosition& pos, Dune::index_constant<1>) const
399 {
400 if(isOutlet_(pos))
401 values.setOutflow(Indices::viscosityTildeIdx);
402 else // walls and inflow
403 values.setDirichlet(Indices::viscosityTildeIdx);
404 }
405
407 void setBcTypes_(BoundaryTypes& values,const GlobalPosition& pos, Dune::index_constant<2>) const
408 {
409 if(isOutlet_(pos))
410 {
411 values.setOutflow(Indices::turbulentKineticEnergyEqIdx);
412 values.setOutflow(Indices::dissipationEqIdx);
413 }
414 else
415 {
416 // walls and inflow
417 values.setDirichlet(Indices::turbulentKineticEnergyIdx);
418 values.setDirichlet(Indices::dissipationIdx);
419 }
420 }
421
423 template<class Element, class FVElementGeometry, class SubControlVolume>
424 bool isDirichletCell_(const Element& element,
425 const FVElementGeometry& fvGeometry,
426 const SubControlVolume& scv,
427 int pvIdx,
428 std::false_type) const
429 {
430 return ParentType::isDirichletCell(element, fvGeometry, scv, pvIdx);
431 }
432
434 template<class Element, class FVElementGeometry, class SubControlVolume>
435 bool isDirichletCell_(const Element& element,
436 const FVElementGeometry& fvGeometry,
437 const SubControlVolume& scv,
438 int pvIdx,
439 std::true_type) const
440 {
441 using SetDirichletCellForBothTurbEq = std::integral_constant<bool, (ModelTraits::turbulenceModel() == TurbulenceModel::kepsilon)>;
442 return isDirichletCellTurbulentTwoEq_(element, fvGeometry, scv, pvIdx, SetDirichletCellForBothTurbEq{});
443 }
444
446 template<class Element>
447 bool isDirichletCellTurbulentTwoEq_(const Element& element,
448 const FVElementGeometry& fvGeometry,
449 const SubControlVolume& scv,
450 int pvIdx,
451 std::true_type) const
452 {
453 const auto eIdx = this->gridGeometry().elementMapper().index(element);
454
455 // set a fixed turbulent kinetic energy and dissipation near the wall
456 if (this->inNearWallRegion(eIdx))
457 return pvIdx == Indices::turbulentKineticEnergyEqIdx || pvIdx == Indices::dissipationEqIdx;
458
459 // set a fixed dissipation at the matching point
460 if (this->isMatchingPoint(eIdx))
461 return pvIdx == Indices::dissipationEqIdx;// set a fixed dissipation (omega) for all cells at the wall
462
463 return false;
464 }
465
467 template<class Element>
468 bool isDirichletCellTurbulentTwoEq_(const Element& element,
469 const FVElementGeometry& fvGeometry,
470 const SubControlVolume& scv,
471 int pvIdx,
472 std::false_type) const
473 {
474 // set a fixed dissipation (omega) for all cells at the wall
475 for (const auto& scvf : scvfs(fvGeometry))
476 if (isOnWallAtPos(scvf.center()) && pvIdx == Indices::dissipationIdx)
477 return true;
478
479 return false;
480
481 }
482
484 template<class Element, class SubControlVolume>
485 PrimaryVariables dirichletTurbulentTwoEq_(const Element& element,
486 const SubControlVolume& scv,
487 std::true_type) const
488 {
489 const auto globalPos = scv.center();
490 PrimaryVariables values(initialAtPos(globalPos));
491 unsigned int elementIdx = this->gridGeometry().elementMapper().index(element);
492
493 // fixed value for the turbulent kinetic energy
494 values[Indices::turbulentKineticEnergyEqIdx] = this->turbulentKineticEnergyWallFunction(elementIdx);
495
496 // fixed value for the dissipation
497 values[Indices::dissipationEqIdx] = this->dissipationWallFunction(elementIdx);
498
499 return values;
500 }
501
503 template<class Element, class SubControlVolume>
504 PrimaryVariables dirichletTurbulentTwoEq_(const Element& element,
505 const SubControlVolume& scv,
506 std::false_type) const
507 {
508 const auto globalPos = scv.center();
509 PrimaryVariables values(initialAtPos(globalPos));
510 unsigned int elementIdx = this->gridGeometry().elementMapper().index(element);
511
512 const auto wallDistance = ParentType::wallDistance_[elementIdx];
513 using std::pow;
514 values[Indices::dissipationEqIdx] = 6.0 * ParentType::kinematicViscosity_[elementIdx]
515 / (ParentType::betaOmega() * pow(wallDistance, 2));
516 return values;
517 }
518
519 Scalar eps_;
520 Scalar inletVelocity_;
521 Scalar inletTemperature_;
522 Scalar wallTemperature_;
523 Scalar sandGrainRoughness_;
524 Scalar initializationTime_;
525 Scalar viscosityTilde_;
526 Scalar turbulentKineticEnergy_;
527 Scalar dissipation_;
528 std::string turbulenceModelName_;
529 TimeLoopPtr timeLoop_;
530};
531} // end namespace Dumux
532
533#endif
The available free flow turbulence models in Dumux.
This file contains different functions for estimating turbulence properties.
A simple class for the air fluid properties.
A gaseous phase consisting of a single component.
Geometry::ctype diameter(const Geometry &geo)
Computes the longest distance between points of a geometry.
Definition: diameter.hh:36
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
std::string turbulenceModelToString(TurbulenceModel turbulenceModel)
return the name of the Turbulence Model
Definition: turbulencemodel.hh:54
constexpr unsigned int numTurbulenceEq(TurbulenceModel model)
Definition: turbulencemodel.hh:41
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 density(int phaseIdx) noexcept
I/O name of density for multiphase systems.
Definition: name.hh:65
The DUNE grid type.
Definition: common/properties.hh:57
Property to specify the type of a problem which has to be solved.
Definition: common/properties.hh:69
Definition: common/properties.hh:169
If disabled, the volume variables are not stored (reduces memory, but is slower)
Definition: common/properties.hh:178
specifies if data on flux vars should be saved (faster, but more memory consuming)
Definition: common/properties.hh:188
The type of the fluid system to use.
Definition: common/properties.hh:223
forward declare
Definition: dumux/freeflow/rans/problem.hh:41
This class contains different functions for estimating turbulence properties.
Definition: turbulenceproperties.hh:39
Scalar viscosityTilde(const Scalar velocity, const Scalar diameter, const Scalar kinematicViscosity, bool print=verbose) const
Estimates the viscosity tilde based on a formula given in in the ANSYS Fluent user guide .
Definition: turbulenceproperties.hh:205
Scalar turbulentKineticEnergy(const Scalar velocity, const Scalar diameter, const Scalar kinematicViscosity, bool print=verbose) const
Estimates the turbulent kinetic energy based on a formula given in the ANSYS Fluent user guide .
Definition: turbulenceproperties.hh:143
Scalar dissipationRate(const Scalar velocity, const Scalar diameter, const Scalar kinematicViscosity, bool print=verbose) const
Estimates the dissipation rate based on a formula given in the ANSYS Fluent user guide .
Definition: turbulenceproperties.hh:183
Scalar dissipation(const Scalar velocity, const Scalar diameter, const Scalar kinematicViscosity, bool print=verbose) const
Estimates the dissipation based on a formula given in the ANSYS Fluent user guide .
Definition: turbulenceproperties.hh:162
A gaseous phase consisting of a single component.
Definition: 1pgas.hh:46
Test problem for the one-phase (Navier-) Stokes problem in a channel.
Definition: test/freeflow/rans/problem.hh:112
BoundaryTypes boundaryTypesAtPos(const GlobalPosition &globalPos) const
Specifies which kind of boundary condition should be used for which equation on a given boundary cont...
Definition: test/freeflow/rans/problem.hh:221
bool shouldWriteRestartFile() const
Definition: test/freeflow/rans/problem.hh:188
Scalar sandGrainRoughnessAtPos(const GlobalPosition &globalPos) const
Definition: test/freeflow/rans/problem.hh:183
PipeLauferProblem(std::shared_ptr< const GridGeometry > gridGeometry)
Definition: test/freeflow/rans/problem.hh:136
PrimaryVariables dirichlet(const Element &element, const SubControlVolume &scv) const
Evaluate the boundary conditions for fixed values at cell centers.
Definition: test/freeflow/rans/problem.hh:297
PrimaryVariables dirichlet(const Element &element, const SubControlVolumeFace &scvf) const
Evaluate the boundary conditions for a dirichlet values at the boundary.
Definition: test/freeflow/rans/problem.hh:275
PrimaryVariables initialAtPos(const GlobalPosition &globalPos) const
Evaluate the initial value for a control volume.
Definition: test/freeflow/rans/problem.hh:326
NumEqVector sourceAtPos(const GlobalPosition &globalPos) const
Returns the sources within the domain.
Definition: test/freeflow/rans/problem.hh:204
void setTimeLoop(TimeLoopPtr timeLoop)
Definition: test/freeflow/rans/problem.hh:350
Scalar temperature() const
Returns the temperature [K] within the domain for the isothermal model.
Definition: test/freeflow/rans/problem.hh:196
bool isOnWallAtPos(const GlobalPosition &globalPos) const
Definition: test/freeflow/rans/problem.hh:177
Scalar time() const
Definition: test/freeflow/rans/problem.hh:355
bool isDirichletCell(const Element &element, const FVElementGeometry &fvGeometry, const SubControlVolume &scv, int pvIdx) const
Returns whether a fixed Dirichlet value shall be used at a given cell.
Definition: test/freeflow/rans/problem.hh:258
Definition: test/freeflow/rans/problem.hh:60
std::tuple< StaggeredFreeFlowModel > InheritsFrom
Definition: test/freeflow/rans/problem.hh:60
Definition: test/freeflow/rans/problem.hh:62
std::tuple< RANSModel, ZeroEq > InheritsFrom
Definition: test/freeflow/rans/problem.hh:62
Definition: test/freeflow/rans/problem.hh:63
std::tuple< RANSModel, OneEq > InheritsFrom
Definition: test/freeflow/rans/problem.hh:63
Definition: test/freeflow/rans/problem.hh:64
std::tuple< RANSModel, KOmega > InheritsFrom
Definition: test/freeflow/rans/problem.hh:64
Definition: test/freeflow/rans/problem.hh:65
std::tuple< RANSModel, LowReKEpsilon > InheritsFrom
Definition: test/freeflow/rans/problem.hh:65
Definition: test/freeflow/rans/problem.hh:66
std::tuple< RANSModel, KEpsilon > InheritsFrom
Definition: test/freeflow/rans/problem.hh:66
Definition: test/freeflow/rans/problem.hh:68
std::tuple< RANSModel, ZeroEqNI > InheritsFrom
Definition: test/freeflow/rans/problem.hh:68
Definition: test/freeflow/rans/problem.hh:69
std::tuple< RANSModel, OneEqNI > InheritsFrom
Definition: test/freeflow/rans/problem.hh:69
Definition: test/freeflow/rans/problem.hh:70
std::tuple< RANSModel, KOmegaNI > InheritsFrom
Definition: test/freeflow/rans/problem.hh:70
Definition: test/freeflow/rans/problem.hh:71
std::tuple< RANSModel, LowReKEpsilonNI > InheritsFrom
Definition: test/freeflow/rans/problem.hh:71
Definition: test/freeflow/rans/problem.hh:72
std::tuple< RANSModel, KEpsilonNI > InheritsFrom
Definition: test/freeflow/rans/problem.hh:72
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: test/freeflow/rans/problem.hh:80
Dune::YaspGrid< 2, Dune::TensorProductCoordinates< GetPropType< TypeTag, Properties::Scalar >, 2 > > type
Definition: test/freeflow/rans/problem.hh:87
Defines a type tag and some properties for ree-flow models using the staggered scheme....
A single-phase, isothermal one-equation turbulence model by Spalart-Allmaras.
A single-phase, isothermal k-epsilon model.
A single-phase, isothermal k-omega 2-Eq. model.
A single-phase, isothermal low-Reynolds k-epsilon model.
A single-phase, isothermal Reynolds-Averaged Navier-Stokes 0-Eq. model.
One-equation turbulence problem base class.
K-epsilon turbulence problem base class.
K-Omega turbulence model problem base class.
Low-Re k-epsilon turbulence problem base class.