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-2024 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 "NoBinaryCollision.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  dsmcInitialiseDict.template lookup<scalar>("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] =
105  numberDensitiesDict.lookup<scalar>(molecules[i]);
106  }
107 
108  numberDensities /= nParticle_;
109 
110  label nLocateBoundaryHits = 0;
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, nLocateBoundaryHits, U, Ei, typeId);
180  }
181  }
182  }
183  }
184 
185  reduce(nLocateBoundaryHits, sumOp<label>());
186  if (nLocateBoundaryHits != 0)
187  {
189  << "Initialisation of cloud " << this->name()
190  << " did not accurately locate " << nLocateBoundaryHits
191  << " particles" << endl;
192  }
193 
194  // Initialise the sigmaTcRMax_ field to the product of the cross section of
195  // the most abundant species and the most probable thermal speed (Bird,
196  // p222-223)
197 
198  label mostAbundantType(findMax(numberDensities));
199 
200  const typename ParcelType::constantProperties& cP = constProps
201  (
202  mostAbundantType
203  );
204 
205  sigmaTcRMax_.primitiveFieldRef() = cP.sigmaT()*maxwellianMostProbableSpeed
206  (
207  temperature,
208  cP.mass()
209  );
210 
211  sigmaTcRMax_.correctBoundaryConditions();
212 }
213 
214 
215 template<class ParcelType>
217 {
218  if (isType<NoBinaryCollision<DSMCCloud<ParcelType>>>(binaryCollision()))
219  {
220  return;
221  }
222 
223  // Temporary storage for subCells
224  List<DynamicList<label>> subCells(8);
225 
226  scalar deltaT = mesh().time().deltaTValue();
227 
228  label collisionCandidates = 0;
229 
230  label collisions = 0;
231 
232  forAll(cellOccupancy_, celli)
233  {
234  const DynamicList<ParcelType*>& cellParcels(cellOccupancy_[celli]);
235 
236  label nC(cellParcels.size());
237 
238  if (nC > 1)
239  {
240  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
241  // Assign particles to one of 8 Cartesian subCells
242 
243  // Clear temporary lists
244  forAll(subCells, i)
245  {
246  subCells[i].clear();
247  }
248 
249  // Inverse addressing specifying which subCell a parcel is in
250  List<label> whichSubCell(cellParcels.size());
251 
252  const point& cC = mesh_.cellCentres()[celli];
253 
254  forAll(cellParcels, i)
255  {
256  const ParcelType& p = *cellParcels[i];
257  vector relPos = p.position(mesh()) - cC;
258 
259  label subCell =
260  pos0(relPos.x()) + 2*pos0(relPos.y()) + 4*pos0(relPos.z());
261 
262  subCells[subCell].append(i);
263  whichSubCell[i] = subCell;
264  }
265 
266  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
267 
268  scalar sigmaTcRMax = sigmaTcRMax_[celli];
269 
270  scalar selectedPairs =
271  collisionSelectionRemainder_[celli]
272  + 0.5*nC*(nC - 1)*nParticle_*sigmaTcRMax*deltaT
273  /mesh_.cellVolumes()[celli];
274 
275  label nCandidates(selectedPairs);
276  collisionSelectionRemainder_[celli] = selectedPairs - nCandidates;
277  collisionCandidates += nCandidates;
278 
279  for (label c = 0; c < nCandidates; c++)
280  {
281  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
282  // subCell candidate selection procedure
283 
284  // Select the first collision candidate
285  label candidateP = rndGen_.sampleAB<label>(0, nC);
286 
287  // Declare the second collision candidate
288  label candidateQ = -1;
289 
290  List<label> subCellPs = subCells[whichSubCell[candidateP]];
291  label nSC = subCellPs.size();
292 
293  if (nSC > 1)
294  {
295  // If there are two or more particle in a subCell, choose
296  // another from the same cell. If the same candidate is
297  // chosen, choose again.
298 
299  do
300  {
301  candidateQ = subCellPs[rndGen_.sampleAB<label>(0, nSC)];
302  } while (candidateP == candidateQ);
303  }
304  else
305  {
306  // Select a possible second collision candidate from the
307  // whole cell. If the same candidate is chosen, choose
308  // again.
309 
310  do
311  {
312  candidateQ = rndGen_.sampleAB<label>(0, nC);
313  } while (candidateP == candidateQ);
314  }
315 
316  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
317  // uniform candidate selection procedure
318 
319  // // Select the first collision candidate
320  // label candidateP = rndGen_.sampleAB<label>(0, nC);
321 
322  // // Select a possible second collision candidate
323  // label candidateQ = rndGen_.sampleAB<label>(0, nC);
324 
325  // // If the same candidate is chosen, choose again
326  // while (candidateP == candidateQ)
327  // {
328  // candidateQ = rndGen_.sampleAB<label>(0, nC);
329  // }
330 
331  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
332 
333  ParcelType& parcelP = *cellParcels[candidateP];
334  ParcelType& parcelQ = *cellParcels[candidateQ];
335 
336  scalar sigmaTcR = binaryCollision().sigmaTcR
337  (
338  parcelP,
339  parcelQ
340  );
341 
342  // Update the maximum value of sigmaTcR stored, but use the
343  // initial value in the acceptance-rejection criteria because
344  // the number of collision candidates selected was based on this
345 
346  if (sigmaTcR > sigmaTcRMax_[celli])
347  {
348  sigmaTcRMax_[celli] = sigmaTcR;
349  }
350 
351  if ((sigmaTcR/sigmaTcRMax) > rndGen_.scalar01())
352  {
353  binaryCollision().collide
354  (
355  parcelP,
356  parcelQ
357  );
358 
359  collisions++;
360  }
361  }
362  }
363  }
364 
365  reduce(collisions, sumOp<label>());
366 
367  reduce(collisionCandidates, sumOp<label>());
368 
369  sigmaTcRMax_.correctBoundaryConditions();
370 
371  if (collisionCandidates)
372  {
373  Info<< " Collisions = "
374  << collisions << nl
375  << " Acceptance rate = "
376  << scalar(collisions)/scalar(collisionCandidates) << nl
377  << endl;
378  }
379  else
380  {
381  Info<< " No collisions" << endl;
382  }
383 }
384 
385 
386 template<class ParcelType>
388 {
389  q_ = dimensionedScalar( dimensionSet(1, 0, -3, 0, 0), 0);
390 
391  fD_ = dimensionedVector
392  (
393  "zero",
394  dimensionSet(1, -1, -2, 0, 0),
395  Zero
396  );
397 
398  rhoN_ = dimensionedScalar( dimensionSet(0, -3, 0, 0, 0), vSmall);
399  rhoM_ = dimensionedScalar( dimensionSet(1, -3, 0, 0, 0), vSmall);
400  dsmcRhoN_ = dimensionedScalar( dimensionSet(0, -3, 0, 0, 0), 0);
401  linearKE_ = dimensionedScalar( dimensionSet(1, -1, -2, 0, 0), 0);
402  internalE_ = dimensionedScalar( dimensionSet(1, -1, -2, 0, 0), 0);
403  iDof_ = dimensionedScalar( dimensionSet(0, -3, 0, 0, 0), vSmall);
404 
405  momentum_ = dimensionedVector
406  (
407  "zero",
408  dimensionSet(1, -2, -1, 0, 0),
409  Zero
410  );
411 }
412 
413 
414 template<class ParcelType>
416 {
417  scalarField& rhoN = rhoN_.primitiveFieldRef();
418  scalarField& rhoM = rhoM_.primitiveFieldRef();
419  scalarField& dsmcRhoN = dsmcRhoN_.primitiveFieldRef();
420  scalarField& linearKE = linearKE_.primitiveFieldRef();
421  scalarField& internalE = internalE_.primitiveFieldRef();
422  scalarField& iDof = iDof_.primitiveFieldRef();
423  vectorField& momentum = momentum_.primitiveFieldRef();
424 
425  forAllConstIter(typename DSMCCloud<ParcelType>, *this, iter)
426  {
427  const ParcelType& p = iter();
428  const label celli = p.cell();
429 
430  rhoN[celli]++;
431  rhoM[celli] += constProps(p.typeId()).mass();
432  dsmcRhoN[celli]++;
433  linearKE[celli] += 0.5*constProps(p.typeId()).mass()*(p.U() & p.U());
434  internalE[celli] += p.Ei();
435  iDof[celli] += constProps(p.typeId()).internalDegreesOfFreedom();
436  momentum[celli] += constProps(p.typeId()).mass()*p.U();
437  }
438 
439  rhoN *= nParticle_/mesh().cellVolumes();
440  rhoN_.correctBoundaryConditions();
441 
442  rhoM *= nParticle_/mesh().cellVolumes();
443  rhoM_.correctBoundaryConditions();
444 
445  dsmcRhoN_.correctBoundaryConditions();
446 
447  linearKE *= nParticle_/mesh().cellVolumes();
448  linearKE_.correctBoundaryConditions();
449 
450  internalE *= nParticle_/mesh().cellVolumes();
451  internalE_.correctBoundaryConditions();
452 
453  iDof *= nParticle_/mesh().cellVolumes();
454  iDof_.correctBoundaryConditions();
455 
456  momentum *= nParticle_/mesh().cellVolumes();
457  momentum_.correctBoundaryConditions();
458 }
459 
460 
461 // * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * * //
462 
463 template<class ParcelType>
465 (
466  const vector& position,
467  const label celli,
468  label& nLocateBoundaryHits,
469  const vector& U,
470  const scalar Ei,
471  const label typeId
472 )
473 {
474  this->addParticle
475  (
476  new ParcelType
477  (
478  mesh_,
479  position,
480  celli,
481  nLocateBoundaryHits,
482  U,
483  Ei,
484  typeId
485  )
486  );
487 }
488 
489 
490 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
491 
492 template<class ParcelType>
494 (
495  const word& cloudName,
496  const fvMesh& mesh,
497  bool readFields
498 )
499 :
500  Cloud<ParcelType>(mesh, cloudName, false),
501  cloudName_(cloudName),
502  mesh_(mesh),
503  particleProperties_
504  (
505  IOobject
506  (
507  cloudName + "Properties",
508  mesh_.time().constant(),
509  mesh_,
510  IOobject::MUST_READ_IF_MODIFIED,
511  IOobject::NO_WRITE
512  )
513  ),
514  typeIdList_(particleProperties_.lookup("typeIdList")),
515  nParticle_
516  (
517  particleProperties_.template lookup<scalar>("nEquivalentParticles")
518  ),
519  cellOccupancy_(mesh_.nCells()),
520  sigmaTcRMax_
521  (
522  IOobject
523  (
524  this->name() + "SigmaTcRMax",
525  mesh_.time().name(),
526  mesh_,
527  IOobject::MUST_READ,
528  IOobject::AUTO_WRITE
529  ),
530  mesh_
531  ),
532  collisionSelectionRemainder_
533  (
534  IOobject
535  (
536  this->name() + ":collisionSelectionRemainder",
537  mesh_.time().name(),
538  mesh_
539  ),
540  mesh_,
542  ),
543  q_
544  (
545  IOobject
546  (
547  "q",
548  mesh_.time().name(),
549  mesh_,
550  IOobject::MUST_READ,
551  IOobject::AUTO_WRITE
552  ),
553  mesh_
554  ),
555  fD_
556  (
557  IOobject
558  (
559  "fD",
560  mesh_.time().name(),
561  mesh_,
562  IOobject::MUST_READ,
563  IOobject::AUTO_WRITE
564  ),
565  mesh_
566  ),
567  rhoN_
568  (
569  IOobject
570  (
571  "rhoN",
572  mesh_.time().name(),
573  mesh_,
574  IOobject::MUST_READ,
575  IOobject::AUTO_WRITE
576  ),
577  mesh_
578  ),
579  rhoM_
580  (
581  IOobject
582  (
583  "rhoM",
584  mesh_.time().name(),
585  mesh_,
586  IOobject::MUST_READ,
587  IOobject::AUTO_WRITE
588  ),
589  mesh_
590  ),
591  dsmcRhoN_
592  (
593  IOobject
594  (
595  "dsmcRhoN",
596  mesh_.time().name(),
597  mesh_,
598  IOobject::MUST_READ,
599  IOobject::AUTO_WRITE
600  ),
601  mesh_
602  ),
603  linearKE_
604  (
605  IOobject
606  (
607  "linearKE",
608  mesh_.time().name(),
609  mesh_,
610  IOobject::MUST_READ,
611  IOobject::AUTO_WRITE
612  ),
613  mesh_
614  ),
615  internalE_
616  (
617  IOobject
618  (
619  "internalE",
620  mesh_.time().name(),
621  mesh_,
622  IOobject::MUST_READ,
623  IOobject::AUTO_WRITE
624  ),
625  mesh_
626  ),
627  iDof_
628  (
629  IOobject
630  (
631  "iDof",
632  mesh_.time().name(),
633  mesh_,
634  IOobject::MUST_READ,
635  IOobject::AUTO_WRITE
636  ),
637  mesh_
638  ),
639  momentum_
640  (
641  IOobject
642  (
643  "momentum",
644  mesh_.time().name(),
645  mesh_,
646  IOobject::MUST_READ,
647  IOobject::AUTO_WRITE
648  ),
649  mesh_
650  ),
651  constProps_(),
652  rndGen_(label(149382906)),
653  stdNormal_(rndGen_.generator()),
654  boundaryT_
655  (
657  (
658  IOobject
659  (
660  "boundaryT",
661  mesh_.time().name(),
662  mesh_,
663  IOobject::MUST_READ,
664  IOobject::AUTO_WRITE
665  ),
666  mesh_
667  )
668  ),
669  boundaryU_
670  (
672  (
673  IOobject
674  (
675  "boundaryU",
676  mesh_.time().name(),
677  mesh_,
678  IOobject::MUST_READ,
679  IOobject::AUTO_WRITE
680  ),
681  mesh_
682  )
683  ),
684  binaryCollisionModel_
685  (
686  BinaryCollisionModel<DSMCCloud<ParcelType>>::New
687  (
688  particleProperties_,
689  *this
690  )
691  ),
692  wallInteractionModel_
693  (
694  WallInteractionModel<DSMCCloud<ParcelType>>::New
695  (
696  particleProperties_,
697  *this
698  )
699  ),
700  inflowBoundaryModel_
701  (
702  InflowBoundaryModel<DSMCCloud<ParcelType>>::New
703  (
704  particleProperties_,
705  *this
706  )
707  )
708 {
709  buildConstProps();
710  buildCellOccupancy();
711 
712  // Initialise the collision selection remainder to a random value between 0
713  // and 1.
714  forAll(collisionSelectionRemainder_, i)
715  {
716  collisionSelectionRemainder_[i] = rndGen_.scalar01();
717  }
718 
719  if (readFields)
720  {
721  ParcelType::readFields(*this);
722  }
723 }
724 
725 
726 template<class ParcelType>
728 (
729  const word& cloudName,
730  const fvMesh& mesh,
731  const IOdictionary& dsmcInitialiseDict
732 )
733 :
734  Cloud<ParcelType>(mesh, cloudName, false),
735  cloudName_(cloudName),
736  mesh_(mesh),
737  particleProperties_
738  (
739  IOobject
740  (
741  cloudName + "Properties",
742  mesh_.time().constant(),
743  mesh_,
744  IOobject::MUST_READ_IF_MODIFIED,
745  IOobject::NO_WRITE
746  )
747  ),
748  typeIdList_(particleProperties_.lookup("typeIdList")),
749  nParticle_
750  (
751  particleProperties_.template lookup<scalar>("nEquivalentParticles")
752  ),
753  cellOccupancy_(),
754  sigmaTcRMax_
755  (
756  IOobject
757  (
758  this->name() + "SigmaTcRMax",
759  mesh_.time().name(),
760  mesh_,
761  IOobject::NO_READ,
762  IOobject::AUTO_WRITE
763  ),
764  mesh_,
765  dimensionedScalar( dimensionSet(0, 3, -1, 0, 0), 0),
766  zeroGradientFvPatchScalarField::typeName
767  ),
768  collisionSelectionRemainder_
769  (
770  IOobject
771  (
772  this->name() + ":collisionSelectionRemainder",
773  mesh_.time().name(),
774  mesh_
775  ),
776  mesh_,
778  ),
779  q_
780  (
781  IOobject
782  (
783  this->name() + "q_",
784  mesh_.time().name(),
785  mesh_,
786  IOobject::NO_READ,
787  IOobject::NO_WRITE
788  ),
789  mesh_,
790  dimensionedScalar( dimensionSet(1, 0, -3, 0, 0), 0)
791  ),
792  fD_
793  (
794  IOobject
795  (
796  this->name() + "fD_",
797  mesh_.time().name(),
798  mesh_,
799  IOobject::NO_READ,
800  IOobject::NO_WRITE
801  ),
802  mesh_,
804  (
805  "zero",
806  dimensionSet(1, -1, -2, 0, 0),
807  Zero
808  )
809  ),
810  rhoN_
811  (
812  IOobject
813  (
814  this->name() + "rhoN_",
815  mesh_.time().name(),
816  mesh_,
817  IOobject::NO_READ,
818  IOobject::NO_WRITE
819  ),
820  mesh_,
821  dimensionedScalar( dimensionSet(0, -3, 0, 0, 0), vSmall)
822  ),
823  rhoM_
824  (
825  IOobject
826  (
827  this->name() + "rhoM_",
828  mesh_.time().name(),
829  mesh_,
830  IOobject::NO_READ,
831  IOobject::NO_WRITE
832  ),
833  mesh_,
834  dimensionedScalar( dimensionSet(1, -3, 0, 0, 0), vSmall)
835  ),
836  dsmcRhoN_
837  (
838  IOobject
839  (
840  this->name() + "dsmcRhoN_",
841  mesh_.time().name(),
842  mesh_,
843  IOobject::NO_READ,
844  IOobject::NO_WRITE
845  ),
846  mesh_,
847  dimensionedScalar( dimensionSet(0, -3, 0, 0, 0), 0)
848  ),
849  linearKE_
850  (
851  IOobject
852  (
853  this->name() + "linearKE_",
854  mesh_.time().name(),
855  mesh_,
856  IOobject::NO_READ,
857  IOobject::NO_WRITE
858  ),
859  mesh_,
860  dimensionedScalar( dimensionSet(1, -1, -2, 0, 0), 0)
861  ),
862  internalE_
863  (
864  IOobject
865  (
866  this->name() + "internalE_",
867  mesh_.time().name(),
868  mesh_,
869  IOobject::NO_READ,
870  IOobject::NO_WRITE
871  ),
872  mesh_,
873  dimensionedScalar( dimensionSet(1, -1, -2, 0, 0), 0)
874  ),
875  iDof_
876  (
877  IOobject
878  (
879  this->name() + "iDof_",
880  mesh_.time().name(),
881  mesh_,
882  IOobject::NO_READ,
883  IOobject::NO_WRITE
884  ),
885  mesh_,
886  dimensionedScalar( dimensionSet(0, -3, 0, 0, 0), vSmall)
887  ),
888  momentum_
889  (
890  IOobject
891  (
892  this->name() + "momentum_",
893  mesh_.time().name(),
894  mesh_,
895  IOobject::NO_READ,
896  IOobject::NO_WRITE
897  ),
898  mesh_,
900  (
901  "zero",
902  dimensionSet(1, -2, -1, 0, 0),
903  Zero
904  )
905  ),
906  constProps_(),
907  rndGen_(label(971501)),
908  stdNormal_(rndGen_.generator()),
909  boundaryT_
910  (
912  (
913  IOobject
914  (
915  "boundaryT",
916  mesh_.time().name(),
917  mesh_,
918  IOobject::NO_READ,
919  IOobject::NO_WRITE
920  ),
921  mesh_,
922  dimensionedScalar( dimensionSet(0, 0, 0, 1, 0), 0)
923  )
924  ),
925  boundaryU_
926  (
928  (
929  IOobject
930  (
931  "boundaryU",
932  mesh_.time().name(),
933  mesh_,
934  IOobject::NO_READ,
935  IOobject::NO_WRITE
936  ),
937  mesh_,
939  (
940  "zero",
941  dimensionSet(0, 1, -1, 0, 0),
942  Zero
943  )
944  )
945  ),
946  binaryCollisionModel_(),
947  wallInteractionModel_(),
948  inflowBoundaryModel_()
949 {
950  clear();
951  buildConstProps();
952  initialise(dsmcInitialiseDict);
953 }
954 
955 
956 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
957 
958 template<class ParcelType>
960 {}
961 
962 
963 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
964 
965 template<class ParcelType>
967 {
968  typename ParcelType::trackingData td(*this);
969 
970  // Reset the data collection fields
971  resetFields();
972 
973  if (debug)
974  {
975  this->dumpParticlePositions();
976  }
977 
978  // Insert new particles from the inflow boundary
979  this->inflowBoundary().inflow();
980 
981  // Move the particles ballistically with their current velocities
982  Cloud<ParcelType>::move(*this, td);
983 
984  // Update cell occupancy
985  buildCellOccupancy();
986 
987  // Calculate new velocities via stochastic collisions
988  collisions();
989 
990  // Calculate the volume field data
991  calculateFields();
992 }
993 
994 
995 template<class ParcelType>
997 {
998  label nDSMCParticles = this->size();
999  reduce(nDSMCParticles, sumOp<label>());
1000 
1001  scalar nMol = nDSMCParticles*nParticle_;
1002 
1003  vector linearMomentum = linearMomentumOfSystem();
1004  reduce(linearMomentum, sumOp<vector>());
1005 
1006  scalar linearKineticEnergy = linearKineticEnergyOfSystem();
1007  reduce(linearKineticEnergy, sumOp<scalar>());
1008 
1009  scalar internalEnergy = internalEnergyOfSystem();
1010  reduce(internalEnergy, sumOp<scalar>());
1011 
1012  Info<< "Cloud name: " << this->name() << nl
1013  << " Number of dsmc particles = "
1014  << nDSMCParticles
1015  << endl;
1016 
1017  if (nDSMCParticles)
1018  {
1019  Info<< " Number of molecules = "
1020  << nMol << nl
1021  << " Mass in system = "
1022  << returnReduce(massInSystem(), sumOp<scalar>()) << nl
1023  << " Average linear momentum = "
1024  << linearMomentum/nMol << nl
1025  << " |Average linear momentum| = "
1026  << mag(linearMomentum)/nMol << nl
1027  << " Average linear kinetic energy = "
1028  << linearKineticEnergy/nMol << nl
1029  << " Average internal energy = "
1030  << internalEnergy/nMol << nl
1031  << " Average total energy = "
1032  << (internalEnergy + linearKineticEnergy)/nMol
1033  << endl;
1034  }
1035 }
1036 
1037 
1038 template<class ParcelType>
1040 (
1041  scalar temperature,
1042  scalar mass
1043 )
1044 {
1045  return
1046  sqrt(physicoChemical::k.value()*temperature/mass)
1047  *stdNormal_.sample<vector>();
1048 }
1049 
1050 
1051 template<class ParcelType>
1053 (
1054  scalar temperature,
1055  direction iDof
1056 )
1057 {
1058  scalar Ei = 0.0;
1059 
1060  if (iDof == 0)
1061  {
1062  return Ei;
1063  }
1064  else if (iDof == 2)
1065  {
1066  // Special case for iDof = 2, i.e. diatomics;
1067  Ei = -log(rndGen_.scalar01())*physicoChemical::k.value()*temperature;
1068  }
1069  else
1070  {
1071  scalar a = 0.5*iDof - 1;
1072  scalar energyRatio;
1073  scalar P = -1;
1074 
1075  do
1076  {
1077  energyRatio = 10*rndGen_.scalar01();
1078  P = pow((energyRatio/a), a)*exp(a - energyRatio);
1079  } while (P < rndGen_.scalar01());
1080 
1081  Ei = energyRatio*physicoChemical::k.value()*temperature;
1082  }
1083 
1084  return Ei;
1085 }
1086 
1087 
1088 template<class ParcelType>
1090 {
1091  OFstream pObj
1092  (
1093  this->db().time().path()/"parcelPositions_"
1094  + this->name() + "_"
1095  + this->db().time().name() + ".obj"
1096  );
1097 
1098  forAllConstIter(typename DSMCCloud<ParcelType>, *this, iter)
1099  {
1100  const point pos = iter().position(mesh());
1101 
1102  pObj<< "v " << pos.x() << " " << pos.y() << " " << pos.z() << nl;
1103  }
1104 
1105  pObj.flush();
1106 }
1107 
1108 
1109 template<class ParcelType>
1111 {
1113 
1114  // Update the cell occupancy field
1115  cellOccupancy_.setSize(mesh_.nCells());
1116  buildCellOccupancy();
1117 
1118  // Update the inflow BCs
1119  this->inflowBoundary().topoChange();
1120 }
1121 
1122 
1123 template<class ParcelType>
1125 {
1127 
1128  // Update the cell occupancy field
1129  cellOccupancy_.setSize(mesh_.nCells());
1130  buildCellOccupancy();
1131 
1132  // Update the inflow BCs
1133  this->inflowBoundary().topoChange();
1134 }
1135 
1136 
1137 template<class ParcelType>
1139 {
1141 
1142  // Update the cell occupancy field
1143  cellOccupancy_.setSize(mesh_.nCells());
1144  buildCellOccupancy();
1145 
1146  // Update the inflow BCs
1147  this->inflowBoundary().topoChange();
1148 }
1149 
1150 
1151 // ************************************************************************* //
label k
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:459
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
Templated DSMC particle collision class.
Base cloud calls templated on particle type.
Definition: Cloud.H:74
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
Definition: Cloud.C:457
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
Definition: Cloud.C:634
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Definition: Cloud.C:504
void move(TrackCloudType &cloud, typename ParticleType::trackingData &td)
Move the particles.
Definition: Cloud.C:281
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:80
void addNewParcel(const vector &position, const label celli, label &nLocateBoundaryHits, const vector &U, const scalar Ei, const label typeId)
Add new parcel.
Definition: DSMCCloud.C:465
DSMCCloud(const word &cloudName, const fvMesh &mesh, bool readFields=true)
Construct given name and mesh, will read Parcels and fields from.
Definition: DSMCCloud.C:494
virtual ~DSMCCloud()
Destructor.
Definition: DSMCCloud.C:959
vector equipartitionLinearVelocity(scalar temperature, scalar mass)
Generate a random velocity sampled from the Maxwellian speed.
Definition: DSMCCloud.C:1040
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
Definition: DSMCCloud.C:1110
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
Definition: DSMCCloud.C:1138
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Definition: DSMCCloud.C:1124
void evolve()
Evolve the cloud (move, collide)
Definition: DSMCCloud.C:966
scalar equipartitionInternalEnergy(scalar temperature, direction internalDegreesOfFreedom)
Generate a random internal energy, sampled from the.
Definition: DSMCCloud.C:1053
void clear()
Clear the Cloud.
Definition: DSMCCloudI.H:488
void dumpParticlePositions() const
Dump particle positions to .obj file.
Definition: DSMCCloud.C:1089
void info() const
Print cloud information.
Definition: DSMCCloud.C:996
Generic GeometricField class.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
Templated inflow boundary model class.
Output to file stream.
Definition: OFstream.H:86
virtual void flush()
Flush stream.
Definition: OSstream.C:223
Templated wall interaction model class.
Dimension set for the base types.
Definition: dimensionSet.H:125
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
static List< tetIndices > cellTetIndices(const polyMesh &mesh, label cI)
Return the tet decomposition of the given cell, see.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
scalar scalar01()
Return a scalar uniformly distributed between zero and one.
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
U
Definition: pEqn.H:72
#define WarningInFunction
Report a warning using Foam::Warning.
const dimensionedScalar c
Speed of light in a vacuum.
Collection of constants.
static const zero Zero
Definition: zero.H:97
dimensionedScalar pos(const dimensionedScalar &ds)
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedScalar pos0(const dimensionedScalar &ds)
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
label findMax(const ListType &, const label start=0)
Find index of max element (and larger than given element).
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
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
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const dimensionSet dimless
bool isType(const Type &t)
Check the typeid.
Definition: typeInfo.H:158
messageStream Info
vector point
Point is a vector.
Definition: point.H:41
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar log(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
dimensionedScalar sqrt(const dimensionedScalar &ds)
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
dimensioned< scalar > mag(const dimensioned< Type > &)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Field< vector > vectorField
Specialisation of Field<T> for vector.
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
error FatalError
static const char nl
Definition: Ostream.H:266
uint8_t direction
Definition: direction.H:45
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
tetrahedron< point, const point & > tetPointRef
Definition: tetPointRef.H:44
volScalarField & p
const word cloudName(propsDict.lookup("cloudName"))