KinematicParcel.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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2016-2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27Class
28 Foam::KinematicParcel
29
30Group
31 grpLagrangianIntermediateParcels
32
33Description
34 Kinematic parcel class with rotational motion (as spherical
35 particles only) and one/two-way coupling with the continuous
36 phase.
37
38 Sub-models include:
39 - drag
40 - turbulent dispersion
41 - wall interactions
42
43SourceFiles
44 KinematicParcelI.H
45 KinematicParcel.C
46 KinematicParcelIO.C
47
48\*---------------------------------------------------------------------------*/
49
50#ifndef KinematicParcel_H
51#define KinematicParcel_H
52
53#include "particle.H"
54#include "IOstream.H"
55#include "autoPtr.H"
56#include "interpolation.H"
57#include "demandDrivenEntry.H"
58#include "labelFieldIOField.H"
59#include "vectorFieldIOField.H"
60
61// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
62
63namespace Foam
64{
65
66template<class ParcelType>
67class KinematicParcel;
68
69template<class Type>
70class AveragingMethod;
71
72// Forward declaration of friend functions
73
74template<class ParcelType>
75Ostream& operator<<
76(
77 Ostream&,
79);
80
81/*---------------------------------------------------------------------------*\
82 Class KinematicParcel Declaration
83\*---------------------------------------------------------------------------*/
84
85template<class ParcelType>
87:
88 public ParcelType
89{
90 // Private Data
91
92 //- Number of particle tracking attempts before we assume that it stalls
93 static label maxTrackAttempts;
94
95
96public:
97
98 //- Size in bytes of the fields
99 static const std::size_t sizeofFields;
100
101
102 //- Class to hold kinematic particle constant properties
104 {
105 protected:
106
107 // Protected Data
108
109 //- Constant properties dictionary
110 const dictionary dict_;
111
112
113 private:
114
115 // Private Data
116
117 //- Parcel type id - used for post-processing to flag the type
118 //- of parcels issued by this cloud
119 demandDrivenEntry<label> parcelTypeId_;
120
121 //- Minimum density [kg/m3]
123
124 //- Particle density [kg/m3] (constant)
126
127 //- Minimum parcel mass [kg]
128 demandDrivenEntry<scalar> minParcelMass_;
129
130
131 public:
132
133 // Constructors
134
135 //- Default construct
136 inline constantProperties();
137
138 //- Copy constructor
140
141 //- Construct from dictionary
142 inline constantProperties(const dictionary& parentDict);
143
144
145 // Member Functions
146
147 //- Return const access to the constant properties dictionary
148 inline const dictionary& dict() const;
149
150 //- Return const access to the parcel type id
151 inline label parcelTypeId() const;
152
153 //- Return const access to the minimum density
154 inline scalar rhoMin() const;
155
156 //- Return const access to the particle density
157 inline scalar rho0() const;
158
159 //- Return const access to the minimum parcel mass
160 inline scalar minParcelMass() const;
161 };
162
164 class trackingData
165 :
167 {
168 public:
170 enum trackPart
175 };
176
177
178 private:
179
180 // Private Data
181
182 // Interpolators for continuous phase fields
183
184 //- Density interpolator
186
187 //- Velocity interpolator
189
190 //- Dynamic viscosity interpolator
192
193
194 // Cached continuous phase properties
195
196 //- Density [kg/m3]
197 scalar rhoc_;
198
199 //- Velocity [m/s]
200 vector Uc_;
201
202 //- Viscosity [Pa.s]
203 scalar muc_;
204
205
206 // MPPIC Averages
207
208 //- Volume average
209 autoPtr<AveragingMethod<scalar>> volumeAverage_;
210
211 //- Radius average [ volume^(1/3) ]
212 autoPtr<AveragingMethod<scalar>> radiusAverage_;
213
214 //- Density average
216
217 //- Velocity average
219
220 //- Magnitude velocity squared average
222
223 //- Frequency average
224 autoPtr<AveragingMethod<scalar>> frequencyAverage_;
225
226 //- Mass average
228
229
230 //- Local gravitational or other body-force acceleration
231 const vector& g_;
232
233 //- label specifying which part of the integration
234 //- algorithm is taking place
235 trackPart part_;
236
237
238 public:
239
240 // Constructors
241
242 //- Construct from components
243 template <class TrackCloudType>
244 inline trackingData
245 (
246 const TrackCloudType& cloud,
248 );
249
250
251 // Member Functions
252
253 //- Return const access to the interpolator for continuous
254 //- phase density field
255 inline const interpolation<scalar>& rhoInterp() const;
256
257 //- Return const access to the interpolator for continuous
258 //- phase velocity field
259 inline const interpolation<vector>& UInterp() const;
260
261 //- Return const access to the interpolator for continuous
262 //- phase dynamic viscosity field
263 inline const interpolation<scalar>& muInterp() const;
264
265 //- Return the continuous phase density
266 inline scalar rhoc() const;
267
268 //- Access the continuous phase density
269 inline scalar& rhoc();
270
271 //- Return the continuous phase velocity
272 inline const vector& Uc() const;
273
274 //- Access the continuous phase velocity
275 inline vector& Uc();
276
277 //- Return the continuous phase viscosity
278 inline scalar muc() const;
279
280 //- Access the continuous phase viscosity
281 inline scalar& muc();
282
283 // Return const access to the gravitational acceleration vector
284 inline const vector& g() const;
285
286 //- Return the part of the tracking operation taking place
287 inline trackPart part() const;
288
289 //- Return access to the part of the tracking operation taking place
290 inline trackPart& part();
291
292 //- Update the MPPIC averages
293 template<class TrackCloudType>
294 inline void updateAverages(const TrackCloudType& cloud);
295
296 };
297
298
299protected:
300
301 // Protected Data
302
303 // Parcel properties
304
305 //- Active flag - tracking inactive when active = false
306 // Store as label for data alignment, but only has bool states
307 label active_;
308
309 //- Parcel type id
310 label typeId_;
311
312 //- Number of particles in Parcel
313 scalar nParticle_;
314
315 //- Diameter [m]
316 scalar d_;
317
318 //- Target diameter [m]
319 scalar dTarget_;
320
321 //- Velocity of Parcel [m/s]
322 vector U_;
323
324 //- Density [kg/m3]
325 scalar rho_;
326
327 //- Age [s]
328 scalar age_;
329
330 //- Time spent in turbulent eddy [s]
331 scalar tTurb_;
332
333 //- Turbulent velocity fluctuation [m/s]
335
336 //- Velocity correction due to collisions MPPIC [m/s]
338
339
340 // Protected Member Functions
341
342 //- Calculate new particle velocity
343 template<class TrackCloudType>
344 const vector calcVelocity
345 (
346 TrackCloudType& cloud,
347 trackingData& td,
348 const scalar dt, // timestep
349 const scalar Re, // Reynolds number
350 const scalar mu, // local carrier viscosity
351 const scalar mass, // mass
352 const vector& Su, // explicit particle momentum source
353 vector& dUTrans, // momentum transfer to carrier
354 scalar& Spu // linearised drag coefficient
355 ) const;
356
357
358public:
359
360 // Static data members
361
362 //- Runtime type information
363 TypeName("KinematicParcel");
364
365 //- String representation of properties
367 (
368 ParcelType,
369 " active"
370 + " typeId"
371 + " nParticle"
372 + " d"
373 + " dTarget"
374 + " (Ux Uy Uz)"
375 + " rho"
376 + " age"
377 + " tTurb"
378 + " (UTurbx UTurby UTurbz)"
379 + " (UCorrectx UCorrecty UCorrectz)"
380 );
381
382
383 // Constructors
384
385 //- Construct from mesh, coordinates and topology
386 // Other properties initialised as null
387 inline KinematicParcel
388 (
389 const polyMesh& mesh,
391 const label celli,
392 const label tetFacei,
393 const label tetPti
394 );
395
396 //- Construct from a position and a cell, searching for the rest of the
397 // required topology. Other properties are initialised as null.
398 inline KinematicParcel
399 (
400 const polyMesh& mesh,
401 const vector& position,
402 const label celli
403 );
404
405 //- Construct from components
406 inline KinematicParcel
407 (
408 const polyMesh& mesh,
410 const label celli,
411 const label tetFacei,
412 const label tetPti,
413 const label typeId,
414 const scalar nParticle0,
415 const scalar d0,
416 const scalar dTarget0,
417 const vector& U0,
418 const constantProperties& constProps
419 );
420
421 //- Construct from Istream
423 (
424 const polyMesh& mesh,
425 Istream& is,
426 bool readFields = true,
427 bool newFormat = true
428 );
429
430 //- Construct as a copy
432
433 //- Construct as a copy
435
436 //- Construct and return a (basic particle) clone
437 virtual autoPtr<particle> clone() const
438 {
439 return autoPtr<particle>(new KinematicParcel(*this));
440 }
441
442 //- Construct and return a (basic particle) clone
443 virtual autoPtr<particle> clone(const polyMesh& mesh) const
444 {
445 return autoPtr<particle>(new KinematicParcel(*this, mesh));
446 }
447
448 //- Factory class to read-construct particles used for
449 // parallel transfer
450 class iNew
451 {
452 const polyMesh& mesh_;
453
454 public:
456 iNew(const polyMesh& mesh)
457 :
458 mesh_(mesh)
459 {}
462 {
464 (
465 new KinematicParcel<ParcelType>(mesh_, is, true)
466 );
467 }
468 };
469
470
471 // Member Functions
472
473 // Access
474
475 //- Return const access to active flag
476 inline bool active() const;
477
478 //- Return const access to type id
479 inline label typeId() const;
480
481 //- Return const access to number of particles
482 inline scalar nParticle() const;
483
484 //- Return const access to diameter
485 inline scalar d() const;
486
487 //- Return const access to target diameter
488 inline scalar dTarget() const;
489
490 //- Return const access to velocity
491 inline const vector& U() const;
492
493 //- Return const access to density
494 inline scalar rho() const;
495
496 //- Return const access to the age
497 inline scalar age() const;
498
499 //- Return const access to time spent in turbulent eddy
500 inline scalar tTurb() const;
501
502 //- Return const access to turbulent velocity fluctuation
503 inline const vector& UTurb() const;
504
505 //- Return const access to correction velocity
506 inline const vector& UCorrect() const;
507
508
509 // Edit
510
511 //- Set active flag to the specified state
512 inline void active(const bool state);
513
514 //- Return access to type id
515 inline label& typeId();
516
517 //- Return access to number of particles
518 inline scalar& nParticle();
519
520 //- Return access to diameter
521 inline scalar& d();
522
523 //- Return access to target diameter
524 inline scalar& dTarget();
525
526 //- Return access to velocity
527 inline vector& U();
528
529 //- Return access to density
530 inline scalar& rho();
531
532 //- Return access to the age
533 inline scalar& age();
534
535 //- Return access to time spent in turbulent eddy
536 inline scalar& tTurb();
537
538 //- Return access to turbulent velocity fluctuation
539 inline vector& UTurb();
540
541 //- Return access to correction velocity
542 inline vector& UCorrect();
543
544
545 // Helper functions
546
547 //- Cell owner mass
548 inline scalar massCell(const trackingData& td) const;
549
550 //- Particle mass
551 inline scalar mass() const;
552
553 //- Particle moment of inertia around diameter axis
554 inline scalar momentOfInertia() const;
555
556 //- Particle volume
557 inline scalar volume() const;
558
559 //- Particle volume for a given diameter
560 inline static scalar volume(const scalar d);
561
562 //- Particle projected area
563 inline scalar areaP() const;
564
565 //- Projected area for given diameter
566 inline static scalar areaP(const scalar d);
567
568 //- Particle surface area
569 inline scalar areaS() const;
570
571 //- Surface area for given diameter
572 inline static scalar areaS(const scalar d);
573
574 //- Reynolds number
575 inline scalar Re(const trackingData& td) const;
576
577 //- Reynolds number for given conditions
578 inline static scalar Re
579 (
580 const scalar rhoc,
581 const vector& U,
582 const vector& Uc,
583 const scalar d,
584 const scalar muc
585 );
586
587 //- Weber number
588 inline scalar We
589 (
590 const trackingData& td,
591 const scalar sigma
592 ) const;
593
594 //- Weber number for given conditions
595 inline static scalar We
596 (
597 const scalar rhoc,
598 const vector& U,
599 const vector& Uc,
600 const scalar d,
601 const scalar sigma
602 );
603
604 //- Eotvos number
605 inline scalar Eo
606 (
607 const trackingData& td,
608 const scalar sigma
609 ) const;
610
611 //- Eotvos number for given conditions
612 inline static scalar Eo
613 (
614 const vector& g,
615 const scalar rho,
616 const scalar rhoc,
617 const vector& U,
618 const scalar d,
619 const scalar sigma
620 );
621
622
623 // Main calculation loop
624
625 //- Set cell values
626 template<class TrackCloudType>
627 void setCellValues(TrackCloudType& cloud, trackingData& td);
628
629 //- Apply dispersion to the carrier phase velocity and update
630 // parcel turbulence parameters
631 template<class TrackCloudType>
632 void calcDispersion
633 (
634 TrackCloudType& cloud,
635 trackingData& td,
636 const scalar dt
637 );
638
639 //- Correct cell values using latest transfer information
640 template<class TrackCloudType>
642 (
643 TrackCloudType& cloud,
644 trackingData& td,
645 const scalar dt
646 );
647
648 //- Update parcel properties over the time interval
649 template<class TrackCloudType>
650 void calc
651 (
652 TrackCloudType& cloud,
653 trackingData& td,
654 const scalar dt
655 );
656
657 //- Correct U following MP-PIC sub-models
658 template<class TrackCloudType>
659 void calcUCorrection
660 (
661 TrackCloudType& cloud,
662 trackingData& td,
663 const scalar dt
664 );
665
666
667 // Tracking
668
669 //- Move the parcel
670 template<class TrackCloudType>
671 bool move
672 (
673 TrackCloudType& cloud,
674 trackingData& td,
675 const scalar trackTime
676 );
677
678
679 // Patch interactions
680
681 //- Overridable function to handle the particle hitting a patch
682 // Executed before other patch-hitting functions
683 template<class TrackCloudType>
684 bool hitPatch(TrackCloudType& cloud, trackingData& td);
685
686 //- Overridable function to handle the particle hitting a
687 // processorPatch
688 template<class TrackCloudType>
689 void hitProcessorPatch(TrackCloudType& cloud, trackingData& td);
690
691 //- Overridable function to handle the particle hitting a wallPatch
692 template<class TrackCloudType>
693 void hitWallPatch(TrackCloudType& cloud, trackingData& td);
694
695 //- Transform the physical properties of the particle
696 // according to the given transformation tensor
697 virtual void transformProperties(const tensor& T);
698
699 //- Transform the physical properties of the particle
700 // according to the given separation vector
701 virtual void transformProperties(const vector& separation);
702
703
704 // I-O
705
706 //- Read
707 template<class TrackCloudType>
708 static void readFields(TrackCloudType& c);
709
710 //- Write
711 template<class TrackCloudType>
712 static void writeFields(const TrackCloudType& c);
713
714 //- Write individual parcel properties to stream
715 void writeProperties
716 (
717 Ostream& os,
718 const wordRes& filters,
719 const word& delim,
720 const bool namesOnly = false
721 ) const;
722
723 //- Read particle fields as objects from the obr registry
724 template<class CloudType>
725 static void readObjects(CloudType& c, const objectRegistry& obr);
726
727 //- Write particle fields as objects into the obr registry
728 template<class CloudType>
729 static void writeObjects(const CloudType& c, objectRegistry& obr);
730
731
732 // Ostream Operator
734 friend Ostream& operator<< <ParcelType>
735 (
736 Ostream&,
738 );
739};
740
741
742// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
743
744} // End namespace Foam
745
746// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
747
748#include "KinematicParcelI.H"
750
751// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
752
753#ifdef NoRepository
754 #include "KinematicParcel.C"
755#endif
756
757// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
758
759#endif
760
761// ************************************************************************* //
const uniformDimensionedVectorField & g
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:75
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:64
Class to hold kinematic particle 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_
Constant properties dictionary.
const dictionary & dict() const
Return const access to the constant properties dictionary.
scalar rhoMin() const
Return const access to the minimum density.
Factory class to read-construct particles used for.
iNew(const polyMesh &mesh)
autoPtr< KinematicParcel< ParcelType > > operator()(Istream &is) const
const vector & Uc() const
Return the continuous phase velocity.
scalar muc() const
Return the continuous phase viscosity.
const interpolation< vector > & UInterp() const
scalar rhoc() const
Return the continuous phase density.
const interpolation< scalar > & muInterp() const
const interpolation< scalar > & rhoInterp() const
trackPart part() const
Return the part of the tracking operation taking place.
void updateAverages(const TrackCloudType &cloud)
Update the MPPIC averages.
Kinematic parcel class with rotational motion (as spherical particles only) and one/two-way coupling ...
const vector & UCorrect() const
Return const access to correction velocity.
virtual void transformProperties(const tensor &T)
Transform the physical properties of the particle.
void calcUCorrection(TrackCloudType &cloud, trackingData &td, const scalar dt)
Correct U following MP-PIC sub-models.
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.
void calcDispersion(TrackCloudType &cloud, trackingData &td, const scalar dt)
Apply dispersion to the carrier phase velocity and update.
scalar tTurb() const
Return const access to time spent in turbulent eddy.
virtual autoPtr< particle > clone(const polyMesh &mesh) const
Construct and return a (basic particle) clone.
scalar rho_
Density [kg/m3].
const vector & U() const
Return const access to velocity.
static void writeObjects(const CloudType &c, objectRegistry &obr)
Write particle fields as objects into the obr registry.
label active_
Active flag - tracking inactive when active = false.
KinematicParcel(const KinematicParcel &p)
Construct as a copy.
static void readFields(TrackCloudType &c)
Read.
scalar d() const
Return const access to diameter.
void cellValueSourceCorrection(TrackCloudType &cloud, trackingData &td, const scalar dt)
Correct cell values using latest transfer information.
scalar Re(const trackingData &td) const
Reynolds number.
scalar tTurb_
Time spent in turbulent eddy [s].
virtual autoPtr< particle > clone() const
Construct and return a (basic particle) clone.
scalar volume() const
Particle volume.
scalar dTarget_
Target diameter [m].
const vector calcVelocity(TrackCloudType &cloud, trackingData &td, const scalar dt, const scalar Re, const scalar mu, const scalar mass, const vector &Su, vector &dUTrans, scalar &Spu) const
Calculate new particle velocity.
static const std::size_t sizeofFields
Size in bytes of the fields.
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.
static void readObjects(CloudType &c, const objectRegistry &obr)
Read particle fields as objects from the obr registry.
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.
AddToPropertyList(ParcelType, " active"+" typeId"+" nParticle"+" d"+" dTarget"+" (Ux Uy Uz)"+" rho"+" age"+" tTurb"+" (UTurbx UTurby UTurbz)"+" (UCorrectx UCorrecty UCorrectz)")
String representation of properties.
static void writeFields(const TrackCloudType &c)
Write.
void hitProcessorPatch(TrackCloudType &cloud, trackingData &td)
Overridable function to handle the particle hitting a.
bool hitPatch(TrackCloudType &cloud, trackingData &td)
Overridable function to handle the particle hitting a patch.
vector UCorrect_
Velocity correction due to collisions MPPIC [m/s].
void setCellValues(TrackCloudType &cloud, trackingData &td)
Set cell values.
scalar areaS() const
Particle surface area.
bool move(TrackCloudType &cloud, trackingData &td, const scalar trackTime)
Move the parcel.
void hitWallPatch(TrackCloudType &cloud, trackingData &td)
Overridable function to handle the particle hitting a wallPatch.
scalar Eo(const trackingData &td, const scalar sigma) const
Eotvos number.
void writeProperties(Ostream &os, const wordRes &filters, const word &delim, const bool namesOnly=false) const
Write individual parcel properties to stream.
TypeName("KinematicParcel")
Runtime type information.
bool active() const
Return const access to active flag.
KinematicParcel(const KinematicParcel &p, const polyMesh &mesh)
Construct as a copy.
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.
void calc(TrackCloudType &cloud, trackingData &td, const scalar dt)
Update parcel properties over the time interval.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A cloud is a registry collection of lagrangian particles.
Definition: cloud.H:60
Class for demand-driven dictionary entries.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Class used to pass data into container.
Abstract base class for volume field interpolation.
Definition: interpolation.H:60
Registry of regIOobjects.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:54
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
const volScalarField & mu
dynamicFvMesh & mesh
PtrList< coordinateSystem > coordinates(solidRegions.size())
OBJstream os(runTime.globalPath()/outputName)
zeroField Su
Definition: alphaSuSp.H:1
Namespace for OpenFOAM.
#define AddToPropertyList(ParcelType, str)
Add to existing static 'propertyList' for particle properties.
const volScalarField & cp
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73