KinematicCloud.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-2020 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 TrackCloudType>
92 (
93  TrackCloudType& cloud,
94  typename parcelType::trackingData& td
95 )
96 {
97  if (solution_.steadyState())
98  {
99  cloud.storeState();
100 
101  cloud.preEvolve();
102 
103  evolveCloud(cloud, td);
104 
105  if (solution_.coupled())
106  {
107  cloud.relaxSources(cloud.cloudCopy());
108  }
109  }
110  else
111  {
112  cloud.preEvolve();
113 
114  evolveCloud(cloud, td);
115 
116  if (solution_.coupled())
117  {
118  cloud.scaleSources();
119  }
120  }
121 
122  cloud.info();
123 
124  cloud.postEvolve();
125 
126  if (solution_.steadyState())
127  {
128  cloud.restoreState();
129  }
130 }
131 
132 
133 template<class CloudType>
135 {
136  if (cellOccupancyPtr_.empty())
137  {
138  cellOccupancyPtr_.reset
139  (
140  new List<DynamicList<parcelType*>>(mesh_.nCells())
141  );
142  }
143  else if (cellOccupancyPtr_().size() != mesh_.nCells())
144  {
145  // If the size of the mesh has changed, reset the
146  // cellOccupancy size
147 
148  cellOccupancyPtr_().setSize(mesh_.nCells());
149  }
150 
151  List<DynamicList<parcelType*>>& cellOccupancy = cellOccupancyPtr_();
152 
153  forAll(cellOccupancy, cO)
154  {
155  cellOccupancy[cO].clear();
156  }
157 
158  forAllIter(typename KinematicCloud<CloudType>, *this, iter)
159  {
160  cellOccupancy[iter().cell()].append(&iter());
161  }
162 }
163 
164 
165 template<class CloudType>
167 {
168  // Only build the cellOccupancy if the pointer is set, i.e. it has
169  // been requested before.
170 
171  if (cellOccupancyPtr_.valid())
172  {
173  buildCellOccupancy();
174  }
175 }
176 
177 
178 template<class CloudType>
179 template<class TrackCloudType>
181 (
182  TrackCloudType& cloud,
183  typename parcelType::trackingData& td
184 )
185 {
186  if (solution_.coupled())
187  {
188  cloud.resetSourceTerms();
189  }
190 
191  if (solution_.transient())
192  {
193  label preInjectionSize = this->size();
194 
195  this->surfaceFilm().inject(cloud);
196 
197  // Update the cellOccupancy if the size of the cloud has changed
198  // during the injection.
199  if (preInjectionSize != this->size())
200  {
201  updateCellOccupancy();
202  preInjectionSize = this->size();
203  }
204 
205  injectors_.inject(cloud, td);
206 
207 
208  // Assume that motion will update the cellOccupancy as necessary
209  // before it is required.
210  cloud.motion(cloud, td);
211 
212  stochasticCollision().update(td, solution_.trackTime());
213  }
214  else
215  {
216 // this->surfaceFilm().injectSteadyState(cloud);
217 
218  injectors_.injectSteadyState(cloud, td, solution_.trackTime());
219 
220  CloudType::move(cloud, td, solution_.trackTime());
221  }
222 }
223 
224 
225 template<class CloudType>
227 {
228  Info<< endl;
229 
230  if (debug)
231  {
232  this->writePositions();
233  }
234 
235  this->dispersion().cacheFields(false);
236 
237  forces_.cacheFields(false);
238 
239  functions_.postEvolve();
240 
241  solution_.nextIter();
242 
243  if (this->db().time().writeTime())
244  {
245  outputProperties_.writeObject
246  (
247  IOstream::ASCII,
248  IOstream::currentVersion,
249  this->db().time().writeCompression(),
250  true
251  );
252  }
253 }
254 
255 
256 template<class CloudType>
258 {
259  CloudType::cloudReset(c);
260 
261  rndGen_ = c.rndGen_;
262 
263  forces_.transfer(c.forces_);
264 
265  functions_.transfer(c.functions_);
266 
267  injectors_.transfer(c.injectors_);
268 
269  dispersionModel_.reset(c.dispersionModel_.ptr());
270  patchInteractionModel_.reset(c.patchInteractionModel_.ptr());
271  stochasticCollisionModel_.reset(c.stochasticCollisionModel_.ptr());
272  surfaceFilmModel_.reset(c.surfaceFilmModel_.ptr());
273 
274  UIntegrator_.reset(c.UIntegrator_.ptr());
275 }
276 
277 
278 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
279 
280 template<class CloudType>
282 (
283  const word& cloudName,
284  const volScalarField& rho,
285  const volVectorField& U,
286  const volScalarField& mu,
287  const dimensionedVector& g,
288  bool readFields
289 )
290 :
291  CloudType(rho.mesh(), cloudName, false),
292  kinematicCloud(),
293  cloudCopyPtr_(nullptr),
294  mesh_(rho.mesh()),
295  particleProperties_
296  (
297  IOobject
298  (
299  cloudName + "Properties",
300  rho.mesh().time().constant(),
301  rho.mesh(),
302  IOobject::MUST_READ_IF_MODIFIED,
303  IOobject::NO_WRITE
304  )
305  ),
306  outputProperties_
307  (
308  IOobject
309  (
310  cloudName + "OutputProperties",
311  mesh_.time().timeName(),
312  "uniform"/cloud::prefix/cloudName,
313  mesh_,
314  IOobject::READ_IF_PRESENT,
315  IOobject::NO_WRITE
316  )
317  ),
318  solution_(mesh_, particleProperties_.subDict("solution")),
319  constProps_(particleProperties_),
320  subModelProperties_
321  (
322  particleProperties_.subOrEmptyDict("subModels", solution_.active())
323  ),
324  rndGen_(0),
325  cellOccupancyPtr_(),
326  cellLengthScale_(mag(cbrt(mesh_.V()))),
327  rho_(rho),
328  U_(U),
329  mu_(mu),
330  g_(g),
331  pAmbient_(0.0),
332  forces_
333  (
334  *this,
335  mesh_,
336  subModelProperties_.subOrEmptyDict
337  (
338  "particleForces",
339  solution_.active()
340  ),
341  solution_.active()
342  ),
343  functions_
344  (
345  *this,
346  particleProperties_.subOrEmptyDict("cloudFunctions"),
347  solution_.active()
348  ),
349  injectors_
350  (
351  subModelProperties_.subOrEmptyDict("injectionModels"),
352  *this
353  ),
354  dispersionModel_(nullptr),
355  patchInteractionModel_(nullptr),
356  stochasticCollisionModel_(nullptr),
357  surfaceFilmModel_(nullptr),
358  UIntegrator_(nullptr),
359  UTrans_
360  (
362  (
363  IOobject
364  (
365  this->name() + ":UTrans",
366  this->db().time().timeName(),
367  this->db(),
368  IOobject::READ_IF_PRESENT,
369  IOobject::AUTO_WRITE
370  ),
371  mesh_,
373  )
374  ),
375  UCoeff_
376  (
378  (
379  IOobject
380  (
381  this->name() + ":UCoeff",
382  this->db().time().timeName(),
383  this->db(),
384  IOobject::READ_IF_PRESENT,
385  IOobject::AUTO_WRITE
386  ),
387  mesh_,
389  )
390  )
391 {
392  if (solution_.active())
393  {
394  setModels();
395 
396  if (readFields)
397  {
398  parcelType::readFields(*this);
399  this->deleteLostParticles();
400  }
401  }
402 
403  if (solution_.resetSourcesOnStartup())
404  {
405  resetSourceTerms();
406  }
407 }
408 
409 
410 template<class CloudType>
412 (
414  const word& name
415 )
416 :
417  CloudType(c.mesh_, name, c),
418  kinematicCloud(),
419  cloudCopyPtr_(nullptr),
420  mesh_(c.mesh_),
421  particleProperties_(c.particleProperties_),
422  outputProperties_(c.outputProperties_),
423  solution_(c.solution_),
424  constProps_(c.constProps_),
425  subModelProperties_(c.subModelProperties_),
426  rndGen_(c.rndGen_),
427  cellOccupancyPtr_(nullptr),
428  cellLengthScale_(c.cellLengthScale_),
429  rho_(c.rho_),
430  U_(c.U_),
431  mu_(c.mu_),
432  g_(c.g_),
433  pAmbient_(c.pAmbient_),
434  forces_(c.forces_),
435  functions_(c.functions_),
436  injectors_(c.injectors_),
437  dispersionModel_(c.dispersionModel_->clone()),
438  patchInteractionModel_(c.patchInteractionModel_->clone()),
439  stochasticCollisionModel_(c.stochasticCollisionModel_->clone()),
440  surfaceFilmModel_(c.surfaceFilmModel_->clone()),
441  UIntegrator_(c.UIntegrator_->clone()),
442  UTrans_
443  (
445  (
446  IOobject
447  (
448  this->name() + ":UTrans",
449  this->db().time().timeName(),
450  this->db(),
451  IOobject::NO_READ,
452  IOobject::NO_WRITE,
453  false
454  ),
455  c.UTrans_()
456  )
457  ),
458  UCoeff_
459  (
461  (
462  IOobject
463  (
464  name + ":UCoeff",
465  this->db().time().timeName(),
466  this->db(),
467  IOobject::NO_READ,
468  IOobject::NO_WRITE,
469  false
470  ),
471  c.UCoeff_()
472  )
473  )
474 {}
475 
476 
477 template<class CloudType>
479 (
480  const fvMesh& mesh,
481  const word& name,
483 )
484 :
485  CloudType(mesh, name, IDLList<parcelType>()),
486  kinematicCloud(),
487  cloudCopyPtr_(nullptr),
488  mesh_(mesh),
489  particleProperties_
490  (
491  IOobject
492  (
493  name + "Properties",
494  mesh.time().constant(),
495  mesh,
496  IOobject::NO_READ,
497  IOobject::NO_WRITE,
498  false
499  )
500  ),
501  outputProperties_
502  (
503  IOobject
504  (
505  name + "OutputProperties",
506  mesh_.time().timeName(),
507  "uniform"/cloud::prefix/name,
508  mesh_,
509  IOobject::NO_READ,
510  IOobject::NO_WRITE,
511  false
512  )
513  ),
514  solution_(mesh),
515  constProps_(),
516  subModelProperties_(dictionary::null),
517  rndGen_(0),
518  cellOccupancyPtr_(nullptr),
519  cellLengthScale_(c.cellLengthScale_),
520  rho_(c.rho_),
521  U_(c.U_),
522  mu_(c.mu_),
523  g_(c.g_),
524  pAmbient_(c.pAmbient_),
525  forces_(*this, mesh),
526  functions_(*this),
527  injectors_(*this),
528  dispersionModel_(nullptr),
529  patchInteractionModel_(nullptr),
530  stochasticCollisionModel_(nullptr),
531  surfaceFilmModel_(nullptr),
532  UIntegrator_(nullptr),
533  UTrans_(nullptr),
534  UCoeff_(nullptr)
535 {}
536 
537 
538 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
539 
540 template<class CloudType>
542 {}
543 
544 
545 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
546 
547 template<class CloudType>
549 (
550  parcelType& parcel,
551  const scalar lagrangianDt
552 )
553 {
554  parcel.rho() = constProps_.rho0();
555 }
556 
557 
558 template<class CloudType>
560 (
561  parcelType& parcel,
562  const scalar lagrangianDt,
563  const bool fullyDescribed
564 )
565 {
566  const scalar carrierDt = mesh_.time().deltaTValue();
567  parcel.stepFraction() = (carrierDt - lagrangianDt)/carrierDt;
568 
569  if (parcel.typeId() == -1)
570  {
571  parcel.typeId() = constProps_.parcelTypeId();
572  }
573 }
574 
575 
576 template<class CloudType>
578 {
579  cloudCopyPtr_.reset
580  (
581  static_cast<KinematicCloud<CloudType>*>
582  (
583  clone(this->name() + "Copy").ptr()
584  )
585  );
586 }
587 
588 
589 template<class CloudType>
591 {
592  cloudReset(cloudCopyPtr_());
593  cloudCopyPtr_.clear();
594 }
595 
596 
597 template<class CloudType>
599 {
600  UTrans().field() = Zero;
601  UCoeff().field() = 0.0;
602 }
603 
604 
605 template<class CloudType>
606 template<class Type>
608 (
610  const DimensionedField<Type, volMesh>& field0,
611  const word& name
612 ) const
613 {
614  const scalar coeff = solution_.relaxCoeff(name);
615  field = field0 + coeff*(field - field0);
616 }
617 
618 
619 template<class CloudType>
620 template<class Type>
622 (
624  const word& name
625 ) const
626 {
627  const scalar coeff = solution_.relaxCoeff(name);
628  field *= coeff;
629 }
630 
631 
632 template<class CloudType>
634 (
635  const KinematicCloud<CloudType>& cloudOldTime
636 )
637 {
638  this->relax(UTrans_(), cloudOldTime.UTrans(), "U");
639  this->relax(UCoeff_(), cloudOldTime.UCoeff(), "U");
640 }
641 
642 
643 template<class CloudType>
645 {
646  this->scale(UTrans_(), "U");
647  this->scale(UCoeff_(), "U");
648 }
649 
650 
651 template<class CloudType>
653 {
654  // force calculation of mesh dimensions - needed for parallel runs
655  // with topology change due to lazy evaluation of valid mesh dimensions
656  label nGeometricD = mesh_.nGeometricD();
657 
658  Info<< "\nSolving " << nGeometricD << "-D cloud " << this->name() << endl;
659 
660  this->dispersion().cacheFields(true);
661  forces_.cacheFields(true);
662  updateCellOccupancy();
663 
664  pAmbient_ = constProps_.dict().template
665  lookupOrDefault<scalar>("pAmbient", pAmbient_);
666 
667  functions_.preEvolve();
668 }
669 
670 
671 template<class CloudType>
673 {
674  if (solution_.canEvolve())
675  {
676  typename parcelType::trackingData td(*this);
677 
678  solve(*this, td);
679  }
680 }
681 
682 
683 template<class CloudType>
684 template<class TrackCloudType>
686 (
687  TrackCloudType& cloud,
688  typename parcelType::trackingData& td
689 )
690 {
691  CloudType::move(cloud, td, solution_.trackTime());
692 
693  updateCellOccupancy();
694 }
695 
696 
697 template<class CloudType>
699 (
700  const parcelType& p,
701  const polyPatch& pp,
702  vector& nw,
703  vector& Up
704 ) const
705 {
706  p.patchData(nw, Up);
707  Up /= p.mesh().time().deltaTValue();
708 
709  // If this is a wall patch, then there may be a non-zero tangential velocity
710  // component; the lid velocity in a lid-driven cavity case, for example. We
711  // want the particle to interact with this velocity, so we look it up in the
712  // velocity field and use it to set the wall-tangential component.
713  if (!mesh_.moving() && isA<wallPolyPatch>(pp))
714  {
715  const label patchi = pp.index();
716  const label patchFacei = pp.whichFace(p.face());
717 
718  // We only want to use the boundary condition value only if it is set
719  // by the boundary condition. If the boundary values are extrapolated
720  // (e.g., slip conditions) then they represent the motion of the fluid
721  // just inside the domain rather than that of the wall itself.
722  if (U_.boundaryField()[patchi].fixesValue())
723  {
724  const vector Uw1 = U_.boundaryField()[patchi][patchFacei];
725  const vector& Uw0 =
726  U_.oldTime().boundaryField()[patchi][patchFacei];
727 
728  const scalar f = p.currentTimeFraction();
729 
730  const vector Uw = Uw0 + f*(Uw1 - Uw0);
731 
732  const tensor nnw = nw*nw;
733 
734  Up = (nnw & Up) + Uw - (nnw & Uw);
735  }
736  }
737 }
738 
739 
740 template<class CloudType>
742 {
743  updateCellOccupancy();
744  injectors_.updateMesh();
745  cellLengthScale_ = mag(cbrt(mesh_.V()));
746 }
747 
748 
749 template<class CloudType>
751 {
753 
754  updateMesh();
755 }
756 
757 
758 template<class CloudType>
760 {
761  vector linearMomentum = linearMomentumOfSystem();
762  reduce(linearMomentum, sumOp<vector>());
763 
764  scalar linearKineticEnergy = linearKineticEnergyOfSystem();
765  reduce(linearKineticEnergy, sumOp<scalar>());
766 
767  Info<< "Cloud: " << this->name() << nl
768  << " Current number of parcels = "
769  << returnReduce(this->size(), sumOp<label>()) << nl
770  << " Current mass in system = "
771  << returnReduce(massInSystem(), sumOp<scalar>()) << nl
772  << " Linear momentum = "
773  << linearMomentum << nl
774  << " |Linear momentum| = "
775  << mag(linearMomentum) << nl
776  << " Linear kinetic energy = "
777  << linearKineticEnergy << nl;
778 
779  injectors_.info(Info);
780  this->surfaceFilm().info(Info);
781  this->patchInteraction().info(Info);
782 }
783 
784 
785 // ************************************************************************* //
Template class for intrusive linked lists.
Definition: ILList.H:50
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
UEqn relax()
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.
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:59
virtual ~KinematicCloud()
Destructor.
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:459
void evolveCloud(TrackCloudType &cloud, typename parcelType::trackingData &td)
Evolve the cloud.
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:251
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.
rhoEqn solve()
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:239
KinematicCloud(const word &cloudName, const volScalarField &rho, const volVectorField &U, const volScalarField &mu, const dimensionedVector &g, bool readFields=true)
Construct given carrier gas fields.
autoPtr< SurfaceFilmModel< KinematicCloud< CloudType > > > surfaceFilmModel_
Surface film model.
Random rndGen_
Random number generator - used by some injection routines.
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/m^3].
Templated patch interaction model class.
void updateMesh()
Update mesh.
const volScalarField & mu_
Dynamic viscosity [Pa.s].
T clone(const T &t)
Definition: List.H:54
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.
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:125
A class for handling words, derived from string.
Definition: word.H:59
void motion(TrackCloudType &cloud, typename parcelType::trackingData &td)
Particle motion.
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:178
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:97
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
regionModels::surfaceFilmModel & surfaceFilm
static const char nl
Definition: Ostream.H:260
const Mesh & mesh() const
Return mesh.
labelList f(nPoints)
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
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.
autoPtr< integrationScheme > UIntegrator_
Velocity integration.
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"))
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(TrackCloudType &cloud, typename parcelType::trackingData &td)
Solve the cloud - calls all evolution functions.
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:383
const dimensionSet dimVelocity
scalarField cellLengthScale_
Cell length scale.