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-2024 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  cloudNames_(dict.lookup("cloudNames"))
61 {}
62 
63 
64 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
65 
67 {}
68 
69 
70 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
71 
74 {
76  (
78  (
79  "a",
80  mesh_,
82  )
83  );
84 
85  forAll(cloudNames_, i)
86  {
87  const parcelCloud& c =
88  mesh_.objectRegistry::lookupObject<parcelCloud>(cloudNames_[i]);
89 
90  ta.ref() += c.ap();
91  }
92 
93  return ta;
94 }
95 
96 
99 (
100  const label bandI
101 ) const
102 {
104  (
106  (
107  "e",
108  mesh_,
110  )
111  );
112 
113  return te;
114 }
115 
116 
119 (
120  const label bandI
121 ) const
122 {
124  (
126  (
127  "E",
128  mesh_,
130  )
131  );
132 
133  forAll(cloudNames_, i)
134  {
135  const parcelCloud& c =
136  mesh_.objectRegistry::lookupObject<parcelCloud>(cloudNames_[i]);
137 
138  tE.ref() += c.Ep();
139  }
140 
141  // Total emission is 4 times the projected emission
142  return 4*tE;
143 }
144 
145 
146 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Macros for easy insertion into run-time selection tables.
static tmp< GeometricField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Internal &, const PtrList< Patch > &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
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:197
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
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
const dimensionSet dimless
const dimensionSet dimLength
const dimensionSet dimTime
void pow3(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
const dimensionSet dimMass
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dictionary dict