fixed writing out entries in advective bc
[OpenFOAM-1.6-ext.git] / src / lagrangian / intermediate / submodels / Kinematic / DispersionModel / StochasticDispersionRAS / StochasticDispersionRAS.C
bloba4036f0d1b07f6663ae5152471498deaedce7bcc
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright held by original author
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 "StochasticDispersionRAS.H"
29 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
31 template<class CloudType>
32 Foam::StochasticDispersionRAS<CloudType>::StochasticDispersionRAS
34     const dictionary& dict,
35     CloudType& owner
38     DispersionRASModel<CloudType>(dict, owner)
42 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
44 template<class CloudType>
45 Foam::StochasticDispersionRAS<CloudType>::~StochasticDispersionRAS()
49 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
51 template<class CloudType>
52 bool Foam::StochasticDispersionRAS<CloudType>::active() const
54     return true;
58 template<class CloudType>
59 Foam::vector Foam::StochasticDispersionRAS<CloudType>::update
61     const scalar dt,
62     const label celli,
63     const vector& U,
64     const vector& Uc,
65     vector& UTurb,
66     scalar& tTurb
69     const scalar cps = 0.16432;
71     const volScalarField& k = *this->kPtr_;
72     const volScalarField& epsilon = *this->epsilonPtr_;
74     const scalar UrelMag = mag(U - Uc - UTurb);
76     const scalar tTurbLoc = min
77     (
78         k[celli]/epsilon[celli],
79         cps*pow(k[celli], 1.5)/epsilon[celli]/(UrelMag + SMALL)
80     );
82     // Parcel is perturbed by the turbulence
83     if (dt < tTurbLoc)
84     {
85         tTurb += dt;
87         if (tTurb > tTurbLoc)
88         {
89             tTurb = 0.0;
91             scalar sigma = sqrt(2.0*k[celli]/3.0);
92             vector dir = 2.0*this->owner().rndGen().vector01() - vector::one;
93             dir /= mag(dir) + SMALL;
95             // Numerical Recipes... Ch. 7. Random Numbers...
96             scalar x1 = 0.0;
97             scalar x2 = 0.0;
98             scalar rsq = 10.0;
99             while ((rsq > 1.0) || (rsq == 0.0))
100             {
101                 x1 = 2.0*this->owner().rndGen().scalar01() - 1.0;
102                 x2 = 2.0*this->owner().rndGen().scalar01() - 1.0;
103                 rsq = x1*x1 + x2*x2;
104             }
106             scalar fac = sqrt(-2.0*log(rsq)/rsq);
108             fac *= mag(x1);
110             UTurb = sigma*fac*dir;
112         }
113     }
114     else
115     {
116         tTurb = GREAT;
117         UTurb = vector::zero;
118     }
120     return Uc + UTurb;
124 // ************************************************************************* //