multiBandZoneAbsorptionEmission.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) 2015-2019 OpenCFD Ltd.
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 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35  namespace radiation
36  {
37  defineTypeNameAndDebug(multiBandZoneAbsorptionEmission, 0);
38 
40  (
41  absorptionEmissionModel,
42  multiBandZoneAbsorptionEmission,
43  dictionary
44  );
45  }
46 }
47 
48 
49 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50 
53 (
54  const dictionary& dict,
55  const fvMesh& mesh
56 )
57 :
59  coeffsDict_(dict.subDict(typeName + "Coeffs")),
60  absCoeffs_(maxBands_),
61  emiCoeffs_(maxBands_),
62  nBands_(0),
63  zoneAbsorptivity_(),
64  zoneEmisivity_(),
65  zoneCells_()
66 {
67  coeffsDict_.readEntry("absorptivity", absCoeffs_);
68  coeffsDict_.readEntry("emissivity", emiCoeffs_);
69  nBands_ = absCoeffs_.size();
70 
71  const dictionary& zoneDict = coeffsDict_.subDict("zones");
72 
73  zoneDict.readEntry("absorptivity", zoneAbsorptivity_);
74  zoneDict.readEntry("emissivity", zoneEmisivity_);
75 
76  zoneCells_.setSize(zoneAbsorptivity_.size(), -1);
77 
78  label i = 0;
79  forAllConstIters(zoneAbsorptivity_, iter)
80  {
81  label zoneID = mesh.cellZones().findZoneID(iter.key());
82  if (zoneID == -1)
83  {
85  << "Cannot find cellZone " << iter.key() << endl
86  << "Valid cellZones are " << mesh.cellZones().names()
87  << exit(FatalError);
88  }
89  zoneCells_[i++] = zoneID;
90  }
91 
92 }
93 
94 
95 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
96 
99 {}
100 
101 
102 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
103 
106 (
107  const label bandI
108 ) const
109 {
111  (
112  new volScalarField
113  (
114  IOobject
115  (
116  "a",
117  mesh().time().timeName(),
118  mesh(),
121  ),
122  mesh(),
123  dimensionedScalar("a", dimless/dimLength, absCoeffs_[bandI])
124  )
125  );
126 
127  volScalarField& a = ta.ref();
128 
129  forAll(zoneCells_, zonei)
130  {
131  const cellZone& cZone = mesh().cellZones()[zoneCells_[zonei]];
132 
133  tmp<volScalarField> tzoneAbs(a*0.0);
134  volScalarField& zoneAbs = tzoneAbs.ref();
135 
136  const scalarList& abs = zoneAbsorptivity_.find(cZone.name())();
137 
138  forAll(cZone, i)
139  {
140  label cellId = cZone[i];
141  zoneAbs[cellId] = abs[bandI] - absCoeffs_[bandI];
142  }
143 
144  a += zoneAbs;
145  }
146 
147  return ta;
148 }
149 
150 
153 (
154  const label bandI
155 ) const
156 {
158  (
159  new volScalarField
160  (
161  IOobject
162  (
163  "e",
164  mesh().time().timeName(),
165  mesh(),
168  ),
169  mesh(),
170  dimensionedScalar("e", dimless/dimLength, emiCoeffs_[bandI])
171  )
172  );
173 
174  volScalarField& e = te.ref();
175 
176  forAll(zoneCells_, zonei)
177  {
178  const cellZone& cZone = mesh().cellZones()[zoneCells_[zonei]];
179 
180  tmp<volScalarField> tzoneEm(e*0.0);
181  volScalarField& zoneEm = tzoneEm.ref();
182 
183  const scalarList& emi = zoneEmisivity_.find(cZone.name())();
184 
185  forAll(cZone, i)
186  {
187  label cellId = cZone[i];
188  zoneEm[cellId] = emi[bandI] - emiCoeffs_[bandI];
189  }
190  e += zoneEm;
191  }
192 
193 
194  return te;
195 }
196 
197 
200 (
201  const label bandI
202 ) const
203 {
205  (
206  new volScalarField
207  (
208  IOobject
209  (
210  "E",
211  mesh().time().timeName(),
212  mesh(),
215  ),
216  mesh(),
218  )
219  );
220 
221  return E;
222 }
223 
224 
225 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::radiation::multiBandZoneAbsorptionEmission::multiBandZoneAbsorptionEmission
multiBandZoneAbsorptionEmission(const dictionary &dict, const fvMesh &mesh)
Construct from components.
Definition: multiBandZoneAbsorptionEmission.C:53
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::radiation::multiBandZoneAbsorptionEmission::~multiBandZoneAbsorptionEmission
virtual ~multiBandZoneAbsorptionEmission()
Destructor.
Definition: multiBandZoneAbsorptionEmission.C:98
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::cellZone
A subset of mesh cells.
Definition: cellZone.H:62
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::radiation::multiBandZoneAbsorptionEmission::ECont
tmp< volScalarField > ECont(const label bandI) const
Emission contribution.
Definition: multiBandZoneAbsorptionEmission.C:200
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:227
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:89
multiBandZoneAbsorptionEmission.H
Foam::dictionary::subDict
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:460
Foam::dictionary::readEntry
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
Definition: dictionaryTemplates.C:302
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
radiation
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
timeName
word timeName
Definition: getTimeIndex.H:3
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:123
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
zoneID
const labelIOList & zoneID
Definition: interpolatedFaces.H:22
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
cellId
label cellId
Definition: interrogateWallPatches.H:67
Foam::radiation::defineTypeNameAndDebug
defineTypeNameAndDebug(cloudAbsorptionEmission, 0)
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
forAllConstIters
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
Foam::zoneIdentifier::name
const word & name() const noexcept
The zone name.
Definition: zoneIdentifier.H:123
Foam::List< scalar >
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::radiation::absorptionEmissionModel
Model to supply absorption and emission coefficients for radiation modelling.
Definition: absorptionEmissionModel.H:54
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:188
Foam::radiation::addToRunTimeSelectionTable
addToRunTimeSelectionTable(absorptionEmissionModel, cloudAbsorptionEmission, dictionary)
Foam::radiation::multiBandZoneAbsorptionEmission::eCont
tmp< volScalarField > eCont(const label bandI) const
Emission coefficient.
Definition: multiBandZoneAbsorptionEmission.C:153
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189
Foam::radiation::multiBandZoneAbsorptionEmission::aCont
tmp< volScalarField > aCont(const label bandI) const
Absorption coefficient.
Definition: multiBandZoneAbsorptionEmission.C:106