BUG: UListIO: byteSize overflowing on really big faceLists
[OpenFOAM-2.0.x.git] / src / lagrangian / intermediate / parcels / Templates / KinematicParcel / KinematicParcel.H
blobc4ef66572fca6fedca1c73d37c7347a3c500d552
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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.
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.
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/>.
24 Class
25     Foam::KinematicParcel
27 Description
28     Kinematic parcel class with rotational motion (as spherical
29     particles only) and one/two-way coupling with the continuous
30     phase.
32     Sub-models include:
33     - drag
34     - turbulent dispersion
35     - wall interactions
37 SourceFiles
38     KinematicParcelI.H
39     KinematicParcel.C
40     KinematicParcelIO.C
42 \*---------------------------------------------------------------------------*/
44 #ifndef KinematicParcel_H
45 #define KinematicParcel_H
47 #include "particle.H"
48 #include "IOstream.H"
49 #include "autoPtr.H"
50 #include "interpolation.H"
52 // #include "ParticleForceList.H" // TODO
54 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
56 namespace Foam
59 template<class ParcelType>
60 class KinematicParcel;
62 // Forward declaration of friend functions
64 template<class ParcelType>
65 Ostream& operator<<
67     Ostream&,
68     const KinematicParcel<ParcelType>&
71 /*---------------------------------------------------------------------------*\
72                          Class KinematicParcel Declaration
73 \*---------------------------------------------------------------------------*/
75 template<class ParcelType>
76 class KinematicParcel
78     public ParcelType
80 public:
82     //- Class to hold kinematic particle constant properties
83     class constantProperties
84     {
85         // Private data
87             //- Constant properties dictionary
88             const dictionary dict_;
90             //- Parcel type id - used for post-processing to flag the type
91             //  of parcels issued by this cloud
92             label parcelTypeId_;
94             //- Minimum density [kg/m3]
95             scalar rhoMin_;
97             //- Particle density [kg/m3] (constant)
98             scalar rho0_;
100             //- Minimum particle mass [kg]
101             scalar minParticleMass_;
103             //- Young's modulus [N/m2]
104             scalar youngsModulus_;
106             //- Poisson's ratio
107             scalar poissonsRatio_;
110     public:
112         // Constructors
114             //- Null constructor
115             constantProperties();
117             //- Copy constructor
118             constantProperties(const constantProperties& cp);
120             //- Constructor from dictionary
121             constantProperties
122             (
123                 const dictionary& parentDict,
124                 const bool readFields = true
125             );
127             //- Construct from components
128             constantProperties
129             (
130                 const label parcelTypeId,
131                 const scalar rhoMin,
132                 const scalar rho0,
133                 const scalar minParticleMass,
134                 const scalar youngsModulus,
135                 const scalar poissonsRatio
136             );
139         // Member functions
141             //- Return const access to the constant properties dictionary
142             inline const dictionary& dict() const;
144             //- Return const access to the parcel type id
145             inline label parcelTypeId() const;
147             //- Return const access to the minimum density
148             inline scalar rhoMin() const;
150             //- Return const access to the particle density
151             inline scalar rho0() const;
153             //- Return const access to the minimum particle mass
154             inline scalar minParticleMass() const;
156             //- Return const access to Young's Modulus
157             inline scalar youngsModulus() const;
159             //- Return const access to Poisson's ratio
160             inline scalar poissonsRatio() const;
161     };
164     template<class CloudType>
165     class TrackingData
166     :
167         public ParcelType::template TrackingData<CloudType>
168     {
169     public:
171         enum trackPart
172         {
173             tpVelocityHalfStep,
174             tpLinearTrack,
175             tpRotationalTrack
176         };
179     private:
181         // Private data
183             // Interpolators for continuous phase fields
185                 //- Density interpolator
186                 autoPtr<interpolation<scalar> > rhoInterp_;
188                 //- Velocity interpolator
189                 autoPtr<interpolation<vector> > UInterp_;
191                 //- Dynamic viscosity interpolator
192                 autoPtr<interpolation<scalar> > muInterp_;
195             //- Local gravitational or other body-force acceleration
196             const vector& g_;
198             // label specifying which part of the integration
199             // algorithm is taking place
200             trackPart part_;
203     public:
205         // Constructors
207             //- Construct from components
208             inline TrackingData
209             (
210                 CloudType& cloud,
211                 trackPart part = tpLinearTrack
212             );
215         // Member functions
217             //- Return conat access to the interpolator for continuous
218             //  phase density field
219             inline const interpolation<scalar>& rhoInterp() const;
221             //- Return conat access to the interpolator for continuous
222             //  phase velocity field
223             inline const interpolation<vector>& UInterp() const;
225             //- Return conat access to the interpolator for continuous
226             //  phase dynamic viscosity field
227             inline const interpolation<scalar>& muInterp() const;
229             // Return const access to the gravitational acceleration vector
230             inline const vector& g() const;
232             //- Return the part of the tracking operation taking place
233             inline trackPart part() const;
235             //- Return access to the part of the tracking operation taking place
236             inline trackPart& part();
237     };
240 protected:
242     // Protected data
244         // Parcel properties
246             //- Active flag - tracking inactive when active = false
247             bool active_;
249             //- Parcel type id
250             label typeId_;
252             //- Number of particles in Parcel
253             scalar nParticle_;
255             //- Diameter [m]
256             scalar d_;
258             //- Target diameter [m]
259             scalar dTarget_;
261             //- Velocity of Parcel [m/s]
262             vector U_;
264             //- Force on particle due to collisions [N]
265             vector f_;
267             //- Angular momentum of Parcel in global reference frame
268             // [kg m2/s]
269             vector angularMomentum_;
271             //- Torque on particle due to collisions in global
272             //  reference frame [Nm]
273             vector torque_;
275             //- Density [kg/m3]
276             scalar rho_;
278             //- Age [s]
279             scalar age_;
281             //- Time spent in turbulent eddy [s]
282             scalar tTurb_;
284             //- Turbulent velocity fluctuation [m/s]
285             vector UTurb_;
288         // Cell-based quantities
290             //- Density [kg/m3]
291             scalar rhoc_;
293             //- Velocity [m/s]
294             vector Uc_;
296             //- Viscosity [Pa.s]
297             scalar muc_;
300     // Protected Member Functions
302         //- Calculate new particle velocity
303         template<class TrackData>
304         const vector calcVelocity
305         (
306             TrackData& td,
307             const scalar dt,           // timestep
308             const label cellI,         // owner cell
309             const scalar Re,           // Reynolds number
310             const scalar mu,           // local carrier viscosity
311             const scalar mass,         // mass
312             const vector& Su,          // explicit particle momentum source
313             vector& dUTrans,           // momentum transfer to carrier
314             scalar& Spu                // linearised drag coefficient
315         ) const;
318 public:
320     // Static data members
322         //- String representation of properties
323         static string propHeader;
325         //- Runtime type information
326         TypeName("KinematicParcel");
329     // Constructors
331         //- Construct from owner, position, and cloud owner
332         //  Other properties initialised as null
333         inline KinematicParcel
334         (
335             const polyMesh& mesh,
336             const vector& position,
337             const label cellI,
338             const label tetFaceI,
339             const label tetPtI
340         );
342         //- Construct from components
343         inline KinematicParcel
344         (
345             const polyMesh& mesh,
346             const vector& position,
347             const label cellI,
348             const label tetFaceI,
349             const label tetPtI,
350             const label typeId,
351             const scalar nParticle0,
352             const scalar d0,
353             const scalar dTarget0,
354             const vector& U0,
355             const vector& f0,
356             const vector& angularMomentum0,
357             const vector& torque0,
358             const constantProperties& constProps
359         );
361         //- Construct from Istream
362         KinematicParcel
363         (
364             const polyMesh& mesh,
365             Istream& is,
366             bool readFields = true
367         );
369         //- Construct as a copy
370         KinematicParcel(const KinematicParcel& p);
372         //- Construct as a copy
373         KinematicParcel(const KinematicParcel& p, const polyMesh& mesh);
375         //- Construct and return a (basic particle) clone
376         virtual autoPtr<particle> clone() const
377         {
378             return autoPtr<particle>(new KinematicParcel(*this));
379         }
381         //- Construct and return a (basic particle) clone
382         virtual autoPtr<particle> clone(const polyMesh& mesh) const
383         {
384             return autoPtr<particle>(new KinematicParcel(*this, mesh));
385         }
387         //- Factory class to read-construct particles used for
388         //  parallel transfer
389         class iNew
390         {
391             const polyMesh& mesh_;
393         public:
395             iNew(const polyMesh& mesh)
396             :
397                 mesh_(mesh)
398             {}
400             autoPtr<KinematicParcel<ParcelType> > operator()(Istream& is) const
401             {
402                 return autoPtr<KinematicParcel<ParcelType> >
403                 (
404                     new KinematicParcel<ParcelType>(mesh_, is, true)
405                 );
406             }
407         };
410     // Member Functions
412         // Access
414             //- Return const access to active flag
415             inline bool active() const;
417             //- Return const access to type id
418             inline label typeId() const;
420             //- Return const access to number of particles
421             inline scalar nParticle() const;
423             //- Return const access to diameter
424             inline scalar d() const;
426             //- Return const access to target diameter
427             inline scalar dTarget() const;
429             //- Return const access to velocity
430             inline const vector& U() const;
432             //- Return const access to force
433             inline const vector& f() const;
435             //- Return const access to angular momentum
436             inline const vector& angularMomentum() const;
438             //- Return const access to torque
439             inline const vector& torque() const;
441             //- Return const access to density
442             inline scalar rho() const;
444             //- Return const access to the age
445             inline scalar age() const;
447             //- Return const access to time spent in turbulent eddy
448             inline scalar tTurb() const;
450             //- Return const access to turbulent velocity fluctuation
451             inline const vector& UTurb() const;
453             //- Return const access to carrier density [kg/m3]
454             inline scalar rhoc() const;
456             //- Return const access to carrier velocity [m/s]
457             inline const vector& Uc() const;
459             //- Return const access to carrier viscosity [Pa.s]
460             inline scalar muc() const;
463         // Edit
465             //- Return const access to active flag
466             inline bool& active();
468             //- Return access to type id
469             inline label& typeId();
471             //- Return access to number of particles
472             inline scalar& nParticle();
474             //- Return access to diameter
475             inline scalar& d();
477             //- Return access to target diameter
478             inline scalar& dTarget();
480             //- Return access to velocity
481             inline vector& U();
483             //- Return access to force
484             inline vector& f();
486             //- Return access to angular momentum
487             inline vector& angularMomentum();
489             //- Return access to torque
490             inline vector& torque();
492             //- Return access to density
493             inline scalar& rho();
495             //- Return access to the age
496             inline scalar& age();
498             //- Return access to time spent in turbulent eddy
499             inline scalar& tTurb();
501             //- Return access to turbulent velocity fluctuation
502             inline vector& UTurb();
505         // Helper functions
507             //- Return the index of the face to be used in the interpolation
508             //  routine
509             inline label faceInterpolation() const;
511             //- Cell owner mass
512             inline scalar massCell(const label cellI) const;
514             //- Particle mass
515             inline scalar mass() const;
517             //- Particle moment of inertia around diameter axis
518             inline scalar momentOfInertia() const;
520             //- Particle angular velocity
521             inline vector omega() const;
523             //- Particle volume
524             inline scalar volume() const;
526             //- Particle volume for a given diameter
527             inline static scalar volume(const scalar d);
529             //- Particle projected area
530             inline scalar areaP() const;
532             //- Projected area for given diameter
533             inline static scalar areaP(const scalar d);
535             //- Particle surface area
536             inline scalar areaS() const;
538             //- Surface area for given diameter
539             inline static scalar areaS(const scalar d);
541             //- Reynolds number
542             inline scalar Re
543             (
544                 const vector& U,        // particle velocity
545                 const scalar d,         // particle diameter
546                 const scalar rhoc,      // carrier density
547                 const scalar muc        // carrier dynamic viscosity
548             ) const;
550             //- Weber number
551             inline scalar We
552             (
553                 const vector& U,        // particle velocity
554                 const scalar d,         // particle diameter
555                 const scalar rhoc,      // carrier density
556                 const scalar sigma      // particle surface tension
557             ) const;
560         // Main calculation loop
562             //- Set cell values
563             template<class TrackData>
564             void setCellValues
565             (
566                 TrackData& td,
567                 const scalar dt,
568                 const label cellI
569             );
571             //- Correct cell values using latest transfer information
572             template<class TrackData>
573             void cellValueSourceCorrection
574             (
575                 TrackData& td,
576                 const scalar dt,
577                 const label cellI
578             );
580             //- Update parcel properties over the time interval
581             template<class TrackData>
582             void calc
583             (
584                 TrackData& td,
585                 const scalar dt,
586                 const label cellI
587             );
590         // Tracking
592             //- Move the parcel
593             template<class TrackData>
594             bool move(TrackData& td, const scalar trackTime);
597         // Patch interactions
599             //- Overridable function to handle the particle hitting a face
600             //  without trackData
601             void hitFace(int& td);
603             //- Overridable function to handle the particle hitting a face
604             template<class TrackData>
605             void hitFace(TrackData& td);
607             //- Overridable function to handle the particle hitting a patch
608             //  Executed before other patch-hitting functions
609             template<class TrackData>
610             bool hitPatch
611             (
612                 const polyPatch& p,
613                 TrackData& td,
614                 const label patchI,
615                 const scalar trackFraction,
616                 const tetIndices& tetIs
617             );
619             //- Overridable function to handle the particle hitting a
620             //  processorPatch
621             template<class TrackData>
622             void hitProcessorPatch
623             (
624                 const processorPolyPatch&,
625                 TrackData& td
626             );
628             //- Overridable function to handle the particle hitting a wallPatch
629             template<class TrackData>
630             void hitWallPatch
631             (
632                 const wallPolyPatch&,
633                 TrackData& td,
634                 const tetIndices&
635             );
637             //- Overridable function to handle the particle hitting a polyPatch
638             template<class TrackData>
639             void hitPatch
640             (
641                 const polyPatch&,
642                 TrackData& td
643             );
645             //- Transform the physical properties of the particle
646             //  according to the given transformation tensor
647             virtual void transformProperties(const tensor& T);
649             //- Transform the physical properties of the particle
650             //  according to the given separation vector
651             virtual void transformProperties(const vector& separation);
653             //- The nearest distance to a wall that the particle can be
654             //  in the n direction
655             virtual scalar wallImpactDistance(const vector& n) const;
658         // I-O
660             //- Read
661             template<class CloudType>
662             static void readFields(CloudType& c);
664             //- Write
665             template<class CloudType>
666             static void writeFields(const CloudType& c);
669     // Ostream Operator
671         friend Ostream& operator<< <ParcelType>
672         (
673             Ostream&,
674             const KinematicParcel<ParcelType>&
675         );
679 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
681 } // End namespace Foam
683 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
685 #include "KinematicParcelI.H"
686 #include "KinematicParcelTrackingDataI.H"
688 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
690 #ifdef NoRepository
691     #include "KinematicParcel.C"
692 #endif
694 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
696 #endif
698 // ************************************************************************* //