ETAB.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-2014 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 "ETAB.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 k1_(0.2),
41 k2_(0.2),
42 WeTransition_(100.0),
43 AWe_(0.0)
44{
45 if (!this->defaultCoeffs(true))
46 {
47 this->coeffDict().readEntry("k1", k1_);
48 this->coeffDict().readEntry("k2", k2_);
49 this->coeffDict().readEntry("WeTransition", WeTransition_);
50 }
51
52 scalar k21 = k2_/k1_;
53 AWe_ = (k21*sqrt(WeTransition_) - 1.0)/pow4(WeTransition_);
54}
55
56
57template<class CloudType>
59:
61 k1_(bum.k1_),
62 k2_(bum.k2_),
63 WeTransition_(bum.WeTransition_),
64 AWe_(bum.AWe_)
65{}
66
67
68// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
69
70template<class CloudType>
72{}
73
74
75// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
76
77template<class CloudType>
79(
80 const scalar dt,
81 const vector& g,
82 scalar& d,
83 scalar& tc,
84 scalar& ms,
85 scalar& nParticle,
86 scalar& KHindex,
87 scalar& y,
88 scalar& yDot,
89 const scalar d0,
90 const scalar rho,
91 const scalar mu,
92 const scalar sigma,
93 const vector& U,
94 const scalar rhoc,
95 const scalar muc,
96 const vector& Urel,
97 const scalar Urmag,
98 const scalar tMom,
99 scalar& dChild,
100 scalar& massChild
101)
102{
103 scalar r = 0.5*d;
104 scalar r2 = r*r;
105 scalar r3 = r*r2;
106
107 scalar semiMass = nParticle*pow3(d);
108
109 // inverse of characteristic viscous damping time
110 scalar rtd = 0.5*this->TABCmu_*mu/(rho*r2);
111
112 // oscillation frequency (squared)
113 scalar omega2 = this->TABComega_*sigma/(rho*r3) - rtd*rtd;
114
115 if (omega2 > 0)
116 {
117 scalar omega = sqrt(omega2);
118 scalar romega = 1.0/omega;
119
120 scalar We = rhoc*sqr(Urmag)*r/sigma;
121 scalar Wetmp = We/this->TABtwoWeCrit_;
122
123 scalar y1 = y - Wetmp;
124 scalar y2 = yDot*romega;
125
126 scalar a = sqrt(y1*y1 + y2*y2);
127
128 // scotty we may have break-up
129 if (a + Wetmp > 1.0)
130 {
131 scalar phic = y1/a;
132
133 // constrain phic within -1 to 1
134 phic = max(min(phic, 1), -1);
135
136 scalar phit = acos(phic);
137 scalar phi = phit;
138 scalar quad = -y2/a;
139 if (quad < 0)
140 {
142 }
143
144 scalar tb = 0;
145
146 if (mag(y) < 1.0)
147 {
148 scalar theta = acos((1.0 - Wetmp)/a);
149
150 if (theta < phi)
151 {
152 if (constant::mathematical::twoPi - theta >= phi)
153 {
154 theta = -theta;
155 }
157 }
158 tb = (theta - phi)*romega;
159
160 // breakup occurs
161 if (dt > tb)
162 {
163 y = 1.0;
164 yDot = -a*omega*sin(omega*tb + phi);
165 }
166 }
167
168 // update droplet size
169 if (dt > tb)
170 {
171 scalar sqrtWe = AWe_*pow4(We) + 1.0;
172 scalar Kbr = k1_*omega*sqrtWe;
173
174 if (We > WeTransition_)
175 {
176 sqrtWe = sqrt(We);
177 Kbr =k2_*omega*sqrtWe;
178 }
179
180 scalar rWetmp = 1.0/Wetmp;
181 scalar cosdtbu = max(-1.0, min(1.0, 1.0 - rWetmp));
182 scalar dtbu = romega*acos(cosdtbu);
183 scalar decay = exp(-Kbr*dtbu);
184
185 scalar rNew = decay*r;
186 if (rNew < r)
187 {
188 d = 2.0*rNew;
189 y = 0.0;
190 yDot = 0.0;
191 }
192 }
193 }
194 }
195 else
196 {
197 // reset droplet distortion parameters
198 y = 0;
199 yDot = 0;
200 }
201
202 // update the nParticle count to conserve mass
203 nParticle = semiMass/pow3(d);
204
205 // Do not add child parcel
206 return false;
207}
208
209
210// ************************************************************************* //
scalar y
surfaceScalarField phic(mixture.cAlpha() *mag(alphaPhic/mesh.magSf()))
const uniformDimensionedVectorField & g
surfaceScalarField & phi
Templated break-up model class.
Definition: BreakupModel.H:61
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:75
The Enhanced TAB model.
Definition: ETAB.H:72
virtual ~ETAB()
Destructor.
Definition: ETAB.C:71
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
virtual bool update()
Update the mesh for both mesh motion and topology change.
const dictionary & coeffDict() const
Return const access to the coefficients dictionary.
Definition: subModelBase.C:131
virtual bool defaultCoeffs(const bool printMsg) const
Returns true if defaultCoeffs is true and outputs on printMsg.
Definition: subModelBase.C:143
U
Definition: pEqn.H:72
const volScalarField & mu
Urel
Definition: pEqn.H:56
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)
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 pow4(const dimensionedScalar &ds)
dimensionedScalar acos(const dimensionedScalar &ds)
dictionary dict