initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / thermophysicalModels / radiation / radiationModel / fvDOM / blackBodyEmission / blackBodyEmission.C
blobbc21f3bfe0bfc027479aef047315895b4db061c4
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2008-2009 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     option) any later version.
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.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "blackBodyEmission.H"
28 #include "dimensionedConstants.H"
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 const Foam::List<Foam::Tuple2<Foam::scalar, Foam::scalar> >
33 Foam::radiation::blackBodyEmission::emissivePowerTable
35     IStringStream
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     )()
140 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
142 Foam::radiation::blackBodyEmission::blackBodyEmission
144     const label nLambda,
145     const volScalarField& T
148     table_
149     (
150         emissivePowerTable,
151         interpolationTable<scalar>::CLAMP,
152         "blackBodyEmissivePower"
153     ),
154     C1_("C1", dimensionSet(1, 4, 3, 0, 0, 0, 0), 3.7419e-16),
155     C2_("C2", dimensionSet(0, 1, 0, 1, 0, 0, 0), 14.388e-6),
156     bLambda_(nLambda),
157     T_(T)
159     forAll(bLambda_, lambdaI)
160     {
161         bLambda_.set
162         (
163             lambdaI,
164             new volScalarField
165             (
166                 IOobject
167                 (
168                     "bLambda_" + Foam::name(lambdaI) ,
169                     T.mesh().time().timeName(),
170                     T.mesh(),
171                     IOobject::NO_READ,
172                     IOobject::NO_WRITE
173                 ),
174                 radiation::sigmaSB*pow4(T)
175             )
176         );
178     }
182 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
184 Foam::radiation::blackBodyEmission::~blackBodyEmission()
188 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
190 Foam::scalar Foam::radiation::blackBodyEmission::fLambdaT
192     const scalar lambdaT
193 ) const
195     return  table_(lambdaT*1.0e6);
199 Foam::tmp<Foam::volScalarField>
200 Foam::radiation::blackBodyEmission::EbDeltaLambdaT
202     const volScalarField& T,
203     const Vector2D<scalar>& band
204 ) const
206     tmp<volScalarField> Eb
207     (
208         new volScalarField
209         (
210             IOobject
211             (
212                 "Eb",
213                 T.mesh().time().timeName(),
214                 T.mesh(),
215                 IOobject::NO_READ,
216                 IOobject::NO_WRITE
217             ),
218             radiation::sigmaSB*pow4(T)
219         )
220     );
223     if (band == Vector2D<scalar>::one)
224     {
225         return Eb;
226     }
227     else
228     {
229         forAll(T, i)
230         {
231             scalar T1 = fLambdaT(band[1]*T[i]);
232             scalar T2 = fLambdaT(band[0]*T[i]);
233             dimensionedScalar fLambdaDelta
234             (
235                 "fLambdaDelta",
236                 dimless,
237                 T1 - T2
238             );
239             Eb()[i] = Eb()[i]*fLambdaDelta.value();
240         }
241         return Eb;
242     }
246 void Foam::radiation::blackBodyEmission::correct
248     const label lambdaI,
249     const Vector2D<scalar>& band
252     bLambda_[lambdaI] = EbDeltaLambdaT(T_, band);
256 // ************************************************************************* //