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-2016 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 vector& position,
77  const label celli,
78  const label tetFacei,
79  const label tetPtI
80 )
81 :
82  ParcelType(owner, position, 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  const label tetFacei,
106  const label tetPtI,
107  const label typeId,
108  const scalar nParticle0,
109  const scalar d0,
110  const scalar dTarget0,
111  const vector& U0,
112  const constantProperties& constProps
113 )
114 :
115  ParcelType(owner, position, celli, tetFacei, tetPtI),
116  active_(true),
117  typeId_(typeId),
118  nParticle_(nParticle0),
119  d_(d0),
120  dTarget_(dTarget0),
121  U_(U0),
122  rho_(constProps.rho0()),
123  age_(0.0),
124  tTurb_(0.0),
125  UTurb_(Zero),
126  rhoc_(0.0),
127  Uc_(Zero),
128  muc_(0.0)
129 {}
130 
131 
132 // * * * * * * * * * constantProperties Member Functions * * * * * * * * * * //
133 
134 template<class ParcelType>
135 inline const Foam::dictionary&
137 {
138  return dict_;
139 }
140 
141 
142 template<class ParcelType>
143 inline Foam::label
145 {
146  return parcelTypeId_.value();
147 }
148 
149 
150 template<class ParcelType>
151 inline Foam::scalar
153 {
154  return rhoMin_.value();
155 }
156 
157 
158 template<class ParcelType>
159 inline Foam::scalar
161 {
162  return rho0_.value();
163 }
164 
165 
166 template<class ParcelType>
167 inline Foam::scalar
169 {
170  return minParcelMass_.value();
171 }
172 
173 
174 // * * * * * * * KinematicParcel Member Functions * * * * * * * //
175 
176 template<class ParcelType>
178 {
179  return active_;
180 }
181 
182 
183 template<class ParcelType>
185 {
186  return typeId_;
187 }
188 
189 
190 template<class ParcelType>
192 {
193  return nParticle_;
194 }
195 
196 
197 template<class ParcelType>
198 inline Foam::scalar Foam::KinematicParcel<ParcelType>::d() const
199 {
200  return d_;
201 }
202 
203 
204 template<class ParcelType>
206 {
207  return dTarget_;
208 }
209 
210 
211 template<class ParcelType>
213 {
214  return U_;
215 }
216 
217 
218 template<class ParcelType>
219 inline Foam::scalar Foam::KinematicParcel<ParcelType>::rho() const
220 {
221  return rho_;
222 }
223 
224 
225 template<class ParcelType>
226 inline Foam::scalar Foam::KinematicParcel<ParcelType>::age() const
227 {
228  return age_;
229 }
230 
231 
232 template<class ParcelType>
233 inline Foam::scalar Foam::KinematicParcel<ParcelType>::tTurb() const
234 {
235  return tTurb_;
236 }
237 
238 
239 template<class ParcelType>
241 {
242  return UTurb_;
243 }
244 
245 
246 template<class ParcelType>
247 inline Foam::scalar Foam::KinematicParcel<ParcelType>::rhoc() const
248 {
249  return rhoc_;
250 }
251 
252 
253 template<class ParcelType>
255 {
256  return Uc_;
257 }
258 
259 
260 template<class ParcelType>
261 inline Foam::scalar Foam::KinematicParcel<ParcelType>::muc() const
262 {
263  return muc_;
264 }
265 
266 
267 template<class ParcelType>
269 {
270  return active_;
271 }
272 
273 
274 template<class ParcelType>
276 {
277  return typeId_;
278 }
279 
280 
281 template<class ParcelType>
283 {
284  return nParticle_;
285 }
286 
287 
288 template<class ParcelType>
290 {
291  return d_;
292 }
293 
294 
295 template<class ParcelType>
297 {
298  return dTarget_;
299 }
300 
301 
302 template<class ParcelType>
304 {
305  return U_;
306 }
307 
308 
309 template<class ParcelType>
311 {
312  return rho_;
313 }
314 
315 
316 template<class ParcelType>
318 {
319  return age_;
320 }
321 
322 
323 template<class ParcelType>
325 {
326  return tTurb_;
327 }
328 
329 
330 template<class ParcelType>
332 {
333  return UTurb_;
334 }
335 
336 
337 template<class ParcelType>
339 {
340  // Use volume-based interpolation if dealing with external faces
341  if (this->cloud().internalFace(this->face()))
342  {
343  return this->face();
344  }
345  else
346  {
347  return -1;
348  }
349 }
350 
351 
352 template<class ParcelType>
354 (
355  const label celli
356 ) const
357 {
358  return rhoc_*this->mesh().cellVolumes()[celli];
359 }
360 
361 
362 template<class ParcelType>
363 inline Foam::scalar Foam::KinematicParcel<ParcelType>::mass() const
364 {
365  return rho_*volume();
366 }
367 
368 
369 template<class ParcelType>
371 {
372  return 0.1*mass()*sqr(d_);
373 }
374 
375 
376 template<class ParcelType>
377 inline Foam::scalar Foam::KinematicParcel<ParcelType>::volume() const
378 {
379  return volume(d_);
380 }
381 
382 
383 template<class ParcelType>
384 inline Foam::scalar Foam::KinematicParcel<ParcelType>::volume(const scalar d)
385 {
386  return pi/6.0*pow3(d);
387 }
388 
389 
390 template<class ParcelType>
391 inline Foam::scalar Foam::KinematicParcel<ParcelType>::areaP() const
392 {
393  return areaP(d_);
394 }
395 
396 
397 template<class ParcelType>
398 inline Foam::scalar Foam::KinematicParcel<ParcelType>::areaP(const scalar d)
399 {
400  return 0.25*areaS(d);
401 }
402 
403 
404 template<class ParcelType>
405 inline Foam::scalar Foam::KinematicParcel<ParcelType>::areaS() const
406 {
407  return areaS(d_);
408 }
409 
410 
411 template<class ParcelType>
412 inline Foam::scalar Foam::KinematicParcel<ParcelType>::areaS(const scalar d)
413 {
414  return pi*d*d;
415 }
416 
417 
418 template<class ParcelType>
419 inline Foam::scalar Foam::KinematicParcel<ParcelType>::Re
420 (
421  const vector& U,
422  const scalar d,
423  const scalar rhoc,
424  const scalar muc
425 ) const
426 {
427  return rhoc*mag(U - Uc_)*d/(muc + ROOTVSMALL);
428 }
429 
430 
431 template<class ParcelType>
432 inline Foam::scalar Foam::KinematicParcel<ParcelType>::We
433 (
434  const vector& U,
435  const scalar d,
436  const scalar rhoc,
437  const scalar sigma
438 ) const
439 {
440  return rhoc*magSqr(U - Uc_)*d/(sigma + ROOTVSMALL);
441 }
442 
443 
444 template<class ParcelType>
445 inline Foam::scalar Foam::KinematicParcel<ParcelType>::Eo
446 (
447  const vector& a,
448  const scalar d,
449  const scalar sigma
450 ) const
451 {
452  vector dir = U_/(mag(U_) + ROOTVSMALL);
453  return mag(a & dir)*(rho_ - rhoc_)*sqr(d)/(sigma + ROOTVSMALL);
454 }
455 
456 
457 // ************************************************************************* //
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 mass() const
Particle mass.
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
const double e
Elementary charge.
Definition: doubleFloat.H:78
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
scalar Eo(const vector &a, const scalar d, const scalar sigma) const
Eotvos number.
scalar Re(const vector &U, const scalar d, const scalar rhoc, const scalar muc) const
Reynolds number.
const dictionary & dict() const
Return const access to the constant properties dictionary.
dictionary subOrEmptyDict(const word &, const bool mustRead=false) const
Find and return a sub-dictionary as a copy, or.
Definition: dictionary.C:668
dimensionedSymmTensor sqr(const dimensionedVector &dv)
scalar rho0() const
Return const access to the particle density.
Class to hold kinematic particle constant properties.
scalar rho() const
Return const access to density.
scalar tTurb_
Time spent in turbulent eddy [s].
scalar age_
Age [s].
const Type & value() const
Return the value.
KinematicParcel(const polyMesh &mesh, const vector &position, const label celli, const label tetFacei, const label tetPtI)
Construct from owner, position, and cloud owner.
label faceInterpolation() const
Return the index of the face used in the interpolation routine.
const vector & Uc() const
Return const access to carrier velocity [m/s].
scalar volume() const
Particle volume.
scalar massCell(const label celli) const
Cell owner mass.
scalar dTarget() const
Return const access to target diameter.
const vector & UTurb() const
Return const access to turbulent velocity fluctuation.
scalar muc() const
Return const access to carrier viscosity [Pa.s].
vector Uc_
Velocity [m/s].
scalar areaS() const
Particle surface area.
scalar dTarget_
Target diameter [m].
scalar rhoc() const
Return const access to carrier density [kg/m3].
scalar tTurb() const
Return const access to time spent in turbulent eddy.
scalar We(const vector &U, const scalar d, const scalar rhoc, const scalar sigma) const
Weber number.
scalar momentOfInertia() const
Particle moment of inertia around diameter axis.
dynamicFvMesh & mesh
scalar rho_
Density [kg/m3].
mathematical constants.
scalar muc_
Viscosity [Pa.s].
bool active_
Active flag - tracking inactive when active = false.
const dictionary dict_
Constant properties dictionary.
A cloud is a collection of lagrangian particles.
Definition: cloud.H:51
label typeId_
Parcel type id.
vector UTurb_
Turbulent velocity fluctuation [m/s].
label parcelTypeId() const
Return const access to the parcel type id.
static const zero Zero
Definition: zero.H:91
scalar age() const
Return const access to the age.
dimensioned< scalar > magSqr(const dimensioned< Type > &)
scalar nParticle_
Number of particles in Parcel.
dimensionedScalar pow3(const dimensionedScalar &ds)
scalar rhoc_
Density [kg/m3].
vector U_
Velocity of Parcel [m/s].
scalar rhoMin() const
Return const access to the minimum density.
bool active() const
Return const access to active flag.
dimensioned< scalar > mag(const dimensioned< Type > &)
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
scalar areaP() const
Particle projected area.
scalar d() const
Return const access to diameter.
scalar d_
Diameter [m].
scalar nParticle() const
Return const access to number of particles.
scalar minParcelMass() const
Return const access to the minimum parcel mass.
const vector & U() const
Return const access to velocity.
label typeId() const
Return const access to type id.