initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / lagrangian / dieselSpray / spraySubModels / dispersionModel / gradientDispersionRAS / gradientDispersionRAS.C
blobf5c846ebf86280b945bf71379e75f881c7d742c7
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
11     OpenFOAM is free software; you can redistribute it and/or modify it
12     under the terms of the GNU General Public License as published by the
13     Free Software Foundation; either version 2 of the License, or (at your
14     option) any later version.
16     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
19     for more details.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "error.H"
29 #include "gradientDispersionRAS.H"
30 #include "addToRunTimeSelectionTable.H"
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 namespace Foam
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 defineTypeNameAndDebug(gradientDispersionRAS, 0);
41 addToRunTimeSelectionTable
43     dispersionModel,
44     gradientDispersionRAS,
45     dictionary
49 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
51 // Construct from components
52 gradientDispersionRAS::gradientDispersionRAS
54     const dictionary& dict,
55     spray& sm
58     dispersionRASModel(dict, sm)
62 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
64 gradientDispersionRAS::~gradientDispersionRAS()
68 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
70 void gradientDispersionRAS::disperseParcels() const
73     const scalar cps = 0.16432;
75     scalar dt = spray_.runTime().deltaT().value();
76     const volScalarField& k = turbulence().k();
77     volVectorField gradk = fvc::grad(k);
78     const volScalarField& epsilon = turbulence().epsilon();
79     const volVectorField& U = spray_.U();
81     for
82     (
83         spray::iterator elmnt = spray_.begin();
84         elmnt != spray_.end();
85         ++elmnt
86     )
87     {
88         label celli = elmnt().cell();
89         scalar UrelMag = mag(elmnt().U() - U[celli] - elmnt().Uturb());
91         scalar Tturb = min
92         (
93             k[celli]/epsilon[celli], 
94             cps*pow(k[celli], 1.5)/epsilon[celli]/(UrelMag + SMALL)
95         );
96         // parcel is perturbed by the turbulence
97         if (dt < Tturb)
98         {
99             elmnt().tTurb() += dt;
101             if (elmnt().tTurb() > Tturb)
102             {
103                 elmnt().tTurb() = 0.0;
105                 scalar sigma = sqrt(2.0*k[celli]/3.0);
106                 vector dir = -gradk[celli]/(mag(gradk[celli]) + SMALL);
108                 // numerical recipes... Ch. 7. Random Numbers...
109                 scalar x1 = 0.0;
110                 scalar x2 = 0.0;
111                 scalar rsq = 10.0;
112                 while((rsq > 1.0) || (rsq == 0.0))
113                 {
114                     x1 = 2.0*spray_.rndGen().scalar01() - 1.0;
115                     x2 = 2.0*spray_.rndGen().scalar01() - 1.0;
116                     rsq = x1*x1 + x2*x2;
117                 }
118                 
119                 scalar fac = sqrt(-2.0*log(rsq)/rsq);
120                 
121                 // in 2D calculations the -grad(k) is always
122                 // away from the axis of symmetry
123                 // This creates a 'hole' in the spray and to
124                 // prevent this we let x1 be both negative/positive
125                 if (spray_.twoD())
126                 {
127                     fac *= x1;
128                 }
129                 else
130                 {
131                     fac *= mag(x1);
132                 }
133                 
134                 elmnt().Uturb() = sigma*fac*dir;
135             }
136         }
137         else
138         {
139             elmnt().tTurb() = GREAT;
140             elmnt().Uturb() = vector::zero;
141         }
142     }
146 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
148 } // End namespace Foam
150 // ************************************************************************* //