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 <incompressibleLESModels/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 incompressible
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& nu = lesModel.nu();
00055     tmp<volScalarField> nuSgs = lesModel.nuSgs();
00056 
00057     volScalarField ystar
00058     (
00059         IOobject
00060         (
00061             "ystar",
00062             mesh_.time().constant(),
00063             mesh_
00064         ),
00065         mesh_,
00066         dimensionedScalar("ystar", dimLength, GREAT)
00067     );
00068 
00069     const fvPatchList& patches = mesh_.boundary();
00070     forAll(patches, patchi)
00071     {
00072         if (isA<wallFvPatch>(patches[patchi]))
00073         {
00074             const fvPatchVectorField& Uw = U.boundaryField()[patchi];
00075             const scalarField& nuw = nu.boundaryField()[patchi];
00076             const scalarField& nuSgsw = nuSgs().boundaryField()[patchi];
00077 
00078             ystar.boundaryField()[patchi] =
00079                 nuw/sqrt((nuw + nuSgsw)*mag(Uw.snGrad()) + VSMALL);
00080         }
00081     }
00082 
00083     wallPointYPlus::yPlusCutOff = 500;
00084     wallDistData<wallPointYPlus> y(mesh_, ystar);
00085 
00086     delta_ = min
00087     (
00088         static_cast<const volScalarField&>(geometricDelta_()),
00089         (kappa_/Cdelta_)*((scalar(1) + SMALL) - exp(-y/ystar/Aplus_))*y
00090     );
00091 }
00092 
00093 
00094 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00095 
00096 vanDriestDelta::vanDriestDelta
00097 (
00098     const word& name,
00099     const fvMesh& mesh,
00100     const dictionary& dd
00101 )
00102 :
00103     LESdelta(name, mesh),
00104     geometricDelta_
00105     (
00106         LESdelta::New("geometricDelta", mesh, dd.subDict(type() + "Coeffs"))
00107     ),
00108     kappa_(dd.lookupOrDefault<scalar>("kappa", 0.41)),
00109     Aplus_
00110     (
00111         dd.subDict(type() + "Coeffs").lookupOrDefault<scalar>("Aplus", 26.0)
00112     ),
00113     Cdelta_
00114     (
00115         dd.subDict(type() + "Coeffs").lookupOrDefault<scalar>("Cdelta", 0.158)
00116     ),
00117     calcInterval_
00118     (
00119         dd.subDict(type() + "Coeffs").lookupOrDefault<label>("calcInterval", 1)
00120     )
00121 {
00122     delta_ = geometricDelta_();
00123 }
00124 
00125 
00126 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00127 
00128 void vanDriestDelta::read(const dictionary& d)
00129 {
00130     const dictionary& dd(d.subDict(type() + "Coeffs"));
00131 
00132     geometricDelta_().read(dd);
00133     d.readIfPresent<scalar>("kappa", kappa_);
00134     dd.readIfPresent<scalar>("Aplus", Aplus_);
00135     dd.readIfPresent<scalar>("Cdelta", Cdelta_);
00136     dd.readIfPresent<label>("calcInterval", calcInterval_);
00137     calcDelta();
00138 }
00139 
00140 
00141 void vanDriestDelta::correct()
00142 {
00143     if (mesh().time().timeIndex() % calcInterval_ == 0)
00144     {
00145         geometricDelta_().correct();
00146         calcDelta();
00147     }
00148 }
00149 
00150 
00151 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00152 
00153 } // End namespace LESModels
00154 } // End namespace incompressible
00155 } // End namespace Foam
00156 
00157 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines