TAB.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2015 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "TAB.H"
29
30// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31
32template<class CloudType>
34(
35 const dictionary& dict,
36 CloudType& owner
37)
38:
39 BreakupModel<CloudType>(dict, owner, typeName, true),
40 SMDCalcMethod_(this->coeffDict().lookup("SMDCalculationMethod"))
41{
42 // calculate the inverse function of the Rossin-Rammler Distribution
43 const scalar xx0 = 12.0;
44 const scalar rrd100 =
45 1.0/(1.0 - exp(-xx0)*(1.0 + xx0 + sqr(xx0)/2.0 + pow3(xx0)/6.0));
46
47 forAll(rrd_, n)
48 {
49 scalar xx = 0.12*(n + 1);
50 rrd_[n] =
51 (1.0 - exp(-xx)*(1.0 + xx + sqr(xx)/2.0 + pow3(xx)/6.0))*rrd100;
52 }
53
54 if (SMDCalcMethod_ == "method1")
55 {
56 SMDMethod_ = method1;
57 }
58 else if (SMDCalcMethod_ == "method2")
59 {
60 SMDMethod_ = method2;
61 }
62 else
63 {
64 SMDMethod_ = method2;
66 << "Unknown SMDCalculationMethod. Valid options are "
67 << "(method1 | method2). Using method2" << endl;
68 }
69}
70
71
72template<class CloudType>
74:
76 SMDCalcMethod_(bum.SMDCalcMethod_)
77{}
78
79
80// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
81
82template<class CloudType>
84{}
85
86
87// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
88
89template<class CloudType>
91(
92 const scalar dt,
93 const vector& g,
94 scalar& d,
95 scalar& tc,
96 scalar& ms,
97 scalar& nParticle,
98 scalar& KHindex,
99 scalar& y,
100 scalar& yDot,
101 const scalar d0,
102 const scalar rho,
103 const scalar mu,
104 const scalar sigma,
105 const vector& U,
106 const scalar rhoc,
107 const scalar muc,
108 const vector& Urel,
109 const scalar Urmag,
110 const scalar tMom,
111 scalar& dChild,
112 scalar& massChild
113)
114{
115 Random& rndGen = this->owner().rndGen();
116
117 scalar r = 0.5*d;
118 scalar r2 = r*r;
119 scalar r3 = r*r2;
120
121 scalar semiMass = nParticle*pow3(d);
122
123 // inverse of characteristic viscous damping time
124 scalar rtd = 0.5*this->TABCmu_*mu/(rho*r2);
125
126 // oscillation frequency (squared)
127 scalar omega2 = this->TABComega_*sigma/(rho*r3) - rtd*rtd;
128
129 if (omega2 > 0)
130 {
131 scalar omega = sqrt(omega2);
132 scalar We = rhoc*sqr(Urmag)*r/sigma;
133 scalar Wetmp = We/this->TABtwoWeCrit_;
134
135 scalar y1 = y - Wetmp;
136 scalar y2 = yDot/omega;
137
138 scalar a = sqrt(y1*y1 + y2*y2);
139
140 // scotty we may have break-up
141 if (a+Wetmp > 1.0)
142 {
143 scalar phic = y1/a;
144
145 // constrain phic within -1 to 1
146 phic = max(min(phic, 1), -1);
147
148 scalar phit = acos(phic);
149 scalar phi = phit;
150 scalar quad = -y2/a;
151 if (quad < 0)
152 {
154 }
155
156 scalar tb = 0;
157
158 if (mag(y) < 1.0)
159 {
160 scalar coste = 1.0;
161 if ((Wetmp - a < -1) && (yDot < 0))
162 {
163 coste = -1.0;
164 }
165
166 scalar theta = acos((coste-Wetmp)/a);
167
168 if (theta < phi)
169 {
170 if (constant::mathematical::twoPi - theta >= phi)
171 {
172 theta = -theta;
173 }
175 }
176 tb = (theta-phi)/omega;
177
178 // breakup occurs
179 if (dt > tb)
180 {
181 y = 1.0;
182 yDot = -a*omega*sin(omega*tb + phi);
183 }
184
185 }
186
187 // update droplet size
188 if (dt > tb)
189 {
190 scalar rs =
191 r/(1.0 + (4.0/3.0)*sqr(y) + rho*r3/(8*sigma)*sqr(yDot));
192
193 label n = 0;
194 scalar rNew = 0.0;
195 switch (SMDMethod_)
196 {
197 case method1:
198 {
199 #include "TABSMDCalcMethod1.H"
200 break;
201 }
202 case method2:
203 {
204 #include "TABSMDCalcMethod2.H"
205 break;
206 }
207 }
208
209 if (rNew < r)
210 {
211 d = 2*rNew;
212 y = 0;
213 yDot = 0;
214 }
215 }
216 }
217 }
218 else
219 {
220 // reset droplet distortion parameters
221 y = 0;
222 yDot = 0;
223 }
224
225 // update the nParticle count to conserve mass
226 nParticle = semiMass/pow3(d);
227
228 // Do not add child parcel
229 return false;
230}
231
232
233// ************************************************************************* //
scalar y
surfaceScalarField phic(mixture.cAlpha() *mag(alphaPhic/mesh.magSf()))
label n
const uniformDimensionedVectorField & g
surfaceScalarField & phi
Templated break-up model class.
Definition: BreakupModel.H:61
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:75
Random number generator.
Definition: Random.H:60
The TAB Method for Numerical Calculation of Spray Droplet Breakup.
Definition: TAB.H:69
virtual ~TAB()
Destructor.
Definition: TAB.C:83
@ method1
Definition: TAB.H:75
@ method2
Definition: TAB.H:76
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
virtual bool update()
Update the mesh for both mesh motion and topology change.
Lookup type of boundary radiation properties.
Definition: lookup.H:66
U
Definition: pEqn.H:72
const volScalarField & mu
Urel
Definition: pEqn.H:56
#define WarningInFunction
Report a warning using Foam::Warning.
constexpr scalar twoPi(2 *M_PI)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensionedScalar sin(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensionedScalar y1(const dimensionedScalar &ds)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
dimensionedScalar acos(const dimensionedScalar &ds)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
Random rndGen
Definition: createFields.H:23