blackBodyEmission.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-2021 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 
26 #include "blackBodyEmission.H"
29 
30 using namespace Foam::constant;
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
36 (
37  {
38  { 1000, 0.00032},
39  { 1100, 0.00091},
40  { 1200, 0.00213},
41  { 1300, 0.00432},
42  { 1400, 0.00779},
43  { 1500, 0.01280},
44  { 1600, 0.01972},
45  { 1700, 0.02853},
46  { 1800, 0.03934},
47  { 1900, 0.05210},
48  { 2000, 0.06672},
49  { 2100, 0.08305},
50  { 2200, 0.10088},
51  { 2300, 0.12002},
52  { 2400, 0.14025},
53  { 2500, 0.16135},
54  { 2600, 0.18311},
55  { 2700, 0.20535},
56  { 2800, 0.22788},
57  { 2900, 0.25055},
58  { 3000, 0.27322},
59  { 3100, 0.29576},
60  { 3200, 0.31809},
61  { 3300, 0.34009},
62  { 3400, 0.36172},
63  { 3500, 0.38290},
64  { 3600, 0.40359},
65  { 3700, 0.42375},
66  { 3800, 0.44336},
67  { 3900, 0.46240},
68  { 4000, 0.48085},
69  { 4100, 0.49872},
70  { 4200, 0.51599},
71  { 4300, 0.53267},
72  { 4400, 0.54877},
73  { 4500, 0.56429},
74  { 4600, 0.57925},
75  { 4700, 0.59366},
76  { 4800, 0.60753},
77  { 4900, 0.62088},
78  { 5000, 0.63372},
79  { 5100, 0.64606},
80  { 5200, 0.65794},
81  { 5300, 0.66935},
82  { 5400, 0.68033},
83  { 5500, 0.69087},
84  { 5600, 0.70101},
85  { 5700, 0.71076},
86  { 5800, 0.72012},
87  { 5900, 0.72913},
88  { 6000, 0.73778},
89  { 6100, 0.74610},
90  { 6200, 0.75410},
91  { 6300, 0.76180},
92  { 6400, 0.76920},
93  { 6500, 0.77631},
94  { 6600, 0.78316},
95  { 6700, 0.78975},
96  { 6800, 0.79609},
97  { 6900, 0.80219},
98  { 7000, 0.80807},
99  { 7100, 0.81373},
100  { 7200, 0.81918},
101  { 7300, 0.82443},
102  { 7400, 0.82949},
103  { 7500, 0.83436},
104  { 7600, 0.83906},
105  { 7700, 0.84359},
106  { 7800, 0.84796},
107  { 7900, 0.85218},
108  { 8000, 0.85625},
109  { 8100, 0.86017},
110  { 8200, 0.86396},
111  { 8300, 0.86762},
112  { 8400, 0.87115},
113  { 8500, 0.87456},
114  { 8600, 0.87786},
115  { 8700, 0.88105},
116  { 8800, 0.88413},
117  { 8900, 0.88711},
118  { 9000, 0.88999},
119  { 9100, 0.89277},
120  { 9200, 0.89547},
121  { 9300, 0.89807},
122  { 9400, 0.90060},
123  { 9500, 0.90304},
124  { 9600, 0.90541},
125  { 9700, 0.90770},
126  { 9800, 0.90992},
127  { 9900, 0.91207},
128  {10000, 0.91415},
129  {12000, 0.94505},
130  {15000, 0.96893},
131  {20000, 0.98555},
132  {30000, 0.99529},
133  {40000, 0.99792},
134  {50000, 0.99890}
135  }
136 );
137 
138 
139 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
140 
142 (
143  const label nLambda,
144  const volScalarField& T
145 )
146 :
147  table_
148  (
149  "blackBodyEmissivePower",
151  linearInterpolationWeights::typeName,
152  autoPtr<TableReader<scalar>>(nullptr),
153  emissivePowerTable
154  ),
155  C1_("C1", dimensionSet(1, 4, 3, 0, 0, 0, 0), 3.7419e-16),
156  C2_("C2", dimensionSet(0, 1, 0, 1, 0, 0, 0), 14.388e-6),
157  bLambda_(nLambda),
158  T_(T)
159 {
160  forAll(bLambda_, lambdaI)
161  {
162  bLambda_.set
163  (
164  lambdaI,
165  new volScalarField
166  (
167  IOobject
168  (
169  "bLambda_" + Foam::name(lambdaI) ,
170  T.mesh().time().timeName(),
171  T.mesh(),
174  ),
176 
177  )
178  );
179  }
180 }
181 
182 
183 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
184 
186 {}
187 
188 
189 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
190 
191 Foam::scalar Foam::radiationModels::blackBodyEmission::fLambdaT
192 (
193  const scalar lambdaT
194 ) const
195 {
196  return table_.value(1e6*lambdaT);
197 }
198 
199 
202 (
203  const volScalarField& T,
204  const Vector2D<scalar>& band
205 ) const
206 {
207  tmp<volScalarField> deltaLambdaT
208  (
210  (
211  "deltaLambdaT",
212  T.mesh(),
214  )
215  );
216 
217  if (band != Vector2D<scalar>::one)
218  {
219  scalarField& deltaLambdaTf = deltaLambdaT.ref();
220 
221  forAll(T, i)
222  {
223  deltaLambdaTf[i] = fLambdaT(band[1]*T[i]) - fLambdaT(band[0]*T[i]);
224  }
225  }
226 
227  return deltaLambdaT;
228 }
229 
230 
233 (
234  const volScalarField& T,
235  const Vector2D<scalar>& band
236 ) const
237 {
239  (
241  (
242  "Eb",
244  )
245  );
246 
247  if (band != Vector2D<scalar>::one)
248  {
249  scalarField& Ebif = Eb.ref();
250 
251  forAll(T, i)
252  {
253  Ebif[i] *= fLambdaT(band[1]*T[i]) - fLambdaT(band[0]*T[i]);
254  }
255 
256  volScalarField::Boundary& EbBf = Eb.ref().boundaryFieldRef();
257 
258  forAll(EbBf, patchi)
259  {
260  fvPatchScalarField& EbPf = EbBf[patchi];
261 
262  if (!EbPf.coupled())
263  {
264  const scalarField& Tpf = T.boundaryField()[patchi];
265 
266  forAll(EbPf, facei)
267  {
268  const scalar T1 = fLambdaT(band[1]*Tpf[facei]);
269  const scalar T2 = fLambdaT(band[0]*Tpf[facei]);
270 
271  EbPf[facei] *= T1 - T2;
272  }
273  }
274  }
275  }
276 
277  return Eb;
278 }
279 
280 
282 (
283  const label lambdaI,
284  const Vector2D<scalar>& band
285 )
286 {
287  bLambda_[lambdaI] = EbDeltaLambdaT(T_, band);
288 }
289 
290 
291 // ************************************************************************* //
Collection of constants.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
tmp< Foam::volScalarField > deltaLambdaT(const volScalarField &T, const Vector2D< scalar > &band) const
Proportion of total energy at T from lambda1 to lambda2.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
const dimensionSet dimless
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
void correct(const label lambdaI, const Vector2D< scalar > &band)
Dimension set for the base types.
Definition: dimensionSet.H:121
blackBodyEmission(const label nLambda, const volScalarField &T)
Construct from components.
virtual bool coupled() const
Return true if this patch field is coupled.
Definition: fvPatchField.H:323
static const List< Tuple2< scalar, scalar > > emissivePowerTable
Static table of black body emissive power.
const Mesh & mesh() const
Return mesh.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
label patchi
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensionedScalar pow4(const dimensionedScalar &ds)
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:105
A class for managing temporary objects.
Definition: PtrList.H:53
tmp< Foam::volScalarField > EbDeltaLambdaT(const volScalarField &T, const Vector2D< scalar > &band) const
Integral energy at T from lambda1 to lambda2.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98