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-2023 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  label& nLocateBoundaryHits
105 )
106 :
107  ParcelType(owner, position, celli, nLocateBoundaryHits),
108  moving_(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 {}
119 
120 
121 // * * * * * * * * * constantProperties Member Functions * * * * * * * * * * //
122 
123 template<class ParcelType>
124 inline const Foam::dictionary&
126 {
127  return dict_;
128 }
129 
130 
131 template<class ParcelType>
132 inline Foam::label
134 {
135  return parcelTypeId_.value();
136 }
137 
138 
139 template<class ParcelType>
140 inline Foam::scalar
142 {
143  return rhoMin_.value();
144 }
145 
146 
147 template<class ParcelType>
148 inline Foam::scalar
150 {
151  return rho0_.value();
152 }
153 
154 
155 template<class ParcelType>
156 inline Foam::scalar
158 {
159  return minParcelMass_.value();
160 }
161 
162 
163 // * * * * * * * MomentumParcel Member Functions * * * * * * * //
164 
165 template<class ParcelType>
167 {
168  return moving_;
169 }
170 
171 
172 template<class ParcelType>
174 {
175  return typeId_;
176 }
177 
178 
179 template<class ParcelType>
181 {
182  return nParticle_;
183 }
184 
185 
186 template<class ParcelType>
187 inline Foam::scalar Foam::MomentumParcel<ParcelType>::d() const
188 {
189  return d_;
190 }
191 
192 
193 template<class ParcelType>
194 inline Foam::scalar Foam::MomentumParcel<ParcelType>::dTarget() const
195 {
196  return dTarget_;
197 }
198 
199 
200 template<class ParcelType>
202 {
203  return U_;
204 }
205 
206 
207 template<class ParcelType>
208 inline Foam::scalar Foam::MomentumParcel<ParcelType>::rho() const
209 {
210  return rho_;
211 }
212 
213 
214 template<class ParcelType>
215 inline Foam::scalar Foam::MomentumParcel<ParcelType>::age() const
216 {
217  return age_;
218 }
219 
220 
221 template<class ParcelType>
222 inline Foam::scalar Foam::MomentumParcel<ParcelType>::tTurb() const
223 {
224  return tTurb_;
225 }
226 
227 
228 template<class ParcelType>
230 {
231  return UTurb_;
232 }
233 
234 
235 template<class ParcelType>
237 {
238  return moving_;
239 }
240 
241 
242 template<class ParcelType>
244 {
245  return typeId_;
246 }
247 
248 
249 template<class ParcelType>
251 {
252  return nParticle_;
253 }
254 
255 
256 template<class ParcelType>
258 {
259  return d_;
260 }
261 
262 
263 template<class ParcelType>
265 {
266  return dTarget_;
267 }
268 
269 
270 template<class ParcelType>
272 {
273  return U_;
274 }
275 
276 
277 template<class ParcelType>
279 {
280  return rho_;
281 }
282 
283 
284 template<class ParcelType>
286 {
287  return age_;
288 }
289 
290 
291 template<class ParcelType>
293 {
294  return tTurb_;
295 }
296 
297 
298 template<class ParcelType>
300 {
301  return UTurb_;
302 }
303 
304 
305 template<class ParcelType>
307 (
308  const trackingData& td
309 ) const
310 {
311  return td.rhoc()*td.mesh.cellVolumes()[this->cell()];
312 }
313 
314 
315 template<class ParcelType>
316 inline Foam::scalar Foam::MomentumParcel<ParcelType>::mass() const
317 {
318  return rho_*volume();
319 }
320 
321 
322 template<class ParcelType>
324 {
325  return 0.1*mass()*sqr(d_);
326 }
327 
328 
329 template<class ParcelType>
330 inline Foam::scalar Foam::MomentumParcel<ParcelType>::volume() const
331 {
332  return volume(d_);
333 }
334 
335 
336 template<class ParcelType>
337 inline Foam::scalar Foam::MomentumParcel<ParcelType>::volume(const scalar d)
338 {
339  return pi/6.0*pow3(d);
340 }
341 
342 
343 template<class ParcelType>
344 inline Foam::scalar Foam::MomentumParcel<ParcelType>::areaP() const
345 {
346  return areaP(d_);
347 }
348 
349 
350 template<class ParcelType>
351 inline Foam::scalar Foam::MomentumParcel<ParcelType>::areaP(const scalar d)
352 {
353  return 0.25*areaS(d);
354 }
355 
356 
357 template<class ParcelType>
358 inline Foam::scalar Foam::MomentumParcel<ParcelType>::areaS() const
359 {
360  return areaS(d_);
361 }
362 
363 
364 template<class ParcelType>
365 inline Foam::scalar Foam::MomentumParcel<ParcelType>::areaS(const scalar d)
366 {
367  return pi*d*d;
368 }
369 
370 
371 template<class ParcelType>
373 (
374  const trackingData& td
375 ) const
376 {
377  return Re(td.rhoc(), U_, td.Uc(), d_, td.muc());
378 }
379 
380 
381 template<class ParcelType>
383 (
384  const scalar rhoc,
385  const vector& U,
386  const vector& Uc,
387  const scalar d,
388  const scalar muc
389 )
390 {
391  return rhoc*mag(U - Uc)*d/max(muc, rootVSmall);
392 }
393 
394 
395 template<class ParcelType>
397 (
398  const trackingData& td,
399  const scalar sigma
400 ) const
401 {
402  return We(td.rhoc(), U_, td.Uc(), d_, sigma);
403 }
404 
405 
406 template<class ParcelType>
408 (
409  const scalar rhoc,
410  const vector& U,
411  const vector& Uc,
412  const scalar d,
413  const scalar sigma
414 )
415 {
416  return rhoc*magSqr(U - Uc)*d/max(sigma, rootVSmall);
417 }
418 
419 
420 template<class ParcelType>
422 (
423  const trackingData& td,
424  const scalar sigma
425 ) const
426 {
427  return Eo(td.g(), rho_, td.rhoc(), U_, d_, sigma);
428 }
429 
430 
431 template<class ParcelType>
433 (
434  const vector& g,
435  const scalar rho,
436  const scalar rhoc,
437  const vector& U,
438  const scalar d,
439  const scalar sigma
440 )
441 {
442  const vector dir = U/max(mag(U), rootVSmall);
443  return mag(g & dir)*mag(rho - rhoc)*sqr(d)/max(sigma, rootVSmall);
444 }
445 
446 
447 // ************************************************************************* //
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:162
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:106
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