Loading...
Searching...
No Matches
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-2024 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_;
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_;
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 (
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
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 //- Return a (basic particle) clone
279 virtual autoPtr<particle> clone() const
280 {
281 return particle::Clone(*this);
282 }
283
284 //- Return a (basic particle) clone
285 virtual autoPtr<particle> clone(const polyMesh& mesh) const
286 {
287 return particle::Clone(*this, mesh);
288 }
289
290 //- Factory class to read-construct particles (for parallel transfer)
291 class iNew
292 {
293 const polyMesh& mesh_;
294
295 public:
296
297 iNew(const polyMesh& mesh)
298 :
299 mesh_(mesh)
300 {}
301
302 autoPtr<SprayParcel<ParcelType>> operator()(Istream& is) const
303 {
304 return autoPtr<SprayParcel<ParcelType>>
305 (
306 new SprayParcel<ParcelType>(mesh_, is, true)
307 );
308 }
309 };
310
311
312 // Member Functions
313
314 // Access
315
316 //- Return const access to initial droplet diameter
317 inline scalar d0() const;
318
319 //- Return const access to initial droplet position
320 inline const vector& position0() const;
321
322 //- Return const access to the liquid surface tension
323 inline scalar sigma() const;
324
325 //- Return const access to the liquid dynamic viscosity
326 inline scalar mu() const;
328 //- Return const access to liquid core
329 inline scalar liquidCore() const;
330
331 //- Return const access to Kelvin-Helmholtz breakup index
332 inline scalar KHindex() const;
333
334 //- Return const access to spherical deviation
335 inline scalar y() const;
337 //- Return const access to rate of change of spherical deviation
338 inline scalar yDot() const;
339
340 //- Return const access to atomization characteristic time
341 inline scalar tc() const;
342
343 //- Return const access to stripped parcel mass
344 inline scalar ms() const;
345
346 //- Return const access to injector id
347 inline scalar injector() const;
348
349 //- Return const access to momentum relaxation time
350 inline scalar tMom() const;
351
352 //- Return const access to passive user scalar
353 inline scalar user() const;
354
355
356 // Edit
358 //- Return access to initial droplet diameter
359 inline scalar& d0();
360
361 //- Return access to initial droplet position
362 inline vector& position0();
364 //- Return access to the liquid surface tension
365 inline scalar& sigma();
366
367 //- Return access to the liquid dynamic viscosity
368 inline scalar& mu();
369
370 //- Return access to liquid core
371 inline scalar& liquidCore();
372
373 //- Return access to Kelvin-Helmholtz breakup index
374 inline scalar& KHindex();
375
376 //- Return access to spherical deviation
377 inline scalar& y();
378
379 //- Return access to rate of change of spherical deviation
380 inline scalar& yDot();
381
382 //- Return access to atomization characteristic time
383 inline scalar& tc();
384
385 //- Return access to stripped parcel mass
386 inline scalar& ms();
387
388 //- Return access to injector id
389 inline scalar& injector();
391 //- Return access to momentum relaxation time
392 inline scalar& tMom();
393
394 //- Return access to passive user scalar
395 inline scalar& user();
396
397
398 // Main calculation loop
399
400 //- Set cell values
401 template<class TrackCloudType>
402 void setCellValues(TrackCloudType& cloud, trackingData& td);
403
404 //- Correct parcel properties according to atomization model
405 template<class TrackCloudType>
406 void calcAtomization
407 (
408 TrackCloudType& cloud,
410 const scalar dt
411 );
412
413 //- Correct parcel properties according to breakup model
414 template<class TrackCloudType>
416 (
417 TrackCloudType& cloud,
419 const scalar dt
420 );
421
422 //- Correct cell values using latest transfer information
423 template<class TrackCloudType>
426 TrackCloudType& cloud,
428 const scalar dt
429 );
431 //- Correct surface values due to emitted species
432 template<class TrackCloudType>
434 (
435 TrackCloudType& cloud,
437 const scalar T,
438 const scalarField& Cs,
439 scalar& rhos,
440 scalar& mus,
441 scalar& Pr,
442 scalar& kappa
443 );
444
445 //- Update parcel properties over the time interval
446 template<class TrackCloudType>
447 void calc
448 (
449 TrackCloudType& cloud,
451 const scalar dt
452 );
454 //- Calculate the chi-factor for flash-boiling for the
455 // atomization model
456 template<class TrackCloudType>
457 scalar chi
459 TrackCloudType& cloud,
461 const scalarField& X
462 ) const;
464 //- Solve the TAB equation
465 template<class TrackCloudType>
466 void solveTABEq
467 (
468 TrackCloudType& cloud,
470 const scalar dt
471 );
472
474 // I-O
475
476 //- Read
477 template<class CloudType, class CompositionType>
478 static void readFields
479 (
480 CloudType& c,
481 const CompositionType& compModel
482 );
484 //- Read - no composition
485 template<class CloudType>
486 static void readFields(CloudType& c);
487
488 //- Write
489 template<class CloudType, class CompositionType>
490 static void writeFields
491 (
492 const CloudType& c,
493 const CompositionType& compModel
494 );
495
496 //- Write - composition supplied
497 template<class CloudType>
498 static void writeFields(const CloudType& c);
499
500 //- Write individual parcel properties to stream
501 void writeProperties
502 (
504 const wordRes& filters,
505 const word& delim,
506 const bool namesOnly
507 ) const;
509 //- Read particle fields as objects from the obr registry
510 // - no composition
511 template<class CloudType>
512 static void readObjects
514 CloudType& c,
515 const objectRegistry& obr
516 );
517
518 //- Read particle fields as objects from the obr registry
519 template<class CloudType, class CompositionType>
520 static void readObjects
521 (
523 const CompositionType& compModel,
524 const objectRegistry& obr
525 );
526
527 //- Write particle fields as objects into the obr registry
528 // - no composition
529 template<class CloudType>
530 static void writeObjects
531 (
532 const CloudType& c,
533 objectRegistry& obr
534 );
535
536 //- Write particle fields as objects into the obr registry
537 template<class CloudType, class CompositionType>
538 static void writeObjects
540 const CloudType& c,
541 const CompositionType& compModel,
542 objectRegistry& obr
543 );
544
545
546 // Ostream Operator
547
548 friend Ostream& operator<< <ParcelType>
549 (
552 );
553};
554
555
556// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
557
558} // End namespace Foam
559
560// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
562#include "SprayParcelI.H"
563
564// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
565
566#ifdef NoRepository
567 #include "SprayParcel.C"
568#endif
569
570
571// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
572
573#endif
574
575// ************************************************************************* //
const dimensionedScalar & pMin
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition Istream.H:60
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
scalar mu0() const
Return const access to the initial dynamic viscosity.
scalar sigma0() const
Return const access to the initial surface tension.
autoPtr< SprayParcel< ParcelType > > operator()(Istream &is) const
iNew(const polyMesh &mesh)
Reacting spray parcel, with added functionality for atomization and breakup.
Definition SprayParcel.H:57
scalar & tc()
Return access to atomization characteristic time.
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
Return a (basic particle) clone.
SprayParcel(const polyMesh &mesh, const vector &position, const label celli)
Construct from a position and a cell, searching for the rest of the.
static void writeObjects(const CloudType &c, const CompositionType &compModel, objectRegistry &obr)
Write particle fields as objects into the obr registry.
scalar & y()
Return access to spherical deviation.
scalar & liquidCore()
Return access to liquid core.
scalar & KHindex()
Return access to Kelvin-Helmholtz breakup index.
static void writeObjects(const CloudType &c, objectRegistry &obr)
Write particle fields as objects into the obr registry.
scalar & yDot()
Return access to rate of change of spherical deviation.
void cellValueSourceCorrection(TrackCloudType &cloud, trackingData &td, const scalar dt)
Correct cell values using latest transfer information.
Definition SprayParcel.C:43
SprayParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti, const label typeId, const scalar nParticle0, const scalar d0, const scalar dTarget0, const vector &U0, const vector &f0, const vector &angularMomentum0, const vector &torque0, const scalarField &Y0, const scalar liquidCore, const scalar KHindex, const scalar y, const scalar yDot, const scalar tc, const scalar ms, const scalar injector, const scalar tMom, const scalar user, const typename ParcelType::constantProperties &constProps)
Construct from components.
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 & tMom()
Return access to momentum relaxation time.
scalar mu() const
Return const access to the liquid dynamic viscosity.
scalar & d0()
Return access to initial droplet diameter.
scalar & user()
Return access to passive user scalar.
virtual autoPtr< particle > clone() const
Return a (basic particle) clone.
void calcAtomization(TrackCloudType &cloud, trackingData &td, const scalar dt)
Correct parcel properties according to atomization model.
scalar & sigma()
Return access to the liquid surface tension.
static void writeFields(const CloudType &c, const CompositionType &compModel)
Write.
TypeName("SprayParcel")
Runtime type information.
AddToPropertyList(ParcelType, " d0"+" position0"+" sigma"+" mu"+" liquidCore"+" KHindex"+" y"+" yDot"+" tc"+" ms"+" injector"+" tMom"+" user")
String representation of properties.
static void readObjects(CloudType &c, const CompositionType &compModel, const objectRegistry &obr)
Read particle fields as objects from the obr registry.
SprayParcel(const polyMesh &mesh, Istream &is, bool readFields=true, bool newFormat=true)
Construct from Istream.
ReactingParcel< ThermoParcel< KinematicParcel< particle > > >::trackingData trackingData
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.
static void writeFields(const CloudType &c)
Write - composition supplied.
vector & position0()
Return access to initial droplet position.
scalar & injector()
Return access to injector id.
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.
SprayParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from mesh, coordinates and topology.
void setCellValues(TrackCloudType &cloud, trackingData &td)
Set cell values.
Definition SprayParcel.C:31
void calcBreakup(TrackCloudType &cloud, trackingData &td, const scalar dt)
Correct parcel properties according to breakup model.
scalar sigma() const
Return const access to the liquid surface tension.
static void readFields(CloudType &c)
Read - no composition.
scalar & ms()
Return access to stripped parcel mass.
SprayParcel(const SprayParcel &p, const polyMesh &mesh)
Construct as a copy.
static void readFields(CloudType &c, const CompositionType &compModel)
const vector & position0() const
Return const access to initial droplet position.
scalar & mu()
Return access to the liquid dynamic viscosity.
void calc(TrackCloudType &cloud, trackingData &td, const scalar dt)
Update parcel properties over the time interval.
Definition SprayParcel.C:56
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A cloud is a registry collection of lagrangian particles.
Definition cloud.H:56
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:133
Registry of regIOobjects.
vector position() const
Return current particle position.
Definition particleI.H:283
static autoPtr< particle > Clone(const Derived &p)
Clone a particle.
Definition particle.H:552
const barycentric & coordinates() const noexcept
Return current particle coordinates.
Definition particleI.H:116
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
A List of wordRe with additional matching capabilities.
Definition wordRes.H:56
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
const dimensionedScalar rhoMin
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
Namespace for OpenFOAM.
DSMCCloud< dsmcParcel > CloudType
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
bool cp(const fileName &src, const fileName &dst, const bool followLink=true)
Copy the source to the destination (recursively if necessary).
Definition POSIX.C:1065
Barycentric< scalar > barycentric
A scalar version of the templated Barycentric.
Definition barycentric.H:45
Vector< scalar > vector
Definition vector.H:57
#define AddToPropertyList(ParcelType, str)
Add to existing static 'propertyList' for particle properties.
scalar rho0
scalarList Y0(nSpecie, Zero)
const scalarField & Cs
scalar T0
dimensionedScalar Pr("Pr", dimless, laminarTransport)
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68