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