KinematicCloud.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-2020 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::KinematicCloud
29 
30 Group
31  grpLagrangianIntermediateClouds
32 
33 Description
34  Templated base class for kinematic cloud
35 
36  - cloud function objects
37 
38  - particle forces, e.g.
39  - buoyancy
40  - drag
41  - pressure gradient
42  - ...
43 
44  - sub-models:
45  - dispersion model
46  - injection model
47  - patch interaction model
48  - stochastic collision model
49  - surface film model
50 
51 SourceFiles
52  KinematicCloudI.H
53  KinematicCloud.C
54 
55 \*---------------------------------------------------------------------------*/
56 
57 #ifndef KinematicCloud_H
58 #define KinematicCloud_H
59 
60 #include "particle.H"
61 #include "Cloud.H"
62 #include "kinematicCloud.H"
63 #include "IOdictionary.H"
64 #include "autoPtr.H"
65 #include "Random.H"
66 #include "fvMesh.H"
67 #include "volFields.H"
68 #include "fvMatrices.H"
69 #include "cloudSolution.H"
70 
71 #include "ParticleForceList.H"
73 
74 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
75 
76 namespace Foam
77 {
78 
79 // Forward declaration of classes
80 
81 class integrationScheme;
82 
83 template<class CloudType>
84 class InjectionModelList;
85 
86 template<class CloudType>
87 class DispersionModel;
88 
89 template<class CloudType>
91 
92 template<class CloudType>
93 class SurfaceFilmModel;
94 
95 template<class CloudType>
97 
98 template<class CloudType>
99 class PackingModel;
101 template<class CloudType>
102 class DampingModel;
103 
104 template<class CloudType>
105 class IsotropyModel;
106 
107 
108 /*---------------------------------------------------------------------------*\
109  Class KinematicCloud Declaration
110 \*---------------------------------------------------------------------------*/
111 
112 template<class CloudType>
113 class KinematicCloud
114 :
115  public CloudType,
116  public kinematicCloud
117 {
118 public:
119 
120  // Public typedefs
121 
122  //- Type of cloud this cloud was instantiated for
123  typedef CloudType cloudType;
124 
125  //- Type of parcel the cloud was instantiated for
126  typedef typename CloudType::particleType parcelType;
127 
128  //- Convenience typedef for this cloud type
131  //- Force models type
133 
134  //- Function object type
136  functionType;
137 
138 
139 private:
140 
141  // Private Data
142 
143  //- Cloud copy pointer
144  autoPtr<KinematicCloud<CloudType>> cloudCopyPtr_;
145 
146 
147  // Private Member Functions
148 
149  //- No copy construct
150  KinematicCloud(const KinematicCloud&) = delete;
151 
152  //- No copy assignment
153  void operator=(const KinematicCloud&) = delete;
154 
155 
156 protected:
157 
158  // Protected Data
159 
160  //- References to the mesh and time databases
161  const fvMesh& mesh_;
162 
163  //- Dictionary of particle properties
165 
166  //- Dictionary of output properties
168 
169  //- Solution properties
171 
172  //- Parcel constant properties
173  typename parcelType::constantProperties constProps_;
175  //- Sub-models dictionary
177 
178  //- Random number generator - used by some injection routines
179  mutable Random rndGen_;
180 
181  //- Cell occupancy information for each parcel, (demand driven)
183 
184  //- Cell length scale
186 
187 
188  // References to the carrier gas fields
190  //- Density [kg/m3]
191  const volScalarField& rho_;
192 
193  //- Velocity [m/s]
195 
196  //- Dynamic viscosity [Pa.s]
197  const volScalarField& mu_;
198 
200  // Environmental properties
201 
202  //- Gravity
203  const dimensionedVector& g_;
205  //- Averaged ambient domain pressure
206  scalar pAmbient_;
207 
208 
209  //- Optional particle forces
211 
212  //- Optional cloud function objects
215 
216  // References to the cloud sub-models
217 
218  //- Injector models
220 
221  //- Dispersion model
224 
225  //- Patch interaction model
228 
229  //- Stochastic collision model
233  //- Surface film model
236 
237  //- Packing model
241  //- Damping model
244 
245  //- Exchange model
248 
249 
250  // Reference to the particle integration schemes
252  //- Velocity integration
254 
255 
256  // Sources
257 
258  //- Momentum
260 
261  //- Coefficient for carrier phase U equation
263 
265  // Initialisation
266 
267  //- Set cloud sub-models
268  void setModels();
269 
271  // Cloud evolution functions
272 
273  //- Solve the cloud - calls all evolution functions
274  template<class TrackCloudType>
275  void solve
276  (
277  TrackCloudType& cloud,
278  typename parcelType::trackingData& td
279  );
280 
281  //- Build the cellOccupancy
283 
284  //- Update (i.e. build) the cellOccupancy if it has
285  // already been used
286  void updateCellOccupancy();
287 
288  //- Evolve the cloud
289  template<class TrackCloudType>
290  void evolveCloud
291  (
292  TrackCloudType& cloud,
293  typename parcelType::trackingData& td
294  );
295 
296  //- Post-evolve
297  void postEvolve(const typename parcelType::trackingData& td);
298 
299  //- Reset state of cloud
301 
302 
303 public:
304 
305  // Public Data
307  //- Flag to write log into Info
308  bool log;
309 
310 
311  // Constructors
312 
313  //- Construct given carrier gas fields
315  (
316  const word& cloudName,
317  const volScalarField& rho,
318  const volVectorField& U,
319  const volScalarField& mu,
320  const dimensionedVector& g,
321  bool readFields = true
322  );
323 
324  //- Copy constructor with new name
326  (
328  const word& name
329  );
330 
331  //- Copy constructor with new name - creates bare cloud
333  (
334  const fvMesh& mesh,
335  const word& name,
337  );
338 
339  //- Construct and return clone based on (this) with new name
340  virtual autoPtr<Cloud<parcelType>> clone(const word& name)
341  {
343  (
344  new KinematicCloud(*this, name)
345  );
346  }
347 
348  //- Construct and return bare clone based on (this) with new name
349  virtual autoPtr<Cloud<parcelType>> cloneBare(const word& name) const
350  {
352  (
353  new KinematicCloud(this->mesh(), name, *this)
354  );
355  }
356 
357 
358  //- Destructor
359  virtual ~KinematicCloud() = default;
360 
361 
362  // Member Functions
363 
364  // Access
365 
366  //- Return a reference to the cloud copy
367  inline const KinematicCloud& cloudCopy() const;
368 
369 
370  // References to the mesh and databases
371 
372  //- Return reference to the mesh
373  inline const fvMesh& mesh() const;
374 
375  //- Return particle properties dictionary
376  inline const IOdictionary& particleProperties() const;
377 
378  //- Return output properties dictionary
379  inline const IOdictionary& outputProperties() const;
380 
381  //- Return non-const access to the output properties dictionary
382  inline IOdictionary& outputProperties();
383 
384  //- Return const access to the solution properties
385  inline const cloudSolution& solution() const;
386 
387  //- Return access to the solution properties
388  inline cloudSolution& solution();
389 
390  //- Return the constant properties
391  inline const typename parcelType::constantProperties&
392  constProps() const;
393 
394  //- Return access to the constant properties
395  inline typename parcelType::constantProperties& constProps();
396 
397  //- Return reference to the sub-models dictionary
398  inline const dictionary& subModelProperties() const;
399 
400 
401  // Cloud data
402 
403  //- Return reference to the random object
404  inline Random& rndGen() const;
405 
406  //- Return the cell occupancy information for each
407  // parcel, non-const access, the caller is
408  // responsible for updating it for its own purposes
409  // if particles are removed or created.
411 
412  //- Return the cell length scale
413  inline const scalarField& cellLengthScale() const;
414 
415 
416  // References to the carrier gas fields
417 
418  //- Return carrier gas velocity
419  inline const volVectorField& U() const;
420 
421  //- Return carrier gas density
422  inline const volScalarField& rho() const;
423 
424  //- Return carrier gas dynamic viscosity
425  inline const volScalarField& mu() const;
426 
427 
428  // Environmental properties
429 
430  //- Gravity
431  inline const dimensionedVector& g() const;
432 
433  //- Return const-access to the ambient pressure
434  inline scalar pAmbient() const;
435 
436  //- Return reference to the ambient pressure
437  inline scalar& pAmbient();
438 
439 
440  //- Optional particle forces
441  inline const forceType& forces() const;
442 
443  //- Return the optional particle forces
444  inline forceType& forces();
445 
446  //- Optional cloud function objects
447  inline functionType& functions();
448 
449 
450  // Sub-models
451 
452  //- Return const access to the injection model
454  injectors() const;
455 
456  //- Return reference to the injection model
458  injectors();
459 
460  //- Return const-access to the dispersion model
462  dispersion() const;
463 
464  //- Return reference to the dispersion model
466  dispersion();
467 
468  //- Return const-access to the patch interaction model
470  patchInteraction() const;
471 
472  //- Return reference to the patch interaction model
475 
476  //- Return const-access to the stochastic collision model
477  inline const
479  stochasticCollision() const;
480 
481  //- Return reference to the stochastic collision model
484 
485  //- Return const-access to the surface film model
487  surfaceFilm() const;
488 
489  //- Return reference to the surface film model
491  surfaceFilm();
492 
493 
494  //- Return const access to the packing model
496  packingModel() const;
497 
498  //- Return a reference to the packing model
500  packingModel();
501 
502  //- Return const access to the damping model
504  dampingModel() const;
505 
506  //- Return a reference to the damping model
508  dampingModel();
509 
510  //- Return const access to the isotropy model
512  isotropyModel() const;
513 
514  //- Return a reference to the isotropy model
516  isotropyModel();
517 
518 
519  // Integration schemes
520 
521  //-Return reference to velocity integration
522  inline const integrationScheme& UIntegrator() const;
523 
524 
525  // Sources
526 
527  // Momentum
528 
529  //- Return reference to momentum source
531 
532  //- Return const reference to momentum source
533  inline const volVectorField::Internal&
534  UTrans() const;
535 
536  //- Return coefficient for carrier phase U equation
538 
539  //- Return const coefficient for carrier phase U equation
540  inline const volScalarField::Internal&
541  UCoeff() const;
542 
543  //- Return tmp momentum source term (compressible)
544  inline tmp<fvVectorMatrix> SU
545  (
546  volVectorField& U,
547  bool incompressible = false
548  ) const;
549 
550 
551  // Check
552 
553  //- Total number of parcels
554  virtual label nParcels() const
555  {
556  return CloudType::nParcels();
557  }
558 
559  //- Total mass in system
560  inline scalar massInSystem() const;
561 
562  //- Total linear momentum of the system
563  inline vector linearMomentumOfSystem() const;
564 
565  //- Average particle per parcel
566  inline scalar totalParticlePerParcel() const;
567 
568  //- Total linear kinetic energy in the system
569  inline scalar linearKineticEnergyOfSystem() const;
570 
571  //- Total rotational kinetic energy in the system
572  inline scalar rotationalKineticEnergyOfSystem() const;
573 
574  //- Mean diameter Dij
575  inline scalar Dij(const label i, const label j) const;
576 
577  //- Max diameter
578  inline scalar Dmax() const;
579 
580 
581  // Fields
582 
583  //- Volume swept rate of parcels per cell
584  inline const tmp<volScalarField> vDotSweep() const;
585 
586  //- Return the particle volume fraction field
587  // Note: for particles belonging to this cloud only
588  inline const tmp<volScalarField> theta() const;
589 
590  //- Return the particle mass fraction field
591  // Note: for particles belonging to this cloud only
592  inline const tmp<volScalarField> alpha() const;
593 
594  //- Return the particle effective density field
595  // Note: for particles belonging to this cloud only
596  inline const tmp<volScalarField> rhoEff() const;
597 
598 
599  // Cloud evolution functions
600 
601  //- Set parcel thermo properties
603  (
604  parcelType& parcel,
605  const scalar lagrangianDt
606  );
607 
608  //- Check parcel properties
610  (
611  parcelType& parcel,
612  const scalar lagrangianDt,
613  const bool fullyDescribed
614  );
615 
616  //- Store the current cloud state
617  void storeState();
618 
619  //- Reset the current cloud to the previously stored state
620  void restoreState();
621 
622  //- Reset the cloud source terms
623  void resetSourceTerms();
624 
625  //- Relax field
626  template<class Type>
627  void relax
628  (
630  const DimensionedField<Type, volMesh>& field0,
631  const word& name
632  ) const;
633 
634  //- Scale field
635  template<class Type>
636  void scale
637  (
639  const word& name
640  ) const;
641 
642  //- Apply relaxation to (steady state) cloud sources
643  void relaxSources(const KinematicCloud<CloudType>& cloudOldTime);
644 
645  //- Apply scaling to (transient) cloud sources
646  void scaleSources();
647 
648  //- Pre-evolve
649  void preEvolve
650  (
651  const typename parcelType::trackingData& td
652  );
653 
654  //- Evolve the cloud
655  void evolve();
656 
657  //- Particle motion
658  template<class TrackCloudType>
659  void motion
660  (
661  TrackCloudType& cloud,
662  typename parcelType::trackingData& td
663  );
664 
665  //- Calculate the patch normal and velocity to interact with,
666  // accounting for patch motion if required.
667  void patchData
668  (
669  const parcelType& p,
670  const polyPatch& pp,
671  vector& normal,
672  vector& Up
673  ) const;
674 
675 
676  // Mapping
677 
678  //- Update mesh
679  void updateMesh();
680 
681  //- Remap the cells of particles corresponding to the
682  // mesh topology change with a default tracking data object
683  virtual void autoMap(const mapPolyMesh&);
684 
685 
686  // I-O
687 
688  //- Print cloud information
689  void info();
690 
691  //- Read particle fields from objects in the obr registry
692  virtual void readObjects(const objectRegistry& obr);
693 
694  //- Write particle fields as objects into the obr registry
695  virtual void writeObjects(objectRegistry& obr) const;
696 };
697 
698 
699 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
700 
701 } // End namespace Foam
702 
703 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
704 
705 #include "KinematicCloudI.H"
706 
707 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
708 
709 #ifdef NoRepository
710  #include "KinematicCloud.C"
711 #endif
712 
713 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
714 
715 #endif
716 
717 // ************************************************************************* //
ParcelType particleType
Definition: Cloud.H:130
tmp< fvVectorMatrix > SU(volVectorField &U, bool incompressible=false) const
Return tmp momentum source term (compressible)
virtual label nParcels() const
Return the number of particles in the cloud.
Definition: Cloud.H:197
autoPtr< IOobject > clone() const
Clone.
Definition: IOobject.H:570
autoPtr< IsotropyModel< KinematicCloud< CloudType > > > isotropyModel_
Exchange model.
const SurfaceFilmModel< KinematicCloud< CloudType > > & surfaceFilm() const
Return const-access to the surface film model.
rDeltaTY field()
void scaleSources()
Apply scaling to (transient) cloud sources.
void setModels()
Set cloud sub-models.
const volVectorField & U_
Velocity [m/s].
void relaxSources(const KinematicCloud< CloudType > &cloudOldTime)
Apply relaxation to (steady state) cloud sources.
const StochasticCollisionModel< KinematicCloud< CloudType > > & stochasticCollision() const
Return const-access to the stochastic collision model.
void storeState()
Store the current cloud state.
scalar totalParticlePerParcel() const
Average particle per parcel.
autoPtr< PatchInteractionModel< KinematicCloud< CloudType > > > patchInteractionModel_
Patch interaction model.
const scalarField & cellLengthScale() const
Return the cell length scale.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:120
scalar rotationalKineticEnergyOfSystem() const
Total rotational kinetic energy in the system.
const tmp< volScalarField > theta() const
Return the particle volume fraction field.
const word & cloudName() const
Return the cloud type.
Definition: DSMCCloudI.H:30
InjectionModelList< KinematicCloud< CloudType > > injectors_
Injector models.
scalar Dmax() const
Max diameter.
CloudType::particleType parcelType
Type of parcel the cloud was instantiated for.
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:150
const InjectionModelList< KinematicCloud< CloudType > > & injectors() const
Return const access to the injection model.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:56
const KinematicCloud & cloudCopy() const
Return a reference to the cloud copy.
void evolveCloud(TrackCloudType &cloud, typename parcelType::trackingData &td)
Evolve the cloud.
const volScalarField & rho() const
Return carrier gas density.
autoPtr< volVectorField::Internal > UTrans_
Momentum.
const fvMesh & mesh() const
Return reference to the mesh.
void setParcelThermoProperties(parcelType &parcel, const scalar lagrangianDt)
Set parcel thermo properties.
void postEvolve(const typename parcelType::trackingData &td)
Post-evolve.
void relax(DimensionedField< Type, volMesh > &field, const DimensionedField< Type, volMesh > &field0, const word &name) const
Relax field.
volVectorField::Internal & UTrans()
Return reference to momentum source.
Base class for packing models.
CloudFunctionObjectList< KinematicCloud< CloudType > > functionType
Function object type.
virtual void writeObjects(objectRegistry &obr) const
Write particle fields as objects into the obr registry.
const IsotropyModel< KinematicCloud< CloudType > > & isotropyModel() const
Return const access to the isotropy model.
const tmp< volScalarField > rhoEff() const
Return the particle effective density field.
const DampingModel< KinematicCloud< CloudType > > & dampingModel() const
Return const access to the damping model.
autoPtr< DampingModel< KinematicCloud< CloudType > > > dampingModel_
Damping model.
List< DynamicList< parcelType * > > & cellOccupancy()
Return the cell occupancy information for each.
autoPtr< SurfaceFilmModel< KinematicCloud< CloudType > > > surfaceFilmModel_
Surface film model.
Random rndGen_
Random number generator - used by some injection routines.
virtual label nParcels() const
Total number of parcels.
void updateCellOccupancy()
Update (i.e. build) the cellOccupancy if it has.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:157
void buildCellOccupancy()
Build the cellOccupancy.
const parcelType::constantProperties & constProps() const
Return the constant properties.
const volScalarField & rho_
Density [kg/m3].
Templated patch interaction model class.
scalar massInSystem() const
Total mass in system.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:50
autoPtr< PackingModel< KinematicCloud< CloudType > > > packingModel_
Packing model.
const forceType & forces() const
Optional particle forces.
void updateMesh()
Update mesh.
const volScalarField & mu_
Dynamic viscosity [Pa.s].
bool log
Flag to write log into Info.
virtual ~KinematicCloud()=default
Destructor.
volScalarField::Internal & UCoeff()
Return coefficient for carrier phase U equation.
void resetSourceTerms()
Reset the cloud source terms.
CloudType cloudType
Type of cloud this cloud was instantiated for.
A class for handling words, derived from Foam::string.
Definition: word.H:63
const PackingModel< KinematicCloud< CloudType > > & packingModel() const
Return const access to the packing model.
void motion(TrackCloudType &cloud, typename parcelType::trackingData &td)
Particle motion.
A cloud is a registry collection of lagrangian particles.
Definition: cloud.H:53
virtual void readObjects(const objectRegistry &obr)
Read particle fields from objects in the obr registry.
scalar pAmbient_
Averaged ambient domain pressure.
List of injection models.
parcelType::constantProperties constProps_
Parcel constant properties.
const fvMesh & mesh_
References to the mesh and time databases.
void checkParcelProperties(parcelType &parcel, const scalar lagrangianDt, const bool fullyDescribed)
Check parcel properties.
const dictionary subModelProperties_
Sub-models dictionary.
Random number generator.
Definition: Random.H:55
scalar linearKineticEnergyOfSystem() const
Total linear kinetic energy in the system.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
const tmp< volScalarField > vDotSweep() const
Volume swept rate of parcels per cell.
const DispersionModel< KinematicCloud< CloudType > > & dispersion() const
Return const-access to the dispersion model.
scalar pAmbient() const
Return const-access to the ambient pressure.
autoPtr< List< DynamicList< parcelType * > > > cellOccupancyPtr_
Cell occupancy information for each parcel, (demand driven)
Templated wall surface film model class.
const cloudSolution & solution() const
Return const access to the solution properties.
Base for a set of schemes which integrate simple ODEs which arise from semi-implcit rate expressions...
Base class for collisional return-to-isotropy models.
const IOdictionary & particleProperties() const
Return particle properties dictionary.
void preEvolve(const typename parcelType::trackingData &td)
Pre-evolve.
Virtual abstract base class for templated KinematicCloud.
cloudSolution solution_
Solution properties.
Templated base class for kinematic cloud.
void patchData(const parcelType &p, const polyPatch &pp, vector &normal, vector &Up) const
Calculate the patch normal and velocity to interact with,.
Stores all relevant solution info for cloud.
Definition: cloudSolution.H:49
KinematicCloud< CloudType > kinematicCloudType
Convenience typedef for this cloud type.
const integrationScheme & UIntegrator() const
Return reference to velocity integration.
void evolve()
Evolve the cloud.
IOdictionary particleProperties_
Dictionary of particle properties.
autoPtr< volScalarField::Internal > UCoeff_
Coefficient for carrier phase U equation.
Templated stochastic collision model class.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:79
const dimensionedScalar c
Speed of light in a vacuum.
Base class for dispersion modelling.
void info()
Print cloud information.
vector linearMomentumOfSystem() const
Total linear momentum of the system.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
forceType forces_
Optional particle forces.
const dimensionedVector & g() const
Gravity.
void restoreState()
Reset the current cloud to the previously stored state.
A special matrix type and solver, designed for finite volume solutions of scalar equations.
const IOdictionary & outputProperties() const
Return output properties dictionary.
const volVectorField & U() const
Return carrier gas velocity.
functionType functions_
Optional cloud function objects.
const PatchInteractionModel< KinematicCloud< CloudType > > & patchInteraction() const
Return const-access to the patch interaction model.
virtual autoPtr< Cloud< parcelType > > cloneBare(const word &name) const
Construct and return bare clone based on (this) with new name.
autoPtr< integrationScheme > UIntegrator_
Velocity integration.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
void scale(DimensionedField< Type, volMesh > &field, const word &name) const
Scale field.
List of cloud function objects.
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject *> &storedObjects)
Read the selected GeometricFields of the templated type.
void cloudReset(KinematicCloud< CloudType > &c)
Reset state of cloud.
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Registry of regIOobjects.
virtual void autoMap(const mapPolyMesh &)
Remap the cells of particles corresponding to the.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:69
scalar Dij(const label i, const label j) const
Mean diameter Dij.
autoPtr< DispersionModel< KinematicCloud< CloudType > > > dispersionModel_
Dispersion model.
const volScalarField & mu() const
Return carrier gas dynamic viscosity.
IOdictionary outputProperties_
Dictionary of output properties.
const dimensionedVector & g_
Gravity.
const tmp< volScalarField > alpha() const
Return the particle mass fraction field.
functionType & functions()
Optional cloud function objects.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:67
Base class for collisional damping models.
void solve(TrackCloudType &cloud, typename parcelType::trackingData &td)
Solve the cloud - calls all evolution functions.
Random & rndGen() const
Return reference to the random object.
List of particle forces.
autoPtr< StochasticCollisionModel< KinematicCloud< CloudType > > > stochasticCollisionModel_
Stochastic collision model.
Namespace for OpenFOAM.
const dictionary & subModelProperties() const
Return reference to the sub-models dictionary.
ParticleForceList< KinematicCloud< CloudType > > forceType
Force models type.
scalarField cellLengthScale_
Cell length scale.