parcelCloudList.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) 2020-2022 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 "parcelCloudList.H"
28 #include "wordIOList.H"
29 #include "fvMatrices.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
34 
36 
38 
39 
40 // * * * * * * * * * * Private Static Member Functions * * * * * * * * * * * //
41 
42 Foam::wordList Foam::parcelCloudList::cloudNames(const objectRegistry& db)
43 {
44  typeIOobject<wordGlobalIOList> cloudNamesIO
45  (
47  db.time().constant(),
48  db,
51  );
52 
53  if (cloudNamesIO.headerOk())
54  {
55  return wordGlobalIOList(cloudNamesIO);
56  }
57 
58  typeIOobject<IOdictionary> cloudIO
59  (
60  defaultCloudName + "Properties",
61  db.time().constant(),
62  db,
65  );
66 
67  if (cloudIO.headerOk())
68  {
69  return wordList(1, defaultCloudName);
70  }
71 
73  << "Cloud properties were not found in either "
74  << cloudNamesIO.relativeObjectPath() << " or "
75  << cloudIO.relativeObjectPath() << exit(FatalError);
76 
77  return wordList::null();
78 }
79 
80 
81 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
82 
84 (
85  const wordList& cloudNames,
86  const volScalarField& rho,
87  const volVectorField& U,
88  const volScalarField& mu,
89  const dimensionedVector& g
90 )
91 :
93  mesh_(rho.mesh())
94 {
95  this->setSize(cloudNames.size());
96 
97  forAll(cloudNames, i)
98  {
99  this->set(i, parcelCloud::New(cloudNames[i], rho, U, mu, g));
100  }
101 }
102 
103 
105 (
106  const wordList& cloudNames,
107  const volScalarField& rho,
108  const volVectorField& U,
109  const dimensionedVector& g,
110  const fluidThermo& carrierThermo
111 )
112 :
113  PtrList<parcelCloud>(),
114  mesh_(rho.mesh())
115 {
116  this->setSize(cloudNames.size());
117 
118  forAll(cloudNames, i)
119  {
120  this->set(i, parcelCloud::New(cloudNames[i], rho, U, g, carrierThermo));
121  }
122 }
123 
124 
126 (
127  const volScalarField& rho,
128  const volVectorField& U,
129  const volScalarField& mu,
130  const dimensionedVector& g
131 )
132 :
133  parcelCloudList(cloudNames(rho.mesh()), rho, U, mu, g)
134 {}
135 
136 
138 (
139  const volScalarField& rho,
140  const volVectorField& U,
141  const dimensionedVector& g,
142  const fluidThermo& carrierThermo
143 )
144 :
145  parcelCloudList(cloudNames(rho.mesh()), rho, U, g, carrierThermo)
146 {}
147 
148 
149 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
150 
152 {}
153 
154 
155 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
156 
158 {
159  tmp<volScalarField> ttheta
160  (
162  (
163  cloudNamesName + ":theta",
164  mesh_,
166  extrapolatedCalculatedFvPatchScalarField::typeName
167  )
168  );
169  forAll(*this, i)
170  {
171  ttheta.ref() += operator[](i).theta();
172  }
173  return ttheta;
174 }
175 
176 
178 (
179  const volVectorField& U
180 ) const
181 {
183  forAll(*this, i)
184  {
185  tSU.ref() += operator[](i).SU(U);
186  }
187  return tSU;
188 }
189 
190 
192 {
194  (
196  (
197  cloudNamesName + ":UTrans",
198  mesh_,
200  )
201  );
202  forAll(*this, i)
203  {
204  tUTrans.ref() += operator[](i).UTrans();
205  }
206  return tUTrans;
207 }
208 
209 
211 {
213  (
215  (
216  cloudNamesName + ":UCoeff",
217  mesh_,
219  )
220  );
221  forAll(*this, i)
222  {
223  tUCoeff.ref() += operator[](i).UCoeff();
224  }
225  return tUCoeff;
226 }
227 
228 
230 (
231  const volScalarField& hs
232 ) const
233 {
235  forAll(*this, i)
236  {
237  tSh.ref() += operator[](i).Sh(hs);
238  }
239  return tSh;
240 }
241 
242 
244 {
246  (
248  (
249  cloudNamesName + ":hsTrans",
250  mesh_,
252  )
253  );
254  forAll(*this, i)
255  {
256  thsTrans.ref() += operator[](i).hsTrans();
257  }
258  return thsTrans;
259 }
260 
261 
263 {
265  (
267  (
268  cloudNamesName + ":hsCoeff",
269  mesh_,
271  )
272  );
273  forAll(*this, i)
274  {
275  thsCoeff.ref() += operator[](i).hsCoeff();
276  }
277  return thsCoeff;
278 }
279 
280 
282 {
284  (
286  (
287  cloudNamesName + ":radiation:Ep",
288  mesh_,
290  )
291  );
292  forAll(*this, i)
293  {
294  tEp.ref() += operator[](i).Ep();
295  }
296  return tEp;
297 }
298 
299 
301 {
303  (
305  (
306  cloudNamesName + ":radiation:ap",
307  mesh_,
309  )
310  );
311  forAll(*this, i)
312  {
313  tap.ref() += operator[](i).ap();
314  }
315  return tap;
316 
317 }
318 
319 
321 {
322  tmp<volScalarField> tsigmap
323  (
325  (
326  cloudNamesName + ":radiation:sigmap",
327  mesh_,
329  )
330  );
331  forAll(*this, i)
332  {
333  tsigmap.ref() += operator[](i).sigmap();
334  }
335  return tsigmap;
336 }
337 
338 
340 (
341  const label speciei,
342  const volScalarField& Yi
343 ) const
344 {
346  forAll(*this, i)
347  {
348  tSYi.ref() += operator[](i).SYi(speciei, Yi);
349  }
350  return tSYi;
351 }
352 
353 
355 (
356  const volScalarField& rho
357 ) const
358 {
360  forAll(*this, i)
361  {
362  tSrho.ref() += operator[](i).Srho(rho);
363  }
364  return tSrho;
365 }
366 
367 
369 {
371  (
373  (
374  cloudNamesName + ":Srho",
375  mesh_,
377  )
378  );
379  forAll(*this, i)
380  {
381  tSrho.ref() += operator[](i).Srho();
382  }
383  return tSrho;
384 }
385 
386 
388 {
389  forAll(*this, i)
390  {
391  operator[](i).info();
392  }
393 }
394 
395 
397 {
398  forAll(*this, i)
399  {
400  operator[](i).evolve();
401  }
402 }
403 
404 
406 {
407  forAll(*this, i)
408  {
409  operator[](i).storeGlobalPositions();
410  }
411 }
412 
413 
415 {
416  forAll(*this, i)
417  {
418  operator[](i).topoChange(map);
419  }
420 }
421 
422 
424 {
425  forAll(*this, i)
426  {
427  operator[](i).mapMesh(map);
428  }
429 }
430 
431 
433 {
434  forAll(*this, i)
435  {
436  operator[](i).distribute(map);
437  }
438 }
439 
440 
441 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
static tmp< DimensionedField< Type, GeoMesh > > New(const word &name, const Mesh &mesh, const dimensionSet &, const Field< Type > &)
Return a temporary field constructed from name, mesh,.
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &)
Return a temporary field constructed from name,.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
static const List< word > & null()
Return a null List.
Definition: ListI.H:118
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
bool set(const label) const
Is element set.
Definition: PtrListI.H:65
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
Definition: PtrList.C:131
Base-class for fluid thermodynamic properties.
Definition: fluidThermo.H:57
List of parcel clouds, with the same interface as an individual parcel cloud. This is the object that...
static const wordList defaultCloudNames
The default cloud names (i.e., a list of length one with the.
tmp< fvScalarMatrix > Sh(const volScalarField &hs) const
Return sensible enthalpy source term [J/s].
tmp< volScalarField::Internal > UCoeff() const
Momentum transfer coefficient [kg].
tmp< volVectorField::Internal > UTrans() const
Momentum transfer [kg m/s].
const tmp< volScalarField > theta() const
Return the particle volume fraction field.
~parcelCloudList()
Destructor.
void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
tmp< fvScalarMatrix > SYi(const label speciei, const volScalarField &Yi) const
Return mass source term for specie [kg/s].
void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
tmp< volScalarField > sigmap() const
Return equivalent particulate scattering factor [1/m].
void evolve()
Evolve the cloud.
parcelCloudList(const wordList &cloudNames, const volScalarField &rho, const volVectorField &U, const volScalarField &mu, const dimensionedVector &g)
Construct specified clouds with given carrier fields.
tmp< volScalarField::Internal > hsCoeff() const
Sensible enthalpy transfer coefficient [J/K].
static const word defaultCloudName
The default cloud name.
void info()
Print cloud information.
tmp< fvVectorMatrix > SU(const volVectorField &U) const
Return momentum source term [kg m/s^2].
tmp< volScalarField::Internal > Srho() const
Return total mass source [kg/m^3/s].
void storeGlobalPositions()
Call this before a topology change. Stores the particles global.
tmp< volScalarField > Ep() const
Return equivalent particulate emission [kg/m/s^3].
tmp< volScalarField::Internal > hsTrans() const
Sensible enthalpy transfer [J].
static const word cloudNamesName
The name of the clouds file in which multiple cloud names are.
tmp< volScalarField > ap() const
Return equivalent particulate absorption [1/m].
Virtual abstract base class for parcel clouds. As parcelCloudBase but with additional virtualisation ...
Definition: parcelCloud.H:57
static autoPtr< parcelCloud > New(const word &name, const volScalarField &rho, const volVectorField &U, const volScalarField &mu, const dimensionedVector &g)
Selectors.
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
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
A special matrix type and solver, designed for finite volume solutions of scalar equations.
U
Definition: pEqn.H:72
const dimensionedScalar mu
Atomic mass unit.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
static const zero Zero
Definition: zero.H:97
List< word > wordList
A List of words.
Definition: fileName.H:54
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
const dimensionSet dimEnergy
dimensionedScalar pow3(const dimensionedScalar &ds)
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
const dimensionSet dimless
GlobalIOList< word > wordGlobalIOList
Definition: wordIOList.H:44
const dimensionSet dimLength
const dimensionSet dimTemperature
const dimensionSet dimAcceleration
const dimensionSet dimTime
const dimensionSet dimDensity
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:45
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
const dimensionSet dimMass
const dimensionSet dimVelocity
error FatalError