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 | Copyright (C) 2011-2016 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 "cachedRandom.H"
60 #include "fvMesh.H"
61 #include "volFields.H"
62 #include "fvMatrices.H"
63 #include "IntegrationSchemesFwd.H"
64 #include "cloudSolution.H"
65 
66 #include "ParticleForceList.H"
68 
69 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
70 
71 namespace Foam
72 {
73 
74 // Forward declaration of classes
75 
76 template<class CloudType>
77 class InjectionModelList;
78 
79 template<class CloudType>
80 class DispersionModel;
81 
82 template<class CloudType>
84 
85 template<class CloudType>
86 class SurfaceFilmModel;
87 
88 template<class CloudType>
90 
91 
92 /*---------------------------------------------------------------------------*\
93  Class KinematicCloud Declaration
94 \*---------------------------------------------------------------------------*/
95 
96 template<class CloudType>
97 class KinematicCloud
98 :
99  public CloudType,
100  public kinematicCloud
101 {
102 public:
103 
104  // Public typedefs
105 
106  //- Type of cloud this cloud was instantiated for
107  typedef CloudType cloudType;
108 
109  //- Type of parcel the cloud was instantiated for
110  typedef typename CloudType::particleType parcelType;
111 
112  //- Convenience typedef for this cloud type
114 
115  //- Force models type
117 
118  //- Function object type
120  functionType;
121 
122 
123 private:
124 
125  // Private data
126 
127  //- Cloud copy pointer
128  autoPtr<KinematicCloud<CloudType>> cloudCopyPtr_;
129 
130 
131  // Private Member Functions
132 
133  //- Disallow default bitwise copy construct
135 
136  //- Disallow default bitwise assignment
137  void operator=(const KinematicCloud&);
138 
139 
140 protected:
141 
142  // Protected data
143 
144  //- References to the mesh and time databases
145  const fvMesh& mesh_;
146 
147  //- Dictionary of particle properties
149 
150  //- Dictionary of output properties
152 
153  //- Solution properties
155 
156  //- Parcel constant properties
157  typename parcelType::constantProperties constProps_;
158 
159  //- Sub-models dictionary
161 
162  //- Random number generator - used by some injection routines
164 
165  //- Cell occupancy information for each parcel, (demand driven)
167 
168  //- Cell length scale
170 
171 
172  // References to the carrier gas fields
173 
174  //- Density [kg/m3]
175  const volScalarField& rho_;
176 
177  //- Velocity [m/s]
178  const volVectorField& U_;
179 
180  //- Dynamic viscosity [Pa.s]
181  const volScalarField& mu_;
182 
183 
184  // Environmental properties
185 
186  //- Gravity
187  const dimensionedVector& g_;
188 
189  //- Averaged ambient domain pressure
190  scalar pAmbient_;
191 
192 
193  //- Optional particle forces
194  forceType forces_;
195 
196  //- Optional cloud function objects
197  functionType functions_;
198 
199 
200  // References to the cloud sub-models
201 
202  //- Injector models
204 
205  //- Dispersion model
208 
209  //- Patch interaction model
212 
213  //- Stochastic collision model
216 
217  //- Surface film model
220 
221 
222  // Reference to the particle integration schemes
223 
224  //- Velocity integration
226 
227 
228  // Sources
229 
230  //- Momentum
232 
233  //- Coefficient for carrier phase U equation
235 
236 
237  // Initialisation
238 
239  //- Set cloud sub-models
240  void setModels();
241 
242 
243  // Cloud evolution functions
244 
245  //- Solve the cloud - calls all evolution functions
246  template<class TrackData>
247  void solve(TrackData& td);
248 
249  //- Build the cellOccupancy
250  void buildCellOccupancy();
251 
252  //- Update (i.e. build) the cellOccupancy if it has
253  // already been used
254  void updateCellOccupancy();
255 
256  //- Evolve the cloud
257  template<class TrackData>
258  void evolveCloud(TrackData& td);
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  //- Construct and return clone based on (this) with new name
298  virtual autoPtr<Cloud<parcelType>> clone(const word& name)
299  {
301  (
302  new KinematicCloud(*this, name)
303  );
304  }
305 
306  //- Construct and return bare clone based on (this) with new name
307  virtual autoPtr<Cloud<parcelType>> cloneBare(const word& name) const
308  {
310  (
311  new KinematicCloud(this->mesh(), name, *this)
312  );
313  }
314 
315 
316  //- Destructor
317  virtual ~KinematicCloud();
318 
319 
320  // Member Functions
321 
322  // Access
323 
324  //- Return a reference to the cloud copy
325  inline const KinematicCloud& cloudCopy() const;
326 
327  //- Switch to specify if particles of the cloud can return
328  // non-zero wall distance values - true for kinematic parcels
329  virtual bool hasWallImpactDistance() const;
330 
331 
332  // References to the mesh and databases
333 
334  //- Return reference to the mesh
335  inline const fvMesh& mesh() const;
336 
337  //- Return particle properties dictionary
338  inline const IOdictionary& particleProperties() const;
339 
340  //- Return output properties dictionary
341  inline const IOdictionary& outputProperties() const;
342 
343  //- Return non-const access to the output properties dictionary
344  inline IOdictionary& outputProperties();
345 
346  //- Return const access to the solution properties
347  inline const cloudSolution& solution() const;
348 
349  //- Return access to the solution properties
350  inline cloudSolution& solution();
351 
352  //- Return the constant properties
353  inline const typename parcelType::constantProperties&
354  constProps() const;
355 
356  //- Return access to the constant properties
357  inline typename parcelType::constantProperties& constProps();
358 
359  //- Return reference to the sub-models dictionary
360  inline const dictionary& subModelProperties() const;
361 
362 
363  // Cloud data
364 
365  //- Return reference to the random object
366  inline cachedRandom& rndGen();
367 
368  //- Return the cell occupancy information for each
369  // parcel, non-const access, the caller is
370  // responsible for updating it for its own purposes
371  // if particles are removed or created.
373 
374  //- Return the cell length scale
375  inline const scalarField& cellLengthScale() const;
376 
377 
378  // References to the carrier gas fields
379 
380  //- Return carrier gas velocity
381  inline const volVectorField& U() const;
382 
383  //- Return carrier gas density
384  inline const volScalarField& rho() const;
385 
386  //- Return carrier gas dynamic viscosity
387  inline const volScalarField& mu() const;
388 
389 
390  // Environmental properties
391 
392  //- Gravity
393  inline const dimensionedVector& g() const;
394 
395  //- Return const-access to the ambient pressure
396  inline scalar pAmbient() const;
397 
398  //- Return reference to the ambient pressure
399  inline scalar& pAmbient();
400 
401 
402  //- Optional particle forces
403  inline const forceType& forces() const;
404 
405  //- Return the optional particle forces
406  inline forceType& forces();
407 
408  //- Optional cloud function objects
409  inline functionType& functions();
410 
411 
412  // Sub-models
413 
414  //- Return const access to the injection model
416  injectors() const;
417 
418  //- Return reference to the injection model
420  injectors();
421 
422  //- Return const-access to the dispersion model
424  dispersion() const;
425 
426  //- Return reference to the dispersion model
428  dispersion();
429 
430  //- Return const-access to the patch interaction model
432  patchInteraction() const;
433 
434  //- Return reference to the patch interaction model
437 
438  //- Return const-access to the stochastic collision model
439  inline const
441  stochasticCollision() const;
442 
443  //- Return reference to the stochastic collision model
446 
447  //- Return const-access to the surface film model
449  surfaceFilm() const;
450 
451  //- Return reference to the surface film model
453  surfaceFilm();
454 
455 
456  // Integration schemes
457 
458  //-Return reference to velocity integration
459  inline const vectorIntegrationScheme& UIntegrator() const;
460 
461 
462  // Sources
463 
464  // Momentum
465 
466  //- Return reference to momentum source
468 
469  //- Return const reference to momentum source
471  UTrans() const;
472 
473  //- Return coefficient for carrier phase U equation
475 
476  //- Return const coefficient for carrier phase U equation
478  UCoeff() const;
479 
480  //- Return tmp momentum source term
481  inline tmp<fvVectorMatrix> SU(volVectorField& U) const;
482 
483 
484  // Check
485 
486  //- Total number of parcels
487  inline label nParcels() const;
488 
489  //- Total mass in system
490  inline scalar massInSystem() const;
491 
492  //- Total linear momentum of the system
493  inline vector linearMomentumOfSystem() const;
494 
495  //- Total linear kinetic energy in the system
496  inline scalar linearKineticEnergyOfSystem() const;
497 
498  //- Total rotational kinetic energy in the system
499  inline scalar rotationalKineticEnergyOfSystem() const;
500 
501  //- Mean diameter Dij
502  inline scalar Dij(const label i, const label j) const;
503 
504  //- Max diameter
505  inline scalar Dmax() const;
506 
507 
508  // Fields
509 
510  //- Volume swept rate of parcels per cell
511  inline const tmp<volScalarField> vDotSweep() const;
512 
513  //- Return the particle volume fraction field
514  // Note: for particles belonging to this cloud only
515  inline const tmp<volScalarField> theta() const;
516 
517  //- Return the particle mass fraction field
518  // Note: for particles belonging to this cloud only
519  inline const tmp<volScalarField> alpha() const;
520 
521  //- Return the particle effective density field
522  // Note: for particles belonging to this cloud only
523  inline const tmp<volScalarField> rhoEff() const;
524 
525 
526  // Cloud evolution functions
527 
528  //- Set parcel thermo properties
530  (
531  parcelType& parcel,
532  const scalar lagrangianDt
533  );
534 
535  //- Check parcel properties
537  (
538  parcelType& parcel,
539  const scalar lagrangianDt,
540  const bool fullyDescribed
541  );
542 
543  //- Store the current cloud state
544  void storeState();
545 
546  //- Reset the current cloud to the previously stored state
547  void restoreState();
548 
549  //- Reset the cloud source terms
550  void resetSourceTerms();
551 
552  //- Relax field
553  template<class Type>
554  void relax
555  (
557  const DimensionedField<Type, volMesh>& field0,
558  const word& name
559  ) const;
560 
561  //- Scale field
562  template<class Type>
563  void scale
564  (
566  const word& name
567  ) const;
568 
569  //- Apply relaxation to (steady state) cloud sources
570  void relaxSources(const KinematicCloud<CloudType>& cloudOldTime);
571 
572  //- Apply scaling to (transient) cloud sources
573  void scaleSources();
574 
575  //- Pre-evolve
576  void preEvolve();
577 
578  //- Evolve the cloud
579  void evolve();
580 
581  //- Particle motion
582  template<class TrackData>
583  void motion(TrackData& td);
584 
585  //- Calculate the patch normal and velocity to interact with,
586  // accounting for patch motion if required.
587  void patchData
588  (
589  const parcelType& p,
590  const polyPatch& pp,
591  const scalar trackFraction,
592  const tetIndices& tetIs,
593  vector& normal,
594  vector& Up
595  ) const;
596 
597 
598  // Mapping
599 
600  //- Update mesh
601  void updateMesh();
602 
603  //- Remap the cells of particles corresponding to the
604  // mesh topology change with a default tracking data object
605  virtual void autoMap(const mapPolyMesh&);
606 
607 
608  // I-O
609 
610  //- Print cloud information
611  void info();
612 };
613 
614 
615 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
616 
617 } // End namespace Foam
618 
619 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
620 
621 #include "KinematicCloudI.H"
622 
623 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
624 
625 #ifdef NoRepository
626  #include "KinematicCloud.C"
627 #endif
628 
629 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
630 
631 #endif
632 
633 // ************************************************************************* //
ParcelType particleType
Definition: Cloud.H:114
autoPtr< DimensionedField< vector, volMesh > > UTrans_
Momentum.
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
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.
void scale(DimensionedField< Type, volMesh > &field, const word &name) const
Scale field.
const DispersionModel< KinematicCloud< CloudType > > & dispersion() const
Return const-access to the dispersion model.
void storeState()
Store the current cloud state.
scalar Dij(const label i, const label j) const
Mean diameter Dij.
const volScalarField & mu() const
Return carrier gas dynamic viscosity.
autoPtr< PatchInteractionModel< KinematicCloud< CloudType > > > patchInteractionModel_
Patch interaction model.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
InjectionModelList< KinematicCloud< CloudType > > injectors_
Injector models.
CloudType::particleType parcelType
Type of parcel the cloud was instantiated for.
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.
void setParcelThermoProperties(parcelType &parcel, const scalar lagrangianDt)
Set parcel thermo properties.
const tmp< volScalarField > theta() const
Return the particle volume fraction field.
const tmp< volScalarField > rhoEff() const
Return the particle effective density field.
DimensionedField< vector, volMesh > & UTrans()
Return reference to momentum source.
void postEvolve()
Post-evolve.
CloudFunctionObjectList< KinematicCloud< CloudType > > functionType
Function object type.
scalar rotationalKineticEnergyOfSystem() const
Total rotational kinetic energy in the system.
Random number generator.
Definition: cachedRandom.H:63
const forceType & forces() const
Optional particle forces.
List< DynamicList< parcelType * > > & cellOccupancy()
Return the cell occupancy information for each.
autoPtr< SurfaceFilmModel< KinematicCloud< CloudType > > > surfaceFilmModel_
Surface film model.
void relax(DimensionedField< Type, volMesh > &field, const DimensionedField< Type, volMesh > &field0, const word &name) const
Relax field.
void updateCellOccupancy()
Update (i.e. build) the cellOccupancy if it has.
virtual autoPtr< Cloud< parcelType > > cloneBare(const word &name) const
Construct and return bare clone based on (this) with new name.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
void buildCellOccupancy()
Build the cellOccupancy.
const scalarField & cellLengthScale() const
Return the cell length scale.
const volScalarField & rho_
Density [kg/m3].
Templated patch interaction model class.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
const cloudSolution & solution() const
Return const access to the solution properties.
void updateMesh()
Update mesh.
const volScalarField & mu_
Dynamic viscosity [Pa.s].
const vectorIntegrationScheme & UIntegrator() const
Return reference to velocity integration.
DimensionedField< scalar, volMesh > & UCoeff()
Return coefficient for carrier phase U equation.
scalar pAmbient() const
Return const-access to the ambient pressure.
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
const InjectionModelList< KinematicCloud< CloudType > > & injectors() const
Return const access to the injection model.
const dimensionedVector & g() const
Gravity.
scalar pAmbient_
Averaged ambient domain pressure.
List of injection models.
Foam::autoPtr< IOobject > clone() const
Clone.
Definition: IOobject.H:239
parcelType::constantProperties constProps_
Parcel constant properties.
const fvMesh & mesh_
References to the mesh and time databases.
const fvMesh & mesh() const
Return reference to the mesh.
Storage and named access for the indices of a tet which is part of the decomposition of a cell...
Definition: tetIndices.H:81
const StochasticCollisionModel< KinematicCloud< CloudType > > & stochasticCollision() const
Return const-access to the stochastic collision model.
void checkParcelProperties(parcelType &parcel, const scalar lagrangianDt, const bool fullyDescribed)
Check parcel properties.
const dictionary subModelProperties_
Sub-models dictionary.
scalar linearKineticEnergyOfSystem() const
Total linear kinetic energy in the system.
cachedRandom & rndGen()
Return reference to the random object.
const IOdictionary & outputProperties() const
Return output properties dictionary.
const KinematicCloud & cloudCopy() const
Return a reference to the cloud copy.
autoPtr< List< DynamicList< parcelType * > > > cellOccupancyPtr_
Cell occupancy information for each parcel, (demand driven)
Templated wall surface film model class.
void evolveCloud(TrackData &td)
Evolve the cloud.
tmp< fvVectorMatrix > SU(volVectorField &U) const
Return tmp momentum source term.
virtual void readFields()
Read the field data for the cloud of particles. Dummy at.
const IOdictionary & particleProperties() const
Return particle properties dictionary.
label nParcels() const
Total number of parcels.
A normal distribution model.
const volVectorField & U() const
Return carrier gas velocity.
void motion(TrackData &td)
Particle motion.
Virtual abstract base class for templated KinematicCloud.
cloudSolution solution_
Solution properties.
Templated base class for kinematic cloud.
Stores all relevant solution info for cloud.
Definition: cloudSolution.H:51
KinematicCloud< CloudType > kinematicCloudType
Convenience typedef for this cloud type.
void evolve()
Evolve the cloud.
IOdictionary particleProperties_
Dictionary of particle properties.
void patchData(const parcelType &p, const polyPatch &pp, const scalar trackFraction, const tetIndices &tetIs, vector &normal, vector &Up) const
Calculate the patch normal and velocity to interact with,.
const PatchInteractionModel< KinematicCloud< CloudType > > & patchInteraction() const
Return const-access to the patch interaction model.
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.
Top level model for Integration schemes.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
virtual bool hasWallImpactDistance() const
Switch to specify if particles of the cloud can return.
forceType forces_
Optional particle forces.
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.
scalar massInSystem() const
Total mass in system.
const tmp< volScalarField > vDotSweep() const
Volume swept rate of parcels per cell.
const word & cloudName() const
Return the cloud type.
Definition: DSMCCloudI.H:34
functionType functions_
Optional cloud function objects.
const parcelType::constantProperties & constProps() const
Return the constant properties.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:53
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:54
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
autoPtr< DispersionModel< KinematicCloud< CloudType > > > dispersionModel_
Dispersion model.
IOdictionary outputProperties_
Dictionary of output properties.
const tmp< volScalarField > alpha() const
Return the particle mass fraction field.
const dimensionedVector & g_
Gravity.
functionType & functions()
Optional cloud function objects.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:68
void preEvolve()
Pre-evolve.
void solve(TrackData &td)
Solve the cloud - calls all evolution functions.
vector linearMomentumOfSystem() const
Total linear momentum of the system.
const SurfaceFilmModel< KinematicCloud< CloudType > > & surfaceFilm() const
Return const-access to the surface film model.
scalar Dmax() const
Max diameter.
List of particle forces.
cachedRandom rndGen_
Random number generator - used by some injection routines.
const word & name() const
Return name.
Definition: IOobject.H:260
const volScalarField & rho() const
Return carrier gas density.
autoPtr< StochasticCollisionModel< KinematicCloud< CloudType > > > stochasticCollisionModel_
Stochastic collision model.
autoPtr< DimensionedField< scalar, volMesh > > UCoeff_
Coefficient for carrier phase U equation.
Namespace for OpenFOAM.
const dictionary & subModelProperties() const
Return reference to the sub-models dictionary.
ParticleForceList< KinematicCloud< CloudType > > forceType
Force models type.
autoPtr< vectorIntegrationScheme > UIntegrator_
Velocity integration.
scalarField cellLengthScale_
Cell length scale.