00001 /*---------------------------------------------------------------------------*\ 00002 ========= | 00003 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox 00004 \\ / O peration | 00005 \\ / A nd | Copyright (C) 1991-2010 OpenCFD Ltd. 00006 \\/ M anipulation | 00007 ------------------------------------------------------------------------------- 00008 License 00009 This file is part of OpenFOAM. 00010 00011 OpenFOAM is free software: you can redistribute it and/or modify it 00012 under the terms of the GNU General Public License as published by 00013 the Free Software Foundation, either version 3 of the License, or 00014 (at your option) any later version. 00015 00016 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT 00017 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or 00018 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License 00019 for more details. 00020 00021 You should have received a copy of the GNU General Public License 00022 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>. 00023 00024 Class 00025 Foam::compressible::LESModels::SpalartAllmaras 00026 00027 Description 00028 SpalartAllmaras for compressible flows 00029 00030 SourceFiles 00031 SpalartAllmaras.C 00032 00033 \*---------------------------------------------------------------------------*/ 00034 00035 #ifndef compressibleSpalartAllmaras_H 00036 #define compressibleSpalartAllmaras_H 00037 00038 #include <compressibleLESModels/LESModel.H> 00039 #include <finiteVolume/volFields.H> 00040 00041 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // 00042 00043 namespace Foam 00044 { 00045 namespace compressible 00046 { 00047 namespace LESModels 00048 { 00049 00050 /*---------------------------------------------------------------------------*\ 00051 Class SpalartAllmaras Declaration 00052 \*---------------------------------------------------------------------------*/ 00053 00054 class SpalartAllmaras 00055 : 00056 public LESModel 00057 { 00058 // Private data 00059 00060 // Model coefficients 00061 00062 dimensionedScalar sigmaNut_; 00063 dimensionedScalar Prt_; 00064 00065 dimensionedScalar Cb1_; 00066 dimensionedScalar Cb2_; 00067 dimensionedScalar Cv1_; 00068 dimensionedScalar Cv2_; 00069 dimensionedScalar CDES_; 00070 dimensionedScalar ck_; 00071 dimensionedScalar kappa_; 00072 dimensionedScalar Cw1_; 00073 dimensionedScalar Cw2_; 00074 dimensionedScalar Cw3_; 00075 00076 00077 // Fields 00078 00079 volScalarField nuTilda_; 00080 volScalarField dTilda_; 00081 volScalarField muSgs_; 00082 volScalarField alphaSgs_; 00083 00084 00085 // Private member functions 00086 00087 //- Update sub-grid scale fields 00088 void updateSubGridScaleFields(); 00089 00090 tmp<volScalarField> fv1() const; 00091 tmp<volScalarField> fv2() const; 00092 tmp<volScalarField> fv3() const; 00093 tmp<volScalarField> fw(const volScalarField& Stilda) const; 00094 00095 // Disallow default bitwise copy construct and assignment 00096 SpalartAllmaras(const SpalartAllmaras&); 00097 SpalartAllmaras& operator=(const SpalartAllmaras&); 00098 00099 00100 public: 00101 00102 //- Runtime type information 00103 TypeName("SpalartAllmaras"); 00104 00105 00106 // Constructors 00107 00108 //- Constructor from components 00109 SpalartAllmaras 00110 ( 00111 const volScalarField& rho, 00112 const volVectorField& U, 00113 const surfaceScalarField& phi, 00114 const basicThermo& thermoPhysicalModel 00115 ); 00116 00117 00118 //- Destructor 00119 virtual ~SpalartAllmaras() 00120 {} 00121 00122 00123 // Member Functions 00124 00125 tmp<volScalarField> nuTilda() const 00126 { 00127 return nuTilda_; 00128 } 00129 00130 //- Return SGS kinetic energy 00131 virtual tmp<volScalarField> k() const 00132 { 00133 return sqr(muSgs()/rho()/ck_/dTilda_); 00134 } 00135 00136 //- Return sub-grid disipation rate 00137 virtual tmp<volScalarField> epsilon() const; 00138 00139 //- Return SGS viscosity 00140 virtual tmp<volScalarField> muSgs() const 00141 { 00142 return muSgs_; 00143 } 00144 00145 //- Return SGS thermal diffusivity 00146 virtual tmp<volScalarField> alphaSgs() const 00147 { 00148 return alphaSgs_; 00149 } 00150 00151 //- Return thermal conductivity 00152 virtual tmp<volScalarField> alphaEff() const 00153 { 00154 return tmp<volScalarField> 00155 ( 00156 new volScalarField("alphaEff", alphaSgs_ + alpha()) 00157 ); 00158 } 00159 00160 //- Return the sub-grid stress tensor. 00161 virtual tmp<volSymmTensorField> B() const; 00162 00163 //- Return the deviatoric part of the effective sub-grid 00164 // turbulence stress tensor including the laminar stress 00165 virtual tmp<volSymmTensorField> devRhoBeff() const; 00166 00167 //- Returns div(rho*dev(B)). 00168 // This is the additional term due to the filtering of the NSE. 00169 virtual tmp<fvVectorMatrix> divDevRhoBeff(volVectorField& U) const; 00170 00171 //- Correct nuTilda and related properties 00172 virtual void correct(const tmp<volTensorField>& gradU); 00173 00174 //- Read LESProperties dictionary 00175 virtual bool read(); 00176 }; 00177 00178 00179 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // 00180 00181 } // End namespace LESModels 00182 } // End namespace compressible 00183 } // End namespace Foam 00184 00185 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // 00186 00187 #endif 00188 00189 // ************************ vim: set sw=4 sts=4 et: ************************ //