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-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 
26 #include "blackBodyEmission.H"
28 
29 using namespace Foam::constant;
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
35 (
36  IStringStream
37  (
38  "("
39  "( 1000 0.00032)"
40  "( 1100 0.00091)"
41  "( 1200 0.00213)"
42  "( 1300 0.00432)"
43  "( 1400 0.00779)"
44  "( 1500 0.01280)"
45  "( 1600 0.01972)"
46  "( 1700 0.02853)"
47  "( 1800 0.03934)"
48  "( 1900 0.05210)"
49  "( 2000 0.06672)"
50  "( 2100 0.08305)"
51  "( 2200 0.10088)"
52  "( 2300 0.12002)"
53  "( 2400 0.14025)"
54  "( 2500 0.16135)"
55  "( 2600 0.18311)"
56  "( 2700 0.20535)"
57  "( 2800 0.22788)"
58  "( 2900 0.25055)"
59  "( 3000 0.27322)"
60  "( 3100 0.29576)"
61  "( 3200 0.31809)"
62  "( 3300 0.34009)"
63  "( 3400 0.36172)"
64  "( 3500 0.38290)"
65  "( 3600 0.40359)"
66  "( 3700 0.42375)"
67  "( 3800 0.44336)"
68  "( 3900 0.46240)"
69  "( 4000 0.48085)"
70  "( 4100 0.49872)"
71  "( 4200 0.51599)"
72  "( 4300 0.53267)"
73  "( 4400 0.54877)"
74  "( 4500 0.56429)"
75  "( 4600 0.57925)"
76  "( 4700 0.59366)"
77  "( 4800 0.60753)"
78  "( 4900 0.62088)"
79  "( 5000 0.63372)"
80  "( 5100 0.64606)"
81  "( 5200 0.65794)"
82  "( 5300 0.66935)"
83  "( 5400 0.68033)"
84  "( 5500 0.69087)"
85  "( 5600 0.70101)"
86  "( 5700 0.71076)"
87  "( 5800 0.72012)"
88  "( 5900 0.72913)"
89  "( 6000 0.73778)"
90  "( 6100 0.74610)"
91  "( 6200 0.75410)"
92  "( 6300 0.76180)"
93  "( 6400 0.76920)"
94  "( 6500 0.77631)"
95  "( 6600 0.78316)"
96  "( 6700 0.78975)"
97  "( 6800 0.79609)"
98  "( 6900 0.80219)"
99  "( 7000 0.80807)"
100  "( 7100 0.81373)"
101  "( 7200 0.81918)"
102  "( 7300 0.82443)"
103  "( 7400 0.82949)"
104  "( 7500 0.83436)"
105  "( 7600 0.83906)"
106  "( 7700 0.84359)"
107  "( 7800 0.84796)"
108  "( 7900 0.85218)"
109  "( 8000 0.85625)"
110  "( 8100 0.86017)"
111  "( 8200 0.86396)"
112  "( 8300 0.86762)"
113  "( 8400 0.87115)"
114  "( 8500 0.87456)"
115  "( 8600 0.87786)"
116  "( 8700 0.88105)"
117  "( 8800 0.88413)"
118  "( 8900 0.88711)"
119  "( 9000 0.88999)"
120  "( 9100 0.89277)"
121  "( 9200 0.89547)"
122  "( 9300 0.89807)"
123  "( 9400 0.90060)"
124  "( 9500 0.90304)"
125  "( 9600 0.90541)"
126  "( 9700 0.90770)"
127  "( 9800 0.90992)"
128  "( 9900 0.91207)"
129  "(10000 0.91415)"
130  "(12000 0.94505)"
131  "(15000 0.96893)"
132  "(20000 0.98555)"
133  "(30000 0.99529)"
134  "(40000 0.99792)"
135  "(50000 0.99890)"
136  ")"
137  )()
138 );
139 
140 
141 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
142 
144 (
145  const label nLambda,
146  const volScalarField& T
147 )
148 :
149  table_
150  (
151  emissivePowerTable,
153  "blackBodyEmissivePower"
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::radiation::blackBodyEmission::fLambdaT
192 (
193  const scalar lambdaT
194 ) const
195 {
196  return table_(1e6*lambdaT);
197 }
198 
199 
202 (
203  const volScalarField& T,
204  const Vector2D<scalar>& band
205 ) const
206 {
207  tmp<volScalarField> deltaLambdaT
208  (
209  new volScalarField
210  (
211  IOobject
212  (
213  "deltaLambdaT",
214  T.mesh().time().timeName(),
215  T.mesh(),
218  ),
219  T.mesh(),
220  dimensionedScalar("deltaLambdaT", dimless, 1.0)
221  )
222  );
223 
224  if (band != Vector2D<scalar>::one)
225  {
226  scalarField& deltaLambdaTf = deltaLambdaT.ref();
227 
228  forAll(T, i)
229  {
230  deltaLambdaTf[i] = fLambdaT(band[1]*T[i]) - fLambdaT(band[0]*T[i]);
231  }
232  }
233 
234  return deltaLambdaT;
235 }
236 
237 
240 (
241  const volScalarField& T,
242  const Vector2D<scalar>& band
243 ) const
244 {
246  (
247  new volScalarField
248  (
249  IOobject
250  (
251  "Eb",
252  T.mesh().time().timeName(),
253  T.mesh(),
256  ),
258  )
259  );
260 
261  if (band != Vector2D<scalar>::one)
262  {
263  scalarField& Ebif = Eb.ref();
264 
265  forAll(T, i)
266  {
267  Ebif[i] *= fLambdaT(band[1]*T[i]) - fLambdaT(band[0]*T[i]);
268  }
269 
270  volScalarField::Boundary& EbBf = Eb.ref().boundaryFieldRef();
271 
272  forAll(EbBf, patchi)
273  {
274  fvPatchScalarField& EbPf = EbBf[patchi];
275 
276  if (!EbPf.coupled())
277  {
278  const scalarField& Tpf = T.boundaryField()[patchi];
279 
280  forAll(EbPf, facei)
281  {
282  const scalar T1 = fLambdaT(band[1]*Tpf[facei]);
283  const scalar T2 = fLambdaT(band[0]*Tpf[facei]);
284 
285  EbPf[facei] *= T1 - T2;
286  }
287  }
288  }
289  }
290 
291  return Eb;
292 }
293 
294 
296 (
297  const label lambdaI,
298  const Vector2D<scalar>& band
299 )
300 {
301  bLambda_[lambdaI] = EbDeltaLambdaT(T_, band);
302 }
303 
304 
305 // ************************************************************************* //
Collection of constants.
static const List< Tuple2< scalar, scalar > > emissivePowerTable
Static table of black body emissive power.
#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
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:60
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
void correct(const label lambdaI, const Vector2D< scalar > &band)
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
tmp< Foam::volScalarField > deltaLambdaT(const volScalarField &T, const Vector2D< scalar > &band) const
Proportion of total energy at T from lambda1 to lambda2.
Dimension set for the base types.
Definition: dimensionSet.H:120
blackBodyEmission(const label nLambda, const volScalarField &T)
Construct from components.
An interpolation/look-up table of scalar vs <Type> values. The reference scalar values must be monoto...
virtual bool coupled() const
Return true if this patch field is coupled.
Definition: fvPatchField.H:328
const Mesh & mesh() const
Return mesh.
tmp< Foam::volScalarField > EbDeltaLambdaT(const volScalarField &T, const Vector2D< scalar > &band) const
Integral energy at T from lambda1 to lambda2.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
label patchi
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensionedScalar pow4(const dimensionedScalar &ds)
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:98
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