Initial commit for version 2.0.x patch release
[OpenFOAM-2.0.x.git] / src / lagrangian / intermediate / parcels / Templates / ThermoParcel / ThermoParcel.H
blobfa3c05e80925069adf8fbd3c4f5694fca068891c
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2008-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::ThermoParcel
27 Description
28     Thermodynamic parcel class with one/two-way coupling with the continuous
29     phase. Includes Kinematic parcel sub-models, plus:
30     - heat transfer
32 SourceFiles
33     ThermoParcelI.H
34     ThermoParcel.C
35     ThermoParcelIO.C
37 \*---------------------------------------------------------------------------*/
39 #ifndef ThermoParcel_H
40 #define ThermoParcel_H
42 #include "particle.H"
43 #include "SLGThermo.H"
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 namespace Foam
50 template<class ParcelType>
51 class ThermoParcel;
53 template<class ParcelType>
54 Ostream& operator<<
56     Ostream&,
57     const ThermoParcel<ParcelType>&
60 /*---------------------------------------------------------------------------*\
61                        Class ThermoParcel Declaration
62 \*---------------------------------------------------------------------------*/
64 template<class ParcelType>
65 class ThermoParcel
67     public ParcelType
69 public:
71     //- Class to hold thermo particle constant properties
72     class constantProperties
73     :
74         public ParcelType::constantProperties
75     {
77         // Private data
79             //- Particle initial temperature [K]
80             scalar T0_;
82             //- Minimum temperature [K]
83             scalar TMin_;
85             //- Particle specific heat capacity [J/(kg.K)]
86             scalar Cp0_;
88             //- Particle emissivity [] (radiation)
89             scalar epsilon0_;
91             //- Particle scattering factor [] (radiation)
92             scalar f0_;
94             //- Default carrier Prandtl number []
95             scalar Pr_;
98     public:
100         // Constructors
102             //- Null constructor
103             constantProperties();
105             //- Copy constructor
106             constantProperties(const constantProperties& cp);
108             //- Constructor from dictionary
109             constantProperties
110             (
111                 const dictionary& parentDict,
112                 const bool readFields = true
113             );
115             //- Construct from components
116             constantProperties
117             (
118                 const label parcelTypeId,
119                 const scalar rhoMin,
120                 const scalar rho0,
121                 const scalar minParticleMass,
122                 const scalar youngsModulus,
123                 const scalar poissonsRatio,
124                 const scalar T0,
125                 const scalar TMin,
126                 const scalar Cp0,
127                 const scalar epsilon0,
128                 const scalar f0,
129                 const scalar Pr
130             );
133         // Member functions
135             // Access
137                 //- Return const access to the particle initial temperature [K]
138                 inline scalar T0() const;
140                 //- Return const access to minimum temperature [K]
141                 inline scalar TMin() const;
143                 //- Return const access to the particle specific heat capacity
144                 //  [J/(kg.K)]
145                 inline scalar Cp0() const;
147                 //- Return const access to the particle emissivity []
148                 //  Active for radiation only
149                 inline scalar epsilon0() const;
151                 //- Return const access to the particle scattering factor []
152                 //  Active for radiation only
153                 inline scalar f0() const;
155                 //- Return const access to the default carrier Prandtl number []
156                 inline scalar Pr() const;
157     };
160     template<class CloudType>
161     class TrackingData
162     :
163         public ParcelType::template TrackingData<CloudType>
164     {
165     private:
167         // Private data
169             //- Local copy of carrier specific heat field
170             //  Cp not stored on carrier thermo, but returned as tmp<...>
171             const volScalarField Cp_;
174             // Interpolators for continuous phase fields
176                 //- Temperature field interpolator
177                 autoPtr<interpolation<scalar> > TInterp_;
179                 //- Specific heat capacity field interpolator
180                 autoPtr<interpolation<scalar> > CpInterp_;
182                 //- Radiation field interpolator
183                 autoPtr<interpolation<scalar> > GInterp_;
187     public:
189         typedef typename ParcelType::template TrackingData<CloudType>::trackPart
190             trackPart;
192         // Constructors
194             //- Construct from components
195             inline TrackingData
196             (
197                 CloudType& cloud,
198                 trackPart part = ParcelType::template
199                     TrackingData<CloudType>::tpLinearTrack
200             );
203         // Member functions
205             //- Return access to the locally stored carrier Cp field
206             inline const volScalarField& Cp() const;
208             //- Return const access to the interpolator for continuous
209             //  phase temperature field
210             inline const interpolation<scalar>& TInterp() const;
212             //- Return const access to the interpolator for continuous
213             //  phase specific heat capacity field
214             inline const interpolation<scalar>& CpInterp() const;
216             //- Return const access to the interpolator for continuous
217             //  radiation field
218             inline const interpolation<scalar>& GInterp() const;
219     };
222 protected:
224     // Protected data
226         // Parcel properties
228             //- Temperature [K]
229             scalar T_;
231             //- Specific heat capacity [J/(kg.K)]
232             scalar Cp_;
235         // Cell-based quantities
237             //- Temperature [K]
238             scalar Tc_;
240             //- Specific heat capacity [J/(kg.K)]
241             scalar Cpc_;
244     // Protected Member Functions
246         //- Calculate new particle temperature
247         template<class TrackData>
248         scalar calcHeatTransfer
249         (
250             TrackData& td,
251             const scalar dt,           // timestep
252             const label cellI,         // owner cell
253             const scalar Re,           // Reynolds number
254             const scalar Pr,           // Prandtl number - surface
255             const scalar kappa,        // Thermal conductivity - surface
256             const scalar NCpW,         // Sum of N*Cp*W of emission species
257             const scalar Sh,           // explicit particle enthalpy source
258             scalar& dhsTrans,          // sensible enthalpy transfer to carrier
259             scalar& Sph                // linearised heat transfer coefficient
260         );
263 public:
265     // Static data members
267         //- String representation of properties
268         static string propHeader;
270         //- Runtime type information
271         TypeName("ThermoParcel");
274     // Constructors
276         //- Construct from owner, position, and cloud owner
277         //  Other properties initialised as null
278         inline ThermoParcel
279         (
280             const polyMesh& mesh,
281             const vector& position,
282             const label cellI,
283             const label tetFaceI,
284             const label tetPtI
285         );
287         //- Construct from components
288         inline ThermoParcel
289         (
290             const polyMesh& mesh,
291             const vector& position,
292             const label cellI,
293             const label tetFaceI,
294             const label tetPtI,
295             const label typeId,
296             const scalar nParticle0,
297             const scalar d0,
298             const scalar dTarget0,
299             const vector& U0,
300             const vector& f0,
301             const vector& angularMomentum0,
302             const vector& torque0,
303             const constantProperties& constProps
304         );
306         //- Construct from Istream
307         ThermoParcel
308         (
309             const polyMesh& mesh,
310             Istream& is,
311             bool readFields = true
312         );
314         //- Construct as a copy
315         ThermoParcel(const ThermoParcel& p);
317         //- Construct as a copy
318         ThermoParcel(const ThermoParcel& p, const polyMesh& mesh);
320         //- Construct and return a (basic particle) clone
321         virtual autoPtr<particle> clone() const
322         {
323             return autoPtr<particle>(new ThermoParcel(*this));
324         }
326         //- Construct and return a (basic particle) clone
327         virtual autoPtr<particle> clone(const polyMesh& mesh) const
328         {
329             return autoPtr<particle>(new ThermoParcel(*this, mesh));
330         }
332         //- Factory class to read-construct particles used for
333         //  parallel transfer
334         class iNew
335         {
336             const polyMesh& mesh_;
338         public:
340             iNew(const polyMesh& mesh)
341             :
342                 mesh_(mesh)
343             {}
345             autoPtr<ThermoParcel<ParcelType> > operator()(Istream& is) const
346             {
347                 return autoPtr<ThermoParcel<ParcelType> >
348                 (
349                     new ThermoParcel<ParcelType>(mesh_, is, true)
350                 );
351             }
352         };
355     // Member Functions
357         // Access
359             //- Return const access to temperature
360             inline scalar T() const;
362             //- Return const access to specific heat capacity
363             inline scalar Cp() const;
365             //- Return the parcel sensible enthalpy
366             inline scalar hs() const;
368             //- Return const access to carrier temperature
369             inline scalar Tc() const;
371             //- Return const access to carrier specific heat capacity
372             inline scalar Cpc() const;
375         // Edit
377             //- Return access to temperature
378             inline scalar& T();
380             //- Return access to specific heat capacity
381             inline scalar& Cp();
384         // Main calculation loop
386             //- Set cell values
387             template<class TrackData>
388             void setCellValues
389             (
390                 TrackData& td,
391                 const scalar dt,
392                 const label cellI
393             );
395             //- Correct cell values using latest transfer information
396             template<class TrackData>
397             void cellValueSourceCorrection
398             (
399                 TrackData& td,
400                 const scalar dt,
401                 const label cellI
402             );
404             //- Calculate surface thermo properties
405             template<class TrackData>
406             void calcSurfaceValues
407             (
408                 TrackData& td,
409                 const label cellI,
410                 const scalar T,
411                 scalar& Ts,
412                 scalar& rhos,
413                 scalar& mus,
414                 scalar& Pr,
415                 scalar& kappas
416             ) const;
418             //- Update parcel properties over the time interval
419             template<class TrackData>
420             void calc
421             (
422                 TrackData& td,
423                 const scalar dt,
424                 const label cellI
425             );
428         // I-O
430             //- Read
431             template<class CloudType>
432             static void readFields(CloudType& c);
434             //- Write
435             template<class CloudType>
436             static void writeFields(const CloudType& c);
439     // Ostream Operator
441         friend Ostream& operator<< <ParcelType>
442         (
443             Ostream&,
444             const ThermoParcel<ParcelType>&
445         );
449 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
451 } // End namespace Foam
453 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
455 #include "ThermoParcelI.H"
456 #include "ThermoParcelTrackingDataI.H"
458 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
460 #ifdef NoRepository
461     #include "ThermoParcel.C"
462 #endif
464 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
466 #endif
468 // ************************************************************************* //