All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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-2018 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 q_.boundaryFieldRef();
132 }
133 
134 
135 template<class ParcelType>
138 {
139  return fD_.boundaryFieldRef();
140 }
141 
142 
143 template<class ParcelType>
146 {
147  return rhoN_.boundaryFieldRef();
148 }
149 
150 
151 template<class ParcelType>
154 {
155  return rhoM_.boundaryFieldRef();
156 }
157 
158 
159 template<class ParcelType>
162 {
163  return linearKE_.boundaryFieldRef();
164 }
165 
166 
167 template<class ParcelType>
170 {
171  return internalE_.boundaryFieldRef();
172 }
173 
174 
175 template<class ParcelType>
178 {
179  return iDof_.boundaryFieldRef();
180 }
181 
182 
183 template<class ParcelType>
186 {
187  return momentum_.boundaryFieldRef();
188 }
189 
190 
191 template<class ParcelType>
192 inline const Foam::volScalarField&
194 {
195  return boundaryT_;
196 }
197 
198 
199 template<class ParcelType>
200 inline const Foam::volVectorField&
202 {
203  return boundaryU_;
204 }
205 
206 
207 template<class ParcelType>
210 {
211  return binaryCollisionModel_;
212 }
213 
214 
215 template<class ParcelType>
218 {
219  return binaryCollisionModel_();
220 }
221 
222 
223 template<class ParcelType>
226 {
227  return wallInteractionModel_;
228 }
229 
230 
231 template<class ParcelType>
234 {
235  return wallInteractionModel_();
236 }
237 
238 
239 template<class ParcelType>
242 {
243  return inflowBoundaryModel_;
244 }
245 
246 
247 template<class ParcelType>
250 {
251  return inflowBoundaryModel_();
252 }
253 
254 
255 template<class ParcelType>
256 inline Foam::scalar Foam::DSMCCloud<ParcelType>::massInSystem() const
257 {
258  scalar sysMass = 0.0;
259 
260  forAllConstIter(typename DSMCCloud<ParcelType>, *this, iter)
261  {
262  const ParcelType& p = iter();
263 
264  const typename ParcelType::constantProperties& cP = constProps
265  (
266  p.typeId()
267  );
268 
269  sysMass += cP.mass();
270  }
271 
272  return nParticle_*sysMass;
273 }
274 
275 
276 template<class ParcelType>
278 {
279  vector linearMomentum(Zero);
280 
281  forAllConstIter(typename DSMCCloud<ParcelType>, *this, iter)
282  {
283  const ParcelType& p = iter();
284 
285  const typename ParcelType::constantProperties& cP = constProps
286  (
287  p.typeId()
288  );
289 
290  linearMomentum += cP.mass()*p.U();
291  }
292 
293  return nParticle_*linearMomentum;
294 }
295 
296 
297 template<class ParcelType>
298 inline Foam::scalar
300 {
301  scalar linearKineticEnergy = 0.0;
302 
303  forAllConstIter(typename DSMCCloud<ParcelType>, *this, iter)
304  {
305  const ParcelType& p = iter();
306 
307  const typename ParcelType::constantProperties& cP = constProps
308  (
309  p.typeId()
310  );
311 
312  linearKineticEnergy += 0.5*cP.mass()*(p.U() & p.U());
313  }
314 
315  return nParticle_*linearKineticEnergy;
316 }
317 
318 
319 template<class ParcelType>
320 inline Foam::scalar
322 {
323  scalar internalEnergy = 0.0;
324 
325  forAllConstIter(typename DSMCCloud<ParcelType>, *this, iter)
326  {
327  const ParcelType& p = iter();
328 
329  internalEnergy += p.Ei();
330  }
331 
332  return nParticle_*internalEnergy;
333 }
334 
335 
336 template<class ParcelType>
338 (
339  scalar temperature,
340  scalar mass
341 ) const
342 {
343  return
344  2.0*sqrt(2.0*physicoChemical::k.value()*temperature/(pi*mass));
345 }
346 
347 
348 template<class ParcelType>
350 (
351  scalarField temperature,
352  scalar mass
353 ) const
354 {
355  tmp<scalarField> tfld =
356  2.0*sqrt(2.0*physicoChemical::k.value()*temperature/(pi*mass));
357  return tfld();
358 }
359 
360 
361 template<class ParcelType>
363 (
364  scalar temperature,
365  scalar mass
366 ) const
367 {
368  return sqrt(3.0*physicoChemical::k.value()*temperature/mass);
369 }
370 
371 
372 template<class ParcelType>
374 (
375  scalarField temperature,
376  scalar mass
377 ) const
378 {
379  tmp<scalarField> tfld =
380  sqrt(3.0*physicoChemical::k.value()*temperature/mass);
381  return tfld();
382 }
383 
384 
385 template<class ParcelType>
386 inline Foam::scalar
388 (
389  scalar temperature,
390  scalar mass
391 ) const
392 {
393  return sqrt(2.0*physicoChemical::k.value()*temperature/mass);
394 }
395 
396 
397 template<class ParcelType>
398 inline Foam::scalarField
400 (
401  scalarField temperature,
402  scalar mass
403 ) const
404 {
405  tmp<scalarField> tfld =
406  sqrt(2.0*physicoChemical::k.value()*temperature/mass);
407  return tfld();
408 }
409 
410 
411 template<class ParcelType>
413 {
414  return q_;
415 }
416 
417 
418 template<class ParcelType>
420 {
421  return fD_;
422 }
423 
424 
425 template<class ParcelType>
426 inline const Foam::volScalarField&
428 {
429  return rhoN_;
430 }
431 
432 
433 template<class ParcelType>
435 {
436  return rhoM_;
437 }
438 
439 
440 template<class ParcelType>
441 inline const Foam::volScalarField&
443 {
444  return dsmcRhoN_;
445 }
446 
447 
448 template<class ParcelType>
449 inline const Foam::volScalarField&
451 {
452  return linearKE_;
453 }
454 
455 
456 template<class ParcelType>
457 inline const Foam::volScalarField&
459 {
460  return internalE_;
461 }
462 
463 
464 template<class ParcelType>
465 inline const Foam::volScalarField&
467 {
468  return iDof_;
469 }
470 
471 
472 template<class ParcelType>
474 {
475  return momentum_;
476 }
477 
478 
479 template<class ParcelType>
481 {
483 }
484 
485 
486 // ************************************************************************* //
Collection of constants.
const volVectorField & momentum() const
Return the momentum density field.
Definition: DSMCCloudI.H:473
const volScalarField & boundaryT() const
Return macroscopic temperature.
Definition: DSMCCloudI.H:193
const volScalarField & linearKE() const
Return the total linear kinetic energy (translational and.
Definition: DSMCCloudI.H:450
const volScalarField & dsmcRhoN() const
Return the field of number of DSMC particles.
Definition: DSMCCloudI.H:442
const List< word > & typeIdList() const
Return the idList.
Definition: DSMCCloudI.H:57
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
vector linearMomentumOfSystem() const
Total linear momentum of the system.
Definition: DSMCCloudI.H:277
const word & cloudName() const
Return the cloud type.
Definition: DSMCCloudI.H:34
Templated inflow boundary model class.
Definition: DSMCCloud.H:61
const BinaryCollisionModel< DSMCCloud< ParcelType > > & binaryCollision() const
Return reference to binary elastic collision model.
Definition: DSMCCloudI.H:209
volScalarField::Boundary & rhoMBF()
Return non-const mass density boundary field reference.
Definition: DSMCCloudI.H:153
Templated DSMC particle collision class.
Definition: DSMCCloud.H:55
dimensionedScalar sqrt(const dimensionedScalar &ds)
scalar massInSystem() const
Total mass in system.
Definition: DSMCCloudI.H:256
volScalarField::Boundary & iDofBF()
Return non-const internal degree of freedom density boundary.
Definition: DSMCCloudI.H:177
const InflowBoundaryModel< DSMCCloud< ParcelType > > & inflowBoundary() const
Return reference to wall interaction model.
Definition: DSMCCloudI.H:241
scalar maxwellianRMSSpeed(scalar temperature, scalar mass) const
RMS particle speed.
Definition: DSMCCloudI.H:363
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
const volScalarField & internalE() const
Return the internal energy density field.
Definition: DSMCCloudI.H:458
const volVectorField & boundaryU() const
Return macroscopic velocity.
Definition: DSMCCloudI.H:201
mathematical constants.
A class for handling words, derived from string.
Definition: word.H:59
Random & rndGen()
Return references to the random object.
Definition: DSMCCloudI.H:121
const volScalarField & q() const
Return heat flux at surface field.
Definition: DSMCCloudI.H:412
void clear()
Clear the contents of the list.
Definition: ILList.C:121
volScalarField::Boundary & linearKEBF()
Return non-const linear kinetic energy density boundary.
Definition: DSMCCloudI.H:161
volVectorField::Boundary & fDBF()
Return non-const force density at boundary field reference.
Definition: DSMCCloudI.H:137
static const zero Zero
Definition: zero.H:97
errorManip< error > abort(error &err)
Definition: errorManip.H:131
scalar maxwellianMostProbableSpeed(scalar temperature, scalar mass) const
Most probable speed.
Definition: DSMCCloudI.H:388
volScalarField::Boundary & internalEBF()
Return non-const internal energy density boundary field.
Definition: DSMCCloudI.H:169
Random number generator.
Definition: Random.H:57
const volScalarField & rhoM() const
Return the particle mass density field.
Definition: DSMCCloudI.H:434
void clear()
Clear the Cloud.
Definition: DSMCCloudI.H:480
const volScalarField & rhoN() const
Return the real particle number density field.
Definition: DSMCCloudI.H:427
static const char nl
Definition: Ostream.H:260
scalar linearKineticEnergyOfSystem() const
Total linear kinetic energy in the system.
Definition: DSMCCloudI.H:299
scalar nParticle() const
Return the number of real particles represented by one.
Definition: DSMCCloudI.H:64
const volScalarField & iDof() const
Return the average internal degrees of freedom field.
Definition: DSMCCloudI.H:466
const fvMesh & mesh() const
Return references to the mesh.
Definition: DSMCCloudI.H:41
const volVectorField & fD() const
Return force density at surface field.
Definition: DSMCCloudI.H:419
const dimensionedScalar k
Boltzmann constant.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
scalar internalEnergyOfSystem() const
Total internal energy in the system.
Definition: DSMCCloudI.H:321
const List< DynamicList< ParcelType * > > & cellOccupancy() const
Return the cell occupancy addressing.
Definition: DSMCCloudI.H:72
volScalarField::Boundary & rhoNBF()
Return non-const number density boundary field reference.
Definition: DSMCCloudI.H:145
const WallInteractionModel< DSMCCloud< ParcelType > > & wallInteraction() const
Return reference to wall interaction model.
Definition: DSMCCloudI.H:225
const IOdictionary & particleProperties() const
Return particle properties dictionary.
Definition: DSMCCloudI.H:49
volScalarField & sigmaTcRMax()
Return the sigmaTcRMax field. non-const access to allow.
Definition: DSMCCloudI.H:79
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:53
Templated wall interaction model class.
Definition: DSMCCloud.H:58
scalarField & collisionSelectionRemainder()
Return the collision selection remainder field. non-const.
Definition: DSMCCloudI.H:87
scalar maxwellianAverageSpeed(scalar temperature, scalar mass) const
Average particle speed.
Definition: DSMCCloudI.H:338
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:75
const List< typename ParcelType::constantProperties > & constProps() const
Return all of the constant properties.
Definition: DSMCCloudI.H:95
volVectorField::Boundary & momentumBF()
Return non-const momentum density boundary field reference.
Definition: DSMCCloudI.H:185