kOmegaSSTSAS.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 "kOmegaSSTSAS.H"
27 
28 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
29 
30 namespace Foam
31 {
32 namespace RASModels
33 {
34 
35 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
36 
37 template<class BasicTurbulenceModel>
39 (
40  const volScalarField& S2,
41  const volScalarField& gamma,
42  const volScalarField& beta
43 ) const
44 {
46  (
47  sqrt(this->k_)/(pow025(this->betaStar_)*this->omega_)
48  );
49 
50  volScalarField Lvk
51  (
52  max
53  (
54  kappa_*sqrt(S2)
55  /(
56  mag(fvc::laplacian(this->U_))
58  (
59  "ROOTVSMALL",
60  dimensionSet(0, -1, -1, 0, 0),
61  ROOTVSMALL
62  )
63  ),
64  Cs_*sqrt(kappa_*zeta2_/(beta/this->betaStar_ - gamma))*delta()
65  )
66  );
67 
68  return fvm::Su
69  (
70  this->alpha_*this->rho_
71  *min
72  (
73  max
74  (
75  zeta2_*kappa_*S2*sqr(L/Lvk)
76  - (2*C_/sigmaPhi_)*this->k_
77  *max
78  (
79  magSqr(fvc::grad(this->omega_))/sqr(this->omega_),
80  magSqr(fvc::grad(this->k_))/sqr(this->k_)
81  ),
82  dimensionedScalar("0", dimensionSet(0, 0, -2, 0, 0), 0)
83  ),
84  // Limit SAS production of omega for numerical stability,
85  // particularly during start-up
86  this->omega_/(0.1*this->omega_.time().deltaT())
87  ),
88  this->omega_
89  );
90 }
91 
92 
93 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
94 
95 template<class BasicTurbulenceModel>
97 (
98  const alphaField& alpha,
99  const rhoField& rho,
100  const volVectorField& U,
101  const surfaceScalarField& alphaRhoPhi,
102  const surfaceScalarField& phi,
103  const transportModel& transport,
104  const word& propertiesName,
105  const word& type
106 )
107 :
109  (
110  alpha,
111  rho,
112  U,
113  alphaRhoPhi,
114  phi,
115  transport,
116  propertiesName,
117  type
118  ),
119 
120  Cs_
121  (
123  (
124  "Cs",
125  this->coeffDict_,
126  0.11
127  )
128  ),
129  kappa_
130  (
132  (
133  "kappa",
134  this->coeffDict_,
135  0.41
136  )
137  ),
138  zeta2_
139  (
141  (
142  "zeta2",
143  this->coeffDict_,
144  3.51
145  )
146  ),
147  sigmaPhi_
148  (
150  (
151  "sigmaPhi",
152  this->coeffDict_,
153  2.0/3.0
154  )
155  ),
156  C_
157  (
159  (
160  "C",
161  this->coeffDict_,
162  2
163  )
164  ),
165 
166  delta_
167  (
169  (
170  IOobject::groupName("delta", U.group()),
171  *this,
172  this->coeffDict_
173  )
174  )
175 {
176  if (type == typeName)
177  {
178  this->correctNut();
179  this->printCoeffs(type);
180  }
181 }
182 
183 
184 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
185 
186 template<class BasicTurbulenceModel>
188 {
190  {
191  Cs_.readIfPresent(this->coeffDict());
192  kappa_.readIfPresent(this->coeffDict());
193  sigmaPhi_.readIfPresent(this->coeffDict());
194  zeta2_.readIfPresent(this->coeffDict());
195  C_.readIfPresent(this->coeffDict());
196 
197  return true;
198  }
199  else
200  {
201  return false;
202  }
203 }
204 
205 
206 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
207 
208 } // End namespace RASModels
209 } // End namespace Foam
210 
211 // ************************************************************************* //
beta
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
Foam::IOobject::groupName
static word groupName(Name name, const word &group)
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::RASModels::kOmegaSST
Implementation of the k-omega-SST turbulence model for incompressible and compressible flows.
Definition: kOmegaSST.H:126
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::RASModels::kOmegaSSTSAS::Qsas
virtual tmp< fvScalarMatrix > Qsas(const volScalarField &S2, const volScalarField &gamma, const volScalarField &beta) const
SAS omega source.
Definition: kOmegaSSTSAS.C:39
Foam::dimensionSet
Dimension set for the base types.
Definition: dimensionSet.H:116
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::fvm::Su
tmp< fvMatrix< Type > > Su(const DimensionedField< Type, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
Foam::RASModels::kOmegaSSTSAS::read
virtual bool read()
Re-read model coefficients if they have changed.
Definition: kOmegaSSTSAS.C:187
Foam::RASModels::kOmegaSST::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: kOmegaSST.H:141
U
U
Definition: pEqn.H:46
Foam::pow025
dimensionedScalar pow025(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:131
Foam::fvc::laplacian
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
Foam::RASModels::kOmegaSST::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: kOmegaSST.H:142
delta
scalar delta
Definition: LISASMDCalcMethod2.H:8
Foam::RASModels::kOmegaSSTSAS::kOmegaSSTSAS
kOmegaSSTSAS(const kOmegaSSTSAS &)
kOmegaSSTSAS.H
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::RASModels::kOmegaSST::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: kOmegaSST.H:143
rho
rho
Definition: pEqn.H:3
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
Foam::LESdelta::New
static autoPtr< LESdelta > New(const word &name, const turbulenceModel &turbulence, const dictionary &dict)
Return a reference to the selected LES delta.
Definition: LESdelta.C:66
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
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::magSqr
dimensioned< scalar > magSqr(const dimensioned< Type > &)