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-2023 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 :
465  Cloud<ParcelType>(mesh, cloudName, false),
466  DSMCBaseCloud(),
467  cloudName_(cloudName),
468  mesh_(mesh),
469  particleProperties_
470  (
471  IOobject
472  (
473  cloudName + "Properties",
474  mesh_.time().constant(),
475  mesh_,
476  IOobject::READ_MODIFIED,
477  IOobject::NO_WRITE,
478  IOobject::REGISTER
479  )
480  ),
481  typeIdList_(particleProperties_.lookup("typeIdList")),
482  nParticle_(particleProperties_.get<scalar>("nEquivalentParticles")),
483  cellOccupancy_(mesh_.nCells()),
484  sigmaTcRMax_
485  (
486  IOobject
487  (
488  this->name() + "SigmaTcRMax",
489  mesh_.time().timeName(),
490  mesh_,
491  IOobject::MUST_READ,
492  IOobject::AUTO_WRITE,
493  IOobject::REGISTER
494  ),
495  mesh_
496  ),
497  collisionSelectionRemainder_
498  (
499  IOobject
500  (
501  IOobject::scopedName(this->name(), "collisionSelectionRemainder"),
502  mesh_.time().timeName(),
503  mesh_
504  ),
505  mesh_,
507  ),
508  q_
509  (
510  IOobject
511  (
512  "q",
513  mesh_.time().timeName(),
514  mesh_,
515  IOobject::MUST_READ,
516  IOobject::AUTO_WRITE,
517  IOobject::REGISTER
518  ),
519  mesh_
520  ),
521  fD_
522  (
523  IOobject
524  (
525  "fD",
526  mesh_.time().timeName(),
527  mesh_,
528  IOobject::MUST_READ,
529  IOobject::AUTO_WRITE,
530  IOobject::REGISTER
531  ),
532  mesh_
533  ),
534  rhoN_
535  (
536  IOobject
537  (
538  "rhoN",
539  mesh_.time().timeName(),
540  mesh_,
541  IOobject::MUST_READ,
542  IOobject::AUTO_WRITE,
543  IOobject::REGISTER
544  ),
545  mesh_
546  ),
547  rhoM_
548  (
549  IOobject
550  (
551  "rhoM",
552  mesh_.time().timeName(),
553  mesh_,
554  IOobject::MUST_READ,
555  IOobject::AUTO_WRITE,
556  IOobject::REGISTER
557  ),
558  mesh_
559  ),
560  dsmcRhoN_
561  (
562  IOobject
563  (
564  "dsmcRhoN",
565  mesh_.time().timeName(),
566  mesh_,
567  IOobject::MUST_READ,
568  IOobject::AUTO_WRITE,
569  IOobject::REGISTER
570  ),
571  mesh_
572  ),
573  linearKE_
574  (
575  IOobject
576  (
577  "linearKE",
578  mesh_.time().timeName(),
579  mesh_,
580  IOobject::MUST_READ,
581  IOobject::AUTO_WRITE,
582  IOobject::REGISTER
583  ),
584  mesh_
585  ),
586  internalE_
587  (
588  IOobject
589  (
590  "internalE",
591  mesh_.time().timeName(),
592  mesh_,
593  IOobject::MUST_READ,
594  IOobject::AUTO_WRITE,
595  IOobject::REGISTER
596  ),
597  mesh_
598  ),
599  iDof_
600  (
601  IOobject
602  (
603  "iDof",
604  mesh_.time().timeName(),
605  mesh_,
606  IOobject::MUST_READ,
607  IOobject::AUTO_WRITE,
608  IOobject::REGISTER
609  ),
610  mesh_
611  ),
612  momentum_
613  (
614  IOobject
615  (
616  "momentum",
617  mesh_.time().timeName(),
618  mesh_,
619  IOobject::MUST_READ,
620  IOobject::AUTO_WRITE,
621  IOobject::REGISTER
622  ),
623  mesh_
624  ),
625  constProps_(),
626  rndGen_(Pstream::myProcNo()),
627  boundaryT_
628  (
629  IOobject
630  (
631  "boundaryT",
632  mesh_.time().timeName(),
633  mesh_,
634  IOobject::MUST_READ,
635  IOobject::AUTO_WRITE,
636  IOobject::REGISTER
637  ),
638  mesh_
639  ),
640  boundaryU_
641  (
642  IOobject
643  (
644  "boundaryU",
645  mesh_.time().timeName(),
646  mesh_,
647  IOobject::MUST_READ,
648  IOobject::AUTO_WRITE,
649  IOobject::REGISTER
650  ),
651  mesh_
652  ),
653  binaryCollisionModel_
654  (
655  BinaryCollisionModel<DSMCCloud<ParcelType>>::New
656  (
657  particleProperties_,
658  *this
659  )
660  ),
661  wallInteractionModel_
662  (
663  WallInteractionModel<DSMCCloud<ParcelType>>::New
664  (
665  particleProperties_,
666  *this
667  )
668  ),
669  inflowBoundaryModel_
670  (
671  InflowBoundaryModel<DSMCCloud<ParcelType>>::New
672  (
673  particleProperties_,
674  *this
675  )
676  )
677 {
678  buildConstProps();
679  buildCellOccupancy();
680 
681  // Initialise the collision selection remainder to a random value between 0
682  // and 1.
683  forAll(collisionSelectionRemainder_, i)
684  {
685  collisionSelectionRemainder_[i] = rndGen_.sample01<scalar>();
686  }
687 
688  if (readFields)
689  {
691  }
692 }
693 
694 
695 template<class ParcelType>
697 (
698  const word& cloudName,
699  const fvMesh& mesh,
700  const IOdictionary& dsmcInitialiseDict
701 )
702 :
703  Cloud<ParcelType>(mesh, cloudName, false),
704  DSMCBaseCloud(),
705  cloudName_(cloudName),
706  mesh_(mesh),
707  particleProperties_
708  (
709  IOobject
710  (
711  cloudName + "Properties",
712  mesh_.time().constant(),
713  mesh_,
714  IOobject::READ_MODIFIED,
715  IOobject::NO_WRITE,
716  IOobject::REGISTER
717  )
718  ),
719  typeIdList_(particleProperties_.lookup("typeIdList")),
720  nParticle_(particleProperties_.get<scalar>("nEquivalentParticles")),
721  cellOccupancy_(),
722  sigmaTcRMax_
723  (
724  IOobject
725  (
726  this->name() + "SigmaTcRMax",
727  mesh_.time().timeName(),
728  mesh_,
729  IOobject::NO_READ,
730  IOobject::AUTO_WRITE,
731  IOobject::REGISTER
732  ),
733  mesh_,
734  dimensionedScalar(dimensionSet(0, 3, -1, 0, 0), Zero),
735  fvPatchFieldBase::zeroGradientType()
736  ),
737  collisionSelectionRemainder_
738  (
739  IOobject
740  (
741  IOobject::scopedName(this->name(), "collisionSelectionRemainder"),
742  mesh_.time().timeName(),
743  mesh_
744  ),
745  mesh_,
747  ),
748  q_
749  (
750  IOobject
751  (
752  this->name() + "q_",
753  mesh_.time().timeName(),
754  mesh_,
755  IOobject::NO_READ,
756  IOobject::NO_WRITE
757  ),
758  mesh_,
759  dimensionedScalar(dimensionSet(1, 0, -3, 0, 0), Zero)
760  ),
761  fD_
762  (
763  IOobject
764  (
765  this->name() + "fD_",
766  mesh_.time().timeName(),
767  mesh_,
768  IOobject::NO_READ,
769  IOobject::NO_WRITE
770  ),
771  mesh_,
772  dimensionedVector(dimensionSet(1, -1, -2, 0, 0), Zero)
773  ),
774  rhoN_
775  (
776  IOobject
777  (
778  this->name() + "rhoN_",
779  mesh_.time().timeName(),
780  mesh_,
781  IOobject::NO_READ,
782  IOobject::NO_WRITE
783  ),
784  mesh_,
785  dimensionedScalar("0", dimensionSet(0, -3, 0, 0, 0), VSMALL)
786  ),
787  rhoM_
788  (
789  IOobject
790  (
791  this->name() + "rhoM_",
792  mesh_.time().timeName(),
793  mesh_,
794  IOobject::NO_READ,
795  IOobject::NO_WRITE
796  ),
797  mesh_,
798  dimensionedScalar("0", dimensionSet(1, -3, 0, 0, 0), VSMALL)
799  ),
800  dsmcRhoN_
801  (
802  IOobject
803  (
804  this->name() + "dsmcRhoN_",
805  mesh_.time().timeName(),
806  mesh_,
807  IOobject::NO_READ,
808  IOobject::NO_WRITE
809  ),
810  mesh_,
811  dimensionedScalar(dimensionSet(0, -3, 0, 0, 0), Zero)
812  ),
813  linearKE_
814  (
815  IOobject
816  (
817  this->name() + "linearKE_",
818  mesh_.time().timeName(),
819  mesh_,
820  IOobject::NO_READ,
821  IOobject::NO_WRITE
822  ),
823  mesh_,
824  dimensionedScalar(dimensionSet(1, -1, -2, 0, 0), Zero)
825  ),
826  internalE_
827  (
828  IOobject
829  (
830  this->name() + "internalE_",
831  mesh_.time().timeName(),
832  mesh_,
833  IOobject::NO_READ,
834  IOobject::NO_WRITE
835  ),
836  mesh_,
837  dimensionedScalar(dimensionSet(1, -1, -2, 0, 0), Zero)
838  ),
839  iDof_
840  (
841  IOobject
842  (
843  this->name() + "iDof_",
844  mesh_.time().timeName(),
845  mesh_,
846  IOobject::NO_READ,
847  IOobject::NO_WRITE
848  ),
849  mesh_,
850  dimensionedScalar("0", dimensionSet(0, -3, 0, 0, 0), VSMALL)
851  ),
852  momentum_
853  (
854  IOobject
855  (
856  this->name() + "momentum_",
857  mesh_.time().timeName(),
858  mesh_,
859  IOobject::NO_READ,
860  IOobject::NO_WRITE
861  ),
862  mesh_,
863  dimensionedVector(dimensionSet(1, -2, -1, 0, 0), Zero)
864  ),
865  constProps_(),
866  rndGen_(Pstream::myProcNo()),
867  boundaryT_
868  (
870  (
871  IOobject
872  (
873  "boundaryT",
874  mesh_.time().timeName(),
875  mesh_,
876  IOobject::NO_READ,
877  IOobject::NO_WRITE
878  ),
879  mesh_,
880  dimensionedScalar(dimensionSet(0, 0, 0, 1, 0), Zero)
881  )
882  ),
883  boundaryU_
884  (
886  (
887  IOobject
888  (
889  "boundaryU",
890  mesh_.time().timeName(),
891  mesh_,
892  IOobject::NO_READ,
893  IOobject::NO_WRITE
894  ),
895  mesh_,
896  dimensionedVector(dimensionSet(0, 1, -1, 0, 0), Zero)
897  )
898  ),
899  binaryCollisionModel_(),
900  wallInteractionModel_(),
901  inflowBoundaryModel_()
902 {
903  clear();
904  buildConstProps();
905  initialise(dsmcInitialiseDict);
906 }
907 
908 
909 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
910 
911 template<class ParcelType>
913 {}
914 
915 
916 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
917 
918 template<class ParcelType>
920 {
921  typename ParcelType::trackingData td(*this);
922 
923  // Reset the data collection fields
924  resetFields();
925 
926  if (debug)
927  {
928  this->dumpParticlePositions();
929  }
930 
931  // Insert new particles from the inflow boundary
932  this->inflowBoundary().inflow();
933 
934  // Move the particles ballistically with their current velocities
935  Cloud<ParcelType>::move(*this, td, mesh_.time().deltaTValue());
936 
937  // Update cell occupancy
938  buildCellOccupancy();
939 
940  // Calculate new velocities via stochastic collisions
941  collisions();
943  // Calculate the volume field data
944  calculateFields();
945 }
946 
947 
948 template<class ParcelType>
950 {
951  label nDSMCParticles = this->size();
952  reduce(nDSMCParticles, sumOp<label>());
953 
954  scalar nMol = nDSMCParticles*nParticle_;
955 
956  vector linearMomentum = linearMomentumOfSystem();
957  reduce(linearMomentum, sumOp<vector>());
958 
959  scalar linearKineticEnergy = linearKineticEnergyOfSystem();
960  reduce(linearKineticEnergy, sumOp<scalar>());
961 
962  scalar internalEnergy = internalEnergyOfSystem();
963  reduce(internalEnergy, sumOp<scalar>());
964 
965  Info<< "Cloud name: " << this->name() << nl
966  << " Number of dsmc particles = "
967  << nDSMCParticles
968  << endl;
969 
970  if (nDSMCParticles)
971  {
972  Info<< " Number of molecules = "
973  << nMol << nl
974  << " Mass in system = "
975  << returnReduce(massInSystem(), sumOp<scalar>()) << nl
976  << " Average linear momentum = "
977  << linearMomentum/nMol << nl
978  << " |Average linear momentum| = "
979  << mag(linearMomentum)/nMol << nl
980  << " Average linear kinetic energy = "
981  << linearKineticEnergy/nMol << nl
982  << " Average internal energy = "
983  << internalEnergy/nMol << nl
984  << " Average total energy = "
985  << (internalEnergy + linearKineticEnergy)/nMol
986  << endl;
987  }
988 }
989 
990 
991 template<class ParcelType>
993 (
994  scalar temperature,
995  scalar mass
996 )
997 {
998  return
999  sqrt(physicoChemical::k.value()*temperature/mass)
1000  *rndGen_.GaussNormal<vector>();
1001 }
1002 
1003 
1004 template<class ParcelType>
1006 (
1007  scalar temperature,
1008  direction iDof
1009 )
1010 {
1011  if (iDof == 0)
1012  {
1013  return 0;
1014  }
1015  else if (iDof == 2)
1016  {
1017  // Special case for iDof = 2, i.e. diatomics;
1018  return
1019  (
1020  -log(rndGen_.sample01<scalar>())
1021  *physicoChemical::k.value()*temperature
1022  );
1023  }
1024 
1025 
1026  const scalar a = 0.5*iDof - 1;
1027  scalar energyRatio = 0;
1028  scalar P = -1;
1029 
1030  do
1031  {
1032  energyRatio = 10*rndGen_.sample01<scalar>();
1033  P = pow((energyRatio/a), a)*exp(a - energyRatio);
1034  } while (P < rndGen_.sample01<scalar>());
1035 
1036  return energyRatio*physicoChemical::k.value()*temperature;
1037 }
1038 
1039 
1040 template<class ParcelType>
1042 {
1043  OFstream pObj
1044  (
1045  this->db().time().path()/"parcelPositions_"
1046  + this->name() + "_"
1047  + this->db().time().timeName() + ".obj"
1048  );
1049 
1050  for (const ParcelType& p : *this)
1051  {
1052  pObj<< "v " << p.position().x()
1053  << ' ' << p.position().y()
1054  << ' ' << p.position().z()
1055  << nl;
1056  }
1057 
1058  pObj.flush();
1059 }
1060 
1061 
1062 template<class ParcelType>
1063 void Foam::DSMCCloud<ParcelType>::autoMap(const mapPolyMesh& mapper)
1064 {
1066 
1067  // Update the cell occupancy field
1068  cellOccupancy_.setSize(mesh_.nCells());
1069  buildCellOccupancy();
1070 
1071  // Update the inflow BCs
1072  this->inflowBoundary().autoMap(mapper);
1073 }
1074 
1075 
1076 // ************************************************************************* //
Different types of constants.
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
label findMax(const ListType &input, label start=0)
Linear search for the index of the max element, similar to std::max_element but for lists and returns...
virtual void autoMap(const mapPolyMesh &)
Remap the particles to the correct cells following mesh change.
Definition: DSMCCloud.C:1056
uint8_t direction
Definition: direction.H:46
const word zeroGradientType
A zeroGradient patch field type.
dimensionedScalar log(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:608
void info() const
Print cloud information.
Definition: DSMCCloud.C:942
tetrahedron< point, const point & > tetPointRef
A tetrahedron using referred points.
Definition: tetrahedron.H:72
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
Templated inflow boundary model class.
Definition: DSMCCloud.H:60
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Templated DSMC particle collision class.
Definition: DSMCCloud.H:54
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:531
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
label k
Boltzmann constant.
const dimensionSet dimless
Dimensionless.
Lookup type of boundary radiation properties.
Definition: lookup.H:57
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:76
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
void autoMap(const mapPolyMesh &)
Remap the cells of particles corresponding to the.
Definition: Cloud.C:297
static List< tetIndices > cellTetIndices(const polyMesh &mesh, label cI)
Return the tet decomposition of the given cell, see.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:72
word timeName
Definition: getTimeIndex.H:3
Type sample01()
Return a sample whose components lie in the range [0,1].
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
dimensionedScalar exp(const dimensionedScalar &ds)
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const NameMatchPredicate &selectedFields, DynamicList< regIOobject *> &storedObjects)
Read the selected GeometricFields of the templated type and store on the objectRegistry.
const word cloudName(propsDict.get< word >("cloud"))
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Inter-processor communications stream.
Definition: Pstream.H:54
Vector< scalar > vector
Definition: vector.H:57
errorManip< error > abort(error &err)
Definition: errorManip.H:139
Base cloud calls templated on particle type.
Definition: Cloud.H:51
scalar equipartitionInternalEnergy(scalar temperature, direction internalDegreesOfFreedom)
Generate a random internal energy, sampled from the.
Definition: DSMCCloud.C:999
void clear()
Clear the Cloud.
Definition: DSMCCloudI.H:468
dimensionedScalar pos0(const dimensionedScalar &ds)
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
int debug
Static debugging option.
fileName path(UMean.rootPath()/UMean.caseName()/"graphs"/UMean.instance())
void addNewParcel(const vector &position, const label celli, const vector &U, const scalar Ei, const label typeId)
Add new parcel.
Definition: DSMCCloud.C:436
virtual ~DSMCCloud()
Destructor.
Definition: DSMCCloud.C:905
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
void move(TrackCloudType &cloud, typename ParticleType::trackingData &td, const scalar trackTime)
Move the particles.
Definition: Cloud.C:137
U
Definition: pEqn.H:72
vector point
Point is a vector.
Definition: point.H:37
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:912
messageStream Info
Information stream (stdout output on master, null elsewhere)
Field< vector > vectorField
Specialisation of Field<T> for vector.
Virtual abstract base class for templated DSMCCloud.
Definition: DSMCBaseCloud.H:46
volScalarField & p
Templated wall interaction model class.
Definition: DSMCCloud.H:57
void dumpParticlePositions() const
Dump particle positions to .obj file.
Definition: DSMCCloud.C:1034
vector equipartitionLinearVelocity(scalar temperature, scalar mass)
Generate a random velocity sampled from the Maxwellian speed.
Definition: DSMCCloud.C:986
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:67
void reduce(T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce) using linear/tree communication schedule.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127