1 /*---------------------------------------------------------------------------*\
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 -------------------------------------------------------------------------------
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 * * * * * * * * * * * * * //
35 defineTypeNameAndDebug(greyMeanAbsorptionEmission, 0);
37 addToRunTimeSelectionTable
39 absorptionEmissionModel,
40 greyMeanAbsorptionEmission,
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49 Foam::radiation::greyMeanAbsorptionEmission::greyMeanAbsorptionEmission
51 const dictionary& dict,
55 absorptionEmissionModel(dict, mesh),
56 coeffsDict_((dict.subDict(typeName + "Coeffs"))),
61 fileName(coeffsDict_.lookup("lookUpTableFileName")),
62 mesh.time().constant(),
65 thermo_(mesh.lookupObject<basicThermo>("thermophysicalProperties")),
66 EhrrCoeff_(readScalar(coeffsDict_.lookup("EhrrCoeff"))),
70 const dictionary& functionDicts = dict.subDict(typeName + "Coeffs");
72 forAllConstIter(dictionary, functionDicts, iter)
79 const word& key = iter().keyword();
80 speciesNames_.insert(key, nFunc);
81 const dictionary& dict = iter().dict();
82 coeffs_[nFunc].initialise(dict);
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
90 forAllConstIter(HashTable<label>, speciesNames_, iter)
92 if (mesh.foundObject<volScalarField>("ft"))
94 if (lookUpTable_.found(iter.key()))
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;
103 else if (mesh.foundObject<volScalarField>(iter.key()))
106 const_cast<volScalarField&>
108 mesh.lookupObject<volScalarField>(iter.key())
111 specieIndex_[iter()] = 0;
113 Info<< "specie: " << iter.key() << " is being solved" << endl;
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
132 "Foam::radiation::greyMeanAbsorptionEmission(const"
133 "dictionary& dict, const fvMesh& mesh)"
134 ) << "specie ft is not present " << nl
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
165 mesh().time().timeName(),
166 mesh(), // HR 101116 -- Todo: Need to bring object registry here!
171 dimensionedScalar("a", dimless/dimLength, 0.0)
175 scalarField& a = ta().internalField();
179 const List<scalar>& species = lookUpTable_.lookUp(ft[i]);
181 for (label n=0; n<nSpecies; n++)
185 if (specieIndex_[n] != 0)
187 //moles x pressure [atm]
188 Yipi = species[specieIndex_[n]]*p[i]*9.869231e-6;
197 const absorptionCoeffs::coeffArray& b = coeffs_[n].coeffs(T[i]);
200 // negative temperature exponents
201 if (coeffs_[n].invTemp())
208 ((((b[5]*Ti + b[4])*Ti + b[3])*Ti + b[2])*Ti + b[1])*Ti
217 Foam::tmp<Foam::volScalarField>
218 Foam::radiation::greyMeanAbsorptionEmission::eCont(const label bandI) const
220 tmp<volScalarField> e
227 mesh().time().timeName(),
228 mesh(), // HR 101116 -- Todo: Need to bring object registry here!
233 dimensionedScalar("e", dimless/dimLength, 0.0)
241 Foam::tmp<Foam::volScalarField>
242 Foam::radiation::greyMeanAbsorptionEmission::ECont(const label bandI) const
244 tmp<volScalarField> E
251 mesh_.time().timeName(),
252 mesh_, // HR 101116 -- Todo: Need to bring object registry here!
257 dimensionedScalar("E", dimMass/dimLength/pow3(dimTime), 0.0)
261 if (mesh_.foundObject<volScalarField>("dQ"))
263 const volScalarField& dQ =
264 mesh_.lookupObject<volScalarField>("dQ");
265 E().internalField() = EhrrCoeff_*dQ;
272 // ************************************************************************* //