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) 2015 OpenCFD Ltd
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"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34  namespace radiation
35  {
36  defineTypeNameAndDebug(greyMeanSolidAbsorptionEmission, 0);
37 
39  (
40  absorptionEmissionModel,
41  greyMeanSolidAbsorptionEmission,
42  dictionary
43  );
44  }
45 }
46 
47 // * * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * //
48 
51 {
52  const volScalarField& T = thermo_.T();
53  const volScalarField& p = thermo_.p();
54 
55  tmp<scalarField> tXj(new scalarField(T.internalField().size(), 0.0));
56  scalarField& Xj = tXj();
57 
58  tmp<scalarField> tRhoInv(new scalarField(T.internalField().size(), 0.0));
59  scalarField& rhoInv = tRhoInv();
60 
61  forAll(mixture_.Y(), specieI)
62  {
63  const scalarField& Yi = mixture_.Y()[specieI];
64 
65  forAll (rhoInv, iCell)
66  {
67  rhoInv[iCell] +=
68  Yi[iCell]/mixture_.rho(specieI, p[iCell], T[iCell]);
69  }
70  }
71  const scalarField& Yj = mixture_.Y(specie);
72  const label mySpecieI = mixture_.species()[specie];
73  forAll(Xj, iCell)
74  {
75  Xj[iCell] = Yj[iCell]/mixture_.rho(mySpecieI, p[iCell], T[iCell]);
76  }
77 
78  return (Xj/rhoInv);
79 }
80 
81 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
82 
85 (
86  const dictionary& dict,
87  const fvMesh& mesh
88 )
89 :
91  coeffsDict_((dict.subDict(typeName + "Coeffs"))),
92  thermo_(mesh.lookupObject<solidThermo>(basicThermo::dictName)),
93  speciesNames_(0),
94  mixture_(dynamic_cast<const basicSpecieMixture&>(thermo_)),
95  solidData_(mixture_.Y().size())
96 {
97  if (!isA<basicSpecieMixture>(thermo_))
98  {
100  << "Model requires a multi-component thermo package"
101  << abort(FatalError);
102  }
103 
104  label nFunc = 0;
105  const dictionary& functionDicts = dict.subDict(typeName + "Coeffs");
106 
107  forAllConstIter(dictionary, functionDicts, iter)
108  {
109  // safety:
110  if (!iter().isDict())
111  {
112  continue;
113  }
114  const word& key = iter().keyword();
115  if (!mixture_.contains(key))
116  {
118  << " specie: " << key << " is not found in the solid mixture"
119  << nl
120  << " specie is the mixture are:" << mixture_.species() << nl
121  << nl << endl;
122  }
123  speciesNames_.insert(key, nFunc);
124  const dictionary& dict = iter().dict();
125  dict.lookup("absorptivity") >> solidData_[nFunc][absorptivity];
126  dict.lookup("emissivity") >> solidData_[nFunc][emissivity];
127 
128  nFunc++;
129  }
130 }
131 
132 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
133 
136 {}
137 
138 
139 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
140 
143 calc(const label propertyId) const
144 {
146  (
147  new volScalarField
148  (
149  IOobject
150  (
151  "a",
152  mesh().time().timeName(),
153  mesh(),
156  ),
157  mesh(),
159  zeroGradientFvPatchVectorField::typeName
160  )
161  );
162 
163  scalarField& a = ta().internalField();
164 
165  forAllConstIter(HashTable<label>, speciesNames_, iter)
166  {
167  if (mixture_.contains(iter.key()))
168  {
169  a += solidData_[iter()][propertyId]*X(iter.key());
170  }
171  }
172 
173  ta().correctBoundaryConditions();
174  return ta;
175 }
176 
177 
180 (
181  const label bandI
182 ) const
183 {
184  return calc(emissivity);
185 }
186 
187 
190 (
191  const label bandI
192 ) const
193 {
194  return calc(absorptivity);
195 }
196 
197 // ************************************************************************* //
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
p
p
Definition: pEqn.H:62
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
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::calc
void calc(const argList &args, const Time &runTime, const fvMesh &mesh)
Definition: Lambda2.C:38
Foam::radiation::greyMeanSolidAbsorptionEmission::eCont
tmp< volScalarField > eCont(const label bandI=0) const
Emission coefficient for continuous phase.
Definition: greyMeanSolidAbsorptionEmission.C:180
unitConversion.H
Unit conversion functions.
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::radiation::greyMeanSolidAbsorptionEmission::calc
tmp< volScalarField > calc(const label) const
Calculate the property mixing.
Definition: greyMeanSolidAbsorptionEmission.C:143
greyMeanSolidAbsorptionEmission.H
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
dictName
const word dictName("particleTrackDict")
Foam::solidThermo
Fundamental solid thermodynamic properties.
Definition: solidThermo.H:54
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::radiation::greyMeanSolidAbsorptionEmission::X
tmp< scalarField > X(const word specie) const
Calculate the volumetric fraction of Yj.
Definition: greyMeanSolidAbsorptionEmission.C:50
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
Foam::specie
Base class of the thermophysical property types.
Definition: specie.H:57
radiation
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::radiation::greyMeanSolidAbsorptionEmission::aCont
tmp< volScalarField > aCont(const label bandI=0) const
Absorption coefficient for continuous phase.
Definition: greyMeanSolidAbsorptionEmission.C:190
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
Foam::radiation::greyMeanSolidAbsorptionEmission::greyMeanSolidAbsorptionEmission
greyMeanSolidAbsorptionEmission(const dictionary &dict, const fvMesh &mesh)
Construct from components.
Definition: greyMeanSolidAbsorptionEmission.C:85
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
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:78
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::radiation::defineTypeNameAndDebug
defineTypeNameAndDebug(cloudAbsorptionEmission, 0)
Foam::HashTable< label >
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::radiation::greyMeanSolidAbsorptionEmission::~greyMeanSolidAbsorptionEmission
virtual ~greyMeanSolidAbsorptionEmission()
Destructor.
Definition: greyMeanSolidAbsorptionEmission.C:135
T
const volScalarField & T
Definition: createFields.H:25
scalarField
volScalarField scalarField(fieldObject, mesh)
timeName
word timeName
Definition: getTimeIndex.H:3
Foam::radiation::absorptionEmissionModel
Model to supply absorption and emission coefficients for radiation modelling.
Definition: absorptionEmissionModel.H:52
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:259
Foam::radiation::addToRunTimeSelectionTable
addToRunTimeSelectionTable(absorptionEmissionModel, cloudAbsorptionEmission, dictionary)