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

gradientDispersionRAS.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 <OpenFOAM/error.H>
00027 
00028 #include "gradientDispersionRAS.H"
00029 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00030 
00031 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00032 
00033 namespace Foam
00034 {
00035 
00036 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
00037 
00038 defineTypeNameAndDebug(gradientDispersionRAS, 0);
00039 
00040 addToRunTimeSelectionTable
00041 (
00042     dispersionModel,
00043     gradientDispersionRAS,
00044     dictionary
00045 );
00046 
00047 
00048 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00049 
00050 // Construct from components
00051 gradientDispersionRAS::gradientDispersionRAS
00052 (
00053     const dictionary& dict,
00054     spray& sm
00055 )
00056 :
00057     dispersionRASModel(dict, sm)
00058 {}
00059 
00060 
00061 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
00062 
00063 gradientDispersionRAS::~gradientDispersionRAS()
00064 {}
00065 
00066 
00067 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00068 
00069 void gradientDispersionRAS::disperseParcels() const
00070 {
00071 
00072     const scalar cps = 0.16432;
00073 
00074     scalar dt = spray_.runTime().deltaT().value();
00075     const volScalarField& k = turbulence().k();
00076     volVectorField gradk = fvc::grad(k);
00077     const volScalarField& epsilon = turbulence().epsilon();
00078     const volVectorField& U = spray_.U();
00079 
00080     for
00081     (
00082         spray::iterator elmnt = spray_.begin();
00083         elmnt != spray_.end();
00084         ++elmnt
00085     )
00086     {
00087         label celli = elmnt().cell();
00088         scalar UrelMag = mag(elmnt().U() - U[celli] - elmnt().Uturb());
00089 
00090         scalar Tturb = min
00091         (
00092             k[celli]/epsilon[celli], 
00093             cps*pow(k[celli], 1.5)/epsilon[celli]/(UrelMag + SMALL)
00094         );
00095         // parcel is perturbed by the turbulence
00096         if (dt < Tturb)
00097         {
00098             elmnt().tTurb() += dt;
00099 
00100             if (elmnt().tTurb() > Tturb)
00101             {
00102                 elmnt().tTurb() = 0.0;
00103 
00104                 scalar sigma = sqrt(2.0*k[celli]/3.0);
00105                 vector dir = -gradk[celli]/(mag(gradk[celli]) + SMALL);
00106 
00107                 // numerical recipes... Ch. 7. Random Numbers...
00108                 scalar x1 = 0.0;
00109                 scalar x2 = 0.0;
00110                 scalar rsq = 10.0;
00111                 while((rsq > 1.0) || (rsq == 0.0))
00112                 {
00113                     x1 = 2.0*spray_.rndGen().scalar01() - 1.0;
00114                     x2 = 2.0*spray_.rndGen().scalar01() - 1.0;
00115                     rsq = x1*x1 + x2*x2;
00116                 }
00117                 
00118                 scalar fac = sqrt(-2.0*log(rsq)/rsq);
00119                 
00120                 // in 2D calculations the -grad(k) is always
00121                 // away from the axis of symmetry
00122                 // This creates a 'hole' in the spray and to
00123                 // prevent this we let x1 be both negative/positive
00124                 if (spray_.twoD())
00125                 {
00126                     fac *= x1;
00127                 }
00128                 else
00129                 {
00130                     fac *= mag(x1);
00131                 }
00132                 
00133                 elmnt().Uturb() = sigma*fac*dir;
00134             }
00135         }
00136         else
00137         {
00138             elmnt().tTurb() = GREAT;
00139             elmnt().Uturb() = vector::zero;
00140         }
00141     }
00142 }
00143 
00144 
00145 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00146 
00147 } // End namespace Foam
00148 
00149 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines