DSMCCloud.C
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 "DSMCCloud.H"
27 #include "BinaryCollisionModel.H"
28 #include "WallInteractionModel.H"
29 #include "InflowBoundaryModel.H"
30 #include "constants.H"
33 
34 using namespace Foam::constant;
35 
36 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
37 
38 template<class ParcelType>
40 {
41  Info<< nl << "Constructing constant properties for" << endl;
42  constProps_.setSize(typeIdList_.size());
43 
44  dictionary moleculeProperties
45  (
46  particleProperties_.subDict("moleculeProperties")
47  );
48 
49  forAll(typeIdList_, i)
50  {
51  const word& id(typeIdList_[i]);
52 
53  Info<< " " << id << endl;
54 
55  const dictionary& molDict(moleculeProperties.subDict(id));
56 
57  constProps_[i] =
58  typename ParcelType::constantProperties(molDict);
59  }
60 }
61 
62 
63 template<class ParcelType>
65 {
66  forAll(cellOccupancy_, cO)
67  {
68  cellOccupancy_[cO].clear();
69  }
70 
71  forAllIter(typename DSMCCloud<ParcelType>, *this, iter)
72  {
73  cellOccupancy_[iter().cell()].append(&iter());
74  }
75 }
76 
77 
78 template<class ParcelType>
80 (
81  const IOdictionary& dsmcInitialiseDict
82 )
83 {
84  Info<< nl << "Initialising particles" << endl;
85 
86  const scalar temperature
87  (
88  readScalar(dsmcInitialiseDict.lookup("temperature"))
89  );
90 
91  const vector velocity(dsmcInitialiseDict.lookup("velocity"));
92 
93  const dictionary& numberDensitiesDict
94  (
95  dsmcInitialiseDict.subDict("numberDensities")
96  );
97 
98  List<word> molecules(numberDensitiesDict.toc());
99 
100  Field<scalar> numberDensities(molecules.size());
101 
102  forAll(molecules, i)
103  {
104  numberDensities[i] = readScalar
105  (
106  numberDensitiesDict.lookup(molecules[i])
107  );
108  }
109 
110  numberDensities /= nParticle_;
111 
112  forAll(mesh_.cells(), celli)
113  {
114  List<tetIndices> cellTets = polyMeshTetDecomposition::cellTetIndices
115  (
116  mesh_,
117  celli
118  );
119 
120  forAll(cellTets, tetI)
121  {
122  const tetIndices& cellTetIs = cellTets[tetI];
123  tetPointRef tet = cellTetIs.tet(mesh_);
124  scalar tetVolume = tet.mag();
125 
126  forAll(molecules, i)
127  {
128  const word& moleculeName(molecules[i]);
129 
130  label typeId(findIndex(typeIdList_, moleculeName));
131 
132  if (typeId == -1)
133  {
135  << "typeId " << moleculeName << "not defined." << nl
136  << abort(FatalError);
137  }
138 
139  const typename ParcelType::constantProperties& cP =
140  constProps(typeId);
141 
142  scalar numberDensity = numberDensities[i];
143 
144  // Calculate the number of particles required
145  scalar particlesRequired = numberDensity*tetVolume;
146 
147  // Only integer numbers of particles can be inserted
148  label nParticlesToInsert = label(particlesRequired);
149 
150  // Add another particle with a probability proportional to the
151  // remainder of taking the integer part of particlesRequired
152  if
153  (
154  (particlesRequired - nParticlesToInsert)
155  > rndGen_.scalar01()
156  )
157  {
158  nParticlesToInsert++;
159  }
160 
161  for (label pI = 0; pI < nParticlesToInsert; pI++)
162  {
163  point p = tet.randomPoint(rndGen_);
164 
165  vector U = equipartitionLinearVelocity
166  (
167  temperature,
168  cP.mass()
169  );
170 
171  scalar Ei = equipartitionInternalEnergy
172  (
173  temperature,
174  cP.internalDegreesOfFreedom()
175  );
176 
177  U += velocity;
178 
179  addNewParcel(p, celli, U, Ei, typeId);
180  }
181  }
182  }
183  }
184 
185  // Initialise the sigmaTcRMax_ field to the product of the cross section of
186  // the most abundant species and the most probable thermal speed (Bird,
187  // p222-223)
188 
189  label mostAbundantType(findMax(numberDensities));
190 
191  const typename ParcelType::constantProperties& cP = constProps
192  (
193  mostAbundantType
194  );
195 
196  sigmaTcRMax_.primitiveFieldRef() = cP.sigmaT()*maxwellianMostProbableSpeed
197  (
198  temperature,
199  cP.mass()
200  );
201 
202  sigmaTcRMax_.correctBoundaryConditions();
203 }
204 
205 
206 template<class ParcelType>
208 {
209  if (!binaryCollision().active())
210  {
211  return;
212  }
213 
214  // Temporary storage for subCells
215  List<DynamicList<label>> subCells(8);
216 
217  scalar deltaT = mesh().time().deltaTValue();
218 
219  label collisionCandidates = 0;
220 
221  label collisions = 0;
222 
223  forAll(cellOccupancy_, celli)
224  {
225  const DynamicList<ParcelType*>& cellParcels(cellOccupancy_[celli]);
226 
227  label nC(cellParcels.size());
228 
229  if (nC > 1)
230  {
231  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
232  // Assign particles to one of 8 Cartesian subCells
233 
234  // Clear temporary lists
235  forAll(subCells, i)
236  {
237  subCells[i].clear();
238  }
239 
240  // Inverse addressing specifying which subCell a parcel is in
241  List<label> whichSubCell(cellParcels.size());
242 
243  const point& cC = mesh_.cellCentres()[celli];
244 
245  forAll(cellParcels, i)
246  {
247  const ParcelType& p = *cellParcels[i];
248  vector relPos = p.position() - cC;
249 
250  label subCell =
251  pos0(relPos.x()) + 2*pos0(relPos.y()) + 4*pos0(relPos.z());
252 
253  subCells[subCell].append(i);
254  whichSubCell[i] = subCell;
255  }
256 
257  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
258 
259  scalar sigmaTcRMax = sigmaTcRMax_[celli];
260 
261  scalar selectedPairs =
262  collisionSelectionRemainder_[celli]
263  + 0.5*nC*(nC - 1)*nParticle_*sigmaTcRMax*deltaT
264  /mesh_.cellVolumes()[celli];
265 
266  label nCandidates(selectedPairs);
267  collisionSelectionRemainder_[celli] = selectedPairs - nCandidates;
268  collisionCandidates += nCandidates;
269 
270  for (label c = 0; c < nCandidates; c++)
271  {
272  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
273  // subCell candidate selection procedure
274 
275  // Select the first collision candidate
276  label candidateP = rndGen_.sampleAB<label>(0, nC);
277 
278  // Declare the second collision candidate
279  label candidateQ = -1;
280 
281  List<label> subCellPs = subCells[whichSubCell[candidateP]];
282  label nSC = subCellPs.size();
283 
284  if (nSC > 1)
285  {
286  // If there are two or more particle in a subCell, choose
287  // another from the same cell. If the same candidate is
288  // chosen, choose again.
289 
290  do
291  {
292  candidateQ = subCellPs[rndGen_.sampleAB<label>(0, nSC)];
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_.sampleAB<label>(0, nC);
304  } while (candidateP == candidateQ);
305  }
306 
307  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
308  // uniform candidate selection procedure
309 
310  // // Select the first collision candidate
311  // label candidateP = rndGen_.sampleAB<label>(0, nC);
312 
313  // // Select a possible second collision candidate
314  // label candidateQ = rndGen_.sampleAB<label>(0, nC);
315 
316  // // If the same candidate is chosen, choose again
317  // while (candidateP == candidateQ)
318  // {
319  // candidateQ = rndGen_.sampleAB<label>(0, nC);
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_.scalar01())
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("zero", dimensionSet(1, 0, -3, 0, 0), 0.0);
381 
382  fD_ = dimensionedVector
383  (
384  "zero",
385  dimensionSet(1, -1, -2, 0, 0),
386  Zero
387  );
388 
389  rhoN_ = dimensionedScalar("zero", dimensionSet(0, -3, 0, 0, 0), vSmall);
390  rhoM_ = dimensionedScalar("zero", dimensionSet(1, -3, 0, 0, 0), vSmall);
391  dsmcRhoN_ = dimensionedScalar("zero", dimensionSet(0, -3, 0, 0, 0), 0.0);
392  linearKE_ = dimensionedScalar("zero", dimensionSet(1, -1, -2, 0, 0), 0.0);
393  internalE_ = dimensionedScalar("zero", dimensionSet(1, -1, -2, 0, 0), 0.0);
394  iDof_ = dimensionedScalar("zero", dimensionSet(0, -3, 0, 0, 0), vSmall);
395 
396  momentum_ = dimensionedVector
397  (
398  "zero",
399  dimensionSet(1, -2, -1, 0, 0),
400  Zero
401  );
402 }
403 
404 
405 template<class ParcelType>
407 {
408  scalarField& rhoN = rhoN_.primitiveFieldRef();
409  scalarField& rhoM = rhoM_.primitiveFieldRef();
410  scalarField& dsmcRhoN = dsmcRhoN_.primitiveFieldRef();
411  scalarField& linearKE = linearKE_.primitiveFieldRef();
412  scalarField& internalE = internalE_.primitiveFieldRef();
413  scalarField& iDof = iDof_.primitiveFieldRef();
414  vectorField& momentum = momentum_.primitiveFieldRef();
415 
416  forAllConstIter(typename DSMCCloud<ParcelType>, *this, iter)
417  {
418  const ParcelType& p = iter();
419  const label celli = p.cell();
420 
421  rhoN[celli]++;
422  rhoM[celli] += constProps(p.typeId()).mass();
423  dsmcRhoN[celli]++;
424  linearKE[celli] += 0.5*constProps(p.typeId()).mass()*(p.U() & p.U());
425  internalE[celli] += p.Ei();
426  iDof[celli] += constProps(p.typeId()).internalDegreesOfFreedom();
427  momentum[celli] += constProps(p.typeId()).mass()*p.U();
428  }
429 
430  rhoN *= nParticle_/mesh().cellVolumes();
431  rhoN_.correctBoundaryConditions();
432 
433  rhoM *= nParticle_/mesh().cellVolumes();
434  rhoM_.correctBoundaryConditions();
435 
436  dsmcRhoN_.correctBoundaryConditions();
437 
438  linearKE *= nParticle_/mesh().cellVolumes();
439  linearKE_.correctBoundaryConditions();
440 
441  internalE *= nParticle_/mesh().cellVolumes();
442  internalE_.correctBoundaryConditions();
443 
444  iDof *= nParticle_/mesh().cellVolumes();
445  iDof_.correctBoundaryConditions();
446 
447  momentum *= nParticle_/mesh().cellVolumes();
448  momentum_.correctBoundaryConditions();
449 }
450 
451 
452 // * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * * //
453 
454 template<class ParcelType>
456 (
457  const vector& position,
458  const label celli,
459  const vector& U,
460  const scalar Ei,
461  const label typeId
462 )
463 {
464  this->addParticle(new ParcelType(mesh_, position, celli, U, Ei, typeId));
465 }
466 
467 
468 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
469 
470 template<class ParcelType>
472 (
473  const word& cloudName,
474  const fvMesh& mesh,
475  bool readFields
476 )
477 :
479  DSMCBaseCloud(),
480  cloudName_(cloudName),
481  mesh_(mesh),
482  particleProperties_
483  (
484  IOobject
485  (
486  cloudName + "Properties",
487  mesh_.time().constant(),
488  mesh_,
491  )
492  ),
493  typeIdList_(particleProperties_.lookup("typeIdList")),
494  nParticle_(readScalar(particleProperties_.lookup("nEquivalentParticles"))),
495  cellOccupancy_(mesh_.nCells()),
496  sigmaTcRMax_
497  (
498  IOobject
499  (
500  this->name() + "SigmaTcRMax",
501  mesh_.time().timeName(),
502  mesh_,
505  ),
506  mesh_
507  ),
508  collisionSelectionRemainder_
509  (
510  IOobject
511  (
512  this->name() + ":collisionSelectionRemainder",
513  mesh_.time().timeName(),
514  mesh_
515  ),
516  mesh_,
517  dimensionedScalar("collisionSelectionRemainder", dimless, 0)
518  ),
519  q_
520  (
521  IOobject
522  (
523  "q",
524  mesh_.time().timeName(),
525  mesh_,
528  ),
529  mesh_
530  ),
531  fD_
532  (
533  IOobject
534  (
535  "fD",
536  mesh_.time().timeName(),
537  mesh_,
540  ),
541  mesh_
542  ),
543  rhoN_
544  (
545  IOobject
546  (
547  "rhoN",
548  mesh_.time().timeName(),
549  mesh_,
552  ),
553  mesh_
554  ),
555  rhoM_
556  (
557  IOobject
558  (
559  "rhoM",
560  mesh_.time().timeName(),
561  mesh_,
564  ),
565  mesh_
566  ),
567  dsmcRhoN_
568  (
569  IOobject
570  (
571  "dsmcRhoN",
572  mesh_.time().timeName(),
573  mesh_,
576  ),
577  mesh_
578  ),
579  linearKE_
580  (
581  IOobject
582  (
583  "linearKE",
584  mesh_.time().timeName(),
585  mesh_,
588  ),
589  mesh_
590  ),
591  internalE_
592  (
593  IOobject
594  (
595  "internalE",
596  mesh_.time().timeName(),
597  mesh_,
600  ),
601  mesh_
602  ),
603  iDof_
604  (
605  IOobject
606  (
607  "iDof",
608  mesh_.time().timeName(),
609  mesh_,
612  ),
613  mesh_
614  ),
615  momentum_
616  (
617  IOobject
618  (
619  "momentum",
620  mesh_.time().timeName(),
621  mesh_,
624  ),
625  mesh_
626  ),
627  constProps_(),
628  rndGen_(label(149382906) + 7183*Pstream::myProcNo()),
629  boundaryT_
630  (
632  (
633  IOobject
634  (
635  "boundaryT",
636  mesh_.time().timeName(),
637  mesh_,
640  ),
641  mesh_
642  )
643  ),
644  boundaryU_
645  (
647  (
648  IOobject
649  (
650  "boundaryU",
651  mesh_.time().timeName(),
652  mesh_,
655  ),
656  mesh_
657  )
658  ),
659  binaryCollisionModel_
660  (
662  (
663  particleProperties_,
664  *this
665  )
666  ),
667  wallInteractionModel_
668  (
670  (
671  particleProperties_,
672  *this
673  )
674  ),
675  inflowBoundaryModel_
676  (
678  (
679  particleProperties_,
680  *this
681  )
682  )
683 {
684  buildConstProps();
685  buildCellOccupancy();
686 
687  // Initialise the collision selection remainder to a random value between 0
688  // and 1.
689  forAll(collisionSelectionRemainder_, i)
690  {
691  collisionSelectionRemainder_[i] = rndGen_.scalar01();
692  }
693 
694  if (readFields)
695  {
696  ParcelType::readFields(*this);
697  }
698 }
699 
700 
701 template<class ParcelType>
703 (
704  const word& cloudName,
705  const fvMesh& mesh,
706  const IOdictionary& dsmcInitialiseDict
707 )
708  :
710  DSMCBaseCloud(),
711  cloudName_(cloudName),
712  mesh_(mesh),
713  particleProperties_
714  (
715  IOobject
716  (
717  cloudName + "Properties",
718  mesh_.time().constant(),
719  mesh_,
722  )
723  ),
724  typeIdList_(particleProperties_.lookup("typeIdList")),
725  nParticle_(readScalar(particleProperties_.lookup("nEquivalentParticles"))),
726  cellOccupancy_(),
727  sigmaTcRMax_
728  (
729  IOobject
730  (
731  this->name() + "SigmaTcRMax",
732  mesh_.time().timeName(),
733  mesh_,
736  ),
737  mesh_,
738  dimensionedScalar("zero", dimensionSet(0, 3, -1, 0, 0), 0.0),
739  zeroGradientFvPatchScalarField::typeName
740  ),
741  collisionSelectionRemainder_
742  (
743  IOobject
744  (
745  this->name() + ":collisionSelectionRemainder",
746  mesh_.time().timeName(),
747  mesh_
748  ),
749  mesh_,
750  dimensionedScalar("collisionSelectionRemainder", dimless, 0)
751  ),
752  q_
753  (
754  IOobject
755  (
756  this->name() + "q_",
757  mesh_.time().timeName(),
758  mesh_,
761  ),
762  mesh_,
763  dimensionedScalar("zero", dimensionSet(1, 0, -3, 0, 0), 0.0)
764  ),
765  fD_
766  (
767  IOobject
768  (
769  this->name() + "fD_",
770  mesh_.time().timeName(),
771  mesh_,
774  ),
775  mesh_,
777  (
778  "zero",
779  dimensionSet(1, -1, -2, 0, 0),
780  Zero
781  )
782  ),
783  rhoN_
784  (
785  IOobject
786  (
787  this->name() + "rhoN_",
788  mesh_.time().timeName(),
789  mesh_,
792  ),
793  mesh_,
794  dimensionedScalar("zero", dimensionSet(0, -3, 0, 0, 0), vSmall)
795  ),
796  rhoM_
797  (
798  IOobject
799  (
800  this->name() + "rhoM_",
801  mesh_.time().timeName(),
802  mesh_,
805  ),
806  mesh_,
807  dimensionedScalar("zero", dimensionSet(1, -3, 0, 0, 0), vSmall)
808  ),
809  dsmcRhoN_
810  (
811  IOobject
812  (
813  this->name() + "dsmcRhoN_",
814  mesh_.time().timeName(),
815  mesh_,
818  ),
819  mesh_,
820  dimensionedScalar("zero", dimensionSet(0, -3, 0, 0, 0), 0.0)
821  ),
822  linearKE_
823  (
824  IOobject
825  (
826  this->name() + "linearKE_",
827  mesh_.time().timeName(),
828  mesh_,
831  ),
832  mesh_,
833  dimensionedScalar("zero", dimensionSet(1, -1, -2, 0, 0), 0.0)
834  ),
835  internalE_
836  (
837  IOobject
838  (
839  this->name() + "internalE_",
840  mesh_.time().timeName(),
841  mesh_,
844  ),
845  mesh_,
846  dimensionedScalar("zero", dimensionSet(1, -1, -2, 0, 0), 0.0)
847  ),
848  iDof_
849  (
850  IOobject
851  (
852  this->name() + "iDof_",
853  mesh_.time().timeName(),
854  mesh_,
857  ),
858  mesh_,
859  dimensionedScalar("zero", dimensionSet(0, -3, 0, 0, 0), vSmall)
860  ),
861  momentum_
862  (
863  IOobject
864  (
865  this->name() + "momentum_",
866  mesh_.time().timeName(),
867  mesh_,
870  ),
871  mesh_,
873  (
874  "zero",
875  dimensionSet(1, -2, -1, 0, 0),
876  Zero
877  )
878  ),
879  constProps_(),
880  rndGen_(label(971501) + 1526*Pstream::myProcNo()),
881  boundaryT_
882  (
884  (
885  IOobject
886  (
887  "boundaryT",
888  mesh_.time().timeName(),
889  mesh_,
892  ),
893  mesh_,
894  dimensionedScalar("zero", dimensionSet(0, 0, 0, 1, 0), 0.0)
895  )
896  ),
897  boundaryU_
898  (
900  (
901  IOobject
902  (
903  "boundaryU",
904  mesh_.time().timeName(),
905  mesh_,
908  ),
909  mesh_,
911  (
912  "zero",
913  dimensionSet(0, 1, -1, 0, 0),
914  Zero
915  )
916  )
917  ),
918  binaryCollisionModel_(),
919  wallInteractionModel_(),
920  inflowBoundaryModel_()
921 {
922  clear();
923  buildConstProps();
924  initialise(dsmcInitialiseDict);
925 }
926 
927 
928 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
929 
930 template<class ParcelType>
932 {}
933 
934 
935 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
936 
937 template<class ParcelType>
939 {
940  typename ParcelType::trackingData td(*this);
941 
942  // Reset the data collection fields
943  resetFields();
944 
945  if (debug)
946  {
947  this->dumpParticlePositions();
948  }
949 
950  // Insert new particles from the inflow boundary
951  this->inflowBoundary().inflow();
952 
953  // Move the particles ballistically with their current velocities
954  Cloud<ParcelType>::move(*this, td, mesh_.time().deltaTValue());
955 
956  // Update cell occupancy
957  buildCellOccupancy();
958 
959  // Calculate new velocities via stochastic collisions
960  collisions();
961 
962  // Calculate the volume field data
963  calculateFields();
964 }
965 
966 
967 template<class ParcelType>
969 {
970  label nDSMCParticles = this->size();
971  reduce(nDSMCParticles, sumOp<label>());
972 
973  scalar nMol = nDSMCParticles*nParticle_;
974 
975  vector linearMomentum = linearMomentumOfSystem();
976  reduce(linearMomentum, sumOp<vector>());
977 
978  scalar linearKineticEnergy = linearKineticEnergyOfSystem();
979  reduce(linearKineticEnergy, sumOp<scalar>());
980 
981  scalar internalEnergy = internalEnergyOfSystem();
982  reduce(internalEnergy, sumOp<scalar>());
983 
984  Info<< "Cloud name: " << this->name() << nl
985  << " Number of dsmc particles = "
986  << nDSMCParticles
987  << endl;
988 
989  if (nDSMCParticles)
990  {
991  Info<< " Number of molecules = "
992  << nMol << nl
993  << " Mass in system = "
994  << returnReduce(massInSystem(), sumOp<scalar>()) << nl
995  << " Average linear momentum = "
996  << linearMomentum/nMol << nl
997  << " |Average linear momentum| = "
998  << mag(linearMomentum)/nMol << nl
999  << " Average linear kinetic energy = "
1000  << linearKineticEnergy/nMol << nl
1001  << " Average internal energy = "
1002  << internalEnergy/nMol << nl
1003  << " Average total energy = "
1004  << (internalEnergy + linearKineticEnergy)/nMol
1005  << endl;
1006  }
1007 }
1008 
1009 
1010 template<class ParcelType>
1013  scalar temperature,
1014  scalar mass
1015 )
1016 {
1017  return
1018  sqrt(physicoChemical::k.value()*temperature/mass)
1019  *rndGen_.sampleNormal<vector>();
1020 }
1021 
1022 
1023 template<class ParcelType>
1026  scalar temperature,
1027  direction iDof
1028 )
1029 {
1030  scalar Ei = 0.0;
1031 
1032  if (iDof == 0)
1033  {
1034  return Ei;
1035  }
1036  else if (iDof == 2)
1037  {
1038  // Special case for iDof = 2, i.e. diatomics;
1039  Ei = -log(rndGen_.scalar01())*physicoChemical::k.value()*temperature;
1040  }
1041  else
1042  {
1043  scalar a = 0.5*iDof - 1;
1044  scalar energyRatio;
1045  scalar P = -1;
1046 
1047  do
1048  {
1049  energyRatio = 10*rndGen_.scalar01();
1050  P = pow((energyRatio/a), a)*exp(a - energyRatio);
1051  } while (P < rndGen_.scalar01());
1052 
1053  Ei = energyRatio*physicoChemical::k.value()*temperature;
1054  }
1055 
1056  return Ei;
1057 }
1058 
1059 
1060 template<class ParcelType>
1062 {
1063  OFstream pObj
1064  (
1065  this->db().time().path()/"parcelPositions_"
1066  + this->name() + "_"
1067  + this->db().time().timeName() + ".obj"
1068  );
1069 
1070  forAllConstIter(typename DSMCCloud<ParcelType>, *this, iter)
1071  {
1072  const ParcelType& p = iter();
1073 
1074  pObj<< "v " << p.position().x()
1075  << " " << p.position().y()
1076  << " " << p.position().z()
1077  << nl;
1078  }
1079 
1080  pObj.flush();
1081 }
1082 
1083 
1084 template<class ParcelType>
1086 {
1088 
1089  // Update the cell occupancy field
1090  cellOccupancy_.setSize(mesh_.nCells());
1091  buildCellOccupancy();
1092 
1093  // Update the inflow BCs
1094  this->inflowBoundary().autoMap(mapper);
1095 }
1096 
1097 
1098 // ************************************************************************* //
Collection of constants.
tetrahedron< point, const point & > tetPointRef
Definition: tetPointRef.H:44
virtual void autoMap(const mapPolyMesh &)
Remap the particles to the correct cells following mesh change.
Definition: DSMCCloud.C:1085
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
dimensionedScalar log(const dimensionedScalar &ds)
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
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
void info() const
Print cloud information.
Definition: DSMCCloud.C:968
tUEqn clear()
Templated inflow boundary model class.
Definition: DSMCCloud.H:62
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:453
Output to file stream.
Definition: OFstream.H:82
uint8_t direction
Definition: direction.H:45
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:427
Templated DSMC particle collision class.
Definition: DSMCCloud.H:56
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
void autoMap(const mapPolyMesh &)
Remap the cells of particles corresponding to the.
Definition: Cloud.C:345
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
static List< tetIndices > cellTetIndices(const polyMesh &mesh, label cI)
Return the tet decomposition of the given cell, see.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:52
Dimension set for the base types.
Definition: dimensionSet.H:120
dynamicFvMesh & mesh
dimensionedScalar exp(const dimensionedScalar &ds)
A class for handling words, derived from string.
Definition: word.H:59
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual void flush()
Flush stream.
Definition: OSstream.C:254
const word & constant() const
Return constant name.
Definition: TimePaths.H:124
const Type & value() const
Return const reference to value.
word timeName
Definition: getTimeIndex.H:3
static const zero Zero
Definition: zero.H:97
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if successful.
Definition: doubleScalar.H:68
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
errorManip< error > abort(error &err)
Definition: errorManip.H:131
scalar equipartitionInternalEnergy(scalar temperature, direction internalDegreesOfFreedom)
Generate a random internal energy, sampled from the.
Definition: DSMCCloud.C:1025
label findMax(const ListType &, const label start=0)
Find index of max element (and larger than given element).
dimensionedScalar pos0(const dimensionedScalar &ds)
static const char nl
Definition: Ostream.H:265
const Time & time() const
Return time.
fileName path(UMean.rootPath()/UMean.caseName()/"graphs"/UMean.instance())
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
void addNewParcel(const vector &position, const label celli, const vector &U, const scalar Ei, const label typeId)
Add new parcel.
Definition: DSMCCloud.C:456
virtual ~DSMCCloud()
Destructor.
Definition: DSMCCloud.C:931
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const dimensionedScalar k
Boltzmann constant.
void move(TrackCloudType &cloud, typename ParticleType::trackingData &td, const scalar trackTime)
Move the particles.
Definition: Cloud.C:141
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
vector point
Point is a vector.
Definition: point.H:41
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
const dimensionedScalar c
Speed of light in a vacuum.
void evolve()
Evolve the cloud (move, collide)
Definition: DSMCCloud.C:938
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
Field< vector > vectorField
Specialisation of Field<T> for vector.
Virtual abstract base class for templated DSMCCloud.
Definition: DSMCBaseCloud.H:48
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
const word cloudName(propsDict.lookup("cloudName"))
Templated wall interaction model class.
Definition: DSMCCloud.H:59
void dumpParticlePositions() const
Dump particle positions to .obj file.
Definition: DSMCCloud.C:1061
vector equipartitionLinearVelocity(scalar temperature, scalar mass)
Generate a random velocity sampled from the Maxwellian speed.
Definition: DSMCCloud.C:1012
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:69