greyMeanSolidAbsorptionEmission.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-2017 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 
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  namespace radiation
37  {
38  defineTypeNameAndDebug(greyMeanSolidAbsorptionEmission, 0);
39 
41  (
42  absorptionEmissionModel,
43  greyMeanSolidAbsorptionEmission,
44  dictionary
45  );
46  }
47 }
48 
49 // * * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * //
50 
51 Foam::tmp<Foam::scalarField> Foam::radiation::
52 greyMeanSolidAbsorptionEmission::X(const word specie) const
53 {
54  const volScalarField& T = thermo_.T();
55  const volScalarField& p = thermo_.p();
56 
57  tmp<scalarField> tXj(new scalarField(T.primitiveField().size(), 0.0));
58  scalarField& Xj = tXj.ref();
59 
60  tmp<scalarField> tRhoInv(new scalarField(T.primitiveField().size(), 0.0));
61  scalarField& rhoInv = tRhoInv.ref();
62 
63  forAll(mixture_.Y(), specieI)
64  {
65  const scalarField& Yi = mixture_.Y()[specieI];
66 
67  forAll(rhoInv, iCell)
68  {
69  rhoInv[iCell] +=
70  Yi[iCell]/mixture_.rho(specieI, p[iCell], T[iCell]);
71  }
72  }
73  const scalarField& Yj = mixture_.Y(specie);
74  const label mySpecieI = mixture_.species()[specie];
75  forAll(Xj, iCell)
76  {
77  Xj[iCell] = Yj[iCell]/mixture_.rho(mySpecieI, p[iCell], T[iCell]);
78  }
79 
80  return (Xj/rhoInv);
81 }
82 
83 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
84 
87 (
88  const dictionary& dict,
89  const fvMesh& mesh
90 )
91 :
92  absorptionEmissionModel(dict, mesh),
93  coeffsDict_((dict.optionalSubDict(typeName + "Coeffs"))),
95  speciesNames_(0),
96  mixture_(dynamic_cast<const basicSpecieMixture&>(thermo_)),
97  solidData_(mixture_.Y().size())
98 {
99  if (!isA<basicSpecieMixture>(thermo_))
100  {
102  << "Model requires a multi-component thermo package"
103  << abort(FatalError);
104  }
105 
106  label nFunc = 0;
107  const dictionary& functionDicts = dict.optionalSubDict(typeName + "Coeffs");
108 
109  forAllConstIter(dictionary, functionDicts, iter)
110  {
111  // safety:
112  if (!iter().isDict())
113  {
114  continue;
115  }
116  const word& key = iter().keyword();
117  if (!mixture_.contains(key))
118  {
120  << " specie: " << key << " is not found in the solid mixture"
121  << nl
122  << " specie is the mixture are:" << mixture_.species() << nl
123  << nl << endl;
124  }
125  speciesNames_.insert(key, nFunc);
126  const dictionary& dict = iter().dict();
127  dict.lookup("absorptivity") >> solidData_[nFunc][absorptivity];
128  dict.lookup("emissivity") >> solidData_[nFunc][emissivity];
129 
130  nFunc++;
131  }
132 }
133 
134 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
135 
138 {}
139 
140 
141 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
142 
144 Foam::radiation::greyMeanSolidAbsorptionEmission::
145 calc(const label propertyId) const
146 {
148  (
149  new volScalarField
150  (
151  IOobject
152  (
153  "a",
154  mesh().time().timeName(),
155  mesh(),
158  ),
159  mesh(),
161  extrapolatedCalculatedFvPatchVectorField::typeName
162  )
163  );
164 
165  scalarField& a = ta.ref().primitiveFieldRef();
166 
167  forAllConstIter(HashTable<label>, speciesNames_, iter)
168  {
169  if (mixture_.contains(iter.key()))
170  {
171  a += solidData_[iter()][propertyId]*X(iter.key());
172  }
173  }
174 
175  ta.ref().correctBoundaryConditions();
176  return ta;
177 }
178 
179 
182 (
183  const label bandI
184 ) const
185 {
186  return calc(emissivity);
187 }
188 
189 
192 (
193  const label bandI
194 ) const
195 {
196  return calc(absorptivity);
197 }
198 
199 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
greyMeanSolidAbsorptionEmission(const dictionary &dict, const fvMesh &mesh)
Construct from components.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
Unit conversion functions.
defineTypeNameAndDebug(cloudAbsorptionEmission, 0)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
Model to supply absorption and emission coefficients for radiation modelling.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Macros for easy insertion into run-time selection tables.
const word dictName() const
Return the local dictionary name (final part of scoped name)
Definition: dictionary.H:115
tmp< volScalarField > eCont(const label bandI=0) const
Emission coefficient for continuous phase.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
dynamicFvMesh & mesh
const dictionary & optionalSubDict(const word &) const
Find and return a sub-dictionary if found.
Definition: dictionary.C:759
A class for handling words, derived from string.
Definition: word.H:59
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
word timeName
Definition: getTimeIndex.H:3
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
errorManip< error > abort(error &err)
Definition: errorManip.H:131
tmp< GeometricField< Type, PatchField, GeoMesh > > T() const
Return transpose (only if it is a tensor field)
static const char nl
Definition: Ostream.H:262
Fundamental solid thermodynamic properties.
Definition: solidThermo.H:54
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
#define WarningInFunction
Report a warning using Foam::Warning.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
addToRunTimeSelectionTable(absorptionEmissionModel, cloudAbsorptionEmission, dictionary)
tmp< volScalarField > aCont(const label bandI=0) const
Absorption coefficient for continuous phase.
A class for managing temporary objects.
Definition: PtrList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:576