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-2020 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 
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  public:
195 
196  // Constructors
197 
198  //- Construct from components
199  template <class TrackCloudType>
200  inline trackingData(const TrackCloudType& cloud);
201 
202 
203  // Member Functions
204 
205  //- Return const access to the interpolator for continuous
206  // phase density field
207  inline const interpolation<scalar>& rhoInterp() const;
208 
209  //- Return const access to the interpolator for continuous
210  // phase velocity field
211  inline const interpolation<vector>& UInterp() const;
212 
213  //- Return const access to the interpolator for continuous
214  // phase dynamic viscosity field
215  inline const interpolation<scalar>& muInterp() const;
216 
217  //- Return the continuous phase density
218  inline scalar rhoc() const;
219 
220  //- Access the continuous phase density
221  inline scalar& rhoc();
222 
223  //- Return the continuous phase velocity
224  inline const vector& Uc() const;
225 
226  //- Access the continuous phase velocity
227  inline vector& Uc();
228 
229  //- Return the continuous phase viscosity
230  inline scalar muc() const;
231 
232  //- Access the continuous phase viscosity
233  inline scalar& muc();
234 
235  // Return const access to the gravitational acceleration vector
236  inline const vector& g() const;
237  };
238 
239 
240 protected:
241 
242  // Protected data
243 
244  // Parcel properties
245 
246  //- Moving flag - tracking stopped when moving = false
247  bool moving_;
248 
249  //- Parcel type id
250  label typeId_;
251 
252  //- Number of particles in Parcel
253  scalar nParticle_;
254 
255  //- Diameter [m]
256  scalar d_;
257 
258  //- Target diameter [m]
259  scalar dTarget_;
260 
261  //- Velocity of Parcel [m/s]
262  vector U_;
263 
264  //- Density [kg/m^3]
265  scalar rho_;
266 
267  //- Age [s]
268  scalar age_;
269 
270  //- Time spent in turbulent eddy [s]
271  scalar tTurb_;
272 
273  //- Turbulent velocity fluctuation [m/s]
274  vector UTurb_;
275 
276 
277  // Protected Member Functions
278 
279  //- Calculate new particle velocity
280  template<class TrackCloudType>
281  const vector calcVelocity
282  (
283  TrackCloudType& cloud,
284  trackingData& td,
285  const scalar dt, // timestep
286  const scalar Re, // Reynolds number
287  const scalar mu, // local carrier viscosity
288  const scalar mass, // mass
289  const vector& Su, // explicit particle momentum source
290  vector& dUTrans, // momentum transfer to carrier
291  scalar& Spu // linearised drag coefficient
292  ) const;
293 
294 
295 public:
296 
297  // Static Data Members
298 
299  //- String representation of properties
301  (
302  ParcelType,
303  " moving"
304  + " typeId"
305  + " nParticle"
306  + " d"
307  + " dTarget "
308  + " (Ux Uy Uz)"
309  + " rho"
310  + " age"
311  + " tTurb"
312  + " (UTurbx UTurby UTurbz)"
313  );
314 
315 
316  // Constructors
317 
318  //- Construct from mesh, coordinates and topology
319  // Other properties initialised as null
320  inline MomentumParcel
321  (
322  const polyMesh& mesh,
323  const barycentric& coordinates,
324  const label celli,
325  const label tetFacei,
326  const label tetPti
327  );
328 
329  //- Construct from a position and a cell, searching for the rest of the
330  // required topology. Other properties are initialised as null.
331  inline MomentumParcel
332  (
333  const polyMesh& mesh,
334  const vector& position,
335  const label celli
336  );
337 
338  //- Construct from Istream
340  (
341  const polyMesh& mesh,
342  Istream& is,
343  bool readFields = true
344  );
345 
346  //- Construct as a copy
348 
349  //- Construct as a copy
350  MomentumParcel(const MomentumParcel& p, const polyMesh& mesh);
351 
352  //- Construct and return a (basic particle) clone
353  virtual autoPtr<particle> clone() const
354  {
355  return autoPtr<particle>(new MomentumParcel(*this));
356  }
357 
358  //- Construct and return a (basic particle) clone
359  virtual autoPtr<particle> clone(const polyMesh& mesh) const
360  {
361  return autoPtr<particle>(new MomentumParcel(*this, mesh));
362  }
363 
364  //- Factory class to read-construct particles used for
365  // parallel transfer
366  class iNew
367  {
368  const polyMesh& mesh_;
369 
370  public:
372  iNew(const polyMesh& mesh)
373  :
374  mesh_(mesh)
375  {}
377  autoPtr<MomentumParcel<ParcelType>> operator()(Istream& is) const
378  {
380  (
381  new MomentumParcel<ParcelType>(mesh_, is, true)
382  );
383  }
384  };
385 
386 
387  // Member Functions
388 
389  // Access
390 
391  //- Return const access to moving flag
392  inline bool moving() const;
393 
394  //- Return const access to type id
395  inline label typeId() const;
396 
397  //- Return const access to number of particles
398  inline scalar nParticle() const;
399 
400  //- Return const access to diameter
401  inline scalar d() const;
402 
403  //- Return const access to target diameter
404  inline scalar dTarget() const;
405 
406  //- Return const access to velocity
407  inline const vector& U() const;
408 
409  //- Return const access to density
410  inline scalar rho() const;
411 
412  //- Return const access to the age
413  inline scalar age() const;
414 
415  //- Return const access to time spent in turbulent eddy
416  inline scalar tTurb() const;
417 
418  //- Return const access to turbulent velocity fluctuation
419  inline const vector& UTurb() const;
420 
421 
422  // Edit
423 
424  //- Return const access to moving flag
425  inline bool& moving();
426 
427  //- Return access to type id
428  inline label& typeId();
429 
430  //- Return access to number of particles
431  inline scalar& nParticle();
432 
433  //- Return access to diameter
434  inline scalar& d();
435 
436  //- Return access to target diameter
437  inline scalar& dTarget();
438 
439  //- Return access to velocity
440  inline vector& U();
441 
442  //- Return access to density
443  inline scalar& rho();
444 
445  //- Return access to the age
446  inline scalar& age();
447 
448  //- Return access to time spent in turbulent eddy
449  inline scalar& tTurb();
450 
451  //- Return access to turbulent velocity fluctuation
452  inline vector& UTurb();
453 
454 
455  // Helper functions
456 
457  //- Cell owner mass
458  inline scalar massCell(const trackingData& td) const;
459 
460  //- Particle mass
461  inline scalar mass() const;
462 
463  //- Particle moment of inertia around diameter axis
464  inline scalar momentOfInertia() const;
465 
466  //- Particle volume
467  inline scalar volume() const;
468 
469  //- Particle volume for a given diameter
470  inline static scalar volume(const scalar d);
471 
472  //- Particle projected area
473  inline scalar areaP() const;
474 
475  //- Projected area for given diameter
476  inline static scalar areaP(const scalar d);
477 
478  //- Particle surface area
479  inline scalar areaS() const;
480 
481  //- Surface area for given diameter
482  inline static scalar areaS(const scalar d);
483 
484  //- Reynolds number
485  inline scalar Re(const trackingData& td) const;
486 
487  //- Reynolds number for given conditions
488  inline static scalar Re
489  (
490  const scalar rhoc,
491  const vector& U,
492  const vector& Uc,
493  const scalar d,
494  const scalar muc
495  );
496 
497  //- Weber number
498  inline scalar We
499  (
500  const trackingData& td,
501  const scalar sigma
502  ) const;
503 
504  //- Weber number for given conditions
505  inline static scalar We
506  (
507  const scalar rhoc,
508  const vector& U,
509  const vector& Uc,
510  const scalar d,
511  const scalar sigma
512  );
513 
514  //- Eotvos number
515  inline scalar Eo
516  (
517  const trackingData& td,
518  const scalar sigma
519  ) const;
520 
521  //- Eotvos number for given conditions
522  inline static scalar Eo
523  (
524  const vector& g,
525  const scalar rho,
526  const scalar rhoc,
527  const vector& U,
528  const scalar d,
529  const scalar sigma
530  );
531 
532 
533  // Main calculation loop
534 
535  //- Set cell values
536  template<class TrackCloudType>
537  void setCellValues(TrackCloudType& cloud, trackingData& td);
538 
539  //- Apply dispersion to the carrier phase velocity and update
540  // parcel turbulence parameters
541  template<class TrackCloudType>
542  void calcDispersion
543  (
544  TrackCloudType& cloud,
545  trackingData& td,
546  const scalar dt
547  );
548 
549  //- Correct cell values using latest transfer information
550  template<class TrackCloudType>
552  (
553  TrackCloudType& cloud,
554  trackingData& td,
555  const scalar dt
556  );
557 
558  //- Update parcel properties over the time interval
559  template<class TrackCloudType>
560  void calc
561  (
562  TrackCloudType& cloud,
563  trackingData& td,
564  const scalar dt
565  );
566 
567 
568  // Tracking
569 
570  //- Move the parcel
571  template<class TrackCloudType>
572  bool move
573  (
574  TrackCloudType& cloud,
575  trackingData& td,
576  const scalar trackTime
577  );
578 
579 
580  // Patch interactions
581 
582  //- Overridable function to handle the particle hitting a patch
583  // Executed before other patch-hitting functions
584  template<class TrackCloudType>
585  bool hitPatch(TrackCloudType& cloud, trackingData& td);
586 
587  //- Overridable function to handle the particle hitting a
588  // processorPatch
589  template<class TrackCloudType>
590  void hitProcessorPatch(TrackCloudType& cloud, trackingData& td);
591 
592  //- Overridable function to handle the particle hitting a wallPatch
593  template<class TrackCloudType>
594  void hitWallPatch(TrackCloudType& cloud, trackingData& td);
595 
596  //- Transform the physical properties of the particle
597  // according to the given transformation
598  virtual void transformProperties(const transformer&);
599 
600 
601  // I-O
602 
603  //- Read
604  template<class TrackCloudType>
605  static void readFields(TrackCloudType& c);
606 
607  //- Write
608  template<class TrackCloudType>
609  static void writeFields(const TrackCloudType& c);
610 
611 
612  // Ostream Operator
613 
614  friend Ostream& operator<< <ParcelType>
615  (
616  Ostream&,
618  );
619 };
620 
621 
622 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
623 
624 } // End namespace Foam
625 
626 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
627 
628 #include "MomentumParcelI.H"
630 
631 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
632 
633 #ifdef NoRepository
634  #include "MomentumParcel.C"
635 #endif
636 
637 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
638 
639 #endif
640 
641 // ************************************************************************* //
scalar volume() const
Particle volume.
const vector & UTurb() const
Return const access to turbulent velocity fluctuation.
TemplateName(blendedSchemeBase)
bool moving_
Moving flag - tracking stopped when moving = false.
scalar d_
Diameter [m].
void hitWallPatch(TrackCloudType &cloud, trackingData &td)
Overridable function to handle the particle hitting a wallPatch.
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 areaS() const
Particle surface area.
Vector-tensor class used to perform translations, rotations and scaling operations in 3D space...
Definition: transformer.H:83
scalar dTarget_
Target diameter [m].
void setCellValues(TrackCloudType &cloud, trackingData &td)
Set cell values.
const dictionary dict_
Constant properties dictionary.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
bool hitPatch(TrackCloudType &cloud, trackingData &td)
Overridable function to handle the particle hitting a patch.
scalar age() const
Return const access to the age.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
scalar rho0() const
Return const access to the particle 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:753
label typeId_
Parcel type id.
scalar minParcelMass() const
Return const access to the minimum parcel mass.
scalar rhoMin() const
Return const access to the minimum density.
scalar We(const trackingData &td, const scalar sigma) const
Weber number.
scalar dTarget() const
Return const access to target diameter.
bool moving() const
Return const access to moving flag.
const dimensionedScalar c
Speed of light in a vacuum.
scalar d() const
Return const access to diameter.
void cellValueSourceCorrection(TrackCloudType &cloud, trackingData &td, const scalar dt)
Correct cell values using latest transfer information.
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 rho_
Density [kg/m^3].
scalar nParticle_
Number of particles in Parcel.
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
void calc(TrackCloudType &cloud, trackingData &td, const scalar dt)
Update parcel properties over the time interval.
void hitProcessorPatch(TrackCloudType &cloud, trackingData &td)
Overridable function to handle the particle hitting a.
static void readFields(TrackCloudType &c)
Read.
scalar massCell(const trackingData &td) const
Cell owner mass.
scalar nParticle() const
Return const access to number of particles.
dynamicFvMesh & mesh
scalar age_
Age [s].
scalar tTurb_
Time spent in turbulent eddy [s].
A cloud is a collection of lagrangian particles.
Definition: cloud.H:51
label parcelTypeId() const
Return const access to the parcel type id.
scalar Eo(const trackingData &td, const scalar sigma) const
Eotvos number.
scalar momentOfInertia() const
Particle moment of inertia around diameter axis.
scalar Re(const trackingData &td) const
Reynolds number.
vector U_
Velocity of Parcel [m/s].
Momentum parcel class with rotational motion (as spherical particles only) and one/two-way coupling w...
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
const dimensionedScalar mu
Atomic mass unit.
scalar rho() const
Return const access to density.
virtual void transformProperties(const transformer &)
Transform the physical properties of the particle.
void calcDispersion(TrackCloudType &cloud, trackingData &td, const scalar dt)
Apply dispersion to the carrier phase velocity and update.
MomentumParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from mesh, coordinates and topology.
scalar mass() const
Particle mass.
static void writeFields(const TrackCloudType &c)
Write.
Factory class to read-construct particles used for.
AddToPropertyList(ParcelType, " moving"+" typeId"+" nParticle"+" d"+" dTarget "+" (Ux Uy Uz)"+" rho"+" age"+" tTurb"+" (UTurbx UTurby UTurbz)")
String representation of properties.
const vector & U() const
Return const access to velocity.
const dictionary & dict() const
Return const access to the constant properties dictionary.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
label typeId() const
Return const access to type id.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
const tmp< volScalarField::Internal > & Su
Definition: alphaSuSp.H:6
bool move(TrackCloudType &cloud, trackingData &td, const scalar trackTime)
Move the parcel.
volScalarField & p
const dimensionedVector & g
scalar areaP() const
Particle projected area.
scalar tTurb() const
Return const access to time spent in turbulent eddy.
virtual autoPtr< particle > clone() const
Construct and return a (basic particle) clone.
vector UTurb_
Turbulent velocity fluctuation [m/s].
Class to hold momentum parcel constant properties.
Namespace for OpenFOAM.