greyMeanAbsorptionEmission.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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.
15 
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.
20 
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/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
28 #include "unitConversion.H"
30 #include "basicSpecieMixture.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  namespace radiation
37  {
38  defineTypeNameAndDebug(greyMeanAbsorptionEmission, 0);
39 
41  (
42  absorptionEmissionModel,
43  greyMeanAbsorptionEmission,
44  dictionary
45  );
46  }
47 }
48 
49 
50 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51 
53 (
54  const dictionary& dict,
55  const fvMesh& mesh
56 )
57 :
59  coeffsDict_((dict.subDict(typeName + "Coeffs"))),
60  speciesNames_(0),
61  specieIndex_(label(0)),
62  lookUpTablePtr_(),
63  thermo_(mesh.lookupObject<fluidThermo>(basicThermo::dictName)),
64  EhrrCoeff_(readScalar(coeffsDict_.lookup("EhrrCoeff"))),
65  Yj_(nSpecies_)
66 {
67  if (!isA<basicSpecieMixture>(thermo_))
68  {
70  << "Model requires a multi-component thermo package"
71  << abort(FatalError);
72  }
73 
74 
75  label nFunc = 0;
76  const dictionary& functionDicts = dict.subDict(typeName + "Coeffs");
77 
78  forAllConstIter(dictionary, functionDicts, iter)
79  {
80  // safety:
81  if (!iter().isDict())
82  {
83  continue;
84  }
85  const word& key = iter().keyword();
86  speciesNames_.insert(key, nFunc);
87  const dictionary& dict = iter().dict();
88  coeffs_[nFunc].initialise(dict);
89  nFunc++;
90  }
91 
92  if (coeffsDict_.found("lookUpTableFileName"))
93  {
94  const word name = coeffsDict_.lookup("lookUpTableFileName");
95  if (name != "none")
96  {
97  lookUpTablePtr_.set
98  (
100  (
101  fileName(coeffsDict_.lookup("lookUpTableFileName")),
102  mesh.time().constant(),
103  mesh
104  )
105  );
106 
107  if (!mesh.foundObject<volScalarField>("ft"))
108  {
110  << "specie ft is not present to use with "
111  << "lookUpTableFileName " << nl
112  << exit(FatalError);
113  }
114  }
115  }
116 
117  // Check that all the species on the dictionary are present in the
118  // look-up table and save the corresponding indices of the look-up table
119 
120  label j = 0;
121  forAllConstIter(HashTable<label>, speciesNames_, iter)
122  {
123  if (!lookUpTablePtr_.empty())
124  {
125  if (lookUpTablePtr_().found(iter.key()))
126  {
127  label index = lookUpTablePtr_().findFieldIndex(iter.key());
128 
129  Info<< "specie: " << iter.key() << " found on look-up table "
130  << " with index: " << index << endl;
131 
132  specieIndex_[iter()] = index;
133  }
134  else if (mesh.foundObject<volScalarField>(iter.key()))
135  {
136  volScalarField& Y =
137  const_cast<volScalarField&>
138  (
139  mesh.lookupObject<volScalarField>(iter.key())
140  );
141  Yj_.set(j, &Y);
142  specieIndex_[iter()] = 0;
143  j++;
144  Info<< "specie: " << iter.key() << " is being solved" << endl;
145  }
146  else
147  {
149  << "specie: " << iter.key()
150  << " is neither in look-up table: "
151  << lookUpTablePtr_().tableName()
152  << " nor is being solved" << nl
153  << exit(FatalError);
154  }
155  }
156  else if (mesh.foundObject<volScalarField>(iter.key()))
157  {
158  volScalarField& Y =
159  const_cast<volScalarField&>
160  (
161  mesh.lookupObject<volScalarField>(iter.key())
162  );
163 
164  Yj_.set(j, &Y);
165  specieIndex_[iter()] = 0;
166  j++;
167  }
168  else
169  {
171  << " there is not lookup table and the specie" << nl
172  << iter.key() << nl
173  << " is not found " << nl
174  << exit(FatalError);
175  }
176  }
177 }
178 
179 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
180 
182 {}
183 
184 
185 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
186 
189 {
190  const basicSpecieMixture& mixture =
191  dynamic_cast<const basicSpecieMixture&>(thermo_);
192 
193  const volScalarField& T = thermo_.T();
194  const volScalarField& p = thermo_.p();
195 
196 
198  (
199  new volScalarField
200  (
201  IOobject
202  (
203  "aCont" + name(bandI),
204  mesh().time().timeName(),
205  mesh(),
208  ),
209  mesh(),
211  zeroGradientFvPatchVectorField::typeName
212  )
213  );
214 
215  scalarField& a = ta().internalField();
216 
217  forAll(a, cellI)
218  {
219  forAllConstIter(HashTable<label>, speciesNames_, iter)
220  {
221  label n = iter();
222  scalar Xipi = 0.0;
223  if (specieIndex_[n] != 0)
224  {
225  //Specie found in the lookUpTable.
226  const volScalarField& ft =
227  mesh_.lookupObject<volScalarField>("ft");
228 
229  const List<scalar>& Ynft = lookUpTablePtr_().lookUp(ft[cellI]);
230  //moles x pressure [atm]
231  Xipi = Ynft[specieIndex_[n]]*paToAtm(p[cellI]);
232  }
233  else
234  {
235  scalar invWt = 0.0;
236  forAll (mixture.Y(), s)
237  {
238  invWt += mixture.Y(s)[cellI]/mixture.W(s);
239  }
240 
241  label index = mixture.species()[iter.key()];
242  scalar Xk = mixture.Y(index)[cellI]/(mixture.W(index)*invWt);
243 
244  Xipi = Xk*paToAtm(p[cellI]);
245  }
246 
247  const absorptionCoeffs::coeffArray& b = coeffs_[n].coeffs(T[cellI]);
248 
249  scalar Ti = T[cellI];
250  // negative temperature exponents
251  if (coeffs_[n].invTemp())
252  {
253  Ti = 1.0/T[cellI];
254  }
255  a[cellI] +=
256  Xipi
257  *(
258  ((((b[5]*Ti + b[4])*Ti + b[3])*Ti + b[2])*Ti + b[1])*Ti
259  + b[0]
260  );
261  }
262  }
263  ta().correctBoundaryConditions();
264  return ta;
265 }
266 
267 
270 {
271  return aCont(bandI);
272 }
273 
274 
277 {
279  (
280  new volScalarField
281  (
282  IOobject
283  (
284  "ECont" + name(bandI),
285  mesh_.time().timeName(),
286  mesh_,
289  ),
290  mesh_,
292  )
293  );
294 
295  if (mesh_.foundObject<volScalarField>("dQ"))
296  {
297  const volScalarField& dQ =
298  mesh_.lookupObject<volScalarField>("dQ");
299 
300  if (dQ.dimensions() == dimEnergy/dimTime)
301  {
302  E().internalField() = EhrrCoeff_*dQ/mesh_.V();
303  }
304  else if (dQ.dimensions() == dimEnergy/dimTime/dimVolume)
305  {
306  E().internalField() = EhrrCoeff_*dQ;
307  }
308  else
309  {
310  if (debug)
311  {
313  << "Incompatible dimensions for dQ field" << endl;
314  }
315  }
316  }
317  else
318  {
320  << "dQ field not found in mesh" << endl;
321  }
322 
323  return E;
324 }
325 
326 
327 // ************************************************************************* //
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
Foam::radiation::greyMeanAbsorptionEmission::~greyMeanAbsorptionEmission
virtual ~greyMeanAbsorptionEmission()
Destructor.
Definition: greyMeanAbsorptionEmission.C:181
p
p
Definition: pEqn.H:62
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::fileName
A class for handling file names.
Definition: fileName.H:69
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
Foam::basicSpecieMixture
Specialization of basicMultiComponentMixture for a mixture consisting of a number for molecular speci...
Definition: basicSpecieMixture.H:49
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::dimEnergy
const dimensionSet dimEnergy
dQ
dQ
Definition: YEEqn.H:14
Foam::fluidThermo
Fundamental fluid thermodynamic properties.
Definition: fluidThermo.H:49
unitConversion.H
Unit conversion functions.
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::paToAtm
scalar paToAtm(const scalar pa)
Conversion from atm to Pa.
Definition: unitConversion.H:63
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
Foam::radiation::greyMeanAbsorptionEmission::ECont
tmp< volScalarField > ECont(const label bandI=0) const
Emission contribution for continuous phase.
Definition: greyMeanAbsorptionEmission.C:276
greyMeanAbsorptionEmission.H
dictName
const word dictName("particleTrackDict")
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:28
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:87
Foam::radiation::greyMeanAbsorptionEmission::aCont
tmp< volScalarField > aCont(const label bandI=0) const
Absorption coefficient for continuous phase.
Definition: greyMeanAbsorptionEmission.C:188
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:55
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
radiation
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
s
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Foam::radiation::defineTypeNameAndDebug
defineTypeNameAndDebug(cloudAbsorptionEmission, 0)
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
readScalar
#define readScalar
Definition: doubleScalar.C:38
Foam::HashTable< label >
found
bool found
Definition: TABSMDCalcMethod2.H:32
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::List< scalar >
Foam::FixedList< scalar, nCoeffs_ >
Foam::radiation::greyMeanAbsorptionEmission::eCont
tmp< volScalarField > eCont(const label bandI=0) const
Emission coefficient for continuous phase.
Definition: greyMeanAbsorptionEmission.C:269
timeName
word timeName
Definition: getTimeIndex.H:3
Foam::radiation::greyMeanAbsorptionEmission::greyMeanAbsorptionEmission
greyMeanAbsorptionEmission(const dictionary &dict, const fvMesh &mesh)
Construct from components.
Definition: greyMeanAbsorptionEmission.C:53
Foam::radiation::absorptionEmissionModel
Model to supply absorption and emission coefficients for radiation modelling.
Definition: absorptionEmissionModel.H:52
Foam::interpolationLookUpTable< scalar >
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
zeroGradientFvPatchFields.H
basicSpecieMixture.H
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:259
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Foam::mixture
Definition: mixture.H:52
Foam::radiation::addToRunTimeSelectionTable
addToRunTimeSelectionTable(absorptionEmissionModel, cloudAbsorptionEmission, dictionary)
Y
PtrList< volScalarField > & Y
Definition: createFields.H:36