MomentumParcelI.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-2020 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  moving_(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  moving_(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 // * * * * * * * * * constantProperties Member Functions * * * * * * * * * * //
119 
120 template<class ParcelType>
121 inline const Foam::dictionary&
123 {
124  return dict_;
125 }
126 
127 
128 template<class ParcelType>
129 inline Foam::label
131 {
132  return parcelTypeId_.value();
133 }
134 
135 
136 template<class ParcelType>
137 inline Foam::scalar
139 {
140  return rhoMin_.value();
141 }
142 
143 
144 template<class ParcelType>
145 inline Foam::scalar
147 {
148  return rho0_.value();
149 }
150 
151 
152 template<class ParcelType>
153 inline Foam::scalar
155 {
156  return minParcelMass_.value();
157 }
158 
159 
160 // * * * * * * * MomentumParcel Member Functions * * * * * * * //
161 
162 template<class ParcelType>
164 {
165  return moving_;
166 }
167 
168 
169 template<class ParcelType>
171 {
172  return typeId_;
173 }
174 
175 
176 template<class ParcelType>
178 {
179  return nParticle_;
180 }
181 
182 
183 template<class ParcelType>
184 inline Foam::scalar Foam::MomentumParcel<ParcelType>::d() const
185 {
186  return d_;
187 }
188 
189 
190 template<class ParcelType>
191 inline Foam::scalar Foam::MomentumParcel<ParcelType>::dTarget() const
192 {
193  return dTarget_;
194 }
195 
196 
197 template<class ParcelType>
199 {
200  return U_;
201 }
202 
203 
204 template<class ParcelType>
205 inline Foam::scalar Foam::MomentumParcel<ParcelType>::rho() const
206 {
207  return rho_;
208 }
209 
210 
211 template<class ParcelType>
212 inline Foam::scalar Foam::MomentumParcel<ParcelType>::age() const
213 {
214  return age_;
215 }
216 
217 
218 template<class ParcelType>
219 inline Foam::scalar Foam::MomentumParcel<ParcelType>::tTurb() const
220 {
221  return tTurb_;
222 }
223 
224 
225 template<class ParcelType>
227 {
228  return UTurb_;
229 }
230 
231 
232 template<class ParcelType>
234 {
235  return moving_;
236 }
237 
238 
239 template<class ParcelType>
241 {
242  return typeId_;
243 }
244 
245 
246 template<class ParcelType>
248 {
249  return nParticle_;
250 }
251 
252 
253 template<class ParcelType>
255 {
256  return d_;
257 }
258 
259 
260 template<class ParcelType>
262 {
263  return dTarget_;
264 }
265 
266 
267 template<class ParcelType>
269 {
270  return U_;
271 }
272 
273 
274 template<class ParcelType>
276 {
277  return rho_;
278 }
279 
280 
281 template<class ParcelType>
283 {
284  return age_;
285 }
286 
287 
288 template<class ParcelType>
290 {
291  return tTurb_;
292 }
293 
294 
295 template<class ParcelType>
297 {
298  return UTurb_;
299 }
300 
301 
302 template<class ParcelType>
304 (
305  const trackingData& td
306 ) const
307 {
308  return td.rhoc()*this->mesh().cellVolumes()[this->cell()];
309 }
310 
311 
312 template<class ParcelType>
313 inline Foam::scalar Foam::MomentumParcel<ParcelType>::mass() const
314 {
315  return rho_*volume();
316 }
317 
318 
319 template<class ParcelType>
321 {
322  return 0.1*mass()*sqr(d_);
323 }
324 
325 
326 template<class ParcelType>
327 inline Foam::scalar Foam::MomentumParcel<ParcelType>::volume() const
328 {
329  return volume(d_);
330 }
331 
332 
333 template<class ParcelType>
334 inline Foam::scalar Foam::MomentumParcel<ParcelType>::volume(const scalar d)
335 {
336  return pi/6.0*pow3(d);
337 }
338 
339 
340 template<class ParcelType>
341 inline Foam::scalar Foam::MomentumParcel<ParcelType>::areaP() const
342 {
343  return areaP(d_);
344 }
345 
346 
347 template<class ParcelType>
348 inline Foam::scalar Foam::MomentumParcel<ParcelType>::areaP(const scalar d)
349 {
350  return 0.25*areaS(d);
351 }
352 
353 
354 template<class ParcelType>
355 inline Foam::scalar Foam::MomentumParcel<ParcelType>::areaS() const
356 {
357  return areaS(d_);
358 }
359 
360 
361 template<class ParcelType>
362 inline Foam::scalar Foam::MomentumParcel<ParcelType>::areaS(const scalar d)
363 {
364  return pi*d*d;
365 }
366 
367 
368 template<class ParcelType>
369 inline Foam::scalar Foam::MomentumParcel<ParcelType>::Re
370 (
371  const trackingData& td
372 ) const
373 {
374  return Re(td.rhoc(), U_, td.Uc(), d_, td.muc());
375 }
376 
377 
378 template<class ParcelType>
379 inline Foam::scalar Foam::MomentumParcel<ParcelType>::Re
380 (
381  const scalar rhoc,
382  const vector& U,
383  const vector& Uc,
384  const scalar d,
385  const scalar muc
386 )
387 {
388  return rhoc*mag(U - Uc)*d/max(muc, rootVSmall);
389 }
390 
391 
392 template<class ParcelType>
393 inline Foam::scalar Foam::MomentumParcel<ParcelType>::We
394 (
395  const trackingData& td,
396  const scalar sigma
397 ) const
398 {
399  return We(td.rhoc(), U_, td.Uc(), d_, sigma);
400 }
401 
402 
403 template<class ParcelType>
404 inline Foam::scalar Foam::MomentumParcel<ParcelType>::We
405 (
406  const scalar rhoc,
407  const vector& U,
408  const vector& Uc,
409  const scalar d,
410  const scalar sigma
411 )
412 {
413  return rhoc*magSqr(U - Uc)*d/max(sigma, rootVSmall);
414 }
415 
416 
417 template<class ParcelType>
418 inline Foam::scalar Foam::MomentumParcel<ParcelType>::Eo
419 (
420  const trackingData& td,
421  const scalar sigma
422 ) const
423 {
424  return Eo(td.g(), rho_, td.rhoc(), U_, d_, sigma);
425 }
426 
427 
428 template<class ParcelType>
429 inline Foam::scalar Foam::MomentumParcel<ParcelType>::Eo
430 (
431  const vector& g,
432  const scalar rho,
433  const scalar rhoc,
434  const vector& U,
435  const scalar d,
436  const scalar sigma
437 )
438 {
439  const vector dir = U/max(mag(U), rootVSmall);
440  return mag(g & dir)*mag(rho - rhoc)*sqr(d)/max(sigma, rootVSmall);
441 }
442 
443 
444 // ************************************************************************* //
scalar volume() const
Particle volume.
const vector & UTurb() const
Return const access to turbulent velocity fluctuation.
bool moving_
Moving flag - tracking stopped when moving = false.
scalar d_
Diameter [m].
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
scalar areaS() const
Particle surface area.
scalar dTarget_
Target diameter [m].
const dictionary dict_
Constant properties dictionary.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
dimensionedSymmTensor sqr(const dimensionedVector &dv)
scalar age() const
Return const access to the age.
scalar rho0() const
Return const access to the particle density.
label typeId_
Parcel type id.
scalar minParcelMass() const
Return const access to the minimum parcel mass.
scalar rhoMin() const
Return const access to the minimum density.
scalar We(const trackingData &td, const scalar sigma) const
Weber number.
fvMesh & mesh
scalar dTarget() const
Return const access to target diameter.
bool moving() const
Return const access to moving flag.
scalar d() const
Return const access to diameter.
const Type & value() const
Return the value.
scalar rho_
Density [kg/m^3].
scalar nParticle_
Number of particles in Parcel.
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
scalar massCell(const trackingData &td) const
Cell owner mass.
scalar nParticle() const
Return const access to number of particles.
mathematical constants.
scalar age_
Age [s].
scalar tTurb_
Time spent in turbulent eddy [s].
label parcelTypeId() const
Return const access to the parcel type id.
scalar Eo(const trackingData &td, const scalar sigma) const
Eotvos number.
static const zero Zero
Definition: zero.H:97
scalar momentOfInertia() const
Particle moment of inertia around diameter axis.
scalar muc() const
Return the continuous phase viscosity.
scalar Re(const trackingData &td) const
Reynolds number.
vector U_
Velocity of Parcel [m/s].
dimensioned< scalar > magSqr(const dimensioned< Type > &)
scalar rho() const
Return const access to density.
MomentumParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from mesh, coordinates and topology.
scalar mass() const
Particle mass.
scalar rhoc() const
Return the continuous phase density.
dimensionedScalar pow3(const dimensionedScalar &ds)
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:57
const vector & U() const
Return const access to velocity.
const dictionary & dict() const
Return const access to the constant properties dictionary.
dimensioned< scalar > mag(const dimensioned< Type > &)
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:105
label typeId() const
Return const access to type id.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:76
const vector & Uc() const
Return the continuous phase velocity.
scalar areaP() const
Particle projected area.
scalar tTurb() const
Return const access to time spent in turbulent eddy.
vector UTurb_
Turbulent velocity fluctuation [m/s].
dictionary subOrEmptyDict(const word &, const bool mustRead=false) const
Find and return a sub-dictionary as a copy, or.
Definition: dictionary.C:1033
Class to hold momentum parcel constant properties.