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-2018 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  // Private Member Functions
133 
134  //- Disallow default bitwise copy construct
136 
137  //- Disallow default bitwise assignment
138  void operator=(const KinematicCloud&);
139 
140 
141 protected:
142 
143  // Protected data
144 
145  //- References to the mesh and time databases
146  const fvMesh& mesh_;
147 
148  //- Dictionary of particle properties
150 
151  //- Dictionary of output properties
153 
154  //- Solution properties
156 
157  //- Parcel constant properties
158  typename parcelType::constantProperties constProps_;
159 
160  //- Sub-models dictionary
162 
163  //- Random number generator - used by some injection routines
164  mutable Random rndGen_;
165 
166  //- Cell occupancy information for each parcel, (demand driven)
168 
169  //- Cell length scale
171 
172 
173  // References to the carrier gas fields
174 
175  //- Density [kg/m3]
176  const volScalarField& rho_;
177 
178  //- Velocity [m/s]
179  const volVectorField& U_;
180 
181  //- Dynamic viscosity [Pa.s]
182  const volScalarField& mu_;
183 
184 
185  // Environmental properties
186 
187  //- Gravity
188  const dimensionedVector& g_;
189 
190  //- Averaged ambient domain pressure
191  scalar pAmbient_;
192 
193 
194  //- Optional particle forces
195  forceType forces_;
196 
197  //- Optional cloud function objects
198  functionType functions_;
199 
200 
201  // References to the cloud sub-models
202 
203  //- Injector models
205 
206  //- Dispersion model
209 
210  //- Patch interaction model
213 
214  //- Stochastic collision model
217 
218  //- Surface film model
221 
222 
223  // Reference to the particle integration schemes
224 
225  //- Velocity integration
227 
228 
229  // Sources
230 
231  //- Momentum
233 
234  //- Coefficient for carrier phase U equation
236 
237 
238  // Initialisation
239 
240  //- Set cloud sub-models
241  void setModels();
242 
243 
244  // Cloud evolution functions
245 
246  //- Solve the cloud - calls all evolution functions
247  template<class TrackCloudType>
248  void solve
249  (
250  TrackCloudType& cloud,
251  typename parcelType::trackingData& td
252  );
253 
254  //- Build the cellOccupancy
255  void buildCellOccupancy();
256 
257  //- Update (i.e. build) the cellOccupancy if it has
258  // already been used
259  void updateCellOccupancy();
260 
261  //- Evolve the cloud
262  template<class TrackCloudType>
263  void evolveCloud
264  (
265  TrackCloudType& cloud,
266  typename parcelType::trackingData& td
267  );
268 
269  //- Post-evolve
270  void postEvolve();
271 
272  //- Reset state of cloud
274 
275 
276 public:
277 
278  // Constructors
279 
280  //- Construct given carrier gas fields
282  (
283  const word& cloudName,
284  const volScalarField& rho,
285  const volVectorField& U,
286  const volScalarField& mu,
287  const dimensionedVector& g,
288  bool readFields = true
289  );
290 
291  //- Copy constructor with new name
293  (
295  const word& name
296  );
297 
298  //- Copy constructor with new name - creates bare cloud
300  (
301  const fvMesh& mesh,
302  const word& name,
304  );
305 
306  //- Construct and return clone based on (this) with new name
307  virtual autoPtr<Cloud<parcelType>> clone(const word& name)
308  {
310  (
311  new KinematicCloud(*this, name)
312  );
313  }
314 
315  //- Construct and return bare clone based on (this) with new name
316  virtual autoPtr<Cloud<parcelType>> cloneBare(const word& name) const
317  {
319  (
320  new KinematicCloud(this->mesh(), name, *this)
321  );
322  }
323 
324 
325  //- Destructor
326  virtual ~KinematicCloud();
327 
328 
329  // Member Functions
330 
331  // Access
332 
333  //- Return a reference to the cloud copy
334  inline const KinematicCloud& cloudCopy() const;
335 
336 
337  // References to the mesh and databases
338 
339  //- Return reference to the mesh
340  inline const fvMesh& mesh() const;
341 
342  //- Return particle properties dictionary
343  inline const IOdictionary& particleProperties() const;
344 
345  //- Return output properties dictionary
346  inline const IOdictionary& outputProperties() const;
347 
348  //- Return non-const access to the output properties dictionary
349  inline IOdictionary& outputProperties();
350 
351  //- Return const access to the solution properties
352  inline const cloudSolution& solution() const;
353 
354  //- Return access to the solution properties
355  inline cloudSolution& solution();
356 
357  //- Return the constant properties
358  inline const typename parcelType::constantProperties&
359  constProps() const;
360 
361  //- Return access to the constant properties
362  inline typename parcelType::constantProperties& constProps();
363 
364  //- Return reference to the sub-models dictionary
365  inline const dictionary& subModelProperties() const;
366 
367 
368  // Cloud data
369 
370  //- Return reference to the random object
371  inline Random& rndGen() const;
372 
373  //- Return the cell occupancy information for each
374  // parcel, non-const access, the caller is
375  // responsible for updating it for its own purposes
376  // if particles are removed or created.
378 
379  //- Return the cell length scale
380  inline const scalarField& cellLengthScale() const;
381 
382 
383  // References to the carrier gas fields
384 
385  //- Return carrier gas velocity
386  inline const volVectorField& U() const;
387 
388  //- Return carrier gas density
389  inline const volScalarField& rho() const;
390 
391  //- Return carrier gas dynamic viscosity
392  inline const volScalarField& mu() const;
393 
394 
395  // Environmental properties
396 
397  //- Gravity
398  inline const dimensionedVector& g() const;
399 
400  //- Return const-access to the ambient pressure
401  inline scalar pAmbient() const;
402 
403  //- Return reference to the ambient pressure
404  inline scalar& pAmbient();
405 
406 
407  //- Optional particle forces
408  inline const forceType& forces() const;
409 
410  //- Return the optional particle forces
411  inline forceType& forces();
412 
413  //- Optional cloud function objects
414  inline functionType& functions();
415 
416 
417  // Sub-models
418 
419  //- Return const access to the injection model
421  injectors() const;
422 
423  //- Return reference to the injection model
425  injectors();
426 
427  //- Return const-access to the dispersion model
429  dispersion() const;
430 
431  //- Return reference to the dispersion model
433  dispersion();
434 
435  //- Return const-access to the patch interaction model
437  patchInteraction() const;
438 
439  //- Return reference to the patch interaction model
442 
443  //- Return const-access to the stochastic collision model
444  inline const
446  stochasticCollision() const;
447 
448  //- Return reference to the stochastic collision model
451 
452  //- Return const-access to the surface film model
454  surfaceFilm() const;
455 
456  //- Return reference to the surface film model
458  surfaceFilm();
459 
460 
461  // Integration schemes
462 
463  //-Return reference to velocity integration
464  inline const integrationScheme& UIntegrator() const;
465 
466 
467  // Sources
468 
469  // Momentum
470 
471  //- Return reference to momentum source
473 
474  //- Return const reference to momentum source
475  inline const volVectorField::Internal&
476  UTrans() const;
477 
478  //- Return coefficient for carrier phase U equation
480 
481  //- Return const coefficient for carrier phase U equation
482  inline const volScalarField::Internal&
483  UCoeff() const;
484 
485  //- Return tmp momentum source term
486  inline tmp<fvVectorMatrix> SU(volVectorField& U) const;
487 
488 
489  // Check
490 
491  //- Total number of parcels
492  inline label nParcels() const;
493 
494  //- Total mass in system
495  inline scalar massInSystem() const;
496 
497  //- Total linear momentum of the system
498  inline vector linearMomentumOfSystem() const;
499 
500  //- Total linear kinetic energy in the system
501  inline scalar linearKineticEnergyOfSystem() const;
502 
503  //- Total rotational kinetic energy in the system
504  inline scalar rotationalKineticEnergyOfSystem() const;
505 
506  //- Mean diameter Dij
507  inline scalar Dij(const label i, const label j) const;
508 
509  //- Max diameter
510  inline scalar Dmax() const;
511 
512 
513  // Fields
514 
515  //- Volume swept rate of parcels per cell
516  inline const tmp<volScalarField> vDotSweep() const;
517 
518  //- Return the particle volume fraction field
519  // Note: for particles belonging to this cloud only
520  inline const tmp<volScalarField> theta() const;
521 
522  //- Return the particle mass fraction field
523  // Note: for particles belonging to this cloud only
524  inline const tmp<volScalarField> alpha() const;
525 
526  //- Return the particle effective density field
527  // Note: for particles belonging to this cloud only
528  inline const tmp<volScalarField> rhoEff() const;
529 
530 
531  // Cloud evolution functions
532 
533  //- Set parcel thermo properties
535  (
536  parcelType& parcel,
537  const scalar lagrangianDt
538  );
539 
540  //- Check parcel properties
542  (
543  parcelType& parcel,
544  const scalar lagrangianDt,
545  const bool fullyDescribed
546  );
547 
548  //- Store the current cloud state
549  void storeState();
550 
551  //- Reset the current cloud to the previously stored state
552  void restoreState();
553 
554  //- Reset the cloud source terms
555  void resetSourceTerms();
556 
557  //- Relax field
558  template<class Type>
559  void relax
560  (
562  const DimensionedField<Type, volMesh>& field0,
563  const word& name
564  ) const;
565 
566  //- Scale field
567  template<class Type>
568  void scale
569  (
571  const word& name
572  ) const;
573 
574  //- Apply relaxation to (steady state) cloud sources
575  void relaxSources(const KinematicCloud<CloudType>& cloudOldTime);
576 
577  //- Apply scaling to (transient) cloud sources
578  void scaleSources();
579 
580  //- Pre-evolve
581  void preEvolve();
582 
583  //- Evolve the cloud
584  void evolve();
585 
586  //- Particle motion
587  template<class TrackCloudType>
588  void motion
589  (
590  TrackCloudType& cloud,
591  typename parcelType::trackingData& td
592  );
593 
594  //- Calculate the patch normal and velocity to interact with,
595  // accounting for patch motion if required.
596  void patchData
597  (
598  const parcelType& p,
599  const polyPatch& pp,
600  vector& normal,
601  vector& Up
602  ) const;
603 
604 
605  // Mapping
606 
607  //- Update mesh
608  void updateMesh();
609 
610  //- Remap the cells of particles corresponding to the
611  // mesh topology change with a default tracking data object
612  virtual void autoMap(const mapPolyMesh&);
613 
614 
615  // I-O
616 
617  //- Print cloud information
618  void info();
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:270
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:297
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:137
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:60
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.
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/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: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.