SprayParcel.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) 2016-2019 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
27Class
28 Foam::SprayParcel
29
30Description
31 Reacting spray parcel, with added functionality for atomization and breakup
32
33\*---------------------------------------------------------------------------*/
34
35#ifndef SprayParcel_H
36#define SprayParcel_H
37
38#include "particle.H"
39#include "demandDrivenEntry.H"
40
41// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42
43namespace Foam
44{
45
46template<class ParcelType>
47class SprayParcel;
48
49template<class ParcelType>
50Ostream& operator<<
51(
52 Ostream&,
54);
55
56/*---------------------------------------------------------------------------*\
57 Class SprayParcel Declaration
58\*---------------------------------------------------------------------------*/
59
60template<class ParcelType>
61class SprayParcel
62:
63 public ParcelType
64{
65public:
66
67 //- Class to hold reacting particle constant properties
69 :
70 public ParcelType::constantProperties
71 {
72 // Private data
73
74 //- Particle initial surface tension [N/m]
76
77 //- Particle initial dynamic viscosity [Pa.s]
79
80
81 public:
82
83 // Constructors
84
85 //- Null constructor
87
88 //- Copy constructor
90
91 //- Construct from dictionary
92 constantProperties(const dictionary& parentDict);
93
94 //- Construct from components
96 (
97 const label parcelTypeId,
98 const scalar rhoMin,
99 const scalar rho0,
100 const scalar minParcelMass,
101 const scalar youngsModulus,
102 const scalar poissonsRatio,
103 const scalar T0,
104 const scalar TMin,
105 const scalar TMax,
106 const scalar Cp0,
107 const scalar epsilon0,
108 const scalar f0,
109 const scalar Pr,
110 const scalar pMin,
111 const bool constantVolume,
112 const scalar sigma0,
113 const scalar mu0
114 );
115
116
117 // Access
118
119 //- Return const access to the initial surface tension
120 inline scalar sigma0() const;
121
122 //- Return const access to the initial dynamic viscosity
123 inline scalar mu0() const;
124 };
125
126
127 //- Use base tracking data
128 typedef typename ParcelType::trackingData trackingData;
129
130
131protected:
132
133 // Protected data
134
135 // Spray parcel properties
136
137 //- Initial droplet diameter
138 scalar d0_;
139
140 //- Injection position
142
143 //- Liquid surface tension [N/m]
144 scalar sigma_;
145
146 //- Liquid dynamic viscosity [Pa.s]
147 scalar mu_;
148
149 //- Part of liquid core ( >0.5=liquid, <0.5=droplet )
150 scalar liquidCore_;
151
152 //- Index for KH Breakup
153 scalar KHindex_;
154
155 //- Spherical deviation
156 scalar y_;
157
158 //- Rate of change of spherical deviation
159 scalar yDot_;
160
161 //- Characteristic time (used in atomization and/or breakup model)
162 scalar tc_;
163
164 //- Stripped parcel mass due to breakup
165 scalar ms_;
166
167 //- Injected from injector (needed e.g. for calculating distance
168 // from injector)
169 scalar injector_;
170
171 //- Momentum relaxation time (needed for calculating parcel acc.)
172 scalar tMom_;
173
174 //- Passive scalar (extra variable to be defined by user)
175 scalar user_;
176
177
178public:
179
180 // Static data members
181
182 //- Size in bytes of the fields
183 static const std::size_t sizeofFields;
184
185 //- Runtime type information
186 TypeName("SprayParcel");
187
188 //- String representation of properties
190 (
191 ParcelType,
192 " d0"
193 + " position0"
194 + " sigma"
195 + " mu"
196 + " liquidCore"
197 + " KHindex"
198 + " y"
199 + " yDot"
200 + " tc"
201 + " ms"
202 + " injector"
203 + " tMom"
204 + " user"
205 );
206
207
208 // Constructors
209
210 //- Construct from mesh, coordinates and topology
211 // Other properties initialised as null
212 inline SprayParcel
213 (
214 const polyMesh& mesh,
216 const label celli,
217 const label tetFacei,
218 const label tetPti
219 );
220
221 //- Construct from a position and a cell, searching for the rest of the
222 // required topology. Other properties are initialised as null.
223 inline SprayParcel
224 (
225 const polyMesh& mesh,
226 const vector& position,
227 const label celli
228 );
229
230 //- Construct from components
231 inline SprayParcel
232 (
233 const polyMesh& mesh,
235 const label celli,
236 const label tetFacei,
237 const label tetPti,
238 const label typeId,
239 const scalar nParticle0,
240 const scalar d0,
241 const scalar dTarget0,
242 const vector& U0,
243 const vector& f0,
244 const vector& angularMomentum0,
245 const vector& torque0,
246 const scalarField& Y0,
247 const scalar liquidCore,
248 const scalar KHindex,
249 const scalar y,
250 const scalar yDot,
251 const scalar tc,
252 const scalar ms,
253 const scalar injector,
254 const scalar tMom,
255 const scalar user,
256 const typename ParcelType::constantProperties& constProps
257 );
258
259 //- Construct from Istream
261 (
262 const polyMesh& mesh,
263 Istream& is,
264 bool readFields = true,
265 bool newFormat = true
266 );
267
268 //- Construct as a copy
270 (
271 const SprayParcel& p,
272 const polyMesh& mesh
273 );
274
275 //- Construct as a copy
276 SprayParcel(const SprayParcel& p);
277
278 //- Construct and return a (basic particle) clone
279 virtual autoPtr<particle> clone() const
280 {
282 }
283
284 //- Construct and return a (basic particle) clone
285 virtual autoPtr<particle> clone(const polyMesh& mesh) const
286 {
287 return autoPtr<particle>
288 (
289 new SprayParcel<ParcelType>(*this, mesh)
290 );
291 }
292
293 //- Factory class to read-construct particles used for
294 // parallel transfer
295 class iNew
296 {
297 const polyMesh& mesh_;
298
299 public:
301 iNew(const polyMesh& mesh)
302 :
303 mesh_(mesh)
304 {}
307 {
309 (
310 new SprayParcel<ParcelType>(mesh_, is, true)
311 );
312 }
313 };
314
315
316 // Member Functions
317
318 // Access
319
320 //- Return const access to initial droplet diameter
321 inline scalar d0() const;
322
323 //- Return const access to initial droplet position
324 inline const vector& position0() const;
325
326 //- Return const access to the liquid surface tension
327 inline scalar sigma() const;
328
329 //- Return const access to the liquid dynamic viscosity
330 inline scalar mu() const;
331
332 //- Return const access to liquid core
333 inline scalar liquidCore() const;
334
335 //- Return const access to Kelvin-Helmholtz breakup index
336 inline scalar KHindex() const;
337
338 //- Return const access to spherical deviation
339 inline scalar y() const;
340
341 //- Return const access to rate of change of spherical deviation
342 inline scalar yDot() const;
343
344 //- Return const access to atomization characteristic time
345 inline scalar tc() const;
346
347 //- Return const access to stripped parcel mass
348 inline scalar ms() const;
349
350 //- Return const access to injector id
351 inline scalar injector() const;
352
353 //- Return const access to momentum relaxation time
354 inline scalar tMom() const;
355
356 //- Return const access to passive user scalar
357 inline scalar user() const;
358
359
360 // Edit
361
362 //- Return access to initial droplet diameter
363 inline scalar& d0();
364
365 //- Return access to initial droplet position
366 inline vector& position0();
367
368 //- Return access to the liquid surface tension
369 inline scalar& sigma();
370
371 //- Return access to the liquid dynamic viscosity
372 inline scalar& mu();
373
374 //- Return access to liquid core
375 inline scalar& liquidCore();
376
377 //- Return access to Kelvin-Helmholtz breakup index
378 inline scalar& KHindex();
379
380 //- Return access to spherical deviation
381 inline scalar& y();
382
383 //- Return access to rate of change of spherical deviation
384 inline scalar& yDot();
385
386 //- Return access to atomization characteristic time
387 inline scalar& tc();
388
389 //- Return access to stripped parcel mass
390 inline scalar& ms();
391
392 //- Return access to injector id
393 inline scalar& injector();
394
395 //- Return access to momentum relaxation time
396 inline scalar& tMom();
397
398 //- Return access to passive user scalar
399 inline scalar& user();
400
401
402 // Main calculation loop
403
404 //- Set cell values
405 template<class TrackCloudType>
406 void setCellValues(TrackCloudType& cloud, trackingData& td);
407
408 //- Correct parcel properties according to atomization model
409 template<class TrackCloudType>
410 void calcAtomization
411 (
412 TrackCloudType& cloud,
413 trackingData& td,
414 const scalar dt
415 );
416
417 //- Correct parcel properties according to breakup model
418 template<class TrackCloudType>
419 void calcBreakup
420 (
421 TrackCloudType& cloud,
422 trackingData& td,
423 const scalar dt
424 );
425
426 //- Correct cell values using latest transfer information
427 template<class TrackCloudType>
429 (
430 TrackCloudType& cloud,
431 trackingData& td,
432 const scalar dt
433 );
434
435 //- Correct surface values due to emitted species
436 template<class TrackCloudType>
438 (
439 TrackCloudType& cloud,
440 trackingData& td,
441 const scalar T,
442 const scalarField& Cs,
443 scalar& rhos,
444 scalar& mus,
445 scalar& Pr,
446 scalar& kappa
447 );
448
449 //- Update parcel properties over the time interval
450 template<class TrackCloudType>
451 void calc
452 (
453 TrackCloudType& cloud,
454 trackingData& td,
455 const scalar dt
456 );
457
458 //- Calculate the chi-factor for flash-boiling for the
459 // atomization model
460 template<class TrackCloudType>
461 scalar chi
462 (
463 TrackCloudType& cloud,
464 trackingData& td,
465 const scalarField& X
466 ) const;
467
468 //- Solve the TAB equation
469 template<class TrackCloudType>
470 void solveTABEq
471 (
472 TrackCloudType& cloud,
473 trackingData& td,
474 const scalar dt
475 );
476
477
478 // I-O
479
480 //- Read
481 template<class CloudType, class CompositionType>
482 static void readFields
483 (
484 CloudType& c,
485 const CompositionType& compModel
486 );
487
488 //- Read - no composition
489 template<class CloudType>
490 static void readFields(CloudType& c);
491
492 //- Write
493 template<class CloudType, class CompositionType>
494 static void writeFields
495 (
496 const CloudType& c,
497 const CompositionType& compModel
498 );
499
500 //- Write - composition supplied
501 template<class CloudType>
502 static void writeFields(const CloudType& c);
503
504 //- Write individual parcel properties to stream
505 void writeProperties
506 (
507 Ostream& os,
508 const wordRes& filters,
509 const word& delim,
510 const bool namesOnly
511 ) const;
512
513 //- Read particle fields as objects from the obr registry
514 // - no composition
515 template<class CloudType>
516 static void readObjects
517 (
518 CloudType& c,
519 const objectRegistry& obr
520 );
521
522 //- Read particle fields as objects from the obr registry
523 template<class CloudType, class CompositionType>
524 static void readObjects
525 (
526 CloudType& c,
527 const CompositionType& compModel,
528 const objectRegistry& obr
529 );
530
531 //- Write particle fields as objects into the obr registry
532 // - no composition
533 template<class CloudType>
534 static void writeObjects
535 (
536 const CloudType& c,
537 objectRegistry& obr
538 );
539
540 //- Write particle fields as objects into the obr registry
541 template<class CloudType, class CompositionType>
542 static void writeObjects
543 (
544 const CloudType& c,
545 const CompositionType& compModel,
546 objectRegistry& obr
547 );
548
549
550 // Ostream Operator
552 friend Ostream& operator<< <ParcelType>
553 (
554 Ostream&,
556 );
557};
558
559
560// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
561
562} // End namespace Foam
563
564// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
565
566#include "SprayParcelI.H"
567
568// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
569
570#ifdef NoRepository
571 #include "SprayParcel.C"
572#endif
573
574
575// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
576
577#endif
578
579// ************************************************************************* //
const dimensionedScalar & pMin
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:75
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:64
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
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
Factory class to read-construct particles used for.
Definition: SprayParcel.H:295
autoPtr< SprayParcel< ParcelType > > operator()(Istream &is) const
Definition: SprayParcel.H:305
iNew(const polyMesh &mesh)
Definition: SprayParcel.H:300
Reacting spray parcel, with added functionality for atomization and breakup.
Definition: SprayParcel.H:63
scalar chi(TrackCloudType &cloud, trackingData &td, const scalarField &X) const
Calculate the chi-factor for flash-boiling for the.
virtual autoPtr< particle > clone(const polyMesh &mesh) const
Construct and return a (basic particle) clone.
Definition: SprayParcel.H:284
scalar injector() const
Return const access to injector id.
Definition: SprayParcelI.H:315
static void writeObjects(const CloudType &c, objectRegistry &obr)
Write particle fields as objects into the obr registry.
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
void cellValueSourceCorrection(TrackCloudType &cloud, trackingData &td, const scalar dt)
Correct cell values using latest transfer information.
Definition: SprayParcel.C:50
void correctSurfaceValues(TrackCloudType &cloud, trackingData &td, const scalar T, const scalarField &Cs, scalar &rhos, scalar &mus, scalar &Pr, scalar &kappa)
Correct surface values due to emitted species.
scalar mu() const
Return const access to the liquid dynamic viscosity.
Definition: SprayParcelI.H:266
virtual autoPtr< particle > clone() const
Construct and return a (basic particle) clone.
Definition: SprayParcel.H:278
void calcAtomization(TrackCloudType &cloud, trackingData &td, const scalar dt)
Correct parcel properties according to atomization model.
Definition: SprayParcel.C:151
scalar mu_
Liquid dynamic viscosity [Pa.s].
Definition: SprayParcel.H:146
static const std::size_t sizeofFields
Size in bytes of the fields.
Definition: SprayParcel.H:182
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
static void writeFields(const CloudType &c, const CompositionType &compModel)
Write.
TypeName("SprayParcel")
Runtime type information.
scalar yDot() const
Return const access to rate of change of spherical deviation.
Definition: SprayParcelI.H:294
AddToPropertyList(ParcelType, " d0"+" position0"+" sigma"+" mu"+" liquidCore"+" KHindex"+" y"+" yDot"+" tc"+" ms"+" injector"+" tMom"+" user")
String representation of properties.
scalar d0() const
Return const access to initial droplet diameter.
Definition: SprayParcelI.H:245
ParcelType::trackingData trackingData
Use base tracking data.
Definition: SprayParcel.H:127
static void readObjects(CloudType &c, const objectRegistry &obr)
Read particle fields as objects from the obr registry.
void solveTABEq(TrackCloudType &cloud, trackingData &td, const scalar dt)
Solve the TAB equation.
Definition: SprayParcel.C:373
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
void writeProperties(Ostream &os, const wordRes &filters, const word &delim, const bool namesOnly) const
Write individual parcel properties to stream.
SprayParcel(const SprayParcel &p)
Construct as a copy.
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
void setCellValues(TrackCloudType &cloud, trackingData &td)
Set cell values.
Definition: SprayParcel.C:38
void calcBreakup(TrackCloudType &cloud, trackingData &td, const scalar dt)
Correct parcel properties according to breakup model.
Definition: SprayParcel.C:220
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
SprayParcel(const SprayParcel &p, const polyMesh &mesh)
Construct as a copy.
static void readFields(CloudType &c, const CompositionType &compModel)
Read.
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
void calc(TrackCloudType &cloud, trackingData &td, const scalar dt)
Update parcel properties over the time interval.
Definition: SprayParcel.C:63
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A cloud is a registry collection of lagrangian particles.
Definition: cloud.H:60
Class for demand-driven dictionary entries.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Class used to pass data into container.
Registry of regIOobjects.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:54
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
const volScalarField & T
dynamicFvMesh & mesh
PtrList< coordinateSystem > coordinates(solidRegions.size())
OBJstream os(runTime.globalPath()/outputName)
const dimensionedScalar rhoMin
Namespace for OpenFOAM.
#define AddToPropertyList(ParcelType, str)
Add to existing static 'propertyList' for particle properties.
scalar rho0
scalarList Y0(nSpecie, Zero)
const volScalarField & cp
const scalarField & Cs
scalar T0
Definition: createFields.H:22
dimensionedScalar Pr("Pr", dimless, laminarTransport)
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73