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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 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 "unitConversion.H"
32 #include "basicSpecieMixture.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38  namespace radiation
39  {
40  defineTypeNameAndDebug(greyMeanAbsorptionEmission, 0);
41 
43  (
44  absorptionEmissionModel,
45  greyMeanAbsorptionEmission,
46  dictionary
47  );
48  }
49 }
50 
51 
52 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53 
55 (
56  const dictionary& dict,
57  const fvMesh& mesh
58 )
59 :
61  coeffsDict_((dict.optionalSubDict(typeName + "Coeffs"))),
62  speciesNames_(0),
63  specieIndex_(Zero),
64  lookUpTablePtr_(),
65  thermo_(mesh.lookupObject<fluidThermo>(basicThermo::dictName)),
66  EhrrCoeff_(coeffsDict_.get<scalar>("EhrrCoeff")),
67  Yj_(nSpecies_)
68 {
69  if (!isA<basicSpecieMixture>(thermo_))
70  {
72  << "Model requires a multi-component thermo package"
73  << abort(FatalError);
74  }
75 
76 
77  label nFunc = 0;
78  const dictionary& functionDicts = dict.optionalSubDict(typeName + "Coeffs");
79 
80  for (const entry& dEntry : functionDicts)
81  {
82  if (!dEntry.isDict()) // safety
83  {
84  continue;
85  }
86 
87  const word& key = dEntry.keyword();
88  const dictionary& dict = dEntry.dict();
89 
90  speciesNames_.insert(key, nFunc);
91 
92  coeffs_[nFunc].initialise(dict);
93  nFunc++;
94  }
95 
96  if
97  (
98  coeffsDict_.found("lookUpTableFileName")
99  && "none" != coeffsDict_.get<word>("lookUpTableFileName")
100  )
101  {
102  lookUpTablePtr_.reset
103  (
105  (
106  coeffsDict_.get<fileName>("lookUpTableFileName"),
107  mesh.time().constant(),
108  mesh
109  )
110  );
111 
112  if (!mesh.foundObject<volScalarField>("ft"))
113  {
115  << "specie ft is not present to use with "
116  << "lookUpTableFileName " << nl
117  << exit(FatalError);
118  }
119  }
120 
121  // Check that all the species on the dictionary are present in the
122  // look-up table and save the corresponding indices of the look-up table
123 
124  label j = 0;
125  forAllConstIters(speciesNames_, iter)
126  {
127  const word& specieName = iter.key();
128  const label index = iter.val();
129 
130  volScalarField* fldPtr = mesh.getObjectPtr<volScalarField>(specieName);
131 
132  if (!lookUpTablePtr_.empty())
133  {
134  if (lookUpTablePtr_().found(specieName))
135  {
136  const label fieldIndex =
137  lookUpTablePtr_().findFieldIndex(specieName);
138 
139  Info<< "specie: " << specieName << " found on look-up table "
140  << " with index: " << fieldIndex << endl;
141 
142  specieIndex_[index] = fieldIndex;
143  }
144  else if (fldPtr)
145  {
146  Yj_.set(j, fldPtr);
147  specieIndex_[index] = 0;
148  j++;
149  Info<< "specie: " << iter.key() << " is being solved" << endl;
150  }
151  else
152  {
154  << "specie: " << specieName
155  << " is neither in look-up table: "
156  << lookUpTablePtr_().tableName()
157  << " nor is being solved" << nl
158  << exit(FatalError);
159  }
160  }
161  else if (fldPtr)
162  {
163  Yj_.set(j, fldPtr);
164  specieIndex_[index] = 0;
165  j++;
166  }
167  else
168  {
170  << "There is no lookup table and the specie" << nl
171  << specieName << nl
172  << " is not found " << nl
173  << exit(FatalError);
174  }
175  }
176 }
177 
178 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
179 
181 {}
182 
183 
184 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
185 
188 {
189  const basicSpecieMixture& mixture =
190  dynamic_cast<const basicSpecieMixture&>(thermo_);
191 
192  const volScalarField& T = thermo_.T();
193  const volScalarField& p = thermo_.p();
194 
195 
197  (
198  new volScalarField
199  (
200  IOobject
201  (
202  "aCont" + name(bandI),
203  mesh().time().timeName(),
204  mesh(),
207  ),
208  mesh(),
210  extrapolatedCalculatedFvPatchVectorField::typeName
211  )
212  );
213 
214  scalarField& a = ta.ref().primitiveFieldRef();
215 
216  forAll(a, celli)
217  {
218  forAllConstIters(speciesNames_, iter)
219  {
220  label n = iter();
221  scalar Xipi = 0.0;
222  if (specieIndex_[n] != 0)
223  {
224  //Specie found in the lookUpTable.
225  const volScalarField& ft =
226  mesh_.lookupObject<volScalarField>("ft");
227 
228  const List<scalar>& Ynft = lookUpTablePtr_().lookUp(ft[celli]);
229  //moles x pressure [atm]
230  Xipi = Ynft[specieIndex_[n]]*paToAtm(p[celli]);
231  }
232  else
233  {
234  scalar invWt = 0.0;
235  forAll(mixture.Y(), s)
236  {
237  invWt += mixture.Y(s)[celli]/mixture.W(s);
238  }
239 
240  label index = mixture.species()[iter.key()];
241  scalar Xk = mixture.Y(index)[celli]/(mixture.W(index)*invWt);
242 
243  Xipi = Xk*paToAtm(p[celli]);
244  }
245 
246  const absorptionCoeffs::coeffArray& b = coeffs_[n].coeffs(T[celli]);
247 
248  scalar Ti = T[celli];
249  // negative temperature exponents
250  if (coeffs_[n].invTemp())
251  {
252  Ti = 1.0/T[celli];
253  }
254  a[celli] +=
255  Xipi
256  *(
257  ((((b[5]*Ti + b[4])*Ti + b[3])*Ti + b[2])*Ti + b[1])*Ti
258  + b[0]
259  );
260  }
261  }
263  return ta;
264 }
265 
266 
269 {
270  return aCont(bandI);
271 }
272 
273 
276 {
278  (
279  new volScalarField
280  (
281  IOobject
282  (
283  "ECont" + name(bandI),
284  mesh_.time().timeName(),
285  mesh_,
288  ),
289  mesh_,
291  )
292  );
293 
294  const volScalarField* QdotPtr = mesh_.findObject<volScalarField>("Qdot");
295 
296  if (QdotPtr)
297  {
298  const volScalarField& Qdot = *QdotPtr;
299 
300  if (Qdot.dimensions() == dimEnergy/dimTime)
301  {
302  E.ref().primitiveFieldRef() = EhrrCoeff_*Qdot/mesh_.V();
303  }
304  else if (Qdot.dimensions() == dimEnergy/dimTime/dimVolume)
305  {
306  E.ref().primitiveFieldRef() = EhrrCoeff_*Qdot;
307  }
308  else
309  {
310  if (debug)
311  {
313  << "Incompatible dimensions for Qdot field" << endl;
314  }
315  }
316  }
317  else
318  {
320  << "Qdot field not found in mesh" << endl;
321  }
322 
323  return E;
324 }
325 
326 
327 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
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::radiation::greyMeanAbsorptionEmission::~greyMeanAbsorptionEmission
virtual ~greyMeanAbsorptionEmission()
Destructor.
Definition: greyMeanAbsorptionEmission.C:180
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::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
Foam::radiation::greyMeanAbsorptionEmission::ECont
tmp< volScalarField > ECont(const label bandI=0) const
Emission contribution for continuous phase.
Definition: greyMeanAbsorptionEmission.C:275
greyMeanAbsorptionEmission.H
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::radiation::greyMeanAbsorptionEmission::aCont
tmp< volScalarField > aCont(const label bandI=0) const
Absorption coefficient for continuous phase.
Definition: greyMeanAbsorptionEmission.C:187
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::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
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
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
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:137
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::GeometricField::correctBoundaryConditions
void correctBoundaryConditions()
Correct boundary field.
Definition: GeometricField.C:909
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
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::greyMeanAbsorptionEmission::eCont
tmp< volScalarField > eCont(const label bandI=0) const
Emission coefficient for continuous phase.
Definition: greyMeanAbsorptionEmission.C:268
Foam::radiation::greyMeanAbsorptionEmission::greyMeanAbsorptionEmission
greyMeanAbsorptionEmission(const dictionary &dict, const fvMesh &mesh)
Construct from components.
Definition: greyMeanAbsorptionEmission.C:55
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)
extrapolatedCalculatedFvPatchFields.H