Merge /u/wyldckat/foam-extend32/ branch master into master
[foam-extend-3.2.git] / src / lagrangian / intermediate / parcels / Templates / KinematicParcel / KinematicParcel.H
blob9f032ee4f4a417027d09df36f7318bb128ddbdd5
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | foam-extend: Open Source CFD
4    \\    /   O peration     | Version:     3.2
5     \\  /    A nd           | Web:         http://www.foam-extend.org
6      \\/     M anipulation  | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
8 License
9     This file is part of foam-extend.
11     foam-extend is free software: you can redistribute it and/or modify it
12     under the terms of the GNU General Public License as published by the
13     Free Software Foundation, either version 3 of the License, or (at your
14     option) any later version.
16     foam-extend is distributed in the hope that it will be useful, but
17     WITHOUT ANY WARRANTY; without even the implied warranty of
18     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19     General Public License for more details.
21     You should have received a copy of the GNU General Public License
22     along with foam-extend.  If not, see <http://www.gnu.org/licenses/>.
24 Class
25     Foam::KinematicParcel
27 Description
28     Kinematic parcel class with one/two-way coupling with the continuous
29     phase.
31     Sub-models include:
32     - drag
33     - turbulent dispersion
34     - wall interactions
36 SourceFiles
37     KinematicParcelI.H
38     KinematicParcel.C
39     KinematicParcelIO.C
41 \*---------------------------------------------------------------------------*/
43 #ifndef KinematicParcel_H
44 #define KinematicParcel_H
46 #include "Particle.H"
47 #include "IOstream.H"
48 #include "autoPtr.H"
49 #include "interpolationCellPoint.H"
50 #include "contiguous.H"
52 #include "KinematicCloudTemplate.H"
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 Particle<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             //- Minimum density [kg/m3]
91             const scalar rhoMin_;
93             //- Particle density [kg/m3] (constant)
94             const scalar rho0_;
96             //- Minimum particle mass [kg]
97             const scalar minParticleMass_;
100     public:
102         //- Constructor
103         constantProperties(const dictionary& parentDict);
105         // Member functions
107             //- Return const access to the constant properties dictionary
108             inline const dictionary& dict() const;
110             //- Return const access to the minimum density
111             inline scalar rhoMin() const;
113             //- Return const access to the particle density
114             inline scalar rho0() const;
116             //- Return const access to the minimum particle mass
117             inline scalar minParticleMass() const;
118     };
121     //- Class used to pass kinematic tracking data to the trackToFace function
122     class trackData
123     :
124         public Particle<ParcelType>::trackData
125     {
126         // Private data
128             //- Reference to the cloud containing this particle
129             KinematicCloud<ParcelType>& cloud_;
131             //- Particle constant properties
132             const constantProperties& constProps_;
135             // Interpolators for continuous phase fields
137                 //- Density interpolator
138                 const interpolation<scalar>& rhoInterp_;
140                 //- Velocity interpolator
141                 const interpolation<vector>& UInterp_;
143                 //- Dynamic viscosity interpolator
144                 const interpolation<scalar>& muInterp_;
146             //- Local gravitational or other body-force acceleration
147             const vector& g_;
150    public:
152         // Constructors
154            //- Construct from components
155            inline trackData
156             (
157                 KinematicCloud<ParcelType>& cloud,
158                 const constantProperties& constProps,
159                 const interpolation<scalar>& rhoInterp,
160                 const interpolation<vector>& UInterp,
161                 const interpolation<scalar>& muInterp,
162                 const vector& g
163             );
166         // Member functions
168             //- Return access to the owner cloud
169             inline KinematicCloud<ParcelType>& cloud();
171             //- Return const access to the constant properties
172             inline const constantProperties& constProps() const;
174             //- Return conat access to the interpolator for continuous
175             //  phase density field
176             inline const interpolation<scalar>& rhoInterp() const;
178             //- Return conat access to the interpolator for continuous
179             //  phase velocity field
180             inline const interpolation<vector>& UInterp() const;
182             //- Return conat access to the interpolator for continuous
183             //  phase dynamic viscosity field
184             inline const interpolation<scalar>& muInterp() const;
186             // Return const access to the gravitational acceleration vector
187             inline const vector& g() const;
188     };
191 protected:
193     // Protected data
195         // Parcel properties
197             //- Parcel type id
198             label typeId_;
200             //- Number of particles in Parcel
201             scalar nParticle_;
203             //- Diameter [m]
204             scalar d_;
206             //- Velocity of Parcel [m/s]
207             vector U_;
209             //- Density [kg/m3]
210             scalar rho_;
212             //- Time spent in turbulent eddy [s]
213             scalar tTurb_;
215             //- Turbulent velocity fluctuation [m/s]
216             vector UTurb_;
219         // Cell-based quantities
221             //- Density [kg/m3]
222             scalar rhoc_;
224             //- Velocity [m/s]
225             vector Uc_;
227             //- Viscosity [Pa.s]
228             scalar muc_;
231     // Protected member functions
233         //- Calculate new particle velocity
234         template<class TrackData>
235         const vector calcVelocity
236         (
237             TrackData& td,
238             const scalar dt,           // timestep
239             const label cellI,         // owner cell
240             const scalar Re,           // Reynolds number
241             const scalar mu,           // local carrier viscosity
242             const scalar d,            // diameter
243             const vector& U,           // velocity
244             const scalar rho,          // density
245             const scalar mass,         // mass
246             const vector& Su,          // explicit particle momentum source
247             vector& dUTrans            // momentum transfer to carrier
248         ) const;
251 public:
253     // Static data members
255         //- String representation of properties
256         static string propHeader;
258         //- Runtime type information
259         TypeName("KinematicParcel");
262     friend class Cloud<ParcelType>;
265     // Constructors
267         //- Construct from owner, position, and cloud owner
268         //  Other properties initialised as null
269         inline KinematicParcel
270         (
271             KinematicCloud<ParcelType>& owner,
272             const vector& position,
273             const label cellI
274         );
276         //- Construct from components
277         inline KinematicParcel
278         (
279             KinematicCloud<ParcelType>& owner,
280             const vector& position,
281             const label cellI,
282             const label typeId,
283             const scalar nParticle0,
284             const scalar d0,
285             const vector& U0,
286             const constantProperties& constProps
287         );
289         //- Construct from Istream
290         KinematicParcel
291         (
292             const Cloud<ParcelType>& c,
293             Istream& is,
294             bool readFields = true
295         );
297         //- Construct as a copy
298         KinematicParcel(const KinematicParcel& p);
300         //- Construct and return a clone
301         autoPtr<KinematicParcel> clone() const
302         {
303             return autoPtr<KinematicParcel>(new KinematicParcel(*this));
304         }
307     // Member Functions
309         // Access
311             //- Return const access to type id
312             inline label typeId() const;
314             //- Return const access to number of particles
315             inline scalar nParticle() const;
317             //- Return const access to diameter
318             inline scalar d() const;
320             //- Return const access to velocity
321             inline const vector& U() const;
323             //- Return const access to density
324             inline scalar rho() const;
326             //- Return const access to time spent in turbulent eddy
327             inline scalar tTurb() const;
329             //- Return const access to turbulent velocity fluctuation
330             inline const vector& UTurb() const;
333         // Edit
335             //- Return access to type id
336             inline label typeId();
338             //- Return access to number of particles
339             inline scalar& nParticle();
341             //- Return access to diameter
342             inline scalar& d();
344             //- Return access to velocity
345             inline vector& U();
347             //- Return access to density
348             inline scalar& rho();
350             //- Return access to time spent in turbulent eddy
351             inline scalar& tTurb();
353             //- Return access to turbulent velocity fluctuation
354             inline vector& UTurb();
357         // Helper functions
359             //- The nearest distance to a wall that
360             //  the particle can be in the n direction
361             inline scalar wallImpactDistance(const vector& n) const;
363             //- Return the index of the face to be used in the interpolation
364             //  routine
365             inline label faceInterpolation() const;
367             //- Cell owner mass
368             inline scalar massCell(const label cellI) const;
370             //- Particle mass
371             inline scalar mass() const;
373             //- Particle volume
374             inline scalar volume() const;
376             //- Particle volume for a given diameter
377             inline scalar volume(const scalar d) const;
379             //- Particle projected area
380             inline scalar areaP() const;
382             //- Projected area for given diameter
383             inline scalar areaP(const scalar d) const;
385             //- Particle surface area
386             inline scalar areaS() const;
388             //- Surface area for given diameter
389             inline scalar areaS(const scalar d) const;
391             //- Reynolds number
392             inline scalar Re
393             (
394                 const vector& U,        // particle velocity
395                 const scalar d,         // particle diameter
396                 const scalar rhoc,      // carrier density
397                 const scalar muc        // carrier dynamic viscosity
398             ) const;
401         // Main calculation loop
403             //- Set cell values
404             template<class TrackData>
405             void setCellValues
406             (
407                 TrackData& td,
408                 const scalar dt,
409                 const label cellI
410             );
412             //- Correct cell values using latest transfer information
413             template<class TrackData>
414             void cellValueSourceCorrection
415             (
416                 TrackData& td,
417                 const scalar dt,
418                 const label cellI
419             );
421             //- Update parcel properties over the time interval
422             template<class TrackData>
423             void calc
424             (
425                 TrackData& td,
426                 const scalar dt,
427                 const label cellI
428             );
431         // Tracking
433             //- Move the parcel
434             template<class TrackData>
435             bool move(TrackData& td);
438         // Patch interactions
440             //- Overridable function to handle the particle hitting a patch
441             //  Executed before other patch-hitting functions
442             template<class TrackData>
443             bool hitPatch
444             (
445                 const polyPatch& p,
446                 TrackData& td,
447                 const label patchI
448             );
451             //- Overridable function to handle the particle hitting a patch
452             //  Executed before other patch-hitting functions without trackData
453             bool hitPatch
454             (
455                 const polyPatch& p,
456                 int& td,
457                 const label patchI
458             );
461             //- Overridable function to handle the particle hitting a
462             //  processorPatch
463             template<class TrackData>
464             void hitProcessorPatch
465             (
466                 const processorPolyPatch&,
467                 TrackData& td
468             );
470             //- Overridable function to handle the particle hitting a
471             //  processorPatch without trackData
472             void hitProcessorPatch
473             (
474                 const processorPolyPatch&,
475                 int&
476             );
478             //- Overridable function to handle the particle hitting a wallPatch
479             template<class TrackData>
480             void hitWallPatch
481             (
482                 const wallPolyPatch&,
483                 TrackData& td
484             );
486             //- Overridable function to handle the particle hitting a wallPatch
487             //  without trackData
488             void hitWallPatch
489             (
490                 const wallPolyPatch&,
491                 int&
492             );
494             //- Overridable function to handle the particle hitting a polyPatch
495             template<class TrackData>
496             void hitPatch
497             (
498                 const polyPatch&,
499                 TrackData& td
500             );
502             //- Overridable function to handle the particle hitting a polyPatch
503             //- without trackData
504             void hitPatch
505             (
506                 const polyPatch&,
507                 int&
508             );
510             //- Transform the physical properties of the particle
511             //  according to the given transformation tensor
512             void transformProperties(const tensor& T);
514             //- Transform the physical properties of the particle
515             //  according to the given separation vector
516             void transformProperties(const vector& separation);
519         // I-O
521             //- Read
522             static void readFields(Cloud<ParcelType>& c);
524             //- Write
525             static void writeFields(const Cloud<ParcelType>& c);
528     // Ostream Operator
530         friend Ostream& operator<< <ParcelType>
531         (
532             Ostream&,
533             const KinematicParcel<ParcelType>&
534         );
538 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
540 } // End namespace Foam
542 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
544 #include "KinematicParcelI.H"
546 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
548 #define defineParcelTypeNameAndDebug(Type, DebugSwitch)                       \
549     template<>                                                                \
550     const Foam::word KinematicParcel<Type>::typeName(#Type);                  \
551     template<>                                                                \
552     Foam::debug::debugSwitch                                                  \
553     KinematicParcel<Type>::debug                                              \
554     (                                                                         \
555         std::string(#Type), DebugSwitch                                       \
556     );
558 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
560 #ifdef NoRepository
561     #include "KinematicParcel.C"
562 #endif
564 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
566 #endif
568 // ************************************************************************* //