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