initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / thermophysicalModels / radiation / submodels / absorptionEmissionModel / wideBandAbsorptionEmission / wideBandAbsorptionEmission.C
blob1665519e3e711fdaa8c80493c77f91c1be972fdb
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2008-2009 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 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 \*---------------------------------------------------------------------------*/
27 #include "wideBandAbsorptionEmission.H"
28 #include "addToRunTimeSelectionTable.H"
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 namespace Foam
34     namespace radiation
35     {
36         defineTypeNameAndDebug(wideBandAbsorptionEmission, 0);
38         addToRunTimeSelectionTable
39         (
40             absorptionEmissionModel,
41             wideBandAbsorptionEmission,
42             dictionary
43         );
44     }
48 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
50 Foam::radiation::wideBandAbsorptionEmission::wideBandAbsorptionEmission
52     const dictionary& dict,
53     const fvMesh& mesh
56     absorptionEmissionModel(dict, mesh),
57     coeffsDict_((dict.subDict(typeName + "Coeffs"))),
58     speciesNames_(0),
59     specieIndex_(0),
60     lookUpTable_
61     (
62         fileName(coeffsDict_.lookup("lookUpTableFileName")),
63         mesh.time().constant(),
64         mesh
65     ),
66     thermo_(mesh.lookupObject<basicThermo>("thermophysicalProperties")),
67     Yj_(nSpecies_),
68     totalWaveLength_(0)
70     label nBand = 0;
71     const dictionary& functionDicts = dict.subDict(typeName +"Coeffs");
72     forAllConstIter(dictionary, functionDicts, iter)
73     {
74         // safety:
75         if (!iter().isDict())
76         {
77             continue;
78         }
80         const dictionary& dict = iter().dict();
81         dict.lookup("bandLimits") >> iBands_[nBand];
82         dict.lookup("EhrrCoeff") >> iEhrrCoeffs_[nBand];
83         totalWaveLength_ += iBands_[nBand][1] - iBands_[nBand][0];
85         label nSpec = 0;
86         const dictionary& specDicts = dict.subDict("species");
87         forAllConstIter(dictionary, specDicts, iter)
88         {
89             const word& key = iter().keyword();
90             if (nBand == 0)
91             {
92                 speciesNames_.insert(key, nSpec);
93             }
94             else
95             {
96                 if (!speciesNames_.found(key))
97                 {
98                     FatalErrorIn
99                     (
100                         "Foam::radiation::wideBandAbsorptionEmission(const"
101                         "dictionary& dict, const fvMesh& mesh)"
102                     )   << "specie: " << key << "is not in all the bands"
103                         << nl << exit(FatalError);
104                 }
105             }
106             coeffs_[nSpec][nBand].initialise(specDicts.subDict(key));
107             nSpec++;
108         }
109         nBand++;
110     }
111     nBands_ = nBand;
113     // Check that all the species on the dictionary are present in the
114     // look-up table  and save the corresponding indices of the look-up table
116     label j = 0;
117     forAllConstIter(HashTable<label>, speciesNames_, iter)
118     {
119         if (lookUpTable_.found(iter.key()))
120         {
121             label index = lookUpTable_.findFieldIndex(iter.key());
122             Info<< "specie: " << iter.key() << " found in look-up table "
123                 << " with index: " << index << endl;
124             specieIndex_[iter()] = index;
125         }
126         else if (mesh.foundObject<volScalarField>(iter.key()))
127         {
128             volScalarField& Y = const_cast<volScalarField&>
129                 (mesh.lookupObject<volScalarField>(iter.key()));
131             Yj_.set(j, &Y);
133             specieIndex_[iter()] = 0.0;
134             j++;
135             Info<< "species: " << iter.key() << " is being solved" << endl;
136         }
137         else
138         {
139             FatalErrorIn
140             (
141                 "radiation::wideBandAbsorptionEmission(const"
142                 "dictionary& dict, const fvMesh& mesh)"
143             )   << "specie: " << iter.key()
144                 << " is neither in look-up table : "
145                 << lookUpTable_.tableName() << " nor is being solved"
146                 << exit(FatalError);
147         }
148     }
153 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
155 Foam::radiation::wideBandAbsorptionEmission::~wideBandAbsorptionEmission()
159 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
161 Foam::tmp<Foam::volScalarField>
162 Foam::radiation::wideBandAbsorptionEmission::aCont(const label bandI) const
164     const volScalarField& T = thermo_.T();
165     const volScalarField& p = thermo_.p();
166     const volScalarField& ft = mesh_.lookupObject<volScalarField>("ft");
168     label nSpecies = speciesNames_.size();
170     tmp<volScalarField> ta
171     (
172         new volScalarField
173         (
174             IOobject
175             (
176                 "a",
177                 mesh().time().timeName(),
178                 mesh(),
179                 IOobject::NO_READ,
180                 IOobject::NO_WRITE
181             ),
182             mesh(),
183             dimensionedScalar("a", dimless/dimLength, 0.0)
184         )
185     );
187     scalarField& a = ta().internalField();
189     forAll(a, i)
190     {
191         const List<scalar>& species = lookUpTable_.lookUp(ft[i]);
193         for (label n=0; n<nSpecies; n++)
194         {
195             label l = 0;
196             scalar Yipi = 0.0;
197             if (specieIndex_[n] != 0)
198             {
199                 // moles x pressure [atm]
200                 Yipi = species[specieIndex_[n]]*p[i]*9.869231e-6;
201             }
202             else
203             {
204                 // mass fraction from species being solved
205                 Yipi = Yj_[l][i];
206                 l++;
207             }
209             scalar Ti = T[i];
211             const absorptionCoeffs::coeffArray& b =
212                 coeffs_[n][bandI].coeffs(T[i]);
214             if (coeffs_[n][bandI].invTemp())
215             {
216                 Ti = 1.0/T[i];
217             }
219             a[i]+=
220                 Yipi
221                *(
222                     ((((b[5]*Ti + b[4])*Ti + b[3])*Ti + b[2])*Ti + b[1])*Ti
223                   + b[0]
224                 );
225         }
226     }
228     return ta;
232 Foam::tmp<Foam::volScalarField>
233 Foam::radiation::wideBandAbsorptionEmission::eCont(const label bandI) const
235     tmp<volScalarField> e
236     (
237         new volScalarField
238         (
239             IOobject
240             (
241                 "e",
242                 mesh().time().timeName(),
243                 mesh(),
244                 IOobject::NO_READ,
245                 IOobject::NO_WRITE
246             ),
247             mesh(),
248             dimensionedScalar("e", dimless/dimLength, 0.0)
249         )
250     );
252     return e;
256 Foam::tmp<Foam::volScalarField>
257 Foam::radiation::wideBandAbsorptionEmission::ECont(const label bandI) const
259     tmp<volScalarField> E
260     (
261         new volScalarField
262         (
263             IOobject
264             (
265                 "E",
266                 mesh().time().timeName(),
267                 mesh(),
268                 IOobject::NO_READ,
269                 IOobject::NO_WRITE
270             ),
271             mesh(),
272             dimensionedScalar("E", dimMass/dimLength/pow3(dimTime), 0.0)
273         )
274     );
276     if (mesh().foundObject<volScalarField>("hrr"))
277     {
278         const volScalarField& hrr = mesh().lookupObject<volScalarField>("hrr");
279         E().internalField() =
280             iEhrrCoeffs_[bandI]
281            *hrr.internalField()
282            *(iBands_[bandI][1] - iBands_[bandI][0])
283            /totalWaveLength_;
284     }
286     return E;
289 Foam::tmp<Foam::volScalarField>
290 Foam::radiation::wideBandAbsorptionEmission::addIntensity
292     const label i,
293     const volScalarField& ILambda
294 ) const
296     return ILambda*(iBands_[i][1] - iBands_[i][0])/totalWaveLength_;
300 void Foam::radiation::wideBandAbsorptionEmission::correct
302     volScalarField& a,
303     PtrList<volScalarField>& aLambda
304 ) const
306     a = dimensionedScalar("zero", dimless/dimLength, 0.0);
308     for (label j=0; j<nBands_; j++)
309     {
310         Info<< "Calculating absorption in band: " << j << endl;
311         aLambda[j].internalField() = this->a(j);
312         Info<< "Calculated absorption in band: " << j << endl;
313         a.internalField() +=
314             aLambda[j].internalField()
315            *(iBands_[j][1] - iBands_[j][0])
316            /totalWaveLength_;
317     }
322 // ************************************************************************* //