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-2017 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().writeTime())
237  {
238  outputProperties_.writeObject
239  (
240  IOstream::ASCII,
241  IOstream::currentVersion,
242  this->db().time().writeCompression(),
243  true
244  );
245  }
246 }
247 
248 
249 template<class CloudType>
251 {
252  CloudType::cloudReset(c);
253 
254  rndGen_ = c.rndGen_;
255 
256  forces_.transfer(c.forces_);
257 
258  functions_.transfer(c.functions_);
259 
260  injectors_.transfer(c.injectors_);
261 
262  dispersionModel_.reset(c.dispersionModel_.ptr());
263  patchInteractionModel_.reset(c.patchInteractionModel_.ptr());
264  stochasticCollisionModel_.reset(c.stochasticCollisionModel_.ptr());
265  surfaceFilmModel_.reset(c.surfaceFilmModel_.ptr());
266 
267  UIntegrator_.reset(c.UIntegrator_.ptr());
268 }
269 
270 
271 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
272 
273 template<class CloudType>
275 (
276  const word& cloudName,
277  const volScalarField& rho,
278  const volVectorField& U,
279  const volScalarField& mu,
280  const dimensionedVector& g,
281  bool readFields
282 )
283 :
284  CloudType(rho.mesh(), cloudName, false),
285  kinematicCloud(),
286  cloudCopyPtr_(nullptr),
287  mesh_(rho.mesh()),
288  particleProperties_
289  (
290  IOobject
291  (
292  cloudName + "Properties",
293  rho.mesh().time().constant(),
294  rho.mesh(),
295  IOobject::MUST_READ_IF_MODIFIED,
296  IOobject::NO_WRITE
297  )
298  ),
299  outputProperties_
300  (
301  IOobject
302  (
303  cloudName + "OutputProperties",
304  mesh_.time().timeName(),
305  "uniform"/cloud::prefix/cloudName,
306  mesh_,
307  IOobject::READ_IF_PRESENT,
308  IOobject::NO_WRITE
309  )
310  ),
311  solution_(mesh_, particleProperties_.subDict("solution")),
312  constProps_(particleProperties_),
313  subModelProperties_
314  (
315  particleProperties_.subOrEmptyDict("subModels", solution_.active())
316  ),
317  rndGen_
318  (
319  label(0),
320  solution_.steadyState() ?
321  particleProperties_.lookupOrDefault<label>("randomSampleSize", 100000)
322  : -1
323  ),
324  cellOccupancyPtr_(),
325  cellLengthScale_(mag(cbrt(mesh_.V()))),
326  rho_(rho),
327  U_(U),
328  mu_(mu),
329  g_(g),
330  pAmbient_(0.0),
331  forces_
332  (
333  *this,
334  mesh_,
335  subModelProperties_.subOrEmptyDict
336  (
337  "particleForces",
338  solution_.active()
339  ),
340  solution_.active()
341  ),
342  functions_
343  (
344  *this,
345  particleProperties_.subOrEmptyDict("cloudFunctions"),
346  solution_.active()
347  ),
348  injectors_
349  (
350  subModelProperties_.subOrEmptyDict("injectionModels"),
351  *this
352  ),
353  dispersionModel_(nullptr),
354  patchInteractionModel_(nullptr),
355  stochasticCollisionModel_(nullptr),
356  surfaceFilmModel_(nullptr),
357  UIntegrator_(nullptr),
358  UTrans_
359  (
361  (
362  IOobject
363  (
364  this->name() + ":UTrans",
365  this->db().time().timeName(),
366  this->db(),
367  IOobject::READ_IF_PRESENT,
368  IOobject::AUTO_WRITE
369  ),
370  mesh_,
372  )
373  ),
374  UCoeff_
375  (
377  (
378  IOobject
379  (
380  this->name() + ":UCoeff",
381  this->db().time().timeName(),
382  this->db(),
383  IOobject::READ_IF_PRESENT,
384  IOobject::AUTO_WRITE
385  ),
386  mesh_,
387  dimensionedScalar("zero", dimMass, 0.0)
388  )
389  )
390 {
391  if (solution_.active())
392  {
393  setModels();
394 
395  if (readFields)
396  {
397  parcelType::readFields(*this);
398  this->deleteLostParticles();
399  }
400  }
401 
402  if (solution_.resetSourcesOnStartup())
403  {
404  resetSourceTerms();
405  }
406 }
407 
408 
409 template<class CloudType>
411 (
413  const word& name
414 )
415 :
416  CloudType(c.mesh_, name, c),
417  kinematicCloud(),
418  cloudCopyPtr_(nullptr),
419  mesh_(c.mesh_),
420  particleProperties_(c.particleProperties_),
421  outputProperties_(c.outputProperties_),
422  solution_(c.solution_),
423  constProps_(c.constProps_),
424  subModelProperties_(c.subModelProperties_),
425  rndGen_(c.rndGen_, true),
426  cellOccupancyPtr_(nullptr),
427  cellLengthScale_(c.cellLengthScale_),
428  rho_(c.rho_),
429  U_(c.U_),
430  mu_(c.mu_),
431  g_(c.g_),
432  pAmbient_(c.pAmbient_),
433  forces_(c.forces_),
434  functions_(c.functions_),
435  injectors_(c.injectors_),
436  dispersionModel_(c.dispersionModel_->clone()),
437  patchInteractionModel_(c.patchInteractionModel_->clone()),
438  stochasticCollisionModel_(c.stochasticCollisionModel_->clone()),
439  surfaceFilmModel_(c.surfaceFilmModel_->clone()),
440  UIntegrator_(c.UIntegrator_->clone()),
441  UTrans_
442  (
444  (
445  IOobject
446  (
447  this->name() + ":UTrans",
448  this->db().time().timeName(),
449  this->db(),
450  IOobject::NO_READ,
451  IOobject::NO_WRITE,
452  false
453  ),
454  c.UTrans_()
455  )
456  ),
457  UCoeff_
458  (
460  (
461  IOobject
462  (
463  name + ":UCoeff",
464  this->db().time().timeName(),
465  this->db(),
466  IOobject::NO_READ,
467  IOobject::NO_WRITE,
468  false
469  ),
470  c.UCoeff_()
471  )
472  )
473 {}
474 
475 
476 template<class CloudType>
478 (
479  const fvMesh& mesh,
480  const word& name,
482 )
483 :
484  CloudType(mesh, name, IDLList<parcelType>()),
485  kinematicCloud(),
486  cloudCopyPtr_(nullptr),
487  mesh_(mesh),
488  particleProperties_
489  (
490  IOobject
491  (
492  name + "Properties",
493  mesh.time().constant(),
494  mesh,
495  IOobject::NO_READ,
496  IOobject::NO_WRITE,
497  false
498  )
499  ),
500  outputProperties_
501  (
502  IOobject
503  (
504  name + "OutputProperties",
505  mesh_.time().timeName(),
506  "uniform"/cloud::prefix/name,
507  mesh_,
508  IOobject::NO_READ,
509  IOobject::NO_WRITE,
510  false
511  )
512  ),
513  solution_(mesh),
514  constProps_(),
515  subModelProperties_(dictionary::null),
516  rndGen_(0, 0),
517  cellOccupancyPtr_(nullptr),
518  cellLengthScale_(c.cellLengthScale_),
519  rho_(c.rho_),
520  U_(c.U_),
521  mu_(c.mu_),
522  g_(c.g_),
523  pAmbient_(c.pAmbient_),
524  forces_(*this, mesh),
525  functions_(*this),
526  injectors_(*this),
527  dispersionModel_(nullptr),
528  patchInteractionModel_(nullptr),
529  stochasticCollisionModel_(nullptr),
530  surfaceFilmModel_(nullptr),
531  UIntegrator_(nullptr),
532  UTrans_(nullptr),
533  UCoeff_(nullptr)
534 {}
535 
536 
537 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
538 
539 template<class CloudType>
541 {}
542 
543 
544 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
545 
546 template<class CloudType>
548 {
549  return true;
550 }
551 
552 
553 template<class CloudType>
555 (
556  parcelType& parcel,
557  const scalar lagrangianDt
558 )
559 {
560  parcel.rho() = constProps_.rho0();
561 }
562 
563 
564 template<class CloudType>
566 (
567  parcelType& parcel,
568  const scalar lagrangianDt,
569  const bool fullyDescribed
570 )
571 {
572  const scalar carrierDt = mesh_.time().deltaTValue();
573  parcel.stepFraction() = (carrierDt - lagrangianDt)/carrierDt;
574 
575  if (parcel.typeId() == -1)
576  {
577  parcel.typeId() = constProps_.parcelTypeId();
578  }
579 }
580 
581 
582 template<class CloudType>
584 {
585  cloudCopyPtr_.reset
586  (
587  static_cast<KinematicCloud<CloudType>*>
588  (
589  clone(this->name() + "Copy").ptr()
590  )
591  );
592 }
593 
594 
595 template<class CloudType>
597 {
598  cloudReset(cloudCopyPtr_());
599  cloudCopyPtr_.clear();
600 }
601 
602 
603 template<class CloudType>
605 {
606  UTrans().field() = Zero;
607  UCoeff().field() = 0.0;
608 }
609 
610 
611 template<class CloudType>
612 template<class Type>
614 (
616  const DimensionedField<Type, volMesh>& field0,
617  const word& name
618 ) const
619 {
620  const scalar coeff = solution_.relaxCoeff(name);
621  field = field0 + coeff*(field - field0);
622 }
623 
624 
625 template<class CloudType>
626 template<class Type>
628 (
630  const word& name
631 ) const
632 {
633  const scalar coeff = solution_.relaxCoeff(name);
634  field *= coeff;
635 }
636 
637 
638 template<class CloudType>
640 (
641  const KinematicCloud<CloudType>& cloudOldTime
642 )
643 {
644  this->relax(UTrans_(), cloudOldTime.UTrans(), "U");
645  this->relax(UCoeff_(), cloudOldTime.UCoeff(), "U");
646 }
647 
648 
649 template<class CloudType>
651 {
652  this->scale(UTrans_(), "U");
653  this->scale(UCoeff_(), "U");
654 }
655 
656 
657 template<class CloudType>
659 {
660  // force calculaion of mesh dimensions - needed for parallel runs
661  // with topology change due to lazy evaluation of valid mesh dimensions
662  label nGeometricD = mesh_.nGeometricD();
663 
664  Info<< "\nSolving " << nGeometricD << "-D cloud " << this->name() << endl;
665 
666  this->dispersion().cacheFields(true);
667  forces_.cacheFields(true);
668  updateCellOccupancy();
669 
670  pAmbient_ = constProps_.dict().template
671  lookupOrDefault<scalar>("pAmbient", pAmbient_);
672 
673  functions_.preEvolve();
674 }
675 
676 
677 template<class CloudType>
679 {
680  if (solution_.canEvolve())
681  {
682  typename parcelType::template
683  TrackingData<KinematicCloud<CloudType>> td(*this);
684 
685  solve(td);
686  }
687 }
688 
689 
690 template<class CloudType>
691 template<class TrackData>
693 {
694  td.part() = TrackData::tpLinearTrack;
695  CloudType::move(td, solution_.trackTime());
696 
697  updateCellOccupancy();
698 }
699 
700 
701 template<class CloudType>
703 (
704  const parcelType& p,
705  const polyPatch& pp,
706  vector& nw,
707  vector& Up
708 ) const
709 {
710  p.patchData(nw, Up);
711 
712  // If this is a wall patch, then there may be a non-zero tangential velocity
713  // component; the lid velocity in a lid-driven cavity case, for example. We
714  // want the particle to interact with this velocity, so we look it up in the
715  // velocity field and use it to set the wall-tangential component.
716  if (isA<wallPolyPatch>(pp))
717  {
718  const label patchi = pp.index();
719  const label patchFacei = pp.whichFace(p.face());
720 
721  // We only want to use the boundary condition value onlyif it is set
722  // by the boundary condition. If the boundary values are extrapolated
723  // (e.g., slip conditions) then they represent the motion of the fluid
724  // just inside the domain rather than that of the wall itself.
725  if (U_.boundaryField()[patchi].fixesValue())
726  {
727  const vector Uw1 = U_.boundaryField()[patchi][patchFacei];
728  const vector& Uw0 =
729  U_.oldTime().boundaryField()[patchi][patchFacei];
730 
731  const scalar f = p.currentTimeFraction();
732 
733  const vector Uw = Uw0 + f*(Uw1 - Uw0);
734 
735  const tensor nnw = nw*nw;
736 
737  Up = (nnw & Up) + Uw - (nnw & Uw);
738  }
739  }
740 }
741 
742 
743 template<class CloudType>
745 {
746  updateCellOccupancy();
747  injectors_.updateMesh();
748  cellLengthScale_ = mag(cbrt(mesh_.V()));
749 }
750 
751 
752 template<class CloudType>
754 {
756 
757  updateMesh();
758 }
759 
760 
761 template<class CloudType>
763 {
764  vector linearMomentum = linearMomentumOfSystem();
765  reduce(linearMomentum, sumOp<vector>());
766 
767  scalar linearKineticEnergy = linearKineticEnergyOfSystem();
768  reduce(linearKineticEnergy, sumOp<scalar>());
769 
770  Info<< "Cloud: " << this->name() << nl
771  << " Current number of parcels = "
772  << returnReduce(this->size(), sumOp<label>()) << nl
773  << " Current mass in system = "
774  << returnReduce(massInSystem(), sumOp<scalar>()) << nl
775  << " Linear momentum = "
776  << linearMomentum << nl
777  << " |Linear momentum| = "
778  << mag(linearMomentum) << nl
779  << " Linear kinetic energy = "
780  << linearKineticEnergy << nl;
781 
782  injectors_.info(Info);
783  this->surfaceFilm().info(Info);
784  this->patchInteraction().info(Info);
785 }
786 
787 
788 // ************************************************************************* //
Template class for intrusive linked lists.
Definition: ILList.H:50
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
DSMCCloud< dsmcParcel > CloudType
void scaleSources()
Apply scaling to (transient) cloud sources.
void setModels()
Set cloud sub-models.
const volVectorField & U_
Velocity [m/s].
void relaxSources(const KinematicCloud< CloudType > &cloudOldTime)
Apply relaxation to (steady state) cloud sources.
void storeState()
Store the current cloud state.
UEqn relax()
autoPtr< PatchInteractionModel< KinematicCloud< CloudType > > > patchInteractionModel_
Patch interaction model.
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
InjectionModelList< KinematicCloud< CloudType > > injectors_
Injector models.
Cloud< basicKinematicCollidingParcel > ::particleType parcelType
Type of parcel the cloud was instantiated for.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:60
virtual ~KinematicCloud()
Destructor.
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:453
autoPtr< volVectorField::Internal > UTrans_
Momentum.
void setParcelThermoProperties(parcelType &parcel, const scalar lagrangianDt)
Set parcel thermo properties.
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
void relax(DimensionedField< Type, volMesh > &field, const DimensionedField< Type, volMesh > &field0, const word &name) const
Relax field.
volVectorField::Internal & UTrans()
Return reference to momentum source.
void postEvolve()
Post-evolve.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
autoPtr< SurfaceFilmModel< KinematicCloud< CloudType > > > surfaceFilmModel_
Surface film model.
void updateCellOccupancy()
Update (i.e. build) the cellOccupancy if it has.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
void buildCellOccupancy()
Build the cellOccupancy.
const volScalarField & rho_
Density [kg/m3].
Templated patch interaction model class.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
void updateMesh()
Update mesh.
const volScalarField & mu_
Dynamic viscosity [Pa.s].
dynamicFvMesh & mesh
volScalarField::Internal & UCoeff()
Return coefficient for carrier phase U equation.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
void resetSourceTerms()
Reset the cloud source terms.
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:124
A class for handling words, derived from string.
Definition: word.H:59
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:184
rhoEqn solve()
dimensionedScalar cbrt(const dimensionedScalar &ds)
scalar pAmbient_
Averaged ambient domain pressure.
const word & constant() const
Return constant name.
Definition: TimePaths.H:124
parcelType::constantProperties constProps_
Parcel constant properties.
const fvMesh & mesh_
References to the mesh and time databases.
word timeName
Definition: getTimeIndex.H:3
static const zero Zero
Definition: zero.H:91
void checkParcelProperties(parcelType &parcel, const scalar lagrangianDt, const bool fullyDescribed)
Check parcel properties.
const dictionary subModelProperties_
Sub-models dictionary.
Base cloud calls templated on particle type.
Definition: Cloud.H:52
bool hasWallImpactDistance() const
Switch to specify if particles of the cloud can return.
static const char nl
Definition: Ostream.H:262
const Mesh & mesh() const
Return mesh.
labelList f(nPoints)
Templated wall surface film model class.
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
void evolveCloud(TrackData &td)
Evolve the cloud.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
void motion(TrackData &td)
Particle motion.
Virtual abstract base class for templated KinematicCloud.
cloudSolution solution_
Solution properties.
Templated base class for kinematic cloud.
void patchData(const parcelType &p, const polyPatch &pp, vector &normal, vector &Up) const
Calculate the patch normal and velocity to interact with,.
label patchi
const List< DynamicList< molecule * > > & cellOccupancy
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
void evolve()
Evolve the cloud.
IOdictionary particleProperties_
Dictionary of particle properties.
autoPtr< volScalarField::Internal > UCoeff_
Coefficient for carrier phase U equation.
Templated stochastic collision model class.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
label index() const
Return the index of this patch in the boundaryMesh.
void info()
Print cloud information.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
forceType forces_
Optional particle forces.
void restoreState()
Reset the current cloud to the previously stored state.
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
functionType functions_
Optional cloud function objects.
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
void scale(DimensionedField< Type, volMesh > &field, const word &name) const
Scale field.
const word cloudName(propsDict.lookup("cloudName"))
filmModelType & surfaceFilm
void cloudReset(KinematicCloud< CloudType > &c)
Reset state of cloud.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
virtual void autoMap(const mapPolyMesh &)
Remap the cells of particles corresponding to the.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
autoPtr< DispersionModel< KinematicCloud< CloudType > > > dispersionModel_
Dispersion model.
IOdictionary outputProperties_
Dictionary of output properties.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
const dimensionedVector & g_
Gravity.
void preEvolve()
Pre-evolve.
void solve(TrackData &td)
Solve the cloud - calls all evolution functions.
cachedRandom rndGen_
Random number generator - used by some injection routines.
label nw
Definition: createFields.H:25
autoPtr< StochasticCollisionModel< KinematicCloud< CloudType > > > stochasticCollisionModel_
Stochastic collision model.
label whichFace(const label l) const
Return label of face in patch from global face label.
Definition: polyPatch.H:377
autoPtr< vectorIntegrationScheme > UIntegrator_
Velocity integration.
const dimensionSet dimVelocity
scalarField cellLengthScale_
Cell length scale.