1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | foam-extend: Open Source CFD
4 \\ / O peration | Version: 3.2
5 \\ / A nd | Web: http://www.foam-extend.org
6 \\/ M anipulation | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
9 This file is part of foam-extend.
11 foam-extend 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 3 of the License, or (at your
14 option) any later version.
16 foam-extend is distributed in the hope that it will be useful, but
17 WITHOUT ANY WARRANTY; without even the implied warranty of
18 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 General Public License for more details.
21 You should have received a copy of the GNU General Public License
22 along with foam-extend. If not, see <http://www.gnu.org/licenses/>.
24 \*---------------------------------------------------------------------------*/
28 #include "gradientDispersionRAS.H"
29 #include "addToRunTimeSelectionTable.H"
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 defineTypeNameAndDebug(gradientDispersionRAS, 0);
40 addToRunTimeSelectionTable
43 gradientDispersionRAS,
48 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50 // Construct from components
51 gradientDispersionRAS::gradientDispersionRAS
53 const dictionary& dict,
57 dispersionRASModel(dict, sm)
61 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
63 gradientDispersionRAS::~gradientDispersionRAS()
67 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
69 void gradientDispersionRAS::disperseParcels() const
72 const scalar cps = 0.16432;
74 scalar dt = spray_.runTime().deltaT().value();
75 const volScalarField& k = turbulence().k();
76 volVectorField gradk = fvc::grad(k);
77 const volScalarField& epsilon = turbulence().epsilon();
78 const volVectorField& U = spray_.U();
82 spray::iterator elmnt = spray_.begin();
83 elmnt != spray_.end();
87 label celli = elmnt().cell();
88 scalar UrelMag = mag(elmnt().U() - U[celli] - elmnt().Uturb());
92 k[celli]/epsilon[celli],
93 cps*pow(k[celli], 1.5)/epsilon[celli]/(UrelMag + SMALL)
95 // parcel is perturbed by the turbulence
98 elmnt().tTurb() += dt;
100 if (elmnt().tTurb() > Tturb)
102 elmnt().tTurb() = 0.0;
104 scalar sigma = sqrt(2.0*k[celli]/3.0);
105 vector dir = -gradk[celli]/(mag(gradk[celli]) + SMALL);
107 // numerical recipes... Ch. 7. Random Numbers...
111 while((rsq > 1.0) || (rsq == 0.0))
113 x1 = 2.0*spray_.rndGen().scalar01() - 1.0;
114 x2 = 2.0*spray_.rndGen().scalar01() - 1.0;
118 scalar fac = sqrt(-2.0*log(rsq)/rsq);
120 // in 2D calculations the -grad(k) is always
121 // away from the axis of symmetry
122 // This creates a 'hole' in the spray and to
123 // prevent this we let x1 be both negative/positive
133 elmnt().Uturb() = sigma*fac*dir;
138 elmnt().tTurb() = GREAT;
139 elmnt().Uturb() = vector::zero;
145 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
147 } // End namespace Foam
149 // ************************************************************************* //