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-2020 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 "parcelCloud.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace radiationModels
35 {
36 namespace absorptionEmissionModels
37 {
39 
41  (
43  cloud,
45  );
46 }
47 }
48 }
49 
50 
51 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52 
54 (
55  const dictionary& dict,
56  const fvMesh& mesh
57 )
58 :
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 parcelCloud& c =
89  mesh_.objectRegistry::lookupObject<parcelCloud>(cloudNames_[i]);
90 
91  ta.ref() += c.ap();
92  }
93 
94  return ta;
95 }
96 
97 
100 (
101  const label bandI
102 ) const
103 {
105  (
107  (
108  "e",
109  mesh_,
111  )
112  );
113 
114  return te;
115 }
116 
117 
120 (
121  const label bandI
122 ) const
123 {
125  (
127  (
128  "E",
129  mesh_,
131  )
132  );
133 
134  forAll(cloudNames_, i)
135  {
136  const parcelCloud& c =
137  mesh_.objectRegistry::lookupObject<parcelCloud>(cloudNames_[i]);
138 
139  tE.ref() += c.Ep();
140  }
141 
142  // Total emission is 4 times the projected emission
143  return 4*tE;
144 }
145 
146 
147 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
Virtual abstract base class for parcel clouds. As parcelCloudBase but with additional virtualisation ...
Definition: parcelCloud.H:57
Model to supply absorption and emission coefficients for radiation modelling.
Retrieves absorption/emission data from a cloud object.
virtual tmp< volScalarField > aDisp(const label bandI=0) const
Absorption coefficient for dispersed phase.
virtual tmp< volScalarField > EDisp(const label bandI=0) const
Return emission contribution for dispersed phase.
virtual tmp< volScalarField > eDisp(const label bandI=0) const
Emission coefficient for dispersed phase.
cloud(const dictionary &dict, const fvMesh &mesh)
Construct from components.
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
const dimensionedScalar c
Speed of light in a vacuum.
addToRunTimeSelectionTable(absorptionEmissionModel, greyMeanCombustion, dictionary)
Namespace for OpenFOAM.
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
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimless
const dimensionSet dimLength
const dimensionSet dimTime
const dimensionSet dimMass
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dictionary dict