wideBandAbsorptionEmission.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2018 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
30 #include "basicSpecieMixture.H"
31 #include "unitConversion.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37  namespace radiation
38  {
39  defineTypeNameAndDebug(wideBandAbsorptionEmission, 0);
40 
42  (
43  absorptionEmissionModel,
44  wideBandAbsorptionEmission,
45  dictionary
46  );
47  }
48 }
49 
50 
51 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52 
54 (
55  const dictionary& dict,
56  const fvMesh& mesh
57 )
58 :
60  coeffsDict_((dict.optionalSubDict(typeName + "Coeffs"))),
61  speciesNames_(0),
62  specieIndex_(Zero),
63  lookUpTablePtr_(),
64  thermo_(mesh.lookupObject<fluidThermo>(basicThermo::dictName)),
65  Yj_(nSpecies_),
66  totalWaveLength_(0)
67 {
68  label nBand = 0;
69  const dictionary& functionDicts = dict.optionalSubDict(typeName +"Coeffs");
70  for (const entry& dEntry : functionDicts)
71  {
72  if (!dEntry.isDict()) // safety
73  {
74  continue;
75  }
76 
77  const dictionary& dict = dEntry.dict();
78 
79  dict.readEntry("bandLimits", iBands_[nBand]);
80  dict.readEntry("EhrrCoeff", iEhrrCoeffs_[nBand]);
81  totalWaveLength_ += iBands_[nBand][1] - iBands_[nBand][0];
82 
83  label nSpec = 0;
84 
85  const dictionary& specDicts = dict.subDict("species");
86  for (const entry& dEntry : specDicts)
87  {
88  const word& key = dEntry.keyword();
89 
90  if (nBand == 0)
91  {
92  speciesNames_.insert(key, nSpec);
93  }
94  else if (!speciesNames_.found(key))
95  {
97  << "specie: " << key << " is not in all the bands"
98  << nl << exit(FatalError);
99  }
100  coeffs_[nBand][nSpec].initialise(specDicts.subDict(key));
101  nSpec++;
102  }
103  nBand++;
104  }
105  nBands_ = nBand;
106 
107  if
108  (
109  coeffsDict_.found("lookUpTableFileName")
110  && "none" != coeffsDict_.get<word>("lookUpTableFileName")
111  )
112  {
113  lookUpTablePtr_.set
114  (
116  (
117  coeffsDict_.get<fileName>("lookUpTableFileName"),
118  mesh.time().constant(),
119  mesh
120  )
121  );
122 
123  if (!mesh.foundObject<volScalarField>("ft"))
124  {
126  << "specie ft is not present to use with "
127  << "lookUpTableFileName " << nl
128  << exit(FatalError);
129  }
130  }
131 
132  // Check that all the species on the dictionary are present in the
133  // look-up table and save the corresponding indices of the look-up table
134 
135  label j = 0;
136  forAllConstIters(speciesNames_, iter)
137  {
138  const word& specieName = iter.key();
139  const label index = iter.val();
140 
141  volScalarField* fldPtr = mesh.getObjectPtr<volScalarField>(specieName);
142 
143  if (!lookUpTablePtr_.empty())
144  {
145  if (lookUpTablePtr_().found(specieName))
146  {
147  const label fieldIndex =
148  lookUpTablePtr_().findFieldIndex(specieName);
149 
150  Info<< "specie: " << specieName << " found on look-up table "
151  << " with index: " << fieldIndex << endl;
152 
153  specieIndex_[index] = fieldIndex;
154  }
155  else if (fldPtr)
156  {
157  Yj_.set(j, fldPtr);
158  specieIndex_[index] = 0;
159  j++;
160  Info<< "specie: " << specieName << " is being solved" << endl;
161  }
162  else
163  {
165  << "specie: " << specieName
166  << " is neither in look-up table: "
167  << lookUpTablePtr_().tableName()
168  << " nor is being solved" << nl
169  << exit(FatalError);
170  }
171  }
172  else if (fldPtr)
173  {
174  Yj_.set(j, fldPtr);
175  specieIndex_[index] = 0;
176  j++;
177  }
178  else
179  {
181  << "There is no lookup table and the specie" << nl
182  << specieName << nl
183  << " is not found " << nl
184  << exit(FatalError);
185 
186  }
187  }
188 }
189 
190 
191 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
192 
194 {}
195 
196 
197 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
198 
201 {
202  const basicSpecieMixture& mixture =
203  dynamic_cast<const basicSpecieMixture&>(thermo_);
204 
205  const volScalarField& T = thermo_.T();
206  const volScalarField& p = thermo_.p();
207 
209  (
210  new volScalarField
211  (
212  IOobject
213  (
214  "a",
215  mesh().time().timeName(),
216  mesh(),
219  ),
220  mesh(),
222  )
223  );
224 
225  scalarField& a = ta.ref().primitiveFieldRef();
226 
227  forAll(a, celli)
228  {
229  forAllConstIters(speciesNames_, iter)
230  {
231  const label n = iter();
232  scalar Xipi = 0;
233  if (specieIndex_[n] != 0)
234  {
235  const volScalarField& ft =
236  mesh_.lookupObject<volScalarField>("ft");
237 
238  const List<scalar>& Ynft = lookUpTablePtr_().lookUp(ft[celli]);
239 
240  // moles*pressure [atm]
241  Xipi = Ynft[specieIndex_[n]]*paToAtm(p[celli]);
242  }
243  else
244  {
245  scalar invWt = 0;
246  forAll(mixture.Y(), s)
247  {
248  invWt += mixture.Y(s)[celli]/mixture.W(s);
249  }
250 
251  const label index = mixture.species()[iter.key()];
252 
253  const scalar Xk =
254  mixture.Y(index)[celli]/(mixture.W(index)*invWt);
255 
256  Xipi = Xk*paToAtm(p[celli]);
257  }
258 
259  scalar Ti = T[celli];
260 
262  coeffs_[bandi][n].coeffs(T[celli]);
263 
264  if (coeffs_[bandi][n].invTemp())
265  {
266  Ti = 1.0/T[celli];
267  }
268 
269  a[celli]+=
270  Xipi
271  *(
272  ((((b[5]*Ti + b[4])*Ti + b[3])*Ti + b[2])*Ti + b[1])*Ti
273  + b[0]
274  );
275  }
276  }
277 
278  return ta;
279 }
280 
281 
284 {
285  return aCont(bandi);
286 }
287 
288 
291 {
293  (
294  new volScalarField
295  (
296  IOobject
297  (
298  "E",
299  mesh().time().timeName(),
300  mesh(),
303  ),
304  mesh(),
306  )
307  );
308 
309  const volScalarField* QdotPtr = mesh().findObject<volScalarField>("Qdot");
310 
311  if (QdotPtr)
312  {
313  const volScalarField& Qdot = *QdotPtr;
314 
315  if (Qdot.dimensions() == dimEnergy/dimTime)
316  {
317  E.ref().primitiveFieldRef() =
318  iEhrrCoeffs_[bandi]
319  *Qdot.primitiveField()
320  *(iBands_[bandi][1] - iBands_[bandi][0])
321  /totalWaveLength_
322  /mesh_.V();
323  }
324  else if (Qdot.dimensions() == dimEnergy/dimTime/dimVolume)
325  {
326  E.ref().primitiveFieldRef() =
327  iEhrrCoeffs_[bandi]
328  *Qdot.primitiveField()
329  *(iBands_[bandi][1] - iBands_[bandi][0])
330  /totalWaveLength_;
331  }
332  else
333  {
335  << "Incompatible dimensions for Qdot field" << endl;
336  }
337  }
338 
339  return E;
340 }
341 
342 
344 (
345  volScalarField& a,
346  PtrList<volScalarField>& aLambda
347 ) const
348 {
350 
351  for (label j=0; j<nBands_; j++)
352  {
353  aLambda[j].primitiveFieldRef() = this->a(j);
354 
355  a.primitiveFieldRef() +=
356  aLambda[j].primitiveField()
357  *(iBands_[j][1] - iBands_[j][0])
358  /totalWaveLength_;
359  }
360 
361 }
362 
363 
364 // ************************************************************************* //
Foam::entry
A keyword and a list of tokens is an 'entry'.
Definition: entry.H:67
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
Definition: dimensionSets.H:50
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::fileName
A class for handling file names.
Definition: fileName.H:69
s
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.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))
Definition: gmvOutputSpray.H:25
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::basicSpecieMixture
Specialization of basicMultiComponentMixture for a mixture consisting of a number for molecular speci...
Definition: basicSpecieMixture.H:54
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::dimEnergy
const dimensionSet dimEnergy
Foam::radiation::wideBandAbsorptionEmission::ECont
tmp< volScalarField > ECont(const label bandi=0) const
Emission contribution for continuous phase.
Definition: wideBandAbsorptionEmission.C:290
Foam::fluidThermo
Fundamental fluid thermodynamic properties.
Definition: fluidThermo.H:52
dictName
const word dictName("blockMeshDict")
unitConversion.H
Unit conversion functions.
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:54
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:258
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:62
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
Foam::Field< scalar >
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:89
Foam::paToAtm
constexpr scalar paToAtm(const scalar pa) noexcept
Conversion from Pa to atm.
Definition: unitConversion.H:111
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
Foam::radiation::wideBandAbsorptionEmission::~wideBandAbsorptionEmission
virtual ~wideBandAbsorptionEmission()
Destructor.
Definition: wideBandAbsorptionEmission.C:193
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:65
radiation
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
timeName
word timeName
Definition: getTimeIndex.H:3
Qdot
scalar Qdot
Definition: solveChemistry.H:2
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:84
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
wideBandAbsorptionEmission.H
Foam::radiation::defineTypeNameAndDebug
defineTypeNameAndDebug(cloudAbsorptionEmission, 0)
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::GeometricField::primitiveFieldRef
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
Definition: GeometricField.C:735
Foam::radiation::wideBandAbsorptionEmission::wideBandAbsorptionEmission
wideBandAbsorptionEmission(const dictionary &dict, const fvMesh &mesh)
Construct from components.
Definition: wideBandAbsorptionEmission.C:54
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::radiation::wideBandAbsorptionEmission::eCont
tmp< volScalarField > eCont(const label bandi=0) const
Emission coefficient for continuous phase.
Definition: wideBandAbsorptionEmission.C:283
Foam::nl
constexpr char nl
Definition: Ostream.H:372
forAllConstIters
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
Foam::List< scalar >
Foam::FixedList< scalar, nCoeffs_ >
Foam::radiation::wideBandAbsorptionEmission::aCont
tmp< volScalarField > aCont(const label bandi=0) const
Absorption coefficient for continuous phase.
Definition: wideBandAbsorptionEmission.C:200
Foam::radiation::absorptionEmissionModel
Model to supply absorption and emission coefficients for radiation modelling.
Definition: absorptionEmissionModel.H:54
Foam::interpolationLookUpTable< scalar >
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:61
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:123
basicSpecieMixture.H
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:294
Foam::mixture
Definition: mixture.H:54
Foam::radiation::addToRunTimeSelectionTable
addToRunTimeSelectionTable(absorptionEmissionModel, cloudAbsorptionEmission, dictionary)
Foam::radiation::wideBandAbsorptionEmission::correct
void correct(volScalarField &a, PtrList< volScalarField > &aLambda) const
Correct rays.
Definition: wideBandAbsorptionEmission.C:344