Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
10License
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
31using 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,
149 bounds::repeatableBounding::CLAMP,
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,
163 (
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
188Foam::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 auto deltaLambdaT = volScalarField::New
205 (
206 "deltaLambdaT",
208 T.mesh(),
209 dimensionedScalar("deltaLambdaT", dimless, 1.0)
210 );
211
212 if (band != Vector2D<scalar>::one)
213 {
214 scalarField& deltaLambdaTf = deltaLambdaT.ref();
215
216 forAll(T, i)
217 {
218 deltaLambdaTf[i] = fLambdaT(band[1]*T[i]) - fLambdaT(band[0]*T[i]);
219 }
221
222 return deltaLambdaT;
223}
224
225
228(
229 const volScalarField& T,
230 const Vector2D<scalar>& band
231) const
232{
233 auto Eb = volScalarField::New
234 (
235 "Eb",
238 );
239
240 if (band != Vector2D<scalar>::one)
241 {
242 scalarField& Ebif = Eb.ref();
243
244 forAll(T, i)
245 {
246 Ebif[i] *= fLambdaT(band[1]*T[i]) - fLambdaT(band[0]*T[i]);
247 }
248
249 volScalarField::Boundary& EbBf = Eb.ref().boundaryFieldRef();
250
251 forAll(EbBf, patchi)
252 {
253 fvPatchScalarField& EbPf = EbBf[patchi];
254
255 if (!EbPf.coupled())
256 {
257 const scalarField& Tpf = T.boundaryField()[patchi];
258
259 forAll(EbPf, facei)
260 {
261 const scalar T1 = fLambdaT(band[1]*Tpf[facei]);
262 const scalar T2 = fLambdaT(band[0]*Tpf[facei]);
263
264 EbPf[facei] *= T1 - T2;
265 }
266 }
268 }
269
270 return Eb;
271}
272
273
275(
276 const label lambdaI,
277 const Vector2D<scalar>& band
278)
279{
280 bLambda_[lambdaI] = EbDeltaLambdaT(T_, band);
281}
282
283
284// ************************************************************************* //
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
Templated 2D Vector derived from VectorSpace adding construction from 2 components,...
Definition Vector2D.H:54
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
virtual bool coupled() const
True if the patch field is coupled.
blackBodyEmission(const label nLambda, const volScalarField &T)
Construct from components.
tmp< Foam::volScalarField > EbDeltaLambdaT(const volScalarField &T, const Vector2D< scalar > &band) const
Integral energy at T from lambda1 to lambda2.
tmp< Foam::volScalarField > deltaLambdaT(const volScalarField &T, const Vector2D< scalar > &band) const
Proportion of total energy at T from lambda1 to lambda2.
void correct(const label lambdaI, const Vector2D< scalar > &band)
static const List< Tuple2< scalar, scalar > > emissivePowerTable
Static table of black body emissive power.
A class for managing temporary objects.
Definition tmp.H:75
const volScalarField & T
Namespace for bounding specifications. At the moment, mostly for tables.
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
Different types of constants.
const dimensionSet dimless
Dimensionless.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar pow4(const dimensionedScalar &ds)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
fvPatchField< scalar > fvPatchScalarField
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299