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-2023 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 #include "fvcLaplacian.H"
28 
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace RASModels
34 {
35 
36 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
37 
38 template<class BasicMomentumTransportModel>
40 (
41  const volScalarField::Internal& S2,
42  const volScalarField::Internal& gamma,
43  const volScalarField::Internal& beta
44 ) const
45 {
47  (
48  sqrt(this->k_())/(pow025(this->betaStar_)*this->omega_())
49  );
50 
52  (
53  max
54  (
55  kappa_*sqrt(S2)
56  /(
57  mag(fvc::laplacian(this->U_))()()
59  (
60  "rootVSmall",
61  dimensionSet(0, -1, -1, 0, 0),
62  rootVSmall
63  )
64  ),
65  Cs_*sqrt(kappa_*zeta2_/(beta/this->betaStar_ - gamma))*delta()()
66  )
67  );
68 
69  return fvm::Su
70  (
71  this->alpha_()*this->rho_()
72  *min
73  (
74  max
75  (
76  zeta2_*kappa_*S2*sqr(L/Lvk)
77  - (2*C_/sigmaPhi_)*this->k_()
78  *max
79  (
80  magSqr(fvc::grad(this->omega_)()())/sqr(this->omega_()),
81  magSqr(fvc::grad(this->k_)()())/sqr(this->k_())
82  ),
83  dimensionedScalar(dimensionSet(0, 0, -2, 0, 0), 0)
84  ),
85  // Limit SAS production of omega for numerical stability,
86  // particularly during start-up
87  this->omega_()/(0.1*this->omega_.time().deltaT())
88  ),
89  this->omega_
90  );
91 }
92 
93 
94 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
95 
96 template<class BasicMomentumTransportModel>
98 (
99  const alphaField& alpha,
100  const rhoField& rho,
101  const volVectorField& U,
102  const surfaceScalarField& alphaRhoPhi,
103  const surfaceScalarField& phi,
104  const viscosity& viscosity,
105  const word& type
106 )
107 :
108  kOmegaSST<BasicMomentumTransportModel>
109  (
110  alpha,
111  rho,
112  U,
113  alphaRhoPhi,
114  phi,
115  viscosity
116  ),
117 
118  Cs_
119  (
120  dimensioned<scalar>::lookupOrAddToDict
121  (
122  "Cs",
123  this->coeffDict_,
124  0.11
125  )
126  ),
127  kappa_
128  (
129  dimensioned<scalar>::lookupOrAddToDict
130  (
131  "kappa",
132  this->coeffDict_,
133  0.41
134  )
135  ),
136  zeta2_
137  (
138  dimensioned<scalar>::lookupOrAddToDict
139  (
140  "zeta2",
141  this->coeffDict_,
142  3.51
143  )
144  ),
145  sigmaPhi_
146  (
147  dimensioned<scalar>::lookupOrAddToDict
148  (
149  "sigmaPhi",
150  this->coeffDict_,
151  2.0/3.0
152  )
153  ),
154  C_
155  (
156  dimensioned<scalar>::lookupOrAddToDict
157  (
158  "C",
159  this->coeffDict_,
160  2
161  )
162  ),
163 
164  delta_
165  (
166  LESdelta::New
167  (
168  this->groupName("delta"),
169  *this,
170  this->coeffDict_
171  )
172  )
173 {}
174 
175 
176 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
177 
178 template<class BasicMomentumTransportModel>
180 {
182  {
183  Cs_.readIfPresent(this->coeffDict());
184  kappa_.readIfPresent(this->coeffDict());
185  sigmaPhi_.readIfPresent(this->coeffDict());
186  zeta2_.readIfPresent(this->coeffDict());
187  C_.readIfPresent(this->coeffDict());
188 
189  return true;
190  }
191  else
192  {
193  return false;
194  }
195 }
196 
197 
198 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
199 
200 } // End namespace RASModels
201 } // End namespace Foam
202 
203 // ************************************************************************* //
scalar delta
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricField class.
Abstract base class for LES deltas.
Definition: LESdelta.H:51
kOmegaSSTSAS(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity, const word &type=typeName)
Construct from components.
Definition: kOmegaSSTSAS.C:98
virtual tmp< fvScalarMatrix > Qsas(const volScalarField::Internal &S2, const volScalarField::Internal &gamma, const volScalarField::Internal &beta) const
SAS omega source.
Definition: kOmegaSSTSAS.C:40
virtual bool read()
Re-read model coefficients if they have changed.
Definition: kOmegaSSTSAS.C:179
Specialisation for RAS of the generic kOmegaSSTBase base class. For more information,...
Definition: kOmegaSST.H:64
Dimension set for the base types.
Definition: dimensionSet.H:125
BasicMomentumTransportModel::alphaField alphaField
BasicMomentumTransportModel::rhoField rhoField
A class for managing temporary objects.
Definition: tmp.H:55
Abstract base class for all fluid physical properties.
Definition: viscosity.H:50
A class for handling words, derived from string.
Definition: word.H:62
Calculate the laplacian of the given field.
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
tmp< VolField< Type > > laplacian(const VolField< Type > &vf, const word &name)
Definition: fvcLaplacian.C:45
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
tmp< fvMatrix< Type > > Su(const DimensionedField< Type, volMesh > &, const VolField< Type > &)
Namespace for OpenFOAM.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< scalar > mag(const dimensioned< Type > &)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
dimensioned< scalar > magSqr(const dimensioned< Type > &)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
dimensionedScalar pow025(const dimensionedScalar &ds)