All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
KinematicParcelI.H
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-2019 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 "mathematicalConstants.H"
27 
28 using namespace Foam::constant::mathematical;
29 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
32 template<class ParcelType>
33 inline
35 :
36  dict_(dictionary::null),
37  parcelTypeId_(dict_, -1),
38  rhoMin_(dict_, 0.0),
39  rho0_(dict_, 0.0),
40  minParcelMass_(dict_, 0.0)
41 {}
42 
43 
44 template<class ParcelType>
46 (
47  const constantProperties& cp
48 )
49 :
50  dict_(cp.dict_),
51  parcelTypeId_(cp.parcelTypeId_),
52  rhoMin_(cp.rhoMin_),
53  rho0_(cp.rho0_),
54  minParcelMass_(cp.minParcelMass_)
55 {}
56 
57 
58 template<class ParcelType>
60 (
61  const dictionary& parentDict
62 )
63 :
64  dict_(parentDict.subOrEmptyDict("constantProperties")),
65  parcelTypeId_(dict_, "parcelTypeId", -1),
66  rhoMin_(dict_, "rhoMin", 1e-15),
67  rho0_(dict_, "rho0"),
68  minParcelMass_(dict_, "minParcelMass", 1e-15)
69 {}
70 
71 
72 template<class ParcelType>
74 (
75  const polyMesh& owner,
76  const barycentric& coordinates,
77  const label celli,
78  const label tetFacei,
79  const label tetPti
80 )
81 :
82  ParcelType(owner, coordinates, celli, tetFacei, tetPti),
83  active_(true),
84  typeId_(-1),
85  nParticle_(0),
86  d_(0.0),
87  dTarget_(0.0),
88  U_(Zero),
89  rho_(0.0),
90  age_(0.0),
91  tTurb_(0.0),
92  UTurb_(Zero)
93 {}
94 
95 
96 template<class ParcelType>
98 (
99  const polyMesh& owner,
100  const vector& position,
101  const label celli
102 )
103 :
104  ParcelType(owner, position, celli),
105  active_(true),
106  typeId_(-1),
107  nParticle_(0),
108  d_(0.0),
109  dTarget_(0.0),
110  U_(Zero),
111  rho_(0.0),
112  age_(0.0),
113  tTurb_(0.0),
114  UTurb_(Zero)
115 {}
116 
117 
118 template<class ParcelType>
120 (
121  const polyMesh& owner,
122  const barycentric& coordinates,
123  const label celli,
124  const label tetFacei,
125  const label tetPti,
126  const label typeId,
127  const scalar nParticle0,
128  const scalar d0,
129  const scalar dTarget0,
130  const vector& U0,
131  const constantProperties& constProps
132 )
133 :
134  ParcelType(owner, coordinates, celli, tetFacei, tetPti),
135  active_(true),
136  typeId_(typeId),
137  nParticle_(nParticle0),
138  d_(d0),
139  dTarget_(dTarget0),
140  U_(U0),
141  rho_(constProps.rho0()),
142  age_(0.0),
143  tTurb_(0.0),
144  UTurb_(Zero)
145 {}
146 
147 
148 // * * * * * * * * * constantProperties Member Functions * * * * * * * * * * //
149 
150 template<class ParcelType>
151 inline const Foam::dictionary&
153 {
154  return dict_;
155 }
156 
157 
158 template<class ParcelType>
159 inline Foam::label
161 {
162  return parcelTypeId_.value();
163 }
164 
165 
166 template<class ParcelType>
167 inline Foam::scalar
169 {
170  return rhoMin_.value();
171 }
172 
173 
174 template<class ParcelType>
175 inline Foam::scalar
177 {
178  return rho0_.value();
179 }
180 
181 
182 template<class ParcelType>
183 inline Foam::scalar
185 {
186  return minParcelMass_.value();
187 }
188 
189 
190 // * * * * * * * KinematicParcel Member Functions * * * * * * * //
191 
192 template<class ParcelType>
194 {
195  return active_;
196 }
197 
198 
199 template<class ParcelType>
201 {
202  return typeId_;
203 }
204 
205 
206 template<class ParcelType>
208 {
209  return nParticle_;
210 }
211 
212 
213 template<class ParcelType>
214 inline Foam::scalar Foam::KinematicParcel<ParcelType>::d() const
215 {
216  return d_;
217 }
218 
219 
220 template<class ParcelType>
222 {
223  return dTarget_;
224 }
225 
226 
227 template<class ParcelType>
229 {
230  return U_;
231 }
232 
233 
234 template<class ParcelType>
235 inline Foam::scalar Foam::KinematicParcel<ParcelType>::rho() const
236 {
237  return rho_;
238 }
239 
240 
241 template<class ParcelType>
242 inline Foam::scalar Foam::KinematicParcel<ParcelType>::age() const
243 {
244  return age_;
245 }
246 
247 
248 template<class ParcelType>
249 inline Foam::scalar Foam::KinematicParcel<ParcelType>::tTurb() const
250 {
251  return tTurb_;
252 }
253 
254 
255 template<class ParcelType>
257 {
258  return UTurb_;
259 }
260 
261 
262 template<class ParcelType>
264 {
265  return active_;
266 }
267 
268 
269 template<class ParcelType>
271 {
272  return typeId_;
273 }
274 
275 
276 template<class ParcelType>
278 {
279  return nParticle_;
280 }
281 
282 
283 template<class ParcelType>
285 {
286  return d_;
287 }
288 
289 
290 template<class ParcelType>
292 {
293  return dTarget_;
294 }
295 
296 
297 template<class ParcelType>
299 {
300  return U_;
301 }
302 
303 
304 template<class ParcelType>
306 {
307  return rho_;
308 }
309 
310 
311 template<class ParcelType>
313 {
314  return age_;
315 }
316 
317 
318 template<class ParcelType>
320 {
321  return tTurb_;
322 }
323 
324 
325 template<class ParcelType>
327 {
328  return UTurb_;
329 }
330 
331 
332 template<class ParcelType>
334 (
335  const trackingData& td
336 ) const
337 {
338  return td.rhoc()*this->mesh().cellVolumes()[this->cell()];
339 }
340 
341 
342 template<class ParcelType>
343 inline Foam::scalar Foam::KinematicParcel<ParcelType>::mass() const
344 {
345  return rho_*volume();
346 }
347 
348 
349 template<class ParcelType>
351 {
352  return 0.1*mass()*sqr(d_);
353 }
354 
355 
356 template<class ParcelType>
357 inline Foam::scalar Foam::KinematicParcel<ParcelType>::volume() const
358 {
359  return volume(d_);
360 }
361 
362 
363 template<class ParcelType>
364 inline Foam::scalar Foam::KinematicParcel<ParcelType>::volume(const scalar d)
365 {
366  return pi/6.0*pow3(d);
367 }
368 
369 
370 template<class ParcelType>
371 inline Foam::scalar Foam::KinematicParcel<ParcelType>::areaP() const
372 {
373  return areaP(d_);
374 }
375 
376 
377 template<class ParcelType>
378 inline Foam::scalar Foam::KinematicParcel<ParcelType>::areaP(const scalar d)
379 {
380  return 0.25*areaS(d);
381 }
382 
383 
384 template<class ParcelType>
385 inline Foam::scalar Foam::KinematicParcel<ParcelType>::areaS() const
386 {
387  return areaS(d_);
388 }
389 
390 
391 template<class ParcelType>
392 inline Foam::scalar Foam::KinematicParcel<ParcelType>::areaS(const scalar d)
393 {
394  return pi*d*d;
395 }
396 
397 
398 template<class ParcelType>
399 inline Foam::scalar Foam::KinematicParcel<ParcelType>::Re
400 (
401  const trackingData& td
402 ) const
403 {
404  return Re(td.rhoc(), U_, td.Uc(), d_, td.muc());
405 }
406 
407 
408 template<class ParcelType>
409 inline Foam::scalar Foam::KinematicParcel<ParcelType>::Re
410 (
411  const scalar rhoc,
412  const vector& U,
413  const vector& Uc,
414  const scalar d,
415  const scalar muc
416 )
417 {
418  return rhoc*mag(U - Uc)*d/max(muc, rootVSmall);
419 }
420 
421 
422 template<class ParcelType>
423 inline Foam::scalar Foam::KinematicParcel<ParcelType>::We
424 (
425  const trackingData& td,
426  const scalar sigma
427 ) const
428 {
429  return We(td.rhoc(), U_, td.Uc(), d_, sigma);
430 }
431 
432 
433 template<class ParcelType>
434 inline Foam::scalar Foam::KinematicParcel<ParcelType>::We
435 (
436  const scalar rhoc,
437  const vector& U,
438  const vector& Uc,
439  const scalar d,
440  const scalar sigma
441 )
442 {
443  return rhoc*magSqr(U - Uc)*d/max(sigma, rootVSmall);
444 }
445 
446 
447 template<class ParcelType>
448 inline Foam::scalar Foam::KinematicParcel<ParcelType>::Eo
449 (
450  const trackingData& td,
451  const scalar sigma
452 ) const
453 {
454  return Eo(td.g(), rho_, td.rhoc(), U_, d_, sigma);
455 }
456 
457 
458 template<class ParcelType>
459 inline Foam::scalar Foam::KinematicParcel<ParcelType>::Eo
460 (
461  const vector& g,
462  const scalar rho,
463  const scalar rhoc,
464  const vector& U,
465  const scalar d,
466  const scalar sigma
467 )
468 {
469  const vector dir = U/max(mag(U), rootVSmall);
470  return mag(g & dir)*mag(rho - rhoc)*sqr(d)/max(sigma, rootVSmall);
471 }
472 
473 
474 // ************************************************************************* //
scalar massCell(const trackingData &td) const
Cell owner mass.
const vector & U() const
Return const access to velocity.
scalar tTurb() const
Return const access to time spent in turbulent eddy.
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
scalar rhoMin() const
Return const access to the minimum density.
const vector & UTurb() const
Return const access to turbulent velocity fluctuation.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
scalar nParticle() const
Return const access to number of particles.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionedScalar & sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
Class to hold kinematic particle constant properties.
scalar mass() const
Particle mass.
scalar d() const
Return const access to diameter.
scalar tTurb_
Time spent in turbulent eddy [s].
scalar age_
Age [s].
scalar rho() const
Return const access to density.
scalar rho0() const
Return const access to the particle density.
scalar areaP() const
Particle projected area.
KinematicParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from mesh, coordinates and topology.
scalar dTarget() const
Return const access to target diameter.
const Type & value() const
Return the value.
label parcelTypeId() const
Return const access to the parcel type id.
scalar dTarget_
Target diameter [m].
scalar age() const
Return const access to the age.
dynamicFvMesh & mesh
scalar Eo(const trackingData &td, const scalar sigma) const
Eotvos number.
scalar rho_
Density [kg/m^3].
scalar momentOfInertia() const
Particle moment of inertia around diameter axis.
mathematical constants.
bool active_
Active flag - tracking inactive when active = false.
const dictionary dict_
Constant properties dictionary.
label typeId_
Parcel type id.
vector UTurb_
Turbulent velocity fluctuation [m/s].
static const zero Zero
Definition: zero.H:97
bool active() const
Return const access to active flag.
dimensioned< scalar > magSqr(const dimensioned< Type > &)
scalar rhoc() const
Return the continuous phase density.
scalar muc() const
Return the continuous phase viscosity.
scalar nParticle_
Number of particles in Parcel.
scalar minParcelMass() const
Return const access to the minimum parcel mass.
dimensionedScalar pow3(const dimensionedScalar &ds)
vector U_
Velocity of Parcel [m/s].
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:56
scalar volume() const
Particle volume.
dimensioned< scalar > mag(const dimensioned< Type > &)
scalar Re(const trackingData &td) const
Reynolds number.
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:105
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
label typeId() const
Return const access to type id.
scalar areaS() const
Particle surface area.
const vector & Uc() const
Return the continuous phase velocity.
scalar d_
Diameter [m].
scalar We(const trackingData &td, const scalar sigma) const
Weber number.
dictionary subOrEmptyDict(const word &, const bool mustRead=false) const
Find and return a sub-dictionary as a copy, or.
Definition: dictionary.C:969
const dictionary & dict() const
Return const access to the constant properties dictionary.