OpenFOAM logo
Open Source CFD Toolkit

C8H10.H

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 1991-2005 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software; you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by the
00013     Free Software Foundation; either version 2 of the License, or (at your
00014     option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM; if not, write to the Free Software Foundation,
00023     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
00024 
00025 Class
00026     C8H10
00027 
00028 Description
00029     ethylBenzene
00030 
00031 SourceFiles
00032     C8H10.C
00033 
00034 \*---------------------------------------------------------------------------*/
00035 
00036 #ifndef C8H10_H
00037 #define C8H10_H
00038 
00039 #include "liquid.H"
00040 #include "NSRDSfunc0.H"
00041 #include "NSRDSfunc1.H"
00042 #include "NSRDSfunc2.H"
00043 #include "NSRDSfunc3.H"
00044 #include "NSRDSfunc4.H"
00045 #include "NSRDSfunc5.H"
00046 #include "NSRDSfunc6.H"
00047 #include "NSRDSfunc7.H"
00048 #include "APIdiffCoefFunc.H"
00049 
00050 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00051 
00052 namespace Foam
00053 {
00054 
00055 /*---------------------------------------------------------------------------*\
00056                            Class C8H10 Declaration
00057 \*---------------------------------------------------------------------------*/
00058 
00059 class C8H10
00060 :
00061     public liquid
00062 {
00063     // Private data
00064 
00065         NSRDSfunc5 rho_;
00066         NSRDSfunc1 pv_;
00067         NSRDSfunc6 hl_;
00068         NSRDSfunc0 cp_;
00069         NSRDSfunc0 h_;
00070         NSRDSfunc7 cpg_;
00071         NSRDSfunc4 B_;
00072         NSRDSfunc1 mu_;
00073         NSRDSfunc2 mug_;
00074         NSRDSfunc0 K_;
00075         NSRDSfunc2 Kg_;
00076         NSRDSfunc6 sigma_;
00077         APIdiffCoefFunc D_;
00078 
00079 
00080 public:
00081 
00082     //- Runtime type information
00083     TypeName("C8H10");
00084 
00085 
00086     // Constructors
00087 
00088         //- Construct null
00089         C8H10()
00090         :
00091             liquid(106.167, 617.17, 3.6094e+6, 0.37381, 0.263, 178.15, 4.038e-3, 409.35, 1.9680e-30, 0.3036, 1.8043e+4),
00092             rho_(76.3765398, 0.26438, 617.17, 0.2921),
00093             pv_(88.246, -7691.1, -9.797, 5.931e-06, 2),
00094             hl_(617.17, 516167.924119547, 0.3882, 0, 0, 0),
00095             cp_(818.521762883005, 6.66873887366131, -0.0248005500767658, 4.23860521631015e-05, 0, 0),
00096             // NN: enthalpy, h_, is not used in the sprayModel.
00097             // For consistency, the enthalpy is derived from hlat and hl.
00098             // It is, however, convenient to have it available.
00099             h_(-524002.612929508, 818.521762883005, 3.33436943683065, -0.00826685002558862, 1.05965130407754e-05, 0),
00100             cpg_(738.835984816374, 3201.5598067196, 1559, 2285.07916772632, -702),
00101             B_(0.00165776559571242, -2.77958310962917, -388067.855359952, -5.86905535618412e+18, 1.58052878954854e+21),
00102             mu_(-10.452, 1048.4, -0.0715, 0, 0),
00103             mug_(1.2e-06, 0.4518, 439, 0),
00104             K_(0.20149, -0.00023988, 0, 0, 0, 0),
00105             Kg_(1.708e-05, 1.319, 565.6, 0),
00106             sigma_(617.17, 0.066, 1.268, 0, 0, 0),
00107             D_(147.18, 20.1, 106.167, 28) // NN: Same as nHeptane
00108         {}
00109         C8H10
00110         (
00111             const liquid& l,
00112             const NSRDSfunc5& density,
00113             const NSRDSfunc1& vapourPressure,
00114             const NSRDSfunc6& heatOfVapourisation,
00115             const NSRDSfunc0& heatCapacity,
00116             const NSRDSfunc0& enthalpy,
00117             const NSRDSfunc7& idealGasHeatCapacity,
00118             const NSRDSfunc4& secondVirialCoeff,
00119             const NSRDSfunc1& dynamicViscosity,
00120             const NSRDSfunc2& vapourDynamicViscosity,
00121             const NSRDSfunc0& thermalConductivity,
00122             const NSRDSfunc2& vapourThermalConductivity,
00123             const NSRDSfunc6& surfaceTension,
00124             const APIdiffCoefFunc& vapourDiffussivity
00125         )
00126         :
00127             liquid(l),
00128             rho_(density),
00129             pv_(vapourPressure),
00130             hl_(heatOfVapourisation),
00131             cp_(heatCapacity),
00132             h_(enthalpy),
00133             cpg_(idealGasHeatCapacity),
00134             B_(secondVirialCoeff),
00135             mu_(dynamicViscosity),
00136             mug_(vapourDynamicViscosity),
00137             K_(thermalConductivity),
00138             Kg_(vapourThermalConductivity),
00139             sigma_(surfaceTension),
00140             D_(vapourDiffussivity)
00141         {}
00142 
00143         //- Construct from Istream
00144         C8H10(Istream& is)
00145         :
00146             liquid(is),
00147             rho_(is),
00148             pv_(is),
00149             hl_(is),
00150             cp_(is),
00151             h_(is),
00152             cpg_(is),
00153             B_(is),
00154             mu_(is),
00155             mug_(is),
00156             K_(is),
00157             Kg_(is),
00158             sigma_(is),
00159             D_(is)
00160         {}
00161 
00162 
00163     // Member Functions
00164 
00165         //- Liquid density [kg/m^3]
00166         scalar rho(scalar p, scalar T) const
00167         {
00168             return rho_.f(p, T);
00169         }
00170 
00171         //- Vapour pressure [Pa]
00172         scalar pv(scalar p, scalar T) const
00173         {
00174             return pv_.f(p, T);
00175         }
00176 
00177         //- Heat of vapourisation [J/kg]
00178         scalar hl(scalar p, scalar T) const
00179         {
00180             return hl_.f(p, T);
00181         }
00182 
00183         //- Liquid heat capacity [J/(kg K)]
00184         scalar cp(scalar p, scalar T) const
00185         {
00186             return cp_.f(p, T);
00187         }
00188 
00189         //- Liquid Enthalpy [J/(kg)]
00190         scalar h(scalar p, scalar T) const
00191         {
00192             return h_.f(p, T);
00193         }
00194 
00195         //- Ideal gas heat capacity [J/(kg K)]
00196         scalar cpg(scalar p, scalar T) const
00197         {
00198             return cpg_.f(p, T);
00199         }
00200 
00201         //- Second Virial Coefficient [m^3/kg]
00202         scalar B(scalar p, scalar T) const
00203         {
00204             return B_.f(p, T);
00205         }
00206 
00207         //- Liquid viscosity [Pa s]
00208         scalar mu(scalar p, scalar T) const
00209         {
00210             return mu_.f(p, T);
00211         }
00212 
00213         //- Vapour viscosity [Pa s]
00214         scalar mug(scalar p, scalar T) const
00215         {
00216             return mug_.f(p, T);
00217         }
00218 
00219         //- Liquid thermal conductivity  [W/(m K)]
00220         scalar K(scalar p, scalar T) const
00221         {
00222             return K_.f(p, T);
00223         }
00224 
00225         //- Vapour thermal conductivity  [W/(m K)]
00226         scalar Kg(scalar p, scalar T) const
00227         {
00228             return Kg_.f(p, T);
00229         }
00230 
00231         //- Surface tension [N/m]
00232         scalar sigma(scalar p, scalar T) const
00233         {
00234             return sigma_.f(p, T);
00235         }
00236 
00237         //- Vapour diffussivity [m2/s]
00238         scalar D(scalar p, scalar T) const
00239         {
00240             return D_.f(p, T);
00241         }
00242 
00243 
00244         //- Write the function coefficients
00245         void writeData(Ostream& os) const
00246         {
00247             liquid::writeData(os); os << nl;
00248             rho_.writeData(os); os << nl;
00249             pv_.writeData(os); os << nl;
00250             hl_.writeData(os); os << nl;
00251             cp_.writeData(os); os << nl;
00252             cpg_.writeData(os); os << nl;
00253             B_.writeData(os); os << nl;
00254             mu_.writeData(os); os << nl;
00255             mug_.writeData(os); os << nl;
00256             K_.writeData(os); os << nl;
00257             Kg_.writeData(os); os << nl;
00258             sigma_.writeData(os); os << nl;
00259             D_.writeData(os); os << endl;
00260         }
00261 
00262 
00263     // Ostream Operator
00264 
00265         friend Ostream& operator<<(Ostream& os, const C8H10& l)
00266         {
00267             l.writeData(os);
00268             return os;
00269         }
00270 };
00271 
00272 
00273 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00274 
00275 } // End namespace Foam
00276 
00277 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00278 
00279 #endif
00280 
00281 // ************************************************************************* //
For further information go to www.openfoam.org