ENH: autoLayerDriver: better layering information message
[OpenFOAM-2.0.x.git] / src / thermophysicalModels / radiationModels / submodels / absorptionEmissionModel / greyMeanAbsorptionEmission / greyMeanAbsorptionEmission.C
blob74f682d6458d2ca63e6e62af0bd9989fccd22ba6
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
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 "greyMeanAbsorptionEmission.H"
27 #include "addToRunTimeSelectionTable.H"
28 #include "unitConversion.H"
29 #include "zeroGradientFvPatchFields.H"
30 #include "basicMultiComponentMixture.H"
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 namespace Foam
36     namespace radiation
37     {
38         defineTypeNameAndDebug(greyMeanAbsorptionEmission, 0);
40         addToRunTimeSelectionTable
41         (
42             absorptionEmissionModel,
43             greyMeanAbsorptionEmission,
44             dictionary
45         );
46     }
50 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
52 Foam::radiation::greyMeanAbsorptionEmission::greyMeanAbsorptionEmission
54     const dictionary& dict,
55     const fvMesh& mesh
58     absorptionEmissionModel(dict, mesh),
59     coeffsDict_((dict.subDict(typeName + "Coeffs"))),
60     speciesNames_(0),
61     specieIndex_(0),
62     lookUpTablePtr_(),
63     thermo_(mesh.lookupObject<basicThermo>("thermophysicalProperties")),
64     EhrrCoeff_(readScalar(coeffsDict_.lookup("EhrrCoeff"))),
65     Yj_(nSpecies_)
67     if (!isA<basicMultiComponentMixture>(thermo_))
68     {
69         FatalErrorIn
70         (
71             "radiation::greyMeanAbsorptionEmission::greyMeanAbsorptionEmission"
72             "("
73                 "const dictionary&, "
74                 "const fvMesh&"
75             ")"
76         )   << "Model requires a multi-component thermo package"
77             << abort(FatalError);
78     }
81     label nFunc = 0;
82     const dictionary& functionDicts = dict.subDict(typeName + "Coeffs");
84     forAllConstIter(dictionary, functionDicts, iter)
85     {
86         // safety:
87         if (!iter().isDict())
88         {
89             continue;
90         }
91         const word& key = iter().keyword();
92         speciesNames_.insert(key, nFunc);
93         const dictionary& dict = iter().dict();
94         coeffs_[nFunc].initialise(dict);
95         nFunc++;
96     }
98     if (coeffsDict_.found("lookUpTableFileName"))
99     {
100         const word name = coeffsDict_.lookup("lookUpTableFileName");
101         if (name != "none")
102         {
103             lookUpTablePtr_.set
104             (
105                 new interpolationLookUpTable<scalar>
106                 (
107                     fileName(coeffsDict_.lookup("lookUpTableFileName")),
108                     mesh.time().constant(),
109                     mesh
110                 )
111             );
113             if (!mesh.foundObject<volScalarField>("ft"))
114             {
115                 FatalErrorIn
116                 (
117                     "Foam::radiation::greyMeanAbsorptionEmission(const"
118                     "dictionary& dict, const fvMesh& mesh)"
119                 )   << "specie ft is not present to use with "
120                     << "lookUpTableFileName " << nl
121                     << exit(FatalError);
122             }
123         }
124     }
126     // Check that all the species on the dictionary are present in the
127     // look-up table and save the corresponding indices of the look-up table
129     label j = 0;
130     forAllConstIter(HashTable<label>, speciesNames_, iter)
131     {
132         if (!lookUpTablePtr_.empty())
133         {
134             if (lookUpTablePtr_().found(iter.key()))
135             {
136                 label index = lookUpTablePtr_().findFieldIndex(iter.key());
138                 Info<< "specie: " << iter.key() << " found on look-up table "
139                     << " with index: " << index << endl;
141                 specieIndex_[iter()] = index;
142             }
143             else if (mesh.foundObject<volScalarField>(iter.key()))
144             {
145                 volScalarField& Y =
146                     const_cast<volScalarField&>
147                     (
148                         mesh.lookupObject<volScalarField>(iter.key())
149                     );
150                 Yj_.set(j, &Y);
151                 specieIndex_[iter()] = 0;
152                 j++;
153                 Info<< "specie: " << iter.key() << " is being solved" << endl;
154             }
155             else
156             {
157                 FatalErrorIn
158                 (
159                     "Foam::radiation::greyMeanAbsorptionEmission(const"
160                     "dictionary& dict, const fvMesh& mesh)"
161                 )   << "specie: " << iter.key()
162                     << " is neither in look-up table: "
163                     << lookUpTablePtr_().tableName()
164                     << " nor is being solved" << nl
165                     << exit(FatalError);
166             }
167         }
168         else if (mesh.foundObject<volScalarField>(iter.key()))
169         {
170             volScalarField& Y =
171                 const_cast<volScalarField&>
172                 (
173                     mesh.lookupObject<volScalarField>(iter.key())
174                 );
176             Yj_.set(j, &Y);
177             specieIndex_[iter()] = 0;
178             j++;
179         }
180         else
181         {
182             FatalErrorIn
183             (
184                 "Foam::radiation::greyMeanAbsorptionEmission(const"
185                 "dictionary& dict, const fvMesh& mesh)"
186             )   << " there is not lookup table and the specie" << nl
187                 << iter.key() << nl
188                 << " is not found " << nl
189                 << exit(FatalError);
191         }
192     }
195 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
197 Foam::radiation::greyMeanAbsorptionEmission::~greyMeanAbsorptionEmission()
201 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
203 Foam::tmp<Foam::volScalarField>
204 Foam::radiation::greyMeanAbsorptionEmission::aCont(const label bandI) const
206     const basicMultiComponentMixture& mixture =
207         dynamic_cast<const basicMultiComponentMixture&>(thermo_);
209     const volScalarField& T = thermo_.T();
210     const volScalarField& p = thermo_.p();
213     tmp<volScalarField> ta
214     (
215         new volScalarField
216         (
217             IOobject
218             (
219                 "aCont" + name(bandI),
220                 mesh().time().timeName(),
221                 mesh(),
222                 IOobject::NO_READ,
223                 IOobject::NO_WRITE
224             ),
225             mesh(),
226             dimensionedScalar("a", dimless/dimLength, 0.0),
227             zeroGradientFvPatchVectorField::typeName
228         )
229     );
231     scalarField& a = ta().internalField();
233     forAll(a, cellI)
234     {
235         forAllConstIter(HashTable<label>, speciesNames_, iter)
236         {
237             label n = iter();
238             scalar Xipi = 0.0;
239             if (specieIndex_[n] != 0)
240             {
241                 //Specie found in the lookUpTable.
242                 const volScalarField& ft =
243                     mesh_.lookupObject<volScalarField>("ft");
245                 const List<scalar>& Ynft = lookUpTablePtr_().lookUp(ft[cellI]);
246                 //moles x pressure [atm]
247                 Xipi = Ynft[specieIndex_[n]]*paToAtm(p[cellI]);
248             }
249             else
250             {
251                 scalar invWt = 0.0;
252                 forAll (mixture.Y(), s)
253                 {
254                     invWt += mixture.Y(s)[cellI]/mixture.W(s);
255                 }
257                 label index = mixture.species()[iter.key()];
258                 scalar Xk = mixture.Y(index)[cellI]/(mixture.W(index)*invWt);
260                 Xipi = Xk*paToAtm(p[cellI]);
261             }
263             const absorptionCoeffs::coeffArray& b = coeffs_[n].coeffs(T[cellI]);
265             scalar Ti = T[cellI];
266             // negative temperature exponents
267             if (coeffs_[n].invTemp())
268             {
269                 Ti = 1.0/T[cellI];
270             }
271             a[cellI] +=
272                 Xipi
273                *(
274                     ((((b[5]*Ti + b[4])*Ti + b[3])*Ti + b[2])*Ti + b[1])*Ti
275                   + b[0]
276                 );
277         }
278     }
279     ta().correctBoundaryConditions();
280     return ta;
284 Foam::tmp<Foam::volScalarField>
285 Foam::radiation::greyMeanAbsorptionEmission::eCont(const label bandI) const
287     tmp<volScalarField> e
288     (
289         new volScalarField
290         (
291             IOobject
292             (
293                 "eCont" + name(bandI),
294                 mesh().time().timeName(),
295                 mesh(),
296                 IOobject::NO_READ,
297                 IOobject::NO_WRITE
298             ),
299             mesh(),
300             dimensionedScalar("e", dimless/dimLength, 0.0)
301         )
302     );
304     return e;
308 Foam::tmp<Foam::volScalarField>
309 Foam::radiation::greyMeanAbsorptionEmission::ECont(const label bandI) const
311     tmp<volScalarField> E
312     (
313         new volScalarField
314         (
315             IOobject
316             (
317                 "ECont" + name(bandI),
318                 mesh_.time().timeName(),
319                 mesh_,
320                 IOobject::NO_READ,
321                 IOobject::NO_WRITE
322             ),
323             mesh_,
324             dimensionedScalar("E", dimMass/dimLength/pow3(dimTime), 0.0)
325         )
326     );
328     if (mesh_.foundObject<volScalarField>("dQ"))
329     {
330         const volScalarField& dQ =
331             mesh_.lookupObject<volScalarField>("dQ");
332         E().internalField() = EhrrCoeff_*dQ;
333     }
335     return E;
339 // ************************************************************************* //