KinematicParcel.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 | Copyright (C) 2011-2017 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::KinematicParcel
26 
27 Description
28  Kinematic parcel class with rotational motion (as spherical
29  particles only) and one/two-way coupling with the continuous
30  phase.
31 
32  Sub-models include:
33  - drag
34  - turbulent dispersion
35  - wall interactions
36 
37 SourceFiles
38  KinematicParcelI.H
39  KinematicParcel.C
40  KinematicParcelIO.C
41 
42 \*---------------------------------------------------------------------------*/
43 
44 #ifndef KinematicParcel_H
45 #define KinematicParcel_H
46 
47 #include "particle.H"
48 #include "IOstream.H"
49 #include "autoPtr.H"
50 #include "interpolation.H"
51 #include "demandDrivenEntry.H"
52 
53 // #include "ParticleForceList.H" // TODO
54 
55 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
56 
57 namespace Foam
58 {
59 
60 template<class ParcelType>
61 class KinematicParcel;
62 
63 // Forward declaration of friend functions
64 
65 template<class ParcelType>
66 Ostream& operator<<
67 (
68  Ostream&,
70 );
71 
72 /*---------------------------------------------------------------------------*\
73  Class KinematicParcel Declaration
74 \*---------------------------------------------------------------------------*/
75 
76 template<class ParcelType>
77 class KinematicParcel
78 :
79  public ParcelType
80 {
81  // Private data
82 
83  //- Size in bytes of the fields
84  static const std::size_t sizeofFields_;
85 
86  //- Number of particle tracking attempts before we assume that it stalls
87  static label maxTrackAttempts;
88 
89 
90 public:
91 
92  //- Class to hold kinematic particle constant properties
93  class constantProperties
94  {
95  protected:
96 
97  // Protected data
98 
99  //- Constant properties dictionary
100  const dictionary dict_;
101 
102 
103  private:
104 
105  // Private data
106 
107  //- Parcel type id - used for post-processing to flag the type
108  // of parcels issued by this cloud
109  demandDrivenEntry<label> parcelTypeId_;
110 
111  //- Minimum density [kg/m3]
113 
114  //- Particle density [kg/m3] (constant)
116 
117  //- Minimum parcel mass [kg]
118  demandDrivenEntry<scalar> minParcelMass_;
119 
120 
121  public:
122 
123  // Constructors
124 
125  //- Null constructor
127 
128  //- Copy constructor
130 
131  //- Construct from dictionary
132  constantProperties(const dictionary& parentDict);
133 
134 
135  // Member functions
136 
137  //- Return const access to the constant properties dictionary
138  inline const dictionary& dict() const;
139 
140  //- Return const access to the parcel type id
141  inline label parcelTypeId() const;
142 
143  //- Return const access to the minimum density
144  inline scalar rhoMin() const;
145 
146  //- Return const access to the particle density
147  inline scalar rho0() const;
148 
149  //- Return const access to the minimum parcel mass
150  inline scalar minParcelMass() const;
151  };
152 
153 
154  template<class CloudType>
155  class TrackingData
156  :
157  public ParcelType::template TrackingData<CloudType>
158  {
159  public:
161  enum trackPart
162  {
163  tpVelocityHalfStep,
164  tpLinearTrack,
165  tpRotationalTrack
166  };
167 
168 
169  private:
170 
171  // Private data
172 
173  // Interpolators for continuous phase fields
174 
175  //- Density interpolator
176  autoPtr<interpolation<scalar>> rhoInterp_;
177 
178  //- Velocity interpolator
180 
181  //- Dynamic viscosity interpolator
183 
184 
185  //- Local gravitational or other body-force acceleration
186  const vector& g_;
187 
188  // label specifying which part of the integration
189  // algorithm is taking place
190  trackPart part_;
191 
192 
193  public:
194 
195  // Constructors
196 
197  //- Construct from components
198  inline TrackingData
199  (
200  CloudType& cloud,
201  trackPart part = tpLinearTrack
202  );
203 
204 
205  // Member functions
206 
207  //- Return conat access to the interpolator for continuous
208  // phase density field
209  inline const interpolation<scalar>& rhoInterp() const;
210 
211  //- Return conat access to the interpolator for continuous
212  // phase velocity field
213  inline const interpolation<vector>& UInterp() const;
214 
215  //- Return conat access to the interpolator for continuous
216  // phase dynamic viscosity field
217  inline const interpolation<scalar>& muInterp() const;
218 
219  // Return const access to the gravitational acceleration vector
220  inline const vector& g() const;
221 
222  //- Return the part of the tracking operation taking place
223  inline trackPart part() const;
224 
225  //- Return access to the part of the tracking operation taking place
226  inline trackPart& part();
227  };
228 
229 
230 protected:
231 
232  // Protected data
233 
234  // Parcel properties
235 
236  //- Active flag - tracking inactive when active = false
237  bool active_;
238 
239  //- Parcel type id
240  label typeId_;
241 
242  //- Number of particles in Parcel
243  scalar nParticle_;
244 
245  //- Diameter [m]
246  scalar d_;
247 
248  //- Target diameter [m]
249  scalar dTarget_;
250 
251  //- Velocity of Parcel [m/s]
252  vector U_;
253 
254  //- Density [kg/m3]
255  scalar rho_;
256 
257  //- Age [s]
258  scalar age_;
259 
260  //- Time spent in turbulent eddy [s]
261  scalar tTurb_;
262 
263  //- Turbulent velocity fluctuation [m/s]
264  vector UTurb_;
265 
266 
267  // Cell-based quantities
268 
269  //- Density [kg/m3]
270  scalar rhoc_;
271 
272  //- Velocity [m/s]
273  vector Uc_;
274 
275  //- Viscosity [Pa.s]
276  scalar muc_;
277 
278 
279  // Protected Member Functions
280 
281  //- Calculate new particle velocity
282  template<class TrackData>
283  const vector calcVelocity
284  (
285  TrackData& td,
286  const scalar dt, // timestep
287  const label celli, // owner cell
288  const scalar Re, // Reynolds number
289  const scalar mu, // local carrier viscosity
290  const scalar mass, // mass
291  const vector& Su, // explicit particle momentum source
292  vector& dUTrans, // momentum transfer to carrier
293  scalar& Spu // linearised drag coefficient
294  ) const;
295 
296 
297 public:
298 
299  // Static data members
300 
301  //- Runtime type information
302  TypeName("KinematicParcel");
303 
304  //- String representation of properties
306  (
307  ParcelType,
308  " active"
309  + " typeId"
310  + " nParticle"
311  + " d"
312  + " dTarget "
313  + " (Ux Uy Uz)"
314  + " rho"
315  + " age"
316  + " tTurb"
317  + " (UTurbx UTurby UTurbz)"
318  );
319 
320 
321  // Constructors
322 
323  //- Construct from mesh, coordinates and topology
324  // Other properties initialised as null
325  inline KinematicParcel
326  (
327  const polyMesh& mesh,
328  const barycentric& coordinates,
329  const label celli,
330  const label tetFacei,
331  const label tetPti
332  );
333 
334  //- Construct from a position and a cell, searching for the rest of the
335  // required topology. Other properties are initialised as null.
336  inline KinematicParcel
337  (
338  const polyMesh& mesh,
339  const vector& position,
340  const label celli
341  );
342 
343  //- Construct from components
344  inline KinematicParcel
345  (
346  const polyMesh& mesh,
347  const barycentric& coordinates,
348  const label celli,
349  const label tetFacei,
350  const label tetPti,
351  const label typeId,
352  const scalar nParticle0,
353  const scalar d0,
354  const scalar dTarget0,
355  const vector& torque0,
356  const constantProperties& constProps
357  );
358 
359  //- Construct from Istream
361  (
362  const polyMesh& mesh,
363  Istream& is,
364  bool readFields = true
365  );
366 
367  //- Construct as a copy
369 
370  //- Construct as a copy
372 
373  //- Construct and return a (basic particle) clone
374  virtual autoPtr<particle> clone() const
375  {
376  return autoPtr<particle>(new KinematicParcel(*this));
377  }
378 
379  //- Construct and return a (basic particle) clone
380  virtual autoPtr<particle> clone(const polyMesh& mesh) const
381  {
382  return autoPtr<particle>(new KinematicParcel(*this, mesh));
383  }
384 
385  //- Factory class to read-construct particles used for
386  // parallel transfer
387  class iNew
388  {
389  const polyMesh& mesh_;
390 
391  public:
393  iNew(const polyMesh& mesh)
394  :
395  mesh_(mesh)
396  {}
398  autoPtr<KinematicParcel<ParcelType>> operator()(Istream& is) const
399  {
401  (
402  new KinematicParcel<ParcelType>(mesh_, is, true)
403  );
404  }
405  };
406 
407 
408  // Member Functions
409 
410  // Access
411 
412  //- Return const access to active flag
413  inline bool active() const;
414 
415  //- Return const access to type id
416  inline label typeId() const;
417 
418  //- Return const access to number of particles
419  inline scalar nParticle() const;
420 
421  //- Return const access to diameter
422  inline scalar d() const;
423 
424  //- Return const access to target diameter
425  inline scalar dTarget() const;
426 
427  //- Return const access to velocity
428  inline const vector& U() const;
429 
430  //- Return const access to density
431  inline scalar rho() const;
432 
433  //- Return const access to the age
434  inline scalar age() const;
435 
436  //- Return const access to time spent in turbulent eddy
437  inline scalar tTurb() const;
438 
439  //- Return const access to turbulent velocity fluctuation
440  inline const vector& UTurb() const;
441 
442  //- Return const access to carrier density [kg/m3]
443  inline scalar rhoc() const;
444 
445  //- Return const access to carrier velocity [m/s]
446  inline const vector& Uc() const;
447 
448  //- Return const access to carrier viscosity [Pa.s]
449  inline scalar muc() const;
450 
451 
452  // Edit
453 
454  //- Return const access to active flag
455  inline bool& active();
456 
457  //- Return access to type id
458  inline label& typeId();
459 
460  //- Return access to number of particles
461  inline scalar& nParticle();
462 
463  //- Return access to diameter
464  inline scalar& d();
465 
466  //- Return access to target diameter
467  inline scalar& dTarget();
468 
469  //- Return access to velocity
470  inline vector& U();
471 
472  //- Return access to density
473  inline scalar& rho();
474 
475  //- Return access to the age
476  inline scalar& age();
477 
478  //- Return access to time spent in turbulent eddy
479  inline scalar& tTurb();
480 
481  //- Return access to turbulent velocity fluctuation
482  inline vector& UTurb();
483 
484 
485  // Helper functions
486 
487  //- Cell owner mass
488  inline scalar massCell(const label celli) const;
489 
490  //- Particle mass
491  inline scalar mass() const;
492 
493  //- Particle moment of inertia around diameter axis
494  inline scalar momentOfInertia() const;
495 
496  //- Particle volume
497  inline scalar volume() const;
498 
499  //- Particle volume for a given diameter
500  inline static scalar volume(const scalar d);
501 
502  //- Particle projected area
503  inline scalar areaP() const;
504 
505  //- Projected area for given diameter
506  inline static scalar areaP(const scalar d);
507 
508  //- Particle surface area
509  inline scalar areaS() const;
510 
511  //- Surface area for given diameter
512  inline static scalar areaS(const scalar d);
513 
514  //- Reynolds number
515  inline scalar Re
516  (
517  const vector& U, // particle velocity
518  const scalar d, // particle diameter
519  const scalar rhoc, // carrier density
520  const scalar muc // carrier dynamic viscosity
521  ) const;
522 
523  //- Weber number
524  inline scalar We
525  (
526  const vector& U, // particle velocity
527  const scalar d, // particle diameter
528  const scalar rhoc, // carrier density
529  const scalar sigma // particle surface tension
530  ) const;
531 
532  //- Eotvos number
533  inline scalar Eo
534  (
535  const vector& a, // acceleration
536  const scalar d, // particle diameter
537  const scalar sigma // particle surface tension
538  ) const;
539 
540 
541  // Main calculation loop
542 
543  //- Set cell values
544  template<class TrackData>
545  void setCellValues
546  (
547  TrackData& td,
548  const scalar dt,
549  const label celli
550  );
551 
552  //- Correct cell values using latest transfer information
553  template<class TrackData>
555  (
556  TrackData& td,
557  const scalar dt,
558  const label celli
559  );
560 
561  //- Update parcel properties over the time interval
562  template<class TrackData>
563  void calc
564  (
565  TrackData& td,
566  const scalar dt,
567  const label celli
568  );
569 
570 
571  // Tracking
572 
573  //- Move the parcel
574  template<class TrackData>
575  bool move(TrackData& td, const scalar trackTime);
576 
577 
578  // Patch interactions
579 
580  //- Overridable function to handle the particle hitting a face
581  // without trackData
582  void hitFace(int& td);
583 
584  //- Overridable function to handle the particle hitting a face
585  template<class TrackData>
586  void hitFace(TrackData& td);
587 
588  //- Overridable function to handle the particle hitting a patch
589  // Executed before other patch-hitting functions
590  template<class TrackData>
591  bool hitPatch
592  (
593  const polyPatch& p,
594  TrackData& td,
595  const label patchi,
596  const scalar trackFraction,
597  const tetIndices& tetIs
598  );
599 
600  //- Overridable function to handle the particle hitting a
601  // processorPatch
602  template<class TrackData>
603  void hitProcessorPatch
604  (
605  const processorPolyPatch&,
606  TrackData& td
607  );
608 
609  //- Overridable function to handle the particle hitting a wallPatch
610  template<class TrackData>
611  void hitWallPatch
612  (
613  const wallPolyPatch&,
614  TrackData& td,
615  const tetIndices&
616  );
617 
618  //- Overridable function to handle the particle hitting a polyPatch
619  template<class TrackData>
620  void hitPatch
621  (
622  const polyPatch&,
623  TrackData& td
624  );
625 
626  //- Transform the physical properties of the particle
627  // according to the given transformation tensor
628  virtual void transformProperties(const tensor& T);
629 
630  //- Transform the physical properties of the particle
631  // according to the given separation vector
632  virtual void transformProperties(const vector& separation);
633 
634  //- The nearest distance to a wall that the particle can be
635  // in the n direction
636  virtual scalar wallImpactDistance(const vector& n) const;
637 
638 
639  // I-O
640 
641  //- Read
642  template<class CloudType>
643  static void readFields(CloudType& c);
644 
645  //- Write
646  template<class CloudType>
647  static void writeFields(const CloudType& c);
648 
649 
650  // Ostream Operator
651 
652  friend Ostream& operator<< <ParcelType>
653  (
654  Ostream&,
656  );
657 };
658 
659 
660 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
661 
662 } // End namespace Foam
663 
664 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
665 
666 #include "KinematicParcelI.H"
668 
669 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
670 
671 #ifdef NoRepository
672  #include "KinematicParcel.C"
673 #endif
674 
675 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
676 
677 #endif
678 
679 // ************************************************************************* //
const vector & U() const
Return const access to velocity.
scalar tTurb() const
Return const access to time spent in turbulent eddy.
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
scalar rhoMin() const
Return const access to the minimum density.
static void readFields(CloudType &c)
Read.
bool hitPatch(const polyPatch &p, TrackData &td, const label patchi, const scalar trackFraction, const tetIndices &tetIs)
Overridable function to handle the particle hitting a patch.
zeroField Su
Definition: alphaSuSp.H:1
const vector & UTurb() const
Return const access to turbulent velocity fluctuation.
virtual scalar wallImpactDistance(const vector &n) const
The nearest distance to a wall that the particle can be.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
Kinematic parcel class with rotational motion (as spherical particles only) and one/two-way coupling ...
scalar muc() const
Return const access to carrier viscosity [Pa.s].
TypeName("KinematicParcel")
Runtime type information.
scalar nParticle() const
Return const access to number of particles.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
Class to hold kinematic particle constant properties.
scalar mass() const
Particle mass.
scalar d() const
Return const access to diameter.
scalar tTurb_
Time spent in turbulent eddy [s].
scalar age_
Age [s].
scalar rho() const
Return const access to density.
bool cp(const fileName &src, const fileName &dst, const bool followLink=true)
Copy, recursively if necessary, the source to the destination.
Definition: POSIX.C:739
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
scalar rho0() const
Return const access to the particle density.
scalar areaP() const
Particle projected area.
Factory class to read-construct particles used for.
KinematicParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from mesh, coordinates and topology.
scalar dTarget() const
Return const access to target diameter.
vector Uc_
Velocity [m/s].
Neighbour processor patch.
bool move(TrackData &td, const scalar trackTime)
Move the parcel.
label parcelTypeId() const
Return const access to the parcel type id.
scalar dTarget_
Target diameter [m].
scalar age() const
Return const access to the age.
dynamicFvMesh & mesh
scalar Eo(const vector &a, const scalar d, const scalar sigma) const
Eotvos number.
scalar rho_
Density [kg/m3].
virtual void transformProperties(const tensor &T)
Transform the physical properties of the particle.
scalar momentOfInertia() const
Particle moment of inertia around diameter axis.
scalar muc_
Viscosity [Pa.s].
void hitFace(int &td)
Overridable function to handle the particle hitting a face.
bool active_
Active flag - tracking inactive when active = false.
const dictionary dict_
Constant properties dictionary.
scalar massCell(const label celli) const
Cell owner mass.
A cloud is a collection of lagrangian particles.
Definition: cloud.H:51
Foam::wallPolyPatch.
Definition: wallPolyPatch.H:48
label typeId_
Parcel type id.
vector UTurb_
Turbulent velocity fluctuation [m/s].
Storage and named access for the indices of a tet which is part of the decomposition of a cell...
Definition: tetIndices.H:81
bool active() const
Return const access to active flag.
const vector & Uc() const
Return const access to carrier velocity [m/s].
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
const dimensionedVector & g
scalar nParticle_
Number of particles in Parcel.
void cellValueSourceCorrection(TrackData &td, const scalar dt, const label celli)
Correct cell values using latest transfer information.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
scalar minParcelMass() const
Return const access to the minimum parcel mass.
void calc(TrackData &td, const scalar dt, const label celli)
Update parcel properties over the time interval.
scalar rhoc() const
Return const access to carrier density [kg/m3].
const dimensionedScalar mu
Atomic mass unit.
scalar rhoc_
Density [kg/m3].
vector U_
Velocity of Parcel [m/s].
label patchi
virtual autoPtr< particle > clone() const
Construct and return a (basic particle) clone.
PtrList< coordinateSystem > coordinates(solidRegions.size())
const dimensionedScalar c
Speed of light in a vacuum.
scalar volume() const
Particle volume.
void setCellValues(TrackData &td, const scalar dt, const label celli)
Set cell values.
static void writeFields(const CloudType &c)
Write.
void hitWallPatch(const wallPolyPatch &, TrackData &td, const tetIndices &)
Overridable function to handle the particle hitting a wallPatch.
label n
scalar We(const vector &U, const scalar d, const scalar rhoc, const scalar sigma) const
Weber number.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
void hitProcessorPatch(const processorPolyPatch &, TrackData &td)
Overridable function to handle the particle hitting a.
volScalarField & p
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
const vector calcVelocity(TrackData &td, const scalar dt, const label celli, const scalar Re, const scalar mu, const scalar mass, const vector &Su, vector &dUTrans, scalar &Spu) const
Calculate new particle velocity.
label typeId() const
Return const access to type id.
AddToPropertyList(ParcelType, " active"+" typeId"+" nParticle"+" d"+" dTarget "+" (Ux Uy Uz)"+" rho"+" age"+" tTurb"+" (UTurbx UTurby UTurbz)")
String representation of properties.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:69
scalar Re(const vector &U, const scalar d, const scalar rhoc, const scalar muc) const
Reynolds number.
scalar areaS() const
Particle surface area.
scalar d_
Diameter [m].
Namespace for OpenFOAM.
const dictionary & dict() const
Return const access to the constant properties dictionary.