multiBandAbsorptionEmission.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(multiBandAbsorptionEmission, 0);
38 
40  (
41  absorptionEmissionModel,
42  multiBandAbsorptionEmission,
43  dictionary
44  );
45  }
46 }
47 
48 
49 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50 
52 (
53  const dictionary& dict,
54  const fvMesh& mesh
55 )
56 :
58  coeffsDict_(dict.subDict(typeName + "Coeffs")),
59  absCoeffs_(maxBands_),
60  emiCoeffs_(maxBands_),
61  nBands_(0)
62 {
63  coeffsDict_.readEntry("absorptivity", absCoeffs_);
64  coeffsDict_.readEntry("emissivity", emiCoeffs_);
65  nBands_ = absCoeffs_.size();
66 }
67 
68 
69 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
70 
73 {}
74 
75 
76 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
77 
80 (
81  const label bandI
82 ) const
83 {
84  return volScalarField::New
85  (
86  "a",
88  mesh(),
89  dimensionedScalar("a", dimless/dimLength, absCoeffs_[bandI])
90  );
91 }
92 
93 
96 (
97  const label bandI
98 ) const
99 {
100  return volScalarField::New
101  (
102  "e",
104  mesh(),
105  dimensionedScalar("e", dimless/dimLength, emiCoeffs_[bandI])
106  );
107 }
108 
109 
112 (
113  const label bandI
114 ) const
115 {
116  return volScalarField::New
117  (
118  "E",
120  mesh(),
122  );
123 }
124 
125 
126 // ************************************************************************* //
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
tmp< volScalarField > ECont(const label bandI) const
Emission contribution.
defineTypeNameAndDebug(cloudAbsorptionEmission, 0)
const dimensionSet dimless
Dimensionless.
Model to supply absorption and emission coefficients for radiation modelling.
Macros for easy insertion into run-time selection tables.
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...
dynamicFvMesh & mesh
tmp< volScalarField > eCont(const label bandI) const
Emission coefficient.
multiBandAbsorptionEmission(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...
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
addToRunTimeSelectionTable(absorptionEmissionModel, cloudAbsorptionEmission, dictionary)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
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.
tmp< volScalarField > aCont(const label bandI) const
Absorption coefficient.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127