LISAAtomization.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-2016 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
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 "LISAAtomization.H"
29 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
32 template<class CloudType>
34 (
35  const dictionary& dict,
36  CloudType& owner
37 )
38 :
39  AtomizationModel<CloudType>(dict, owner, typeName),
40  Cl_(this->coeffDict().getScalar("Cl")),
41  cTau_(this->coeffDict().getScalar("cTau")),
42  Q_(this->coeffDict().getScalar("Q")),
43  lisaExp_(this->coeffDict().getScalar("lisaExp")),
44  injectorDirection_
45  (
46  this->coeffDict().template get<vector>("injectorDirection")
47  ),
48  SMDCalcMethod_(this->coeffDict().getWord("SMDCalculationMethod"))
49 {
50  // Note: Would be good if this could be picked up from the injector
51  injectorDirection_.normalise();
52 
53  if (SMDCalcMethod_ == "method1")
54  {
55  SMDMethod_ = method1;
56  }
57  else if (SMDCalcMethod_ == "method2")
58  {
59  SMDMethod_ = method2;
60  }
61  else
62  {
63  SMDMethod_ = method2;
64  Info<< "Warning: SMDCalculationMethod " << SMDCalcMethod_
65  << " unknown. Options are (method1 | method2). Using method2"
66  << endl;
67  }
68 }
69 
70 template<class CloudType>
72 (
74 )
75 :
77  Cl_(am.Cl_),
78  cTau_(am.cTau_),
79  Q_(am.Q_),
80  lisaExp_(am.lisaExp_),
81  injectorDirection_(am.injectorDirection_),
82  SMDCalcMethod_(am.SMDCalcMethod_)
83 {}
84 
85 
86 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
87 
88 template<class CloudType>
90 {}
91 
92 
93 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
94 
95 template<class CloudType>
97 {
98  return 1.0;
99 }
100 
101 
102 template<class CloudType>
104 {
105  return true;
106 }
107 
108 
109 template<class CloudType>
111 (
112  const scalar dt,
113  scalar& d,
114  scalar& liquidCore,
115  scalar& tc,
116  const scalar rho,
117  const scalar mu,
118  const scalar sigma,
119  const scalar volFlowRate,
120  const scalar rhoAv,
121  const scalar Urel,
122  const vector& pos,
123  const vector& injectionPos,
124  const scalar pAmbient,
125  const scalar chi,
126  Random& rndGen
127 ) const
128 {
129  if (volFlowRate < SMALL)
130  {
131  return;
132  }
133 
134  scalar tau = 0.0;
135  scalar dL = 0.0;
136  scalar k = 0.0;
137 
138  // update atomization characteristic time
139  tc += dt;
140 
141  scalar We = 0.5*rhoAv*sqr(Urel)*d/sigma;
142  scalar nu = mu/rho;
143 
144  scalar Q = rhoAv/rho;
145 
146  vector diff = pos - injectionPos;
147  scalar pWalk = mag(diff);
148  scalar traveledTime = pWalk/Urel;
149 
150  scalar h = diff & injectorDirection_;
151  scalar delta = sqrt(sqr(pWalk) - sqr(h));
152 
153  scalar hSheet = volFlowRate/(constant::mathematical::pi*delta*Urel);
154 
155  // update drop diameter
156  d = min(d, hSheet);
157 
158  if (We > 27.0/16.0)
159  {
160  scalar kPos = 0.0;
161  scalar kNeg = Q*sqr(Urel)*rho/sigma;
162 
163  scalar derivPos = sqrt(Q*sqr(Urel));
164 
165  scalar derivNeg =
166  (
167  8.0*sqr(nu)*pow3(kNeg)
168  + Q*sqr(Urel)*kNeg
169  - 3.0*sigma/2.0/rho*sqr(kNeg)
170  )
171  / sqrt
172  (
173  4.0*sqr(nu)*pow4(kNeg)
174  + Q*sqr(Urel)*sqr(kNeg)
175  - sigma*pow3(kNeg)/rho
176  )
177  - 4.0*nu*kNeg;
178 
179  scalar kOld = 0.0;
180 
181  for (label i=0; i<40; i++)
182  {
183  k = kPos - (derivPos/((derivNeg - derivPos)/(kNeg - kPos)));
184 
185  scalar derivk =
186  (
187  8.0*sqr(nu)*pow3(k)
188  + Q*sqr(Urel)*k
189  - 3.0*sigma/2.0/rho*sqr(k)
190  )
191  / sqrt
192  (
193  4.0*sqr(nu)*pow4(k)
194  + Q*sqr(Urel)*sqr(k)
195  - sigma*pow3(k)/rho
196  )
197  - 4.0*nu*k;
198 
199  if (derivk > 0)
200  {
201  derivPos = derivk;
202  kPos = k;
203  }
204  else
205  {
206  derivNeg = derivk;
207  kNeg = k;
208  }
209 
210  if (mag(k - kOld)/k < 1e-4)
211  {
212  break;
213  }
214 
215  kOld = k;
216  }
217 
218  scalar omegaS =
219  - 2.0*nu*sqr(k)
220  + sqrt
221  (
222  4.0*sqr(nu)*pow4(k)
223  + Q*sqr(Urel)*sqr(k)
224  - sigma*pow3(k)/rho
225  );
226 
227  tau = cTau_/omegaS;
228 
229  dL = sqrt(8.0*d/k);
230  }
231  else
232  {
233  k = rhoAv*sqr(Urel)/(2.0*sigma);
234 
235  scalar J = 0.5*traveledTime*hSheet;
236 
237  tau = pow(3.0*cTau_, 2.0/3.0)*cbrt(J*sigma/(sqr(Q)*pow4(Urel)*rho));
238 
239  dL = sqrt(4.0*d/k);
240  }
241 
242  scalar kL = 1.0/(dL*sqrt(0.5 + 1.5 * mu/sqrt(rho*sigma*dL)));
243 
244  scalar dD = cbrt(3.0*constant::mathematical::pi*sqr(dL)/kL);
245 
246  scalar atmPressure = 1.0e+5;
247 
248  scalar pRatio = pAmbient/atmPressure;
249 
250  dD = dD*pow(pRatio, lisaExp_);
251 
252  scalar pExp = 0.135;
253 
254  // modifying dD to take account of flash boiling
255  dD = dD*(1.0 - chi*pow(pRatio, -pExp));
256  scalar lBU = Cl_ * mag(Urel)*tau;
257 
258  if (pWalk > lBU)
259  {
260  scalar x = 0;
261 
262  switch (SMDMethod_)
263  {
264  case method1:
265  {
266  #include "LISASMDCalcMethod1.H"
267  break;
268  }
269  case method2:
270  {
271  #include "LISASMDCalcMethod2.H"
272  break;
273  }
274  }
275 
276  // New droplet properties
277  liquidCore = 0.0;
278  d = x;
279  tc = 0.0;
280  }
281 }
282 
283 
284 // ************************************************************************* //
Foam::Random
Random number generator.
Definition: Random.H:59
Foam::LISAAtomization::calcChi
virtual bool calcChi() const
Flag to indicate if chi needs to be calculated.
Definition: LISAAtomization.C:103
Foam::constant::physicoChemical::mu
const dimensionedScalar mu
Atomic mass unit.
Definition: createFieldRefs.H:4
Urel
Urel
Definition: pEqn.H:56
Foam::AtomizationModel
Templated atomization model class.
Definition: SprayCloud.H:49
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
LISASMDCalcMethod1.H
rho
rho
Definition: readInitialConditions.H:88
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::diff
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
Definition: triad.C:378
nu
volScalarField & nu
Definition: readMechanicalProperties.H:176
LISASMDCalcMethod2.H
Foam::pow4
dimensionedScalar pow4(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:100
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:89
Foam::constant::universal::h
const dimensionedScalar h
Planck constant.
Definition: setRegionSolidFields.H:33
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
delta
scalar delta
Definition: LISASMDCalcMethod2.H:8
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
Foam::DSMCCloud
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:71
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::LISAAtomization::update
virtual void update(const scalar dt, scalar &d, scalar &liquidCore, scalar &tc, const scalar rho, const scalar mu, const scalar sigma, const scalar volFlowRate, const scalar rhoAv, const scalar Urel, const vector &pos, const vector &injectionPos, const scalar pAmbient, const scalar chi, Random &rndGen) const
Definition: LISAAtomization.C:111
Foam::LISAAtomization::LISAAtomization
LISAAtomization(const dictionary &, CloudType &)
Construct from dictionary.
Definition: LISAAtomization.C:34
Foam::LISAAtomization
Primary Breakup Model for pressure swirl atomizers.
Definition: LISAAtomization.H:66
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::constant::mathematical::pi
constexpr scalar pi(M_PI)
Foam::Vector< scalar >
LISAAtomization.H
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::LISAAtomization::initLiquidCore
virtual scalar initLiquidCore() const
Initial value of liquidCore.
Definition: LISAAtomization.C:96
x
x
Definition: LISASMDCalcMethod2.H:52
rndGen
Random rndGen
Definition: createFields.H:23
sigma
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
Foam::cbrt
dimensionedScalar cbrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:155
Foam::LISAAtomization::~LISAAtomization
virtual ~LISAAtomization()
Destructor.
Definition: LISAAtomization.C:89
Foam::pos
dimensionedScalar pos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:177