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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2016-2024 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Class
28  Foam::KinematicParcel
29 
30 Group
31  grpLagrangianIntermediateParcels
32 
33 Description
34  Kinematic parcel class with rotational motion (as spherical
35  particles only) and one/two-way coupling with the continuous
36  phase.
37 
38  Sub-models include:
39  - drag
40  - turbulent dispersion
41  - wall interactions
42 
43 SourceFiles
44  KinematicParcelI.H
45  KinematicParcel.C
46  KinematicParcelIO.C
47 
48 \*---------------------------------------------------------------------------*/
49 
50 #ifndef KinematicParcel_H
51 #define KinematicParcel_H
52 
53 #include "particle.H"
54 #include "IOstream.H"
55 #include "interpolation.H"
56 #include "demandDrivenEntry.H"
57 #include "labelFieldIOField.H"
58 #include "vectorFieldIOField.H"
59 
60 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
61 
62 namespace Foam
63 {
64 
65 template<class ParcelType>
66 class KinematicParcel;
67 
68 template<class Type>
69 class AveragingMethod;
70 
71 // Forward declaration of friend functions
72 
73 template<class ParcelType>
74 Ostream& operator<<
75 (
76  Ostream&,
78 );
79 
80 /*---------------------------------------------------------------------------*\
81  Class KinematicParcel Declaration
82 \*---------------------------------------------------------------------------*/
83 
84 template<class ParcelType>
85 class KinematicParcel
86 :
87  public ParcelType
88 {
89  // Private Data
90 
91  //- Number of particle tracking attempts before we assume that it stalls
92  static label maxTrackAttempts;
93 
94 
95 public:
96 
97  //- Size in bytes of the fields
98  static const std::size_t sizeofFields;
99 
100 
101  //- Class to hold kinematic particle constant properties
102  class constantProperties
103  {
104  protected:
105 
106  // Protected Data
107 
108  //- Constant properties dictionary
109  const dictionary dict_;
110 
111 
112  private:
113 
114  // Private Data
115 
116  //- Parcel type id - used for post-processing to flag the type
117  //- of parcels issued by this cloud
118  demandDrivenEntry<label> parcelTypeId_;
119 
120  //- Minimum density [kg/m3]
122 
123  //- Particle density [kg/m3] (constant)
125 
126  //- Minimum parcel mass [kg]
127  demandDrivenEntry<scalar> minParcelMass_;
128 
129 
130  public:
131 
132  // Constructors
133 
134  //- Default construct
135  inline constantProperties();
136 
137  //- Copy constructor
139 
140  //- Construct from dictionary
141  inline constantProperties(const dictionary& parentDict);
142 
143 
144  // Member Functions
145 
146  //- Return const access to the constant properties dictionary
147  inline const dictionary& dict() const;
148 
149  //- Return const access to the parcel type id
150  inline label parcelTypeId() const;
151 
152  //- Return const access to the minimum density
153  inline scalar rhoMin() const;
154 
155  //- Return const access to the particle density
156  inline scalar rho0() const;
157 
158  //- Return const access to the minimum parcel mass
159  inline scalar minParcelMass() const;
160  };
161 
162 
163  class trackingData
164  :
165  public ParcelType::trackingData
166  {
167  public:
168 
169  enum trackPart
170  {
174  };
175 
176 
177  private:
178 
179  // Private Data
180 
181  // Interpolators for continuous phase fields
182 
183  //- Density interpolator
184  autoPtr<interpolation<scalar>> rhoInterp_;
185 
186  //- Velocity interpolator
188 
189  //- Dynamic viscosity interpolator
191 
192 
193  // Cached continuous phase properties
194 
195  //- Density [kg/m3]
196  scalar rhoc_;
197 
198  //- Velocity [m/s]
199  vector Uc_;
201  //- Viscosity [Pa.s]
202  scalar muc_;
203 
204 
205  // MPPIC Averages
206 
207  //- Volume average
208  autoPtr<AveragingMethod<scalar>> volumeAverage_;
209 
210  //- Radius average [ volume^(1/3) ]
211  autoPtr<AveragingMethod<scalar>> radiusAverage_;
212 
213  //- Density average
214  autoPtr<AveragingMethod<scalar>> rhoAverage_;
215 
216  //- Velocity average
218 
219  //- Magnitude velocity squared average
220  autoPtr<AveragingMethod<scalar>> uSqrAverage_;
221 
222  //- Frequency average
223  autoPtr<AveragingMethod<scalar>> frequencyAverage_;
224 
225  //- Mass average
226  autoPtr<AveragingMethod<scalar>> massAverage_;
227 
228 
229  //- Local gravitational or other body-force acceleration
230  const vector& g_;
231 
232  //- label specifying which part of the integration
233  //- algorithm is taking place
234  trackPart part_;
235 
236 
237  public:
238 
239  // Constructors
240 
241  //- Construct from components
242  template <class TrackCloudType>
243  inline trackingData
244  (
245  const TrackCloudType& cloud,
247  );
248 
249 
250  // Member Functions
251 
252  //- Return const access to the interpolator for continuous
253  //- phase density field
254  inline const interpolation<scalar>& rhoInterp() const;
255 
256  //- Return const access to the interpolator for continuous
257  //- phase velocity field
258  inline const interpolation<vector>& UInterp() const;
259 
260  //- Return const access to the interpolator for continuous
261  //- phase dynamic viscosity field
262  inline const interpolation<scalar>& muInterp() const;
263 
264  //- Return the continuous phase density
265  inline scalar rhoc() const;
266 
267  //- Access the continuous phase density
268  inline scalar& rhoc();
269 
270  //- Return the continuous phase velocity
271  inline const vector& Uc() const;
272 
273  //- Access the continuous phase velocity
274  inline vector& Uc();
275 
276  //- Return the continuous phase viscosity
277  inline scalar muc() const;
278 
279  //- Access the continuous phase viscosity
280  inline scalar& muc();
281 
282  // Return const access to the gravitational acceleration vector
283  inline const vector& g() const;
284 
285  //- Return the part of the tracking operation taking place
286  inline trackPart part() const;
287 
288  //- Return access to the part of the tracking operation taking place
289  inline trackPart& part();
290 
291  //- Update the MPPIC averages
292  template<class TrackCloudType>
293  inline void updateAverages(const TrackCloudType& cloud);
294 
295  };
296 
297 
298 protected:
299 
300  // Protected Data
301 
302  // Parcel properties
303 
304  //- Active flag - tracking inactive when active = false
305  // Store as label for data alignment, but only has bool states
306  label active_;
307 
308  //- Parcel type id
309  label typeId_;
310 
311  //- Number of particles in Parcel
312  scalar nParticle_;
313 
314  //- Diameter [m]
315  scalar d_;
316 
317  //- Target diameter [m]
318  scalar dTarget_;
319 
320  //- Velocity of Parcel [m/s]
321  vector U_;
322 
323  //- Density [kg/m3]
324  scalar rho_;
325 
326  //- Age [s]
327  scalar age_;
328 
329  //- Time spent in turbulent eddy [s]
330  scalar tTurb_;
331 
332  //- Turbulent velocity fluctuation [m/s]
333  vector UTurb_;
334 
335  //- Velocity correction due to collisions MPPIC [m/s]
337 
338 
339  // Protected Member Functions
340 
341  //- Calculate new particle velocity
342  template<class TrackCloudType>
343  const vector calcVelocity
344  (
345  TrackCloudType& cloud,
346  trackingData& td,
347  const scalar dt, // timestep
348  const scalar Re, // Reynolds number
349  const scalar mu, // local carrier viscosity
350  const scalar mass, // mass
351  const vector& Su, // explicit particle momentum source
352  vector& dUTrans, // momentum transfer to carrier
353  scalar& Spu // linearised drag coefficient
354  ) const;
355 
356 
357 public:
358 
359  // Static data members
360 
361  //- Runtime type information
362  TypeName("KinematicParcel");
363 
364  //- String representation of properties
366  (
367  ParcelType,
368  " active"
369  + " typeId"
370  + " nParticle"
371  + " d"
372  + " dTarget"
373  + " (Ux Uy Uz)"
374  + " rho"
375  + " age"
376  + " tTurb"
377  + " (UTurbx UTurby UTurbz)"
378  + " (UCorrectx UCorrecty UCorrectz)"
379  );
380 
381 
382  // Constructors
383 
384  //- Construct from mesh, coordinates and topology
385  // Other properties initialised as null
386  inline KinematicParcel
387  (
388  const polyMesh& mesh,
389  const barycentric& coordinates,
390  const label celli,
391  const label tetFacei,
392  const label tetPti
393  );
394 
395  //- Construct from a position and a cell, searching for the rest of the
396  // required topology. Other properties are initialised as null.
398  (
399  const polyMesh& mesh,
400  const vector& position,
401  const label celli
402  );
403 
404  //- Construct from components
405  inline KinematicParcel
406  (
407  const polyMesh& mesh,
408  const barycentric& coordinates,
409  const label celli,
410  const label tetFacei,
411  const label tetPti,
412  const label typeId,
413  const scalar nParticle0,
414  const scalar d0,
415  const scalar dTarget0,
416  const vector& U0,
417  const constantProperties& constProps
418  );
419 
420  //- Construct from Istream
422  (
423  const polyMesh& mesh,
424  Istream& is,
425  bool readFields = true,
426  bool newFormat = true
427  );
428 
429  //- Construct as a copy
431 
432  //- Construct as a copy
434 
435  //- Return a (basic particle) clone
436  virtual autoPtr<particle> clone() const
437  {
438  return particle::Clone(*this);
439  }
440 
441  //- Return a (basic particle) clone
442  virtual autoPtr<particle> clone(const polyMesh& mesh) const
443  {
444  return particle::Clone(*this, mesh);
445  }
446 
447  //- Factory class to read-construct particles (for parallel transfer)
448  class iNew
449  {
450  const polyMesh& mesh_;
451 
452  public:
453 
454  iNew(const polyMesh& mesh)
455  :
456  mesh_(mesh)
457  {}
458 
459  autoPtr<KinematicParcel<ParcelType>> operator()(Istream& is) const
460  {
461  return autoPtr<KinematicParcel<ParcelType>>
462  (
463  new KinematicParcel<ParcelType>(mesh_, is, true)
464  );
465  }
466  };
467 
468 
469  // Member Functions
470 
471  // Access
472 
473  //- Return const access to active flag
474  inline bool active() const;
475 
476  //- Return const access to type id
477  inline label typeId() const;
478 
479  //- Return const access to number of particles
480  inline scalar nParticle() const;
481 
482  //- Return const access to diameter
483  inline scalar d() const;
484 
485  //- Return const access to target diameter
486  inline scalar dTarget() const;
487 
488  //- Return const access to velocity
489  inline const vector& U() const;
490 
491  //- Return const access to density
492  inline scalar rho() const;
493 
494  //- Return const access to the age
495  inline scalar age() const;
496 
497  //- Return const access to time spent in turbulent eddy
498  inline scalar tTurb() const;
499 
500  //- Return const access to turbulent velocity fluctuation
501  inline const vector& UTurb() const;
502 
503  //- Return const access to correction velocity
504  inline const vector& UCorrect() const;
505 
506 
507  // Edit
508 
509  //- Set active flag to the specified state
510  inline void active(const bool state);
511 
512  //- Return access to type id
513  inline label& typeId();
514 
515  //- Return access to number of particles
516  inline scalar& nParticle();
517 
518  //- Return access to diameter
519  inline scalar& d();
520 
521  //- Return access to target diameter
522  inline scalar& dTarget();
523 
524  //- Return access to velocity
525  inline vector& U();
526 
527  //- Return access to density
528  inline scalar& rho();
529 
530  //- Return access to the age
531  inline scalar& age();
532 
533  //- Return access to time spent in turbulent eddy
534  inline scalar& tTurb();
535 
536  //- Return access to turbulent velocity fluctuation
537  inline vector& UTurb();
538 
539  //- Return access to correction velocity
540  inline vector& UCorrect();
541 
542 
543  // Helper functions
544 
545  //- Cell owner mass
546  inline scalar massCell(const trackingData& td) const;
547 
548  //- Particle mass
549  inline scalar mass() const;
550 
551  //- Particle moment of inertia around diameter axis
552  inline scalar momentOfInertia() const;
553 
554  //- Particle volume
555  inline scalar volume() const;
556 
557  //- Particle volume for a given diameter
558  inline static scalar volume(const scalar d);
559 
560  //- Particle projected area
561  inline scalar areaP() const;
562 
563  //- Projected area for given diameter
564  inline static scalar areaP(const scalar d);
565 
566  //- Particle surface area
567  inline scalar areaS() const;
568 
569  //- Surface area for given diameter
570  inline static scalar areaS(const scalar d);
571 
572  //- Reynolds number
573  inline scalar Re(const trackingData& td) const;
574 
575  //- Reynolds number for given conditions
576  inline static scalar Re
577  (
578  const scalar rhoc,
579  const vector& U,
580  const vector& Uc,
581  const scalar d,
582  const scalar muc
583  );
584 
585  //- Weber number
586  inline scalar We
587  (
588  const trackingData& td,
589  const scalar sigma
590  ) const;
592  //- Weber number for given conditions
593  inline static scalar We
594  (
595  const scalar rhoc,
596  const vector& U,
597  const vector& Uc,
598  const scalar d,
599  const scalar sigma
600  );
601 
602  //- Eotvos number
603  inline scalar Eo
604  (
605  const trackingData& td,
606  const scalar sigma
607  ) const;
608 
609  //- Eotvos number for given conditions
610  inline static scalar Eo
611  (
612  const vector& g,
613  const scalar rho,
614  const scalar rhoc,
615  const vector& U,
616  const scalar d,
617  const scalar sigma
618  );
619 
620 
621  // Main calculation loop
622 
623  //- Set cell values
624  template<class TrackCloudType>
625  void setCellValues(TrackCloudType& cloud, trackingData& td);
626 
627  //- Apply dispersion to the carrier phase velocity and update
628  // parcel turbulence parameters
629  template<class TrackCloudType>
630  void calcDispersion
631  (
632  TrackCloudType& cloud,
633  trackingData& td,
634  const scalar dt
635  );
636 
637  //- Correct cell values using latest transfer information
638  template<class TrackCloudType>
640  (
641  TrackCloudType& cloud,
642  trackingData& td,
643  const scalar dt
644  );
645 
646  //- Update parcel properties over the time interval
647  template<class TrackCloudType>
648  void calc
649  (
650  TrackCloudType& cloud,
651  trackingData& td,
652  const scalar dt
653  );
654 
655  //- Correct U following MP-PIC sub-models
656  template<class TrackCloudType>
657  void calcUCorrection
658  (
659  TrackCloudType& cloud,
660  trackingData& td,
661  const scalar dt
662  );
663 
664 
665  // Tracking
666 
667  //- Move the parcel
668  template<class TrackCloudType>
669  bool move
670  (
671  TrackCloudType& cloud,
672  trackingData& td,
673  const scalar trackTime
674  );
675 
676 
677  // Patch interactions
678 
679  //- Overridable function to handle the particle hitting a patch
680  // Executed before other patch-hitting functions
681  template<class TrackCloudType>
682  bool hitPatch(TrackCloudType& cloud, trackingData& td);
683 
684  //- Overridable function to handle the particle hitting a
685  // processorPatch
686  template<class TrackCloudType>
687  void hitProcessorPatch(TrackCloudType& cloud, trackingData& td);
688 
689  //- Overridable function to handle the particle hitting a wallPatch
690  template<class TrackCloudType>
691  void hitWallPatch(TrackCloudType& cloud, trackingData& td);
692 
693  //- Transform the physical properties of the particle
694  // according to the given transformation tensor
695  virtual void transformProperties(const tensor& T);
696 
697  //- Transform the physical properties of the particle
698  // according to the given separation vector
699  virtual void transformProperties(const vector& separation);
700 
701 
702  // I-O
703 
704  //- Read
705  template<class TrackCloudType>
706  static void readFields(TrackCloudType& c);
707 
708  //- Write
709  template<class TrackCloudType>
710  static void writeFields(const TrackCloudType& c);
711 
712  //- Write individual parcel properties to stream
713  void writeProperties
714  (
715  Ostream& os,
716  const wordRes& filters,
717  const word& delim,
718  const bool namesOnly = false
719  ) const;
720 
721  //- Read particle fields as objects from the obr registry
722  template<class CloudType>
723  static void readObjects(CloudType& c, const objectRegistry& obr);
724 
725  //- Write particle fields as objects into the obr registry
726  template<class CloudType>
727  static void writeObjects(const CloudType& c, objectRegistry& obr);
728 
729 
730  // Ostream Operator
731 
732  friend Ostream& operator<< <ParcelType>
733  (
734  Ostream&,
736  );
737 };
738 
739 
740 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
741 
742 } // End namespace Foam
743 
744 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
745 
746 #include "KinematicParcelI.H"
748 
749 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
750 
751 #ifdef NoRepository
752  #include "KinematicParcel.C"
753 #endif
754 
755 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
756 
757 #endif
758 
759 // ************************************************************************* //
scalar massCell(const trackingData &td) const
Cell owner mass.
const vector & U() const
Return const access to velocity.
static void writeFields(const TrackCloudType &c)
Write.
scalar tTurb() const
Return const access to time spent in turbulent eddy.
scalar rhoMin() const
Return const access to the minimum density.
zeroField Su
Definition: alphaSuSp.H:1
const vector & UTurb() const
Return const access to turbulent velocity fluctuation.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
Kinematic parcel class with rotational motion (as spherical particles only) and one/two-way coupling ...
autoPtr< KinematicParcel< ParcelType > > operator()(Istream &is) const
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
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 the source to the destination (recursively if necessary).
Definition: POSIX.C:1063
const interpolation< vector > & UInterp() const
Return const access to the interpolator for continuous phase velocity field.
scalar rho0() const
Return const access to the particle density.
void setCellValues(TrackCloudType &cloud, trackingData &td)
Set cell values.
scalar areaP() const
Particle projected area.
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
KinematicParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from mesh, coordinates and topology.
const vector & UCorrect() const
Return const access to correction velocity.
scalar dTarget() const
Return const access to target diameter.
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.
trackPart part() const
Return the part of the tracking operation taking place.
void cellValueSourceCorrection(TrackCloudType &cloud, trackingData &td, const scalar dt)
Correct cell values using latest transfer information.
dynamicFvMesh & mesh
scalar Eo(const trackingData &td, const scalar sigma) const
Eotvos number.
AddToPropertyList(ParcelType, " active"+" typeId"+" nParticle"+" d"+" dTarget"+" (Ux Uy Uz)"+" rho"+" age"+" tTurb"+" (UTurbx UTurby UTurbz)"+" (UCorrectx UCorrecty UCorrectz)")
String representation of properties.
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.
A class for handling words, derived from Foam::string.
Definition: word.H:63
static void writeObjects(const CloudType &c, objectRegistry &obr)
Write particle fields as objects into the obr registry.
const dictionary dict_
Constant properties dictionary.
A cloud is a registry collection of lagrangian particles.
Definition: cloud.H:53
label typeId_
Parcel type id.
void writeProperties(Ostream &os, const wordRes &filters, const word &delim, const bool namesOnly=false) const
Write individual parcel properties to stream.
vector UTurb_
Turbulent velocity fluctuation [m/s].
Base class for lagrangian averaging methods.
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:53
static void readFields(TrackCloudType &c)
Read.
void calc(TrackCloudType &cloud, trackingData &td, const scalar dt)
Update parcel properties over the time interval.
Vector< scalar > vector
Definition: vector.H:57
bool active() const
Return const access to active flag.
static const std::size_t sizeofFields
Size in bytes of the fields.
static autoPtr< particle > Clone(const Derived &p)
Clone a particle.
Definition: particle.H:552
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.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
const uniformDimensionedVectorField & g
static void readObjects(CloudType &c, const objectRegistry &obr)
Read particle fields as objects from the obr registry.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
scalar rhoc() const
Return the continuous phase density.
OBJstream os(runTime.globalPath()/outputName)
scalar muc() const
Return the continuous phase viscosity.
scalar nParticle_
Number of particles in Parcel.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
scalar minParcelMass() const
Return const access to the minimum parcel mass.
trackingData(const TrackCloudType &cloud, trackPart part=tpLinearTrack)
Construct from components.
const dimensionedScalar mu
Atomic mass unit.
iNew(const polyMesh &mesh)
vector U_
Velocity of Parcel [m/s].
virtual autoPtr< particle > clone() const
Return a (basic particle) clone.
PtrList< coordinateSystem > coordinates(solidRegions.size())
const dimensionedScalar c
Speed of light in a vacuum.
scalar volume() const
Particle volume.
const interpolation< scalar > & rhoInterp() const
Return const access to the interpolator for continuous phase density field.
const interpolation< scalar > & muInterp() const
Return const access to the interpolator for continuous phase dynamic viscosity field.
scalar Re(const trackingData &td) const
Reynolds number.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:75
bool move(TrackCloudType &cloud, trackingData &td, const scalar trackTime)
Move the parcel.
label active_
Active flag - tracking inactive when active = false.
volScalarField & p
Registry of regIOobjects.
label typeId() const
Return const access to type id.
Tensor of scalars, i.e. Tensor<scalar>.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:67
void hitProcessorPatch(TrackCloudType &cloud, trackingData &td)
Overridable function to handle the particle hitting a.
scalar areaS() const
Particle surface area.
bool hitPatch(TrackCloudType &cloud, trackingData &td)
Overridable function to handle the particle hitting a patch.
const vector & Uc() const
Return the continuous phase velocity.
scalar d_
Diameter [m].
vector UCorrect_
Velocity correction due to collisions MPPIC [m/s].
void calcUCorrection(TrackCloudType &cloud, trackingData &td, const scalar dt)
Correct U following MP-PIC sub-models.
scalar We(const trackingData &td, const scalar sigma) const
Weber number.
void calcDispersion(TrackCloudType &cloud, trackingData &td, const scalar dt)
Apply dispersion to the carrier phase velocity and update.
Namespace for OpenFOAM.
const dictionary & dict() const
Return const access to the constant properties dictionary.
void hitWallPatch(TrackCloudType &cloud, trackingData &td)
Overridable function to handle the particle hitting a wallPatch.
void updateAverages(const TrackCloudType &cloud)
Update the MPPIC averages.