kOmegaSSTDES.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2016 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "kOmegaSSTDES.H"
27 
28 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
29 
30 namespace Foam
31 {
32 namespace LESModels
33 {
34 
35 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
36 
37 template<class BasicTurbulenceModel>
39 {
40  return sqrt(this->k_)/(this->betaStar_*this->omega_);
41 }
42 
43 
44 template<class BasicTurbulenceModel>
46 (
47  const volScalarField& F1,
48  const volScalarField& F2
49 ) const
50 {
51  switch (FSST_)
52  {
53  case 0:
54  return max(Lt()/(CDES_*this->delta()), scalar(1));
55  case 1:
56  return max(Lt()*(1 - F1)/(CDES_*this->delta()), scalar(1));
57  case 2:
58  return max(Lt()*(1 - F2)/(CDES_*this->delta()), scalar(1));
59  default:
61  << "Incorrect FSST = " << FSST_ << ", should be 0, 1 or 2"
62  << exit(FatalError);
63  return F1;
64  }
65 }
66 
67 
68 template<class BasicTurbulenceModel>
70 (
71  const volScalarField& F1,
72  const volScalarField& F2
73 ) const
74 {
75  return this->betaStar_*this->omega_*FDES(F1, F2);
76 }
77 
78 
79 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
80 
81 template<class BasicTurbulenceModel>
83 (
84  const alphaField& alpha,
85  const rhoField& rho,
86  const volVectorField& U,
87  const surfaceScalarField& alphaRhoPhi,
88  const surfaceScalarField& phi,
89  const transportModel& transport,
90  const word& propertiesName,
91  const word& type
92 )
93 :
94  kOmegaSST
95  <
96  LESeddyViscosity<BasicTurbulenceModel>,
97  BasicTurbulenceModel
98  >
99  (
100  type,
101  alpha,
102  rho,
103  U,
104  alphaRhoPhi,
105  phi,
106  transport,
107  propertiesName
108  ),
109 
110  CDES_
111  (
113  (
114  "CDES",
115  this->coeffDict_,
116  0.61
117  )
118  ),
119  FSST_(this->coeffDict_.lookupOrDefault("FSST", 2))
120 {
121  if (type == typeName)
122  {
123  this->printCoeffs(type);
124  }
125 }
126 
127 
128 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
129 
130 template<class BasicTurbulenceModel>
132 {
133  if
134  (
135  kOmegaSST<LESeddyViscosity<BasicTurbulenceModel>, BasicTurbulenceModel>
136  ::read()
137  )
138  {
139  CDES_.readIfPresent(this->coeffDict());
140  this->coeffDict().readIfPresent("FSST", FSST_);
141 
142  return true;
143  }
144  else
145  {
146  return false;
147  }
148 }
149 
150 
151 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
152 
153 } // End namespace LESModels
154 } // End namespace Foam
155 
156 // ************************************************************************* //
scalar delta
surfaceScalarField & phi
U
Definition: pEqn.H:83
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define F1(B, C, D)
Definition: SHA1.C:173
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
dimensionedScalar sqrt(const dimensionedScalar &ds)
BasicTurbulenceModel::rhoField rhoField
Definition: kOmegaSSTDES.H:124
Generic dimensioned Type class.
BasicTurbulenceModel::transportModel transportModel
Definition: kOmegaSSTDES.H:125
virtual tmp< volScalarField > epsilonByk(const volScalarField &F1, const volScalarField &F2) const
Return epsilon/k which for standard RAS is betaStar*omega.
Definition: kOmegaSSTDES.C:70
bool read(const char *, int32_t &)
Definition: int32IO.C:85
A class for handling words, derived from string.
Definition: word.H:59
tmp< volScalarField > Lt() const
Return the turbulent length-scale.
Definition: kOmegaSSTDES.C:38
Implementation of the k-omega-SST-DES turbulence model for incompressible and compressible flows...
virtual tmp< volScalarField > FDES(const volScalarField &F1, const volScalarField &F2) const
The DES dissipation-rate multiplier with options zonal filtering.
Definition: kOmegaSSTDES.C:46
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:461
kOmegaSSTDES(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName, const word &type=typeName)
Construct from components.
Definition: kOmegaSSTDES.C:83
A class for managing temporary objects.
Definition: PtrList.H:54
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Namespace for OpenFOAM.
virtual bool read()
Read model coefficients if they have changed.
Definition: kOmegaSSTDES.C:131
BasicTurbulenceModel::alphaField alphaField
Definition: kOmegaSSTDES.H:123