SprayParcelI.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2016 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
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.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
27 
28 template<class ParcelType>
30 :
31  ParcelType::constantProperties(),
32  sigma0_(this->dict_, 0.0),
33  mu0_(this->dict_, 0.0)
34 {}
35 
36 
37 template<class ParcelType>
39 (
40  const constantProperties& cp
41 )
42 :
43  ParcelType::constantProperties(cp),
44  sigma0_(cp.sigma0_),
45  mu0_(cp.mu0_)
46 {}
47 
48 
49 template<class ParcelType>
51 (
52  const dictionary& parentDict
53 )
54 :
55  ParcelType::constantProperties(parentDict),
56  sigma0_(this->dict_, "sigma0"),
57  mu0_(this->dict_, "mu0")
58 {}
59 
60 
61 template<class ParcelType>
63 (
64  const label parcelTypeId,
65  const scalar rhoMin,
66  const scalar rho0,
67  const scalar minParcelMass,
68  const scalar youngsModulus,
69  const scalar poissonsRatio,
70  const scalar T0,
71  const scalar TMin,
72  const scalar TMax,
73  const scalar Cp0,
74  const scalar epsilon0,
75  const scalar f0,
76  const scalar Pr,
77  const scalar pMin,
78  const Switch& constantVolume,
79  const scalar sigma0,
80  const scalar mu0
81 )
82 :
83  ParcelType::constantProperties
84  (
85  parcelTypeId,
86  rhoMin,
87  rho0,
88  minParcelMass,
89  youngsModulus,
90  poissonsRatio,
91  T0,
92  TMin,
93  TMax,
94  Cp0,
95  epsilon0,
96  f0,
97  Pr,
98  pMin,
99  constantVolume
100  ),
101  sigma0_(this->dict_, sigma0),
102  mu0_(this->dict_, mu0)
103 {}
104 
105 
106 template<class ParcelType>
108 (
109  const polyMesh& mesh,
110  const vector& position,
111  const label celli,
112  const label tetFacei,
113  const label tetPtI
114 )
115 :
116  ParcelType(mesh, position, celli, tetFacei, tetPtI),
117  d0_(this->d()),
118  position0_(position),
119  sigma_(0.0),
120  mu_(0.0),
121  liquidCore_(0.0),
122  KHindex_(0.0),
123  y_(0.0),
124  yDot_(0.0),
125  tc_(0.0),
126  ms_(0.0),
127  injector_(1.0),
128  tMom_(GREAT),
129  user_(0.0)
130 {}
131 
132 
133 template<class ParcelType>
135 (
136  const polyMesh& mesh,
137  const vector& position,
138  const label celli,
139  const label tetFacei,
140  const label tetPtI,
141  const label typeId,
142  const scalar nParticle0,
143  const scalar d0,
144  const scalar dTarget0,
145  const vector& U0,
146  const vector& f0,
147  const vector& angularMomentum0,
148  const vector& torque0,
149  const scalarField& Y0,
150  const scalar liquidCore,
151  const scalar KHindex,
152  const scalar y,
153  const scalar yDot,
154  const scalar tc,
155  const scalar ms,
156  const scalar injector,
157  const scalar tMom,
158  const scalar user,
159  const typename ParcelType::constantProperties& constProps
160 )
161 :
162  ParcelType
163  (
164  mesh,
165  position,
166  celli,
167  tetFacei,
168  tetPtI,
169  typeId,
170  nParticle0,
171  d0,
172  dTarget0,
173  U0,
174  f0,
175  angularMomentum0,
176  torque0,
177  Y0,
178  constProps
179  ),
180  d0_(d0),
181  position0_(position),
182  sigma_(constProps.sigma0()),
183  mu_(constProps.mu0()),
184  liquidCore_(liquidCore),
185  KHindex_(KHindex),
186  y_(y),
187  yDot_(yDot),
188  tc_(tc),
189  ms_(ms),
190  injector_(injector),
191  tMom_(tMom),
192  user_(user)
193 {}
194 
195 
196 // * * * * * * * * * constantProperties Member Functions * * * * * * * * * * //
197 
198 template<class ParcelType>
199 inline Foam::scalar
201 {
202  return sigma0_.value();
203 }
204 
205 
206 template<class ParcelType>
207 inline Foam::scalar
209 {
210  return mu0_.value();
211 }
212 
213 
214 // * * * * * * * * * * SprayParcel Member Functions * * * * * * * * * * * * //
215 
216 template<class ParcelType>
217 inline Foam::scalar Foam::SprayParcel<ParcelType>::d0() const
218 {
219  return d0_;
220 }
221 
222 
223 template<class ParcelType>
225 {
226  return position0_;
227 }
228 
229 
230 template<class ParcelType>
231 inline Foam::scalar Foam::SprayParcel<ParcelType>::sigma() const
232 {
233  return sigma_;
234 }
235 
236 
237 template<class ParcelType>
238 inline Foam::scalar Foam::SprayParcel<ParcelType>::mu() const
239 {
240  return mu_;
241 }
242 
243 
244 template<class ParcelType>
245 inline Foam::scalar Foam::SprayParcel<ParcelType>::liquidCore() const
246 {
247  return liquidCore_;
248 }
249 
250 
251 template<class ParcelType>
252 inline Foam::scalar Foam::SprayParcel<ParcelType>::KHindex() const
253 {
254  return KHindex_;
255 }
256 
257 
258 template<class ParcelType>
259 inline Foam::scalar Foam::SprayParcel<ParcelType>::y() const
260 {
261  return y_;
262 }
263 
264 
265 template<class ParcelType>
266 inline Foam::scalar Foam::SprayParcel<ParcelType>::yDot() const
267 {
268  return yDot_;
269 }
270 
271 
272 template<class ParcelType>
273 inline Foam::scalar Foam::SprayParcel<ParcelType>::tc() const
274 {
275  return tc_;
276 }
277 
278 
279 template<class ParcelType>
280 inline Foam::scalar Foam::SprayParcel<ParcelType>::ms() const
281 {
282  return ms_;
283 }
284 
285 
286 template<class ParcelType>
287 inline Foam::scalar Foam::SprayParcel<ParcelType>::injector() const
288 {
289  return injector_;
290 }
291 
292 
293 template<class ParcelType>
294 inline Foam::scalar Foam::SprayParcel<ParcelType>::tMom() const
295 {
296  return tMom_;
297 }
298 
299 
300 template<class ParcelType>
301 inline Foam::scalar Foam::SprayParcel<ParcelType>::user() const
302 {
303  return user_;
304 }
305 
306 
307 template<class ParcelType>
308 inline Foam::scalar& Foam::SprayParcel<ParcelType>::d0()
309 {
310  return d0_;
311 }
312 
313 
314 template<class ParcelType>
316 {
317  return position0_;
318 }
319 
320 
321 template<class ParcelType>
323 {
324  return sigma_;
325 }
326 
327 
328 template<class ParcelType>
329 inline Foam::scalar& Foam::SprayParcel<ParcelType>::mu()
330 {
331  return mu_;
332 }
333 
334 
335 template<class ParcelType>
337 {
338  return liquidCore_;
339 }
340 
341 
342 template<class ParcelType>
344 {
345  return KHindex_;
346 }
347 
348 
349 template<class ParcelType>
350 inline Foam::scalar& Foam::SprayParcel<ParcelType>::y()
351 {
352  return y_;
353 }
354 
355 
356 template<class ParcelType>
358 {
359  return yDot_;
360 }
361 
362 
363 template<class ParcelType>
364 inline Foam::scalar& Foam::SprayParcel<ParcelType>::tc()
365 {
366  return tc_;
367 }
368 
369 
370 template<class ParcelType>
371 inline Foam::scalar& Foam::SprayParcel<ParcelType>::ms()
372 {
373  return ms_;
374 }
375 
376 
377 template<class ParcelType>
379 {
380  return injector_;
381 }
382 
383 
384 template<class ParcelType>
386 {
387  return tMom_;
388 }
389 
390 
391 template<class ParcelType>
393 {
394  return user_;
395 }
396 
397 
398 // ************************************************************************* //
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
scalar d0() const
Return const access to initial droplet diameter.
Definition: SprayParcelI.H:217
scalar tMom_
Momentum relaxation time (needed for calculating parcel acc.)
Definition: SprayParcel.H:170
scalar mu0() const
Return const access to the initial dynamic viscosity.
Definition: SprayParcelI.H:208
scalar yDot() const
Return const access to rate of change of spherical deviation.
Definition: SprayParcelI.H:266
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, y/n, t/f, or none.
Definition: Switch.H:60
const Type & value() const
Return the value.
scalar y() const
Return const access to spherical deviation.
Definition: SprayParcelI.H:259
scalar liquidCore_
Part of liquid core ( >0.5=liquid, <0.5=droplet )
Definition: SprayParcel.H:148
scalar tMom() const
Return const access to momentum relaxation time.
Definition: SprayParcelI.H:294
scalar liquidCore() const
Return const access to liquid core.
Definition: SprayParcelI.H:245
scalar injector() const
Return const access to injector id.
Definition: SprayParcelI.H:287
scalar sigma0() const
Return const access to the initial surface tension.
Definition: SprayParcelI.H:200
SprayParcel(const polyMesh &mesh, const vector &position, const label celli, const label tetFacei, const label tetPtI)
Construct from owner, position, and cloud owner.
Definition: SprayParcelI.H:108
const vector & position0() const
Return const access to initial droplet position.
Definition: SprayParcelI.H:224
scalar ms() const
Return const access to stripped parcel mass.
Definition: SprayParcelI.H:280
scalar injector_
Injected from injector (needed e.g. for calculating distance.
Definition: SprayParcel.H:167
Class to hold reacting particle constant properties.
Definition: SprayParcel.H:70
scalar user_
Passive scalar (extra variable to be defined by user)
Definition: SprayParcel.H:173
scalar KHindex() const
Return const access to Kelvin-Helmholtz breakup index.
Definition: SprayParcelI.H:252
constantProperties()
Null constructor.
Definition: SprayParcelI.H:29
scalar sigma_
Liquid surface tension [N/m].
Definition: SprayParcel.H:142
scalar sigma() const
Return const access to the liquid surface tension.
Definition: SprayParcelI.H:231
vector position0_
Injection position.
Definition: SprayParcel.H:139
scalar yDot_
Rate of change of spherical deviation.
Definition: SprayParcel.H:157
scalar ms_
Stripped parcel mass due to breakup.
Definition: SprayParcel.H:163
scalar d0_
Initial droplet diameter.
Definition: SprayParcel.H:136
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
scalar mu() const
Return const access to the liquid dynamic viscosity.
Definition: SprayParcelI.H:238
scalar tc() const
Return const access to atomization characteristic time.
Definition: SprayParcelI.H:273
scalar tc_
Characteristic time (used in atomization and/or breakup model)
Definition: SprayParcel.H:160
scalar mu_
Liquid dynamic viscosity [Pa.s].
Definition: SprayParcel.H:145
scalar y_
Spherical deviation.
Definition: SprayParcel.H:154
scalar KHindex_
Index for KH Breakup.
Definition: SprayParcel.H:151
scalar user() const
Return const access to passive user scalar.
Definition: SprayParcelI.H:301