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 "GlobalIOLists.H"
29 #include "fvMatrices.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
34 
35 
36 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
37 
38 template<class ... Args>
39 void Foam::parcelCloudList::initialise
40 (
41  const Args& ... args
42 )
43 {
44  typeIOobject<wordGlobalIOList> cloudsIO
45  (
46  cloudsName,
47  mesh_.time().constant(),
48  mesh_,
51  );
52 
53  if (cloudsIO.headerOk())
54  {
55  wordGlobalIOList cloudNames(cloudsIO);
56 
57  this->setSize(cloudNames.size());
58 
59  forAll(cloudNames, i)
60  {
61  this->set(i, parcelCloud::New(cloudNames[i], args ...));
62  }
63  }
64  else
65  {
66  typeIOobject<IOdictionary> cloudIO
67  (
68  "cloudProperties",
69  mesh_.time().constant(),
70  mesh_,
73  );
74 
75  if (cloudIO.headerOk())
76  {
77  this->setSize(1);
78 
79  this->set(0, parcelCloud::New("cloud", args ...));
80  }
81  else
82  {
83  Info<< "Clouds not active: Neither "
84  << cloudsIO.relativeObjectPath()
85  << " nor " << cloudIO.relativeObjectPath() << " found" << endl;
86  }
87  }
88 }
89 
90 
91 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
92 
94 (
95  const volScalarField& rho,
96  const volVectorField& U,
97  const volScalarField& mu,
98  const dimensionedVector& g
99 )
100 :
102  mesh_(rho.mesh())
103 {
104  initialise(rho, U, mu, g);
105 }
106 
107 
109 (
110  const volScalarField& rho,
111  const volVectorField& U,
112  const dimensionedVector& g,
113  const fluidThermo& carrierThermo
114 )
115 :
117  mesh_(rho.mesh())
118 {
119  initialise(rho, U, g, carrierThermo);
120 }
121 
122 
123 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
124 
126 {}
127 
128 
129 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
130 
132 {
133  tmp<volScalarField> ttheta
134  (
136  (
137  cloudsName + ":theta",
138  mesh_,
140  extrapolatedCalculatedFvPatchScalarField::typeName
141  )
142  );
143  forAll(*this, i)
144  {
145  ttheta.ref() += operator[](i).theta();
146  }
147  return ttheta;
148 }
149 
150 
152 (
153  const volVectorField& U
154 ) const
155 {
157  forAll(*this, i)
158  {
159  tSU.ref() += operator[](i).SU(U);
160  }
161  return tSU;
162 }
163 
164 
166 {
168  (
170  (
171  cloudsName + ":UTrans",
172  mesh_,
174  )
175  );
176  forAll(*this, i)
177  {
178  tUTrans.ref() += operator[](i).UTrans();
179  }
180  return tUTrans;
181 }
182 
183 
185 {
187  (
189  (
190  cloudsName + ":UCoeff",
191  mesh_,
193  )
194  );
195  forAll(*this, i)
196  {
197  tUCoeff.ref() += operator[](i).UCoeff();
198  }
199  return tUCoeff;
200 }
201 
202 
204 (
205  const volScalarField& hs
206 ) const
207 {
209  forAll(*this, i)
210  {
211  tSh.ref() += operator[](i).Sh(hs);
212  }
213  return tSh;
214 }
215 
216 
218 {
220  (
222  (
223  cloudsName + ":hsTrans",
224  mesh_,
226  )
227  );
228  forAll(*this, i)
229  {
230  thsTrans.ref() += operator[](i).hsTrans();
231  }
232  return thsTrans;
233 }
234 
235 
237 {
239  (
241  (
242  cloudsName + ":hsCoeff",
243  mesh_,
245  )
246  );
247  forAll(*this, i)
248  {
249  thsCoeff.ref() += operator[](i).hsCoeff();
250  }
251  return thsCoeff;
252 }
253 
254 
256 {
258  (
260  (
261  cloudsName + ":radiation:Ep",
262  mesh_,
264  )
265  );
266  forAll(*this, i)
267  {
268  tEp.ref() += operator[](i).Ep();
269  }
270  return tEp;
271 }
272 
273 
275 {
277  (
279  (
280  cloudsName + ":radiation:ap",
281  mesh_,
283  )
284  );
285  forAll(*this, i)
286  {
287  tap.ref() += operator[](i).ap();
288  }
289  return tap;
290 
291 }
292 
293 
295 {
296  tmp<volScalarField> tsigmap
297  (
299  (
300  cloudsName + ":radiation:sigmap",
301  mesh_,
303  )
304  );
305  forAll(*this, i)
306  {
307  tsigmap.ref() += operator[](i).sigmap();
308  }
309  return tsigmap;
310 }
311 
312 
314 (
315  const label speciei,
316  const volScalarField& Yi
317 ) const
318 {
320  forAll(*this, i)
321  {
322  tSYi.ref() += operator[](i).SYi(speciei, Yi);
323  }
324  return tSYi;
325 }
326 
327 
329 (
330  const volScalarField& rho
331 ) const
332 {
334  forAll(*this, i)
335  {
336  tSrho.ref() += operator[](i).Srho(rho);
337  }
338  return tSrho;
339 }
340 
341 
343 {
345  (
347  (
348  cloudsName + ":Srho",
349  mesh_,
351  )
352  );
353  forAll(*this, i)
354  {
355  tSrho.ref() += operator[](i).Srho();
356  }
357  return tSrho;
358 }
359 
360 
362 {
363  forAll(*this, i)
364  {
365  operator[](i).info();
366  }
367 }
368 
369 
371 {
372  forAll(*this, i)
373  {
374  operator[](i).evolve();
375  }
376 }
377 
378 
380 {
381  forAll(*this, i)
382  {
383  operator[](i).storeGlobalPositions();
384  }
385 }
386 
387 
389 {
390  forAll(*this, i)
391  {
392  operator[](i).distribute(map);
393  }
394 }
395 
396 
397 // ************************************************************************* //
const tmp< volScalarField > theta() const
Return the particle volume fraction field.
tmp< fvScalarMatrix > SYi(const label speciei, const volScalarField &Yi) const
Return mass source term for specie [kg/s].
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
static tmp< DimensionedField< Type, GeoMesh > > New(const word &name, const Mesh &mesh, const dimensionSet &)
Return a temporary field constructed from name, mesh.
tmp< volScalarField::Internal > hsCoeff() const
Sensible enthalpy transfer coefficient [J/K].
static autoPtr< parcelCloud > New(const word &name, const volScalarField &rho, const volVectorField &U, const volScalarField &mu, const dimensionedVector &g)
Selectors.
parcelCloudList(const volScalarField &rho, const volVectorField &U, const volScalarField &mu, const dimensionedVector &g)
Construct with given carrier fields.
tmp< volScalarField > sigmap() const
Return equivalent particulate scattering factor [1/m].
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
~parcelCloudList()
Destructor.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
const T & operator[](const label) const
Return element const reference.
Definition: UPtrListI.H:96
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
const dimensionSet dimless
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:372
tmp< fvVectorMatrix > SU(const volVectorField &U) const
Return momentum source term [kg m/s^2].
const dimensionSet dimLength
tmp< volVectorField::Internal > UTrans() const
Momentum transfer [kg m/s].
const dimensionSet dimTime
void info()
Print cloud information.
void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
const dimensionSet dimAcceleration
GlobalIOList< word > wordGlobalIOList
Definition: GlobalIOLists.H:49
A class for handling words, derived from string.
Definition: word.H:59
tmp< volScalarField::Internal > UCoeff() const
Momentum transfer coefficient [kg].
Base-class for fluid thermodynamic properties.
Definition: fluidThermo.H:53
const dimensionSet dimDensity
const word & constant() const
Return constant name.
Definition: TimePaths.H:123
void storeGlobalPositions()
Call this before a topology change. Stores the particles global.
static const zero Zero
Definition: zero.H:97
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
Definition: PtrList.C:131
tmp< volScalarField::Internal > Srho() const
Return total mass source [kg/m^3/s].
const dimensionSet dimVelocity
const Mesh & mesh() const
Return mesh.
tmp< volScalarField > ap() const
Return equivalent particulate absorption [1/m].
const dimensionSet dimEnergy
const dimensionSet dimMass
tmp< volScalarField::Internal > hsTrans() const
Sensible enthalpy transfer [J].
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:45
A special matrix type and solver, designed for finite volume solutions of scalar equations.
messageStream Info
void evolve()
Evolve the cloud.
static const word cloudsName
The name of the clouds.
A class for managing temporary objects.
Definition: PtrList.H:53
tmp< volScalarField > Ep() const
Return equivalent particulate emission [kg/m/s^3].
tmp< fvScalarMatrix > Sh(const volScalarField &hs) const
Return sensible enthalpy source term [J/s].
const dimensionSet dimTemperature