MomentumParcel.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-2022 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::MomentumParcel
26 
27 Description
28  Momentum 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  MomentumParcelI.H
39  MomentumParcel.C
40  MomentumParcelIO.C
41 
42 \*---------------------------------------------------------------------------*/
43 
44 #ifndef MomentumParcel_H
45 #define MomentumParcel_H
46 
47 #include "particle.H"
48 #include "interpolation.H"
49 #include "demandDrivenEntry.H"
50 
51 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52 
53 namespace Foam
54 {
55 
56 template<class ParcelType>
57 class MomentumParcel;
58 
59 // Forward declaration of friend functions
60 
61 template<class ParcelType>
62 Ostream& operator<<
63 (
64  Ostream&,
66 );
67 
68 /*---------------------------------------------------------------------------*\
69  Class MomentumParcelName Declaration
70 \*---------------------------------------------------------------------------*/
71 
73 
74 
75 /*---------------------------------------------------------------------------*\
76  Class MomentumParcel Declaration
77 \*---------------------------------------------------------------------------*/
78 
79 template<class ParcelType>
80 class MomentumParcel
81 :
82  public ParcelType,
83  public MomentumParcelName
84 {
85  // Private Data
86 
87  //- Size in bytes of the fields
88  static const std::size_t sizeofFields_;
89 
90  //- Number of particle tracking attempts before we assume that it stalls
91  static label maxTrackAttempts;
92 
93 
94 public:
95 
96  //- Class to hold momentum parcel constant properties
97  class constantProperties
98  {
99  protected:
100 
101  // Protected data
102 
103  //- Constant properties dictionary
104  const dictionary dict_;
105 
106 
107  private:
108 
109  // Private Data
110 
111  //- Parcel type id - used for post-processing to flag the type
112  // of parcels issued by this cloud
113  demandDrivenEntry<label> parcelTypeId_;
114 
115  //- Minimum density [kg/m^3]
117 
118  //- Particle density [kg/m^3] (constant)
120 
121  //- Minimum parcel mass [kg]
122  demandDrivenEntry<scalar> minParcelMass_;
123 
124 
125  public:
126 
127  // Constructors
128 
129  //- Null constructor
131 
132  //- Copy constructor
134 
135  //- Construct from dictionary
136  constantProperties(const dictionary& parentDict);
137 
138 
139  // Member Functions
140 
141  //- Return const access to the constant properties dictionary
142  inline const dictionary& dict() const;
143 
144  //- Return const access to the parcel type id
145  inline label parcelTypeId() const;
146 
147  //- Return const access to the minimum density
148  inline scalar rhoMin() const;
149 
150  //- Return const access to the particle density
151  inline scalar rho0() const;
152 
153  //- Return const access to the minimum parcel mass
154  inline scalar minParcelMass() const;
155  };
156 
157 
158  class trackingData
159  :
160  public ParcelType::trackingData
161  {
162  private:
163 
164  // Private Data
165 
166  // Interpolators for continuous phase fields
167 
168  //- Density interpolator
169  autoPtr<interpolation<scalar>> rhoInterp_;
170 
171  //- Velocity interpolator
173 
174  //- Dynamic viscosity interpolator
176 
177 
178  // Cached continuous phase properties
179 
180  //- Density [kg/m^3]
181  scalar rhoc_;
182 
183  //- Velocity [m/s]
184  vector Uc_;
185 
186  //- Viscosity [Pa.s]
187  scalar muc_;
188 
189 
190  //- Local gravitational or other body-force acceleration
191  const vector& g_;
192 
193 
194  // Tracking Parameters
195 
196  //- Track time. Total across the entire step.
197  scalar trackTime_;
198 
199  //- Step fraction range to track between. Allows for the
200  // creation of sub-steps.
201  Pair<scalar> stepFractionRange_;
202 
203 
204  public:
205 
206  // Constructors
207 
208  //- Construct from components
209  template <class TrackCloudType>
210  inline trackingData(const TrackCloudType& cloud);
211 
212 
213  // Member Functions
214 
215  //- Return const access to the interpolator for continuous
216  // phase density field
217  inline const interpolation<scalar>& rhoInterp() const;
218 
219  //- Return const access to the interpolator for continuous
220  // phase velocity field
221  inline const interpolation<vector>& UInterp() const;
222 
223  //- Return const access to the interpolator for continuous
224  // phase dynamic viscosity field
225  inline const interpolation<scalar>& muInterp() const;
226 
227  //- Return the continuous phase density
228  inline scalar rhoc() const;
229 
230  //- Access the continuous phase density
231  inline scalar& rhoc();
232 
233  //- Return the continuous phase velocity
234  inline const vector& Uc() const;
235 
236  //- Access the continuous phase velocity
237  inline vector& Uc();
238 
239  //- Return the continuous phase viscosity
240  inline scalar muc() const;
241 
242  //- Access the continuous phase viscosity
243  inline scalar& muc();
244 
245  // Return the gravitational acceleration vector
246  inline const vector& g() const;
247 
248  //- Return the tracking time
249  inline scalar trackTime() const;
250 
251  //- Access the tracking time
252  inline scalar& trackTime();
253 
254  //- Return the step fraction range to track between
255  inline Pair<scalar> stepFractionRange() const;
256 
257  //- Access the step fraction range to track between
259  };
260 
261 
262 protected:
263 
264  // Protected data
265 
266  // Parcel properties
267 
268  //- Moving flag - tracking stopped when moving = false
269  bool moving_;
270 
271  //- Parcel type id
272  label typeId_;
273 
274  //- Number of particles in Parcel
275  scalar nParticle_;
276 
277  //- Diameter [m]
278  scalar d_;
279 
280  //- Target diameter [m]
281  scalar dTarget_;
282 
283  //- Velocity of Parcel [m/s]
284  vector U_;
285 
286  //- Density [kg/m^3]
287  scalar rho_;
288 
289  //- Age [s]
290  scalar age_;
291 
292  //- Time spent in turbulent eddy [s]
293  scalar tTurb_;
294 
295  //- Turbulent velocity fluctuation [m/s]
296  vector UTurb_;
297 
298 
299  // Protected Member Functions
300 
301  //- Set cell values
302  template<class TrackCloudType>
303  void setCellValues(TrackCloudType& cloud, trackingData& td);
304 
305  //- Apply dispersion to the carrier phase velocity and update
306  // parcel turbulence parameters
307  template<class TrackCloudType>
308  void calcDispersion
309  (
310  TrackCloudType& cloud,
311  trackingData& td,
312  const scalar dt
313  );
314 
315  //- Correct cell values using latest transfer information
316  template<class TrackCloudType>
318  (
319  TrackCloudType& cloud,
320  trackingData& td,
321  const scalar dt
322  );
323 
324  //- Update parcel properties over the time interval
325  template<class TrackCloudType>
326  void calc
327  (
328  TrackCloudType& cloud,
329  trackingData& td,
330  const scalar dt
331  );
332 
333  //- Calculate new particle velocity
334  template<class TrackCloudType>
335  const vector calcVelocity
336  (
337  TrackCloudType& cloud,
338  trackingData& td,
339  const scalar dt, // timestep
340  const scalar Re, // Reynolds number
341  const scalar mu, // local carrier viscosity
342  const scalar mass, // mass
343  const vector& Su, // explicit particle momentum source
344  vector& dUTrans, // momentum transfer to carrier
345  scalar& Spu // linearised drag coefficient
346  ) const;
347 
348 
349 public:
350 
351  // Static Data Members
352 
353  //- String representation of properties
355  (
356  ParcelType,
357  " moving"
358  + " typeId"
359  + " nParticle"
360  + " d"
361  + " dTarget "
362  + " (Ux Uy Uz)"
363  + " rho"
364  + " age"
365  + " tTurb"
366  + " (UTurbx UTurby UTurbz)"
367  );
368 
369 
370  // Constructors
371 
372  //- Construct from mesh, coordinates and topology
373  // Other properties initialised as null
374  inline MomentumParcel
375  (
376  const polyMesh& mesh,
377  const barycentric& coordinates,
378  const label celli,
379  const label tetFacei,
380  const label tetPti,
381  const label facei
382  );
383 
384  //- Construct from a position and a cell, searching for the rest of the
385  // required topology. Other properties are initialised as null.
386  inline MomentumParcel
387  (
388  const polyMesh& mesh,
389  const vector& position,
390  const label celli
391  );
392 
393  //- Construct from Istream
394  MomentumParcel(Istream& is, bool readFields = true);
395 
396  //- Construct as a copy
398 
399  //- Construct and return a clone
400  virtual autoPtr<particle> clone() const
401  {
402  return autoPtr<particle>(new MomentumParcel(*this));
403  }
404 
405  //- Construct from Istream and return
407  {
408  return autoPtr<MomentumParcel>(new MomentumParcel(is));
409  }
410 
411 
412  // Member Functions
413 
414  // Access
415 
416  //- Return const access to moving flag
417  inline bool moving() const;
418 
419  //- Return const access to type id
420  inline label typeId() const;
421 
422  //- Return const access to number of particles
423  inline scalar nParticle() const;
424 
425  //- Return const access to diameter
426  inline scalar d() const;
427 
428  //- Return const access to target diameter
429  inline scalar dTarget() const;
430 
431  //- Return const access to velocity
432  inline const vector& U() const;
433 
434  //- Return const access to density
435  inline scalar rho() const;
436 
437  //- Return const access to the age
438  inline scalar age() const;
439 
440  //- Return const access to time spent in turbulent eddy
441  inline scalar tTurb() const;
442 
443  //- Return const access to turbulent velocity fluctuation
444  inline const vector& UTurb() const;
445 
446 
447  // Edit
448 
449  //- Return const access to moving flag
450  inline bool& moving();
451 
452  //- Return access to type id
453  inline label& typeId();
454 
455  //- Return access to number of particles
456  inline scalar& nParticle();
457 
458  //- Return access to diameter
459  inline scalar& d();
460 
461  //- Return access to target diameter
462  inline scalar& dTarget();
463 
464  //- Return access to velocity
465  inline vector& U();
466 
467  //- Return access to density
468  inline scalar& rho();
469 
470  //- Return access to the age
471  inline scalar& age();
472 
473  //- Return access to time spent in turbulent eddy
474  inline scalar& tTurb();
475 
476  //- Return access to turbulent velocity fluctuation
477  inline vector& UTurb();
478 
479 
480  // Helper functions
481 
482  //- Cell owner mass
483  inline scalar massCell(const trackingData& td) const;
484 
485  //- Particle mass
486  inline scalar mass() const;
487 
488  //- Particle moment of inertia around diameter axis
489  inline scalar momentOfInertia() const;
490 
491  //- Particle volume
492  inline scalar volume() const;
493 
494  //- Particle volume for a given diameter
495  inline static scalar volume(const scalar d);
496 
497  //- Particle projected area
498  inline scalar areaP() const;
499 
500  //- Projected area for given diameter
501  inline static scalar areaP(const scalar d);
502 
503  //- Particle surface area
504  inline scalar areaS() const;
505 
506  //- Surface area for given diameter
507  inline static scalar areaS(const scalar d);
508 
509  //- Reynolds number
510  inline scalar Re(const trackingData& td) const;
511 
512  //- Reynolds number for given conditions
513  inline static scalar Re
514  (
515  const scalar rhoc,
516  const vector& U,
517  const vector& Uc,
518  const scalar d,
519  const scalar muc
520  );
521 
522  //- Weber number
523  inline scalar We
524  (
525  const trackingData& td,
526  const scalar sigma
527  ) const;
528 
529  //- Weber number for given conditions
530  inline static scalar We
531  (
532  const scalar rhoc,
533  const vector& U,
534  const vector& Uc,
535  const scalar d,
536  const scalar sigma
537  );
538 
539  //- Eotvos number
540  inline scalar Eo
541  (
542  const trackingData& td,
543  const scalar sigma
544  ) const;
545 
546  //- Eotvos number for given conditions
547  inline static scalar Eo
548  (
549  const vector& g,
550  const scalar rho,
551  const scalar rhoc,
552  const vector& U,
553  const scalar d,
554  const scalar sigma
555  );
556 
557 
558  // Tracking
559 
560  //- Move the parcel
561  template<class TrackCloudType>
562  bool move(TrackCloudType& cloud, trackingData& td);
563 
564 
565  // Transformations
566 
567  //- Transform the physical properties of the particle
568  // according to the given transformation
569  virtual void transformProperties(const transformer&);
570 
571 
572  // Transfers
573 
574  //- Make changes following a parallel transfer
575  template<class TrackCloudType>
576  void correctAfterParallelTransfer(TrackCloudType&, trackingData&);
577 
578 
579  // Patch interactions
580 
581  //- Overridable function to handle the particle hitting a patch
582  // Executed before other patch-hitting functions
583  template<class TrackCloudType>
584  bool hitPatch(TrackCloudType& cloud, trackingData& td);
585 
586  //- Overridable function to handle the particle hitting a wedgePatch
587  template<class TrackCloudType>
588  void hitWedgePatch(TrackCloudType&, trackingData&);
589 
590  //- Overridable function to handle the particle hitting a
591  // symmetryPlanePatch
592  template<class TrackCloudType>
593  void hitSymmetryPlanePatch(TrackCloudType&, trackingData&);
594 
595  //- Overridable function to handle the particle hitting a
596  // symmetryPatch
597  template<class TrackCloudType>
598  void hitSymmetryPatch(TrackCloudType&, trackingData&);
599 
600  //- Overridable function to handle the particle hitting a
601  // cyclicPatch
602  template<class TrackCloudType>
603  void hitCyclicPatch(TrackCloudType&, trackingData&);
604 
605  //- Overridable function to handle the particle hitting an
606  // nonConformalCyclicPolyPatch
607  template<class TrackCloudType>
609  (
610  const vector& displacement,
611  const scalar fraction,
612  const label patchi,
613  TrackCloudType& cloud,
614  trackingData& td
615  );
616 
617  //- Overridable function to handle the particle hitting a
618  // processorPatch
619  using ParcelType::hitProcessorPatch;
620 
621  //- Overridable function to handle the particle hitting a wallPatch
622  template<class TrackCloudType>
623  void hitWallPatch(TrackCloudType& cloud, trackingData& td);
624 
625  //- Overridable function to handle the particle hitting a basic
626  // patch. Fall-through for the above.
627  template<class TrackCloudType>
628  void hitBasicPatch(TrackCloudType& cloud, trackingData& td);
629 
630 
631  // I-O
632 
633  //- Read
634  template<class TrackCloudType>
635  static void readFields(TrackCloudType& c);
636 
637  //- Write
638  template<class TrackCloudType>
639  static void writeFields(const TrackCloudType& c);
640 
641 
642  // Ostream Operator
643 
644  friend Ostream& operator<< <ParcelType>
645  (
646  Ostream&,
648  );
649 };
650 
651 
652 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
653 
654 } // End namespace Foam
655 
656 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
657 
658 #include "MomentumParcelI.H"
660 
661 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
662 
663 #ifdef NoRepository
664  #include "MomentumParcel.C"
665 #endif
666 
667 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
668 
669 #endif
670 
671 // ************************************************************************* //
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
Class to hold momentum parcel constant properties.
scalar minParcelMass() const
Return const access to the minimum parcel mass.
scalar rho0() const
Return const access to the particle density.
label parcelTypeId() const
Return const access to the parcel type id.
const dictionary dict_
Constant properties dictionary.
const dictionary & dict() const
Return const access to the constant properties dictionary.
scalar rhoMin() const
Return const access to the minimum density.
const vector & Uc() const
Return the continuous phase velocity.
scalar muc() const
Return the continuous phase viscosity.
const interpolation< vector > & UInterp() const
Return const access to the interpolator for continuous.
scalar trackTime() const
Return the tracking time.
trackingData(const TrackCloudType &cloud)
Construct from components.
scalar rhoc() const
Return the continuous phase density.
Pair< scalar > stepFractionRange() const
Return the step fraction range to track between.
const interpolation< scalar > & muInterp() const
Return const access to the interpolator for continuous.
const interpolation< scalar > & rhoInterp() const
Return const access to the interpolator for continuous.
Momentum parcel class with rotational motion (as spherical particles only) and one/two-way coupling w...
label typeId_
Parcel type id.
label typeId() const
Return const access to type id.
vector UTurb_
Turbulent velocity fluctuation [m/s].
scalar momentOfInertia() const
Particle moment of inertia around diameter axis.
void calcDispersion(TrackCloudType &cloud, trackingData &td, const scalar dt)
Apply dispersion to the carrier phase velocity and update.
scalar tTurb() const
Return const access to time spent in turbulent eddy.
scalar rho_
Density [kg/m^3].
const vector & U() const
Return const access to velocity.
bool moving_
Moving flag - tracking stopped when moving = false.
bool move(TrackCloudType &cloud, trackingData &td)
Move the parcel.
static void readFields(TrackCloudType &c)
Read.
scalar d() const
Return const access to diameter.
void cellValueSourceCorrection(TrackCloudType &cloud, trackingData &td, const scalar dt)
Correct cell values using latest transfer information.
void hitCyclicPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a.
static autoPtr< MomentumParcel > New(Istream &is)
Construct from Istream and return.
scalar Re(const trackingData &td) const
Reynolds number.
scalar tTurb_
Time spent in turbulent eddy [s].
scalar volume() const
Particle volume.
scalar dTarget_
Target diameter [m].
const vector calcVelocity(TrackCloudType &cloud, trackingData &td, const scalar dt, const scalar Re, const scalar mu, const scalar mass, const vector &Su, vector &dUTrans, scalar &Spu) const
Calculate new particle velocity.
scalar nParticle() const
Return const access to number of particles.
scalar massCell(const trackingData &td) const
Cell owner mass.
scalar d_
Diameter [m].
scalar dTarget() const
Return const access to target diameter.
void hitSymmetryPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a.
virtual autoPtr< particle > clone() const
Construct and return a clone.
void hitSymmetryPlanePatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a.
scalar We(const trackingData &td, const scalar sigma) const
Weber number.
scalar areaP() const
Particle projected area.
scalar rho() const
Return const access to density.
scalar mass() const
Particle mass.
MomentumParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti, const label facei)
Construct from mesh, coordinates and topology.
static void writeFields(const TrackCloudType &c)
Write.
AddToPropertyList(ParcelType, " moving"+" typeId"+" nParticle"+" d"+" dTarget "+" (Ux Uy Uz)"+" rho"+" age"+" tTurb"+" (UTurbx UTurby UTurbz)")
String representation of properties.
bool hitPatch(TrackCloudType &cloud, trackingData &td)
Overridable function to handle the particle hitting a patch.
void hitWedgePatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a wedgePatch.
void setCellValues(TrackCloudType &cloud, trackingData &td)
Set cell values.
scalar areaS() const
Particle surface area.
void correctAfterParallelTransfer(TrackCloudType &, trackingData &)
Make changes following a parallel transfer.
void hitWallPatch(TrackCloudType &cloud, trackingData &td)
Overridable function to handle the particle hitting a wallPatch.
scalar Eo(const trackingData &td, const scalar sigma) const
Eotvos number.
bool hitNonConformalCyclicPatch(const vector &displacement, const scalar fraction, const label patchi, TrackCloudType &cloud, trackingData &td)
Overridable function to handle the particle hitting an.
scalar age_
Age [s].
virtual void transformProperties(const transformer &)
Transform the physical properties of the particle.
scalar age() const
Return const access to the age.
vector U_
Velocity of Parcel [m/s].
const vector & UTurb() const
Return const access to turbulent velocity fluctuation.
scalar nParticle_
Number of particles in Parcel.
void hitBasicPatch(TrackCloudType &cloud, trackingData &td)
Overridable function to handle the particle hitting a basic.
bool moving() const
Return const access to moving flag.
void calc(TrackCloudType &cloud, trackingData &td, const scalar dt)
Update parcel properties over the time interval.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
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
Abstract base class for interpolation.
Definition: interpolation.H:55
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
Vector-tensor class used to perform translations, rotations and scaling operations in 3D space.
Definition: transformer.H:84
label patchi
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
const dimensionedScalar mu
Atomic mass unit.
const dimensionedScalar c
Speed of light in a vacuum.
tmp< VolField< Type > > Su(const VolField< Type > &su, const VolField< Type > &vf)
Definition: fvcSup.C:44
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)
volScalarField & p