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-2023 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  rndGen_ = c.rndGen_;
243 
244  forces_.transfer(c.forces_);
245 
246  functions_.transfer(c.functions_);
247 
248  injectors_.transfer(c.injectors_);
249 
250  dispersionModel_.reset(c.dispersionModel_.ptr());
251  patchInteractionModel_.reset(c.patchInteractionModel_.ptr());
252  stochasticCollisionModel_.reset(c.stochasticCollisionModel_.ptr());
253  filmModel_.reset(c.filmModel_.ptr());
254 
255  UIntegrator_.reset(c.UIntegrator_.ptr());
256 }
257 
258 
259 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
260 
261 template<class CloudType>
263 (
264  const word& cloudName,
265  const volScalarField& rho,
266  const volVectorField& U,
267  const volScalarField& mu,
268  const dimensionedVector& g,
269  const bool readFields
270 )
271 :
272  CloudType(cloudName, rho, U, mu, g, false),
273  cloudCopyPtr_(nullptr),
274  particleProperties_
275  (
276  IOobject
277  (
278  cloudName + "Properties",
279  this->mesh().time().constant(),
280  this->mesh(),
281  IOobject::MUST_READ_IF_MODIFIED,
282  IOobject::NO_WRITE
283  )
284  ),
285  outputProperties_
286  (
287  IOobject
288  (
289  cloudName + "OutputProperties",
290  this->mesh().time().name(),
291  "uniform"/cloud::prefix/cloudName,
292  this->mesh(),
293  IOobject::READ_IF_PRESENT,
294  IOobject::NO_WRITE
295  )
296  ),
297  solution_(this->mesh(), particleProperties_.subDict("solution")),
298  constProps_(particleProperties_),
299  subModelProperties_
300  (
301  particleProperties_.subOrEmptyDict("subModels", true)
302  ),
303  rndGen_(0),
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  true
327  ),
328  injectors_
329  (
330  subModelProperties_.subOrEmptyDict("injectionModels"),
331  *this
332  ),
333  dispersionModel_(nullptr),
334  patchInteractionModel_(nullptr),
335  stochasticCollisionModel_(nullptr),
336  filmModel_(nullptr),
337  UIntegrator_(nullptr),
338  UTrans_
339  (
340  new volVectorField::Internal
341  (
342  IOobject
343  (
344  this->name() + ":UTrans",
345  this->db().time().name(),
346  this->db(),
347  IOobject::READ_IF_PRESENT,
348  IOobject::AUTO_WRITE
349  ),
350  this->mesh(),
352  )
353  ),
354  UCoeff_
355  (
356  new volScalarField::Internal
357  (
358  IOobject
359  (
360  this->name() + ":UCoeff",
361  this->db().time().name(),
362  this->db(),
363  IOobject::READ_IF_PRESENT,
364  IOobject::AUTO_WRITE
365  ),
366  this->mesh(),
368  )
369  )
370 {
371  setModels();
372 
373  if (readFields)
374  {
375  parcelType::readFields(*this);
376  this->deleteLostParticles();
377  }
378 
380  {
382  }
383 }
384 
385 
386 template<class CloudType>
388 (
389  const word& cloudName,
390  const volScalarField& rho,
391  const volVectorField& U,
392  const dimensionedVector& g,
393  const fluidThermo& carrierThermo,
394  const bool readFields
395 )
396 :
397  MomentumCloud(cloudName, rho, U, carrierThermo.mu(), g, readFields)
398 {}
399 
400 
401 template<class CloudType>
403 (
405  const word& name
406 )
407 :
408  CloudType(c, name),
409  cloudCopyPtr_(nullptr),
410  particleProperties_(c.particleProperties_),
411  outputProperties_(c.outputProperties_),
412  solution_(c.solution_),
413  constProps_(c.constProps_),
414  subModelProperties_(c.subModelProperties_),
415  rndGen_(c.rndGen_),
416  cellOccupancyPtr_(nullptr),
417  cellLengthScale_(c.cellLengthScale_),
418  rho_(c.rho_),
419  U_(c.U_),
420  mu_(c.mu_),
421  g_(c.g_),
422  pAmbient_(c.pAmbient_),
423  forces_(c.forces_),
424  functions_(c.functions_),
425  injectors_(c.injectors_),
426  dispersionModel_(c.dispersionModel_->clone()),
427  patchInteractionModel_(c.patchInteractionModel_->clone()),
428  stochasticCollisionModel_(c.stochasticCollisionModel_->clone()),
429  filmModel_(c.filmModel_->clone()),
430  UIntegrator_(c.UIntegrator_->clone()),
431  UTrans_
432  (
433  new volVectorField::Internal
434  (
435  IOobject
436  (
437  this->name() + ":UTrans",
438  this->db().time().name(),
439  this->db(),
440  IOobject::NO_READ,
441  IOobject::NO_WRITE,
442  false
443  ),
444  c.UTrans_()
445  )
446  ),
447  UCoeff_
448  (
449  new volScalarField::Internal
450  (
451  IOobject
452  (
453  name + ":UCoeff",
454  this->db().time().name(),
455  this->db(),
456  IOobject::NO_READ,
457  IOobject::NO_WRITE,
458  false
459  ),
460  c.UCoeff_()
461  )
462  )
463 {}
464 
465 
466 template<class CloudType>
468 (
469  const fvMesh& mesh,
470  const word& name,
472 )
473 :
474  CloudType(mesh, name, c),
475  cloudCopyPtr_(nullptr),
476  particleProperties_
477  (
478  IOobject
479  (
480  name + "Properties",
481  mesh.time().constant(),
482  mesh,
483  IOobject::NO_READ,
484  IOobject::NO_WRITE,
485  false
486  )
487  ),
488  outputProperties_
489  (
490  IOobject
491  (
492  name + "OutputProperties",
493  this->mesh().time().name(),
494  "uniform"/cloud::prefix/name,
495  this->mesh(),
496  IOobject::NO_READ,
497  IOobject::NO_WRITE,
498  false
499  )
500  ),
501  solution_(mesh),
502  constProps_(),
503  subModelProperties_(dictionary::null),
504  rndGen_(0),
505  cellOccupancyPtr_(nullptr),
506  cellLengthScale_(c.cellLengthScale_),
507  rho_(c.rho_),
508  U_(c.U_),
509  mu_(c.mu_),
510  g_(c.g_),
511  pAmbient_(c.pAmbient_),
512  forces_(*this, mesh),
513  functions_(*this),
514  injectors_(*this),
515  dispersionModel_(nullptr),
516  patchInteractionModel_(nullptr),
517  stochasticCollisionModel_(nullptr),
518  filmModel_(nullptr),
519  UIntegrator_(nullptr),
520  UTrans_(nullptr),
521  UCoeff_(nullptr)
522 {}
523 
524 
525 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
526 
527 template<class CloudType>
529 {}
530 
531 
532 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
533 
534 template<class CloudType>
536 (
537  parcelType& parcel
538 )
539 {
540  parcel.rho() = constProps_.rho0();
541 }
542 
543 
544 template<class CloudType>
546 (
547  parcelType& parcel,
548  const label injectori
549 )
550 {
551  if (parcel.typeId() == -1)
552  {
553  parcel.typeId() = constProps_.parcelTypeId();
554  }
555 }
556 
557 
558 template<class CloudType>
560 {
561  cloudCopyPtr_.reset
562  (
563  static_cast<MomentumCloud<CloudType>*>
564  (
565  clone(this->name() + "Copy").ptr()
566  )
567  );
568 }
569 
570 
571 template<class CloudType>
573 {
574  cloudReset(cloudCopyPtr_());
575  cloudCopyPtr_.clear();
576 }
577 
578 
579 template<class CloudType>
581 {
582  UTransRef().field() = Zero;
583  UCoeffRef().field() = 0.0;
584 }
585 
586 
587 template<class CloudType>
588 template<class Type>
590 (
592  const DimensionedField<Type, volMesh>& field0,
593  const word& name
594 ) const
595 {
596  const scalar coeff = solution_.relaxCoeff(name);
597  field = field0 + coeff*(field - field0);
598 }
599 
600 
601 template<class CloudType>
602 template<class Type>
604 (
606  const word& name
607 ) const
608 {
609  const scalar coeff = solution_.relaxCoeff(name);
610  field *= coeff;
611 }
612 
613 
614 template<class CloudType>
616 (
617  const MomentumCloud<CloudType>& cloudOldTime
618 )
619 {
620  this->relax(UTrans_(), cloudOldTime.UTrans_(), "U");
621  this->relax(UCoeff_(), cloudOldTime.UCoeff_(), "U");
622 }
623 
624 
625 template<class CloudType>
627 {
628  this->scale(UTrans_(), "U");
629  this->scale(UCoeff_(), "U");
630 }
631 
632 
633 template<class CloudType>
635 {
636  // force calculation of mesh dimensions - needed for parallel runs
637  // with topology change due to lazy evaluation of valid mesh dimensions
638  label nGeometricD = this->mesh().nGeometricD();
639 
640  Info<< nl << "Solving " << nGeometricD << "-D cloud " << this->name()
641  << endl;
642 
643  this->dispersion().cacheFields(true);
644  forces_.cacheFields(true);
645  updateCellOccupancy();
646 
647  pAmbient_ = constProps_.dict().template
648  lookupOrDefault<scalar>("pAmbient", pAmbient_);
649 
650  functions_.preEvolve();
651 }
652 
653 
654 template<class CloudType>
656 {
657  if (solution_.canEvolve())
658  {
659  typename parcelType::trackingData td(*this);
660 
661  solve(*this, td);
662  }
663 }
664 
665 
666 template<class CloudType>
667 template<class TrackCloudType>
669 (
670  TrackCloudType& cloud,
671  typename parcelType::trackingData& td
672 )
673 {
674  CloudType::move(cloud, td);
675 
676  updateCellOccupancy();
677 }
678 
679 
680 template<class CloudType>
682 (
683  const parcelType& p,
684  const polyPatch& pp,
685  vector& nw,
686  vector& Up
687 ) const
688 {
689  p.patchData(this->mesh(), nw, Up);
690  Up /= this->mesh().time().deltaTValue();
691 
692  // If this is a wall patch, then there may be a non-zero tangential velocity
693  // component; the lid velocity in a lid-driven cavity case, for example. We
694  // want the particle to interact with this velocity, so we look it up in the
695  // velocity field and use it to set the wall-tangential component.
696  if (!this->mesh().moving() && isA<wallPolyPatch>(pp))
697  {
698  const label patchi = pp.index();
699  const label patchFacei = pp.whichFace(p.face());
700 
701  // We only want to use the boundary condition value only if it is set
702  // by the boundary condition. If the boundary values are extrapolated
703  // (e.g., slip conditions) then they represent the motion of the fluid
704  // just inside the domain rather than that of the wall itself.
705  if (U_.boundaryField()[patchi].fixesValue())
706  {
707  const vector Uw1 = U_.boundaryField()[patchi][patchFacei];
708  const vector& Uw0 =
709  U_.oldTime().boundaryField()[patchi][patchFacei];
710 
711  const vector Uw = Uw0 + p.stepFraction()*(Uw1 - Uw0);
712 
713  const tensor nnw = nw*nw;
714 
715  Up = (nnw & Up) + Uw - (nnw & Uw);
716  }
717  }
718 }
719 
720 
721 template<class CloudType>
723 {
725 
726  updateCellOccupancy();
727 
728  injectors_.topoChange();
729 
730  cellLengthScale_ = mag(cbrt(this->mesh().V()));
731 }
732 
733 
734 template<class CloudType>
736 {
738 
739  updateCellOccupancy();
740 
741  injectors_.topoChange();
742 
743  cellLengthScale_ = mag(cbrt(this->mesh().V()));
744 }
745 
746 
747 template<class CloudType>
749 {
751 
752  updateCellOccupancy();
753 
754  injectors_.topoChange();
755 
756  cellLengthScale_ = mag(cbrt(this->mesh().V()));
757 }
758 
759 
760 template<class CloudType>
762 {
763  vector linearMomentum = linearMomentumOfSystem();
764  reduce(linearMomentum, sumOp<vector>());
765 
766  scalar linearKineticEnergy = linearKineticEnergyOfSystem();
767  reduce(linearKineticEnergy, sumOp<scalar>());
768 
769  Info<< "Cloud: " << this->name() << nl
770  << " Current number of parcels = "
771  << returnReduce(this->size(), sumOp<label>()) << nl
772  << " Current mass in system = "
773  << returnReduce(massInSystem(), sumOp<scalar>()) << nl
774  << " Linear momentum = "
775  << linearMomentum << nl
776  << " |Linear momentum| = "
777  << mag(linearMomentum) << nl
778  << " Linear kinetic energy = "
779  << linearKineticEnergy << nl;
780 
781  injectors_.info(Info);
782  this->surfaceFilm().info(Info);
783  this->patchInteraction().info(Info);
784 }
785 
786 
787 // ************************************************************************* //
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:419
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
Definition: Cloud.C:558
void deleteLostParticles()
Remove lost particles from cloud and delete.
Definition: Cloud.C:232
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Definition: Cloud.C:449
void move(TrackCloudType &cloud, typename ParticleType::trackingData &td)
Move the particles.
Definition: Cloud.C:275
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:79
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:221
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:160
Base-class for fluid thermodynamic properties.
Definition: fluidThermo.H:57
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
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:251
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)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
static const char nl
Definition: Ostream.H:260
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.
volScalarField & p
const word cloudName(propsDict.lookup("cloudName"))