fvDOM updates
[OpenFOAM-1.6.x.git] / src / thermophysicalModels / radiation / submodels / absorptionEmissionModel / greyMeanAbsorptionEmission / greyMeanAbsorptionEmission.C
blob6204188f1aa045c098c7516b435a63cb548fbb76
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 "greyMeanAbsorptionEmission.H"
28 #include "addToRunTimeSelectionTable.H"
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 namespace Foam
34     namespace radiation
35     {
36         defineTypeNameAndDebug(greyMeanAbsorptionEmission, 0);
38         addToRunTimeSelectionTable
39         (
40             absorptionEmissionModel,
41             greyMeanAbsorptionEmission,
42             dictionary
43         );
44     }
48 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
50 Foam::radiation::greyMeanAbsorptionEmission::greyMeanAbsorptionEmission
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     EhrrCoeff_(readScalar(coeffsDict_.lookup("EhrrCoeff"))),
68     Yj_(nSpecies_)
70     label nFunc = 0;
71     const dictionary& functionDicts = dict.subDict(typeName + "Coeffs");
73     forAllConstIter(dictionary, functionDicts, iter)
74     {
75         // safety:
76         if (!iter().isDict())
77         {
78             continue;
79         }
80         const word& key = iter().keyword();
81         speciesNames_.insert(key, nFunc);
82         const dictionary& dict = iter().dict();
83         coeffs_[nFunc].initialise(dict);
84         nFunc++;
85     }
87     // Check that all the species on the dictionary are present in the
88     // look-up table and save the corresponding indices of the look-up table
90     label j = 0;
91     forAllConstIter(HashTable<label>, speciesNames_, iter)
92     {
93         if (mesh.foundObject<volScalarField>("ft"))
94         {
95             if (lookUpTable_.found(iter.key()))
96             {
97                 label index = lookUpTable_.findFieldIndex(iter.key());
99                 Info<< "specie: " << iter.key() << " found on look-up table "
100                     << " with index: " << index << endl;
102                 specieIndex_[iter()] = index;
103             }
104             else if (mesh.foundObject<volScalarField>(iter.key()))
105             {
106                 volScalarField& Y =
107                     const_cast<volScalarField&>
108                     (
109                         mesh.lookupObject<volScalarField>(iter.key())
110                     );
111                 Yj_.set(j, &Y);
112                 specieIndex_[iter()] = 0;
113                 j++;
114                 Info<< "specie: " << iter.key() << " is being solved" << endl;
115             }
116             else
117             {
118                 FatalErrorIn
119                 (
120                     "Foam::radiation::greyMeanAbsorptionEmission(const"
121                     "dictionary& dict, const fvMesh& mesh)"
122                 )   << "specie: " << iter.key()
123                     << " is neither in look-up table: "
124                     << lookUpTable_.tableName()
125                     << " nor is being solved" << nl
126                     << exit(FatalError);
127             }
128         }
129         else
130         {
131             FatalErrorIn
132             (
133                 "Foam::radiation::greyMeanAbsorptionEmission(const"
134                 "dictionary& dict, const fvMesh& mesh)"
135             )   << "specie ft is not present " << nl
136                 << exit(FatalError);
138         }
139     }
142 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
144 Foam::radiation::greyMeanAbsorptionEmission::~greyMeanAbsorptionEmission()
148 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
150 Foam::tmp<Foam::volScalarField>
151 Foam::radiation::greyMeanAbsorptionEmission::aCont(const label bandI) const
153     const volScalarField& T = thermo_.T();
154     const volScalarField& p = thermo_.p();
155     const volScalarField& ft = mesh_.lookupObject<volScalarField>("ft");
157     label nSpecies = speciesNames_.size();
159     tmp<volScalarField> ta
160     (
161         new volScalarField
162         (
163             IOobject
164             (
165                 "a",
166                 mesh().time().timeName(),
167                 mesh(),
168                 IOobject::NO_READ,
169                 IOobject::NO_WRITE
170             ),
171             mesh(),
172             dimensionedScalar("a", dimless/dimLength, 0.0)
173         )
174     );
176     scalarField& a = ta().internalField();
178     forAll(a, i)
179     {
180         const List<scalar>& species = lookUpTable_.lookUp(ft[i]);
182         for (label n=0; n<nSpecies; n++)
183         {
184             label l = 0;
185             scalar Yipi = 0;
186             if (specieIndex_[n] != 0)
187             {
188                 //moles x pressure [atm]
189                 Yipi = species[specieIndex_[n]]*p[i]*9.869231e-6;
190             }
191             else
192             {
193                 // mass fraction
194                 Yipi = Yj_[l][i];
195                 l++;
196             }
198             const absorptionCoeffs::coeffArray& b = coeffs_[n].coeffs(T[i]);
200             scalar Ti = T[i];
201             // negative temperature exponents
202             if (coeffs_[n].invTemp())
203             {
204                 Ti = 1./T[i];
205             }
206             a[i] +=
207                 Yipi
208                *(
209                     ((((b[5]*Ti + b[4])*Ti + b[3])*Ti + b[2])*Ti + b[1])*Ti
210                   + b[0]
211                 );
212         }
213     }
214     return ta;
218 Foam::tmp<Foam::volScalarField>
219 Foam::radiation::greyMeanAbsorptionEmission::eCont(const label bandI) const
221     tmp<volScalarField> e
222     (
223         new volScalarField
224         (
225             IOobject
226             (
227                 "e",
228                 mesh().time().timeName(),
229                 mesh(),
230                 IOobject::NO_READ,
231                 IOobject::NO_WRITE
232             ),
233             mesh(),
234             dimensionedScalar("e", dimless/dimLength, 0.0)
235         )
236     );
238     return e;
242 Foam::tmp<Foam::volScalarField>
243 Foam::radiation::greyMeanAbsorptionEmission::ECont(const label bandI) const
245     tmp<volScalarField> E
246     (
247         new volScalarField
248         (
249             IOobject
250             (
251                 "E",
252                 mesh_.time().timeName(),
253                 mesh_,
254                 IOobject::NO_READ,
255                 IOobject::NO_WRITE
256             ),
257             mesh_,
258             dimensionedScalar("E", dimMass/dimLength/pow3(dimTime), 0.0)
259         )
260     );
262     if (mesh_.foundObject<volScalarField>("dQ"))
263     {
264         const volScalarField& dQ =
265             mesh_.lookupObject<volScalarField>("dQ");
266         E().internalField() = EhrrCoeff_*dQ;
267     }
269     return E;
273 // ************************************************************************* //