kOmegaSSTSAS.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2015-2020 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 BasicMomentumTransportModel>
39 (
40  const volScalarField::Internal& S2,
41  const volScalarField::Internal& gamma,
42  const volScalarField::Internal& beta
43 ) const
44 {
46  (
47  sqrt(this->k_())/(pow025(this->betaStar_)*this->omega_())
48  );
49 
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(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 BasicMomentumTransportModel>
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& type
105 )
106 :
108  (
109  alpha,
110  rho,
111  U,
112  alphaRhoPhi,
113  phi,
114  transport
115  ),
116 
117  Cs_
118  (
120  (
121  "Cs",
122  this->coeffDict_,
123  0.11
124  )
125  ),
126  kappa_
127  (
129  (
130  "kappa",
131  this->coeffDict_,
132  0.41
133  )
134  ),
135  zeta2_
136  (
138  (
139  "zeta2",
140  this->coeffDict_,
141  3.51
142  )
143  ),
144  sigmaPhi_
145  (
147  (
148  "sigmaPhi",
149  this->coeffDict_,
150  2.0/3.0
151  )
152  ),
153  C_
154  (
156  (
157  "C",
158  this->coeffDict_,
159  2
160  )
161  ),
162 
163  delta_
164  (
166  (
167  IOobject::groupName("delta", alphaRhoPhi.group()),
168  *this,
169  this->coeffDict_
170  )
171  )
172 {}
173 
174 
175 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
176 
177 template<class BasicMomentumTransportModel>
179 {
181  {
182  Cs_.readIfPresent(this->coeffDict());
183  kappa_.readIfPresent(this->coeffDict());
184  sigmaPhi_.readIfPresent(this->coeffDict());
185  zeta2_.readIfPresent(this->coeffDict());
186  C_.readIfPresent(this->coeffDict());
187 
188  return true;
189  }
190  else
191  {
192  return false;
193  }
194 }
195 
196 
197 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
198 
199 } // End namespace RASModels
200 } // End namespace Foam
201 
202 // ************************************************************************* //
scalar delta
BasicMomentumTransportModel::alphaField alphaField
Definition: kOmegaSST.H:68
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
static autoPtr< LESdelta > New(const word &name, const momentumTransportModel &turbulence, const dictionary &dict)
Return a reference to the selected LES delta.
Definition: LESdelta.C:66
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
tmp< fvMatrix< Type > > Su(const DimensionedField< Type, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensionedScalar pow025(const dimensionedScalar &ds)
BasicMomentumTransportModel::transportModel transportModel
Definition: kOmegaSST.H:70
DimensionedField< scalar, volMesh > Internal
Type of the internal field from which this GeometricField is derived.
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
phi
Definition: pEqn.H:104
Dimension set for the base types.
Definition: dimensionSet.H:120
A class for handling words, derived from string.
Definition: word.H:59
static word groupName(Name name, const word &group)
Specialisation for RAS of the generic kOmegaSSTBase base class. For more information, see Description of kOmegaSSTBase.H.
Definition: kOmegaSST.H:57
kOmegaSSTSAS(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &type=typeName)
Construct from components.
Definition: kOmegaSSTSAS.C:97
dimensioned< scalar > magSqr(const dimensioned< Type > &)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
virtual bool read()
Re-read model coefficients if they have changed.
Definition: kOmegaSSTSAS.C:178
virtual tmp< fvScalarMatrix > Qsas(const volScalarField::Internal &S2, const volScalarField::Internal &gamma, const volScalarField::Internal &beta) const
SAS omega source.
Definition: kOmegaSSTSAS.C:39
U
Definition: pEqn.H:72
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensioned< scalar > mag(const dimensioned< Type > &)
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
Namespace for OpenFOAM.
BasicMomentumTransportModel::rhoField rhoField
Definition: kOmegaSST.H:69