SprayParcel.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2019 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 Class
25  Foam::SprayParcel
26 
27 Description
28  Reaching spray parcel, with added functionality for atomization and breakup
29 
30 \*---------------------------------------------------------------------------*/
31 
32 #ifndef SprayParcel_H
33 #define SprayParcel_H
34 
35 #include "particle.H"
36 #include "demandDrivenEntry.H"
37 
38 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42 
43 template<class ParcelType>
44 class SprayParcel;
45 
46 template<class ParcelType>
47 Ostream& operator<<
48 (
49  Ostream&,
51 );
52 
53 /*---------------------------------------------------------------------------*\
54  Class SprayParcel Declaration
55 \*---------------------------------------------------------------------------*/
56 
57 template<class ParcelType>
58 class SprayParcel
59 :
60  public ParcelType
61 {
62  // Private Data
63 
64  //- Size in bytes of the fields
65  static const std::size_t sizeofFields_;
66 
67 
68 public:
69 
70  //- Class to hold reacting particle constant properties
71  class constantProperties
72  :
73  public ParcelType::constantProperties
74  {
75  // Private Data
76 
77  //- Particle initial surface tension [N/m]
79 
80  //- Particle initial dynamic viscosity [Pa.s]
82 
83 
84  public:
85 
86  // Constructors
87 
88  //- Null constructor
90 
91  //- Copy constructor
93 
94  //- Construct from dictionary
95  constantProperties(const dictionary& parentDict);
96 
97  //- Construct from components
99  (
100  const label parcelTypeId,
101  const scalar rhoMin,
102  const scalar rho0,
103  const scalar minParcelMass,
104  const scalar youngsModulus,
105  const scalar poissonsRatio,
106  const scalar T0,
107  const scalar TMin,
108  const scalar TMax,
109  const scalar Cp0,
110  const scalar epsilon0,
111  const scalar f0,
112  const scalar Pr,
113  const scalar pMin,
114  const Switch& constantVolume,
115  const scalar sigma0,
116  const scalar mu0
117  );
118 
119 
120  // Access
121 
122  //- Return const access to the initial surface tension
123  inline scalar sigma0() const;
124 
125  //- Return const access to the initial dynamic viscosity
126  inline scalar mu0() const;
127  };
128 
129 
130  //- Use base tracking data
131  typedef typename ParcelType::trackingData trackingData;
132 
133 
134 protected:
135 
136  // Protected data
137 
138  // Spray parcel properties
139 
140  //- Initial droplet diameter
141  scalar d0_;
142 
143  //- Injection position
145 
146  //- Liquid surface tension [N/m]
147  scalar sigma_;
148 
149  //- Liquid dynamic viscosity [Pa.s]
150  scalar mu_;
151 
152  //- Part of liquid core ( >0.5=liquid, <0.5=droplet )
153  scalar liquidCore_;
154 
155  //- Index for KH Breakup
156  scalar KHindex_;
157 
158  //- Spherical deviation
159  scalar y_;
160 
161  //- Rate of change of spherical deviation
162  scalar yDot_;
163 
164  //- Characteristic time (used in atomization and/or breakup model)
165  scalar tc_;
166 
167  //- Stripped parcel mass due to breakup
168  scalar ms_;
169 
170  //- Injected from injector (needed e.g. for calculating distance
171  // from injector)
172  scalar injector_;
173 
174  //- Momentum relaxation time (needed for calculating parcel acc.)
175  scalar tMom_;
176 
177  //- Passive scalar (extra variable to be defined by user)
178  scalar user_;
179 
180 
181 public:
182 
183  // Static Data Members
184 
185  //- Runtime type information
186  TypeName("SprayParcel");
187 
188 
189  // Constructors
190 
191  //- Construct from mesh, coordinates and topology
192  // Other properties initialised as null
193  inline SprayParcel
194  (
195  const polyMesh& mesh,
196  const barycentric& coordinates,
197  const label celli,
198  const label tetFacei,
199  const label tetPti
200  );
201 
202  //- Construct from a position and a cell, searching for the rest of the
203  // required topology. Other properties are initialised as null.
204  inline SprayParcel
205  (
206  const polyMesh& mesh,
207  const vector& position,
208  const label celli
209  );
210 
211  //- Construct from components
212  inline SprayParcel
213  (
214  const polyMesh& mesh,
215  const barycentric& coordinates,
216  const label celli,
217  const label tetFacei,
218  const label tetPti,
219  const label typeId,
220  const scalar nParticle0,
221  const scalar d0,
222  const scalar dTarget0,
223  const vector& U0,
224  const vector& f0,
225  const vector& angularMomentum0,
226  const vector& torque0,
227  const scalarField& Y0,
228  const scalar liquidCore,
229  const scalar KHindex,
230  const scalar y,
231  const scalar yDot,
232  const scalar tc,
233  const scalar ms,
234  const scalar injector,
235  const scalar tMom,
236  const scalar user,
237  const typename ParcelType::constantProperties& constProps
238  );
239 
240  //- Construct from Istream
242  (
243  const polyMesh& mesh,
244  Istream& is,
245  bool readFields = true
246  );
247 
248  //- Construct as a copy
250  (
251  const SprayParcel& p,
252  const polyMesh& mesh
253  );
254 
255  //- Construct as a copy
256  SprayParcel(const SprayParcel& p);
257 
258  //- Construct and return a (basic particle) clone
259  virtual autoPtr<particle> clone() const
260  {
261  return autoPtr<particle>(new SprayParcel<ParcelType>(*this));
262  }
263 
264  //- Construct and return a (basic particle) clone
265  virtual autoPtr<particle> clone(const polyMesh& mesh) const
266  {
267  return autoPtr<particle>
268  (
269  new SprayParcel<ParcelType>(*this, mesh)
270  );
271  }
272 
273  //- Factory class to read-construct particles used for
274  // parallel transfer
275  class iNew
276  {
277  const polyMesh& mesh_;
278 
279  public:
281  iNew(const polyMesh& mesh)
282  :
283  mesh_(mesh)
284  {}
286  autoPtr<SprayParcel<ParcelType>> operator()(Istream& is) const
287  {
289  (
290  new SprayParcel<ParcelType>(mesh_, is, true)
291  );
292  }
293  };
294 
295 
296  // Member Functions
297 
298  // Access
299 
300  //- Return const access to initial droplet diameter
301  inline scalar d0() const;
302 
303  //- Return const access to initial droplet position
304  inline const vector& position0() const;
305 
306  //- Return const access to the liquid surface tension
307  inline scalar sigma() const;
308 
309  //- Return const access to the liquid dynamic viscosity
310  inline scalar mu() const;
311 
312  //- Return const access to liquid core
313  inline scalar liquidCore() const;
314 
315  //- Return const access to Kelvin-Helmholtz breakup index
316  inline scalar KHindex() const;
317 
318  //- Return const access to spherical deviation
319  inline scalar y() const;
320 
321  //- Return const access to rate of change of spherical deviation
322  inline scalar yDot() const;
323 
324  //- Return const access to atomization characteristic time
325  inline scalar tc() const;
326 
327  //- Return const access to stripped parcel mass
328  inline scalar ms() const;
329 
330  //- Return const access to injector id
331  inline scalar injector() const;
332 
333  //- Return const access to momentum relaxation time
334  inline scalar tMom() const;
335 
336  //- Return const access to passive user scalar
337  inline scalar user() const;
338 
339 
340  // Edit
341 
342  //- Return access to initial droplet diameter
343  inline scalar& d0();
344 
345  //- Return access to initial droplet position
346  inline vector& position0();
347 
348  //- Return access to the liquid surface tension
349  inline scalar& sigma();
350 
351  //- Return access to the liquid dynamic viscosity
352  inline scalar& mu();
353 
354  //- Return access to liquid core
355  inline scalar& liquidCore();
356 
357  //- Return access to Kelvin-Helmholtz breakup index
358  inline scalar& KHindex();
359 
360  //- Return access to spherical deviation
361  inline scalar& y();
362 
363  //- Return access to rate of change of spherical deviation
364  inline scalar& yDot();
365 
366  //- Return access to atomization characteristic time
367  inline scalar& tc();
368 
369  //- Return access to stripped parcel mass
370  inline scalar& ms();
371 
372  //- Return access to injector id
373  inline scalar& injector();
374 
375  //- Return access to momentum relaxation time
376  inline scalar& tMom();
377 
378  //- Return access to passive user scalar
379  inline scalar& user();
380 
381 
382  // Main calculation loop
383 
384  //- Set cell values
385  template<class TrackCloudType>
386  void setCellValues(TrackCloudType& cloud, trackingData& td);
387 
388  //- Correct parcel properties according to atomization model
389  template<class TrackCloudType>
390  void calcAtomization
391  (
392  TrackCloudType& cloud,
393  trackingData& td,
394  const scalar dt
395  );
396 
397  //- Correct parcel properties according to breakup model
398  template<class TrackCloudType>
399  void calcBreakup
400  (
401  TrackCloudType& cloud,
402  trackingData& td,
403  const scalar dt
404  );
405 
406  //- Correct cell values using latest transfer information
407  template<class TrackCloudType>
409  (
410  TrackCloudType& cloud,
411  trackingData& td,
412  const scalar dt
413  );
414 
415  //- Correct surface values due to emitted species
416  template<class TrackCloudType>
418  (
419  TrackCloudType& cloud,
420  trackingData& td,
421  const scalar T,
422  const scalarField& Cs,
423  scalar& rhos,
424  scalar& mus,
425  scalar& Pr,
426  scalar& kappa
427  );
428 
429  //- Update parcel properties over the time interval
430  template<class TrackCloudType>
431  void calc
432  (
433  TrackCloudType& cloud,
434  trackingData& td,
435  const scalar dt
436  );
437 
438  //- Calculate the chi-factor for flash-boiling for the
439  // atomization model
440  template<class TrackCloudType>
441  scalar chi
442  (
443  TrackCloudType& cloud,
444  trackingData& td,
445  const scalarField& X
446  ) const;
447 
448  //- Solve the TAB equation
449  template<class TrackCloudType>
450  void solveTABEq
451  (
452  TrackCloudType& cloud,
453  trackingData& td,
454  const scalar dt
455  );
456 
457 
458  // I-O
459 
460  //- Read
461  template<class CloudType, class CompositionType>
462  static void readFields
463  (
464  CloudType& c,
465  const CompositionType& compModel
466  );
467 
468  //- Read - no composition
469  template<class CloudType>
470  static void readFields(CloudType& c);
471 
472  //- Write
473  template<class CloudType, class CompositionType>
474  static void writeFields
475  (
476  const CloudType& c,
477  const CompositionType& compModel
478  );
479 
480  //- Write - composition supplied
481  template<class CloudType>
482  static void writeFields(const CloudType& c);
483 
484 
485  // Ostream Operator
486 
487  friend Ostream& operator<< <ParcelType>
488  (
489  Ostream&,
491  );
492 };
493 
494 
495 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
496 
497 } // End namespace Foam
498 
499 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
500 
501 #include "SprayParcelI.H"
502 
503 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
504 
505 #ifdef NoRepository
506  #include "SprayParcel.C"
507 #endif
508 
509 
510 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
511 
512 #endif
513 
514 // ************************************************************************* //
scalar mu0() const
Return const access to the initial dynamic viscosity.
Definition: SprayParcelI.H:233
Factory class to read-construct particles used for.
Definition: SprayParcel.H:274
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.
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
const dimensionedScalar & kappa
Coulomb constant: default SI units: [N.m2/C2].
const dimensionedScalar & epsilon0
Electric constant: default SI units: [F/m].
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
scalar tMom_
Momentum relaxation time (needed for calculating parcel acc.)
Definition: SprayParcel.H:174
scalar KHindex() const
Return const access to Kelvin-Helmholtz breakup index.
Definition: SprayParcelI.H:277
scalar tMom() const
Return const access to momentum relaxation time.
Definition: SprayParcelI.H:319
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
scalar sigma0() const
Return const access to the initial surface tension.
Definition: SprayParcelI.H:225
scalar liquidCore() const
Return const access to liquid core.
Definition: SprayParcelI.H:270
bool cp(const fileName &src, const fileName &dst, const bool followLink=true)
Copy, recursively if necessary, the source to the destination.
Definition: POSIX.C:753
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/any.
Definition: Switch.H:60
void calcAtomization(TrackCloudType &cloud, trackingData &td, const scalar dt)
Correct parcel properties according to atomization model.
Definition: SprayParcel.C:149
scalar rho0
void calcBreakup(TrackCloudType &cloud, trackingData &td, const scalar dt)
Correct parcel properties according to breakup model.
Definition: SprayParcel.C:216
const dimensionedScalar & c
Speed of light in a vacuum.
scalar ms() const
Return const access to stripped parcel mass.
Definition: SprayParcelI.H:305
scalar liquidCore_
Part of liquid core ( >0.5=liquid, <0.5=droplet )
Definition: SprayParcel.H:152
SprayParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from mesh, coordinates and topology.
Definition: SprayParcelI.H:108
scalar y() const
Return const access to spherical deviation.
Definition: SprayParcelI.H:284
scalarList Y0(nSpecie, 0.0)
scalar chi(TrackCloudType &cloud, trackingData &td, const scalarField &X) const
Calculate the chi-factor for flash-boiling for the.
scalar d0() const
Return const access to initial droplet diameter.
Definition: SprayParcelI.H:242
dynamicFvMesh & mesh
A cloud is a collection of lagrangian particles.
Definition: cloud.H:51
scalar mu() const
Return const access to the liquid dynamic viscosity.
Definition: SprayParcelI.H:263
scalar tc() const
Return const access to atomization characteristic time.
Definition: SprayParcelI.H:298
scalar injector_
Injected from injector (needed e.g. for calculating distance.
Definition: SprayParcel.H:171
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:177
static void readFields(CloudType &c, const CompositionType &compModel)
Read.
dimensionedScalar pMin("pMin", dimPressure, fluid)
scalar yDot() const
Return const access to rate of change of spherical deviation.
Definition: SprayParcelI.H:291
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
void calc(TrackCloudType &cloud, trackingData &td, const scalar dt)
Update parcel properties over the time interval.
Definition: SprayParcel.C:60
scalar user() const
Return const access to passive user scalar.
Definition: SprayParcelI.H:326
constantProperties()
Null constructor.
Definition: SprayParcelI.H:29
scalar sigma_
Liquid surface tension [N/m].
Definition: SprayParcel.H:146
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
static void writeFields(const CloudType &c, const CompositionType &compModel)
Write.
vector position0_
Injection position.
Definition: SprayParcel.H:143
PtrList< coordinateSystem > coordinates(solidRegions.size())
scalar injector() const
Return const access to injector id.
Definition: SprayParcelI.H:312
const vector & position0() const
Return const access to initial droplet position.
Definition: SprayParcelI.H:249
scalar yDot_
Rate of change of spherical deviation.
Definition: SprayParcel.H:161
void solveTABEq(TrackCloudType &cloud, trackingData &td, const scalar dt)
Solve the TAB equation.
Definition: SprayParcel.C:364
void setCellValues(TrackCloudType &cloud, trackingData &td)
Set cell values.
Definition: SprayParcel.C:35
scalar ms_
Stripped parcel mass due to breakup.
Definition: SprayParcel.H:167
void cellValueSourceCorrection(TrackCloudType &cloud, trackingData &td, const scalar dt)
Correct cell values using latest transfer information.
Definition: SprayParcel.C:47
scalar d0_
Initial droplet diameter.
Definition: SprayParcel.H:140
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
virtual autoPtr< particle > clone() const
Construct and return a (basic particle) clone.
Definition: SprayParcel.H:258
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
scalar tc_
Characteristic time (used in atomization and/or breakup model)
Definition: SprayParcel.H:164
volScalarField & p
scalar mu_
Liquid dynamic viscosity [Pa.s].
Definition: SprayParcel.H:149
TypeName("SprayParcel")
Runtime type information.
scalar y_
Spherical deviation.
Definition: SprayParcel.H:158
ParcelType::trackingData trackingData
Use base tracking data.
Definition: SprayParcel.H:130
scalar KHindex_
Index for KH Breakup.
Definition: SprayParcel.H:155
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:69
scalar sigma() const
Return const access to the liquid surface tension.
Definition: SprayParcelI.H:256
Namespace for OpenFOAM.
Reaching spray parcel, with added functionality for atomization and breakup.
Definition: SprayParcel.H:43
scalar T0
Definition: createFields.H:22