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 :
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::MUST_READ_IF_MODIFIED,
477  IOobject::NO_WRITE
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_,
490  IOobject::MUST_READ,
491  IOobject::AUTO_WRITE
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_,
513  IOobject::MUST_READ,
514  IOobject::AUTO_WRITE
515  ),
516  mesh_
517  ),
518  fD_
519  (
520  IOobject
521  (
522  "fD",
523  mesh_.time().timeName(),
524  mesh_,
525  IOobject::MUST_READ,
526  IOobject::AUTO_WRITE
527  ),
528  mesh_
529  ),
530  rhoN_
531  (
532  IOobject
533  (
534  "rhoN",
535  mesh_.time().timeName(),
536  mesh_,
537  IOobject::MUST_READ,
538  IOobject::AUTO_WRITE
539  ),
540  mesh_
541  ),
542  rhoM_
543  (
544  IOobject
545  (
546  "rhoM",
547  mesh_.time().timeName(),
548  mesh_,
549  IOobject::MUST_READ,
550  IOobject::AUTO_WRITE
551  ),
552  mesh_
553  ),
554  dsmcRhoN_
555  (
556  IOobject
557  (
558  "dsmcRhoN",
559  mesh_.time().timeName(),
560  mesh_,
561  IOobject::MUST_READ,
562  IOobject::AUTO_WRITE
563  ),
564  mesh_
565  ),
566  linearKE_
567  (
568  IOobject
569  (
570  "linearKE",
571  mesh_.time().timeName(),
572  mesh_,
573  IOobject::MUST_READ,
574  IOobject::AUTO_WRITE
575  ),
576  mesh_
577  ),
578  internalE_
579  (
580  IOobject
581  (
582  "internalE",
583  mesh_.time().timeName(),
584  mesh_,
585  IOobject::MUST_READ,
586  IOobject::AUTO_WRITE
587  ),
588  mesh_
589  ),
590  iDof_
591  (
592  IOobject
593  (
594  "iDof",
595  mesh_.time().timeName(),
596  mesh_,
597  IOobject::MUST_READ,
598  IOobject::AUTO_WRITE
599  ),
600  mesh_
601  ),
602  momentum_
603  (
604  IOobject
605  (
606  "momentum",
607  mesh_.time().timeName(),
608  mesh_,
609  IOobject::MUST_READ,
610  IOobject::AUTO_WRITE
611  ),
612  mesh_
613  ),
614  constProps_(),
615  rndGen_(Pstream::myProcNo()),
616  boundaryT_
617  (
618  IOobject
619  (
620  "boundaryT",
621  mesh_.time().timeName(),
622  mesh_,
623  IOobject::MUST_READ,
624  IOobject::AUTO_WRITE
625  ),
626  mesh_
627  ),
628  boundaryU_
629  (
630  IOobject
631  (
632  "boundaryU",
633  mesh_.time().timeName(),
634  mesh_,
635  IOobject::MUST_READ,
636  IOobject::AUTO_WRITE
637  ),
638  mesh_
639  ),
640  binaryCollisionModel_
641  (
642  BinaryCollisionModel<DSMCCloud<ParcelType>>::New
643  (
644  particleProperties_,
645  *this
646  )
647  ),
648  wallInteractionModel_
649  (
650  WallInteractionModel<DSMCCloud<ParcelType>>::New
651  (
652  particleProperties_,
653  *this
654  )
655  ),
656  inflowBoundaryModel_
657  (
658  InflowBoundaryModel<DSMCCloud<ParcelType>>::New
659  (
660  particleProperties_,
661  *this
662  )
663  )
664 {
665  buildConstProps();
666  buildCellOccupancy();
667 
668  // Initialise the collision selection remainder to a random value between 0
669  // and 1.
670  forAll(collisionSelectionRemainder_, i)
671  {
672  collisionSelectionRemainder_[i] = rndGen_.sample01<scalar>();
673  }
674 
675  if (readFields)
676  {
678  }
679 }
680 
681 
682 template<class ParcelType>
684 (
685  const word& cloudName,
686  const fvMesh& mesh,
687  const IOdictionary& dsmcInitialiseDict
688 )
689 :
690  Cloud<ParcelType>(mesh, cloudName, false),
691  DSMCBaseCloud(),
692  cloudName_(cloudName),
693  mesh_(mesh),
694  particleProperties_
695  (
696  IOobject
697  (
698  cloudName + "Properties",
699  mesh_.time().constant(),
700  mesh_,
701  IOobject::MUST_READ_IF_MODIFIED,
702  IOobject::NO_WRITE
703  )
704  ),
705  typeIdList_(particleProperties_.lookup("typeIdList")),
706  nParticle_(particleProperties_.get<scalar>("nEquivalentParticles")),
707  cellOccupancy_(),
708  sigmaTcRMax_
709  (
710  IOobject
711  (
712  this->name() + "SigmaTcRMax",
713  mesh_.time().timeName(),
714  mesh_,
715  IOobject::NO_READ,
716  IOobject::AUTO_WRITE
717  ),
718  mesh_,
719  dimensionedScalar(dimensionSet(0, 3, -1, 0, 0), Zero),
720  fvPatchFieldBase::zeroGradientType()
721  ),
722  collisionSelectionRemainder_
723  (
724  IOobject
725  (
726  this->name() + ":collisionSelectionRemainder",
727  mesh_.time().timeName(),
728  mesh_
729  ),
730  mesh_,
732  ),
733  q_
734  (
735  IOobject
736  (
737  this->name() + "q_",
738  mesh_.time().timeName(),
739  mesh_,
740  IOobject::NO_READ,
741  IOobject::NO_WRITE
742  ),
743  mesh_,
744  dimensionedScalar(dimensionSet(1, 0, -3, 0, 0), Zero)
745  ),
746  fD_
747  (
748  IOobject
749  (
750  this->name() + "fD_",
751  mesh_.time().timeName(),
752  mesh_,
753  IOobject::NO_READ,
754  IOobject::NO_WRITE
755  ),
756  mesh_,
757  dimensionedVector(dimensionSet(1, -1, -2, 0, 0), Zero)
758  ),
759  rhoN_
760  (
761  IOobject
762  (
763  this->name() + "rhoN_",
764  mesh_.time().timeName(),
765  mesh_,
766  IOobject::NO_READ,
767  IOobject::NO_WRITE
768  ),
769  mesh_,
770  dimensionedScalar("0", dimensionSet(0, -3, 0, 0, 0), VSMALL)
771  ),
772  rhoM_
773  (
774  IOobject
775  (
776  this->name() + "rhoM_",
777  mesh_.time().timeName(),
778  mesh_,
779  IOobject::NO_READ,
780  IOobject::NO_WRITE
781  ),
782  mesh_,
783  dimensionedScalar("0", dimensionSet(1, -3, 0, 0, 0), VSMALL)
784  ),
785  dsmcRhoN_
786  (
787  IOobject
788  (
789  this->name() + "dsmcRhoN_",
790  mesh_.time().timeName(),
791  mesh_,
792  IOobject::NO_READ,
793  IOobject::NO_WRITE
794  ),
795  mesh_,
796  dimensionedScalar(dimensionSet(0, -3, 0, 0, 0), Zero)
797  ),
798  linearKE_
799  (
800  IOobject
801  (
802  this->name() + "linearKE_",
803  mesh_.time().timeName(),
804  mesh_,
805  IOobject::NO_READ,
806  IOobject::NO_WRITE
807  ),
808  mesh_,
809  dimensionedScalar(dimensionSet(1, -1, -2, 0, 0), Zero)
810  ),
811  internalE_
812  (
813  IOobject
814  (
815  this->name() + "internalE_",
816  mesh_.time().timeName(),
817  mesh_,
818  IOobject::NO_READ,
819  IOobject::NO_WRITE
820  ),
821  mesh_,
822  dimensionedScalar(dimensionSet(1, -1, -2, 0, 0), Zero)
823  ),
824  iDof_
825  (
826  IOobject
827  (
828  this->name() + "iDof_",
829  mesh_.time().timeName(),
830  mesh_,
831  IOobject::NO_READ,
832  IOobject::NO_WRITE
833  ),
834  mesh_,
835  dimensionedScalar("0", dimensionSet(0, -3, 0, 0, 0), VSMALL)
836  ),
837  momentum_
838  (
839  IOobject
840  (
841  this->name() + "momentum_",
842  mesh_.time().timeName(),
843  mesh_,
844  IOobject::NO_READ,
845  IOobject::NO_WRITE
846  ),
847  mesh_,
848  dimensionedVector(dimensionSet(1, -2, -1, 0, 0), Zero)
849  ),
850  constProps_(),
851  rndGen_(Pstream::myProcNo()),
852  boundaryT_
853  (
855  (
856  IOobject
857  (
858  "boundaryT",
859  mesh_.time().timeName(),
860  mesh_,
861  IOobject::NO_READ,
862  IOobject::NO_WRITE
863  ),
864  mesh_,
865  dimensionedScalar(dimensionSet(0, 0, 0, 1, 0), Zero)
866  )
867  ),
868  boundaryU_
869  (
871  (
872  IOobject
873  (
874  "boundaryU",
875  mesh_.time().timeName(),
876  mesh_,
877  IOobject::NO_READ,
878  IOobject::NO_WRITE
879  ),
880  mesh_,
881  dimensionedVector(dimensionSet(0, 1, -1, 0, 0), Zero)
882  )
883  ),
884  binaryCollisionModel_(),
885  wallInteractionModel_(),
886  inflowBoundaryModel_()
887 {
888  clear();
889  buildConstProps();
890  initialise(dsmcInitialiseDict);
891 }
892 
893 
894 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
895 
896 template<class ParcelType>
898 {}
899 
900 
901 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
902 
903 template<class ParcelType>
905 {
906  typename ParcelType::trackingData td(*this);
907 
908  // Reset the data collection fields
909  resetFields();
910 
911  if (debug)
912  {
913  this->dumpParticlePositions();
914  }
915 
916  // Insert new particles from the inflow boundary
917  this->inflowBoundary().inflow();
918 
919  // Move the particles ballistically with their current velocities
920  Cloud<ParcelType>::move(*this, td, mesh_.time().deltaTValue());
921 
922  // Update cell occupancy
923  buildCellOccupancy();
924 
925  // Calculate new velocities via stochastic collisions
926  collisions();
928  // Calculate the volume field data
929  calculateFields();
930 }
931 
932 
933 template<class ParcelType>
935 {
936  label nDSMCParticles = this->size();
937  reduce(nDSMCParticles, sumOp<label>());
938 
939  scalar nMol = nDSMCParticles*nParticle_;
940 
941  vector linearMomentum = linearMomentumOfSystem();
942  reduce(linearMomentum, sumOp<vector>());
943 
944  scalar linearKineticEnergy = linearKineticEnergyOfSystem();
945  reduce(linearKineticEnergy, sumOp<scalar>());
946 
947  scalar internalEnergy = internalEnergyOfSystem();
948  reduce(internalEnergy, sumOp<scalar>());
949 
950  Info<< "Cloud name: " << this->name() << nl
951  << " Number of dsmc particles = "
952  << nDSMCParticles
953  << endl;
954 
955  if (nDSMCParticles)
956  {
957  Info<< " Number of molecules = "
958  << nMol << nl
959  << " Mass in system = "
960  << returnReduce(massInSystem(), sumOp<scalar>()) << nl
961  << " Average linear momentum = "
962  << linearMomentum/nMol << nl
963  << " |Average linear momentum| = "
964  << mag(linearMomentum)/nMol << nl
965  << " Average linear kinetic energy = "
966  << linearKineticEnergy/nMol << nl
967  << " Average internal energy = "
968  << internalEnergy/nMol << nl
969  << " Average total energy = "
970  << (internalEnergy + linearKineticEnergy)/nMol
971  << endl;
972  }
973 }
974 
975 
976 template<class ParcelType>
978 (
979  scalar temperature,
980  scalar mass
981 )
982 {
983  return
984  sqrt(physicoChemical::k.value()*temperature/mass)
985  *rndGen_.GaussNormal<vector>();
986 }
987 
988 
989 template<class ParcelType>
991 (
992  scalar temperature,
993  direction iDof
994 )
995 {
996  if (iDof == 0)
997  {
998  return 0;
999  }
1000  else if (iDof == 2)
1001  {
1002  // Special case for iDof = 2, i.e. diatomics;
1003  return
1004  (
1005  -log(rndGen_.sample01<scalar>())
1006  *physicoChemical::k.value()*temperature
1007  );
1008  }
1009 
1010 
1011  const scalar a = 0.5*iDof - 1;
1012  scalar energyRatio = 0;
1013  scalar P = -1;
1014 
1015  do
1016  {
1017  energyRatio = 10*rndGen_.sample01<scalar>();
1018  P = pow((energyRatio/a), a)*exp(a - energyRatio);
1019  } while (P < rndGen_.sample01<scalar>());
1020 
1021  return energyRatio*physicoChemical::k.value()*temperature;
1022 }
1023 
1024 
1025 template<class ParcelType>
1027 {
1028  OFstream pObj
1029  (
1030  this->db().time().path()/"parcelPositions_"
1031  + this->name() + "_"
1032  + this->db().time().timeName() + ".obj"
1033  );
1034 
1035  for (const ParcelType& p : *this)
1036  {
1037  pObj<< "v " << p.position().x()
1038  << ' ' << p.position().y()
1039  << ' ' << p.position().z()
1040  << nl;
1041  }
1042 
1043  pObj.flush();
1044 }
1045 
1046 
1047 template<class ParcelType>
1048 void Foam::DSMCCloud<ParcelType>::autoMap(const mapPolyMesh& mapper)
1049 {
1051 
1052  // Update the cell occupancy field
1053  cellOccupancy_.setSize(mesh_.nCells());
1054  buildCellOccupancy();
1055 
1056  // Update the inflow BCs
1057  this->inflowBoundary().autoMap(mapper);
1058 }
1059 
1060 
1061 // ************************************************************************* //
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:1041
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:598
void info() const
Print cloud information.
Definition: DSMCCloud.C:927
tetrahedron< point, const point & > tetPointRef
A tetrahedron using referred points.
Definition: tetrahedron.H:72
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:82
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:81
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:57
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:984
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:890
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:897
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
Reduce inplace (cf. MPI Allreduce) using specified communication schedule.
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:1019
vector equipartitionLinearVelocity(scalar temperature, scalar mass)
Generate a random velocity sampled from the Maxwellian speed.
Definition: DSMCCloud.C:971
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:67
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127