KinematicCloud.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 | Copyright (C) 2011-2013 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 "KinematicCloud.H"
27 #include "IntegrationScheme.H"
28 #include "interpolation.H"
29 #include "subCycleTime.H"
30 
31 #include "InjectionModelList.H"
32 #include "DispersionModel.H"
33 #include "PatchInteractionModel.H"
35 #include "SurfaceFilmModel.H"
36 
37 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
38 
39 template<class CloudType>
41 {
42  dispersionModel_.reset
43  (
45  (
46  subModelProperties_,
47  *this
48  ).ptr()
49  );
50 
51  patchInteractionModel_.reset
52  (
54  (
55  subModelProperties_,
56  *this
57  ).ptr()
58  );
59 
60  stochasticCollisionModel_.reset
61  (
63  (
64  subModelProperties_,
65  *this
66  ).ptr()
67  );
68 
69  surfaceFilmModel_.reset
70  (
72  (
73  subModelProperties_,
74  *this
75  ).ptr()
76  );
77 
78  UIntegrator_.reset
79  (
81  (
82  "U",
83  solution_.integrationSchemes()
84  ).ptr()
85  );
86 }
87 
88 
89 template<class CloudType>
90 template<class TrackData>
92 {
93  if (solution_.steadyState())
94  {
95  td.cloud().storeState();
96 
97  td.cloud().preEvolve();
98 
99  evolveCloud(td);
100 
101  if (solution_.coupled())
102  {
103  td.cloud().relaxSources(td.cloud().cloudCopy());
104  }
105  }
106  else
107  {
108  td.cloud().preEvolve();
109 
110  evolveCloud(td);
111 
112  if (solution_.coupled())
113  {
114  td.cloud().scaleSources();
115  }
116  }
117 
118  td.cloud().info();
119 
120  td.cloud().postEvolve();
121 
122  if (solution_.steadyState())
123  {
124  td.cloud().restoreState();
125  }
126 }
127 
128 
129 template<class CloudType>
131 {
132  if (cellOccupancyPtr_.empty())
133  {
134  cellOccupancyPtr_.reset
135  (
136  new List<DynamicList<parcelType*> >(mesh_.nCells())
137  );
138  }
139  else if (cellOccupancyPtr_().size() != mesh_.nCells())
140  {
141  // If the size of the mesh has changed, reset the
142  // cellOccupancy size
143 
144  cellOccupancyPtr_().setSize(mesh_.nCells());
145  }
146 
147  List<DynamicList<parcelType*> >& cellOccupancy = cellOccupancyPtr_();
148 
149  forAll(cellOccupancy, cO)
150  {
151  cellOccupancy[cO].clear();
152  }
153 
154  forAllIter(typename KinematicCloud<CloudType>, *this, iter)
155  {
156  cellOccupancy[iter().cell()].append(&iter());
157  }
158 }
159 
160 
161 template<class CloudType>
163 {
164  // Only build the cellOccupancy if the pointer is set, i.e. it has
165  // been requested before.
166 
167  if (cellOccupancyPtr_.valid())
168  {
169  buildCellOccupancy();
170  }
171 }
172 
173 
174 template<class CloudType>
175 template<class TrackData>
177 {
178  if (solution_.coupled())
179  {
180  td.cloud().resetSourceTerms();
181  }
182 
183  if (solution_.transient())
184  {
185  label preInjectionSize = this->size();
186 
187  this->surfaceFilm().inject(td);
188 
189  // Update the cellOccupancy if the size of the cloud has changed
190  // during the injection.
191  if (preInjectionSize != this->size())
192  {
193  updateCellOccupancy();
194  preInjectionSize = this->size();
195  }
196 
197  injectors_.inject(td);
198 
199 
200  // Assume that motion will update the cellOccupancy as necessary
201  // before it is required.
202  td.cloud().motion(td);
203 
204  stochasticCollision().update(solution_.trackTime());
205  }
206  else
207  {
208 // this->surfaceFilm().injectSteadyState(td);
209 
210  injectors_.injectSteadyState(td, solution_.trackTime());
211 
212  td.part() = TrackData::tpLinearTrack;
213  CloudType::move(td, solution_.trackTime());
214  }
215 }
216 
217 
218 template<class CloudType>
220 {
221  Info<< endl;
222 
223  if (debug)
224  {
225  this->writePositions();
226  }
227 
228  this->dispersion().cacheFields(false);
229 
230  forces_.cacheFields(false);
231 
232  functions_.postEvolve();
233 
234  solution_.nextIter();
235 
236  if (this->db().time().outputTime())
237  {
238  outputProperties_.writeObject
239  (
240  IOstream::ASCII,
241  IOstream::currentVersion,
242  this->db().time().writeCompression()
243  );
244  }
245 }
246 
247 
248 template<class CloudType>
250 {
251  CloudType::cloudReset(c);
252 
253  rndGen_ = c.rndGen_;
254 
255  forces_.transfer(c.forces_);
256 
257  functions_.transfer(c.functions_);
258 
259  injectors_.transfer(c.injectors_);
260 
261  dispersionModel_.reset(c.dispersionModel_.ptr());
262  patchInteractionModel_.reset(c.patchInteractionModel_.ptr());
263  stochasticCollisionModel_.reset(c.stochasticCollisionModel_.ptr());
264  surfaceFilmModel_.reset(c.surfaceFilmModel_.ptr());
265 
266  UIntegrator_.reset(c.UIntegrator_.ptr());
267 }
268 
269 
270 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
271 
272 template<class CloudType>
274 (
275  const word& cloudName,
276  const volScalarField& rho,
277  const volVectorField& U,
278  const volScalarField& mu,
279  const dimensionedVector& g,
280  bool readFields
281 )
282 :
283  CloudType(rho.mesh(), cloudName, false),
284  kinematicCloud(),
285  cloudCopyPtr_(NULL),
286  mesh_(rho.mesh()),
287  particleProperties_
288  (
289  IOobject
290  (
291  cloudName + "Properties",
292  rho.mesh().time().constant(),
293  rho.mesh(),
294  IOobject::MUST_READ_IF_MODIFIED,
295  IOobject::NO_WRITE
296  )
297  ),
298  outputProperties_
299  (
300  IOobject
301  (
302  cloudName + "OutputProperties",
303  mesh_.time().timeName(),
304  "uniform"/cloud::prefix/cloudName,
305  mesh_,
306  IOobject::READ_IF_PRESENT,
307  IOobject::NO_WRITE
308  )
309  ),
310  solution_(mesh_, particleProperties_.subDict("solution")),
311  constProps_(particleProperties_),
312  subModelProperties_
313  (
314  particleProperties_.subOrEmptyDict("subModels", solution_.active())
315  ),
316  rndGen_
317  (
318  label(0),
319  solution_.steadyState() ?
320  particleProperties_.lookupOrDefault<label>("randomSampleSize", 100000)
321  : -1
322  ),
323  cellOccupancyPtr_(),
324  cellLengthScale_(cbrt(mesh_.V())),
325  rho_(rho),
326  U_(U),
327  mu_(mu),
328  g_(g),
329  pAmbient_(0.0),
330  forces_
331  (
332  *this,
333  mesh_,
334  subModelProperties_.subOrEmptyDict
335  (
336  "particleForces",
337  solution_.active()
338  ),
339  solution_.active()
340  ),
341  functions_
342  (
343  *this,
344  particleProperties_.subOrEmptyDict("cloudFunctions"),
345  solution_.active()
346  ),
347  injectors_
348  (
349  subModelProperties_.subOrEmptyDict("injectionModels"),
350  *this
351  ),
352  dispersionModel_(NULL),
353  patchInteractionModel_(NULL),
354  stochasticCollisionModel_(NULL),
355  surfaceFilmModel_(NULL),
356  UIntegrator_(NULL),
357  UTrans_
358  (
360  (
361  IOobject
362  (
363  this->name() + ":UTrans",
364  this->db().time().timeName(),
365  this->db(),
366  IOobject::READ_IF_PRESENT,
367  IOobject::AUTO_WRITE
368  ),
369  mesh_,
370  dimensionedVector("zero", dimMass*dimVelocity, vector::zero)
371  )
372  ),
373  UCoeff_
374  (
376  (
377  IOobject
378  (
379  this->name() + ":UCoeff",
380  this->db().time().timeName(),
381  this->db(),
382  IOobject::READ_IF_PRESENT,
383  IOobject::AUTO_WRITE
384  ),
385  mesh_,
386  dimensionedScalar("zero", dimMass, 0.0)
387  )
388  )
389 {
390  if (solution_.active())
391  {
392  setModels();
393 
394  if (readFields)
395  {
396  parcelType::readFields(*this);
397  }
398  }
399 
400  if (solution_.resetSourcesOnStartup())
401  {
402  resetSourceTerms();
403  }
404 }
405 
406 
407 template<class CloudType>
409 (
411  const word& name
412 )
413 :
414  CloudType(c.mesh_, name, c),
415  kinematicCloud(),
416  cloudCopyPtr_(NULL),
417  mesh_(c.mesh_),
418  particleProperties_(c.particleProperties_),
419  outputProperties_(c.outputProperties_),
420  solution_(c.solution_),
421  constProps_(c.constProps_),
422  subModelProperties_(c.subModelProperties_),
423  rndGen_(c.rndGen_, true),
424  cellOccupancyPtr_(NULL),
425  cellLengthScale_(c.cellLengthScale_),
426  rho_(c.rho_),
427  U_(c.U_),
428  mu_(c.mu_),
429  g_(c.g_),
430  pAmbient_(c.pAmbient_),
431  forces_(c.forces_),
432  functions_(c.functions_),
433  injectors_(c.injectors_),
434  dispersionModel_(c.dispersionModel_->clone()),
435  patchInteractionModel_(c.patchInteractionModel_->clone()),
436  stochasticCollisionModel_(c.stochasticCollisionModel_->clone()),
437  surfaceFilmModel_(c.surfaceFilmModel_->clone()),
438  UIntegrator_(c.UIntegrator_->clone()),
439  UTrans_
440  (
442  (
443  IOobject
444  (
445  this->name() + ":UTrans",
446  this->db().time().timeName(),
447  this->db(),
448  IOobject::NO_READ,
449  IOobject::NO_WRITE,
450  false
451  ),
452  c.UTrans_()
453  )
454  ),
455  UCoeff_
456  (
458  (
459  IOobject
460  (
461  name + ":UCoeff",
462  this->db().time().timeName(),
463  this->db(),
464  IOobject::NO_READ,
465  IOobject::NO_WRITE,
466  false
467  ),
468  c.UCoeff_()
469  )
470  )
471 {}
472 
473 
474 template<class CloudType>
476 (
477  const fvMesh& mesh,
478  const word& name,
480 )
481 :
483  kinematicCloud(),
484  cloudCopyPtr_(NULL),
485  mesh_(mesh),
486  particleProperties_
487  (
488  IOobject
489  (
490  name + "Properties",
491  mesh.time().constant(),
492  mesh,
493  IOobject::NO_READ,
494  IOobject::NO_WRITE,
495  false
496  )
497  ),
498  outputProperties_
499  (
500  IOobject
501  (
502  name + "OutputProperties",
503  mesh_.time().timeName(),
504  "uniform"/cloud::prefix/name,
505  mesh_,
506  IOobject::NO_READ,
507  IOobject::NO_WRITE,
508  false
509  )
510  ),
511  solution_(mesh),
512  constProps_(),
513  subModelProperties_(dictionary::null),
514  rndGen_(0, 0),
515  cellOccupancyPtr_(NULL),
516  cellLengthScale_(c.cellLengthScale_),
517  rho_(c.rho_),
518  U_(c.U_),
519  mu_(c.mu_),
520  g_(c.g_),
521  pAmbient_(c.pAmbient_),
522  forces_(*this, mesh),
523  functions_(*this),
524  injectors_(*this),
525  dispersionModel_(NULL),
526  patchInteractionModel_(NULL),
527  stochasticCollisionModel_(NULL),
528  surfaceFilmModel_(NULL),
529  UIntegrator_(NULL),
530  UTrans_(NULL),
531  UCoeff_(NULL)
532 {}
533 
534 
535 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
536 
537 template<class CloudType>
539 {}
540 
541 
542 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
543 
544 template<class CloudType>
546 {
547  return true;
548 }
549 
550 
551 template<class CloudType>
553 (
554  parcelType& parcel,
555  const scalar lagrangianDt
556 )
557 {
558  parcel.rho() = constProps_.rho0();
559 }
560 
561 
562 template<class CloudType>
564 (
565  parcelType& parcel,
566  const scalar lagrangianDt,
567  const bool fullyDescribed
568 )
569 {
570  const scalar carrierDt = mesh_.time().deltaTValue();
571  parcel.stepFraction() = (carrierDt - lagrangianDt)/carrierDt;
572 
573  if (parcel.typeId() == -1)
574  {
575  parcel.typeId() = constProps_.parcelTypeId();
576  }
577 }
578 
579 
580 template<class CloudType>
582 {
583  cloudCopyPtr_.reset
584  (
585  static_cast<KinematicCloud<CloudType>*>
586  (
587  clone(this->name() + "Copy").ptr()
588  )
589  );
590 }
591 
592 
593 template<class CloudType>
595 {
596  cloudReset(cloudCopyPtr_());
597  cloudCopyPtr_.clear();
598 }
599 
600 
601 template<class CloudType>
603 {
604  UTrans().field() = vector::zero;
605  UCoeff().field() = 0.0;
606 }
607 
608 
609 template<class CloudType>
610 template<class Type>
612 (
614  const DimensionedField<Type, volMesh>& field0,
615  const word& name
616 ) const
617 {
618  const scalar coeff = solution_.relaxCoeff(name);
619  field = field0 + coeff*(field - field0);
620 }
621 
622 
623 template<class CloudType>
624 template<class Type>
626 (
628  const word& name
629 ) const
630 {
631  const scalar coeff = solution_.relaxCoeff(name);
632  field *= coeff;
633 }
634 
635 
636 template<class CloudType>
638 (
639  const KinematicCloud<CloudType>& cloudOldTime
640 )
641 {
642  this->relax(UTrans_(), cloudOldTime.UTrans(), "U");
643  this->relax(UCoeff_(), cloudOldTime.UCoeff(), "U");
644 }
645 
646 
647 template<class CloudType>
649 {
650  this->scale(UTrans_(), "U");
651  this->scale(UCoeff_(), "U");
652 }
653 
654 
655 template<class CloudType>
657 {
658  // force calculaion of mesh dimensions - needed for parallel runs
659  // with topology change due to lazy evaluation of valid mesh dimensions
660  label nGeometricD = mesh_.nGeometricD();
661 
662  Info<< "\nSolving " << nGeometricD << "-D cloud " << this->name() << endl;
663 
664  this->dispersion().cacheFields(true);
665  forces_.cacheFields(true);
666  updateCellOccupancy();
667 
668  pAmbient_ = constProps_.dict().template
669  lookupOrDefault<scalar>("pAmbient", pAmbient_);
670 
671  functions_.preEvolve();
672 }
673 
674 
675 template<class CloudType>
677 {
678  if (solution_.canEvolve())
679  {
680  typename parcelType::template
681  TrackingData<KinematicCloud<CloudType> > td(*this);
682 
683  solve(td);
684  }
685 }
686 
687 
688 template<class CloudType>
689 template<class TrackData>
691 {
692  td.part() = TrackData::tpLinearTrack;
693  CloudType::move(td, solution_.trackTime());
694 
695  updateCellOccupancy();
696 }
697 
698 
699 template<class CloudType>
701 (
702  const parcelType& p,
703  const polyPatch& pp,
704  const scalar trackFraction,
705  const tetIndices& tetIs,
706  vector& nw,
707  vector& Up
708 ) const
709 {
710  label patchI = pp.index();
711  label patchFaceI = pp.whichFace(p.face());
712 
713  vector n = tetIs.faceTri(mesh_).normal();
714  n /= mag(n);
715 
716  vector U = U_.boundaryField()[patchI][patchFaceI];
717 
718  // Unless the face is rotating, the required normal is n;
719  nw = n;
720 
721  if (!mesh_.moving())
722  {
723  // Only wall patches may have a non-zero wall velocity from
724  // the velocity field when the mesh is not moving.
725 
726  if (isA<wallPolyPatch>(pp))
727  {
728  Up = U;
729  }
730  else
731  {
732  Up = vector::zero;
733  }
734  }
735  else
736  {
737  vector U00 = U_.oldTime().boundaryField()[patchI][patchFaceI];
738 
739  vector n00 = tetIs.oldFaceTri(mesh_).normal();
740 
741  // Difference in normal over timestep
742  vector dn = vector::zero;
743 
744  if (mag(n00) > SMALL)
745  {
746  // If the old normal is zero (for example in layer
747  // addition) then use the current normal, meaning that the
748  // motion can only be translational, and dn remains zero,
749  // otherwise, calculate dn:
750 
751  n00 /= mag(n00);
752 
753  dn = n - n00;
754  }
755 
756  // Total fraction through the timestep of the motion,
757  // including stepFraction before the current tracking step
758  // and the current trackFraction
759  // i.e.
760  // let s = stepFraction, t = trackFraction
761  // Motion of x in time:
762  // |-----------------|---------|---------|
763  // x00 x0 xi x
764  //
765  // where xi is the correct value of x at the required
766  // tracking instant.
767  //
768  // x0 = x00 + s*(x - x00) = s*x + (1 - s)*x00
769  //
770  // i.e. the motion covered by previous tracking portions
771  // within this timestep, and
772  //
773  // xi = x0 + t*(x - x0)
774  // = t*x + (1 - t)*x0
775  // = t*x + (1 - t)*(s*x + (1 - s)*x00)
776  // = (s + t - s*t)*x + (1 - (s + t - s*t))*x00
777  //
778  // let m = (s + t - s*t)
779  //
780  // xi = m*x + (1 - m)*x00 = x00 + m*(x - x00);
781  //
782  // In the same form as before.
783 
784  scalar m =
785  p.stepFraction()
786  + trackFraction
787  - (p.stepFraction()*trackFraction);
788 
789  // When the mesh is moving, the velocity field on wall patches
790  // will contain the velocity associated with the motion of the
791  // mesh, in which case it is interpolated in time using m.
792  // For other patches the face velocity will need to be
793  // reconstructed from the face centre motion.
794 
795  const vector& Cf = mesh_.faceCentres()[p.face()];
796 
797  vector Cf00 = mesh_.faces()[p.face()].centre(mesh_.oldPoints());
798 
799  if (isA<wallPolyPatch>(pp))
800  {
801  Up = U00 + m*(U - U00);
802  }
803  else
804  {
805  Up = (Cf - Cf00)/mesh_.time().deltaTValue();
806  }
807 
808  if (mag(dn) > SMALL)
809  {
810  // Rotational motion, nw requires interpolation and a
811  // rotational velocity around face centre correction to Up
812  // is required.
813 
814  nw = n00 + m*dn;
815 
816  // Cf at tracking instant
817  vector Cfi = Cf00 + m*(Cf - Cf00);
818 
819  // Normal vector cross product
820  vector omega = (n00 ^ n);
821 
822  scalar magOmega = mag(omega);
823 
824  // magOmega = sin(angle between unit normals)
825  // Normalise omega vector by magOmega, then multiply by
826  // angle/dt to give the correct angular velocity vector.
827  omega *= Foam::asin(magOmega)/(magOmega*mesh_.time().deltaTValue());
828 
829  // Project position onto face and calculate this position
830  // relative to the face centre.
831  vector facePos =
832  p.position()
833  - ((p.position() - Cfi) & nw)*nw
834  - Cfi;
835 
836  Up += (omega ^ facePos);
837  }
838 
839  // No further action is required if the motion is
840  // translational only, nw and Up have already been set.
841  }
842 }
843 
844 
845 template<class CloudType>
847 {
848  updateCellOccupancy();
849  injectors_.updateMesh();
850  cellLengthScale_ = cbrt(mesh_.V());
851 }
852 
853 
854 template<class CloudType>
856 {
857  typedef typename particle::TrackingData<KinematicCloud<CloudType> > tdType;
858 
859  tdType td(*this);
860 
861  Cloud<parcelType>::template autoMap<tdType>(td, mapper);
862 
863  updateMesh();
864 }
865 
866 
867 template<class CloudType>
869 {
870  vector linearMomentum = linearMomentumOfSystem();
871  reduce(linearMomentum, sumOp<vector>());
872 
873  scalar linearKineticEnergy = linearKineticEnergyOfSystem();
874  reduce(linearKineticEnergy, sumOp<scalar>());
875 
876  Info<< "Cloud: " << this->name() << nl
877  << " Current number of parcels = "
878  << returnReduce(this->size(), sumOp<label>()) << nl
879  << " Current mass in system = "
880  << returnReduce(massInSystem(), sumOp<scalar>()) << nl
881  << " Linear momentum = "
882  << linearMomentum << nl
883  << " |Linear momentum| = "
884  << mag(linearMomentum) << nl
885  << " Linear kinetic energy = "
886  << linearKineticEnergy << nl;
887 
888  injectors_.info(Info);
889  this->surfaceFilm().info(Info);
890  this->patchInteraction().info(Info);
891 }
892 
893 
894 // ************************************************************************* //
Foam::KinematicCloud::UCoeff
DimensionedField< scalar, volMesh > & UCoeff()
Return coefficient for carrier phase U equation.
Definition: KinematicCloudI.H:540
PatchInteractionModel.H
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::KinematicCloud::scaleSources
void scaleSources()
Apply scaling to (transient) cloud sources.
Definition: KinematicCloud.C:648
relax
p relax()
p
p
Definition: pEqn.H:62
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::compressible::New
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Definition: turbulentFluidThermoModel.C:36
forAllIter
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:431
Foam::KinematicCloud::setModels
void setModels()
Set cloud sub-models.
Definition: KinematicCloud.C:40
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:86
Foam::constant::physicoChemical::mu
const dimensionedScalar mu
Atomic mass unit.
Definition: createFields.H:13
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::DispersionModel
Definition: KinematicCloud.H:79
Foam::dimVelocity
const dimensionSet dimVelocity
Foam::KinematicCloud::patchData
void patchData(const parcelType &p, const polyPatch &pp, const scalar trackFraction, const tetIndices &tetIs, vector &normal, vector &Up) const
Calculate the patch normal and velocity to interact with,.
Definition: KinematicCloud.C:701
Foam::DynamicList
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:56
Foam::KinematicCloud::relaxSources
void relaxSources(const KinematicCloud< CloudType > &cloudOldTime)
Apply relaxation to (steady state) cloud sources.
Definition: KinematicCloud.C:638
Foam::KinematicCloud::storeState
void storeState()
Store the current cloud state.
Definition: KinematicCloud.C:581
Foam::tetIndices::oldFaceTri
triPointRef oldFaceTri(const polyMesh &mesh) const
Return the geometry corresponding to the tri on the.
Definition: tetIndicesI.H:140
nw
label nw
Definition: createFields.H:25
g
const dimensionedVector & g
Definition: setRegionFluidFields.H:33
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
interpolation.H
Foam::KinematicCloud::solve
void solve(TrackData &td)
Solve the cloud - calls all evolution functions.
Definition: KinematicCloud.C:91
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::readFields
This function object reads fields from the time directories and adds them to the mesh database for fu...
Definition: readFields.H:104
U
U
Definition: pEqn.H:46
Foam::dimensionedVector
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Definition: dimensionedVector.H:48
Foam::KinematicCloud::motion
void motion(TrackData &td)
Particle motion.
Definition: KinematicCloud.C:690
n
label n
Definition: TABSMDCalcMethod2.H:31
solve
rhoEqn solve()
Foam::kinematicCloud
Virtual abstract base class for templated KinematicCloud.
Definition: kinematicCloud.H:49
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:43
Foam::KinematicCloud::evolve
void evolve()
Evolve the cloud.
Definition: KinematicCloud.C:676
Foam::KinematicCloud
Templated base class for kinematic cloud.
Definition: KinematicCloud.H:96
Foam::PatchInteractionModel
Templated patch interaction model class.
Definition: KinematicCloud.H:82
Foam::KinematicCloud::hasWallImpactDistance
virtual bool hasWallImpactDistance() const
Switch to specify if particles of the cloud can return.
Definition: KinematicCloud.C:545
Foam::StochasticCollisionModel
Templated stochastic collision model class.
Definition: KinematicCloud.H:88
Foam::label
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
Foam::IDLList
Intrusive doubly-linked list.
Definition: IDLList.H:47
Foam::triangle::normal
vector normal() const
Return vector normal.
Definition: triangleI.H:116
cellOccupancy
const List< DynamicList< molecule * > > & cellOccupancy
Definition: calculateMDFields.H:1
Foam::KinematicCloud::~KinematicCloud
virtual ~KinematicCloud()
Destructor.
Definition: KinematicCloud.C:538
Foam::tetIndices::faceTri
triPointRef faceTri(const polyMesh &mesh) const
Return the geometry corresponding to the tri on the.
Definition: tetIndicesI.H:109
Foam::KinematicCloud::checkParcelProperties
void checkParcelProperties(parcelType &parcel, const scalar lagrangianDt, const bool fullyDescribed)
Check parcel properties.
Definition: KinematicCloud.C:564
Foam::SurfaceFilmModel
Templated wall surface film model class.
Definition: KinematicCloud.H:85
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Foam::KinematicCloud::UTrans
DimensionedField< vector, volMesh > & UTrans()
Return reference to momentum source.
Definition: KinematicCloudI.H:524
subCycleTime.H
Foam::Vector::centre
const Vector< Cmpt > & centre(const Foam::List< Vector< Cmpt > > &) const
Return *this (used for point which is a typedef to Vector<scalar>.
Definition: VectorI.H:106
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
Foam::CloudType
DSMCCloud< dsmcParcel > CloudType
Definition: makeDSMCParcelBinaryCollisionModels.C:36
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
KinematicCloud.H
Foam::KinematicCloud::evolveCloud
void evolveCloud(TrackData &td)
Evolve the cloud.
Definition: KinematicCloud.C:176
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
coeff
const scalar coeff[]
Definition: Test-Polynomial.C:38
Foam::KinematicCloud::autoMap
virtual void autoMap(const mapPolyMesh &)
Remap the cells of particles corresponding to the.
Definition: KinematicCloud.C:855
Foam::KinematicCloud::setParcelThermoProperties
void setParcelThermoProperties(parcelType &parcel, const scalar lagrangianDt)
Set parcel thermo properties.
Definition: KinematicCloud.C:553
Foam::KinematicCloud::resetSourceTerms
void resetSourceTerms()
Reset the cloud source terms.
Definition: KinematicCloud.C:602
rho
rho
Definition: pEqn.H:3
DispersionModel.H
Foam::tetIndices
Storage and named access for the indices of a tet which is part of the decomposition of a cell.
Definition: tetIndices.H:73
Foam::polyPatch::whichFace
label whichFace(const label l) const
Return label of face in patch from global face label.
Definition: polyPatch.H:389
Foam::sumOp
Definition: ops.H:162
Foam::Vector< scalar >
Foam::KinematicCloud::cloudReset
void cloudReset(KinematicCloud< CloudType > &c)
Reset state of cloud.
Definition: KinematicCloud.C:249
Foam::KinematicCloud::postEvolve
void postEvolve()
Post-evolve.
Definition: KinematicCloud.C:219
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::KinematicCloud::updateMesh
void updateMesh()
Update mesh.
Definition: KinematicCloud.C:846
readFields
void readFields(const boolList &haveMesh, const fvMesh &mesh, const autoPtr< fvMeshSubset > &subsetterPtr, IOobjectList &allObjects, PtrList< GeoField > &fields)
Definition: redistributePar.C:589
Foam::KinematicCloud::relax
void relax(DimensionedField< Type, volMesh > &field, const DimensionedField< Type, volMesh > &field0, const word &name) const
Relax field.
Definition: KinematicCloud.C:612
Foam::Cloud
Base cloud calls templated on particle type.
Definition: Cloud.H:52
surfaceFilm
filmModelType & surfaceFilm
Definition: createSurfaceFilmModel.H:6
timeName
word timeName
Definition: getTimeIndex.H:3
StochasticCollisionModel.H
Foam::KinematicCloud::KinematicCloud
KinematicCloud(const KinematicCloud &)
Disallow default bitwise copy construct.
cloudName
const word cloudName(propsDict.lookup("cloudName"))
IntegrationScheme.H
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::KinematicCloud::restoreState
void restoreState()
Reset the current cloud to the previously stored state.
Definition: KinematicCloud.C:594
Foam::DSMCCloud::parcelType
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:217
SurfaceFilmModel.H
Foam::KinematicCloud::updateCellOccupancy
void updateCellOccupancy()
Update (i.e. build) the cellOccupancy if it has.
Definition: KinematicCloud.C:162
Foam::KinematicCloud::info
void info()
Print cloud information.
Definition: KinematicCloud.C:868
Foam::KinematicCloud::buildCellOccupancy
void buildCellOccupancy()
Build the cellOccupancy.
Definition: KinematicCloud.C:130
InjectionModelList.H
Foam::cbrt
dimensionedScalar cbrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:153
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::patchIdentifier::index
label index() const
Return the index of this patch in the boundaryMesh.
Definition: patchIdentifier.H:133
Foam::particle::TrackingData
Definition: particle.H:94
Foam::asin
dimensionedScalar asin(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:258
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Foam::KinematicCloud::preEvolve
void preEvolve()
Pre-evolve.
Definition: KinematicCloud.C:656
Foam::KinematicCloud::scale
void scale(DimensionedField< Type, volMesh > &field, const word &name) const
Scale field.
Definition: KinematicCloud.C:626
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51