Forward compatibility: flex
[foam-extend-3.2.git] / src / thermophysicalModels / radiation / submodels / absorptionEmissionModel / greyMeanAbsorptionEmission / greyMeanAbsorptionEmission.C
blob7e9041a44e9a72dafd646b092a7813782731e53a
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | foam-extend: Open Source CFD
4    \\    /   O peration     | Version:     3.2
5     \\  /    A nd           | Web:         http://www.foam-extend.org
6      \\/     M anipulation  | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
8 License
9     This file is part of foam-extend.
11     foam-extend 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 3 of the License, or (at your
14     option) any later version.
16     foam-extend is distributed in the hope that it will be useful, but
17     WITHOUT ANY WARRANTY; without even the implied warranty of
18     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19     General Public License for more details.
21     You should have received a copy of the GNU General Public License
22     along with foam-extend.  If not, see <http://www.gnu.org/licenses/>.
24 \*---------------------------------------------------------------------------*/
26 #include "greyMeanAbsorptionEmission.H"
27 #include "addToRunTimeSelectionTable.H"
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 namespace Foam
33     namespace radiation
34     {
35         defineTypeNameAndDebug(greyMeanAbsorptionEmission, 0);
37         addToRunTimeSelectionTable
38         (
39             absorptionEmissionModel,
40             greyMeanAbsorptionEmission,
41             dictionary
42         );
43     }
47 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
49 Foam::radiation::greyMeanAbsorptionEmission::greyMeanAbsorptionEmission
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     EhrrCoeff_(readScalar(coeffsDict_.lookup("EhrrCoeff"))),
67     Yj_(nSpecies_)
69     label nFunc = 0;
70     const dictionary& functionDicts = dict.subDict(typeName + "Coeffs");
72     forAllConstIter(dictionary, functionDicts, iter)
73     {
74         // safety:
75         if (!iter().isDict())
76         {
77             continue;
78         }
79         const word& key = iter().keyword();
80         speciesNames_.insert(key, nFunc);
81         const dictionary& dict = iter().dict();
82         coeffs_[nFunc].initialise(dict);
83         nFunc++;
84     }
86     // Check that all the species on the dictionary are present in the
87     // look-up table and save the corresponding indices of the look-up table
89     label j = 0;
90     forAllConstIter(HashTable<label>, speciesNames_, iter)
91     {
92         if (mesh.foundObject<volScalarField>("ft"))
93         {
94             if (lookUpTable_.found(iter.key()))
95             {
96                 label index = lookUpTable_.findFieldIndex(iter.key());
98                 Info<< "specie: " << iter.key() << " found on look-up table "
99                     << " with index: " << index << endl;
101                 specieIndex_[iter()] = index;
102             }
103             else if (mesh.foundObject<volScalarField>(iter.key()))
104             {
105                 volScalarField& Y =
106                     const_cast<volScalarField&>
107                     (
108                         mesh.lookupObject<volScalarField>(iter.key())
109                     );
110                 Yj_.set(j, &Y);
111                 specieIndex_[iter()] = 0;
112                 j++;
113                 Info<< "specie: " << iter.key() << " is being solved" << endl;
114             }
115             else
116             {
117                 FatalErrorIn
118                 (
119                     "Foam::radiation::greyMeanAbsorptionEmission(const"
120                     "dictionary& dict, const fvMesh& mesh)"
121                 )   << "specie: " << iter.key()
122                     << " is neither in look-up table: "
123                     << lookUpTable_.tableName()
124                     << " nor is being solved" << nl
125                     << exit(FatalError);
126             }
127         }
128         else
129         {
130             FatalErrorIn
131             (
132                 "Foam::radiation::greyMeanAbsorptionEmission(const"
133                 "dictionary& dict, const fvMesh& mesh)"
134             )   << "specie ft is not present " << nl
135                 << exit(FatalError);
137         }
138     }
141 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
143 Foam::radiation::greyMeanAbsorptionEmission::~greyMeanAbsorptionEmission()
147 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
149 Foam::tmp<Foam::volScalarField>
150 Foam::radiation::greyMeanAbsorptionEmission::aCont(const label bandI) const
152     const volScalarField& T = thermo_.T();
153     const volScalarField& p = thermo_.p();
154     const volScalarField& ft = mesh_.lookupObject<volScalarField>("ft");
156     label nSpecies = speciesNames_.size();
158     tmp<volScalarField> ta
159     (
160         new volScalarField
161         (
162             IOobject
163             (
164                 "a",
165                 mesh().time().timeName(),
166                 mesh(), // HR 101116 -- Todo: Need to bring object registry here!
167                 IOobject::NO_READ,
168                 IOobject::NO_WRITE
169             ),
170             mesh(),
171             dimensionedScalar("a", dimless/dimLength, 0.0)
172         )
173     );
175     scalarField& a = ta().internalField();
177     forAll(a, i)
178     {
179         const List<scalar>& species = lookUpTable_.lookUp(ft[i]);
181         for (label n=0; n<nSpecies; n++)
182         {
183             label l = 0;
184             scalar Yipi = 0;
185             if (specieIndex_[n] != 0)
186             {
187                 //moles x pressure [atm]
188                 Yipi = species[specieIndex_[n]]*p[i]*9.869231e-6;
189             }
190             else
191             {
192                 // mass fraction
193                 Yipi = Yj_[l][i];
194                 l++;
195             }
197             const absorptionCoeffs::coeffArray& b = coeffs_[n].coeffs(T[i]);
199             scalar Ti = T[i];
200             // negative temperature exponents
201             if (coeffs_[n].invTemp())
202             {
203                 Ti = 1./T[i];
204             }
205             a[i] +=
206                 Yipi
207                *(
208                     ((((b[5]*Ti + b[4])*Ti + b[3])*Ti + b[2])*Ti + b[1])*Ti
209                   + b[0]
210                 );
211         }
212     }
213     return ta;
217 Foam::tmp<Foam::volScalarField>
218 Foam::radiation::greyMeanAbsorptionEmission::eCont(const label bandI) const
220     tmp<volScalarField> e
221     (
222         new volScalarField
223         (
224             IOobject
225             (
226                 "e",
227                 mesh().time().timeName(),
228                 mesh(), // HR 101116 -- Todo: Need to bring object registry here!
229                 IOobject::NO_READ,
230                 IOobject::NO_WRITE
231             ),
232             mesh(),
233             dimensionedScalar("e", dimless/dimLength, 0.0)
234         )
235     );
237     return e;
241 Foam::tmp<Foam::volScalarField>
242 Foam::radiation::greyMeanAbsorptionEmission::ECont(const label bandI) const
244     tmp<volScalarField> E
245     (
246         new volScalarField
247         (
248             IOobject
249             (
250                 "E",
251                 mesh_.time().timeName(),
252                 mesh_, // HR 101116 -- Todo: Need to bring object registry here!
253                 IOobject::NO_READ,
254                 IOobject::NO_WRITE
255             ),
256             mesh_,
257             dimensionedScalar("E", dimMass/dimLength/pow3(dimTime), 0.0)
258         )
259     );
261     if (mesh_.foundObject<volScalarField>("dQ"))
262     {
263         const volScalarField& dQ =
264             mesh_.lookupObject<volScalarField>("dQ");
265         E().internalField() = EhrrCoeff_*dQ;
266     }
268     return E;
272 // ************************************************************************* //