MomentumCloud.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 "MomentumCloud.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  filmModel_.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  this->changeTimeStep();
98 
99  if (solution_.steadyState())
100  {
101  cloud.storeState();
102  }
103 
104  cloud.preEvolve();
105 
106  if (solution_.coupled())
107  {
108  cloud.resetSourceTerms();
109  }
110 
111  if (solution_.transient())
112  {
113  label preInjectionSize = this->size();
114 
115  this->surfaceFilm().inject(cloud);
116 
117  // Update the cellOccupancy if the size of the cloud has changed
118  // during the injection.
119  if (preInjectionSize != this->size())
120  {
121  updateCellOccupancy();
122  preInjectionSize = this->size();
123  }
124 
125  injectors_.inject(cloud, td);
126 
127  // Assume that motion will update the cellOccupancy as necessary
128  // before it is required.
129  cloud.motion(cloud, td);
130 
131  stochasticCollision().update(td);
132  }
133  else
134  {
135  injectors_.injectSteadyState(cloud, td);
136 
137  CloudType::move(cloud, td);
138  }
139 
140  if (solution_.coupled() && solution_.transient())
141  {
142  cloud.scaleSources();
143  }
144 
145  if (solution_.coupled() && solution_.steadyState())
146  {
147  cloud.relaxSources(cloud.cloudCopy());
148  }
149 
150  cloud.info();
151 
152  cloud.postEvolve();
153 
154  if (solution_.steadyState())
155  {
156  cloud.restoreState();
157  }
158 }
159 
160 
161 template<class CloudType>
163 {
164  if (cellOccupancyPtr_.empty())
165  {
166  cellOccupancyPtr_.reset
167  (
168  new List<DynamicList<parcelType*>>(this->mesh().nCells())
169  );
170  }
171  else if (cellOccupancyPtr_().size() != this->mesh().nCells())
172  {
173  // If the size of the mesh has changed, reset the
174  // cellOccupancy size
175 
176  cellOccupancyPtr_().setSize(this->mesh().nCells());
177  }
178 
179  List<DynamicList<parcelType*>>& cellOccupancy = cellOccupancyPtr_();
180 
181  forAll(cellOccupancy, cO)
182  {
183  cellOccupancy[cO].clear();
184  }
185 
186  forAllIter(typename MomentumCloud<CloudType>, *this, iter)
187  {
188  cellOccupancy[iter().cell()].append(&iter());
189  }
190 }
191 
192 
193 template<class CloudType>
195 {
196  // Only build the cellOccupancy if the pointer is set, i.e. it has
197  // been requested before.
198 
199  if (cellOccupancyPtr_.valid())
200  {
201  buildCellOccupancy();
202  }
203 }
204 
205 
206 template<class CloudType>
208 {
209  Info<< endl;
210 
211  if (debug)
212  {
213  this->writePositions();
214  }
215 
216  this->dispersion().cacheFields(false);
217 
218  forces_.cacheFields(false);
219 
220  functions_.postEvolve();
221 
222  solution_.nextIter();
223 
224  if (this->db().time().writeTime())
225  {
226  outputProperties_.writeObject
227  (
228  IOstream::ASCII,
229  IOstream::currentVersion,
230  this->db().time().writeCompression(),
231  true
232  );
233  }
234 }
235 
236 
237 template<class CloudType>
239 {
240  CloudType::cloudReset(c);
241 
242  forces_.transfer(c.forces_);
243 
244  functions_.transfer(c.functions_);
245 
246  injectors_.transfer(c.injectors_);
247 
248  dispersionModel_.reset(c.dispersionModel_.ptr());
249  patchInteractionModel_.reset(c.patchInteractionModel_.ptr());
250  stochasticCollisionModel_.reset(c.stochasticCollisionModel_.ptr());
251  filmModel_.reset(c.filmModel_.ptr());
252 
253  UIntegrator_.reset(c.UIntegrator_.ptr());
254 }
255 
256 
257 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
258 
259 template<class CloudType>
261 (
262  const word& cloudName,
263  const volScalarField& rho,
264  const volVectorField& U,
265  const volScalarField& mu,
266  const dimensionedVector& g,
267  const bool readFields
268 )
269 :
270  CloudType(cloudName, rho, U, mu, g, false),
271  cloudCopyPtr_(nullptr),
272  particleProperties_
273  (
274  IOobject
275  (
276  cloudName + "Properties",
277  this->mesh().time().constant(),
278  this->mesh(),
279  IOobject::MUST_READ_IF_MODIFIED,
280  IOobject::NO_WRITE
281  )
282  ),
283  outputProperties_
284  (
285  IOobject
286  (
287  cloudName + "OutputProperties",
288  this->mesh().time().name(),
289  "uniform"/cloud::prefix/cloudName,
290  this->mesh(),
291  IOobject::READ_IF_PRESENT,
292  IOobject::NO_WRITE
293  )
294  ),
295  solution_(this->mesh(), particleProperties_.subDict("solution")),
296  constProps_(particleProperties_),
297  subModelProperties_
298  (
299  particleProperties_.subOrEmptyDict("subModels", true)
300  ),
301  cpuLoad_(particleProperties_.lookupOrDefault("cpuLoad", false)),
302  rndGen_(0),
303  stdNormal_(rndGen_.generator()),
304  cellOccupancyPtr_(),
305  cellLengthScale_(mag(cbrt(this->mesh().V()))),
306  rho_(rho),
307  U_(U),
308  mu_(mu),
309  g_(g),
310  pAmbient_(0.0),
311  forces_
312  (
313  *this,
314  this->mesh(),
315  subModelProperties_.subOrEmptyDict
316  (
317  "particleForces",
318  true
319  ),
320  true
321  ),
322  functions_
323  (
324  *this,
325  particleProperties_.subOrEmptyDict("cloudFunctions")
326  ),
327  injectors_
328  (
329  subModelProperties_.subOrEmptyDict("injectionModels"),
330  *this
331  ),
332  dispersionModel_(nullptr),
333  patchInteractionModel_(nullptr),
334  stochasticCollisionModel_(nullptr),
335  filmModel_(nullptr),
336  UIntegrator_(nullptr),
337  UTrans_
338  (
339  new volVectorField::Internal
340  (
341  IOobject
342  (
343  this->name() + ":UTrans",
344  this->db().time().name(),
345  this->db(),
346  IOobject::READ_IF_PRESENT,
347  IOobject::AUTO_WRITE
348  ),
349  this->mesh(),
351  )
352  ),
353  UCoeff_
354  (
355  new volScalarField::Internal
356  (
357  IOobject
358  (
359  this->name() + ":UCoeff",
360  this->db().time().name(),
361  this->db(),
362  IOobject::READ_IF_PRESENT,
363  IOobject::AUTO_WRITE
364  ),
365  this->mesh(),
367  )
368  )
369 {
370  setModels();
371 
372  if (readFields)
373  {
374  parcelType::readFields(*this);
375  this->deleteLostParticles();
376  }
377 
379  {
381  }
382 }
383 
384 
385 template<class CloudType>
387 (
388  const word& cloudName,
389  const volScalarField& rho,
390  const volVectorField& U,
391  const dimensionedVector& g,
392  const fluidThermo& carrierThermo,
393  const bool readFields
394 )
395 :
396  MomentumCloud(cloudName, rho, U, carrierThermo.mu(), g, readFields)
397 {}
398 
399 
400 template<class CloudType>
402 (
404  const word& name
405 )
406 :
407  CloudType(c, name),
408  cloudCopyPtr_(nullptr),
409  particleProperties_(c.particleProperties_),
410  outputProperties_(c.outputProperties_),
411  solution_(c.solution_),
412  constProps_(c.constProps_),
413  subModelProperties_(c.subModelProperties_),
414  cpuLoad_(c.cpuLoad_),
415  rndGen_(c.rndGen_),
416  stdNormal_(c.stdNormal_),
417  cellOccupancyPtr_(nullptr),
418  cellLengthScale_(c.cellLengthScale_),
419  rho_(c.rho_),
420  U_(c.U_),
421  mu_(c.mu_),
422  g_(c.g_),
423  pAmbient_(c.pAmbient_),
424  forces_(c.forces_),
425  functions_(c.functions_),
426  injectors_(c.injectors_),
427  dispersionModel_(c.dispersionModel_->clone()),
428  patchInteractionModel_(c.patchInteractionModel_->clone()),
429  stochasticCollisionModel_(c.stochasticCollisionModel_->clone()),
430  filmModel_(c.filmModel_->clone()),
431  UIntegrator_(c.UIntegrator_->clone()),
432  UTrans_
433  (
434  new volVectorField::Internal
435  (
436  IOobject
437  (
438  this->name() + ":UTrans",
439  this->db().time().name(),
440  this->db(),
441  IOobject::NO_READ,
442  IOobject::NO_WRITE,
443  false
444  ),
445  c.UTrans_()
446  )
447  ),
448  UCoeff_
449  (
450  new volScalarField::Internal
451  (
452  IOobject
453  (
454  name + ":UCoeff",
455  this->db().time().name(),
456  this->db(),
457  IOobject::NO_READ,
458  IOobject::NO_WRITE,
459  false
460  ),
461  c.UCoeff_()
462  )
463  )
464 {}
465 
466 
467 template<class CloudType>
469 (
470  const fvMesh& mesh,
471  const word& name,
473 )
474 :
475  CloudType(mesh, name, c),
476  cloudCopyPtr_(nullptr),
477  particleProperties_
478  (
479  IOobject
480  (
481  name + "Properties",
482  mesh.time().constant(),
483  mesh,
484  IOobject::NO_READ,
485  IOobject::NO_WRITE,
486  false
487  )
488  ),
489  outputProperties_
490  (
491  IOobject
492  (
493  name + "OutputProperties",
494  this->mesh().time().name(),
495  "uniform"/cloud::prefix/name,
496  this->mesh(),
497  IOobject::NO_READ,
498  IOobject::NO_WRITE,
499  false
500  )
501  ),
502  solution_(mesh),
503  constProps_(),
504  subModelProperties_(dictionary::null),
505  cpuLoad_(c.cpuLoad_),
506  rndGen_(0),
507  stdNormal_(rndGen_.generator()),
508  cellOccupancyPtr_(nullptr),
509  cellLengthScale_(c.cellLengthScale_),
510  rho_(c.rho_),
511  U_(c.U_),
512  mu_(c.mu_),
513  g_(c.g_),
514  pAmbient_(c.pAmbient_),
515  forces_(*this, mesh),
516  functions_(*this),
517  injectors_(*this),
518  dispersionModel_(nullptr),
519  patchInteractionModel_(nullptr),
520  stochasticCollisionModel_(nullptr),
521  filmModel_(nullptr),
522  UIntegrator_(nullptr),
523  UTrans_(nullptr),
524  UCoeff_(nullptr)
525 {}
526 
527 
528 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
529 
530 template<class CloudType>
532 {}
533 
534 
535 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
536 
537 template<class CloudType>
539 (
540  parcelType& parcel
541 )
542 {
543  parcel.rho() = constProps_.rho0();
544 }
545 
546 
547 template<class CloudType>
549 (
550  parcelType& parcel,
551  const label injectori
552 )
553 {
554  if (parcel.typeId() == -1)
555  {
556  parcel.typeId() = constProps_.parcelTypeId();
557  }
558 }
559 
560 
561 template<class CloudType>
563 {
564  cloudCopyPtr_.reset
565  (
566  static_cast<MomentumCloud<CloudType>*>
567  (
568  clone(this->name() + "Copy").ptr()
569  )
570  );
571 }
572 
573 
574 template<class CloudType>
576 {
577  cloudReset(cloudCopyPtr_());
578  cloudCopyPtr_.clear();
579 }
580 
581 
582 template<class CloudType>
584 {
585  UTransRef().primitiveFieldRef() = Zero;
586  UCoeffRef().primitiveFieldRef() = 0.0;
587 }
588 
589 
590 template<class CloudType>
591 template<class Type>
593 (
595  const DimensionedField<Type, volMesh>& field0,
596  const word& name
597 ) const
598 {
599  const scalar coeff = solution_.relaxCoeff(name);
600  field = field0 + coeff*(field - field0);
601 }
602 
603 
604 template<class CloudType>
605 template<class Type>
607 (
609  const word& name
610 ) const
611 {
612  const scalar coeff = solution_.relaxCoeff(name);
613  field *= coeff;
614 }
615 
616 
617 template<class CloudType>
619 (
620  const MomentumCloud<CloudType>& cloudOldTime
621 )
622 {
623  this->relax(UTrans_(), cloudOldTime.UTrans_(), "U");
624  this->relax(UCoeff_(), cloudOldTime.UCoeff_(), "U");
625 }
626 
627 
628 template<class CloudType>
630 {
631  this->scale(UTrans_(), "U");
632  this->scale(UCoeff_(), "U");
633 }
634 
635 
636 template<class CloudType>
638 {
639  // force calculation of mesh dimensions - needed for parallel runs
640  // with topology change due to lazy evaluation of valid mesh dimensions
641  label nGeometricD = this->mesh().nGeometricD();
642 
643  Info<< nl << "Solving " << nGeometricD << "-D cloud " << this->name()
644  << endl;
645 
646  this->dispersion().cacheFields(true);
647  forces_.cacheFields(true);
648  updateCellOccupancy();
649 
650  pAmbient_ = constProps_.dict().template
651  lookupOrDefault<scalar>("pAmbient", pAmbient_);
652 
653  functions_.preEvolve();
654 }
655 
656 
657 template<class CloudType>
659 {
660  if (solution_.canEvolve())
661  {
662  typename parcelType::trackingData td(*this);
663 
664  solve(*this, td);
665  }
666 }
667 
668 
669 template<class CloudType>
670 template<class TrackCloudType>
672 (
673  TrackCloudType& cloud,
674  typename parcelType::trackingData& td
675 )
676 {
677  CloudType::move(cloud, td);
678 
679  updateCellOccupancy();
680 }
681 
682 
683 template<class CloudType>
685 (
686  const parcelType& p,
687  const polyPatch& pp,
688  vector& nw,
689  vector& Up
690 ) const
691 {
692  p.patchData(this->mesh(), nw, Up);
693  Up /= this->mesh().time().deltaTValue();
694 
695  // If this is a wall patch, then there may be a non-zero tangential velocity
696  // component; the lid velocity in a lid-driven cavity case, for example. We
697  // want the particle to interact with this velocity, so we look it up in the
698  // velocity field and use it to set the wall-tangential component.
699  if (!this->mesh().moving() && isA<wallPolyPatch>(pp))
700  {
701  const label patchi = pp.index();
702  const label patchFacei = pp.whichFace(p.face());
703 
704  // We only want to use the boundary condition value only if it is set
705  // by the boundary condition. If the boundary values are extrapolated
706  // (e.g., slip conditions) then they represent the motion of the fluid
707  // just inside the domain rather than that of the wall itself.
708  if (U_.boundaryField()[patchi].fixesValue())
709  {
710  const vector Uw1 = U_.boundaryField()[patchi][patchFacei];
711  const vector& Uw0 =
712  U_.oldTime().boundaryField()[patchi][patchFacei];
713 
714  const vector Uw = Uw0 + p.stepFraction()*(Uw1 - Uw0);
715 
716  const tensor nnw = nw*nw;
717 
718  Up = (nnw & Up) + Uw - (nnw & Uw);
719  }
720  }
721 }
722 
723 
724 template<class CloudType>
726 {
728 
729  updateCellOccupancy();
730 
731  injectors_.topoChange();
732 
733  cellLengthScale_ = mag(cbrt(this->mesh().V()));
734 }
735 
736 
737 template<class CloudType>
739 {
741 
742  updateCellOccupancy();
743 
744  injectors_.topoChange();
745 
746  cellLengthScale_ = mag(cbrt(this->mesh().V()));
747 }
748 
749 
750 template<class CloudType>
752 {
754 
755  updateCellOccupancy();
756 
757  injectors_.topoChange();
758 
759  cellLengthScale_ = mag(cbrt(this->mesh().V()));
760 }
761 
762 
763 template<class CloudType>
765 {
766  vector linearMomentum = linearMomentumOfSystem();
767  reduce(linearMomentum, sumOp<vector>());
768 
769  scalar linearKineticEnergy = linearKineticEnergyOfSystem();
770  reduce(linearKineticEnergy, sumOp<scalar>());
771 
772  Info<< "Cloud: " << this->name() << nl
773  << " Current number of parcels = "
774  << returnReduce(this->size(), sumOp<label>()) << nl
775  << " Current mass in system = "
776  << returnReduce(massInSystem(), sumOp<scalar>()) << nl
777  << " Linear momentum = "
778  << linearMomentum << nl
779  << " |Linear momentum| = "
780  << mag(linearMomentum) << nl
781  << " Linear kinetic energy = "
782  << linearKineticEnergy << nl;
783 
784  injectors_.info(Info);
785  this->surfaceFilm().info(Info);
786  this->patchInteraction().info(Info);
787 }
788 
789 
790 // ************************************************************************* //
EaEqn relax()
#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
const List< DynamicList< molecule * > > & cellOccupancy
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
void deleteLostParticles()
Remove lost particles from cloud and delete.
Definition: Cloud.C:233
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
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:225
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:78
Generic GeometricField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
InfoProxy< IOobject > info() const
Return info proxy.
Definition: IOobject.H:486
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:91
Templated base class for momentum cloud.
void postEvolve()
Post-evolve.
cloudSolution solution_
Solution properties.
autoPtr< volVectorField::Internal > UTrans_
Momentum.
void setModels()
Set cloud sub-models.
Definition: MomentumCloud.C:40
void motion(TrackCloudType &cloud, typename parcelType::trackingData &td)
Particle motion.
virtual ~MomentumCloud()
Destructor.
void scale(DimensionedField< Type, volMesh > &field, const word &name) const
Scale field.
void storeState()
Store the current cloud state.
void patchData(const parcelType &p, const polyPatch &pp, vector &normal, vector &Up) const
Calculate the patch normal and velocity to interact with,.
void relax(DimensionedField< Type, volMesh > &field, const DimensionedField< Type, volMesh > &field0, const word &name) const
Relax field.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
void scaleSources()
Apply scaling to (transient) cloud sources.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
MomentumCloud(const word &cloudName, const volScalarField &rho, const volVectorField &U, const volScalarField &mu, const dimensionedVector &g, const bool readFields=true)
Construct given carrier fields.
void updateCellOccupancy()
Update (i.e. build) the cellOccupancy if it has.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
void evolve()
Evolve the cloud.
void cloudReset(MomentumCloud< CloudType > &c)
Reset state of cloud.
void info()
Print cloud information.
void restoreState()
Reset the current cloud to the previously stored state.
void checkParcelProperties(parcelType &parcel, const label injectori)
Check parcel properties.
void preEvolve()
Pre-evolve.
void resetSourceTerms()
Reset the cloud source terms.
void buildCellOccupancy()
Build the cellOccupancy.
autoPtr< volScalarField::Internal > UCoeff_
Coefficient for carrier phase U equation.
void setParcelThermoProperties(parcelType &parcel)
Set parcel thermo properties.
void relaxSources(const MomentumCloud< CloudType > &cloudOldTime)
Apply relaxation to (steady state) cloud sources.
void solve(TrackCloudType &cloud, typename parcelType::trackingData &td)
Solve the cloud - calls all evolution functions.
Definition: MomentumCloud.C:92
Templated patch interaction model class.
Templated stochastic collision model class.
Templated wall surface film model class.
const Switch resetSourcesOnStartup() const
Return const access to the reset sources flag.
A cloud is a collection of lagrangian particles.
Definition: cloud.H:55
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Base-class for fluid thermodynamic properties.
Definition: fluidThermo.H:57
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
label index() const
Return the index of this patch in the boundaryMesh.
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
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
label whichFace(const label l) const
Return label of face in patch from global face label.
Definition: polyPatch.H:360
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A class for handling words, derived from string.
Definition: word.H:62
label patchi
U
Definition: pEqn.H:72
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
const dimensionedScalar mu
Atomic mass unit.
const dimensionedScalar c
Speed of light in a vacuum.
static const zero Zero
Definition: zero.H:97
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
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
messageStream Info
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
dimensioned< scalar > mag(const dimensioned< Type > &)
T clone(const T &t)
Definition: List.H:55
const dimensionSet dimMass
const dimensionSet dimVelocity
dimensionedScalar cbrt(const dimensionedScalar &ds)
static const char nl
Definition: Ostream.H:266
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.
volScalarField & p
const word cloudName(propsDict.lookup("cloudName"))