3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
turbulencemodel.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 2 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_FREEFLOW_TURBLENCEMODEL_HH
25#define DUMUX_FREEFLOW_TURBLENCEMODEL_HH
26
27#include <string>
29
30namespace Dumux {
31
37 enum class TurbulenceModel
38 {
40 };
41
42 constexpr unsigned int numTurbulenceEq(TurbulenceModel model)
43 {
44 if (model == TurbulenceModel::none || model == TurbulenceModel::zeroeq)
45 return 0;
46 else if (model == TurbulenceModel::oneeq)
47 return 1;
48 else
49 return 2;
50 }
51
55 std::string turbulenceModelToString(TurbulenceModel turbulenceModel)
56 {
57 switch (turbulenceModel)
58 {
59 case TurbulenceModel::none: return "No_TurbModel";
60 case TurbulenceModel::zeroeq: return "ZeroEq";
61 case TurbulenceModel::oneeq: return "OneEq";
62 case TurbulenceModel::kepsilon: return "KEpsilon";
63 case TurbulenceModel::lowrekepsilon: return "LowReKEpsilon";
64 case TurbulenceModel::komega: return "KOmega";
65 case TurbulenceModel::sst: return "KOmegaSST";
66 default: return "Invalid"; // should never be reached
67 }
68 }
69
74 enum class SSTModel
75 { BSL, SST };
76
80 std::string sstModelToString(SSTModel sstModel)
81 {
82 switch (sstModel)
83 {
84 case SSTModel::BSL: return "BSL";
85 case SSTModel::SST: return "SST";
86 default: return "Invalid";
87 }
88 }
89
94 SSTModel sstModelFromString(const std::string& sstModel)
95 {
96 if (sstModel == "BSL") return SSTModel::BSL;
97 if (sstModel == "SST") return SSTModel::SST;
98 DUNE_THROW(ParameterException, "\nThis SST Model approach : \"" << sstModel << "\" is not implemented.\n"
99 << "The available SST models are as follows: \n"
100 << sstModelToString(SSTModel::BSL) << ": The Baseline SST Model n\n"
101 << sstModelToString(SSTModel::SST) << ": The full standard SST Model");
102 }
103
104} // end namespace Dumux
105
106#endif
Some exceptions thrown in DuMux
TurbulenceModel
The available free flow turbulence models in Dumux.
Definition: turbulencemodel.hh:38
SSTModel
The available variations of the SST Turbulence Model.
Definition: turbulencemodel.hh:75
Definition: adapt.hh:29
std::string turbulenceModelToString(TurbulenceModel turbulenceModel)
return the name of the Turbulence Model
Definition: turbulencemodel.hh:55
std::string sstModelToString(SSTModel sstModel)
return the name of the sst Model as a string
Definition: turbulencemodel.hh:80
constexpr unsigned int numTurbulenceEq(TurbulenceModel model)
Definition: turbulencemodel.hh:42
SSTModel sstModelFromString(const std::string &sstModel)
Convenience function to convert user input given as std::string to the corresponding enum class used ...
Definition: turbulencemodel.hh:94
Exception thrown if a run-time parameter is not specified correctly.
Definition: exceptions.hh:60