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-2021 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  IOobject cloudsIO
45  (
46  cloudsName,
47  mesh_.time().constant(),
48  mesh_,
51  );
52 
53  if (cloudsIO.typeHeaderOk<wordGlobalIOList>(false))
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  IOobject cloudIO
67  (
68  "cloudProperties",
69  mesh_.time().constant(),
70  mesh_,
73  );
74 
75  if (cloudIO.typeHeaderOk<IOdictionary>(false))
76  {
77  this->setSize(1);
78 
79  this->set(0, parcelCloud::New("cloud", args ...));
80  }
81  else
82  {
83  Info<< "Clouds not active: Neither " << cloudsIO.localObjectPath()
84  << " nor " << cloudIO.localObjectPath() << " found" << endl;
85  }
86  }
87 }
88 
89 
90 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
91 
93 (
94  const volScalarField& rho,
95  const volVectorField& U,
96  const volScalarField& mu,
97  const dimensionedVector& g
98 )
99 :
101  mesh_(rho.mesh())
102 {
103  initialise(rho, U, mu, g);
104 }
105 
106 
108 (
109  const volScalarField& rho,
110  const volVectorField& U,
111  const dimensionedVector& g,
112  const fluidThermo& carrierThermo
113 )
114 :
116  mesh_(rho.mesh())
117 {
118  initialise(rho, U, g, carrierThermo);
119 }
120 
121 
122 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
123 
125 {}
126 
127 
128 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
129 
131 {
132  tmp<volScalarField> ttheta
133  (
135  (
136  cloudsName + ":theta",
137  mesh_,
139  extrapolatedCalculatedFvPatchScalarField::typeName
140  )
141  );
142  forAll(*this, i)
143  {
144  ttheta.ref() += operator[](i).theta();
145  }
146  return ttheta;
147 }
148 
149 
151 (
152  const volVectorField& U
153 ) const
154 {
156  forAll(*this, i)
157  {
158  tSU.ref() += operator[](i).SU(U);
159  }
160  return tSU;
161 }
162 
163 
165 {
167  (
169  (
170  cloudsName + ":UTrans",
171  mesh_,
173  )
174  );
175  forAll(*this, i)
176  {
177  tUTrans.ref() += operator[](i).UTrans();
178  }
179  return tUTrans;
180 }
181 
182 
184 {
186  (
188  (
189  cloudsName + ":UCoeff",
190  mesh_,
192  )
193  );
194  forAll(*this, i)
195  {
196  tUCoeff.ref() += operator[](i).UCoeff();
197  }
198  return tUCoeff;
199 }
200 
201 
203 (
204  const volScalarField& hs
205 ) const
206 {
208  forAll(*this, i)
209  {
210  tSh.ref() += operator[](i).Sh(hs);
211  }
212  return tSh;
213 }
214 
215 
217 {
219  (
221  (
222  cloudsName + ":hsTrans",
223  mesh_,
225  )
226  );
227  forAll(*this, i)
228  {
229  thsTrans.ref() += operator[](i).hsTrans();
230  }
231  return thsTrans;
232 }
233 
234 
236 {
238  (
240  (
241  cloudsName + ":hsCoeff",
242  mesh_,
244  )
245  );
246  forAll(*this, i)
247  {
248  thsCoeff.ref() += operator[](i).hsCoeff();
249  }
250  return thsCoeff;
251 }
252 
253 
255 {
257  (
259  (
260  cloudsName + ":radiation:Ep",
261  mesh_,
263  )
264  );
265  forAll(*this, i)
266  {
267  tEp.ref() += operator[](i).Ep();
268  }
269  return tEp;
270 }
271 
272 
274 {
276  (
278  (
279  cloudsName + ":radiation:ap",
280  mesh_,
282  )
283  );
284  forAll(*this, i)
285  {
286  tap.ref() += operator[](i).ap();
287  }
288  return tap;
289 
290 }
291 
292 
294 {
295  tmp<volScalarField> tsigmap
296  (
298  (
299  cloudsName + ":radiation:sigmap",
300  mesh_,
302  )
303  );
304  forAll(*this, i)
305  {
306  tsigmap.ref() += operator[](i).sigmap();
307  }
308  return tsigmap;
309 }
310 
311 
313 (
314  const label speciei,
315  const volScalarField& Yi
316 ) const
317 {
319  forAll(*this, i)
320  {
321  tSYi.ref() += operator[](i).SYi(speciei, Yi);
322  }
323  return tSYi;
324 }
325 
326 
328 (
329  const volScalarField& rho
330 ) const
331 {
333  forAll(*this, i)
334  {
335  tSrho.ref() += operator[](i).Srho(rho);
336  }
337  return tSrho;
338 }
339 
340 
342 {
344  (
346  (
347  cloudsName + ":Srho",
348  mesh_,
350  )
351  );
352  forAll(*this, i)
353  {
354  tSrho.ref() += operator[](i).Srho();
355  }
356  return tSrho;
357 }
358 
359 
361 {
362  forAll(*this, i)
363  {
364  operator[](i).info();
365  }
366 }
367 
368 
370 {
371  forAll(*this, i)
372  {
373  operator[](i).evolve();
374  }
375 }
376 
377 
379 {
380  forAll(*this, i)
381  {
382  operator[](i).storeGlobalPositions();
383  }
384 }
385 
386 
387 // ************************************************************************* //
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
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
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:239
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.
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