summaryrefslogtreecommitdiff
path: root/libm/double/planck.c
blob: 834c85dffefd66ee948e24dc9a5b88d6016ab6ca (plain)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
/*							planck.c
 *
 *	Integral of Planck's black body radiation formula
 *
 *
 *
 * SYNOPSIS:
 *
 * double lambda, T, y, plancki();
 *
 * y = plancki( lambda, T );
 *
 *
 *
 * DESCRIPTION:
 *
 *  Evaluates the definite integral, from wavelength 0 to lambda,
 *  of Planck's radiation formula
 *                      -5
 *            c1  lambda
 *     E =  ------------------
 *            c2/(lambda T)
 *           e             - 1
 *
 * Physical constants c1 = 3.7417749e-16 and c2 = 0.01438769 are built in
 * to the function program.  They are scaled to provide a result
 * in watts per square meter.  Argument T represents temperature in degrees
 * Kelvin; lambda is wavelength in meters.
 *
 * The integral is expressed in closed form, in terms of polylogarithms
 * (see polylog.c).
 *
 * The total area under the curve is
 *      (-1/8) (42 zeta(4) - 12 pi^2 zeta(2) + pi^4 ) c1 (T/c2)^4
 *       = (pi^4 / 15)  c1 (T/c2)^4
 *       =  5.6705032e-8 T^4
 * where sigma = 5.6705032e-8 W m^2 K^-4 is the Stefan-Boltzmann constant.
 *
 *
 * ACCURACY:
 *
 * The left tail of the function experiences some relative error
 * amplification in computing the dominant term exp(-c2/(lambda T)).
 * For the right-hand tail see planckc, below.
 *
 *                      Relative error.
 *   The domain refers to lambda T / c2.
 * arithmetic   domain     # trials      peak         rms
 *    IEEE      0.1, 10      50000      7.1e-15     5.4e-16
 *
 */


/*
Cephes Math Library Release 2.8:  July, 1999
Copyright 1999 by Stephen L. Moshier
*/

#include <math.h>
#ifdef ANSIPROT
extern double polylog (int, double);
extern double exp (double);
extern double log1p (double); /* log(1+x) */
extern double expm1 (double); /* exp(x) - 1 */
double planckc(double, double);
double plancki(double, double);
#else
double polylog(), exp(), log1p(), expm1();
double planckc(), plancki();
#endif

/*  NIST value (1999): 2 pi h c^2 = 3.741 7749(22)  10-16 W m2  */
double planck_c1 = 3.7417749e-16;
/*  NIST value (1999):  h c / k  = 0.014 387 69 m K */
double planck_c2 = 0.01438769;


double
plancki(w, T)
  double w, T;
{
  double b, h, y, bw;

  b = T / planck_c2;
  bw = b * w;

  if (bw > 0.59375)
    {
      y = b * b;
      h = y * y;
      /* Right tail.  */
      y = planckc (w, T);
      /* pi^4 / 15  */
      y =  6.493939402266829149096 * planck_c1 * h  -  y;
      return y;
    }

  h = exp(-planck_c2/(w*T));
  y =      6. * polylog (4, h)  * bw;
  y = (y + 6. * polylog (3, h)) * bw;
  y = (y + 3. * polylog (2, h)) * bw;
  y = (y          - log1p (-h)) * bw;
  h = w * w;
  h = h * h;
  y = y * (planck_c1 / h);
  return y;
}

/*							planckc
 *
 *	Complemented Planck radiation integral
 *
 *
 *
 * SYNOPSIS:
 *
 * double lambda, T, y, planckc();
 *
 * y = planckc( lambda, T );
 *
 *
 *
 * DESCRIPTION:
 *
 *  Integral from w to infinity (area under right hand tail)
 *  of Planck's radiation formula.
 *
 *  The program for large lambda uses an asymptotic series in inverse
 *  powers of the wavelength.
 *
 * ACCURACY:
 *
 *                      Relative error.
 *   The domain refers to lambda T / c2.
 * arithmetic   domain     # trials      peak         rms
 *    IEEE      0.6, 10      50000      1.1e-15     2.2e-16
 *
 */

double
planckc (w, T)
     double w;
     double T;
{
  double b, d, p, u, y;

  b = T / planck_c2;
  d = b*w;
  if (d <= 0.59375)
    {
      y =  6.493939402266829149096 * planck_c1 * b*b*b*b;
      return (y - plancki(w,T));
    }
  u = 1.0/d;
  p = u * u;
#if 0
  y = 236364091.*p/365866013534056632601804800000.;
  y = (y - 15458917./475677107995483570176000000.)*p;
  y = (y + 174611./123104841613737984000000.)*p;
  y = (y - 43867./643745871363538944000.)*p;
  y = ((y + 3617./1081289781411840000.)*p - 1./5928123801600.)*p;
  y = ((y + 691./78460462080000.)*p - 1./2075673600.)*p;
  y = ((((y + 1./35481600.)*p - 1.0/544320.)*p + 1.0/6720.)*p -  1./40.)*p;
  y = y + log(d * expm1(u));
  y = y - 5.*u/8. + 1./3.;
#else
  y = -236364091.*p/45733251691757079075225600000.;
  y = (y + 77683./352527500984795136000000.)*p;
  y = (y - 174611./18465726242060697600000.)*p;
  y = (y + 43867./107290978560589824000.)*p;
  y = ((y - 3617./202741834014720000.)*p + 1./1270312243200.)*p;
  y = ((y - 691./19615115520000.)*p + 1./622702080.)*p;
  y = ((((y - 1./13305600.)*p + 1./272160.)*p - 1./5040.)*p + 1./60.)*p;
  y = y - 0.125*u + 1./3.;
#endif
  y = y * planck_c1 * b / (w*w*w);
  return y;
}


/*							planckd
 *
 *	Planck's black body radiation formula
 *
 *
 *
 * SYNOPSIS:
 *
 * double lambda, T, y, planckd();
 *
 * y = planckd( lambda, T );
 *
 *
 *
 * DESCRIPTION:
 *
 *  Evaluates Planck's radiation formula
 *                      -5
 *            c1  lambda
 *     E =  ------------------
 *            c2/(lambda T)
 *           e             - 1
 *
 */

double
planckd(w, T)
  double w, T;
{
   return (planck_c2 / ((w*w*w*w*w) * (exp(planck_c2/(w*T)) - 1.0)));
}


/* Wavelength, w, of maximum radiation at given temperature T.
   c2/wT = constant
   Wein displacement law.
  */
double
planckw(T)
  double T;
{
  return (planck_c2 / (4.96511423174427630 * T));
}