kOmegaSSTDDES.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 |
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 "kOmegaSSTDDES.H"
27 
28 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
29 
30 namespace Foam
31 {
32 namespace LESModels
33 {
34 
35 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
36 
37 template<class BasicTurbulenceModel>
39 (
40  const volScalarField& magGradU
41 ) const
42 {
44  (
45  min
46  (
47  this->nuEff()
48  /(
49  max
50  (
51  magGradU,
52  dimensionedScalar("SMALL", magGradU.dimensions(), SMALL)
53  )
54  *sqr(this->kappa_*this->y_)
55  ),
56  scalar(10)
57  )
58  );
59  tr().boundaryField() == 0.0;
60 
61  return tr;
62 }
63 
64 
65 template<class BasicTurbulenceModel>
67 (
68  const volScalarField& magGradU
69 ) const
70 {
71  return 1 - tanh(pow(Cd1_*rd(magGradU), Cd2_));
72 }
73 
74 
75 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
76 
77 template<class BasicTurbulenceModel>
79 (
80  const volScalarField& magGradU,
81  const volScalarField& CDES
82 ) const
83 {
84  const volScalarField& k = this->k_;
85  const volScalarField& omega = this->omega_;
86 
87  const volScalarField lRAS(sqrt(k)/(this->betaStar_*omega));
88  const volScalarField lLES(CDES*this->delta());
89 
90  return max
91  (
92  lRAS
93  - fd(magGradU)
94  *max
95  (
96  lRAS - lLES,
97  dimensionedScalar("zero", dimLength, 0)
98  ),
99  dimensionedScalar("small", dimLength, SMALL)
100  );
101 }
102 
103 
104 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
105 
106 template<class BasicTurbulenceModel>
108 (
109  const alphaField& alpha,
110  const rhoField& rho,
111  const volVectorField& U,
112  const surfaceScalarField& alphaRhoPhi,
113  const surfaceScalarField& phi,
114  const transportModel& transport,
115  const word& propertiesName,
116  const word& type
117 )
118 :
120  (
121  alpha,
122  rho,
123  U,
124  alphaRhoPhi,
125  phi,
126  transport,
127  propertiesName,
128  type
129  ),
130 
131  Cd1_
132  (
134  (
135  "Cd1",
136  this->coeffDict_,
137  20
138  )
139  ),
140  Cd2_
141  (
143  (
144  "Cd2",
145  this->coeffDict_,
146  3
147  )
148  )
149 {
150  if (type == typeName)
151  {
152  this->printCoeffs(type);
153  }
154 }
155 
156 
157 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
158 
159 template<class BasicTurbulenceModel>
161 {
163  {
164  Cd1_.readIfPresent(this->coeffDict());
165  Cd2_.readIfPresent(this->coeffDict());
166 
167  return true;
168  }
169  else
170  {
171  return false;
172  }
173 }
174 
175 
176 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
177 
178 } // End namespace LESModels
179 } // End namespace Foam
180 
181 // ************************************************************************* //
Foam::LESModels::kOmegaSSTDDES::kOmegaSSTDDES
kOmegaSSTDDES(const kOmegaSSTDDES &)
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
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::kOmegaSSTDDES::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: kOmegaSSTDDES.H:100
Foam::LESModels::kOmegaSSTDDES::dTilda
virtual tmp< volScalarField > dTilda(const volScalarField &magGradU, const volScalarField &CDES) const
Length scale.
Definition: kOmegaSSTDDES.C:79
U
U
Definition: pEqn.H:46
Foam::tanh
dimensionedScalar tanh(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:263
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::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::LESModels::kOmegaSSTDDES::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: kOmegaSSTDDES.H:99
rho
rho
Definition: pEqn.H:3
Foam::LESModels::kOmegaSSTDDES::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: kOmegaSSTDDES.H:101
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::LESModels::kOmegaSSTDDES::read
virtual bool read()
Re-read model coefficients if they have changed.
Definition: kOmegaSSTDDES.C:160
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
kOmegaSSTDDES.H
Foam::LESModels::kOmegaSSTDDES::fd
tmp< volScalarField > fd(const volScalarField &magGradU) const
Definition: kOmegaSSTDDES.C:67
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
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::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::LESModels::kOmegaSSTDDES::rd
tmp< volScalarField > rd(const volScalarField &magGradU) const
Definition: kOmegaSSTDDES.C:39