kOmegaSSTIDDES.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) 2015 OpenFOAM Foundation
6  \\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd.
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 "kOmegaSSTIDDES.H"
27 
28 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
29 
30 namespace Foam
31 {
32 namespace LESModels
33 {
34 
35 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
36 
37 template<class BasicTurbulenceModel>
39 {
40  if (!isA<IDDESDelta>(this->delta_()))
41  {
43  << "The delta function must be set to a " << IDDESDelta::typeName
44  << " -based model" << exit(FatalError);
45  }
46 
47  return refCast<const IDDESDelta>(this->delta_());
48 }
49 
50 
51 template<class BasicTurbulenceModel>
53 {
54  return max(0.25 - this->y_/IDDESDelta_.hmax(), scalar(-5));
55 }
56 
57 
58 template<class BasicTurbulenceModel>
60 (
61  const volScalarField& magGradU
62 ) const
63 {
64  return tanh(pow3(sqr(Ct_)*rd(this->nut_, magGradU)));
65 }
66 
67 
68 template<class BasicTurbulenceModel>
70 (
71  const volScalarField& magGradU
72 ) const
73 {
74  return tanh(pow(sqr(Cl_)*rd(this->nu(), magGradU), 10));
75 }
76 
77 
78 template<class BasicTurbulenceModel>
80 (
81  const volScalarField& nur,
82  const volScalarField& magGradU
83 ) const
84 {
86  (
87  min
88  (
89  nur
90  /(
91  max
92  (
93  magGradU,
94  dimensionedScalar("SMALL", magGradU.dimensions(), SMALL)
95  )
96  *sqr(this->kappa_*this->y_)
97  ),
98  scalar(10)
99  )
100  );
101  tr().boundaryField() == 0.0;
102 
103  return tr;
104 }
105 
106 
107 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
108 
109 template<class BasicTurbulenceModel>
111 (
112  const volScalarField& magGradU
113 ) const
114 {
115  return 1 - tanh(pow(Cdt1_*rd(this->nut_, magGradU), Cdt2_));
116 }
117 
118 
119 template<class BasicTurbulenceModel>
121 (
122  const volScalarField& magGradU,
123  const volScalarField& CDES
124 ) const
125 {
126  const volScalarField& k = this->k_;
127  const volScalarField& omega = this->omega_;
128 
129  const volScalarField lRAS(sqrt(k)/(this->betaStar_*omega));
130  const volScalarField lLES(CDES*this->delta());
131 
132  const volScalarField alpha(this->alpha());
133  const volScalarField expTerm(exp(sqr(alpha)));
134 
135  tmp<volScalarField> fB = min(2*pow(expTerm, -9.0), scalar(1));
136  tmp<volScalarField> fe1 =
137  2*(pos(alpha)*pow(expTerm, -11.09) + neg(alpha)*pow(expTerm, -9.0));
138  tmp<volScalarField> fe2 = 1 - max(ft(magGradU), fl(magGradU));
139  tmp<volScalarField> fe = max(fe1 - 1, scalar(0))*fe2;
140 
141  const volScalarField fdTilda(max(1 - fdt(magGradU), fB));
142 
143  // Simplified formulation from Gritskevich et al. paper (2011) where fe = 0
144  // return max
145  // (
146  // fdTilda*lRAS + (1 - fdTilda)*lLES,
147  // dimensionedScalar("SMALL", dimLength, SMALL)
148  // );
149 
150  // Original formulation from Shur et al. paper (2008)
151  return max
152  (
153  fdTilda*(1 + fe)*lRAS + (1 - fdTilda)*lLES,
154  dimensionedScalar("SMALL", dimLength, SMALL)
155  );
156 }
157 
158 
159 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
160 
161 template<class BasicTurbulenceModel>
163 (
164  const alphaField& alpha,
165  const rhoField& rho,
166  const volVectorField& U,
167  const surfaceScalarField& alphaRhoPhi,
168  const surfaceScalarField& phi,
169  const transportModel& transport,
170  const word& propertiesName,
171  const word& type
172 )
173 :
175  (
176  alpha,
177  rho,
178  U,
179  alphaRhoPhi,
180  phi,
181  transport,
182  propertiesName,
183  type
184  ),
185 
186  Cdt1_
187  (
189  (
190  "Cdt1",
191  this->coeffDict_,
192  20
193  )
194  ),
195  Cdt2_
196  (
198  (
199  "Cdt2",
200  this->coeffDict_,
201  3
202  )
203  ),
204  Cl_
205  (
207  (
208  "Cl",
209  this->coeffDict_,
210  5
211  )
212  ),
213  Ct_
214  (
216  (
217  "Ct",
218  this->coeffDict_,
219  1.87
220  )
221  ),
222  IDDESDelta_(setDelta())
223 {
224  if (type == typeName)
225  {
226  this->printCoeffs(type);
227  }
228 }
229 
230 
231 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
232 
233 template<class BasicTurbulenceModel>
235 {
237  {
238  Cdt1_.readIfPresent(this->coeffDict());
239  Cdt2_.readIfPresent(this->coeffDict());
240  Cl_.readIfPresent(this->coeffDict());
241  Ct_.readIfPresent(this->coeffDict());
242 
243  return true;
244  }
245  else
246  {
247  return false;
248  }
249 }
250 
251 
252 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
253 
254 } // End namespace LESModels
255 } // End namespace Foam
256 
257 // ************************************************************************* //
Foam::LESModels::kOmegaSSTIDDES::fl
tmp< volScalarField > fl(const volScalarField &magGradU) const
Definition: kOmegaSSTIDDES.C:70
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::LESModels::kOmegaSSTDES
k-omega-SST DES turbulence model for incompressible and compressible flows
Definition: kOmegaSSTDES.H:65
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::LESModels::kOmegaSSTIDDES::fdt
tmp< volScalarField > fdt(const volScalarField &magGradU) const
Delay function.
Definition: kOmegaSSTIDDES.C:111
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:216
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:252
Foam::LESModels::IDDESDelta
Definition: IDDESDelta.H:52
Foam::LESModels::kOmegaSSTIDDES::rd
tmp< volScalarField > rd(const volScalarField &nur, const volScalarField &magGradU) const
Definition: kOmegaSSTIDDES.C:80
Foam::LESModels::kOmegaSSTIDDES::read
virtual bool read()
Re-read model coefficients if they have changed.
Definition: kOmegaSSTIDDES.C:234
U
U
Definition: pEqn.H:46
nu
volScalarField & nu
Definition: readMechanicalProperties.H:179
Foam::tanh
dimensionedScalar tanh(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:263
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:87
delta
scalar delta
Definition: LISASMDCalcMethod2.H:8
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:73
Foam::LESModels::kOmegaSSTIDDES::ft
tmp< volScalarField > ft(const volScalarField &magGradU) const
Definition: kOmegaSSTIDDES.C:60
Foam::FatalError
error FatalError
kOmegaSSTIDDES.H
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::LESModels::kOmegaSSTDES::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: kOmegaSSTDES.H:110
rho
rho
Definition: pEqn.H:3
Foam::LESModels::kOmegaSSTIDDES::dTilda
virtual tmp< volScalarField > dTilda(const volScalarField &magGradU, const volScalarField &CDES) const
Length scale.
Definition: kOmegaSSTIDDES.C:121
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
Foam::LESModels::kOmegaSSTDES::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: kOmegaSSTDES.H:108
Foam::LESModels::kOmegaSSTIDDES::alpha
tmp< volScalarField > alpha() const
Definition: kOmegaSSTIDDES.C:52
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::LESModels::kOmegaSSTIDDES::setDelta
const IDDESDelta & setDelta() const
Check that the supplied delta is an IDDESDelta.
Definition: kOmegaSSTIDDES.C:38
Foam::tr
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
Definition: dimensionedSphericalTensor.C:49
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::LESModels::kOmegaSSTIDDES::kOmegaSSTIDDES
kOmegaSSTIDDES(const kOmegaSSTIDDES &)
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::LESModels::kOmegaSSTDES::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: kOmegaSSTDES.H:109
Foam::neg
dimensionedScalar neg(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:201
Foam::pos
dimensionedScalar pos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:190