Go to the documentation of this file.00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026 #include <OpenFOAM/error.H>
00027
00028 #include "stochasticDispersionRAS.H"
00029 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00030
00031
00032
00033 namespace Foam
00034 {
00035
00036
00037
00038 defineTypeNameAndDebug(stochasticDispersionRAS, 0);
00039
00040 addToRunTimeSelectionTable
00041 (
00042 dispersionModel,
00043 stochasticDispersionRAS,
00044 dictionary
00045 );
00046
00047
00048
00049
00050
00051 stochasticDispersionRAS::stochasticDispersionRAS
00052 (
00053 const dictionary& dict,
00054 spray& sm
00055 )
00056 :
00057 dispersionRASModel(dict, sm)
00058 {}
00059
00060
00061
00062
00063 stochasticDispersionRAS::~stochasticDispersionRAS()
00064 {}
00065
00066
00067
00068
00069 void stochasticDispersionRAS::disperseParcels() const
00070 {
00071
00072 const scalar cps = 0.16432;
00073 const vector one(1.0, 1.0, 1.0);
00074
00075 scalar dt = spray_.runTime().deltaT().value();
00076 const volScalarField& k = turbulence().k();
00077
00078 const volScalarField& epsilon = turbulence().epsilon();
00079 const volVectorField& U = spray_.U();
00080
00081 for
00082 (
00083 spray::iterator elmnt = spray_.begin();
00084 elmnt != spray_.end();
00085 ++elmnt
00086 )
00087 {
00088 label celli = elmnt().cell();
00089 scalar UrelMag = mag(elmnt().U() - U[celli] - elmnt().Uturb());
00090
00091 scalar Tturb = min
00092 (
00093 k[celli]/epsilon[celli],
00094 cps*pow(k[celli], 1.5)/epsilon[celli]/(UrelMag + SMALL)
00095 );
00096
00097
00098 if (dt < Tturb)
00099 {
00100 elmnt().tTurb() += dt;
00101
00102 if (elmnt().tTurb() > Tturb)
00103 {
00104 elmnt().tTurb() = 0.0;
00105
00106 scalar sigma = sqrt(2.0*k[celli]/3.0);
00107 vector dir = 2.0*spray_.rndGen().vector01() - one;
00108 dir /= mag(dir) + SMALL;
00109
00110
00111 scalar x1,x2;
00112 scalar rsq = 10.0;
00113 while((rsq > 1.0) || (rsq == 0.0))
00114 {
00115 x1 = 2.0*spray_.rndGen().scalar01() - 1.0;
00116 x2 = 2.0*spray_.rndGen().scalar01() - 1.0;
00117 rsq = x1*x1 + x2*x2;
00118 }
00119
00120 scalar fac = sqrt(-2.0*log(rsq)/rsq);
00121
00122 fac *= mag(x1);
00123
00124 elmnt().Uturb() = sigma*fac*dir;
00125
00126 }
00127 }
00128 else
00129 {
00130 elmnt().tTurb() = GREAT;
00131 elmnt().Uturb() = vector::zero;
00132 }
00133 }
00134 }
00135
00136
00137
00138
00139 }
00140
00141