KinematicCloud.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-2019 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::KinematicCloud
26 
27 Description
28  Templated base class for kinematic cloud
29 
30  - cloud function objects
31 
32  - particle forces, e.g.
33  - buoyancy
34  - drag
35  - pressure gradient
36  - ...
37 
38  - sub-models:
39  - dispersion model
40  - injection model
41  - patch interaction model
42  - stochastic collision model
43  - surface film model
44 
45 SourceFiles
46  KinematicCloudI.H
47  KinematicCloud.C
48 
49 \*---------------------------------------------------------------------------*/
50 
51 #ifndef KinematicCloud_H
52 #define KinematicCloud_H
53 
54 #include "particle.H"
55 #include "Cloud.H"
56 #include "kinematicCloud.H"
57 #include "IOdictionary.H"
58 #include "autoPtr.H"
59 #include "Random.H"
60 #include "fvMesh.H"
61 #include "volFields.H"
62 #include "fvMatrices.H"
63 #include "cloudSolution.H"
64 
65 #include "ParticleForceList.H"
67 
68 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
69 
70 namespace Foam
71 {
72 
73 // Forward declaration of classes
74 
75 class integrationScheme;
76 
77 template<class CloudType>
78 class InjectionModelList;
79 
80 template<class CloudType>
81 class DispersionModel;
82 
83 template<class CloudType>
85 
86 template<class CloudType>
87 class SurfaceFilmModel;
88 
89 template<class CloudType>
91 
92 
93 /*---------------------------------------------------------------------------*\
94  Class KinematicCloud Declaration
95 \*---------------------------------------------------------------------------*/
96 
97 template<class CloudType>
98 class KinematicCloud
99 :
100  public CloudType,
101  public kinematicCloud
102 {
103 public:
104 
105  // Public Typedefs
106 
107  //- Type of cloud this cloud was instantiated for
108  typedef CloudType cloudType;
109 
110  //- Type of parcel the cloud was instantiated for
111  typedef typename CloudType::particleType parcelType;
112 
113  //- Convenience typedef for this cloud type
115 
116  //- Force models type
118 
119  //- Function object type
121  functionType;
122 
123 
124 private:
125 
126  // Private Data
127 
128  //- Cloud copy pointer
129  autoPtr<KinematicCloud<CloudType>> cloudCopyPtr_;
130 
131 
132 protected:
133 
134  // Protected data
135 
136  //- References to the mesh and time databases
137  const fvMesh& mesh_;
138 
139  //- Dictionary of particle properties
141 
142  //- Dictionary of output properties
144 
145  //- Solution properties
147 
148  //- Parcel constant properties
149  typename parcelType::constantProperties constProps_;
150 
151  //- Sub-models dictionary
153 
154  //- Random number generator - used by some injection routines
155  mutable Random rndGen_;
156 
157  //- Cell occupancy information for each parcel, (demand driven)
159 
160  //- Cell length scale
162 
163 
164  // References to the carrier gas fields
165 
166  //- Density [kg/m^3]
167  const volScalarField& rho_;
168 
169  //- Velocity [m/s]
170  const volVectorField& U_;
171 
172  //- Dynamic viscosity [Pa.s]
173  const volScalarField& mu_;
174 
175 
176  // Environmental properties
177 
178  //- Gravity
179  const dimensionedVector& g_;
180 
181  //- Averaged ambient domain pressure
182  scalar pAmbient_;
183 
184 
185  //- Optional particle forces
186  forceType forces_;
187 
188  //- Optional cloud function objects
189  functionType functions_;
190 
191 
192  // References to the cloud sub-models
193 
194  //- Injector models
196 
197  //- Dispersion model
200 
201  //- Patch interaction model
204 
205  //- Stochastic collision model
208 
209  //- Surface film model
212 
213 
214  // Reference to the particle integration schemes
215 
216  //- Velocity integration
218 
219 
220  // Sources
221 
222  //- Momentum
224 
225  //- Coefficient for carrier phase U equation
227 
228 
229  // Initialisation
230 
231  //- Set cloud sub-models
232  void setModels();
233 
234 
235  // Cloud evolution functions
236 
237  //- Solve the cloud - calls all evolution functions
238  template<class TrackCloudType>
239  void solve
240  (
241  TrackCloudType& cloud,
242  typename parcelType::trackingData& td
243  );
244 
245  //- Build the cellOccupancy
246  void buildCellOccupancy();
247 
248  //- Update (i.e. build) the cellOccupancy if it has
249  // already been used
250  void updateCellOccupancy();
251 
252  //- Evolve the cloud
253  template<class TrackCloudType>
254  void evolveCloud
255  (
256  TrackCloudType& cloud,
257  typename parcelType::trackingData& td
258  );
259 
260  //- Post-evolve
261  void postEvolve();
262 
263  //- Reset state of cloud
265 
266 
267 public:
268 
269  // Constructors
270 
271  //- Construct given carrier gas fields
273  (
274  const word& cloudName,
275  const volScalarField& rho,
276  const volVectorField& U,
277  const volScalarField& mu,
278  const dimensionedVector& g,
279  bool readFields = true
280  );
281 
282  //- Copy constructor with new name
284  (
286  const word& name
287  );
288 
289  //- Copy constructor with new name - creates bare cloud
291  (
292  const fvMesh& mesh,
293  const word& name,
295  );
296 
297  //- Disallow default bitwise copy construction
298  KinematicCloud(const KinematicCloud&) = delete;
299 
300  //- Construct and return clone based on (this) with new name
301  virtual autoPtr<Cloud<parcelType>> clone(const word& name)
302  {
304  (
305  new KinematicCloud(*this, name)
306  );
307  }
308 
309  //- Construct and return bare clone based on (this) with new name
310  virtual autoPtr<Cloud<parcelType>> cloneBare(const word& name) const
311  {
313  (
314  new KinematicCloud(this->mesh(), name, *this)
315  );
316  }
317 
318 
319  //- Destructor
320  virtual ~KinematicCloud();
321 
322 
323  // Member Functions
324 
325  // Access
326 
327  //- Return a reference to the cloud copy
328  inline const KinematicCloud& cloudCopy() const;
329 
330 
331  // References to the mesh and databases
332 
333  //- Return reference to the mesh
334  inline const fvMesh& mesh() const;
335 
336  //- Return particle properties dictionary
337  inline const IOdictionary& particleProperties() const;
338 
339  //- Return output properties dictionary
340  inline const IOdictionary& outputProperties() const;
341 
342  //- Return non-const access to the output properties dictionary
343  inline IOdictionary& outputProperties();
344 
345  //- Return const access to the solution properties
346  inline const cloudSolution& solution() const;
347 
348  //- Return access to the solution properties
349  inline cloudSolution& solution();
350 
351  //- Return the constant properties
352  inline const typename parcelType::constantProperties&
353  constProps() const;
354 
355  //- Return access to the constant properties
356  inline typename parcelType::constantProperties& constProps();
357 
358  //- Return reference to the sub-models dictionary
359  inline const dictionary& subModelProperties() const;
360 
361 
362  // Cloud data
363 
364  //- Return reference to the random object
365  inline Random& rndGen() const;
366 
367  //- Return the cell occupancy information for each
368  // parcel, non-const access, the caller is
369  // responsible for updating it for its own purposes
370  // if particles are removed or created.
372 
373  //- Return the cell length scale
374  inline const scalarField& cellLengthScale() const;
375 
376 
377  // References to the carrier gas fields
378 
379  //- Return carrier gas velocity
380  inline const volVectorField& U() const;
381 
382  //- Return carrier gas density
383  inline const volScalarField& rho() const;
384 
385  //- Return carrier gas dynamic viscosity
386  inline const volScalarField& mu() const;
387 
388 
389  // Environmental properties
390 
391  //- Gravity
392  inline const dimensionedVector& g() const;
393 
394  //- Return const-access to the ambient pressure
395  inline scalar pAmbient() const;
396 
397  //- Return reference to the ambient pressure
398  inline scalar& pAmbient();
399 
400 
401  //- Optional particle forces
402  inline const forceType& forces() const;
403 
404  //- Return the optional particle forces
405  inline forceType& forces();
406 
407  //- Optional cloud function objects
408  inline functionType& functions();
409 
410 
411  // Sub-models
412 
413  //- Return const access to the injection model
415  injectors() const;
416 
417  //- Return reference to the injection model
419  injectors();
420 
421  //- Return const-access to the dispersion model
423  dispersion() const;
424 
425  //- Return reference to the dispersion model
427  dispersion();
428 
429  //- Return const-access to the patch interaction model
431  patchInteraction() const;
432 
433  //- Return reference to the patch interaction model
436 
437  //- Return const-access to the stochastic collision model
438  inline const
440  stochasticCollision() const;
441 
442  //- Return reference to the stochastic collision model
445 
446  //- Return const-access to the surface film model
448  surfaceFilm() const;
449 
450  //- Return reference to the surface film model
452  surfaceFilm();
453 
454 
455  // Integration schemes
456 
457  //-Return reference to velocity integration
458  inline const integrationScheme& UIntegrator() const;
459 
460 
461  // Sources
462 
463  // Momentum
464 
465  //- Return reference to momentum source
467 
468  //- Return const reference to momentum source
469  inline const volVectorField::Internal&
470  UTrans() const;
471 
472  //- Return coefficient for carrier phase U equation
474 
475  //- Return const coefficient for carrier phase U equation
476  inline const volScalarField::Internal&
477  UCoeff() const;
478 
479  //- Return tmp momentum source term
480  inline tmp<fvVectorMatrix> SU(volVectorField& U) const;
481 
482 
483  // Check
484 
485  //- Total number of parcels
486  inline label nParcels() const;
487 
488  //- Total mass in system
489  inline scalar massInSystem() const;
490 
491  //- Total linear momentum of the system
492  inline vector linearMomentumOfSystem() const;
493 
494  //- Total linear kinetic energy in the system
495  inline scalar linearKineticEnergyOfSystem() const;
496 
497  //- Total rotational kinetic energy in the system
498  inline scalar rotationalKineticEnergyOfSystem() const;
499 
500  //- Mean diameter Dij
501  inline scalar Dij(const label i, const label j) const;
502 
503  //- Max diameter
504  inline scalar Dmax() const;
505 
506 
507  // Fields
508 
509  //- Volume swept rate of parcels per cell
510  inline const tmp<volScalarField> vDotSweep() const;
511 
512  //- Return the particle volume fraction field
513  // Note: for particles belonging to this cloud only
514  inline const tmp<volScalarField> theta() const;
515 
516  //- Return the particle mass fraction field
517  // Note: for particles belonging to this cloud only
518  inline const tmp<volScalarField> alpha() const;
519 
520  //- Return the particle effective density field
521  // Note: for particles belonging to this cloud only
522  inline const tmp<volScalarField> rhoEff() const;
523 
524 
525  // Cloud evolution functions
526 
527  //- Set parcel thermo properties
529  (
530  parcelType& parcel,
531  const scalar lagrangianDt
532  );
533 
534  //- Check parcel properties
536  (
537  parcelType& parcel,
538  const scalar lagrangianDt,
539  const bool fullyDescribed
540  );
541 
542  //- Store the current cloud state
543  void storeState();
544 
545  //- Reset the current cloud to the previously stored state
546  void restoreState();
547 
548  //- Reset the cloud source terms
549  void resetSourceTerms();
550 
551  //- Relax field
552  template<class Type>
553  void relax
554  (
556  const DimensionedField<Type, volMesh>& field0,
557  const word& name
558  ) const;
559 
560  //- Scale field
561  template<class Type>
562  void scale
563  (
565  const word& name
566  ) const;
567 
568  //- Apply relaxation to (steady state) cloud sources
569  void relaxSources(const KinematicCloud<CloudType>& cloudOldTime);
570 
571  //- Apply scaling to (transient) cloud sources
572  void scaleSources();
573 
574  //- Pre-evolve
575  void preEvolve();
576 
577  //- Evolve the cloud
578  void evolve();
579 
580  //- Particle motion
581  template<class TrackCloudType>
582  void motion
583  (
584  TrackCloudType& cloud,
585  typename parcelType::trackingData& td
586  );
587 
588  //- Calculate the patch normal and velocity to interact with,
589  // accounting for patch motion if required.
590  void patchData
591  (
592  const parcelType& p,
593  const polyPatch& pp,
594  vector& normal,
595  vector& Up
596  ) const;
597 
598 
599  // Mapping
600 
601  //- Update mesh
602  void updateMesh();
603 
604  //- Remap the cells of particles corresponding to the
605  // mesh topology change with a default tracking data object
606  virtual void autoMap(const mapPolyMesh&);
607 
608 
609  // I-O
610 
611  //- Print cloud information
612  void info();
613 
614 
615  // Member Operators
616 
617  //- Disallow default bitwise assignment
618  void operator=(const KinematicCloud&) = delete;
619 };
620 
621 
622 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
623 
624 } // End namespace Foam
625 
626 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
627 
628 #include "KinematicCloudI.H"
629 
630 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
631 
632 #ifdef NoRepository
633  #include "KinematicCloud.C"
634 #endif
635 
636 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
637 
638 #endif
639 
640 // ************************************************************************* //
ParcelType particleType
Definition: Cloud.H:105
autoPtr< IOobject > clone() const
Clone.
Definition: IOobject.H:268
const SurfaceFilmModel< KinematicCloud< CloudType > > & surfaceFilm() const
Return const-access to the surface film model.
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
const word & name() const
Return name.
Definition: IOobject.H:295
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.
autoPtr< PatchInteractionModel< KinematicCloud< CloudType > > > patchInteractionModel_
Patch interaction model.
const scalarField & cellLengthScale() const
Return the cell length scale.
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const HashSet< word > &selectedFields, LIFOStack< regIOobject *> &storedObjects)
Read the selected GeometricFields of the specified type.
Definition: ReadFields.C:244
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
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:34
InjectionModelList< KinematicCloud< CloudType > > injectors_
Injector models.
scalar Dmax() const
Max diameter.
CloudType::particleType parcelType
Type of parcel the cloud was instantiated for.
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: HashTable.H:59
virtual ~KinematicCloud()
Destructor.
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 relax(DimensionedField< Type, volMesh > &field, const DimensionedField< Type, volMesh > &field0, const word &name) const
Relax field.
volVectorField::Internal & UTrans()
Return reference to momentum source.
void postEvolve()
Post-evolve.
CloudFunctionObjectList< KinematicCloud< CloudType > > functionType
Function object type.
const tmp< volScalarField > rhoEff() const
Return the particle effective density field.
List< DynamicList< parcelType * > > & cellOccupancy()
Return the cell occupancy information for each.
KinematicCloud(const word &cloudName, const volScalarField &rho, const volVectorField &U, const volScalarField &mu, const dimensionedVector &g, bool readFields=true)
Construct given carrier gas fields.
autoPtr< SurfaceFilmModel< KinematicCloud< CloudType > > > surfaceFilmModel_
Surface film model.
Random rndGen_
Random number generator - used by some injection routines.
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:158
void buildCellOccupancy()
Build the cellOccupancy.
const parcelType::constantProperties & constProps() const
Return the constant properties.
const volScalarField & rho_
Density [kg/m^3].
void operator=(const KinematicCloud &)=delete
Disallow default bitwise assignment.
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:52
const forceType & forces() const
Optional particle forces.
void updateMesh()
Update mesh.
const volScalarField & mu_
Dynamic viscosity [Pa.s].
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 string.
Definition: word.H:59
void motion(TrackCloudType &cloud, typename parcelType::trackingData &td)
Particle motion.
A cloud is a collection of lagrangian particles.
Definition: cloud.H:51
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:57
scalar linearKineticEnergyOfSystem() const
Total linear kinetic energy in the system.
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)
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...
const IOdictionary & particleProperties() const
Return particle properties dictionary.
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:51
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:78
const dimensionedScalar c
Speed of light in a vacuum.
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...
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.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
void scale(DimensionedField< Type, volMesh > &field, const word &name) const
Scale field.
List of cloud function objects.
void cloudReset(KinematicCloud< CloudType > &c)
Reset state of cloud.
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:53
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:66
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.
tmp< fvVectorMatrix > SU(volVectorField &U) const
Return tmp momentum source term.
functionType & functions()
Optional cloud function objects.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:69
void preEvolve()
Pre-evolve.
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.
label nParcels() const
Total number of parcels.
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.