DSMCCloudI.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 \*---------------------------------------------------------------------------*/
25 
26 #include "constants.H"
27 
28 using namespace Foam::constant;
29 using namespace Foam::constant::mathematical;
30 
31 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
32 
33 template<class ParcelType>
35 {
36  return cloudName_;
37 }
38 
39 
40 template<class ParcelType>
42 {
43  return mesh_;
44 }
45 
46 
47 template<class ParcelType>
48 inline const Foam::IOdictionary&
50 {
51  return particleProperties_;
52 }
53 
54 
55 template<class ParcelType>
56 inline const Foam::List<Foam::word>&
58 {
59  return typeIdList_;
60 }
61 
62 
63 template<class ParcelType>
64 inline Foam::scalar Foam::DSMCCloud<ParcelType>::nParticle() const
65 {
66  return nParticle_;
67 }
68 
69 
70 template<class ParcelType>
73 {
74  return cellOccupancy_;
75 }
76 
77 
78 template<class ParcelType>
80 {
81  return sigmaTcRMax_;
82 }
83 
84 
85 template<class ParcelType>
86 inline Foam::scalarField&
88 {
89  return collisionSelectionRemainder_;
90 }
91 
92 
93 template<class ParcelType>
96 {
97  return constProps_;
98 }
99 
100 
101 template<class ParcelType>
102 inline const typename ParcelType::constantProperties&
104 (
105  label typeId
106 ) const
107 {
108  if (typeId < 0 || typeId >= constProps_.size())
109  {
111  << "constantProperties for requested typeId index "
112  << typeId << " do not exist" << nl
113  << abort(FatalError);
114  }
115 
116  return constProps_[typeId];
117 }
118 
119 
120 template<class ParcelType>
122 {
123  return rndGen_;
124 }
125 
126 
127 template<class ParcelType>
130 {
131  return stdNormal_;
132 }
133 
134 
135 template<class ParcelType>
138 {
139  return q_.boundaryFieldRef();
140 }
141 
142 
143 template<class ParcelType>
146 {
147  return fD_.boundaryFieldRef();
148 }
149 
150 
151 template<class ParcelType>
154 {
155  return rhoN_.boundaryFieldRef();
156 }
157 
158 
159 template<class ParcelType>
162 {
163  return rhoM_.boundaryFieldRef();
164 }
165 
166 
167 template<class ParcelType>
170 {
171  return linearKE_.boundaryFieldRef();
172 }
173 
174 
175 template<class ParcelType>
178 {
179  return internalE_.boundaryFieldRef();
180 }
181 
182 
183 template<class ParcelType>
186 {
187  return iDof_.boundaryFieldRef();
188 }
189 
190 
191 template<class ParcelType>
194 {
195  return momentum_.boundaryFieldRef();
196 }
197 
198 
199 template<class ParcelType>
200 inline const Foam::volScalarField&
202 {
203  return boundaryT_;
204 }
205 
206 
207 template<class ParcelType>
208 inline const Foam::volVectorField&
210 {
211  return boundaryU_;
212 }
213 
214 
215 template<class ParcelType>
218 {
219  return binaryCollisionModel_;
220 }
221 
222 
223 template<class ParcelType>
226 {
227  return binaryCollisionModel_();
228 }
229 
230 
231 template<class ParcelType>
234 {
235  return wallInteractionModel_;
236 }
237 
238 
239 template<class ParcelType>
242 {
243  return wallInteractionModel_();
244 }
245 
246 
247 template<class ParcelType>
250 {
251  return inflowBoundaryModel_;
252 }
253 
254 
255 template<class ParcelType>
258 {
259  return inflowBoundaryModel_();
260 }
261 
262 
263 template<class ParcelType>
264 inline Foam::scalar Foam::DSMCCloud<ParcelType>::massInSystem() const
265 {
266  scalar sysMass = 0.0;
267 
268  forAllConstIter(typename DSMCCloud<ParcelType>, *this, iter)
269  {
270  const ParcelType& p = iter();
271 
272  const typename ParcelType::constantProperties& cP = constProps
273  (
274  p.typeId()
275  );
276 
277  sysMass += cP.mass();
278  }
279 
280  return nParticle_*sysMass;
281 }
282 
283 
284 template<class ParcelType>
286 {
287  vector linearMomentum(Zero);
288 
289  forAllConstIter(typename DSMCCloud<ParcelType>, *this, iter)
290  {
291  const ParcelType& p = iter();
292 
293  const typename ParcelType::constantProperties& cP = constProps
294  (
295  p.typeId()
296  );
297 
298  linearMomentum += cP.mass()*p.U();
299  }
300 
301  return nParticle_*linearMomentum;
302 }
303 
304 
305 template<class ParcelType>
306 inline Foam::scalar
308 {
309  scalar linearKineticEnergy = 0.0;
310 
311  forAllConstIter(typename DSMCCloud<ParcelType>, *this, iter)
312  {
313  const ParcelType& p = iter();
314 
315  const typename ParcelType::constantProperties& cP = constProps
316  (
317  p.typeId()
318  );
319 
320  linearKineticEnergy += 0.5*cP.mass()*(p.U() & p.U());
321  }
322 
323  return nParticle_*linearKineticEnergy;
324 }
325 
326 
327 template<class ParcelType>
328 inline Foam::scalar
330 {
331  scalar internalEnergy = 0.0;
332 
333  forAllConstIter(typename DSMCCloud<ParcelType>, *this, iter)
334  {
335  const ParcelType& p = iter();
336 
337  internalEnergy += p.Ei();
338  }
339 
340  return nParticle_*internalEnergy;
341 }
342 
343 
344 template<class ParcelType>
346 (
347  scalar temperature,
348  scalar mass
349 ) const
350 {
351  return
352  2.0*sqrt(2.0*physicoChemical::k.value()*temperature/(pi*mass));
353 }
354 
355 
356 template<class ParcelType>
358 (
359  scalarField temperature,
360  scalar mass
361 ) const
362 {
363  tmp<scalarField> tfld =
364  2.0*sqrt(2.0*physicoChemical::k.value()*temperature/(pi*mass));
365  return tfld();
366 }
367 
368 
369 template<class ParcelType>
371 (
372  scalar temperature,
373  scalar mass
374 ) const
375 {
376  return sqrt(3.0*physicoChemical::k.value()*temperature/mass);
377 }
378 
379 
380 template<class ParcelType>
382 (
383  scalarField temperature,
384  scalar mass
385 ) const
386 {
387  tmp<scalarField> tfld =
388  sqrt(3.0*physicoChemical::k.value()*temperature/mass);
389  return tfld();
390 }
391 
392 
393 template<class ParcelType>
394 inline Foam::scalar
396 (
397  scalar temperature,
398  scalar mass
399 ) const
400 {
401  return sqrt(2.0*physicoChemical::k.value()*temperature/mass);
402 }
403 
404 
405 template<class ParcelType>
406 inline Foam::scalarField
408 (
409  scalarField temperature,
410  scalar mass
411 ) const
412 {
413  tmp<scalarField> tfld =
414  sqrt(2.0*physicoChemical::k.value()*temperature/mass);
415  return tfld();
416 }
417 
418 
419 template<class ParcelType>
421 {
422  return q_;
423 }
424 
425 
426 template<class ParcelType>
428 {
429  return fD_;
430 }
431 
432 
433 template<class ParcelType>
434 inline const Foam::volScalarField&
436 {
437  return rhoN_;
438 }
439 
440 
441 template<class ParcelType>
443 {
444  return rhoM_;
445 }
446 
447 
448 template<class ParcelType>
449 inline const Foam::volScalarField&
451 {
452  return dsmcRhoN_;
453 }
454 
455 
456 template<class ParcelType>
457 inline const Foam::volScalarField&
459 {
460  return linearKE_;
461 }
462 
463 
464 template<class ParcelType>
465 inline const Foam::volScalarField&
467 {
468  return internalE_;
469 }
470 
471 
472 template<class ParcelType>
473 inline const Foam::volScalarField&
475 {
476  return iDof_;
477 }
478 
479 
480 template<class ParcelType>
482 {
483  return momentum_;
484 }
485 
486 
487 template<class ParcelType>
489 {
491 }
492 
493 
494 // ************************************************************************* //
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
Templated DSMC particle collision class.
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
scalar massInSystem() const
Total mass in system.
Definition: DSMCCloudI.H:264
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
volScalarField::Boundary & internalEBF()
Return non-const internal energy density boundary field.
Definition: DSMCCloudI.H:177
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
const volScalarField & iDof() const
Return the average internal degrees of freedom field.
Definition: DSMCCloudI.H:474
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
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
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
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
scalar maxwellianAverageSpeed(scalar temperature, scalar mass) const
Average particle speed.
Definition: DSMCCloudI.H:346
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
const WallInteractionModel< DSMCCloud< ParcelType > > & wallInteraction() const
Return reference to wall interaction model.
Definition: DSMCCloudI.H:233
Generic GeometricBoundaryField class.
Generic GeometricField class.
void clear()
Clear the contents of the list.
Definition: ILList.C:121
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.
Standard normal distribution. Not selectable.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
Random number generator.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
mathematical constants.
const dimensionedScalar k
Boltzmann constant.
Collection of constants.
static const zero Zero
Definition: zero.H:97
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
errorManip< error > abort(error &err)
Definition: errorManip.H:131
dimensionedScalar sqrt(const dimensionedScalar &ds)
error FatalError
static const char nl
Definition: Ostream.H:266
volScalarField & p