cloudAbsorptionEmission.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2019 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 
27 #include "thermoCloud.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace radiationModels
35 {
36 namespace absorptionEmissionModels
37 {
38  defineTypeNameAndDebug(cloud, 0);
39 
41  (
42  absorptionEmissionModel,
43  cloud,
44  dictionary
45  );
46 }
47 }
48 }
49 
50 
51 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52 
54 (
55  const dictionary& dict,
56  const fvMesh& mesh
57 )
58 :
59  absorptionEmissionModel(dict, mesh),
60  coeffsDict_(dict.subDict(typeName + "Coeffs")),
61  cloudNames_(coeffsDict_.lookup("cloudNames"))
62 {}
63 
64 
65 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
66 
68 {}
69 
70 
71 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
72 
75 {
77  (
79  (
80  "a",
81  mesh_,
83  )
84  );
85 
86  forAll(cloudNames_, i)
87  {
88  const thermoCloud& tc
89  (
90  mesh_.objectRegistry::lookupObject<thermoCloud>(cloudNames_[i])
91  );
92 
93  ta.ref() += tc.ap();
94  }
95 
96  return ta;
97 }
98 
99 
102 (
103  const label bandI
104 ) const
105 {
107  (
109  (
110  "e",
111  mesh_,
113  )
114  );
115 
116  return te;
117 }
118 
119 
122 (
123  const label bandI
124 ) const
125 {
127  (
129  (
130  "E",
131  mesh_,
133  )
134  );
135 
136  forAll(cloudNames_, i)
137  {
138  const thermoCloud& tc
139  (
140  mesh_.objectRegistry::lookupObject<thermoCloud>(cloudNames_[i])
141  );
142 
143  tE.ref() += tc.Ep();
144  }
145 
146  // Total emission is 4 times the projected emission
147  return 4*tE;
148 }
149 
150 
151 // ************************************************************************* //
addToRunTimeSelectionTable(absorptionEmissionModel, greyMeanCombustion, dictionary)
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
virtual tmp< volScalarField > ap() const =0
Return tmp equivalent particulate absorption.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
virtual tmp< volScalarField > Ep() const =0
Return tmp equivalent particulate emission.
Macros for easy insertion into run-time selection tables.
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:934
virtual tmp< volScalarField > EDisp(const label bandI=0) const
Return emission contribution for dispersed phase.
cloud(const dictionary &dict, const fvMesh &mesh)
Construct from components.
Virtual abstract base class for templated ThermoCloud.
Definition: thermoCloud.H:48
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
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
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
virtual tmp< volScalarField > eDisp(const label bandI=0) const
Emission coefficient for dispersed phase.
virtual tmp< volScalarField > aDisp(const label bandI=0) const
Absorption coefficient for dispersed phase.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
Model to supply absorption and emission coefficients for radiation modelling.
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.