DSMCCloud.C
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  Copyright (C) 2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "DSMCCloud.H"
30 #include "BinaryCollisionModel.H"
31 #include "WallInteractionModel.H"
32 #include "InflowBoundaryModel.H"
33 #include "constants.H"
36 
37 using namespace Foam::constant;
38 
39 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
40 
41 template<class ParcelType>
43 {
44  Info<< nl << "Constructing constant properties for" << endl;
45  constProps_.setSize(typeIdList_.size());
46 
47  dictionary moleculeProperties
48  (
49  particleProperties_.subDict("moleculeProperties")
50  );
51 
52  forAll(typeIdList_, i)
53  {
54  const word& id = typeIdList_[i];
55 
56  Info<< " " << id << endl;
57 
58  const dictionary& molDict(moleculeProperties.subDict(id));
59 
60  constProps_[i] = typename ParcelType::constantProperties(molDict);
61  }
62 }
63 
64 
65 template<class ParcelType>
67 {
68  for (auto& list : cellOccupancy_)
69  {
70  list.clear();
71  }
72 
73  for (ParcelType& p : *this)
74  {
75  cellOccupancy_[p.cell()].append(&p);
76  }
77 }
78 
79 
80 template<class ParcelType>
82 (
83  const IOdictionary& dsmcInitialiseDict
84 )
85 {
86  Info<< nl << "Initialising particles" << endl;
87 
88  const scalar temperature
89  (
90  dsmcInitialiseDict.get<scalar>("temperature")
91  );
92 
93  const vector velocity(dsmcInitialiseDict.get<vector>("velocity"));
94 
95  const dictionary& numberDensitiesDict
96  (
97  dsmcInitialiseDict.subDict("numberDensities")
98  );
99 
100  List<word> molecules(numberDensitiesDict.toc());
101 
102  Field<scalar> numberDensities(molecules.size());
103 
104  forAll(molecules, i)
105  {
106  numberDensities[i] = numberDensitiesDict.get<scalar>(molecules[i]);
107  }
108 
109  numberDensities /= nParticle_;
110 
111  forAll(mesh_.cells(), celli)
112  {
113  List<tetIndices> cellTets = polyMeshTetDecomposition::cellTetIndices
114  (
115  mesh_,
116  celli
117  );
118 
119  forAll(cellTets, tetI)
120  {
121  const tetIndices& cellTetIs = cellTets[tetI];
122  tetPointRef tet = cellTetIs.tet(mesh_);
123  scalar tetVolume = tet.mag();
124 
125  forAll(molecules, i)
126  {
127  const word& moleculeName(molecules[i]);
128 
129  label typeId = typeIdList_.find(moleculeName);
130 
131  if (typeId == -1)
132  {
134  << "typeId " << moleculeName << "not defined." << nl
135  << abort(FatalError);
136  }
137 
138  const typename ParcelType::constantProperties& cP =
139  constProps(typeId);
140 
141  scalar numberDensity = numberDensities[i];
142 
143  // Calculate the number of particles required
144  scalar particlesRequired = numberDensity*tetVolume;
145 
146  // Only integer numbers of particles can be inserted
147  label nParticlesToInsert = label(particlesRequired);
148 
149  // Add another particle with a probability proportional to the
150  // remainder of taking the integer part of particlesRequired
151  if
152  (
153  (particlesRequired - nParticlesToInsert)
154  > rndGen_.sample01<scalar>()
155  )
156  {
157  nParticlesToInsert++;
158  }
159 
160  for (label pI = 0; pI < nParticlesToInsert; pI++)
161  {
162  point p = tet.randomPoint(rndGen_);
163 
164  vector U = equipartitionLinearVelocity
165  (
166  temperature,
167  cP.mass()
168  );
169 
170  scalar Ei = equipartitionInternalEnergy
171  (
172  temperature,
173  cP.internalDegreesOfFreedom()
174  );
175 
176  U += velocity;
177 
178  addNewParcel(p, celli, U, Ei, typeId);
179  }
180  }
181  }
182  }
183 
184  // Initialise the sigmaTcRMax_ field to the product of the cross section of
185  // the most abundant species and the most probable thermal speed (Bird,
186  // p222-223)
187 
188  label mostAbundantType(findMax(numberDensities));
189 
190  const typename ParcelType::constantProperties& cP = constProps
191  (
192  mostAbundantType
193  );
194 
195  sigmaTcRMax_.primitiveFieldRef() = cP.sigmaT()*maxwellianMostProbableSpeed
196  (
197  temperature,
198  cP.mass()
199  );
200 
201  sigmaTcRMax_.correctBoundaryConditions();
202 }
203 
204 
205 template<class ParcelType>
207 {
208  if (!binaryCollision().active())
209  {
210  return;
211  }
212 
213  // Temporary storage for subCells
214  List<DynamicList<label>> subCells(8);
215 
216  scalar deltaT = mesh().time().deltaTValue();
217 
218  label collisionCandidates = 0;
219 
220  label collisions = 0;
221 
222  forAll(cellOccupancy_, celli)
223  {
224  const DynamicList<ParcelType*>& cellParcels(cellOccupancy_[celli]);
225 
226  label nC(cellParcels.size());
227 
228  if (nC > 1)
229  {
230  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
231  // Assign particles to one of 8 Cartesian subCells
232 
233  // Clear temporary lists
234  forAll(subCells, i)
235  {
236  subCells[i].clear();
237  }
238 
239  // Inverse addressing specifying which subCell a parcel is in
240  List<label> whichSubCell(cellParcels.size());
241 
242  const point& cC = mesh_.cellCentres()[celli];
243 
244  forAll(cellParcels, i)
245  {
246  const ParcelType& p = *cellParcels[i];
247  vector relPos = p.position() - cC;
248 
249  label subCell =
250  pos0(relPos.x()) + 2*pos0(relPos.y()) + 4*pos0(relPos.z());
251 
252  subCells[subCell].append(i);
253  whichSubCell[i] = subCell;
254  }
255 
256  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
257 
258  scalar sigmaTcRMax = sigmaTcRMax_[celli];
259 
260  scalar selectedPairs =
261  collisionSelectionRemainder_[celli]
262  + 0.5*nC*(nC - 1)*nParticle_*sigmaTcRMax*deltaT
263  /mesh_.cellVolumes()[celli];
264 
265  label nCandidates(selectedPairs);
266  collisionSelectionRemainder_[celli] = selectedPairs - nCandidates;
267  collisionCandidates += nCandidates;
268 
269  for (label c = 0; c < nCandidates; c++)
270  {
271  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
272  // subCell candidate selection procedure
273 
274  // Select the first collision candidate
275  label candidateP = rndGen_.position<label>(0, nC - 1);
276 
277  // Declare the second collision candidate
278  label candidateQ = -1;
279 
280  List<label> subCellPs = subCells[whichSubCell[candidateP]];
281  label nSC = subCellPs.size();
282 
283  if (nSC > 1)
284  {
285  // If there are two or more particle in a subCell, choose
286  // another from the same cell. If the same candidate is
287  // chosen, choose again.
288 
289  do
290  {
291  label i = rndGen_.position<label>(0, nSC - 1);
292  candidateQ = subCellPs[i];
293  } while (candidateP == candidateQ);
294  }
295  else
296  {
297  // Select a possible second collision candidate from the
298  // whole cell. If the same candidate is chosen, choose
299  // again.
300 
301  do
302  {
303  candidateQ = rndGen_.position<label>(0, nC - 1);
304  } while (candidateP == candidateQ);
305  }
306 
307  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
308  // uniform candidate selection procedure
309 
310  // // Select the first collision candidate
311  // label candidateP = rndGen_.position<label>(0, nC-1);
312 
313  // // Select a possible second collision candidate
314  // label candidateQ = rndGen_.position<label>(0, nC-1);
315 
316  // // If the same candidate is chosen, choose again
317  // while (candidateP == candidateQ)
318  // {
319  // candidateQ = rndGen_.position<label>(0, nC-1);
320  // }
321 
322  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
323 
324  ParcelType& parcelP = *cellParcels[candidateP];
325  ParcelType& parcelQ = *cellParcels[candidateQ];
326 
327  scalar sigmaTcR = binaryCollision().sigmaTcR
328  (
329  parcelP,
330  parcelQ
331  );
332 
333  // Update the maximum value of sigmaTcR stored, but use the
334  // initial value in the acceptance-rejection criteria because
335  // the number of collision candidates selected was based on this
336 
337  if (sigmaTcR > sigmaTcRMax_[celli])
338  {
339  sigmaTcRMax_[celli] = sigmaTcR;
340  }
341 
342  if ((sigmaTcR/sigmaTcRMax) > rndGen_.sample01<scalar>())
343  {
344  binaryCollision().collide
345  (
346  parcelP,
347  parcelQ
348  );
349 
350  collisions++;
351  }
352  }
353  }
354  }
355 
356  reduce(collisions, sumOp<label>());
357 
358  reduce(collisionCandidates, sumOp<label>());
359 
360  sigmaTcRMax_.correctBoundaryConditions();
361 
362  if (collisionCandidates)
363  {
364  Info<< " Collisions = "
365  << collisions << nl
366  << " Acceptance rate = "
367  << scalar(collisions)/scalar(collisionCandidates) << nl
368  << endl;
369  }
370  else
371  {
372  Info<< " No collisions" << endl;
373  }
374 }
375 
376 
377 template<class ParcelType>
379 {
380  q_ = dimensionedScalar("0", dimensionSet(1, 0, -3, 0, 0), Zero);
381  fD_ = dimensionedVector("0", dimensionSet(1, -1, -2, 0, 0), Zero);
382 
383  rhoN_ = dimensionedScalar("0", dimensionSet(0, -3, 0, 0, 0), VSMALL);
384  rhoM_ = dimensionedScalar("0", dimensionSet(1, -3, 0, 0, 0), VSMALL);
385  dsmcRhoN_ = dimensionedScalar("0", dimensionSet(0, -3, 0, 0, 0), Zero);
386  linearKE_ = dimensionedScalar("0", dimensionSet(1, -1, -2, 0, 0), Zero);
387  internalE_ = dimensionedScalar("0", dimensionSet(1, -1, -2, 0, 0), Zero);
388  iDof_ = dimensionedScalar("0", dimensionSet(0, -3, 0, 0, 0), VSMALL);
389  momentum_ = dimensionedVector("0", dimensionSet(1, -2, -1, 0, 0), Zero);
390 }
391 
392 
393 template<class ParcelType>
395 {
396  scalarField& rhoN = rhoN_.primitiveFieldRef();
397  scalarField& rhoM = rhoM_.primitiveFieldRef();
398  scalarField& dsmcRhoN = dsmcRhoN_.primitiveFieldRef();
399  scalarField& linearKE = linearKE_.primitiveFieldRef();
400  scalarField& internalE = internalE_.primitiveFieldRef();
401  scalarField& iDof = iDof_.primitiveFieldRef();
402  vectorField& momentum = momentum_.primitiveFieldRef();
403 
404  for (const ParcelType& p : *this)
405  {
406  const label celli = p.cell();
407 
408  rhoN[celli]++;
409  rhoM[celli] += constProps(p.typeId()).mass();
410  dsmcRhoN[celli]++;
411  linearKE[celli] += 0.5*constProps(p.typeId()).mass()*(p.U() & p.U());
412  internalE[celli] += p.Ei();
413  iDof[celli] += constProps(p.typeId()).internalDegreesOfFreedom();
414  momentum[celli] += constProps(p.typeId()).mass()*p.U();
415  }
416 
417  rhoN *= nParticle_/mesh().cellVolumes();
418  rhoN_.correctBoundaryConditions();
419 
420  rhoM *= nParticle_/mesh().cellVolumes();
421  rhoM_.correctBoundaryConditions();
422 
423  dsmcRhoN_.correctBoundaryConditions();
424 
425  linearKE *= nParticle_/mesh().cellVolumes();
426  linearKE_.correctBoundaryConditions();
427 
428  internalE *= nParticle_/mesh().cellVolumes();
429  internalE_.correctBoundaryConditions();
430 
431  iDof *= nParticle_/mesh().cellVolumes();
432  iDof_.correctBoundaryConditions();
433 
434  momentum *= nParticle_/mesh().cellVolumes();
435  momentum_.correctBoundaryConditions();
436 }
437 
438 
439 // * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * * //
440 
441 template<class ParcelType>
443 (
444  const vector& position,
445  const label celli,
446  const vector& U,
447  const scalar Ei,
448  const label typeId
449 )
450 {
451  this->addParticle(new ParcelType(mesh_, position, celli, U, Ei, typeId));
452 }
453 
454 
455 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
456 
457 template<class ParcelType>
459 (
460  const word& cloudName,
461  const fvMesh& mesh,
462  bool readFields
463 )
464 :
466  DSMCBaseCloud(),
467  cloudName_(cloudName),
468  mesh_(mesh),
469  particleProperties_
470  (
471  IOobject
472  (
473  cloudName + "Properties",
474  mesh_.time().constant(),
475  mesh_,
478  )
479  ),
480  typeIdList_(particleProperties_.lookup("typeIdList")),
481  nParticle_(particleProperties_.get<scalar>("nEquivalentParticles")),
482  cellOccupancy_(mesh_.nCells()),
483  sigmaTcRMax_
484  (
485  IOobject
486  (
487  this->name() + "SigmaTcRMax",
488  mesh_.time().timeName(),
489  mesh_,
492  ),
493  mesh_
494  ),
495  collisionSelectionRemainder_
496  (
497  IOobject
498  (
499  this->name() + ":collisionSelectionRemainder",
500  mesh_.time().timeName(),
501  mesh_
502  ),
503  mesh_,
505  ),
506  q_
507  (
508  IOobject
509  (
510  "q",
511  mesh_.time().timeName(),
512  mesh_,
515  ),
516  mesh_
517  ),
518  fD_
519  (
520  IOobject
521  (
522  "fD",
523  mesh_.time().timeName(),
524  mesh_,
527  ),
528  mesh_
529  ),
530  rhoN_
531  (
532  IOobject
533  (
534  "rhoN",
535  mesh_.time().timeName(),
536  mesh_,
539  ),
540  mesh_
541  ),
542  rhoM_
543  (
544  IOobject
545  (
546  "rhoM",
547  mesh_.time().timeName(),
548  mesh_,
551  ),
552  mesh_
553  ),
554  dsmcRhoN_
555  (
556  IOobject
557  (
558  "dsmcRhoN",
559  mesh_.time().timeName(),
560  mesh_,
563  ),
564  mesh_
565  ),
566  linearKE_
567  (
568  IOobject
569  (
570  "linearKE",
571  mesh_.time().timeName(),
572  mesh_,
575  ),
576  mesh_
577  ),
578  internalE_
579  (
580  IOobject
581  (
582  "internalE",
583  mesh_.time().timeName(),
584  mesh_,
587  ),
588  mesh_
589  ),
590  iDof_
591  (
592  IOobject
593  (
594  "iDof",
595  mesh_.time().timeName(),
596  mesh_,
599  ),
600  mesh_
601  ),
602  momentum_
603  (
604  IOobject
605  (
606  "momentum",
607  mesh_.time().timeName(),
608  mesh_,
611  ),
612  mesh_
613  ),
614  constProps_(),
615  rndGen_(Pstream::myProcNo()),
616  boundaryT_
617  (
619  (
620  IOobject
621  (
622  "boundaryT",
623  mesh_.time().timeName(),
624  mesh_,
627  ),
628  mesh_
629  )
630  ),
631  boundaryU_
632  (
634  (
635  IOobject
636  (
637  "boundaryU",
638  mesh_.time().timeName(),
639  mesh_,
642  ),
643  mesh_
644  )
645  ),
646  binaryCollisionModel_
647  (
649  (
650  particleProperties_,
651  *this
652  )
653  ),
654  wallInteractionModel_
655  (
657  (
658  particleProperties_,
659  *this
660  )
661  ),
662  inflowBoundaryModel_
663  (
665  (
666  particleProperties_,
667  *this
668  )
669  )
670 {
671  buildConstProps();
672  buildCellOccupancy();
673 
674  // Initialise the collision selection remainder to a random value between 0
675  // and 1.
676  forAll(collisionSelectionRemainder_, i)
677  {
678  collisionSelectionRemainder_[i] = rndGen_.sample01<scalar>();
679  }
680 
681  if (readFields)
682  {
683  ParcelType::readFields(*this);
684  }
685 }
686 
687 
688 template<class ParcelType>
690 (
691  const word& cloudName,
692  const fvMesh& mesh,
693  const IOdictionary& dsmcInitialiseDict
694 )
695 :
697  DSMCBaseCloud(),
698  cloudName_(cloudName),
699  mesh_(mesh),
700  particleProperties_
701  (
702  IOobject
703  (
704  cloudName + "Properties",
705  mesh_.time().constant(),
706  mesh_,
709  )
710  ),
711  typeIdList_(particleProperties_.lookup("typeIdList")),
712  nParticle_(particleProperties_.get<scalar>("nEquivalentParticles")),
713  cellOccupancy_(),
714  sigmaTcRMax_
715  (
716  IOobject
717  (
718  this->name() + "SigmaTcRMax",
719  mesh_.time().timeName(),
720  mesh_,
723  ),
724  mesh_,
725  dimensionedScalar(dimensionSet(0, 3, -1, 0, 0), Zero),
726  zeroGradientFvPatchScalarField::typeName
727  ),
728  collisionSelectionRemainder_
729  (
730  IOobject
731  (
732  this->name() + ":collisionSelectionRemainder",
733  mesh_.time().timeName(),
734  mesh_
735  ),
736  mesh_,
738  ),
739  q_
740  (
741  IOobject
742  (
743  this->name() + "q_",
744  mesh_.time().timeName(),
745  mesh_,
748  ),
749  mesh_,
750  dimensionedScalar(dimensionSet(1, 0, -3, 0, 0), Zero)
751  ),
752  fD_
753  (
754  IOobject
755  (
756  this->name() + "fD_",
757  mesh_.time().timeName(),
758  mesh_,
761  ),
762  mesh_,
763  dimensionedVector(dimensionSet(1, -1, -2, 0, 0), Zero)
764  ),
765  rhoN_
766  (
767  IOobject
768  (
769  this->name() + "rhoN_",
770  mesh_.time().timeName(),
771  mesh_,
774  ),
775  mesh_,
776  dimensionedScalar("0", dimensionSet(0, -3, 0, 0, 0), VSMALL)
777  ),
778  rhoM_
779  (
780  IOobject
781  (
782  this->name() + "rhoM_",
783  mesh_.time().timeName(),
784  mesh_,
787  ),
788  mesh_,
789  dimensionedScalar("0", dimensionSet(1, -3, 0, 0, 0), VSMALL)
790  ),
791  dsmcRhoN_
792  (
793  IOobject
794  (
795  this->name() + "dsmcRhoN_",
796  mesh_.time().timeName(),
797  mesh_,
800  ),
801  mesh_,
802  dimensionedScalar(dimensionSet(0, -3, 0, 0, 0), Zero)
803  ),
804  linearKE_
805  (
806  IOobject
807  (
808  this->name() + "linearKE_",
809  mesh_.time().timeName(),
810  mesh_,
813  ),
814  mesh_,
815  dimensionedScalar(dimensionSet(1, -1, -2, 0, 0), Zero)
816  ),
817  internalE_
818  (
819  IOobject
820  (
821  this->name() + "internalE_",
822  mesh_.time().timeName(),
823  mesh_,
826  ),
827  mesh_,
828  dimensionedScalar(dimensionSet(1, -1, -2, 0, 0), Zero)
829  ),
830  iDof_
831  (
832  IOobject
833  (
834  this->name() + "iDof_",
835  mesh_.time().timeName(),
836  mesh_,
839  ),
840  mesh_,
841  dimensionedScalar("0", dimensionSet(0, -3, 0, 0, 0), VSMALL)
842  ),
843  momentum_
844  (
845  IOobject
846  (
847  this->name() + "momentum_",
848  mesh_.time().timeName(),
849  mesh_,
852  ),
853  mesh_,
854  dimensionedVector(dimensionSet(1, -2, -1, 0, 0), Zero)
855  ),
856  constProps_(),
857  rndGen_(Pstream::myProcNo()),
858  boundaryT_
859  (
861  (
862  IOobject
863  (
864  "boundaryT",
865  mesh_.time().timeName(),
866  mesh_,
869  ),
870  mesh_,
871  dimensionedScalar(dimensionSet(0, 0, 0, 1, 0), Zero)
872  )
873  ),
874  boundaryU_
875  (
877  (
878  IOobject
879  (
880  "boundaryU",
881  mesh_.time().timeName(),
882  mesh_,
885  ),
886  mesh_,
887  dimensionedVector(dimensionSet(0, 1, -1, 0, 0), Zero)
888  )
889  ),
890  binaryCollisionModel_(),
891  wallInteractionModel_(),
892  inflowBoundaryModel_()
893 {
894  clear();
895  buildConstProps();
896  initialise(dsmcInitialiseDict);
897 }
898 
899 
900 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
901 
902 template<class ParcelType>
904 {}
905 
906 
907 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
908 
909 template<class ParcelType>
911 {
912  typename ParcelType::trackingData td(*this);
913 
914  // Reset the data collection fields
915  resetFields();
916 
917  if (debug)
918  {
919  this->dumpParticlePositions();
920  }
921 
922  // Insert new particles from the inflow boundary
923  this->inflowBoundary().inflow();
924 
925  // Move the particles ballistically with their current velocities
926  Cloud<ParcelType>::move(*this, td, mesh_.time().deltaTValue());
927 
928  // Update cell occupancy
929  buildCellOccupancy();
930 
931  // Calculate new velocities via stochastic collisions
932  collisions();
933 
934  // Calculate the volume field data
935  calculateFields();
936 }
937 
938 
939 template<class ParcelType>
941 {
942  label nDSMCParticles = this->size();
943  reduce(nDSMCParticles, sumOp<label>());
944 
945  scalar nMol = nDSMCParticles*nParticle_;
946 
947  vector linearMomentum = linearMomentumOfSystem();
948  reduce(linearMomentum, sumOp<vector>());
949 
950  scalar linearKineticEnergy = linearKineticEnergyOfSystem();
951  reduce(linearKineticEnergy, sumOp<scalar>());
952 
953  scalar internalEnergy = internalEnergyOfSystem();
954  reduce(internalEnergy, sumOp<scalar>());
955 
956  Info<< "Cloud name: " << this->name() << nl
957  << " Number of dsmc particles = "
958  << nDSMCParticles
959  << endl;
960 
961  if (nDSMCParticles)
962  {
963  Info<< " Number of molecules = "
964  << nMol << nl
965  << " Mass in system = "
966  << returnReduce(massInSystem(), sumOp<scalar>()) << nl
967  << " Average linear momentum = "
968  << linearMomentum/nMol << nl
969  << " |Average linear momentum| = "
970  << mag(linearMomentum)/nMol << nl
971  << " Average linear kinetic energy = "
972  << linearKineticEnergy/nMol << nl
973  << " Average internal energy = "
974  << internalEnergy/nMol << nl
975  << " Average total energy = "
976  << (internalEnergy + linearKineticEnergy)/nMol
977  << endl;
978  }
979 }
980 
981 
982 template<class ParcelType>
984 (
985  scalar temperature,
986  scalar mass
987 )
988 {
989  return
990  sqrt(physicoChemical::k.value()*temperature/mass)
991  *rndGen_.GaussNormal<vector>();
992 }
993 
994 
995 template<class ParcelType>
997 (
998  scalar temperature,
999  direction iDof
1000 )
1001 {
1002  if (iDof == 0)
1003  {
1004  return 0;
1005  }
1006  else if (iDof == 2)
1007  {
1008  // Special case for iDof = 2, i.e. diatomics;
1009  return
1010  (
1011  -log(rndGen_.sample01<scalar>())
1012  *physicoChemical::k.value()*temperature
1013  );
1014  }
1015 
1016 
1017  const scalar a = 0.5*iDof - 1;
1018  scalar energyRatio = 0;
1019  scalar P = -1;
1020 
1021  do
1022  {
1023  energyRatio = 10*rndGen_.sample01<scalar>();
1024  P = pow((energyRatio/a), a)*exp(a - energyRatio);
1025  } while (P < rndGen_.sample01<scalar>());
1026 
1027  return energyRatio*physicoChemical::k.value()*temperature;
1028 }
1029 
1030 
1031 template<class ParcelType>
1033 {
1034  OFstream pObj
1035  (
1036  this->db().time().path()/"parcelPositions_"
1037  + this->name() + "_"
1038  + this->db().time().timeName() + ".obj"
1039  );
1040 
1041  for (const ParcelType& p : *this)
1042  {
1043  pObj<< "v " << p.position().x()
1044  << ' ' << p.position().y()
1045  << ' ' << p.position().z()
1046  << nl;
1047  }
1048 
1049  pObj.flush();
1050 }
1051 
1052 
1053 template<class ParcelType>
1055 {
1057 
1058  // Update the cell occupancy field
1059  cellOccupancy_.setSize(mesh_.nCells());
1060  buildCellOccupancy();
1061 
1062  // Update the inflow BCs
1063  this->inflowBoundary().autoMap(mapper);
1064 }
1065 
1066 
1067 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:54
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:194
Foam::DSMCCloud::dumpParticlePositions
void dumpParticlePositions() const
Dump particle positions to .obj file.
Definition: DSMCCloud.C:1032
cloudName
const word cloudName(propsDict.get< word >("cloud"))
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:94
Foam::constant::physicoChemical::k
const dimensionedScalar k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::BinaryCollisionModel
Templated DSMC particle collision class.
Definition: DSMCCloud.H:58
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
BinaryCollisionModel.H
Foam::constant
Different types of constants.
Definition: atomicConstants.C:38
Foam::DSMCCloud::equipartitionLinearVelocity
vector equipartitionLinearVelocity(scalar temperature, scalar mass)
Generate a random velocity sampled from the Maxwellian speed.
Definition: DSMCCloud.C:984
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:434
Foam::pos0
dimensionedScalar pos0(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:188
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:261
Foam::dimensionSet
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:108
polyMeshTetDecomposition.H
Foam::sumOp
Definition: ops.H:213
Foam::tetPointRef
tetrahedron< point, const point & > tetPointRef
A tetrahedron using referred points.
Definition: tetPointRef.H:47
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
Foam::dimensionedVector
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Definition: dimensionedVector.H:50
Foam::DSMCCloud::clear
void clear()
Clear the Cloud.
Definition: DSMCCloudI.H:475
Foam::DSMCCloud::info
void info() const
Print cloud information.
Definition: DSMCCloud.C:940
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:51
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::Cloud::autoMap
void autoMap(const mapPolyMesh &)
Remap the cells of particles corresponding to the.
Definition: Cloud.C:345
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
InflowBoundaryModel.H
Foam::DSMCCloud::~DSMCCloud
virtual ~DSMCCloud()
Destructor.
Definition: DSMCCloud.C:903
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
Foam::DSMCCloud
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:71
Foam::WallInteractionModel
Templated wall interaction model class.
Definition: DSMCCloud.H:61
timeName
word timeName
Definition: getTimeIndex.H:3
Foam::DSMCCloud::addNewParcel
void addNewParcel(const vector &position, const label celli, const vector &U, const scalar Ei, const label typeId)
Add new parcel.
Definition: DSMCCloud.C:443
Foam::Cloud::move
void move(TrackCloudType &cloud, typename ParticleType::trackingData &td, const scalar trackTime)
Move the particles.
Definition: Cloud.C:147
Foam::FatalError
error FatalError
Foam::log
dimensionedScalar log(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:262
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
WallInteractionModel.H
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
constants.H
Foam::OFstream
Output to file stream, using an OSstream.
Definition: OFstream.H:53
Foam::readFields
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type.
Definition: ReadFieldsTemplates.C:312
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
U
U
Definition: pEqn.H:72
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=worldComm)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:463
clear
patchWriters clear()
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::Vector< scalar >
Foam::findMax
label findMax(const ListType &input, label start=0)
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::Cloud< ParcelType >
DSMCCloud.H
path
fileName path(UMean.rootPath()/UMean.caseName()/"graphs"/UMean.instance())
Foam::direction
uint8_t direction
Definition: direction.H:52
Foam::DSMCCloud::autoMap
virtual void autoMap(const mapPolyMesh &)
Remap the particles to the correct cells following mesh change.
Definition: DSMCCloud.C:1054
Foam::IOobject::MUST_READ_IF_MODIFIED
Definition: IOobject.H:186
Foam::polyMeshTetDecomposition::cellTetIndices
static List< tetIndices > cellTetIndices(const polyMesh &mesh, label cI)
Return the tet decomposition of the given cell, see.
Definition: polyMeshTetDecomposition.C:566
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:161
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::DSMCCloud::equipartitionInternalEnergy
scalar equipartitionInternalEnergy(scalar temperature, direction internalDegreesOfFreedom)
Generate a random internal energy, sampled from the.
Definition: DSMCCloud.C:997
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::DSMCBaseCloud
Virtual abstract base class for templated DSMCCloud.
Definition: DSMCBaseCloud.H:50
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:188
Foam::InflowBoundaryModel
Templated inflow boundary model class.
Definition: DSMCCloud.H:64
zeroGradientFvPatchFields.H
Foam::DSMCCloud::evolve
void evolve()
Evolve the cloud (move, collide)
Definition: DSMCCloud.C:910
Foam::OSstream::flush
virtual void flush()
Flush stream.
Definition: OSstream.C:275
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189
Foam::IOobject::MUST_READ
Definition: IOobject.H:185