KinematicParcelI.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2017 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  rhoc_(0.0),
94  Uc_(Zero),
95  muc_(0.0)
96 {}
97 
98 
99 template<class ParcelType>
101 (
102  const polyMesh& owner,
103  const vector& position,
104  const label celli
105 )
106 :
107  ParcelType(owner, position, celli),
108  active_(true),
109  typeId_(-1),
110  nParticle_(0),
111  d_(0.0),
112  dTarget_(0.0),
113  U_(Zero),
114  rho_(0.0),
115  age_(0.0),
116  tTurb_(0.0),
117  UTurb_(Zero),
118  rhoc_(0.0),
119  Uc_(Zero),
120  muc_(0.0)
121 {}
122 
123 
124 template<class ParcelType>
126 (
127  const polyMesh& owner,
128  const barycentric& coordinates,
129  const label celli,
130  const label tetFacei,
131  const label tetPti,
132  const label typeId,
133  const scalar nParticle0,
134  const scalar d0,
135  const scalar dTarget0,
136  const vector& U0,
137  const constantProperties& constProps
138 )
139 :
140  ParcelType(owner, coordinates, celli, tetFacei, tetPti),
141  active_(true),
142  typeId_(typeId),
143  nParticle_(nParticle0),
144  d_(d0),
145  dTarget_(dTarget0),
146  U_(U0),
147  rho_(constProps.rho0()),
148  age_(0.0),
149  tTurb_(0.0),
150  UTurb_(Zero),
151  rhoc_(0.0),
152  Uc_(Zero),
153  muc_(0.0)
154 {}
155 
156 
157 // * * * * * * * * * constantProperties Member Functions * * * * * * * * * * //
158 
159 template<class ParcelType>
160 inline const Foam::dictionary&
162 {
163  return dict_;
164 }
165 
166 
167 template<class ParcelType>
168 inline Foam::label
170 {
171  return parcelTypeId_.value();
172 }
173 
174 
175 template<class ParcelType>
176 inline Foam::scalar
178 {
179  return rhoMin_.value();
180 }
181 
182 
183 template<class ParcelType>
184 inline Foam::scalar
186 {
187  return rho0_.value();
188 }
189 
190 
191 template<class ParcelType>
192 inline Foam::scalar
194 {
195  return minParcelMass_.value();
196 }
197 
198 
199 // * * * * * * * KinematicParcel Member Functions * * * * * * * //
200 
201 template<class ParcelType>
203 {
204  return active_;
205 }
206 
207 
208 template<class ParcelType>
210 {
211  return typeId_;
212 }
213 
214 
215 template<class ParcelType>
217 {
218  return nParticle_;
219 }
220 
221 
222 template<class ParcelType>
223 inline Foam::scalar Foam::KinematicParcel<ParcelType>::d() const
224 {
225  return d_;
226 }
227 
228 
229 template<class ParcelType>
231 {
232  return dTarget_;
233 }
234 
235 
236 template<class ParcelType>
238 {
239  return U_;
240 }
241 
242 
243 template<class ParcelType>
244 inline Foam::scalar Foam::KinematicParcel<ParcelType>::rho() const
245 {
246  return rho_;
247 }
248 
249 
250 template<class ParcelType>
251 inline Foam::scalar Foam::KinematicParcel<ParcelType>::age() const
252 {
253  return age_;
254 }
255 
256 
257 template<class ParcelType>
258 inline Foam::scalar Foam::KinematicParcel<ParcelType>::tTurb() const
259 {
260  return tTurb_;
261 }
262 
263 
264 template<class ParcelType>
266 {
267  return UTurb_;
268 }
269 
270 
271 template<class ParcelType>
272 inline Foam::scalar Foam::KinematicParcel<ParcelType>::rhoc() const
273 {
274  return rhoc_;
275 }
276 
277 
278 template<class ParcelType>
280 {
281  return Uc_;
282 }
283 
284 
285 template<class ParcelType>
286 inline Foam::scalar Foam::KinematicParcel<ParcelType>::muc() const
287 {
288  return muc_;
289 }
290 
291 
292 template<class ParcelType>
294 {
295  return active_;
296 }
297 
298 
299 template<class ParcelType>
301 {
302  return typeId_;
303 }
304 
305 
306 template<class ParcelType>
308 {
309  return nParticle_;
310 }
311 
312 
313 template<class ParcelType>
315 {
316  return d_;
317 }
318 
319 
320 template<class ParcelType>
322 {
323  return dTarget_;
324 }
325 
326 
327 template<class ParcelType>
329 {
330  return U_;
331 }
332 
333 
334 template<class ParcelType>
336 {
337  return rho_;
338 }
339 
340 
341 template<class ParcelType>
343 {
344  return age_;
345 }
346 
347 
348 template<class ParcelType>
350 {
351  return tTurb_;
352 }
353 
354 
355 template<class ParcelType>
357 {
358  return UTurb_;
359 }
360 
361 
362 template<class ParcelType>
364 (
365  const label celli
366 ) const
367 {
368  return rhoc_*this->mesh().cellVolumes()[celli];
369 }
370 
371 
372 template<class ParcelType>
373 inline Foam::scalar Foam::KinematicParcel<ParcelType>::mass() const
374 {
375  return rho_*volume();
376 }
377 
378 
379 template<class ParcelType>
381 {
382  return 0.1*mass()*sqr(d_);
383 }
384 
385 
386 template<class ParcelType>
387 inline Foam::scalar Foam::KinematicParcel<ParcelType>::volume() const
388 {
389  return volume(d_);
390 }
391 
392 
393 template<class ParcelType>
394 inline Foam::scalar Foam::KinematicParcel<ParcelType>::volume(const scalar d)
395 {
396  return pi/6.0*pow3(d);
397 }
398 
399 
400 template<class ParcelType>
401 inline Foam::scalar Foam::KinematicParcel<ParcelType>::areaP() const
402 {
403  return areaP(d_);
404 }
405 
406 
407 template<class ParcelType>
408 inline Foam::scalar Foam::KinematicParcel<ParcelType>::areaP(const scalar d)
409 {
410  return 0.25*areaS(d);
411 }
412 
413 
414 template<class ParcelType>
415 inline Foam::scalar Foam::KinematicParcel<ParcelType>::areaS() const
416 {
417  return areaS(d_);
418 }
419 
420 
421 template<class ParcelType>
422 inline Foam::scalar Foam::KinematicParcel<ParcelType>::areaS(const scalar d)
423 {
424  return pi*d*d;
425 }
426 
427 
428 template<class ParcelType>
429 inline Foam::scalar Foam::KinematicParcel<ParcelType>::Re
430 (
431  const vector& U,
432  const scalar d,
433  const scalar rhoc,
434  const scalar muc
435 ) const
436 {
437  return rhoc*mag(U - Uc_)*d/(muc + ROOTVSMALL);
438 }
439 
440 
441 template<class ParcelType>
442 inline Foam::scalar Foam::KinematicParcel<ParcelType>::We
443 (
444  const vector& U,
445  const scalar d,
446  const scalar rhoc,
447  const scalar sigma
448 ) const
449 {
450  return rhoc*magSqr(U - Uc_)*d/(sigma + ROOTVSMALL);
451 }
452 
453 
454 template<class ParcelType>
455 inline Foam::scalar Foam::KinematicParcel<ParcelType>::Eo
456 (
457  const vector& a,
458  const scalar d,
459  const scalar sigma
460 ) const
461 {
462  vector dir = U_/(mag(U_) + ROOTVSMALL);
463  return mag(a & dir)*(rho_ - rhoc_)*sqr(d)/(sigma + ROOTVSMALL);
464 }
465 
466 
467 // ************************************************************************* //
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 double e
Elementary charge.
Definition: doubleFloat.H:78
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:137
scalar muc() const
Return const access to carrier viscosity [Pa.s].
scalar nParticle() const
Return const access to number of particles.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
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.
vector Uc_
Velocity [m/s].
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 vector &a, const scalar d, const scalar sigma) const
Eotvos number.
scalar rho_
Density [kg/m3].
scalar momentOfInertia() const
Particle moment of inertia around diameter axis.
mathematical constants.
scalar muc_
Viscosity [Pa.s].
bool active_
Active flag - tracking inactive when active = false.
const dictionary dict_
Constant properties dictionary.
scalar massCell(const label celli) const
Cell owner mass.
label typeId_
Parcel type id.
vector UTurb_
Turbulent velocity fluctuation [m/s].
static const zero Zero
Definition: zero.H:91
bool active() const
Return const access to active flag.
const vector & Uc() const
Return const access to carrier velocity [m/s].
dimensioned< scalar > magSqr(const dimensioned< Type > &)
scalar nParticle_
Number of particles in Parcel.
scalar minParcelMass() const
Return const access to the minimum parcel mass.
scalar rhoc() const
Return const access to carrier density [kg/m3].
dimensionedScalar pow3(const dimensionedScalar &ds)
scalar rhoc_
Density [kg/m3].
vector U_
Velocity of Parcel [m/s].
scalar volume() const
Particle volume.
dimensioned< scalar > mag(const dimensioned< Type > &)
scalar We(const vector &U, const scalar d, const scalar rhoc, const scalar sigma) const
Weber number.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
label typeId() const
Return const access to type id.
scalar Re(const vector &U, const scalar d, const scalar rhoc, const scalar muc) const
Reynolds number.
scalar areaS() const
Particle surface area.
scalar d_
Diameter [m].
dictionary subOrEmptyDict(const word &, const bool mustRead=false) const
Find and return a sub-dictionary as a copy, or.
Definition: dictionary.C:727
const dictionary & dict() const
Return const access to the constant properties dictionary.