fixed writing out entries in advective bc
[OpenFOAM-1.6-ext.git] / src / lagrangian / intermediate / parcels / Templates / ThermoParcel / ThermoParcel.H
blob56baeeed4584b680be572407c3168611265e7201
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright held by original author
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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     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, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 Class
26     Foam::ThermoParcel
28 Description
29     Thermodynamic parcel class with one/two-way coupling with the continuous
30     phase. Includes Kinematic parcel sub-models, plus:
31     - heat transfer
33 SourceFiles
34     ThermoParcelI.H
35     ThermoParcel.C
36     ThermoParcelIO.C
38 \*---------------------------------------------------------------------------*/
40 #ifndef ThermoParcel_H
41 #define ThermoParcel_H
43 #include "IOstream.H"
44 #include "autoPtr.H"
45 #include "interpolationCellPoint.H"
46 #include "contiguous.H"
48 #include "KinematicParcel.H"
49 #include "ThermoCloud.H"
51 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53 namespace Foam
56 template<class ParcelType>
57 class ThermoParcel;
59 template<class ParcelType>
60 Ostream& operator<<
62     Ostream&,
63     const ThermoParcel<ParcelType>&
66 /*---------------------------------------------------------------------------*\
67                        Class ThermoParcel Declaration
68 \*---------------------------------------------------------------------------*/
70 template<class ParcelType>
71 class ThermoParcel
73     public KinematicParcel<ParcelType>
75 public:
77     //- Class to hold thermo particle constant properties
78     class constantProperties
79     :
80         public KinematicParcel<ParcelType>::constantProperties
81     {
83         // Private data
85             //- Particle initial temperature [K]
86             const scalar T0_;
88             //- Minimum temperature [K]
89             const scalar TMin_;
91             //- Particle specific heat capacity [J/(kg.K)]
92             const scalar cp0_;
94             //- Particle emissivity [] (radiation)
95             const scalar epsilon0_;
97             //- Particle scattering factor [] (radiation)
98             const scalar f0_;
100             //- Default carrier Prandtl number []
101             const scalar Pr_;
104     public:
106         // Constructors
107         constantProperties(const dictionary& parentDict);
109         // Member functions
111             // Access
113                 //- Return const access to the particle initial temperature [K]
114                 inline scalar T0() const;
116                 //- Return const access to minimum temperature [K]
117                 inline scalar TMin() const;
119                 //- Return const access to the particle specific heat capacity
120                 //  [J/(kg.K)]
121                 inline scalar cp0() const;
123                 //- Return const access to the particle emissivity []
124                 //  Active for radiation only
125                 inline scalar epsilon0() const;
127                 //- Return const access to the particle scattering factor []
128                 //  Active for radiation only
129                 inline scalar f0() const;
131                 //- Return const access to the default carrier Prandtl number []
132                 inline scalar Pr() const;
133     };
136     //- Class used to pass thermo tracking data to the trackToFace function
137     class trackData
138     :
139         public KinematicParcel<ParcelType>::trackData
140     {
142         // Private data
144             //- Reference to the cloud containing this particle
145             ThermoCloud<ParcelType>& cloud_;
147             //- Particle constant properties
148             const constantProperties& constProps_;
150             // Interpolators for continuous phase fields
152                 //- Temperature field interpolator
153                 const interpolation<scalar>& TInterp_;
155                 //- Specific heat capacity field interpolator
156                 const interpolation<scalar>& cpInterp_;
159     public:
161         // Constructors
163             //- Construct from components
164             inline trackData
165             (
166                 ThermoCloud<ParcelType>& cloud,
167                 const constantProperties& constProps,
168                 const interpolation<scalar>& rhoInterp,
169                 const interpolation<vector>& UInterp,
170                 const interpolation<scalar>& muInterp,
171                 const interpolation<scalar>& TInterp,
172                 const interpolation<scalar>& cpInterp,
173                 const vector& g
174             );
177         // Member functions
179             //- Return access to the owner cloud
180             inline ThermoCloud<ParcelType>& cloud();
182             //- Return const access to the owner cloud
183             inline const constantProperties& constProps() const;
185             //- Return const access to the interpolator for continuous
186             //  phase temperature field
187             inline const interpolation<scalar>& TInterp() const;
189             //- Return const access to the interpolator for continuous
190             //  phase specific heat capacity field
191             inline const interpolation<scalar>& cpInterp() const;
192     };
195 protected:
197     // Protected data
199         // Parcel properties
201             //- Temperature [K]
202             scalar T_;
204             //- Specific heat capacity [J/(kg.K)]
205             scalar cp_;
208         // Cell-based quantities
210             //- Temperature [K]
211             scalar Tc_;
213             //- Specific heat capacity [J/(kg.K)]
214             scalar cpc_;
217     // Protected member functions
219         //- Calculate new particle temperature
220         template<class TrackData>
221         scalar calcHeatTransfer
222         (
223             TrackData& td,
224             const scalar dt,           // timestep
225             const label cellI,         // owner cell
226             const scalar Re,           // Reynolds number
227             const scalar Pr,           // Prandtl number - surface
228             const scalar kappa,        // Thermal conductivity - surface
229             const scalar d,            // diameter
230             const scalar rho,          // density
231             const scalar T,            // temperature
232             const scalar cp,           // specific heat capacity
233             const scalar NCpW,         // Sum of N*Cp*W of emission species
234             const scalar Sh,           // explicit particle enthalpy source
235             scalar& dhsTrans           // sensible enthalpy transfer to carrier
236         );
239 public:
241     // Static data members
243         //- String representation of properties
244         static string propHeader;
246         //- Runtime type information
247         TypeName("ThermoParcel");
250     friend class Cloud<ParcelType>;
253     // Constructors
255         //- Construct from owner, position, and cloud owner
256         //  Other properties initialised as null
257         inline ThermoParcel
258         (
259             ThermoCloud<ParcelType>& owner,
260             const vector& position,
261             const label cellI
262         );
264         //- Construct from components
265         inline ThermoParcel
266         (
267             ThermoCloud<ParcelType>& owner,
268             const vector& position,
269             const label cellI,
270             const label typeId,
271             const scalar nParticle0,
272             const scalar d0,
273             const vector& U0,
274             const constantProperties& constProps
275         );
277         //- Construct from Istream
278         ThermoParcel
279         (
280             const Cloud<ParcelType>& c,
281             Istream& is,
282             bool readFields = true
283         );
285         //- Construct as a copy
286         ThermoParcel(const ThermoParcel& p);
288         //- Construct and return a clone
289         autoPtr<ThermoParcel> clone() const
290         {
291             return autoPtr<ThermoParcel>(new ThermoParcel(*this));
292         }
295     // Member Functions
297         // Access
299             //- Return const access to temperature
300             inline scalar T() const;
302             //- Return const access to specific heat capacity
303             inline scalar cp() const;
306         // Edit
308             //- Return access to temperature
309             inline scalar& T();
311             //- Return access to specific heat capacity
312             inline scalar& cp();
315         // Main calculation loop
317             //- Set cell values
318             template<class TrackData>
319             void setCellValues
320             (
321                 TrackData& td,
322                 const scalar dt,
323                 const label cellI
324             );
326             //- Correct cell values using latest transfer information
327             template<class TrackData>
328             void cellValueSourceCorrection
329             (
330                 TrackData& td,
331                 const scalar dt,
332                 const label cellI
333             );
335             //- Calculate surface thermo properties
336             template<class TrackData>
337             void calcSurfaceValues
338             (
339                 TrackData& td,
340                 const label cellI,
341                 const scalar T,
342                 scalar& Ts,
343                 scalar& rhos,
344                 scalar& mus,
345                 scalar& Pr,
346                 scalar& kappa
347             ) const;
349             //- Update parcel properties over the time interval
350             template<class TrackData>
351             void calc
352             (
353                 TrackData& td,
354                 const scalar dt,
355                 const label cellI
356             );
359         // I-O
361             //- Read
362             static void readFields(Cloud<ParcelType>& c);
364             //- Write
365             static void writeFields(const Cloud<ParcelType>& c);
368     // Ostream Operator
370         friend Ostream& operator<< <ParcelType>
371         (
372             Ostream&,
373             const ThermoParcel<ParcelType>&
374         );
378 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
380 } // End namespace Foam
382 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
384 #include "ThermoParcelI.H"
386 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
388 #ifdef NoRepository
389     #include "ThermoParcel.C"
390 #endif
392 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
394 #endif
396 // ************************************************************************* //