1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2008-2011 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
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
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 * * * * * * * * * * * * * //
38 defineTypeNameAndDebug(greyMeanAbsorptionEmission, 0);
40 addToRunTimeSelectionTable
42 absorptionEmissionModel,
43 greyMeanAbsorptionEmission,
50 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52 Foam::radiation::greyMeanAbsorptionEmission::greyMeanAbsorptionEmission
54 const dictionary& dict,
58 absorptionEmissionModel(dict, mesh),
59 coeffsDict_((dict.subDict(typeName + "Coeffs"))),
63 thermo_(mesh.lookupObject<basicThermo>("thermophysicalProperties")),
64 EhrrCoeff_(readScalar(coeffsDict_.lookup("EhrrCoeff"))),
67 if (!isA<basicMultiComponentMixture>(thermo_))
71 "radiation::greyMeanAbsorptionEmission::greyMeanAbsorptionEmission"
76 ) << "Model requires a multi-component thermo package"
82 const dictionary& functionDicts = dict.subDict(typeName + "Coeffs");
84 forAllConstIter(dictionary, functionDicts, iter)
91 const word& key = iter().keyword();
92 speciesNames_.insert(key, nFunc);
93 const dictionary& dict = iter().dict();
94 coeffs_[nFunc].initialise(dict);
98 if (coeffsDict_.found("lookUpTableFileName"))
100 const word name = coeffsDict_.lookup("lookUpTableFileName");
105 new interpolationLookUpTable<scalar>
107 fileName(coeffsDict_.lookup("lookUpTableFileName")),
108 mesh.time().constant(),
113 if (!mesh.foundObject<volScalarField>("ft"))
117 "Foam::radiation::greyMeanAbsorptionEmission(const"
118 "dictionary& dict, const fvMesh& mesh)"
119 ) << "specie ft is not present to use with "
120 << "lookUpTableFileName " << nl
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
130 forAllConstIter(HashTable<label>, speciesNames_, iter)
132 if (!lookUpTablePtr_.empty())
134 if (lookUpTablePtr_().found(iter.key()))
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;
143 else if (mesh.foundObject<volScalarField>(iter.key()))
146 const_cast<volScalarField&>
148 mesh.lookupObject<volScalarField>(iter.key())
151 specieIndex_[iter()] = 0;
153 Info<< "specie: " << iter.key() << " is being solved" << endl;
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
168 else if (mesh.foundObject<volScalarField>(iter.key()))
171 const_cast<volScalarField&>
173 mesh.lookupObject<volScalarField>(iter.key())
177 specieIndex_[iter()] = 0;
184 "Foam::radiation::greyMeanAbsorptionEmission(const"
185 "dictionary& dict, const fvMesh& mesh)"
186 ) << " there is not lookup table and the specie" << nl
188 << " is not found " << nl
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
219 "aCont" + name(bandI),
220 mesh().time().timeName(),
226 dimensionedScalar("a", dimless/dimLength, 0.0),
227 zeroGradientFvPatchVectorField::typeName
231 scalarField& a = ta().internalField();
235 forAllConstIter(HashTable<label>, speciesNames_, iter)
239 if (specieIndex_[n] != 0)
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]);
252 forAll (mixture.Y(), s)
254 invWt += mixture.Y(s)[cellI]/mixture.W(s);
257 label index = mixture.species()[iter.key()];
258 scalar Xk = mixture.Y(index)[cellI]/(mixture.W(index)*invWt);
260 Xipi = Xk*paToAtm(p[cellI]);
263 const absorptionCoeffs::coeffArray& b = coeffs_[n].coeffs(T[cellI]);
265 scalar Ti = T[cellI];
266 // negative temperature exponents
267 if (coeffs_[n].invTemp())
274 ((((b[5]*Ti + b[4])*Ti + b[3])*Ti + b[2])*Ti + b[1])*Ti
279 ta().correctBoundaryConditions();
284 Foam::tmp<Foam::volScalarField>
285 Foam::radiation::greyMeanAbsorptionEmission::eCont(const label bandI) const
287 tmp<volScalarField> e
293 "eCont" + name(bandI),
294 mesh().time().timeName(),
300 dimensionedScalar("e", dimless/dimLength, 0.0)
308 Foam::tmp<Foam::volScalarField>
309 Foam::radiation::greyMeanAbsorptionEmission::ECont(const label bandI) const
311 tmp<volScalarField> E
317 "ECont" + name(bandI),
318 mesh_.time().timeName(),
324 dimensionedScalar("E", dimMass/dimLength/pow3(dimTime), 0.0)
328 if (mesh_.foundObject<volScalarField>("dQ"))
330 const volScalarField& dQ =
331 mesh_.lookupObject<volScalarField>("dQ");
332 E().internalField() = EhrrCoeff_*dQ;
339 // ************************************************************************* //