SpalartAllmarasIDDES.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) 2011-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 "SpalartAllmarasIDDES.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& chi,
123  const volScalarField& fv1,
124  const volTensorField& gradU
125 ) const
126 {
127  const volScalarField magGradU(mag(gradU));
128  const volScalarField psi(this->psi(chi, fv1));
129 
130  const volScalarField& lRAS(this->y_);
131  const volScalarField lLES(psi*this->CDES_*this->delta());
132 
133  const volScalarField alpha(this->alpha());
134  const volScalarField expTerm(exp(sqr(alpha)));
135 
136  tmp<volScalarField> fB = min(2*pow(expTerm, -9.0), scalar(1));
137  tmp<volScalarField> fe1 =
138  2*(pos(alpha)*pow(expTerm, -11.09) + neg(alpha)*pow(expTerm, -9.0));
139  tmp<volScalarField> fe2 = 1 - max(ft(magGradU), fl(magGradU));
140  tmp<volScalarField> fe = max(fe1 - 1, scalar(0))*psi*fe2;
141 
142  const volScalarField fdTilda(max(1 - fdt(magGradU), fB));
143 
144  // Simplified formulation from Gritskevich et al. paper (2011) where fe = 0
145  // return max
146  // (
147  // fdTilda*lRAS + (1 - fdTilda)*lLES,
148  // dimensionedScalar("SMALL", dimLength, SMALL)
149  // );
150 
151  // Original formulation from Shur et al. paper (2008)
152  return max
153  (
154  fdTilda*(1 + fe)*lRAS + (1 - fdTilda)*lLES,
155  dimensionedScalar("SMALL", dimLength, SMALL)
156  );
157 }
158 
159 
160 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
161 
162 template<class BasicTurbulenceModel>
164 (
165  const alphaField& alpha,
166  const rhoField& rho,
167  const volVectorField& U,
168  const surfaceScalarField& alphaRhoPhi,
169  const surfaceScalarField& phi,
170  const transportModel& transport,
171  const word& propertiesName,
172  const word& type
173 )
174 :
176  (
177  alpha,
178  rho,
179  U,
180  alphaRhoPhi,
181  phi,
182  transport,
183  propertiesName,
184  type
185  ),
186 
187  Cdt1_
188  (
190  (
191  "Cdt1",
192  this->coeffDict_,
193  8
194  )
195  ),
196  Cdt2_
197  (
199  (
200  "Cdt2",
201  this->coeffDict_,
202  3
203  )
204  ),
205  Cl_
206  (
208  (
209  "Cl",
210  this->coeffDict_,
211  3.55
212  )
213  ),
214  Ct_
215  (
217  (
218  "Ct",
219  this->coeffDict_,
220  1.63
221  )
222  ),
223  IDDESDelta_(setDelta())
224 {
225  if (type == typeName)
226  {
227  this->printCoeffs(type);
228  }
229 }
230 
231 
232 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
233 
234 template<class BasicTurbulenceModel>
236 {
238  {
239  Cdt1_.readIfPresent(this->coeffDict());
240  Cdt2_.readIfPresent(this->coeffDict());
241  Cl_.readIfPresent(this->coeffDict());
242  Ct_.readIfPresent(this->coeffDict());
243 
244  return true;
245  }
246  else
247  {
248  return false;
249  }
250 }
251 
252 
253 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
254 
255 } // End namespace LESModels
256 } // End namespace Foam
257 
258 // ************************************************************************* //
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
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::SpalartAllmarasIDDES::dTilda
virtual tmp< volScalarField > dTilda(const volScalarField &chi, const volScalarField &fv1, const volTensorField &gradU) const
Length scale.
Definition: SpalartAllmarasIDDES.C:121
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::LESModels::SpalartAllmarasIDDES::alpha
tmp< volScalarField > alpha() const
Definition: SpalartAllmarasIDDES.C:52
Foam::LESModels::SpalartAllmarasDES
SpalartAllmarasDES DES turbulence model for incompressible and compressible flows.
Definition: SpalartAllmarasDES.H:72
Foam::LESModels::DESModel::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: DESModel.H:71
Foam::LESModels::SpalartAllmarasIDDES::setDelta
const IDDESDelta & setDelta() const
Check that the supplied delta is an IDDESDelta.
Definition: SpalartAllmarasIDDES.C:38
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:252
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::LESModels::IDDESDelta
Definition: IDDESDelta.H:52
U
U
Definition: pEqn.H:46
Foam::LESModels::SpalartAllmarasIDDES::read
virtual bool read()
Re-read model coefficients if they have changed.
Definition: SpalartAllmarasIDDES.C:235
nu
volScalarField & nu
Definition: readMechanicalProperties.H:179
Foam::LESModels::SpalartAllmarasIDDES::fl
tmp< volScalarField > fl(const volScalarField &magGradU) const
Definition: SpalartAllmarasIDDES.C:70
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::FatalError
error FatalError
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam::LESModels::DESModel::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: DESModel.H:72
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
rho
rho
Definition: pEqn.H:3
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
psi
const volScalarField & psi
Definition: setRegionFluidFields.H:13
Foam::LESModels::SpalartAllmarasIDDES::SpalartAllmarasIDDES
SpalartAllmarasIDDES(const SpalartAllmarasIDDES &)
Foam::LESModels::SpalartAllmarasIDDES::fdt
tmp< volScalarField > fdt(const volScalarField &magGradU) const
Delay function.
Definition: SpalartAllmarasIDDES.C:111
Foam::tr
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
Definition: dimensionedSphericalTensor.C:49
Foam::LESModels::DESModel::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: DESModel.H:70
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::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::neg
dimensionedScalar neg(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:201
Foam::LESModels::SpalartAllmarasIDDES::ft
tmp< volScalarField > ft(const volScalarField &magGradU) const
Definition: SpalartAllmarasIDDES.C:60
SpalartAllmarasIDDES.H
Foam::LESModels::SpalartAllmarasIDDES::rd
tmp< volScalarField > rd(const volScalarField &nur, const volScalarField &magGradU) const
Definition: SpalartAllmarasIDDES.C:80
Foam::pos
dimensionedScalar pos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:190