DOC: Corrected class names in the file descriptions
[freefoam.git] / src / thermophysicalModels / radiation / submodels / absorptionEmissionModel / wideBandAbsorptionEmission / wideBandAbsorptionEmission.C
blob5556c7c6cfd799268796474d614b9cff3fa6c5b2
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2008-2010 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 \*---------------------------------------------------------------------------*/
26 #include "wideBandAbsorptionEmission.H"
27 #include <OpenFOAM/addToRunTimeSelectionTable.H>
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 namespace Foam
33     namespace radiation
34     {
35         defineTypeNameAndDebug(wideBandAbsorptionEmission, 0);
37         addToRunTimeSelectionTable
38         (
39             absorptionEmissionModel,
40             wideBandAbsorptionEmission,
41             dictionary
42         );
43     }
47 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
49 Foam::radiation::wideBandAbsorptionEmission::wideBandAbsorptionEmission
51     const dictionary& dict,
52     const fvMesh& mesh
55     absorptionEmissionModel(dict, mesh),
56     coeffsDict_((dict.subDict(typeName + "Coeffs"))),
57     speciesNames_(0),
58     specieIndex_(0),
59     lookUpTable_
60     (
61         fileName(coeffsDict_.lookup("lookUpTableFileName")),
62         mesh.time().constant(),
63         mesh
64     ),
65     thermo_(mesh.lookupObject<basicThermo>("thermophysicalProperties")),
66     Yj_(nSpecies_),
67     totalWaveLength_(0)
69     label nBand = 0;
70     const dictionary& functionDicts = dict.subDict(typeName +"Coeffs");
71     forAllConstIter(dictionary, functionDicts, iter)
72     {
73         // safety:
74         if (!iter().isDict())
75         {
76             continue;
77         }
79         const dictionary& dict = iter().dict();
80         dict.lookup("bandLimits") >> iBands_[nBand];
81         dict.lookup("EhrrCoeff") >> iEhrrCoeffs_[nBand];
82         totalWaveLength_ += iBands_[nBand][1] - iBands_[nBand][0];
84         label nSpec = 0;
85         const dictionary& specDicts = dict.subDict("species");
86         forAllConstIter(dictionary, specDicts, iter)
87         {
88             const word& key = iter().keyword();
89             if (nBand == 0)
90             {
91                 speciesNames_.insert(key, nSpec);
92             }
93             else
94             {
95                 if (!speciesNames_.found(key))
96                 {
97                     FatalErrorIn
98                     (
99                         "Foam::radiation::wideBandAbsorptionEmission(const"
100                         "dictionary& dict, const fvMesh& mesh)"
101                     )   << "specie: " << key << "is not in all the bands"
102                         << nl << exit(FatalError);
103                 }
104             }
105             coeffs_[nSpec][nBand].initialise(specDicts.subDict(key));
106             nSpec++;
107         }
108         nBand++;
109     }
110     nBands_ = nBand;
112     // Check that all the species on the dictionary are present in the
113     // look-up table  and save the corresponding indices of the look-up table
115     label j = 0;
116     forAllConstIter(HashTable<label>, speciesNames_, iter)
117     {
118         if (lookUpTable_.found(iter.key()))
119         {
120             label index = lookUpTable_.findFieldIndex(iter.key());
121             Info<< "specie: " << iter.key() << " found in look-up table "
122                 << " with index: " << index << endl;
123             specieIndex_[iter()] = index;
124         }
125         else if (mesh.foundObject<volScalarField>(iter.key()))
126         {
127             volScalarField& Y = const_cast<volScalarField&>
128                 (mesh.lookupObject<volScalarField>(iter.key()));
130             Yj_.set(j, &Y);
132             specieIndex_[iter()] = 0.0;
133             j++;
134             Info<< "species: " << iter.key() << " is being solved" << endl;
135         }
136         else
137         {
138             FatalErrorIn
139             (
140                 "radiation::wideBandAbsorptionEmission(const"
141                 "dictionary& dict, const fvMesh& mesh)"
142             )   << "specie: " << iter.key()
143                 << " is neither in look-up table : "
144                 << lookUpTable_.tableName() << " nor is being solved"
145                 << exit(FatalError);
146         }
147     }
152 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
154 Foam::radiation::wideBandAbsorptionEmission::~wideBandAbsorptionEmission()
158 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
160 Foam::tmp<Foam::volScalarField>
161 Foam::radiation::wideBandAbsorptionEmission::aCont(const label bandI) const
163     const volScalarField& T = thermo_.T();
164     const volScalarField& p = thermo_.p();
165     const volScalarField& ft = mesh_.lookupObject<volScalarField>("ft");
167     label nSpecies = speciesNames_.size();
169     tmp<volScalarField> ta
170     (
171         new volScalarField
172         (
173             IOobject
174             (
175                 "a",
176                 mesh().time().timeName(),
177                 mesh(),
178                 IOobject::NO_READ,
179                 IOobject::NO_WRITE
180             ),
181             mesh(),
182             dimensionedScalar("a", dimless/dimLength, 0.0)
183         )
184     );
186     scalarField& a = ta().internalField();
188     forAll(a, i)
189     {
190         const List<scalar>& species = lookUpTable_.lookUp(ft[i]);
192         for (label n=0; n<nSpecies; n++)
193         {
194             label l = 0;
195             scalar Yipi = 0.0;
196             if (specieIndex_[n] != 0)
197             {
198                 // moles x pressure [atm]
199                 Yipi = species[specieIndex_[n]]*p[i]*9.869231e-6;
200             }
201             else
202             {
203                 // mass fraction from species being solved
204                 Yipi = Yj_[l][i];
205                 l++;
206             }
208             scalar Ti = T[i];
210             const absorptionCoeffs::coeffArray& b =
211                 coeffs_[n][bandI].coeffs(T[i]);
213             if (coeffs_[n][bandI].invTemp())
214             {
215                 Ti = 1.0/T[i];
216             }
218             a[i]+=
219                 Yipi
220                *(
221                     ((((b[5]*Ti + b[4])*Ti + b[3])*Ti + b[2])*Ti + b[1])*Ti
222                   + b[0]
223                 );
224         }
225     }
227     return ta;
231 Foam::tmp<Foam::volScalarField>
232 Foam::radiation::wideBandAbsorptionEmission::eCont(const label bandI) const
234     tmp<volScalarField> e
235     (
236         new volScalarField
237         (
238             IOobject
239             (
240                 "e",
241                 mesh().time().timeName(),
242                 mesh(),
243                 IOobject::NO_READ,
244                 IOobject::NO_WRITE
245             ),
246             mesh(),
247             dimensionedScalar("e", dimless/dimLength, 0.0)
248         )
249     );
251     return e;
255 Foam::tmp<Foam::volScalarField>
256 Foam::radiation::wideBandAbsorptionEmission::ECont(const label bandI) const
258     tmp<volScalarField> E
259     (
260         new volScalarField
261         (
262             IOobject
263             (
264                 "E",
265                 mesh().time().timeName(),
266                 mesh(),
267                 IOobject::NO_READ,
268                 IOobject::NO_WRITE
269             ),
270             mesh(),
271             dimensionedScalar("E", dimMass/dimLength/pow3(dimTime), 0.0)
272         )
273     );
275     if (mesh().foundObject<volScalarField>("hrr"))
276     {
277         const volScalarField& hrr = mesh().lookupObject<volScalarField>("hrr");
278         E().internalField() =
279             iEhrrCoeffs_[bandI]
280            *hrr.internalField()
281            *(iBands_[bandI][1] - iBands_[bandI][0])
282            /totalWaveLength_;
283     }
285     return E;
288 Foam::tmp<Foam::volScalarField>
289 Foam::radiation::wideBandAbsorptionEmission::addIntensity
291     const label i,
292     const volScalarField& ILambda
293 ) const
295     return ILambda*(iBands_[i][1] - iBands_[i][0])/totalWaveLength_;
299 void Foam::radiation::wideBandAbsorptionEmission::correct
301     volScalarField& a,
302     PtrList<volScalarField>& aLambda
303 ) const
305     a = dimensionedScalar("zero", dimless/dimLength, 0.0);
307     for (label j=0; j<nBands_; j++)
308     {
309         Info<< "Calculating absorption in band: " << j << endl;
310         aLambda[j].internalField() = this->a(j);
311         Info<< "Calculated absorption in band: " << j << endl;
312         a.internalField() +=
313             aLambda[j].internalField()
314            *(iBands_[j][1] - iBands_[j][0])
315            /totalWaveLength_;
316     }
321 // ************************ vim: set sw=4 sts=4 et: ************************ //