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