3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
test/freeflow/ransnc/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 *****************************************************************************/
24#ifndef DUMUX_RANS_NC_TEST_PROBLEM_HH
25#define DUMUX_RANS_NC_TEST_PROBLEM_HH
26
27#include <dune/grid/yaspgrid.hh>
28
33
36
39
42
45
48
49
50namespace Dumux {
51
52template <class TypeTag>
53class FlatPlateNCTestProblem;
54
55namespace Properties {
56
57// Create new type tags
58namespace TTag {
59// Base Typetag
60struct RANSNCModel { using InheritsFrom = std::tuple<StaggeredFreeFlowModel>; };
61// Isothermal Typetags
62struct FlatPlateNCZeroEq { using InheritsFrom = std::tuple<RANSNCModel, ZeroEqNC>; };
63struct FlatPlateNCOneEq { using InheritsFrom = std::tuple<RANSNCModel, OneEqNC>; };
64struct FlatPlateNCKOmega { using InheritsFrom = std::tuple<RANSNCModel, KOmegaNC>; };
65struct FlatPlateNCLowReKEpsilon { using InheritsFrom = std::tuple<RANSNCModel, LowReKEpsilonNC>; };
66struct FlatPlateNCKEpsilon { using InheritsFrom = std::tuple<RANSNCModel, KEpsilonNC>; };
67// Isothermal Typetags
68struct FlatPlateNCNIZeroEq { using InheritsFrom = std::tuple<RANSNCModel, ZeroEqNCNI>; };
69struct FlatPlateNCNIOneEq { using InheritsFrom = std::tuple<RANSNCModel, OneEqNCNI>; };
70struct FlatPlateNCNIKOmega { using InheritsFrom = std::tuple<RANSNCModel, KOmegaNCNI>; };
71struct FlatPlateNCNILowReKEpsilon { using InheritsFrom = std::tuple<RANSNCModel, LowReKEpsilonNCNI>; };
72struct FlatPlateNCNIKEpsilon { using InheritsFrom = std::tuple<RANSNCModel, KEpsilonNCNI>; };
73} // end namespace TTag
74
75// The fluid system
76template<class TypeTag>
77struct FluidSystem<TypeTag, TTag::RANSNCModel>
78{
80 static constexpr auto phaseIdx = H2OAir::gasPhaseIdx; // simulate the air phase
82};
83
84// replace the main component balance eq with a total balance eq
85template<class TypeTag>
86struct ReplaceCompEqIdx<TypeTag, TTag::RANSNCModel> { static constexpr int value = 0; };
87
88// Set the grid type
89template<class TypeTag>
90struct Grid<TypeTag, TTag::RANSNCModel>
91{ using type = Dune::YaspGrid<2, Dune::TensorProductCoordinates<GetPropType<TypeTag, Properties::Scalar>, 2> >; };
92
93// Set the problem property
94template<class TypeTag>
95struct Problem<TypeTag, TTag::RANSNCModel> { using type = Dumux::FlatPlateNCTestProblem<TypeTag> ; };
96
97template<class TypeTag>
98struct EnableGridGeometryCache<TypeTag, TTag::RANSNCModel> { static constexpr bool value = true; };
99template<class TypeTag>
100struct EnableGridFluxVariablesCache<TypeTag, TTag::RANSNCModel> { static constexpr bool value = true; };
101template<class TypeTag>
102struct EnableGridVolumeVariablesCache<TypeTag, TTag::RANSNCModel> { static constexpr bool value = true; };
103
104template<class TypeTag>
105struct UseMoles<TypeTag, TTag::RANSNCModel> { static constexpr bool value = true; };
106} // end namespace Properties
107
117template <class TypeTag>
118class FlatPlateNCTestProblem : public RANSProblem<TypeTag>
119{
121
130
132 using Element = typename GridGeometry::GridView::template Codim<0>::Entity;
133 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
134 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
135 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
136 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
137
138 using TimeLoopPtr = std::shared_ptr<CheckPointTimeLoop<Scalar>>;
139
140 static constexpr auto dimWorld = GetPropType<TypeTag, Properties::GridView>::dimensionworld;
141 static constexpr auto transportEqIdx = Indices::conti0EqIdx + 1;
142 static constexpr auto transportCompIdx = Indices::conti0EqIdx + 1;
143
144public:
145 FlatPlateNCTestProblem(std::shared_ptr<const GridGeometry> gridGeometry)
146 : ParentType(gridGeometry), eps_(1e-6)
147 {
148 inletVelocity_ = getParam<Scalar>("Problem.InletVelocity", 0.1);
149 inletTemperature_ = getParam<Scalar>("Problem.InletTemperature", 283.15);
150 wallTemperature_ = getParam<Scalar>("Problem.WallTemperature", 313.15);
151 inletMoleFraction_ = getParam<Scalar>("Problem.InletMoleFraction", 1e-3);
152
153 FluidSystem::init();
155 FluidState fluidState;
156 const auto phaseIdx = 0;
157 fluidState.setPressure(phaseIdx, 1e5);
158 fluidState.setTemperature(temperature());
159 fluidState.setMassFraction(phaseIdx, phaseIdx, 1.0);
160 Scalar density = FluidSystem::density(fluidState, phaseIdx);
161 Scalar kinematicViscosity = FluidSystem::viscosity(fluidState, phaseIdx) / density;
162 Scalar diameter = this->gridGeometry().bBoxMax()[1] - this->gridGeometry().bBoxMin()[1];
163 viscosityTilde_ = 1e-3 * turbulenceProperties.viscosityTilde(inletVelocity_, diameter, kinematicViscosity);
164 turbulentKineticEnergy_ = turbulenceProperties.turbulentKineticEnergy(inletVelocity_, diameter, kinematicViscosity);
165 if (ModelTraits::turbulenceModel() == TurbulenceModel::komega)
166 dissipation_ = turbulenceProperties.dissipationRate(inletVelocity_, diameter, kinematicViscosity);
167 else
168 dissipation_ = turbulenceProperties.dissipation(inletVelocity_, diameter, kinematicViscosity);
169 turbulenceModelName_ = turbulenceModelToString(ModelTraits::turbulenceModel());
170 std::cout << "Using the "<< turbulenceModelName_ << " Turbulence Model. \n";
171 std::cout << std::endl;
172 }
173
177 // \{
178
179 bool isOnWallAtPos(const GlobalPosition& globalPos) const
180 {
181 return globalPos[1] < eps_;
182 }
183
184
186 {
187 return false;
188 }
189
195 Scalar temperature() const
196 { return 273.15 + 10; } // 10C
197
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 // turbulence model-specific boundary types
226 static constexpr auto numEq = numTurbulenceEq(ModelTraits::turbulenceModel());
227 setBcTypes_(values, globalPos, Dune::index_constant<numEq>{});
228
229 if(isInlet_(globalPos))
230 {
231 values.setDirichlet(Indices::velocityXIdx);
232 values.setDirichlet(Indices::velocityYIdx);
233 values.setDirichlet(transportCompIdx);
234#if NONISOTHERMAL
235 values.setDirichlet(Indices::temperatureIdx);
236#endif
237 }
238 else if(isOutlet_(globalPos))
239 {
240 values.setDirichlet(Indices::pressureIdx);
241 values.setOutflow(transportEqIdx);
242#if NONISOTHERMAL
243 values.setOutflow(Indices::energyEqIdx);
244#endif
245 }
246 else if(isOnWallAtPos(globalPos))
247 {
248 values.setDirichlet(Indices::velocityXIdx);
249 values.setDirichlet(Indices::velocityYIdx);
250 values.setNeumann(transportEqIdx);
251#if NONISOTHERMAL
252 values.setDirichlet(Indices::temperatureIdx);
253#endif
254 }
255 else
256 values.setAllSymmetry();
257
258 return values;
259 }
260
269 template<class Element, class FVElementGeometry, class SubControlVolume>
270 bool isDirichletCell(const Element& element,
271 const FVElementGeometry& fvGeometry,
272 const SubControlVolume& scv,
273 int pvIdx) const
274 {
275 using IsKOmegaKEpsilon = std::integral_constant<bool, (ModelTraits::turbulenceModel() == TurbulenceModel::komega
276 || ModelTraits::turbulenceModel() == TurbulenceModel::kepsilon)>;
277 return isDirichletCell_(element, fvGeometry, scv, pvIdx, IsKOmegaKEpsilon{});
278 }
279
287 PrimaryVariables dirichlet(const Element& element, const SubControlVolumeFace& scvf) const
288 {
289 const auto globalPos = scvf.ipGlobal();
290 PrimaryVariables values(initialAtPos(globalPos));
291
292 if (isInlet_(globalPos)
293 && globalPos[1] > 0.4 * this->gridGeometry().bBoxMax()[1]
294 && globalPos[1] < 0.6 * this->gridGeometry().bBoxMax()[1])
295 {
296 values[transportCompIdx] = (time() > 10.0) ? inletMoleFraction_ : 0.0;
297 }
298
299#if NONISOTHERMAL
300 if (time() > 10.0 && isOnWallAtPos(globalPos))
301 {
302 values[Indices::temperatureIdx] = wallTemperature_;
303 }
304#endif
305
306 return values;
307 }
308
316 template<bool enable = (ModelTraits::turbulenceModel() == TurbulenceModel::komega
317 || ModelTraits::turbulenceModel() == TurbulenceModel::kepsilon),
318 std::enable_if_t<!enable, int> = 0>
319 PrimaryVariables dirichlet(const Element& element, const SubControlVolume& scv) const
320 {
321 const auto globalPos = scv.center();
322 PrimaryVariables values(initialAtPos(globalPos));
323 return values;
324 }
325
333 template<bool enable = (ModelTraits::turbulenceModel() == TurbulenceModel::komega
334 || ModelTraits::turbulenceModel() == TurbulenceModel::kepsilon),
335 std::enable_if_t<enable, int> = 0>
336 PrimaryVariables dirichlet(const Element& element, const SubControlVolume& scv) const
337 {
338 using SetDirichletCellForBothTurbEq = std::integral_constant<bool, (ModelTraits::turbulenceModel() == TurbulenceModel::kepsilon)>;
339
340 return dirichletTurbulentTwoEq_(element, scv, SetDirichletCellForBothTurbEq{});
341 }
342
348 PrimaryVariables initialAtPos(const GlobalPosition &globalPos) const
349 {
350 PrimaryVariables values(0.0);
351 values[Indices::pressureIdx] = 1.0e+5;
352 values[transportCompIdx] = 0.0;
353#if NONISOTHERMAL
354 values[Indices::temperatureIdx] = temperature();
355#endif
356
357 // block velocity profile
358 values[Indices::velocityXIdx] = 0.0;
359 if (!isOnWallAtPos(globalPos))
360 values[Indices::velocityXIdx] = inletVelocity_;
361 values[Indices::velocityYIdx] = 0.0;
362
363 // turbulence model-specific initial conditions
364 static constexpr auto numEq = numTurbulenceEq(ModelTraits::turbulenceModel());
365 setInitialAtPos_(values, globalPos, Dune::index_constant<numEq>{});
366
367 return values;
368 }
369
370 // \}
371
372 void setTimeLoop(TimeLoopPtr timeLoop)
373 {
374 timeLoop_ = timeLoop;
375 }
376
377 Scalar time() const
378 {
379 return timeLoop_->time();
380 }
381
382private:
383 bool isInlet_(const GlobalPosition& globalPos) const
384 {
385 return globalPos[0] < eps_;
386 }
387
388 bool isOutlet_(const GlobalPosition& globalPos) const
389 {
390 return globalPos[0] > this->gridGeometry().bBoxMax()[0] - eps_;
391 }
392
394 void setInitialAtPos_(PrimaryVariables& values, const GlobalPosition &globalPos, Dune::index_constant<0>) const {}
395
397 void setInitialAtPos_(PrimaryVariables& values, const GlobalPosition &globalPos, Dune::index_constant<1>) const
398 {
399 values[Indices::viscosityTildeIdx] = viscosityTilde_;
400 if (isOnWallAtPos(globalPos))
401 values[Indices::viscosityTildeIdx] = 0.0;
402 }
403
405 void setInitialAtPos_(PrimaryVariables& values, const GlobalPosition &globalPos, Dune::index_constant<2>) const
406 {
407 values[Indices::turbulentKineticEnergyIdx] = turbulentKineticEnergy_;
408 values[Indices::dissipationIdx] = dissipation_;
409 if (isOnWallAtPos(globalPos))
410 {
411 values[Indices::turbulentKineticEnergyIdx] = 0.0;
412 values[Indices::dissipationIdx] = 0.0;
413 }
414 }
415
417 void setBcTypes_(BoundaryTypes& values, const GlobalPosition& pos, Dune::index_constant<0>) const {}
418
420 void setBcTypes_(BoundaryTypes& values, const GlobalPosition& pos, Dune::index_constant<1>) const
421 {
422 if(isOutlet_(pos))
423 values.setOutflow(Indices::viscosityTildeIdx);
424 else // walls and inflow
425 values.setDirichlet(Indices::viscosityTildeIdx);
426 }
427
429 void setBcTypes_(BoundaryTypes& values,const GlobalPosition& pos, Dune::index_constant<2>) const
430 {
431 if(isOutlet_(pos))
432 {
433 values.setOutflow(Indices::turbulentKineticEnergyEqIdx);
434 values.setOutflow(Indices::dissipationEqIdx);
435 }
436 else
437 {
438 // walls and inflow
439 values.setDirichlet(Indices::turbulentKineticEnergyIdx);
440 values.setDirichlet(Indices::dissipationIdx);
441 }
442 }
443
445 template<class Element, class FVElementGeometry, class SubControlVolume>
446 bool isDirichletCell_(const Element& element,
447 const FVElementGeometry& fvGeometry,
448 const SubControlVolume& scv,
449 int pvIdx,
450 std::false_type) const
451 {
452 return ParentType::isDirichletCell(element, fvGeometry, scv, pvIdx);
453 }
454
456 template<class Element, class FVElementGeometry, class SubControlVolume>
457 bool isDirichletCell_(const Element& element,
458 const FVElementGeometry& fvGeometry,
459 const SubControlVolume& scv,
460 int pvIdx,
461 std::true_type) const
462 {
463 using SetDirichletCellForBothTurbEq = std::integral_constant<bool, (ModelTraits::turbulenceModel() == TurbulenceModel::kepsilon)>;
464 return isDirichletCellTurbulentTwoEq_(element, fvGeometry, scv, pvIdx, SetDirichletCellForBothTurbEq{});
465 }
466
468 template<class Element>
469 bool isDirichletCellTurbulentTwoEq_(const Element& element,
470 const FVElementGeometry& fvGeometry,
471 const SubControlVolume& scv,
472 int pvIdx,
473 std::true_type) const
474 {
475 const auto eIdx = this->gridGeometry().elementMapper().index(element);
476
477 // set a fixed turbulent kinetic energy and dissipation near the wall
478 if (this->inNearWallRegion(eIdx))
479 return pvIdx == Indices::turbulentKineticEnergyEqIdx || pvIdx == Indices::dissipationEqIdx;
480
481 // set a fixed dissipation at the matching point
482 if (this->isMatchingPoint(eIdx))
483 return pvIdx == Indices::dissipationEqIdx;// set a fixed dissipation (omega) for all cells at the wall
484
485 return false;
486 }
487
489 template<class Element>
490 bool isDirichletCellTurbulentTwoEq_(const Element& element,
491 const FVElementGeometry& fvGeometry,
492 const SubControlVolume& scv,
493 int pvIdx,
494 std::false_type) const
495 {
496 // set a fixed dissipation (omega) for all cells at the wall
497 for (const auto& scvf : scvfs(fvGeometry))
498 if (isOnWallAtPos(scvf.center()) && pvIdx == Indices::dissipationIdx)
499 return true;
500
501 return false;
502
503 }
504
506 template<class Element, class SubControlVolume>
507 PrimaryVariables dirichletTurbulentTwoEq_(const Element& element,
508 const SubControlVolume& scv,
509 std::true_type) const
510 {
511 const auto globalPos = scv.center();
512 PrimaryVariables values(initialAtPos(globalPos));
513 unsigned int elementIdx = this->gridGeometry().elementMapper().index(element);
514
515 // fixed value for the turbulent kinetic energy
516 values[Indices::turbulentKineticEnergyEqIdx] = this->turbulentKineticEnergyWallFunction(elementIdx);
517
518 // fixed value for the dissipation
519 values[Indices::dissipationEqIdx] = this->dissipationWallFunction(elementIdx);
520
521 return values;
522 }
523
525 template<class Element, class SubControlVolume>
526 PrimaryVariables dirichletTurbulentTwoEq_(const Element& element,
527 const SubControlVolume& scv,
528 std::false_type) const
529 {
530 const auto globalPos = scv.center();
531 PrimaryVariables values(initialAtPos(globalPos));
532 unsigned int elementIdx = this->gridGeometry().elementMapper().index(element);
533
534 const auto wallDistance = ParentType::wallDistance_[elementIdx];
535 using std::pow;
536 values[Indices::dissipationEqIdx] = 6.0 * ParentType::kinematicViscosity_[elementIdx]
537 / (ParentType::betaOmega() * pow(wallDistance, 2));
538 return values;
539 }
540
541 const Scalar eps_;
542 Scalar inletVelocity_;
543 Scalar inletMoleFraction_;
544 Scalar inletTemperature_;
545 Scalar wallTemperature_;
546 Scalar viscosityTilde_;
547 Scalar turbulentKineticEnergy_;
548 Scalar dissipation_;
549 TimeLoopPtr timeLoop_;
550 std::string turbulenceModelName_;
551};
552} // end namespace Dumux
553
554#endif
A single-phase, multi-component k-epsilon model.
A single-phase, multi-component k-omega model.
A single-phase, multi-component low-Re k-epsilon model.
A single-phase, multi-component one-equation model.
A single-phase, multi-component Reynolds-Averaged Navier-Stokes 0-Eq. model.
This file contains different functions for estimating turbulence properties.
An adapter for multi-phase fluid systems to be used with (compositional) one-phase models.
A compositional two-phase fluid system with water and air as components in both, the liquid and the g...
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
Property whether to use moles or kg as amount unit for balance equations.
Definition: common/properties.hh:102
The component balance index that should be replaced by the total mass/mole balance.
Definition: common/properties.hh:104
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
An adapter for multi-phase fluid systems to be used with (compositional) one-phase models.
Definition: 1padapter.hh:46
A compositional two-phase fluid system with water and air as components in both, the liquid and the g...
Definition: h2oair.hh:75
Test problem for the one-phase model.
Definition: test/freeflow/ransnc/problem.hh:119
Scalar temperature() const
Returns the temperature within the domain in [K].
Definition: test/freeflow/ransnc/problem.hh:195
NumEqVector sourceAtPos(const GlobalPosition &globalPos) const
Returns the sources within the domain.
Definition: test/freeflow/ransnc/problem.hh:204
PrimaryVariables dirichlet(const Element &element, const SubControlVolume &scv) const
Evaluate the boundary conditions for fixed values at cell centers.
Definition: test/freeflow/ransnc/problem.hh:319
PrimaryVariables initialAtPos(const GlobalPosition &globalPos) const
Evaluates the initial value for a control volume.
Definition: test/freeflow/ransnc/problem.hh:348
FlatPlateNCTestProblem(std::shared_ptr< const GridGeometry > gridGeometry)
Definition: test/freeflow/ransnc/problem.hh:145
void setTimeLoop(TimeLoopPtr timeLoop)
Definition: test/freeflow/ransnc/problem.hh:372
bool shouldWriteRestartFile() const
Definition: test/freeflow/ransnc/problem.hh:185
bool isOnWallAtPos(const GlobalPosition &globalPos) const
Definition: test/freeflow/ransnc/problem.hh:179
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/ransnc/problem.hh:270
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/ransnc/problem.hh:221
Scalar time() const
Definition: test/freeflow/ransnc/problem.hh:377
PrimaryVariables dirichlet(const Element &element, const SubControlVolumeFace &scvf) const
Evaluate the boundary conditions for a dirichlet values at the boundary.
Definition: test/freeflow/ransnc/problem.hh:287
Definition: test/freeflow/ransnc/problem.hh:60
std::tuple< StaggeredFreeFlowModel > InheritsFrom
Definition: test/freeflow/ransnc/problem.hh:60
Definition: test/freeflow/ransnc/problem.hh:62
std::tuple< RANSNCModel, ZeroEqNC > InheritsFrom
Definition: test/freeflow/ransnc/problem.hh:62
Definition: test/freeflow/ransnc/problem.hh:63
std::tuple< RANSNCModel, OneEqNC > InheritsFrom
Definition: test/freeflow/ransnc/problem.hh:63
Definition: test/freeflow/ransnc/problem.hh:64
std::tuple< RANSNCModel, KOmegaNC > InheritsFrom
Definition: test/freeflow/ransnc/problem.hh:64
Definition: test/freeflow/ransnc/problem.hh:65
std::tuple< RANSNCModel, LowReKEpsilonNC > InheritsFrom
Definition: test/freeflow/ransnc/problem.hh:65
Definition: test/freeflow/ransnc/problem.hh:66
std::tuple< RANSNCModel, KEpsilonNC > InheritsFrom
Definition: test/freeflow/ransnc/problem.hh:66
Definition: test/freeflow/ransnc/problem.hh:68
std::tuple< RANSNCModel, ZeroEqNCNI > InheritsFrom
Definition: test/freeflow/ransnc/problem.hh:68
Definition: test/freeflow/ransnc/problem.hh:69
std::tuple< RANSNCModel, OneEqNCNI > InheritsFrom
Definition: test/freeflow/ransnc/problem.hh:69
Definition: test/freeflow/ransnc/problem.hh:70
std::tuple< RANSNCModel, KOmegaNCNI > InheritsFrom
Definition: test/freeflow/ransnc/problem.hh:70
Definition: test/freeflow/ransnc/problem.hh:71
std::tuple< RANSNCModel, LowReKEpsilonNCNI > InheritsFrom
Definition: test/freeflow/ransnc/problem.hh:71
Definition: test/freeflow/ransnc/problem.hh:72
std::tuple< RANSNCModel, KEpsilonNCNI > InheritsFrom
Definition: test/freeflow/ransnc/problem.hh:72
Dune::YaspGrid< 2, Dune::TensorProductCoordinates< GetPropType< TypeTag, Properties::Scalar >, 2 > > type
Definition: test/freeflow/ransnc/problem.hh:91
Defines a type tag and some properties for ree-flow models using the staggered scheme....
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.
Zero-equation turbulence problem base class.