DSMCCloud.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::DSMCCloud
26 
27 Description
28  Templated base class for dsmc cloud
29 
30 SourceFiles
31  DSMCCloudI.H
32  DSMCCloud.C
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #ifndef DSMCCloud_H
37 #define DSMCCloud_H
38 
39 #include "Cloud.H"
40 #include "DSMCBaseCloud.H"
41 #include "IOdictionary.H"
42 #include "autoPtr.H"
43 #include "Random.H"
44 #include "fvMesh.H"
45 #include "volFields.H"
46 #include "scalarIOField.H"
47 
48 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49 
50 namespace Foam
51 {
52 
53 // Forward declaration of classes
54 
55 template<class CloudType>
57 
58 template<class CloudType>
60 
61 template<class CloudType>
63 
64 /*---------------------------------------------------------------------------*\
65  Class DSMCCloud Declaration
66 \*---------------------------------------------------------------------------*/
67 
68 template<class ParcelType>
69 class DSMCCloud
70 :
71  public Cloud<ParcelType>,
72  public DSMCBaseCloud
73 {
74  // Private data
75 
76  //- Cloud type - used to set the name of the parcel properties
77  // dictionary by appending "Properties"
78  const word cloudName_;
79 
80  //- References to the mesh and time databases
81  const fvMesh& mesh_;
82 
83  //- Dictionary of particle properties
84  IOdictionary particleProperties_;
85 
86  //- A list of unique instances of molecule types in the
87  // simulation. The position of an entry in the list maps to
88  // the label identifying the typeId, i.e. where typeIdList_ =
89  // (N2 O2 CO2) N2 has typeId label = 0, O2 = 1, CO2 = 2.
90  List<word> typeIdList_;
91 
92  //- Number of real atoms/molecules represented by a parcel
93  scalar nParticle_;
94 
95  //- A data structure holding which particles are in which cell
96  List<DynamicList<ParcelType*>> cellOccupancy_;
97 
98  //- A field holding the value of (sigmaT * cR)max for each
99  // cell (see Bird p220). Initialised with the parcels,
100  // updated as required, and read in on start/restart.
101  volScalarField sigmaTcRMax_;
102 
103  //- A field holding the remainder from the previous collision selections
104  DimensionedField<scalar, volMesh> collisionSelectionRemainder_;
105 
106  //- Heat flux at surface field
107  volScalarField q_;
108 
109  //- Force density at surface field
110  volVectorField fD_;
111 
112  //- Number density field
113  volScalarField rhoN_;
114 
115  //- Mass density field
116  volScalarField rhoM_;
117 
118  //- Dsmc particle density field
119  volScalarField dsmcRhoN_;
120 
121  //- Linear kinetic energy density field
122  volScalarField linearKE_;
123 
124  //- Internal energy density field
125  volScalarField internalE_;
126 
127  // Internal degree of freedom density field
128  volScalarField iDof_;
129 
130  //- Momentum density field
131  volVectorField momentum_;
132 
133  //- Parcel constant properties - one for each type
135 
136  //- Random number generator
137  Random rndGen_;
138 
139 
140  // boundary value fields
141 
142  //- Boundary temperature
143  volScalarField boundaryT_;
144 
145  //- Boundary velocity
146  volVectorField boundaryU_;
147 
148 
149  // References to the cloud sub-models
150 
151  //- Binary collision model
153  binaryCollisionModel_;
154 
155  //- Wall interaction model
157  wallInteractionModel_;
158 
159  //- Inflow boundary model
161  inflowBoundaryModel_;
162 
163 
164  // Private Member Functions
165 
166  //- Build the constant properties for all of the species
167  void buildConstProps();
168 
169  //- Record which particles are in which cell
170  void buildCellOccupancy();
171 
172  //- Initialise the system
173  void initialise(const IOdictionary& dsmcInitialiseDict);
174 
175  //- Calculate collisions between molecules
176  void collisions();
177 
178  //- Reset the data accumulation field values to zero
179  void resetFields();
180 
181  //- Calculate the volume field data
182  void calculateFields();
183 
184  //- Disallow default bitwise copy construct
185  DSMCCloud(const DSMCCloud&);
186 
187  //- Disallow default bitwise assignment
188  void operator=(const DSMCCloud&);
189 
190 
191 public:
192 
193  // Constructors
194 
195  //- Construct given name and mesh, will read Parcels and fields from
196  // file
197  DSMCCloud
198  (
199  const word& cloudName,
200  const fvMesh& mesh,
201  bool readFields = true
202  );
203 
204  //- Construct given name, mesh and initialisation dictionary.
205  DSMCCloud
206  (
207  const word& cloudName,
208  const fvMesh& mesh,
209  const IOdictionary& dsmcInitialiseDict
210  );
211 
212 
213  //- Destructor
214  virtual ~DSMCCloud();
215 
216 
217  //- Type of parcel the cloud was instantiated for
218  typedef ParcelType parcelType;
219 
220 
221  // Member Functions
222 
223  // Access
224 
225  // References to the mesh and databases
226 
227  //- Return the cloud type
228  inline const word& cloudName() const;
229 
230  //- Return refernce to the mesh
231  inline const fvMesh& mesh() const;
232 
233 
234  // References to the dsmc specific data
235 
236  //- Return particle properties dictionary
237  inline const IOdictionary& particleProperties() const;
238 
239  //- Return the idList
240  inline const List<word>& typeIdList() const;
241 
242  //- Return the number of real particles represented by one
243  // parcel
244  inline scalar nParticle() const;
245 
246  //- Return the cell occupancy addressing
247  inline const List<DynamicList<ParcelType*>>&
248  cellOccupancy() const;
249 
250  //- Return the sigmaTcRMax field. non-const access to allow
251  // updating.
252  inline volScalarField& sigmaTcRMax();
253 
254  //- Return the collision selection remainder field. non-const
255  // access to allow updating.
257 
258  //- Return all of the constant properties
260  constProps() const;
261 
262  //- Return the constant properties of the given typeId
263  inline const typename ParcelType::constantProperties&
264  constProps(label typeId) const;
265 
266  //- Return refernce to the random object
267  inline Random& rndGen();
268 
269 
270  // References to the boundary fields for surface data collection
271 
272  //- Return non-const heat flux boundary field reference
273  inline volScalarField::Boundary& qBF();
274 
275  //- Return non-const force density at boundary field reference
276  inline volVectorField::Boundary& fDBF();
277 
278  //- Return non-const number density boundary field reference
279  inline volScalarField::Boundary& rhoNBF();
280 
281  //- Return non-const mass density boundary field reference
282  inline volScalarField::Boundary& rhoMBF();
283 
284  //- Return non-const linear kinetic energy density boundary
285  // field reference
286  inline volScalarField::Boundary& linearKEBF();
287 
288  //- Return non-const internal energy density boundary field
289  // reference
290  inline volScalarField::Boundary& internalEBF();
291 
292  //- Return non-const internal degree of freedom density boundary
293  // field reference
294  inline volScalarField::Boundary& iDofBF();
295 
296  //- Return non-const momentum density boundary field reference
297  inline volVectorField::Boundary& momentumBF();
298 
299 
300  // References to the macroscopic fields
301 
302  //- Return macroscopic temperature
303  inline const volScalarField& boundaryT() const;
304 
305  //- Return macroscopic velocity
306  inline const volVectorField& boundaryU() const;
307 
308  //- Return heat flux at surface field
309  inline const volScalarField& q() const;
310 
311  //- Return force density at surface field
312  inline const volVectorField& fD() const;
313 
314  //- Return the real particle number density field
315  inline const volScalarField& rhoN() const;
316 
317  //- Return the particle mass density field
318  inline const volScalarField& rhoM() const;
319 
320  //- Return the field of number of DSMC particles
321  inline const volScalarField& dsmcRhoN() const;
322 
323  //- Return the total linear kinetic energy (translational and
324  // thermal density field
325  inline const volScalarField& linearKE() const;
326 
327  //- Return the internal energy density field
328  inline const volScalarField& internalE() const;
329 
330  //- Return the average internal degrees of freedom field
331  inline const volScalarField& iDof() const;
332 
333  //- Return the momentum density field
334  inline const volVectorField& momentum() const;
335 
336 
337  // Kinetic theory helper functions
338 
339  //- Generate a random velocity sampled from the Maxwellian speed
340  // distribution
342  (
343  scalar temperature,
344  scalar mass
345  );
346 
347  //- Generate a random internal energy, sampled from the
348  // equilibrium distribution (Bird eqn 11.22 and 11.23 and
349  // adapting code from DSMC3.FOR)
351  (
352  scalar temperature,
353  direction internalDegreesOfFreedom
354  );
355 
356 
357  // From the Maxwellian distribution:
358  //- Average particle speed
359  inline scalar maxwellianAverageSpeed
360  (
361  scalar temperature,
362  scalar mass
363  ) const;
364 
366  (
367  scalarField temperature,
368  scalar mass
369  ) const;
370 
371  //- RMS particle speed
372  inline scalar maxwellianRMSSpeed
373  (
374  scalar temperature,
375  scalar mass
376  ) const;
377 
379  (
380  scalarField temperature,
381  scalar mass
382  ) const;
383 
384  //- Most probable speed
385  inline scalar maxwellianMostProbableSpeed
386  (
387  scalar temperature,
388  scalar mass
389  ) const;
390 
392  (
393  scalarField temperature,
394  scalar mass
395  ) const;
396 
397 
398  // Sub-models
399 
400  //- Return reference to binary elastic collision model
402  binaryCollision() const;
403 
404  //- Return non-const reference to binary elastic collision model
406  binaryCollision();
407 
408  //- Return reference to wall interaction model
410  wallInteraction() const;
411 
412  //- Return non-const reference to wall interaction model
414  wallInteraction();
415 
416  //- Return reference to wall interaction model
418  inflowBoundary() const;
419 
420  //- Return non-const reference to wall interaction model
422  inflowBoundary();
423 
424 
425  // Check
426 
427  //- Total mass in system
428  inline scalar massInSystem() const;
429 
430  //- Total linear momentum of the system
431  inline vector linearMomentumOfSystem() const;
432 
433  //- Total linear kinetic energy in the system
434  inline scalar linearKineticEnergyOfSystem() const;
435 
436  //- Total internal energy in the system
437  inline scalar internalEnergyOfSystem() const;
438 
439  //- Print cloud information
440  void info() const;
441 
442  //- Dump particle positions to .obj file
443  void dumpParticlePositions() const;
444 
445 
446 
447 
448  // Cloud evolution functions
449 
450  //- Add new parcel
451  void addNewParcel
452  (
453  const vector& position,
454  const vector& U,
455  const scalar Ei,
456  const label celli,
457  const label tetFacei,
458  const label tetPtI,
459  const label typeId
460  );
461 
462  //- Evolve the cloud (move, collide)
463  void evolve();
464 
465  //- Clear the Cloud
466  inline void clear();
467 
468 
469  // Mapping
470 
471  //- Remap the particles to the correct cells following mesh change
472  virtual void autoMap(const mapPolyMesh&);
473 };
474 
475 
476 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
477 
478 } // End namespace Foam
479 
480 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
481 
482 #include "DSMCCloudI.H"
483 
484 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
485 
486 #ifdef NoRepository
487  #include "DSMCCloud.C"
488 #endif
489 
490 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
491 
492 #endif
493 
494 // ************************************************************************* //
const volVectorField & fD() const
Return force density at surface field.
Definition: DSMCCloudI.H:419
void addNewParcel(const vector &position, const vector &U, const scalar Ei, const label celli, const label tetFacei, const label tetPtI, const label typeId)
Add new parcel.
Definition: DSMCCloud.C:465
virtual void autoMap(const mapPolyMesh &)
Remap the particles to the correct cells following mesh change.
Definition: DSMCCloud.C:1113
U
Definition: pEqn.H:83
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
uint8_t direction
Definition: direction.H:46
const volVectorField & momentum() const
Return the momentum density field.
Definition: DSMCCloudI.H:473
Templated inflow boundary model class.
Definition: DSMCCloud.H:61
volScalarField::Boundary & rhoMBF()
Return non-const mass density boundary field reference.
Definition: DSMCCloudI.H:153
Templated DSMC particle collision class.
Definition: DSMCCloud.H:55
const List< word > & typeIdList() const
Return the idList.
Definition: DSMCCloudI.H:57
const IOdictionary & particleProperties() const
Return particle properties dictionary.
Definition: DSMCCloudI.H:49
volScalarField::Boundary & iDofBF()
Return non-const internal degree of freedom density boundary.
Definition: DSMCCloudI.H:177
void dumpParticlePositions() const
Dump particle positions to .obj file.
Definition: DSMCCloud.C:1089
const volScalarField & boundaryT() const
Return macroscopic temperature.
Definition: DSMCCloudI.H:193
const InflowBoundaryModel< DSMCCloud< ParcelType > > & inflowBoundary() const
Return reference to wall interaction model.
Definition: DSMCCloudI.H:241
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
volScalarField::Boundary & qBF()
Return non-const heat flux boundary field reference.
Definition: DSMCCloudI.H:129
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
scalar maxwellianMostProbableSpeed(scalar temperature, scalar mass) const
Most probable speed.
Definition: DSMCCloudI.H:388
const volScalarField & dsmcRhoN() const
Return the field of number of DSMC particles.
Definition: DSMCCloudI.H:442
scalar massInSystem() const
Total mass in system.
Definition: DSMCCloudI.H:256
const volScalarField & rhoN() const
Return the real particle number density field.
Definition: DSMCCloudI.H:427
scalar linearKineticEnergyOfSystem() const
Total linear kinetic energy in the system.
Definition: DSMCCloudI.H:299
const volScalarField & internalE() const
Return the internal energy density field.
Definition: DSMCCloudI.H:458
A class for handling words, derived from string.
Definition: word.H:59
const volScalarField & q() const
Return heat flux at surface field.
Definition: DSMCCloudI.H:412
Random & rndGen()
Return refernce to the random object.
Definition: DSMCCloudI.H:121
volScalarField::Boundary & linearKEBF()
Return non-const linear kinetic energy density boundary.
Definition: DSMCCloudI.H:161
const volScalarField & rhoM() const
Return the particle mass density field.
Definition: DSMCCloudI.H:434
volVectorField::Boundary & fDBF()
Return non-const force density at boundary field reference.
Definition: DSMCCloudI.H:137
scalar nParticle() const
Return the number of real particles represented by one.
Definition: DSMCCloudI.H:64
volScalarField::Boundary & internalEBF()
Return non-const internal energy density boundary field.
Definition: DSMCCloudI.H:169
Simple random number generator.
Definition: Random.H:49
Base cloud calls templated on particle type.
Definition: Cloud.H:52
scalar internalEnergyOfSystem() const
Total internal energy in the system.
Definition: DSMCCloudI.H:321
scalar equipartitionInternalEnergy(scalar temperature, direction internalDegreesOfFreedom)
Generate a random internal energy, sampled from the.
Definition: DSMCCloud.C:1053
const volScalarField & iDof() const
Return the average internal degrees of freedom field.
Definition: DSMCCloudI.H:466
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:217
void clear()
Clear the Cloud.
Definition: DSMCCloudI.H:480
scalar maxwellianRMSSpeed(scalar temperature, scalar mass) const
RMS particle speed.
Definition: DSMCCloudI.H:363
const List< typename ParcelType::constantProperties > & constProps() const
Return all of the constant properties.
Definition: DSMCCloudI.H:95
virtual void readFields()
Read the field data for the cloud of particles. Dummy at.
const volScalarField & linearKE() const
Return the total linear kinetic energy (translational and.
Definition: DSMCCloudI.H:450
virtual ~DSMCCloud()
Destructor.
Definition: DSMCCloud.C:954
const fvMesh & mesh() const
Return refernce to the mesh.
Definition: DSMCCloudI.H:41
const List< DynamicList< ParcelType * > > & cellOccupancy() const
Return the cell occupancy addressing.
Definition: DSMCCloudI.H:72
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
const volVectorField & boundaryU() const
Return macroscopic velocity.
Definition: DSMCCloudI.H:201
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
void evolve()
Evolve the cloud (move, collide)
Definition: DSMCCloud.C:961
const word & cloudName() const
Return the cloud type.
Definition: DSMCCloudI.H:34
Virtual abstract base class for templated DSMCCloud.
Definition: DSMCBaseCloud.H:48
volScalarField::Boundary & rhoNBF()
Return non-const number density boundary field reference.
Definition: DSMCCloudI.H:145
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:53
volScalarField & sigmaTcRMax()
Return the sigmaTcRMax field. non-const access to allow.
Definition: DSMCCloudI.H:79
scalar maxwellianAverageSpeed(scalar temperature, scalar mass) const
Average particle speed.
Definition: DSMCCloudI.H:338
Templated wall interaction model class.
Definition: DSMCCloud.H:58
vector linearMomentumOfSystem() const
Total linear momentum of the system.
Definition: DSMCCloudI.H:277
const WallInteractionModel< DSMCCloud< ParcelType > > & wallInteraction() const
Return reference to wall interaction model.
Definition: DSMCCloudI.H:225
scalarField & collisionSelectionRemainder()
Return the collision selection remainder field. non-const.
Definition: DSMCCloudI.H:87
vector equipartitionLinearVelocity(scalar temperature, scalar mass)
Generate a random velocity sampled from the Maxwellian speed.
Definition: DSMCCloud.C:1035
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:68
void info() const
Print cloud information.
Definition: DSMCCloud.C:991
const BinaryCollisionModel< DSMCCloud< ParcelType > > & binaryCollision() const
Return reference to binary elastic collision model.
Definition: DSMCCloudI.H:209
Namespace for OpenFOAM.
volVectorField::Boundary & momentumBF()
Return non-const momentum density boundary field reference.
Definition: DSMCCloudI.H:185