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  zoneEmissivity_(),
65  zoneIds_()
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", zoneEmissivity_);
75 
76  zoneIds_.resize(zoneAbsorptivity_.size(), -1);
77 
78  label numZones = 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  zoneIds_[numZones] = zoneID;
90  ++numZones;
91  }
92  // zoneIds_.resize(numZones);
93 }
94 
95 
96 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
97 
100 {}
101 
102 
103 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
104 
107 (
108  const label bandI
109 ) const
110 {
111  auto ta = volScalarField::New
112  (
113  "a",
115  mesh(),
116  dimensionedScalar("a", dimless/dimLength, absCoeffs_[bandI])
117  );
118  scalarField& a = ta.ref().primitiveFieldRef();
119 
120  for (const label zonei : zoneIds_)
121  {
122  const cellZone& zn = mesh().cellZones()[zonei];
123  const auto iter = zoneAbsorptivity_.cfind(zn.name());
124 
125  if (iter.good()) // Check is redundant (cannot fail)
126  {
127  const scalarList& absorb = iter.val();
128 
129  UIndirectList<scalar>(a, zn) = absorb[bandI];
130  }
131  }
133  return ta;
134 }
135 
136 
139 (
140  const label bandI
141 ) const
142 {
143  auto te = volScalarField::New
144  (
145  "e",
147  mesh(),
148  dimensionedScalar("e", dimless/dimLength, emiCoeffs_[bandI])
149  );
150  scalarField& e = te.ref().primitiveFieldRef();
151 
152  for (const label zonei : zoneIds_)
153  {
154  const cellZone& zn = mesh().cellZones()[zonei];
155  const auto iter = zoneEmissivity_.cfind(zn.name());
156 
157  if (iter.good()) // Check is redundant (cannot fail)
158  {
159  const scalarList& emit = iter.val();
160 
161  UIndirectList<scalar>(e, zn) = emit[bandI];
162  }
163  }
165  return te;
166 }
167 
168 
171 (
172  const label bandI
173 ) const
174 {
175  return volScalarField::New
176  (
177  "E",
179  mesh(),
181  );
182 }
183 
184 
185 // ************************************************************************* //
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:153
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:608
defineTypeNameAndDebug(cloudAbsorptionEmission, 0)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
const fvMesh & mesh() const
Reference to the mesh.
const dimensionSet dimless
Dimensionless.
Model to supply absorption and emission coefficients for radiation modelling.
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:441
Macros for easy insertion into run-time selection tables.
tmp< volScalarField > eCont(const label bandI) const
Emission coefficient.
tmp< volScalarField > aCont(const label bandI) const
Absorption coefficient.
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, IOobjectOption::readOption readOpt=IOobjectOption::MUST_READ) const
Find entry and assign to T val. FatalIOError if it is found and the number of tokens is incorrect...
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
dynamicFvMesh & mesh
label findZoneID(const word &zoneName) const
Find zone index by name, return -1 if not found.
Definition: ZoneMesh.C:718
multiBandZoneAbsorptionEmission(const dictionary &dict, const fvMesh &mesh)
Construct from components.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
tmp< volScalarField > ECont(const label bandI) const
Emission contribution.
A subset of mesh cells.
Definition: cellZone.H:58
dimensionedScalar pow3(const dimensionedScalar &ds)
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
A List with indirect addressing. Like IndirectList but does not store addressing. ...
Definition: faMatrix.H:52
addToRunTimeSelectionTable(absorptionEmissionModel, cloudAbsorptionEmission, dictionary)
const word & name() const noexcept
The zone name.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
Definition: polyMesh.H:679
wordList names() const
A list of the zone names.
Definition: ZoneMesh.C:451
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Do not request registration (bool: false)
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
Namespace for OpenFOAM.
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127