FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

greyMeanAbsorptionEmission.C

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) 2008-2010 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
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your 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, see <http://www.gnu.org/licenses/>.
00023 
00024 \*---------------------------------------------------------------------------*/
00025 
00026 #include "greyMeanAbsorptionEmission.H"
00027 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00028 
00029 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
00030 
00031 namespace Foam
00032 {
00033     namespace radiation
00034     {
00035         defineTypeNameAndDebug(greyMeanAbsorptionEmission, 0);
00036 
00037         addToRunTimeSelectionTable
00038         (
00039             absorptionEmissionModel,
00040             greyMeanAbsorptionEmission,
00041             dictionary
00042         );
00043     }
00044 }
00045 
00046 
00047 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00048 
00049 Foam::radiation::greyMeanAbsorptionEmission::greyMeanAbsorptionEmission
00050 (
00051     const dictionary& dict,
00052     const fvMesh& mesh
00053 )
00054 :
00055     absorptionEmissionModel(dict, mesh),
00056     coeffsDict_((dict.subDict(typeName + "Coeffs"))),
00057     speciesNames_(0),
00058     specieIndex_(0),
00059     lookUpTable_
00060     (
00061         fileName(coeffsDict_.lookup("lookUpTableFileName")),
00062         mesh.time().constant(),
00063         mesh
00064     ),
00065     thermo_(mesh.lookupObject<basicThermo>("thermophysicalProperties")),
00066     EhrrCoeff_(readScalar(coeffsDict_.lookup("EhrrCoeff"))),
00067     Yj_(nSpecies_)
00068 {
00069     label nFunc = 0;
00070     const dictionary& functionDicts = dict.subDict(typeName + "Coeffs");
00071 
00072     forAllConstIter(dictionary, functionDicts, iter)
00073     {
00074         // safety:
00075         if (!iter().isDict())
00076         {
00077             continue;
00078         }
00079         const word& key = iter().keyword();
00080         speciesNames_.insert(key, nFunc);
00081         const dictionary& dict = iter().dict();
00082         coeffs_[nFunc].initialise(dict);
00083         nFunc++;
00084     }
00085 
00086     // Check that all the species on the dictionary are present in the
00087     // look-up table and save the corresponding indices of the look-up table
00088 
00089     label j = 0;
00090     forAllConstIter(HashTable<label>, speciesNames_, iter)
00091     {
00092         if (mesh.foundObject<volScalarField>("ft"))
00093         {
00094             if (lookUpTable_.found(iter.key()))
00095             {
00096                 label index = lookUpTable_.findFieldIndex(iter.key());
00097 
00098                 Info<< "specie: " << iter.key() << " found on look-up table "
00099                     << " with index: " << index << endl;
00100 
00101                 specieIndex_[iter()] = index;
00102             }
00103             else if (mesh.foundObject<volScalarField>(iter.key()))
00104             {
00105                 volScalarField& Y =
00106                     const_cast<volScalarField&>
00107                     (
00108                         mesh.lookupObject<volScalarField>(iter.key())
00109                     );
00110                 Yj_.set(j, &Y);
00111                 specieIndex_[iter()] = 0;
00112                 j++;
00113                 Info<< "specie: " << iter.key() << " is being solved" << endl;
00114             }
00115             else
00116             {
00117                 FatalErrorIn
00118                 (
00119                     "Foam::radiation::greyMeanAbsorptionEmission(const"
00120                     "dictionary& dict, const fvMesh& mesh)"
00121                 )   << "specie: " << iter.key()
00122                     << " is neither in look-up table: "
00123                     << lookUpTable_.tableName()
00124                     << " nor is being solved" << nl
00125                     << exit(FatalError);
00126             }
00127         }
00128         else
00129         {
00130             FatalErrorIn
00131             (
00132                 "Foam::radiation::greyMeanAbsorptionEmission(const"
00133                 "dictionary& dict, const fvMesh& mesh)"
00134             )   << "specie ft is not present " << nl
00135                 << exit(FatalError);
00136 
00137         }
00138     }
00139 }
00140 
00141 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
00142 
00143 Foam::radiation::greyMeanAbsorptionEmission::~greyMeanAbsorptionEmission()
00144 {}
00145 
00146 
00147 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00148 
00149 Foam::tmp<Foam::volScalarField>
00150 Foam::radiation::greyMeanAbsorptionEmission::aCont(const label bandI) const
00151 {
00152     const volScalarField& T = thermo_.T();
00153     const volScalarField& p = thermo_.p();
00154     const volScalarField& ft = mesh_.lookupObject<volScalarField>("ft");
00155 
00156     label nSpecies = speciesNames_.size();
00157 
00158     tmp<volScalarField> ta
00159     (
00160         new volScalarField
00161         (
00162             IOobject
00163             (
00164                 "a",
00165                 mesh().time().timeName(),
00166                 mesh(),
00167                 IOobject::NO_READ,
00168                 IOobject::NO_WRITE
00169             ),
00170             mesh(),
00171             dimensionedScalar("a", dimless/dimLength, 0.0)
00172         )
00173     );
00174 
00175     scalarField& a = ta().internalField();
00176 
00177     forAll(a, i)
00178     {
00179         const List<scalar>& species = lookUpTable_.lookUp(ft[i]);
00180 
00181         for (label n=0; n<nSpecies; n++)
00182         {
00183             label l = 0;
00184             scalar Yipi = 0;
00185             if (specieIndex_[n] != 0)
00186             {
00187                 //moles x pressure [atm]
00188                 Yipi = species[specieIndex_[n]]*p[i]*9.869231e-6;
00189             }
00190             else
00191             {
00192                 // mass fraction
00193                 Yipi = Yj_[l][i];
00194                 l++;
00195             }
00196 
00197             const absorptionCoeffs::coeffArray& b = coeffs_[n].coeffs(T[i]);
00198 
00199             scalar Ti = T[i];
00200             // negative temperature exponents
00201             if (coeffs_[n].invTemp())
00202             {
00203                 Ti = 1./T[i];
00204             }
00205             a[i] +=
00206                 Yipi
00207                *(
00208                     ((((b[5]*Ti + b[4])*Ti + b[3])*Ti + b[2])*Ti + b[1])*Ti
00209                   + b[0]
00210                 );
00211         }
00212     }
00213     return ta;
00214 }
00215 
00216 
00217 Foam::tmp<Foam::volScalarField>
00218 Foam::radiation::greyMeanAbsorptionEmission::eCont(const label bandI) const
00219 {
00220     tmp<volScalarField> e
00221     (
00222         new volScalarField
00223         (
00224             IOobject
00225             (
00226                 "e",
00227                 mesh().time().timeName(),
00228                 mesh(),
00229                 IOobject::NO_READ,
00230                 IOobject::NO_WRITE
00231             ),
00232             mesh(),
00233             dimensionedScalar("e", dimless/dimLength, 0.0)
00234         )
00235     );
00236 
00237     return e;
00238 }
00239 
00240 
00241 Foam::tmp<Foam::volScalarField>
00242 Foam::radiation::greyMeanAbsorptionEmission::ECont(const label bandI) const
00243 {
00244     tmp<volScalarField> E
00245     (
00246         new volScalarField
00247         (
00248             IOobject
00249             (
00250                 "E",
00251                 mesh_.time().timeName(),
00252                 mesh_,
00253                 IOobject::NO_READ,
00254                 IOobject::NO_WRITE
00255             ),
00256             mesh_,
00257             dimensionedScalar("E", dimMass/dimLength/pow3(dimTime), 0.0)
00258         )
00259     );
00260 
00261     if (mesh_.foundObject<volScalarField>("dQ"))
00262     {
00263         const volScalarField& dQ =
00264             mesh_.lookupObject<volScalarField>("dQ");
00265         E().internalField() = EhrrCoeff_*dQ;
00266     }
00267 
00268     return E;
00269 }
00270 
00271 
00272 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines