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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30
31template<class ParcelType>
33:
34 ParcelType::constantProperties(),
35 sigma0_(this->dict_, 0.0),
36 mu0_(this->dict_, 0.0)
37{}
38
39
40template<class ParcelType>
42(
44)
45:
46 ParcelType::constantProperties(cp),
47 sigma0_(cp.sigma0_),
48 mu0_(cp.mu0_)
49{}
50
51
52template<class ParcelType>
54(
55 const dictionary& parentDict
56)
57:
58 ParcelType::constantProperties(parentDict),
59 sigma0_(this->dict_, "sigma0"),
60 mu0_(this->dict_, "mu0")
61{}
62
63
64template<class ParcelType>
66(
67 const label parcelTypeId,
68 const scalar rhoMin,
69 const scalar rho0,
70 const scalar minParcelMass,
71 const scalar youngsModulus,
72 const scalar poissonsRatio,
73 const scalar T0,
74 const scalar TMin,
75 const scalar TMax,
76 const scalar Cp0,
77 const scalar epsilon0,
78 const scalar f0,
79 const scalar Pr,
80 const scalar pMin,
81 const bool constantVolume,
82 const scalar sigma0,
83 const scalar mu0
84)
85:
86 ParcelType::constantProperties
87 (
88 parcelTypeId,
89 rhoMin,
90 rho0,
91 minParcelMass,
92 youngsModulus,
93 poissonsRatio,
94 T0,
95 TMin,
96 TMax,
97 Cp0,
98 epsilon0,
99 f0,
100 Pr,
101 pMin,
102 constantVolume
103 ),
104 sigma0_(this->dict_, sigma0),
105 mu0_(this->dict_, mu0)
106{}
107
108
109template<class ParcelType>
111(
112 const polyMesh& mesh,
114 const label celli,
115 const label tetFacei,
116 const label tetPti
117)
118:
119 ParcelType(mesh, coordinates, celli, tetFacei, tetPti),
120 d0_(this->d()),
121 position0_(this->position()),
122 sigma_(0.0),
123 mu_(0.0),
124 liquidCore_(0.0),
125 KHindex_(0.0),
126 y_(0.0),
127 yDot_(0.0),
128 tc_(0.0),
129 ms_(0.0),
130 injector_(1.0),
131 tMom_(GREAT),
132 user_(0.0)
133{}
134
135
136template<class ParcelType>
138(
139 const polyMesh& mesh,
140 const vector& position,
141 const label celli
142)
143:
144 ParcelType(mesh, position, celli),
145 d0_(this->d()),
146 position0_(this->position()),
147 sigma_(0.0),
148 mu_(0.0),
149 liquidCore_(0.0),
150 KHindex_(0.0),
151 y_(0.0),
152 yDot_(0.0),
153 tc_(0.0),
154 ms_(0.0),
155 injector_(1.0),
156 tMom_(GREAT),
157 user_(0.0)
158{}
159
160
161template<class ParcelType>
163(
164 const polyMesh& mesh,
166 const label celli,
167 const label tetFacei,
168 const label tetPti,
169 const label typeId,
170 const scalar nParticle0,
171 const scalar d0,
172 const scalar dTarget0,
173 const vector& U0,
174 const vector& f0,
175 const vector& angularMomentum0,
176 const vector& torque0,
177 const scalarField& Y0,
178 const scalar liquidCore,
179 const scalar KHindex,
180 const scalar y,
181 const scalar yDot,
182 const scalar tc,
183 const scalar ms,
184 const scalar injector,
185 const scalar tMom,
186 const scalar user,
187 const typename ParcelType::constantProperties& constProps
188)
189:
190 ParcelType
191 (
192 mesh,
194 celli,
195 tetFacei,
196 tetPti,
197 typeId,
198 nParticle0,
199 d0,
200 dTarget0,
201 U0,
202 f0,
203 angularMomentum0,
204 torque0,
205 Y0,
206 constProps
207 ),
208 d0_(d0),
209 position0_(this->position()),
210 sigma_(constProps.sigma0()),
211 mu_(constProps.mu0()),
212 liquidCore_(liquidCore),
213 KHindex_(KHindex),
214 y_(y),
215 yDot_(yDot),
216 tc_(tc),
217 ms_(ms),
218 injector_(injector),
219 tMom_(tMom),
220 user_(user)
221{}
222
223
224// * * * * * * * * * constantProperties Member Functions * * * * * * * * * * //
225
226template<class ParcelType>
227inline Foam::scalar
229{
230 return sigma0_.value();
231}
232
233
234template<class ParcelType>
235inline Foam::scalar
237{
238 return mu0_.value();
239}
240
241
242// * * * * * * * * * * SprayParcel Member Functions * * * * * * * * * * * * //
243
244template<class ParcelType>
245inline Foam::scalar Foam::SprayParcel<ParcelType>::d0() const
246{
247 return d0_;
248}
249
250
251template<class ParcelType>
253{
254 return position0_;
255}
256
257
258template<class ParcelType>
259inline Foam::scalar Foam::SprayParcel<ParcelType>::sigma() const
260{
261 return sigma_;
262}
263
264
265template<class ParcelType>
266inline Foam::scalar Foam::SprayParcel<ParcelType>::mu() const
267{
268 return mu_;
269}
270
271
272template<class ParcelType>
274{
275 return liquidCore_;
276}
277
278
279template<class ParcelType>
280inline Foam::scalar Foam::SprayParcel<ParcelType>::KHindex() const
281{
282 return KHindex_;
283}
284
285
286template<class ParcelType>
287inline Foam::scalar Foam::SprayParcel<ParcelType>::y() const
288{
289 return y_;
290}
291
292
293template<class ParcelType>
294inline Foam::scalar Foam::SprayParcel<ParcelType>::yDot() const
295{
296 return yDot_;
297}
298
299
300template<class ParcelType>
301inline Foam::scalar Foam::SprayParcel<ParcelType>::tc() const
302{
303 return tc_;
304}
305
306
307template<class ParcelType>
308inline Foam::scalar Foam::SprayParcel<ParcelType>::ms() const
309{
310 return ms_;
311}
312
313
314template<class ParcelType>
316{
317 return injector_;
318}
319
320
321template<class ParcelType>
322inline Foam::scalar Foam::SprayParcel<ParcelType>::tMom() const
323{
324 return tMom_;
325}
326
327
328template<class ParcelType>
329inline Foam::scalar Foam::SprayParcel<ParcelType>::user() const
330{
331 return user_;
332}
333
334
335template<class ParcelType>
337{
338 return d0_;
339}
340
341
342template<class ParcelType>
344{
345 return position0_;
346}
347
348
349template<class ParcelType>
351{
352 return sigma_;
353}
354
355
356template<class ParcelType>
358{
359 return mu_;
360}
361
362
363template<class ParcelType>
365{
366 return liquidCore_;
367}
368
369
370template<class ParcelType>
372{
373 return KHindex_;
374}
375
376
377template<class ParcelType>
379{
380 return y_;
381}
382
383
384template<class ParcelType>
386{
387 return yDot_;
388}
389
390
391template<class ParcelType>
393{
394 return tc_;
395}
396
397
398template<class ParcelType>
400{
401 return ms_;
402}
403
404
405template<class ParcelType>
407{
408 return injector_;
409}
410
411
412template<class ParcelType>
414{
415 return tMom_;
416}
417
418
419template<class ParcelType>
421{
422 return user_;
423}
424
425
426// ************************************************************************* //
scalar y
const dimensionedScalar & pMin
Class to hold reacting particle constant properties.
Definition: SprayParcel.H:70
scalar mu0() const
Return const access to the initial dynamic viscosity.
Definition: SprayParcelI.H:236
scalar sigma0() const
Return const access to the initial surface tension.
Definition: SprayParcelI.H:228
Reacting spray parcel, with added functionality for atomization and breakup.
Definition: SprayParcel.H:63
scalar injector() const
Return const access to injector id.
Definition: SprayParcelI.H:315
scalar liquidCore() const
Return const access to liquid core.
Definition: SprayParcelI.H:273
scalar d0_
Initial droplet diameter.
Definition: SprayParcel.H:137
scalar tMom() const
Return const access to momentum relaxation time.
Definition: SprayParcelI.H:322
vector position0_
Injection position.
Definition: SprayParcel.H:140
scalar mu() const
Return const access to the liquid dynamic viscosity.
Definition: SprayParcelI.H:266
scalar mu_
Liquid dynamic viscosity [Pa.s].
Definition: SprayParcel.H:146
scalar sigma_
Liquid surface tension [N/m].
Definition: SprayParcel.H:143
scalar user_
Passive scalar (extra variable to be defined by user)
Definition: SprayParcel.H:174
scalar yDot() const
Return const access to rate of change of spherical deviation.
Definition: SprayParcelI.H:294
scalar d0() const
Return const access to initial droplet diameter.
Definition: SprayParcelI.H:245
scalar tMom_
Momentum relaxation time (needed for calculating parcel acc.)
Definition: SprayParcel.H:171
scalar KHindex() const
Return const access to Kelvin-Helmholtz breakup index.
Definition: SprayParcelI.H:280
scalar tc() const
Return const access to atomization characteristic time.
Definition: SprayParcelI.H:301
scalar y_
Spherical deviation.
Definition: SprayParcel.H:155
scalar tc_
Characteristic time (used in atomization and/or breakup model)
Definition: SprayParcel.H:161
scalar injector_
Injected from injector (needed e.g. for calculating distance.
Definition: SprayParcel.H:168
scalar KHindex_
Index for KH Breakup.
Definition: SprayParcel.H:152
scalar ms() const
Return const access to stripped parcel mass.
Definition: SprayParcelI.H:308
scalar y() const
Return const access to spherical deviation.
Definition: SprayParcelI.H:287
scalar sigma() const
Return const access to the liquid surface tension.
Definition: SprayParcelI.H:259
scalar ms_
Stripped parcel mass due to breakup.
Definition: SprayParcel.H:164
scalar yDot_
Rate of change of spherical deviation.
Definition: SprayParcel.H:158
scalar user() const
Return const access to passive user scalar.
Definition: SprayParcelI.H:329
scalar liquidCore_
Part of liquid core ( >0.5=liquid, <0.5=droplet )
Definition: SprayParcel.H:149
const vector & position0() const
Return const access to initial droplet position.
Definition: SprayParcelI.H:252
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
dynamicFvMesh & mesh
PtrList< coordinateSystem > coordinates(solidRegions.size())
const dimensionedScalar rhoMin
scalar rho0
scalarList Y0(nSpecie, Zero)
const volScalarField & cp
scalar T0
Definition: createFields.H:22
dimensionedScalar Pr("Pr", dimless, laminarTransport)