Initial commit for version 2.0.x patch release
[OpenFOAM-2.0.x.git] / src / lagrangian / spray / parcels / Templates / SprayParcel / SprayParcel.H
blob9dc9b94e114aa1173acdf7fcd717abac55d2a78a
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011-2011 OpenCFD Ltd.
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::SprayParcel
27 Description
28     Reacing spray parcel, with added functionality for atomization and breakup
30 \*---------------------------------------------------------------------------*/
32 #ifndef SprayParcel_H
33 #define SprayParcel_H
35 #include "particle.H"
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 namespace Foam
42 template<class ParcelType>
43 class SprayParcel;
45 template<class ParcelType>
46 Ostream& operator<<
48     Ostream&,
49     const SprayParcel<ParcelType>&
52 /*---------------------------------------------------------------------------*\
53                           Class SprayParcel Declaration
54 \*---------------------------------------------------------------------------*/
56 template<class ParcelType>
57 class SprayParcel
59     public ParcelType
62 protected:
64     // Protected data
66         // Spray parcel properties
68             //- Initial droplet diameter
69             scalar d0_;
71             //- Injection position
72             vector position0_;
74             //- Part of liquid core ( >0.5=liquid, <0.5=droplet )
75             scalar liquidCore_;
77             //- Index for KH Breakup
78             scalar KHindex_;
80             //- Spherical deviation
81             scalar y_;
83             //- Rate of change of spherical deviation
84             scalar yDot_;
86             //- Characteristic time (used in atomization and/or breakup model)
87             scalar tc_;
89             //- Stripped parcel mass due to breakup
90             scalar ms_;
92             //- Injected from injector (needed e.g. for calculating distance
93             //  from injector)
94             scalar injector_;
96             //- Momentum relaxation time (needed for calculating parcel acc.)
97             scalar tMom_;
99             //- Passive scalar (extra variable to be defined by user)
100             scalar user_;
103 public:
105     // Static data members
107         //- String representation of properties
108         static string propHeader;
110         //- Runtime type information
111         TypeName("SprayParcel");
114     // Constructors
116         //- Construct from owner, position, and cloud owner
117         //  Other properties initialised as null
118         inline SprayParcel
119         (
120             const polyMesh& mesh,
121             const vector& position,
122             const label cellI,
123             const label tetFaceI,
124             const label tetPtI
125         );
127         //- Construct from components
128         inline SprayParcel
129         (
130             const polyMesh& mesh,
131             const vector& position,
132             const label cellI,
133             const label tetFaceI,
134             const label tetPtI,
135             const label typeId,
136             const scalar nParticle0,
137             const scalar d0,
138             const scalar dTarget0,
139             const vector& U0,
140             const vector& f0,
141             const vector& angularMomentum0,
142             const vector& torque0,
143             const scalarField& Y0,
144             const scalar liquidCore,
145             const scalar KHindex,
146             const scalar y,
147             const scalar yDot,
148             const scalar tc,
149             const scalar ms,
150             const scalar injector,
151             const scalar tMom,
152             const scalar user,
153             const typename ParcelType::constantProperties& constProps
154         );
156         //- Construct from Istream
157         SprayParcel
158         (
159             const polyMesh& mesh,
160             Istream& is,
161             bool readFields = true
162         );
164         //- Construct as a copy
165         SprayParcel
166         (
167             const SprayParcel& p,
168             const polyMesh& mesh
169         );
171         //- Construct as a copy
172         SprayParcel(const SprayParcel& p);
174         //- Construct and return a (basic particle) clone
175         virtual autoPtr<particle> clone() const
176         {
177             return autoPtr<particle>(new SprayParcel<ParcelType>(*this));
178         }
180         //- Construct and return a (basic particle) clone
181         virtual autoPtr<particle> clone(const polyMesh& mesh) const
182         {
183             return autoPtr<particle>
184             (
185                 new SprayParcel<ParcelType>(*this, mesh)
186             );
187         }
189         //- Factory class to read-construct particles used for
190         //  parallel transfer
191         class iNew
192         {
193             const polyMesh& mesh_;
195         public:
197             iNew(const polyMesh& mesh)
198             :
199                 mesh_(mesh)
200             {}
202             autoPtr<SprayParcel<ParcelType> > operator()(Istream& is) const
203             {
204                 return autoPtr<SprayParcel<ParcelType> >
205                 (
206                     new SprayParcel<ParcelType>(mesh_, is, true)
207                 );
208             }
209         };
212     // Member Functions
214         // Access
216             //- Return const access to initial droplet diameter
217             inline scalar d0() const;
219             //- Return const access to initial droplet position
220             inline const vector& position0() const;
222             //- Return const access to liquid core
223             inline scalar liquidCore() const;
225             //- Return const access to Kelvin-Helmholtz breakup index
226             inline scalar KHindex() const;
228             //- Return const access to spherical deviation
229             inline scalar y() const;
231             //- Return const access to rate of change of spherical deviation
232             inline scalar yDot() const;
234             //- Return const access to atomization characteristic time
235             inline scalar tc() const;
237             //- Return const access to stripped parcel mass
238             inline scalar ms() const;
240             //- Return const access to injector id
241             inline scalar injector() const;
243             //- Return const access to momentum relaxation time
244             inline scalar tMom() const;
246             //- Return const access to passive user scalar
247             inline scalar user() const;
250         // Edit
252             //- Return access to initial droplet diameter
253             inline scalar& d0();
255             //- Return access to initial droplet position
256             inline vector& position0();
258             //- Return access to liquid core
259             inline scalar& liquidCore();
261             //- Return access to Kelvin-Helmholtz breakup index
262             inline scalar& KHindex();
264             //- Return access to spherical deviation
265             inline scalar& y();
267             //- Return access to rate of change of spherical deviation
268             inline scalar& yDot();
270             //- Return access to atomization characteristic time
271             inline scalar& tc();
273             //- Return access to stripped parcel mass
274             inline scalar& ms();
276             //- Return access to injector id
277             inline scalar& injector();
279             //- Return access to momentum relaxation time
280             inline scalar& tMom();
282             //- Return access to passive user scalar
283             inline scalar& user();
286         // Main calculation loop
288             //- Set cell values
289             template<class TrackData>
290             void setCellValues
291             (
292                 TrackData& td,
293                 const scalar dt,
294                 const label cellI
295             );
297             //- Correct parcel properties according to atomization model
298             template<class TrackData>
299             void calcAtomization
300             (
301                 TrackData& td,
302                 const scalar dt,
303                 const label cellI
304             );
306             //- Correct parcel properties according to breakup model
307             template<class TrackData>
308             void calcBreakup
309             (
310                 TrackData& td,
311                 const scalar dt,
312                 const label cellI
313             );
315             //- Correct cell values using latest transfer information
316             template<class TrackData>
317             void cellValueSourceCorrection
318             (
319                 TrackData& td,
320                 const scalar dt,
321                 const label cellI
322             );
324             //- Correct surface values due to emitted species
325             template<class TrackData>
326             void correctSurfaceValues
327             (
328                 TrackData& td,
329                 const label cellI,
330                 const scalar T,
331                 const scalarField& Cs,
332                 scalar& rhos,
333                 scalar& mus,
334                 scalar& Pr,
335                 scalar& kappa
336             );
338             //- Update parcel properties over the time interval
339             template<class TrackData>
340             void calc
341             (
342                 TrackData& td,
343                 const scalar dt,
344                 const label cellI
345             );
347             //- Calculate the chi-factor for flash-boiling for the
348             //  atomization model
349             template<class TrackData>
350             scalar chi
351             (
352                 TrackData& td,
353                 const scalarField& X
354             ) const;
356             //- Solve the TAB equation
357             template<class TrackData>
358             void solveTABEq
359             (
360                 TrackData& td,
361                 const scalar dt
362             );
365         // I-O
367             //- Read
368             template<class CloudType, class CompositionType>
369             static void readFields
370             (
371                 CloudType& c,
372                 const CompositionType& compModel
373             );
375             //- Read - no composition
376             template<class CloudType>
377             static void readFields(CloudType& c);
379             //- Write
380             template<class CloudType, class CompositionType>
381             static void writeFields
382             (
383                 const CloudType& c,
384                 const CompositionType& compModel
385             );
387             //- Write - composition supplied
388             template<class CloudType>
389             static void writeFields(const CloudType& c);
392     // Ostream Operator
394         friend Ostream& operator<< <ParcelType>
395         (
396             Ostream&,
397             const SprayParcel<ParcelType>&
398         );
402 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
404 } // End namespace Foam
406 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
408 #include "SprayParcelI.H"
410 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
412 #ifdef NoRepository
413     #include "SprayParcel.C"
414 #endif
417 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
419 #endif
421 // ************************************************************************* //