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-2015 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
110 
111  //- Minimum density [kg/m3]
113 
114  //- Particle density [kg/m3] (constant)
116 
117  //- Minimum parcel mass [kg]
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:
160 
161  enum trackPart
162  {
166  };
167 
168 
169  private:
170 
171  // Private data
172 
173  // Interpolators for continuous phase fields
174 
175  //- Density interpolator
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
191 
192 
193  public:
194 
195  // Constructors
196 
197  //- Construct from components
198  inline TrackingData
199  (
200  CloudType& cloud,
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 owner, position, and cloud owner
324  // Other properties initialised as null
325  inline KinematicParcel
326  (
327  const polyMesh& mesh,
328  const vector& position,
329  const label cellI,
330  const label tetFaceI,
331  const label tetPtI
332  );
333 
334  //- Construct from components
335  inline KinematicParcel
336  (
337  const polyMesh& mesh,
338  const vector& position,
339  const label cellI,
340  const label tetFaceI,
341  const label tetPtI,
342  const label typeId,
343  const scalar nParticle0,
344  const scalar d0,
345  const scalar dTarget0,
346  const vector& torque0,
347  const constantProperties& constProps
348  );
349 
350  //- Construct from Istream
352  (
353  const polyMesh& mesh,
354  Istream& is,
355  bool readFields = true
356  );
357 
358  //- Construct as a copy
360 
361  //- Construct as a copy
363 
364  //- Construct and return a (basic particle) clone
365  virtual autoPtr<particle> clone() const
366  {
367  return autoPtr<particle>(new KinematicParcel(*this));
368  }
369 
370  //- Construct and return a (basic particle) clone
371  virtual autoPtr<particle> clone(const polyMesh& mesh) const
372  {
373  return autoPtr<particle>(new KinematicParcel(*this, mesh));
374  }
375 
376  //- Factory class to read-construct particles used for
377  // parallel transfer
378  class iNew
379  {
380  const polyMesh& mesh_;
381 
382  public:
383 
384  iNew(const polyMesh& mesh)
385  :
386  mesh_(mesh)
387  {}
388 
390  {
392  (
393  new KinematicParcel<ParcelType>(mesh_, is, true)
394  );
395  }
396  };
397 
398 
399  // Member Functions
400 
401  // Access
402 
403  //- Return const access to active flag
404  inline bool active() const;
405 
406  //- Return const access to type id
407  inline label typeId() const;
408 
409  //- Return const access to number of particles
410  inline scalar nParticle() const;
411 
412  //- Return const access to diameter
413  inline scalar d() const;
414 
415  //- Return const access to target diameter
416  inline scalar dTarget() const;
417 
418  //- Return const access to velocity
419  inline const vector& U() const;
420 
421  //- Return const access to density
422  inline scalar rho() const;
423 
424  //- Return const access to the age
425  inline scalar age() const;
426 
427  //- Return const access to time spent in turbulent eddy
428  inline scalar tTurb() const;
429 
430  //- Return const access to turbulent velocity fluctuation
431  inline const vector& UTurb() const;
432 
433  //- Return const access to carrier density [kg/m3]
434  inline scalar rhoc() const;
435 
436  //- Return const access to carrier velocity [m/s]
437  inline const vector& Uc() const;
438 
439  //- Return const access to carrier viscosity [Pa.s]
440  inline scalar muc() const;
441 
442 
443  // Edit
444 
445  //- Return const access to active flag
446  inline bool& active();
447 
448  //- Return access to type id
449  inline label& typeId();
450 
451  //- Return access to number of particles
452  inline scalar& nParticle();
453 
454  //- Return access to diameter
455  inline scalar& d();
456 
457  //- Return access to target diameter
458  inline scalar& dTarget();
459 
460  //- Return access to velocity
461  inline vector& U();
462 
463  //- Return access to density
464  inline scalar& rho();
465 
466  //- Return access to the age
467  inline scalar& age();
468 
469  //- Return access to time spent in turbulent eddy
470  inline scalar& tTurb();
471 
472  //- Return access to turbulent velocity fluctuation
473  inline vector& UTurb();
474 
475 
476  // Helper functions
477 
478  //- Return the index of the face used in the interpolation routine
479  inline label faceInterpolation() const;
480 
481  //- Cell owner mass
482  inline scalar massCell(const label cellI) const;
483 
484  //- Particle mass
485  inline scalar mass() const;
486 
487  //- Particle moment of inertia around diameter axis
488  inline scalar momentOfInertia() const;
489 
490  //- Particle volume
491  inline scalar volume() const;
492 
493  //- Particle volume for a given diameter
494  inline static scalar volume(const scalar d);
495 
496  //- Particle projected area
497  inline scalar areaP() const;
498 
499  //- Projected area for given diameter
500  inline static scalar areaP(const scalar d);
501 
502  //- Particle surface area
503  inline scalar areaS() const;
504 
505  //- Surface area for given diameter
506  inline static scalar areaS(const scalar d);
507 
508  //- Reynolds number
509  inline scalar Re
510  (
511  const vector& U, // particle velocity
512  const scalar d, // particle diameter
513  const scalar rhoc, // carrier density
514  const scalar muc // carrier dynamic viscosity
515  ) const;
516 
517  //- Weber number
518  inline scalar We
519  (
520  const vector& U, // particle velocity
521  const scalar d, // particle diameter
522  const scalar rhoc, // carrier density
523  const scalar sigma // particle surface tension
524  ) const;
525 
526  //- Eotvos number
527  inline scalar Eo
528  (
529  const vector& a, // acceleration
530  const scalar d, // particle diameter
531  const scalar sigma // particle surface tension
532  ) const;
533 
534 
535  // Main calculation loop
536 
537  //- Set cell values
538  template<class TrackData>
539  void setCellValues
540  (
541  TrackData& td,
542  const scalar dt,
543  const label cellI
544  );
545 
546  //- Correct cell values using latest transfer information
547  template<class TrackData>
549  (
550  TrackData& td,
551  const scalar dt,
552  const label cellI
553  );
554 
555  //- Update parcel properties over the time interval
556  template<class TrackData>
557  void calc
558  (
559  TrackData& td,
560  const scalar dt,
561  const label cellI
562  );
563 
564 
565  // Tracking
566 
567  //- Move the parcel
568  template<class TrackData>
569  bool move(TrackData& td, const scalar trackTime);
570 
571 
572  // Patch interactions
573 
574  //- Overridable function to handle the particle hitting a face
575  // without trackData
576  void hitFace(int& td);
577 
578  //- Overridable function to handle the particle hitting a face
579  template<class TrackData>
580  void hitFace(TrackData& td);
581 
582  //- Overridable function to handle the particle hitting a patch
583  // Executed before other patch-hitting functions
584  template<class TrackData>
585  bool hitPatch
586  (
587  const polyPatch& p,
588  TrackData& td,
589  const label patchI,
590  const scalar trackFraction,
591  const tetIndices& tetIs
592  );
593 
594  //- Overridable function to handle the particle hitting a
595  // processorPatch
596  template<class TrackData>
597  void hitProcessorPatch
598  (
599  const processorPolyPatch&,
600  TrackData& td
601  );
602 
603  //- Overridable function to handle the particle hitting a wallPatch
604  template<class TrackData>
605  void hitWallPatch
606  (
607  const wallPolyPatch&,
608  TrackData& td,
609  const tetIndices&
610  );
611 
612  //- Overridable function to handle the particle hitting a polyPatch
613  template<class TrackData>
614  void hitPatch
615  (
616  const polyPatch&,
617  TrackData& td
618  );
619 
620  //- Transform the physical properties of the particle
621  // according to the given transformation tensor
622  virtual void transformProperties(const tensor& T);
623 
624  //- Transform the physical properties of the particle
625  // according to the given separation vector
626  virtual void transformProperties(const vector& separation);
627 
628  //- The nearest distance to a wall that the particle can be
629  // in the n direction
630  virtual scalar wallImpactDistance(const vector& n) const;
631 
632 
633  // I-O
634 
635  //- Read
636  template<class CloudType>
637  static void readFields(CloudType& c);
638 
639  //- Write
640  template<class CloudType>
641  static void writeFields(const CloudType& c);
642 
643 
644  // Ostream Operator
645 
646  friend Ostream& operator<< <ParcelType>
647  (
648  Ostream&,
650  );
651 };
652 
653 
654 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
655 
656 } // End namespace Foam
657 
658 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
659 
660 #include "KinematicParcelI.H"
662 
663 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
664 
665 #ifdef NoRepository
666  #include "KinematicParcel.C"
667 #endif
668 
669 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
670 
671 #endif
672 
673 // ************************************************************************* //
Foam::KinematicParcel::sizeofFields_
static const std::size_t sizeofFields_
Size in bytes of the fields.
Definition: KinematicParcel.H:83
Foam::KinematicParcel::iNew
Factory class to read-construct particles used for.
Definition: KinematicParcel.H:377
Foam::Tensor
Templated 3D tensor derived from VectorSpace adding construction from 9 components,...
Definition: complexI.H:224
Foam::KinematicParcel::constantProperties::parcelTypeId
label parcelTypeId() const
Return const access to the parcel type id.
Definition: KinematicParcelI.H:144
Foam::KinematicParcel::typeId
label typeId() const
Return const access to type id.
Definition: KinematicParcelI.H:184
Foam::KinematicParcel::tTurb_
scalar tTurb_
Time spent in turbulent eddy [s].
Definition: KinematicParcel.H:260
Foam::KinematicParcel::UTurb
const vector & UTurb() const
Return const access to turbulent velocity fluctuation.
Definition: KinematicParcelI.H:240
Foam::KinematicParcel::AddToPropertyList
AddToPropertyList(ParcelType, " active"+" typeId"+" nParticle"+" d"+" dTarget "+" (Ux Uy Uz)"+" rho"+" age"+" tTurb"+" (UTurbx UTurby UTurbz)")
String representation of properties.
p
p
Definition: pEqn.H:62
Foam::KinematicParcel::cellValueSourceCorrection
void cellValueSourceCorrection(TrackData &td, const scalar dt, const label cellI)
Correct cell values using latest transfer information.
Definition: KinematicParcel.C:85
Foam::constant::physicoChemical::mu
const dimensionedScalar mu
Atomic mass unit.
Definition: createFields.H:13
Foam::KinematicParcel::nParticle
scalar nParticle() const
Return const access to number of particles.
Definition: KinematicParcelI.H:191
Foam::KinematicParcel::muc_
scalar muc_
Viscosity [Pa.s].
Definition: KinematicParcel.H:275
Foam::KinematicParcel::dTarget_
scalar dTarget_
Target diameter [m].
Definition: KinematicParcel.H:248
Foam::KinematicParcel::iNew::mesh_
const polyMesh & mesh_
Definition: KinematicParcel.H:379
Foam::KinematicParcel::constantProperties::minParcelMass_
demandDrivenEntry< scalar > minParcelMass_
Minimum parcel mass [kg].
Definition: KinematicParcel.H:117
Foam::KinematicParcel::U
const vector & U() const
Return const access to velocity.
Definition: KinematicParcelI.H:212
Foam::KinematicParcel::tTurb
scalar tTurb() const
Return const access to time spent in turbulent eddy.
Definition: KinematicParcelI.H:233
Foam::KinematicParcel::wallImpactDistance
virtual scalar wallImpactDistance(const vector &n) const
The nearest distance to a wall that the particle can be.
Definition: KinematicParcel.C:486
demandDrivenEntry.H
Foam::KinematicParcel::TrackingData::trackPart
trackPart
Definition: KinematicParcel.H:160
Foam::KinematicParcel::dTarget
scalar dTarget() const
Return const access to target diameter.
Definition: KinematicParcelI.H:205
Foam::KinematicParcel::TrackingData::TrackingData
TrackingData(CloudType &cloud, trackPart part=tpLinearTrack)
Construct from components.
Definition: KinematicParcelTrackingDataI.H:29
Foam::KinematicParcel::TrackingData::tpLinearTrack
@ tpLinearTrack
Definition: KinematicParcel.H:163
Foam::KinematicParcel::constantProperties::rho0
scalar rho0() const
Return const access to the particle density.
Definition: KinematicParcelI.H:160
Foam::KinematicParcel::TrackingData::muInterp
const interpolation< scalar > & muInterp() const
Return conat access to the interpolator for continuous.
Definition: KinematicParcelTrackingDataI.H:85
interpolation.H
Foam::KinematicParcel::constantProperties::parcelTypeId_
demandDrivenEntry< label > parcelTypeId_
Parcel type id - used for post-processing to flag the type.
Definition: KinematicParcel.H:108
Foam::cp
bool cp(const fileName &src, const fileName &dst)
Copy, recursively if necessary, the source to the destination.
Definition: POSIX.C:755
Foam::KinematicParcel::calcVelocity
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.
Foam::readFields
This function object reads fields from the time directories and adds them to the mesh database for fu...
Definition: readFields.H:104
Foam::KinematicParcel::clone
virtual autoPtr< particle > clone(const polyMesh &mesh) const
Construct and return a (basic particle) clone.
Definition: KinematicParcel.H:370
Foam::KinematicParcel::move
bool move(TrackData &td, const scalar trackTime)
Move the parcel.
Definition: KinematicParcel.C:255
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
KinematicParcelI.H
Foam::KinematicParcel::age_
scalar age_
Age [s].
Definition: KinematicParcel.H:257
Foam::KinematicParcel::TrackingData::g
const vector & g() const
Definition: KinematicParcelTrackingDataI.H:94
Foam::KinematicParcel::age
scalar age() const
Return const access to the age.
Definition: KinematicParcelI.H:226
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::KinematicParcel::iNew::operator()
autoPtr< KinematicParcel< ParcelType > > operator()(Istream &is) const
Definition: KinematicParcel.H:388
Foam::KinematicParcel::TrackingData::tpVelocityHalfStep
@ tpVelocityHalfStep
Definition: KinematicParcel.H:162
Foam::KinematicParcel::muc
scalar muc() const
Return const access to carrier viscosity [Pa.s].
Definition: KinematicParcelI.H:261
Foam::KinematicParcel::TrackingData::part_
trackPart part_
Definition: KinematicParcel.H:189
Foam::label
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
Foam::fvc::Su
tmp< GeometricField< Type, fvPatchField, volMesh > > Su(const GeometricField< Type, fvPatchField, volMesh > &su, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcSup.C:44
Foam::KinematicParcel::rhoc
scalar rhoc() const
Return const access to carrier density [kg/m3].
Definition: KinematicParcelI.H:247
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:57
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Foam::KinematicParcel::TrackingData::UInterp
const interpolation< vector > & UInterp() const
Return conat access to the interpolator for continuous.
Definition: KinematicParcelTrackingDataI.H:76
Foam::KinematicParcel::UTurb_
vector UTurb_
Turbulent velocity fluctuation [m/s].
Definition: KinematicParcel.H:263
IOstream.H
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:55
Foam::KinematicParcel::TrackingData::part
trackPart part() const
Return the part of the tracking operation taking place.
Definition: KinematicParcelTrackingDataI.H:104
Foam::KinematicParcel::Uc_
vector Uc_
Velocity [m/s].
Definition: KinematicParcel.H:272
Foam::KinematicParcel::momentOfInertia
scalar momentOfInertia() const
Particle moment of inertia around diameter axis.
Definition: KinematicParcelI.H:370
Foam::KinematicParcel::hitProcessorPatch
void hitProcessorPatch(const processorPolyPatch &, TrackData &td)
Overridable function to handle the particle hitting a.
Definition: KinematicParcel.C:429
Foam::processorPolyPatch
Neighbour processor patch.
Definition: processorPolyPatch.H:55
Foam::KinematicParcel::volume
scalar volume() const
Particle volume.
Definition: KinematicParcelI.H:377
Foam::KinematicParcel::TrackingData::g_
const vector & g_
Local gravitational or other body-force acceleration.
Definition: KinematicParcel.H:185
Foam::KinematicParcel::constantProperties::minParcelMass
scalar minParcelMass() const
Return const access to the minimum parcel mass.
Definition: KinematicParcelI.H:168
Foam::KinematicParcel::areaP
scalar areaP() const
Particle projected area.
Definition: KinematicParcelI.H:391
Foam::DSMCCloud
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:68
Foam::KinematicParcel::faceInterpolation
label faceInterpolation() const
Return the index of the face used in the interpolation routine.
Definition: KinematicParcelI.H:338
Foam::interpolation< scalar >
Foam::KinematicParcel::U_
vector U_
Velocity of Parcel [m/s].
Definition: KinematicParcel.H:251
Foam::KinematicParcel::d_
scalar d_
Diameter [m].
Definition: KinematicParcel.H:245
Foam::KinematicParcel::Re
scalar Re(const vector &U, const scalar d, const scalar rhoc, const scalar muc) const
Reynolds number.
Definition: KinematicParcelI.H:420
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::KinematicParcel::massCell
scalar massCell(const label cellI) const
Cell owner mass.
Definition: KinematicParcelI.H:354
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::KinematicParcel::nParticle_
scalar nParticle_
Number of particles in Parcel.
Definition: KinematicParcel.H:242
Foam::wallPolyPatch
Foam::wallPolyPatch.
Definition: wallPolyPatch.H:48
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::KinematicParcel::rho_
scalar rho_
Density [kg/m3].
Definition: KinematicParcel.H:254
Foam::KinematicParcel::active_
bool active_
Active flag - tracking inactive when active = false.
Definition: KinematicParcel.H:236
Foam::KinematicParcel::hitFace
void hitFace(int &td)
Overridable function to handle the particle hitting a face.
Foam::KinematicParcel::active
bool active() const
Return const access to active flag.
Definition: KinematicParcelI.H:177
Foam::KinematicParcel::TrackingData::tpRotationalTrack
@ tpRotationalTrack
Definition: KinematicParcel.H:164
Foam::KinematicParcel::TrackingData::rhoInterp_
autoPtr< interpolation< scalar > > rhoInterp_
Density interpolator.
Definition: KinematicParcel.H:175
Foam::KinematicParcel::typeId_
label typeId_
Parcel type id.
Definition: KinematicParcel.H:239
Foam::KinematicParcel::rhoc_
scalar rhoc_
Density [kg/m3].
Definition: KinematicParcel.H:269
Foam::KinematicParcel::areaS
scalar areaS() const
Particle surface area.
Definition: KinematicParcelI.H:405
Foam::KinematicParcel::TrackingData::UInterp_
autoPtr< interpolation< vector > > UInterp_
Velocity interpolator.
Definition: KinematicParcel.H:178
Foam::KinematicParcel::hitPatch
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.
Definition: KinematicParcel.C:384
Foam::KinematicParcel::constantProperties::rhoMin
scalar rhoMin() const
Return const access to the minimum density.
Definition: KinematicParcelI.H:152
Foam::KinematicParcel::We
scalar We(const vector &U, const scalar d, const scalar rhoc, const scalar sigma) const
Weber number.
Definition: KinematicParcelI.H:433
Foam::tetIndices
Storage and named access for the indices of a tet which is part of the decomposition of a cell.
Definition: tetIndices.H:73
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
Foam::cloud
A cloud is a collection of lagrangian particles.
Definition: cloud.H:51
Foam::KinematicParcel::constantProperties::rho0_
demandDrivenEntry< scalar > rho0_
Particle density [kg/m3] (constant)
Definition: KinematicParcel.H:114
Foam::KinematicParcel::Eo
scalar Eo(const vector &a, const scalar d, const scalar sigma) const
Eotvos number.
Definition: KinematicParcelI.H:446
Foam::KinematicParcel::setCellValues
void setCellValues(TrackData &td, const scalar dt, const label cellI)
Set cell values.
Definition: KinematicParcel.C:43
Foam::KinematicParcel::maxTrackAttempts
static label maxTrackAttempts
Number of particle tracking attempts before we assume that it stalls.
Definition: KinematicParcel.H:86
KinematicParcel.C
Foam::Vector< scalar >
Foam::demandDrivenEntry< label >
Foam::KinematicParcel
Kinematic parcel class with rotational motion (as spherical particles only) and one/two-way coupling ...
Definition: KinematicParcel.H:60
Foam::KinematicParcel::TypeName
TypeName("KinematicParcel")
Runtime type information.
Foam::constant::physicoChemical::sigma
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
particle.H
Foam::KinematicParcel::constantProperties::dict
const dictionary & dict() const
Return const access to the constant properties dictionary.
Definition: KinematicParcelI.H:136
Foam::KinematicParcel::transformProperties
virtual void transformProperties(const tensor &T)
Transform the physical properties of the particle.
Definition: KinematicParcel.C:466
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::KinematicParcel::TrackingData::rhoInterp
const interpolation< scalar > & rhoInterp() const
Return conat access to the interpolator for continuous.
Definition: KinematicParcelTrackingDataI.H:67
Foam::KinematicParcel::Uc
const vector & Uc() const
Return const access to carrier velocity [m/s].
Definition: KinematicParcelI.H:254
Foam::KinematicParcel::TrackingData
Definition: KinematicParcel.H:154
Foam::KinematicParcel::d
scalar d() const
Return const access to diameter.
Definition: KinematicParcelI.H:198
Foam::KinematicParcel::rho
scalar rho() const
Return const access to density.
Definition: KinematicParcelI.H:219
Foam::KinematicParcel::KinematicParcel
KinematicParcel(const polyMesh &mesh, const vector &position, const label cellI, const label tetFaceI, const label tetPtI)
Construct from owner, position, and cloud owner.
Foam::KinematicParcel::calc
void calc(TrackData &td, const scalar dt, const label cellI)
Update parcel properties over the time interval.
Definition: KinematicParcel.C:98
Foam::KinematicParcel::writeFields
static void writeFields(const CloudType &c)
Write.
Definition: KinematicParcelIO.C:166
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::KinematicParcel::TrackingData::muInterp_
autoPtr< interpolation< scalar > > muInterp_
Dynamic viscosity interpolator.
Definition: KinematicParcel.H:181
Foam::KinematicParcel::constantProperties::dict_
const dictionary dict_
Constant properties dictionary.
Definition: KinematicParcel.H:99
Foam::KinematicParcel::clone
virtual autoPtr< particle > clone() const
Construct and return a (basic particle) clone.
Definition: KinematicParcel.H:364
Foam::KinematicParcel::constantProperties::rhoMin_
demandDrivenEntry< scalar > rhoMin_
Minimum density [kg/m3].
Definition: KinematicParcel.H:111
Foam::KinematicParcel::constantProperties::constantProperties
constantProperties()
Null constructor.
Definition: KinematicParcelI.H:34
Foam::KinematicParcel::mass
scalar mass() const
Particle mass.
Definition: KinematicParcelI.H:363
Foam::KinematicParcel::readFields
static void readFields(CloudType &c)
Read.
Definition: KinematicParcelIO.C:102
Foam::KinematicParcel::constantProperties
Class to hold kinematic particle constant properties.
Definition: KinematicParcel.H:92
KinematicParcelTrackingDataI.H
Foam::KinematicParcel::hitWallPatch
void hitWallPatch(const wallPolyPatch &, TrackData &td, const tetIndices &)
Overridable function to handle the particle hitting a wallPatch.
Definition: KinematicParcel.C:441
Foam::KinematicParcel::iNew::iNew
iNew(const polyMesh &mesh)
Definition: KinematicParcel.H:383
autoPtr.H