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-2023 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 [m]
150  scalar d0_;
151 
152  //- Initial mass [kg]
153  scalar mass0_;
154 
155  //- Injection position
157 
158  //- Liquid surface tension [N/m]
159  scalar sigma_;
160 
161  //- Liquid dynamic viscosity [Pa.s]
162  scalar mu_;
163 
164  //- Part of liquid core ( >0.5=liquid, <0.5=droplet )
165  scalar liquidCore_;
166 
167  //- Index for KH Breakup
168  scalar KHindex_;
169 
170  //- Spherical deviation
171  scalar y_;
172 
173  //- Rate of change of spherical deviation
174  scalar yDot_;
175 
176  //- Characteristic time (used in atomisation and/or breakup model)
177  scalar tc_;
178 
179  //- Stripped parcel mass due to breakup
180  scalar ms_;
181 
182  //- Injector id
184 
185  //- Momentum relaxation time (needed for calculating parcel acc.)
186  scalar tMom_;
187 
188 
189 public:
190 
191  // Constructors
192 
193  //- Construct from mesh, coordinates and topology
194  // Other properties initialised as null
195  inline SprayParcel
196  (
197  const polyMesh& mesh,
198  const barycentric& coordinates,
199  const label celli,
200  const label tetFacei,
201  const label tetPti,
202  const label facei
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  label& nLocateBoundaryHits
213  );
214 
215  //- Construct from Istream
216  SprayParcel(Istream& is, bool readFields = true);
217 
218  //- Construct as a copy
219  SprayParcel(const SprayParcel& p);
220 
221  //- Construct and return a clone
222  virtual autoPtr<particle> clone() const
223  {
224  return autoPtr<particle>(new SprayParcel<ParcelType>(*this));
225  }
226 
227  //- Construct from Istream and return
228  static autoPtr<SprayParcel> New(Istream& is)
229  {
230  return autoPtr<SprayParcel>(new SprayParcel(is));
231  }
232 
233 
234  // Member Functions
235 
236  // Access
237 
238  //- Return const access to initial droplet diameter
239  inline scalar d0() const;
240 
241  //- Return const access to initial mass [kg]
242  inline scalar mass0() const;
243 
244  //- Return const access to initial droplet position
245  inline const vector& position0() const;
246 
247  //- Return const access to the liquid surface tension
248  inline scalar sigma() const;
249 
250  //- Return const access to the liquid dynamic viscosity
251  inline scalar mu() const;
252 
253  //- Return const access to liquid core
254  inline scalar liquidCore() const;
255 
256  //- Return const access to Kelvin-Helmholtz breakup index
257  inline scalar KHindex() const;
258 
259  //- Return const access to spherical deviation
260  inline scalar y() const;
261 
262  //- Return const access to rate of change of spherical deviation
263  inline scalar yDot() const;
264 
265  //- Return const access to atomisation characteristic time
266  inline scalar tc() const;
267 
268  //- Return const access to stripped parcel mass
269  inline scalar ms() const;
270 
271  //- Return const access to injector id
272  inline label injector() const;
273 
274  //- Return const access to momentum relaxation time
275  inline scalar tMom() const;
276 
277 
278  // Edit
279 
280  //- Return access to initial droplet diameter
281  inline scalar& d0();
282 
283  //- Return access to initial mass [kg]
284  inline scalar& mass0();
285 
286  //- Return access to initial droplet position
287  inline vector& position0();
288 
289  //- Return access to the liquid surface tension
290  inline scalar& sigma();
291 
292  //- Return access to the liquid dynamic viscosity
293  inline scalar& mu();
294 
295  //- Return access to liquid core
296  inline scalar& liquidCore();
297 
298  //- Return access to Kelvin-Helmholtz breakup index
299  inline scalar& KHindex();
300 
301  //- Return access to spherical deviation
302  inline scalar& y();
303 
304  //- Return access to rate of change of spherical deviation
305  inline scalar& yDot();
306 
307  //- Return access to atomisation characteristic time
308  inline scalar& tc();
309 
310  //- Return access to stripped parcel mass
311  inline scalar& ms();
312 
313  //- Return access to injector id
314  inline label& injector();
315 
316  //- Return access to momentum relaxation time
317  inline scalar& tMom();
318 
319 
320  // Main calculation loop
321 
322  //- Set cell values
323  template<class TrackCloudType>
324  void setCellValues(TrackCloudType& cloud, trackingData& td);
325 
326  //- Correct parcel properties according to atomisation model
327  template<class TrackCloudType>
328  void calcAtomisation
329  (
330  TrackCloudType& cloud,
331  trackingData& td,
332  const scalar dt
333  );
334 
335  //- Correct parcel properties according to breakup model
336  template<class TrackCloudType>
337  void calcBreakup
338  (
339  TrackCloudType& cloud,
340  trackingData& td,
341  const scalar dt
342  );
343 
344  //- Correct cell values using latest transfer information
345  template<class TrackCloudType>
347  (
348  TrackCloudType& cloud,
349  trackingData& td,
350  const scalar dt
351  );
352 
353  //- Correct surface values due to emitted species
354  template<class TrackCloudType>
356  (
357  TrackCloudType& cloud,
358  trackingData& td,
359  const scalar T,
360  const scalarField& Cs,
361  scalar& rhos,
362  scalar& mus,
363  scalar& Pr,
364  scalar& kappa
365  );
366 
367  //- Update parcel properties over the time interval
368  template<class TrackCloudType>
369  void calc
370  (
371  TrackCloudType& cloud,
372  trackingData& td,
373  const scalar dt
374  );
375 
376  //- Calculate the chi-factor for flash-boiling for the
377  // atomisation model
378  template<class TrackCloudType>
379  scalar chi
380  (
381  TrackCloudType& cloud,
382  trackingData& td,
383  const scalarField& X
384  ) const;
385 
386  //- Solve the TAB equation
387  template<class TrackCloudType>
388  void solveTABEq
389  (
390  TrackCloudType& cloud,
391  trackingData& td,
392  const scalar dt
393  );
394 
395 
396  // I-O
397 
398  //- Read
399  template<class CloudType, class CompositionType>
400  static void readFields
401  (
402  CloudType& c,
403  const CompositionType& compModel
404  );
405 
406  //- Read - no composition
407  template<class CloudType>
408  static void readFields(CloudType& c);
409 
410  //- Write
411  template<class CloudType, class CompositionType>
412  static void writeFields
413  (
414  const CloudType& c,
415  const CompositionType& compModel
416  );
417 
418  //- Write - composition supplied
419  template<class CloudType>
420  static void writeFields(const CloudType& c);
421 
422 
423  // Ostream Operator
424 
425  friend Ostream& operator<< <ParcelType>
426  (
427  Ostream&,
429  );
430 };
431 
432 
433 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
434 
435 } // End namespace Foam
436 
437 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
438 
439 #include "SprayParcelI.H"
440 
441 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
442 
443 #ifdef NoRepository
444  #include "SprayParcel.C"
445 #endif
446 
447 
448 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
449 
450 #endif
451 
452 // ************************************************************************* //
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:80
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:57
Class to hold reacting particle constant properties.
Definition: SprayParcel.H:82
scalar mu0() const
Return const access to the initial dynamic viscosity.
Definition: SprayParcelI.H:174
scalar sigma0() const
Return const access to the initial surface tension.
Definition: SprayParcelI.H:166
Reaching spray parcel, with added functionality for atomisation and breakup.
Definition: SprayParcel.H:69
scalar chi(TrackCloudType &cloud, trackingData &td, const scalarField &X) const
Calculate the chi-factor for flash-boiling for the.
SprayParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti, const label facei)
Construct from mesh, coordinates and topology.
Definition: SprayParcelI.H:110
scalar liquidCore() const
Return const access to liquid core.
Definition: SprayParcelI.H:218
scalar d0_
Initial droplet diameter [m].
Definition: SprayParcel.H:149
label injector() const
Return const access to injector id.
Definition: SprayParcelI.H:260
scalar tMom() const
Return const access to momentum relaxation time.
Definition: SprayParcelI.H:267
vector position0_
Injection position.
Definition: SprayParcel.H:155
void cellValueSourceCorrection(TrackCloudType &cloud, trackingData &td, const scalar dt)
Correct cell values using latest transfer information.
Definition: SprayParcel.C:48
void calcAtomisation(TrackCloudType &cloud, trackingData &td, const scalar dt)
Correct parcel properties according to atomisation model.
Definition: SprayParcel.C:150
label injector_
Injector id.
Definition: SprayParcel.H:182
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:211
scalar mu_
Liquid dynamic viscosity [Pa.s].
Definition: SprayParcel.H:161
scalar sigma_
Liquid surface tension [N/m].
Definition: SprayParcel.H:158
static void writeFields(const CloudType &c, const CompositionType &compModel)
Write.
scalar yDot() const
Return const access to rate of change of spherical deviation.
Definition: SprayParcelI.H:239
scalar d0() const
Return const access to initial droplet diameter.
Definition: SprayParcelI.H:183
ParcelType::trackingData trackingData
Use base tracking data.
Definition: SprayParcel.H:139
void solveTABEq(TrackCloudType &cloud, trackingData &td, const scalar dt)
Solve the TAB equation.
Definition: SprayParcel.C:369
virtual autoPtr< particle > clone() const
Construct and return a clone.
Definition: SprayParcel.H:221
static autoPtr< SprayParcel > New(Istream &is)
Construct from Istream and return.
Definition: SprayParcel.H:227
scalar tMom_
Momentum relaxation time (needed for calculating parcel acc.)
Definition: SprayParcel.H:185
scalar KHindex() const
Return const access to Kelvin-Helmholtz breakup index.
Definition: SprayParcelI.H:225
scalar tc() const
Return const access to atomisation characteristic time.
Definition: SprayParcelI.H:246
scalar y_
Spherical deviation.
Definition: SprayParcel.H:170
scalar mass0() const
Return const access to initial mass [kg].
Definition: SprayParcelI.H:190
scalar tc_
Characteristic time (used in atomisation and/or breakup model)
Definition: SprayParcel.H:176
scalar KHindex_
Index for KH Breakup.
Definition: SprayParcel.H:167
void setCellValues(TrackCloudType &cloud, trackingData &td)
Set cell values.
Definition: SprayParcel.C:36
void calcBreakup(TrackCloudType &cloud, trackingData &td, const scalar dt)
Correct parcel properties according to breakup model.
Definition: SprayParcel.C:221
scalar ms() const
Return const access to stripped parcel mass.
Definition: SprayParcelI.H:253
scalar y() const
Return const access to spherical deviation.
Definition: SprayParcelI.H:232
scalar mass0_
Initial mass [kg].
Definition: SprayParcel.H:152
scalar sigma() const
Return const access to the liquid surface tension.
Definition: SprayParcelI.H:204
scalar ms_
Stripped parcel mass due to breakup.
Definition: SprayParcel.H:179
scalar yDot_
Rate of change of spherical deviation.
Definition: SprayParcel.H:173
static void readFields(CloudType &c, const CompositionType &compModel)
Read.
scalar liquidCore_
Part of liquid core ( >0.5=liquid, <0.5=droplet )
Definition: SprayParcel.H:164
const vector & position0() const
Return const access to initial droplet position.
Definition: SprayParcelI.H:197
void calc(TrackCloudType &cloud, trackingData &td, const scalar dt)
Update parcel properties over the time interval.
Definition: SprayParcel.C:61
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:61
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A cloud is a collection of lagrangian particles.
Definition: cloud.H:55
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
const dimensionedScalar epsilon0
Electric constant: default SI units: [F/m].
const dimensionedScalar c
Speed of light in a vacuum.
Namespace for OpenFOAM.
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
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
TemplateName(FvFaceCellWave)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
scalar rho0
volScalarField & p
scalar T0
Definition: createFields.H:22