FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

vanDriestDelta.C

Go to the documentation of this file.
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 \*---------------------------------------------------------------------------*/
00025 
00026 #include "vanDriestDelta.H"
00027 #include <compressibleLESModels/LESModel.H>
00028 #include <finiteVolume/wallFvPatch.H>
00029 #include <finiteVolume/wallDistData.H>
00030 #include <finiteVolume/wallPointYPlus.H>
00031 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00032 
00033 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00034 
00035 namespace Foam
00036 {
00037 namespace compressible
00038 {
00039 namespace LESModels
00040 {
00041 
00042 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
00043 
00044 defineTypeNameAndDebug(vanDriestDelta, 0);
00045 addToRunTimeSelectionTable(LESdelta, vanDriestDelta, dictionary);
00046 
00047 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00048 
00049 void vanDriestDelta::calcDelta()
00050 {
00051     const LESModel& lesModel = mesh_.lookupObject<LESModel>("LESProperties");
00052 
00053     const volVectorField& U = lesModel.U();
00054     const volScalarField& rho = lesModel.rho();
00055     const volScalarField& mu = lesModel.mu();
00056     tmp<volScalarField> muSgs = lesModel.muSgs();
00057 
00058     volScalarField ystar
00059     (
00060         IOobject
00061         (
00062             "ystar",
00063             mesh_.time().constant(),
00064             mesh_
00065         ),
00066         mesh_,
00067         dimensionedScalar("ystar", dimLength, GREAT)
00068     );
00069 
00070     const fvPatchList& patches = mesh_.boundary();
00071     forAll(patches, patchi)
00072     {
00073         if (isA<wallFvPatch>(patches[patchi]))
00074         {
00075             const fvPatchVectorField& Uw = U.boundaryField()[patchi];
00076             const scalarField& rhow = rho.boundaryField()[patchi];
00077             const scalarField& muw = mu.boundaryField()[patchi];
00078             const scalarField& muSgsw = muSgs().boundaryField()[patchi];
00079 
00080             ystar.boundaryField()[patchi] =
00081                 muw/(rhow*sqrt((muw + muSgsw)*mag(Uw.snGrad())/rhow + VSMALL));
00082         }
00083     }
00084 
00085     wallPointYPlus::yPlusCutOff = 500;
00086     wallDistData<wallPointYPlus> y(mesh_, ystar);
00087 
00088     delta_ = min
00089     (
00090         static_cast<const volScalarField&>(geometricDelta_()),
00091         (kappa_/Cdelta_)*((scalar(1) + SMALL) - exp(-y/ystar/Aplus_))*y
00092     );
00093 }
00094 
00095 
00096 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00097 
00098 vanDriestDelta::vanDriestDelta
00099 (
00100     const word& name,
00101     const fvMesh& mesh,
00102     const dictionary& dd
00103 )
00104 :
00105     LESdelta(name, mesh),
00106     geometricDelta_
00107     (
00108         LESdelta::New("geometricDelta", mesh, dd.subDict(type() + "Coeffs"))
00109     ),
00110     kappa_(dd.lookupOrDefault<scalar>("kappa", 0.41)),
00111     Aplus_
00112     (
00113         dd.subDict(type() + "Coeffs").lookupOrDefault<scalar>("Aplus", 26.0)
00114     ),
00115     Cdelta_
00116     (
00117         dd.subDict(type() + "Coeffs").lookupOrDefault<scalar>("Cdelta", 0.158)
00118     ),
00119     calcInterval_
00120     (
00121         dd.subDict(type() + "Coeffs").lookupOrDefault<label>("calcInterval", 1)
00122     )
00123 {
00124     delta_ = geometricDelta_();
00125 }
00126 
00127 
00128 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00129 
00130 void vanDriestDelta::read(const dictionary& d)
00131 {
00132     const dictionary& dd(d.subDict(type() + "Coeffs"));
00133 
00134     geometricDelta_().read(dd);
00135     d.readIfPresent<scalar>("kappa", kappa_);
00136     dd.readIfPresent<scalar>("Aplus", Aplus_);
00137     dd.readIfPresent<scalar>("Cdelta", Cdelta_);
00138     dd.readIfPresent<label>("calcInterval", calcInterval_);
00139     calcDelta();
00140 }
00141 
00142 
00143 void vanDriestDelta::correct()
00144 {
00145     if (mesh().time().timeIndex() % calcInterval_ == 0)
00146     {
00147         geometricDelta_().correct();
00148         calcDelta();
00149     }
00150 }
00151 
00152 
00153 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00154 
00155 } // End namespace LESModels
00156 } // End namespace compressible
00157 } // End namespace Foam
00158 
00159 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines