blackBodyEmission.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2018 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "blackBodyEmission.H"
30 
31 using namespace Foam::constant;
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
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 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
139 
141 (
142  const label nLambda,
143  const volScalarField& T
144 )
145 :
146  table_
147  (
148  emissivePowerTable,
150  "blackBodyEmissivePower"
151  ),
152  C1_("C1", dimensionSet(1, 4, 3, 0, 0, 0, 0), 3.7419e-16),
153  C2_("C2", dimensionSet(0, 1, 0, 1, 0, 0, 0), 14.388e-6),
154  bLambda_(nLambda),
155  T_(T)
156 {
157  forAll(bLambda_, lambdaI)
158  {
159  bLambda_.set
160  (
161  lambdaI,
162  new volScalarField
163  (
164  IOobject
165  (
166  "bLambda_" + Foam::name(lambdaI) ,
167  T.mesh().time().timeName(),
168  T.mesh(),
171  ),
173 
174  )
175  );
176  }
177 }
178 
179 
180 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
181 
183 {}
184 
185 
186 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
187 
188 Foam::scalar Foam::radiation::blackBodyEmission::fLambdaT
189 (
190  const scalar lambdaT
191 ) const
192 {
193  return table_(1e6*lambdaT);
194 }
195 
196 
199 (
200  const volScalarField& T,
201  const Vector2D<scalar>& band
202 ) const
203 {
204  tmp<volScalarField> deltaLambdaT
205  (
206  new volScalarField
207  (
208  IOobject
209  (
210  "deltaLambdaT",
211  T.mesh().time().timeName(),
212  T.mesh(),
215  ),
216  T.mesh(),
217  dimensionedScalar("deltaLambdaT", dimless, 1.0)
218  )
219  );
220 
221  if (band != Vector2D<scalar>::one)
222  {
223  scalarField& deltaLambdaTf = deltaLambdaT.ref();
224 
225  forAll(T, i)
226  {
227  deltaLambdaTf[i] = fLambdaT(band[1]*T[i]) - fLambdaT(band[0]*T[i]);
228  }
229  }
230 
231  return deltaLambdaT;
232 }
233 
234 
237 (
238  const volScalarField& T,
239  const Vector2D<scalar>& band
240 ) const
241 {
243  (
244  new volScalarField
245  (
246  IOobject
247  (
248  "Eb",
249  T.mesh().time().timeName(),
250  T.mesh(),
253  ),
255  )
256  );
257 
258  if (band != Vector2D<scalar>::one)
259  {
260  scalarField& Ebif = Eb.ref();
261 
262  forAll(T, i)
263  {
264  Ebif[i] *= fLambdaT(band[1]*T[i]) - fLambdaT(band[0]*T[i]);
265  }
266 
267  volScalarField::Boundary& EbBf = Eb.ref().boundaryFieldRef();
268 
269  forAll(EbBf, patchi)
270  {
271  fvPatchScalarField& EbPf = EbBf[patchi];
272 
273  if (!EbPf.coupled())
274  {
275  const scalarField& Tpf = T.boundaryField()[patchi];
276 
277  forAll(EbPf, facei)
278  {
279  const scalar T1 = fLambdaT(band[1]*Tpf[facei]);
280  const scalar T2 = fLambdaT(band[0]*Tpf[facei]);
281 
282  EbPf[facei] *= T1 - T2;
283  }
284  }
285  }
286  }
287 
288  return Eb;
289 }
290 
291 
293 (
294  const label lambdaI,
295  const Vector2D<scalar>& band
296 )
297 {
298  bLambda_[lambdaI] = EbDeltaLambdaT(T_, band);
299 }
300 
301 
302 // ************************************************************************* //
Foam::fvPatchField< scalar >
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::radiation::blackBodyEmission::EbDeltaLambdaT
tmp< Foam::volScalarField > EbDeltaLambdaT(const volScalarField &T, const Vector2D< scalar > &band) const
Integral energy at T from lambda1 to lambda2.
Definition: blackBodyEmission.C:237
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::constant
Different types of constants.
Definition: atomicConstants.C:38
Foam::dimensionSet
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:108
Foam::Vector2D< scalar >
Foam::radiation::blackBodyEmission::blackBodyEmission
blackBodyEmission(const label nLambda, const volScalarField &T)
Construct from components.
Definition: blackBodyEmission.C:141
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::fvPatchField::coupled
virtual bool coupled() const
Return true if this patch field is coupled.
Definition: fvPatchField.H:345
Foam::pow4
dimensionedScalar pow4(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:100
blackBodyEmission.H
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:227
Foam::Field< scalar >
Foam::radiation::blackBodyEmission::correct
void correct(const label lambdaI, const Vector2D< scalar > &band)
Definition: blackBodyEmission.C:293
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
Foam::bounds::normalBounding::CLAMP
Clamp value to the start/end value.
physicoChemicalConstants.H
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:63
Foam::radiation::blackBodyEmission::emissivePowerTable
static const List< Tuple2< scalar, scalar > > emissivePowerTable
Static table of black body emissive power.
Definition: blackBodyEmission.H:63
Foam::radiation::blackBodyEmission::~blackBodyEmission
~blackBodyEmission()
Destructor.
Definition: blackBodyEmission.C:182
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::constant::physicoChemical::sigma
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::radiation::blackBodyEmission::deltaLambdaT
tmp< Foam::volScalarField > deltaLambdaT(const volScalarField &T, const Vector2D< scalar > &band) const
Proportion of total energy at T from lambda1 to lambda2.
Definition: blackBodyEmission.C:199
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:188
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189