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  );
213 
214  //- Construct from Istream
215  SprayParcel(Istream& is, bool readFields = true);
216 
217  //- Construct as a copy
218  SprayParcel(const SprayParcel& p);
219 
220  //- Construct and return a clone
221  virtual autoPtr<particle> clone() const
222  {
223  return autoPtr<particle>(new SprayParcel<ParcelType>(*this));
224  }
225 
226  //- Construct from Istream and return
227  static autoPtr<SprayParcel> New(Istream& is)
228  {
229  return autoPtr<SprayParcel>(new SprayParcel(is));
230  }
231 
232 
233  // Member Functions
234 
235  // Access
236 
237  //- Return const access to initial droplet diameter
238  inline scalar d0() const;
239 
240  //- Return const access to initial mass [kg]
241  inline scalar mass0() const;
242 
243  //- Return const access to initial droplet position
244  inline const vector& position0() const;
245 
246  //- Return const access to the liquid surface tension
247  inline scalar sigma() const;
248 
249  //- Return const access to the liquid dynamic viscosity
250  inline scalar mu() const;
251 
252  //- Return const access to liquid core
253  inline scalar liquidCore() const;
254 
255  //- Return const access to Kelvin-Helmholtz breakup index
256  inline scalar KHindex() const;
257 
258  //- Return const access to spherical deviation
259  inline scalar y() const;
260 
261  //- Return const access to rate of change of spherical deviation
262  inline scalar yDot() const;
263 
264  //- Return const access to atomisation characteristic time
265  inline scalar tc() const;
266 
267  //- Return const access to stripped parcel mass
268  inline scalar ms() const;
269 
270  //- Return const access to injector id
271  inline label injector() const;
272 
273  //- Return const access to momentum relaxation time
274  inline scalar tMom() const;
275 
276 
277  // Edit
278 
279  //- Return access to initial droplet diameter
280  inline scalar& d0();
281 
282  //- Return access to initial mass [kg]
283  inline scalar& mass0();
284 
285  //- Return access to initial droplet position
286  inline vector& position0();
287 
288  //- Return access to the liquid surface tension
289  inline scalar& sigma();
290 
291  //- Return access to the liquid dynamic viscosity
292  inline scalar& mu();
293 
294  //- Return access to liquid core
295  inline scalar& liquidCore();
296 
297  //- Return access to Kelvin-Helmholtz breakup index
298  inline scalar& KHindex();
299 
300  //- Return access to spherical deviation
301  inline scalar& y();
302 
303  //- Return access to rate of change of spherical deviation
304  inline scalar& yDot();
305 
306  //- Return access to atomisation characteristic time
307  inline scalar& tc();
308 
309  //- Return access to stripped parcel mass
310  inline scalar& ms();
311 
312  //- Return access to injector id
313  inline label& injector();
314 
315  //- Return access to momentum relaxation time
316  inline scalar& tMom();
317 
318 
319  // Main calculation loop
320 
321  //- Set cell values
322  template<class TrackCloudType>
323  void setCellValues(TrackCloudType& cloud, trackingData& td);
324 
325  //- Correct parcel properties according to atomisation model
326  template<class TrackCloudType>
327  void calcAtomisation
328  (
329  TrackCloudType& cloud,
330  trackingData& td,
331  const scalar dt
332  );
333 
334  //- Correct parcel properties according to breakup model
335  template<class TrackCloudType>
336  void calcBreakup
337  (
338  TrackCloudType& cloud,
339  trackingData& td,
340  const scalar dt
341  );
342 
343  //- Correct cell values using latest transfer information
344  template<class TrackCloudType>
346  (
347  TrackCloudType& cloud,
348  trackingData& td,
349  const scalar dt
350  );
351 
352  //- Correct surface values due to emitted species
353  template<class TrackCloudType>
355  (
356  TrackCloudType& cloud,
357  trackingData& td,
358  const scalar T,
359  const scalarField& Cs,
360  scalar& rhos,
361  scalar& mus,
362  scalar& Pr,
363  scalar& kappa
364  );
365 
366  //- Update parcel properties over the time interval
367  template<class TrackCloudType>
368  void calc
369  (
370  TrackCloudType& cloud,
371  trackingData& td,
372  const scalar dt
373  );
374 
375  //- Calculate the chi-factor for flash-boiling for the
376  // atomisation model
377  template<class TrackCloudType>
378  scalar chi
379  (
380  TrackCloudType& cloud,
381  trackingData& td,
382  const scalarField& X
383  ) const;
384 
385  //- Solve the TAB equation
386  template<class TrackCloudType>
387  void solveTABEq
388  (
389  TrackCloudType& cloud,
390  trackingData& td,
391  const scalar dt
392  );
393 
394 
395  // I-O
396 
397  //- Read
398  template<class CloudType, class CompositionType>
399  static void readFields
400  (
401  CloudType& c,
402  const CompositionType& compModel
403  );
404 
405  //- Read - no composition
406  template<class CloudType>
407  static void readFields(CloudType& c);
408 
409  //- Write
410  template<class CloudType, class CompositionType>
411  static void writeFields
412  (
413  const CloudType& c,
414  const CompositionType& compModel
415  );
416 
417  //- Write - composition supplied
418  template<class CloudType>
419  static void writeFields(const CloudType& c);
420 
421 
422  // Ostream Operator
423 
424  friend Ostream& operator<< <ParcelType>
425  (
426  Ostream&,
428  );
429 };
430 
431 
432 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
433 
434 } // End namespace Foam
435 
436 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
437 
438 #include "SprayParcelI.H"
439 
440 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
441 
442 #ifdef NoRepository
443  #include "SprayParcel.C"
444 #endif
445 
446 
447 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
448 
449 #endif
450 
451 // ************************************************************************* //
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:79
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:173
scalar sigma0() const
Return const access to the initial surface tension.
Definition: SprayParcelI.H:165
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:217
scalar d0_
Initial droplet diameter [m].
Definition: SprayParcel.H:149
label injector() const
Return const access to injector id.
Definition: SprayParcelI.H:259
scalar tMom() const
Return const access to momentum relaxation time.
Definition: SprayParcelI.H:266
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:210
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:238
scalar d0() const
Return const access to initial droplet diameter.
Definition: SprayParcelI.H:182
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:220
static autoPtr< SprayParcel > New(Istream &is)
Construct from Istream and return.
Definition: SprayParcel.H:226
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:224
scalar tc() const
Return const access to atomisation characteristic time.
Definition: SprayParcelI.H:245
scalar y_
Spherical deviation.
Definition: SprayParcel.H:170
scalar mass0() const
Return const access to initial mass [kg].
Definition: SprayParcelI.H:189
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:252
scalar y() const
Return const access to spherical deviation.
Definition: SprayParcelI.H:231
scalar mass0_
Initial mass [kg].
Definition: SprayParcel.H:152
scalar sigma() const
Return const access to the liquid surface tension.
Definition: SprayParcelI.H:203
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:196
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:160
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