BUG: UListIO: byteSize overflowing on really big faceLists
[OpenFOAM-2.0.x.git] / src / lagrangian / intermediate / parcels / Templates / ReactingParcel / ReactingParcel.H
blobc69c7e3d3021e96d3f787b979b88582603f94cc7
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::ReactingParcel
27 Description
28     Reacting parcel class with one/two-way coupling with the continuous
29     phase.
31 SourceFiles
32     ReactingParcelI.H
33     ReactingParcel.C
34     ReactingParcelIO.C
36 \*---------------------------------------------------------------------------*/
38 #ifndef ReactingParcel_H
39 #define ReactingParcel_H
41 #include "particle.H"
42 #include "SLGThermo.H"
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 namespace Foam
49 template<class ParcelType>
50 class ReactingParcel;
52 template<class ParcelType>
53 Ostream& operator<<
55     Ostream&,
56     const ReactingParcel<ParcelType>&
59 /*---------------------------------------------------------------------------*\
60                         Class ReactingParcel Declaration
61 \*---------------------------------------------------------------------------*/
63 template<class ParcelType>
64 class ReactingParcel
66     public ParcelType
68 public:
70     //- Class to hold reacting particle constant properties
71     class constantProperties
72     :
73         public ParcelType::constantProperties
74     {
75         // Private data
77             //- Minimum pressure [Pa]
78             scalar pMin_;
80             //- Constant volume flag - e.g. during mass transfer
81             Switch constantVolume_;
83             //- Vaporisation temperature [K]
84             scalar Tvap_;
86             //- Boiling point [K]
87             scalar Tbp_;
90     public:
92         // Constructors
94             //- Null constructor
95             constantProperties();
97             //- Copy constructor
98             constantProperties(const constantProperties& cp);
100             //- Constructor from dictionary
101             constantProperties
102             (
103                 const dictionary& parentDict,
104                 const bool readFields = true
105             );
107             //- Construct from components
108             constantProperties
109             (
110                 const label parcelTypeId,
111                 const scalar rhoMin,
112                 const scalar rho0,
113                 const scalar minParticleMass,
114                 const scalar youngsModulus,
115                 const scalar poissonsRatio,
116                 const scalar T0,
117                 const scalar TMin,
118                 const scalar Cp0,
119                 const scalar epsilon0,
120                 const scalar f0,
121                 const scalar Pr,
122                 const scalar pMin,
123                 const Switch& constantVolume,
124                 const scalar Tvap,
125                 const scalar Tbp
126             );
129         // Access
131             //- Return const access to the minimum pressure
132             inline scalar pMin() const;
134             //- Return const access to the constant volume flag
135             inline Switch constantVolume() const;
137             //- Return const access to the vaporisation temperature
138             inline scalar Tvap() const;
140             //- Return const access to the boiling point
141             inline scalar Tbp() const;
142     };
145     template<class CloudType>
146     class TrackingData
147     :
148         public ParcelType::template TrackingData<CloudType>
149     {
150     private:
152         // Private data
154             // Interpolators for continuous phase fields
156                 //- Interpolator for continuous phase pressure field
157                 autoPtr<interpolation<scalar> > pInterp_;
160     public:
162         typedef typename ParcelType::template TrackingData<CloudType>::trackPart
163             trackPart;
165         // Constructors
167             //- Construct from components
168             inline TrackingData
169             (
170                 CloudType& cloud,
171                 trackPart part = ParcelType::template
172                     TrackingData<CloudType>::tpLinearTrack
173             );
176         // Member functions
178             //- Return const access to the interpolator for continuous phase
179             //  pressure field
180             inline const interpolation<scalar>& pInterp() const;
181     };
184 protected:
186     // Protected data
188         // Parcel properties
190             //- Initial particle mass [kg]
191             scalar mass0_;
193             //- Mass fractions of mixture []
194             scalarField Y_;
197         // Cell-based quantities
199             //- Pressure [Pa]
200             scalar pc_;
203     // Protected Member Functions
205         //- Calculate Phase change
206         template<class TrackData>
207         void calcPhaseChange
208         (
209             TrackData& td,
210             const scalar dt,           // timestep
211             const label cellI,         // owner cell
212             const scalar Re,           // Reynolds number
213             const scalar Ts,           // Surface temperature
214             const scalar nus,          // Surface kinematic viscosity
215             const scalar d,            // diameter
216             const scalar T,            // temperature
217             const scalar mass,         // mass
218             const label idPhase,       // id of phase involved in phase change
219             const scalar YPhase,       // total mass fraction
220             const scalarField& YComponents, // component mass fractions
221             scalarField& dMassPC,      // mass transfer - local to particle
222             scalar& Sh,                // explicit particle enthalpy source
223             scalar& N,                 // flux of species emitted from particle
224             scalar& NCpW,              // sum of N*Cp*W of emission species
225             scalarField& Cs            // carrier conc. of emission species
226         );
228         //- Update mass fraction
229         scalar updateMassFraction
230         (
231             const scalar mass0,
232             const scalarField& dMass,
233             scalarField& Y
234         ) const;
237 public:
239     // Static data members
241         //- String representation of properties
242         static string propHeader;
244         //- Runtime type information
245         TypeName("ReactingParcel");
248     // Constructors
250         //- Construct from owner, position, and cloud owner
251         //  Other properties initialised as null
252         inline ReactingParcel
253         (
254             const polyMesh& mesh,
255             const vector& position,
256             const label cellI,
257             const label tetFaceI,
258             const label tetPtI
259         );
261         //- Construct from components
262         inline ReactingParcel
263         (
264             const polyMesh& mesh,
265             const vector& position,
266             const label cellI,
267             const label tetFaceI,
268             const label tetPtI,
269             const label typeId,
270             const scalar nParticle0,
271             const scalar d0,
272             const scalar dTarget0,
273             const vector& U0,
274             const vector& f0,
275             const vector& angularMomentum0,
276             const vector& torque0,
277             const scalarField& Y0,
278             const constantProperties& constProps
279         );
281         //- Construct from Istream
282         ReactingParcel
283         (
284             const polyMesh& mesh,
285             Istream& is,
286             bool readFields = true
287         );
289         //- Construct as a copy
290         ReactingParcel
291         (
292             const ReactingParcel& p,
293             const polyMesh& mesh
294         );
296         //- Construct as a copy
297         ReactingParcel(const ReactingParcel& p);
299         //- Construct and return a (basic particle) clone
300         virtual autoPtr<particle> clone() const
301         {
302             return autoPtr<particle>(new ReactingParcel<ParcelType>(*this));
303         }
305         //- Construct and return a (basic particle) clone
306         virtual autoPtr<particle> clone(const polyMesh& mesh) const
307         {
308             return autoPtr<particle>
309             (
310                 new ReactingParcel<ParcelType>(*this, mesh)
311             );
312         }
314         //- Factory class to read-construct particles used for
315         //  parallel transfer
316         class iNew
317         {
318             const polyMesh& mesh_;
320         public:
322             iNew(const polyMesh& mesh)
323             :
324                 mesh_(mesh)
325             {}
327             autoPtr<ReactingParcel<ParcelType> > operator()(Istream& is) const
328             {
329                 return autoPtr<ReactingParcel<ParcelType> >
330                 (
331                     new ReactingParcel<ParcelType>(mesh_, is, true)
332                 );
333             }
334         };
337     // Member Functions
339         // Access
341             //- Return const access to initial mass
342             inline scalar mass0() const;
344             //- Return const access to mass fractions of mixture
345             inline const scalarField& Y() const;
347             //- Return the owner cell pressure
348             inline scalar pc() const;
350             //- Return reference to the owner cell pressure
351             inline scalar& pc();
354         // Edit
356             //- Return access to initial mass
357             inline scalar& mass0();
359             //- Return access to mass fractions of mixture
360             inline scalarField& Y();
363         // Main calculation loop
365             //- Set cell values
366             template<class TrackData>
367             void setCellValues
368             (
369                 TrackData& td,
370                 const scalar dt,
371                 const label cellI
372             );
374             //- Correct cell values using latest transfer information
375             template<class TrackData>
376             void cellValueSourceCorrection
377             (
378                 TrackData& td,
379                 const scalar dt,
380                 const label cellI
381             );
383             //- Correct surface values due to emitted species
384             template<class TrackData>
385             void correctSurfaceValues
386             (
387                 TrackData& td,
388                 const label cellI,
389                 const scalar T,
390                 const scalarField& Cs,
391                 scalar& rhos,
392                 scalar& mus,
393                 scalar& Prs,
394                 scalar& kappas
395             );
397             //- Update parcel properties over the time interval
398             template<class TrackData>
399             void calc
400             (
401                 TrackData& td,
402                 const scalar dt,
403                 const label cellI
404             );
407         // I-O
409             //- Read
410             template<class CloudType, class CompositionType>
411             static void readFields
412             (
413                 CloudType& c,
414                 const CompositionType& compModel
415             );
417             //- Read - no composition
418             template<class CloudType>
419             static void readFields(CloudType& c);
421             //- Write
422             template<class CloudType, class CompositionType>
423             static void writeFields
424             (
425                 const CloudType& c,
426                 const CompositionType& compModel
427             );
429             //- Write - composition supplied
430             template<class CloudType>
431             static void writeFields(const CloudType& c);
434     // Ostream Operator
436         friend Ostream& operator<< <ParcelType>
437         (
438             Ostream&,
439             const ReactingParcel<ParcelType>&
440         );
444 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
446 } // End namespace Foam
448 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
450 #include "ReactingParcelI.H"
451 #include "ReactingParcelTrackingDataI.H"
453 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
455 #ifdef NoRepository
456     #include "ReactingParcel.C"
457 #endif
459 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
461 #endif
463 // ************************************************************************* //