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