![]() |
|
00001 /*---------------------------------------------------------------------------*\ 00002 ========= | 00003 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox 00004 \\ / O peration | 00005 \\ / A nd | Copyright (C) 1991-2005 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 the 00013 Free Software Foundation; either version 2 of the License, or (at your 00014 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, write to the Free Software Foundation, 00023 Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA 00024 00025 Class 00026 DeardorffDiffStress 00027 00028 Description 00029 <pre> 00030 Differential SGS Stress Equation Model 00031 ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ 00032 The DSEM uses a model version of the full balance equation for the SGS 00033 stress tensor to simulate the behaviour of B, hence, 00034 00035 d/dt(B) + div(U*B) - div(nuSgs*grad(B)) 00036 = 00037 P - c1*epsilon/k*B - 0.667*(1 - c1)*epsilon*I - c2*(P - 0.333*trP*I) 00038 00039 where 00040 00041 k = 0.5*tr(B), 00042 epsilon = ce*k^3/2/delta, 00043 epsilon/k = ce*k^1/2/delta 00044 P = -(B'L + L'B) 00045 nuSgs = ck*sqrt(k)*delta 00046 nuEff = nuSgs + nu 00047 </pre> 00048 00049 SourceFiles 00050 DeardorffDiffStress.C 00051 00052 \*---------------------------------------------------------------------------*/ 00053 00054 #ifndef DeardorffDiffStress_H 00055 #define DeardorffDiffStress_H 00056 00057 #include "GenSGSStress.H" 00058 00059 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // 00060 00061 namespace Foam 00062 { 00063 namespace LESmodels 00064 { 00065 00066 /*---------------------------------------------------------------------------*\ 00067 Class DeardorffDiffStress Declaration 00068 \*---------------------------------------------------------------------------*/ 00069 00070 class DeardorffDiffStress 00071 : 00072 public GenSGSStress 00073 { 00074 // Private data 00075 00076 dimensionedScalar ck_; 00077 dimensionedScalar cm_; 00078 00079 00080 // Private Member Functions 00081 00082 // Disallow default bitwise copy construct and assignment 00083 DeardorffDiffStress(const DeardorffDiffStress&); 00084 DeardorffDiffStress& operator=(const DeardorffDiffStress&); 00085 00086 00087 public: 00088 00089 //- Runtime type information 00090 TypeName("DeardorffDiffStress"); 00091 00092 // Constructors 00093 00094 //- Constructor from components 00095 DeardorffDiffStress 00096 ( 00097 const volVectorField& U, 00098 const surfaceScalarField& phi, 00099 transportModel& transport 00100 ); 00101 00102 00103 // Destructor 00104 00105 ~DeardorffDiffStress() 00106 {} 00107 00108 00109 // Member Functions 00110 00111 //- Return the effective diffusivity for B 00112 tmp<volScalarField> DBEff() const 00113 { 00114 return tmp<volScalarField> 00115 ( 00116 new volScalarField("DBEff", nuSgs_ + nu()) 00117 ); 00118 } 00119 00120 //- Correct Eddy-Viscosity and related properties 00121 void correct(const tmp<volTensorField>& gradU); 00122 00123 //- Read turbulenceProperties dictionary 00124 bool read(); 00125 }; 00126 00127 00128 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // 00129 00130 } // End namespace LESmodels 00131 } // End namespace Foam 00132 00133 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // 00134 00135 #endif 00136 00137 // ************************************************************************* //