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-2018 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 radiation
35  {
37 
39  (
43  );
44  }
45 }
46 
47 
48 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49 
51 (
52  const dictionary& dict,
53  const fvMesh& mesh
54 )
55 :
56  absorptionEmissionModel(dict, mesh),
57  coeffsDict_(dict.subDict(typeName + "Coeffs")),
58  cloudNames_(coeffsDict_.lookup("cloudNames"))
59 {}
60 
61 
62 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
63 
65 {}
66 
67 
68 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
69 
72 {
74  (
75  new volScalarField
76  (
77  IOobject
78  (
79  "a",
80  mesh_.time().timeName(),
81  mesh_,
84  false
85  ),
86  mesh_,
88  )
89  );
90 
91  forAll(cloudNames_, i)
92  {
93  const thermoCloud& tc
94  (
95  mesh_.objectRegistry::lookupObject<thermoCloud>(cloudNames_[i])
96  );
97 
98  ta.ref() += tc.ap();
99  }
100 
101  return ta;
102 }
103 
104 
107 {
109  (
110  new volScalarField
111  (
112  IOobject
113  (
114  "e",
115  mesh_.time().timeName(),
116  mesh_,
119  false
120  ),
121  mesh_,
123  )
124  );
125 
126  return te;
127 }
128 
129 
132 {
134  (
135  new volScalarField
136  (
137  IOobject
138  (
139  "E",
140  mesh_.time().timeName(),
141  mesh_,
144  false
145  ),
146  mesh_,
148  )
149  );
150 
151  forAll(cloudNames_, i)
152  {
153  const thermoCloud& tc
154  (
155  mesh_.objectRegistry::lookupObject<thermoCloud>(cloudNames_[i])
156  );
157 
158  tE.ref() += tc.Ep();
159  }
160 
161  // Total emission is 4 times the projected emission
162  return 4*tE;
163 }
164 
165 
166 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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.
virtual tmp< volScalarField > eDisp(const label bandI=0) const
Emission coefficient for dispersed phase.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
virtual tmp< volScalarField > EDisp(const label bandI=0) const
Return emission contribution for dispersed phase.
defineTypeNameAndDebug(cloudAbsorptionEmission, 0)
virtual tmp< volScalarField > Ep() const =0
Return tmp equivalent particulate emission.
Model to supply absorption and emission coefficients for radiation modelling.
Macros for easy insertion into run-time selection tables.
virtual tmp< volScalarField > aDisp(const label bandI=0) const
Absorption coefficient for dispersed phase.
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:692
Retrieves absorption/emission data from a cloud object.
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.
cloudAbsorptionEmission(const dictionary &dict, const fvMesh &mesh)
Construct from components.
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: PtrList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
Namespace for OpenFOAM.