ViewVC Help
View File | Revision Log | Show Annotations | Download File | Root Listing
root/radiance/ray/src/gen/gendaymtx.c
Revision: 2.31
Committed: Sat Dec 28 18:05:14 2019 UTC (4 years, 4 months ago) by greg
Content type: text/plain
Branch: MAIN
Changes since 2.30: +1 -3 lines
Log Message:
Removed redundant include files

File Contents

# User Rev Content
1 greg 2.1 #ifndef lint
2 greg 2.31 static const char RCSid[] = "$Id: gendaymtx.c,v 2.30 2019/11/07 23:15:07 greg Exp $";
3 greg 2.1 #endif
4     /*
5     * gendaymtx.c
6     *
7     * Generate a daylight matrix based on Perez Sky Model.
8     *
9     * Most of this code is borrowed (see copyright below) from Ian Ashdown's
10     * excellent re-implementation of Jean-Jacques Delaunay's gendaylit.c
11     *
12     * Created by Greg Ward on 1/16/13.
13     */
14    
15     /*********************************************************************
16     *
17     * H32_gendaylit.CPP - Perez Sky Model Calculation
18     *
19     * Version: 1.00A
20     *
21     * History: 09/10/01 - Created.
22     * 11/10/08 - Modified for Unix compilation.
23     * 11/10/12 - Fixed conditional __max directive.
24     * 1/11/13 - Tweaks and optimizations (G.Ward)
25     *
26     * Compilers: Microsoft Visual C/C++ Professional V10.0
27     *
28     * Author: Ian Ashdown, P.Eng.
29     * byHeart Consultants Limited
30     * 620 Ballantree Road
31     * West Vancouver, B.C.
32     * Canada V7S 1W3
33     * e-mail: [email protected]
34     *
35     * References: Perez, R., P. Ineichen, R. Seals, J. Michalsky, and R.
36     * Stewart. 1990. ìModeling Daylight Availability and
37     * Irradiance Components from Direct and Global
38     * Irradiance,î Solar Energy 44(5):271-289.
39     *
40     * Perez, R., R. Seals, and J. Michalsky. 1993.
41     * ìAll-Weather Model for Sky Luminance Distribution -
42     * Preliminary Configuration and Validation,î Solar Energy
43     * 50(3):235-245.
44     *
45     * Perez, R., R. Seals, and J. Michalsky. 1993. "ERRATUM to
46     * All-Weather Model for Sky Luminance Distribution -
47     * Preliminary Configuration and Validation,î Solar Energy
48     * 51(5):423.
49     *
50     * NOTE: This program is a completely rewritten version of
51     * gendaylit.c written by Jean-Jacques Delaunay (1994).
52     *
53     * Copyright 2009-2012 byHeart Consultants Limited. All rights
54     * reserved.
55     *
56     * Redistribution and use in source and binary forms, with or without
57     * modification, are permitted for personal and commercial purposes
58     * provided that redistribution of source code must retain the above
59     * copyright notice, this list of conditions and the following
60     * disclaimer:
61     *
62     * THIS SOFTWARE IS PROVIDED "AS IS" AND ANY EXPRESSED OR IMPLIED
63     * WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED WARRANTIES
64     * OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE
65     * DISCLAIMED. IN NO EVENT SHALL byHeart Consultants Limited OR
66     * ITS EMPLOYEES BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL,
67     * SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT
68     * LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF
69     * USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED
70     * AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
71     * LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
72     * ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
73     * POSSIBILITY OF SUCH DAMAGE.
74     *
75     *********************************************************************/
76    
77     /* Zenith is along the Z-axis */
78     /* X-axis points east */
79     /* Y-axis points north */
80     /* azimuth is measured as degrees or radians east of North */
81    
82     /* Include files */
83     #define _USE_MATH_DEFINES
84     #include <stdlib.h>
85     #include <ctype.h>
86 greg 2.30 #include "platform.h"
87 greg 2.1 #include "rtmath.h"
88 greg 2.23 #include "rtio.h"
89 greg 2.1 #include "color.h"
90 greg 2.30 #include "sun.h"
91 greg 2.1
92     char *progname; /* Program name */
93     char errmsg[128]; /* Error message buffer */
94     const double DC_SolarConstantE = 1367.0; /* Solar constant W/m^2 */
95     const double DC_SolarConstantL = 127.5; /* Solar constant klux */
96    
97     double altitude; /* Solar altitude (radians) */
98     double azimuth; /* Solar azimuth (radians) */
99     double apwc; /* Atmospheric precipitable water content */
100     double dew_point = 11.0; /* Surface dew point temperature (deg. C) */
101     double diff_illum; /* Diffuse illuminance */
102     double diff_irrad; /* Diffuse irradiance */
103     double dir_illum; /* Direct illuminance */
104     double dir_irrad; /* Direct irradiance */
105     int julian_date; /* Julian date */
106     double perez_param[5]; /* Perez sky model parameters */
107     double sky_brightness; /* Sky brightness */
108     double sky_clearness; /* Sky clearness */
109     double solar_rad; /* Solar radiance */
110     double sun_zenith; /* Sun zenith angle (radians) */
111     int input = 0; /* Input type */
112 greg 2.11 int output = 0; /* Output type */
113 greg 2.1
114     extern double dmax( double, double );
115     extern double CalcAirMass();
116     extern double CalcDiffuseIllumRatio( int );
117     extern double CalcDiffuseIrradiance();
118     extern double CalcDirectIllumRatio( int );
119     extern double CalcDirectIrradiance();
120     extern double CalcEccentricity();
121     extern double CalcPrecipWater( double );
122     extern double CalcRelHorzIllum( float *parr );
123     extern double CalcRelLuminance( double, double );
124     extern double CalcSkyBrightness();
125     extern double CalcSkyClearness();
126     extern int CalcSkyParamFromIllum();
127     extern int GetCategoryIndex();
128     extern void CalcPerezParam( double, double, double, int );
129     extern void CalcSkyPatchLumin( float *parr );
130     extern void ComputeSky( float *parr );
131    
132     /* Degrees into radians */
133     #define DegToRad(deg) ((deg)*(PI/180.))
134    
135     /* Radiuans into degrees */
136     #define RadToDeg(rad) ((rad)*(180./PI))
137    
138    
139     /* Perez sky model coefficients */
140    
141     /* Reference: Perez, R., R. Seals, and J. Michalsky, 1993. "All- */
142     /* Weather Model for Sky Luminance Distribution - */
143     /* Preliminary Configuration and Validation," Solar */
144     /* Energy 50(3):235-245, Table 1. */
145    
146     static const double PerezCoeff[8][20] =
147     {
148     /* Sky clearness (epsilon): 1.000 to 1.065 */
149     { 1.3525, -0.2576, -0.2690, -1.4366, -0.7670,
150     0.0007, 1.2734, -0.1233, 2.8000, 0.6004,
151     1.2375, 1.0000, 1.8734, 0.6297, 0.9738,
152     0.2809, 0.0356, -0.1246, -0.5718, 0.9938 },
153     /* Sky clearness (epsilon): 1.065 to 1.230 */
154     { -1.2219, -0.7730, 1.4148, 1.1016, -0.2054,
155     0.0367, -3.9128, 0.9156, 6.9750, 0.1774,
156     6.4477, -0.1239, -1.5798, -0.5081, -1.7812,
157     0.1080, 0.2624, 0.0672, -0.2190, -0.4285 },
158     /* Sky clearness (epsilon): 1.230 to 1.500 */
159     { -1.1000, -0.2515, 0.8952, 0.0156, 0.2782,
160     -0.1812, - 4.5000, 1.1766, 24.7219, -13.0812,
161     -37.7000, 34.8438, -5.0000, 1.5218, 3.9229,
162     -2.6204, -0.0156, 0.1597, 0.4199, -0.5562 },
163     /* Sky clearness (epsilon): 1.500 to 1.950 */
164     { -0.5484, -0.6654, -0.2672, 0.7117, 0.7234,
165     -0.6219, -5.6812, 2.6297, 33.3389, -18.3000,
166     -62.2500, 52.0781, -3.5000, 0.0016, 1.1477,
167     0.1062, 0.4659, -0.3296, -0.0876, -0.0329 },
168     /* Sky clearness (epsilon): 1.950 to 2.800 */
169     { -0.6000, -0.3566, -2.5000, 2.3250, 0.2937,
170     0.0496, -5.6812, 1.8415, 21.0000, -4.7656 ,
171     -21.5906, 7.2492, -3.5000, -0.1554, 1.4062,
172     0.3988, 0.0032, 0.0766, -0.0656, -0.1294 },
173     /* Sky clearness (epsilon): 2.800 to 4.500 */
174     { -1.0156, -0.3670, 1.0078, 1.4051, 0.2875,
175     -0.5328, -3.8500, 3.3750, 14.0000, -0.9999,
176     -7.1406, 7.5469, -3.4000, -0.1078, -1.0750,
177     1.5702, -0.0672, 0.4016, 0.3017, -0.4844 },
178     /* Sky clearness (epsilon): 4.500 to 6.200 */
179     { -1.0000, 0.0211, 0.5025, -0.5119, -0.3000,
180     0.1922, 0.7023, -1.6317, 19.0000, -5.0000,
181     1.2438, -1.9094, -4.0000, 0.0250, 0.3844,
182     0.2656, 1.0468, -0.3788, -2.4517, 1.4656 },
183     /* Sky clearness (epsilon): 6.200 to ... */
184     { -1.0500, 0.0289, 0.4260, 0.3590, -0.3250,
185     0.1156, 0.7781, 0.0025, 31.0625, -14.5000,
186     -46.1148, 55.3750, -7.2312, 0.4050, 13.3500,
187     0.6234, 1.5000, -0.6426, 1.8564, 0.5636 }
188     };
189    
190     /* Perez irradiance component model coefficients */
191    
192     /* Reference: Perez, R., P. Ineichen, R. Seals, J. Michalsky, and R. */
193     /* Stewart. 1990. ìModeling Daylight Availability and */
194     /* Irradiance Components from Direct and Global */
195     /* Irradiance,î Solar Energy 44(5):271-289. */
196    
197     typedef struct
198     {
199     double lower; /* Lower bound */
200     double upper; /* Upper bound */
201     } CategoryBounds;
202    
203     /* Perez sky clearness (epsilon) categories (Table 1) */
204     static const CategoryBounds SkyClearCat[8] =
205     {
206     { 1.000, 1.065 }, /* Overcast */
207     { 1.065, 1.230 },
208     { 1.230, 1.500 },
209     { 1.500, 1.950 },
210     { 1.950, 2.800 },
211     { 2.800, 4.500 },
212     { 4.500, 6.200 },
213 greg 2.12 { 6.200, 12.01 } /* Clear */
214 greg 2.1 };
215    
216     /* Luminous efficacy model coefficients */
217     typedef struct
218     {
219     double a;
220     double b;
221     double c;
222     double d;
223     } ModelCoeff;
224    
225     /* Diffuse luminous efficacy model coefficients (Table 4, Eqn. 7) */
226     static const ModelCoeff DiffuseLumEff[8] =
227     {
228     { 97.24, -0.46, 12.00, -8.91 },
229     { 107.22, 1.15, 0.59, -3.95 },
230     { 104.97, 2.96, -5.53, -8.77 },
231     { 102.39, 5.59, -13.95, -13.90 },
232     { 100.71, 5.94, -22.75, -23.74 },
233     { 106.42, 3.83, -36.15, -28.83 },
234     { 141.88, 1.90, -53.24, -14.03 },
235     { 152.23, 0.35, -45.27, -7.98 }
236     };
237    
238     /* Direct luminous efficacy model coefficients (Table 4, Eqn. 8) */
239     static const ModelCoeff DirectLumEff[8] =
240     {
241     { 57.20, -4.55, -2.98, 117.12 },
242     { 98.99, -3.46, -1.21, 12.38 },
243     { 109.83, -4.90, -1.71, -8.81 },
244     { 110.34, -5.84, -1.99, -4.56 },
245     { 106.36, -3.97, -1.75, -6.16 },
246     { 107.19, -1.25, -1.51, -26.73 },
247     { 105.75, 0.77, -1.26, -34.44 },
248     { 101.18, 1.58, -1.10, -8.29 }
249     };
250    
251 greg 2.3 #ifndef NSUNPATCH
252 greg 2.10 #define NSUNPATCH 4 /* max. # patches to spread sun into */
253 greg 2.3 #endif
254    
255 greg 2.10 int nsuns = NSUNPATCH; /* number of sun patches to use */
256     double fixed_sun_sa = -1; /* fixed solid angle per sun? */
257    
258 greg 2.1 int verbose = 0; /* progress reports to stderr? */
259    
260     int outfmt = 'a'; /* output format */
261    
262     int rhsubdiv = 1; /* Reinhart sky subdivisions */
263    
264 greg 2.4 COLOR skycolor = {.96, 1.004, 1.118}; /* sky coloration */
265     COLOR suncolor = {1., 1., 1.}; /* sun color */
266     COLOR grefl = {.2, .2, .2}; /* ground reflectance */
267 greg 2.1
268     int nskypatch; /* number of Reinhart patches */
269     float *rh_palt; /* sky patch altitudes (radians) */
270     float *rh_pazi; /* sky patch azimuths (radians) */
271     float *rh_dom; /* sky patch solid angle (sr) */
272    
273     #define vector(v,alt,azi) ( (v)[1] = tcos(alt), \
274     (v)[0] = (v)[1]*tsin(azi), \
275     (v)[1] *= tcos(azi), \
276     (v)[2] = tsin(alt) )
277    
278     #define rh_vector(v,i) vector(v,rh_palt[i],rh_pazi[i])
279    
280     #define rh_cos(i) tsin(rh_palt[i])
281    
282     extern int rh_init(void);
283     extern float * resize_dmatrix(float *mtx_data, int nsteps, int npatch);
284     extern void AddDirect(float *parr);
285    
286 greg 2.14
287     static const char *
288     getfmtname(int fmt)
289     {
290     switch (fmt) {
291     case 'a':
292     return("ascii");
293     case 'f':
294     return("float");
295     case 'd':
296     return("double");
297     }
298     return("unknown");
299     }
300    
301    
302 greg 2.1 int
303     main(int argc, char *argv[])
304     {
305     char buf[256];
306 greg 2.14 int doheader = 1; /* output header? */
307 greg 2.8 double rotation = 0; /* site rotation (degrees) */
308 greg 2.1 double elevation; /* site elevation (meters) */
309     int dir_is_horiz; /* direct is meas. on horizontal? */
310     float *mtx_data = NULL; /* our matrix data */
311 greg 2.28 int avgSky = 0; /* compute average sky r.t. matrix? */
312 greg 2.27 int ntsteps = 0; /* number of time steps */
313 greg 2.28 int tstorage = 0; /* number of allocated time steps */
314     int nstored = 0; /* number of time steps in matrix */
315 greg 2.1 int last_monthly = 0; /* month of last report */
316 greg 2.22 int inconsistent = 0; /* inconsistent options set? */
317 greg 2.1 int mo, da; /* month (1-12) and day (1-31) */
318     double hr; /* hour (local standard time) */
319     double dir, dif; /* direct and diffuse values */
320     int mtx_offset;
321     int i, j;
322    
323     progname = argv[0];
324     /* get options */
325     for (i = 1; i < argc && argv[i][0] == '-'; i++)
326     switch (argv[i][1]) {
327 greg 2.4 case 'g': /* ground reflectance */
328     grefl[0] = atof(argv[++i]);
329     grefl[1] = atof(argv[++i]);
330     grefl[2] = atof(argv[++i]);
331 greg 2.1 break;
332 greg 2.4 case 'v': /* verbose progress reports */
333 greg 2.1 verbose++;
334     break;
335 greg 2.14 case 'h': /* turn off header */
336     doheader = 0;
337     break;
338 greg 2.4 case 'o': /* output format */
339 greg 2.1 switch (argv[i][2]) {
340     case 'f':
341     case 'd':
342     case 'a':
343     outfmt = argv[i][2];
344     break;
345     default:
346     goto userr;
347     }
348     break;
349 greg 2.11 case 'O': /* output type */
350     switch (argv[i][2]) {
351     case '0':
352     output = 0;
353     break;
354     case '1':
355     output = 1;
356     break;
357     default:
358     goto userr;
359     }
360     if (argv[i][3])
361     goto userr;
362     break;
363 greg 2.4 case 'm': /* Reinhart subdivisions */
364 greg 2.1 rhsubdiv = atoi(argv[++i]);
365     break;
366 greg 2.4 case 'c': /* sky color */
367 greg 2.22 inconsistent |= (skycolor[1] <= 1e-4);
368 greg 2.1 skycolor[0] = atof(argv[++i]);
369     skycolor[1] = atof(argv[++i]);
370     skycolor[2] = atof(argv[++i]);
371     break;
372 greg 2.4 case 'd': /* solar (direct) only */
373 greg 2.1 skycolor[0] = skycolor[1] = skycolor[2] = 0;
374 greg 2.22 if (suncolor[1] <= 1e-4) {
375     inconsistent = 1;
376 greg 2.4 suncolor[0] = suncolor[1] = suncolor[2] = 1;
377 greg 2.22 }
378 greg 2.1 break;
379 greg 2.4 case 's': /* sky only (no direct) */
380     suncolor[0] = suncolor[1] = suncolor[2] = 0;
381 greg 2.22 if (skycolor[1] <= 1e-4) {
382     inconsistent = 1;
383 greg 2.1 skycolor[0] = skycolor[1] = skycolor[2] = 1;
384 greg 2.22 }
385 greg 2.1 break;
386 greg 2.8 case 'r': /* rotate distribution */
387     if (argv[i][2] && argv[i][2] != 'z')
388     goto userr;
389     rotation = atof(argv[++i]);
390     break;
391 greg 2.10 case '5': /* 5-phase calculation */
392     nsuns = 1;
393 greg 2.19 fixed_sun_sa = PI/360.*atof(argv[++i]);
394 greg 2.21 if (fixed_sun_sa <= 0) {
395     fprintf(stderr, "%s: missing solar disk size argument for '-5' option\n",
396     argv[0]);
397     exit(1);
398     }
399 greg 2.19 fixed_sun_sa *= fixed_sun_sa*PI;
400 greg 2.10 break;
401 greg 2.27 case 'A': /* compute average sky */
402     avgSky = 1;
403     break;
404 greg 2.1 default:
405     goto userr;
406     }
407     if (i < argc-1)
408     goto userr;
409 greg 2.22 if (inconsistent)
410     fprintf(stderr, "%s: WARNING: inconsistent -s, -d, -c options!\n",
411     progname);
412 greg 2.1 if (i == argc-1 && freopen(argv[i], "r", stdin) == NULL) {
413     fprintf(stderr, "%s: cannot open '%s' for input\n",
414     progname, argv[i]);
415     exit(1);
416     }
417     if (verbose) {
418     if (i == argc-1)
419     fprintf(stderr, "%s: reading weather tape '%s'\n",
420     progname, argv[i]);
421     else
422     fprintf(stderr, "%s: reading weather tape from <stdin>\n",
423     progname);
424     }
425     /* read weather tape header */
426 greg 2.2 if (scanf("place %[^\r\n] ", buf) != 1)
427 greg 2.1 goto fmterr;
428     if (scanf("latitude %lf\n", &s_latitude) != 1)
429     goto fmterr;
430     if (scanf("longitude %lf\n", &s_longitude) != 1)
431     goto fmterr;
432     if (scanf("time_zone %lf\n", &s_meridian) != 1)
433     goto fmterr;
434     if (scanf("site_elevation %lf\n", &elevation) != 1)
435     goto fmterr;
436     if (scanf("weather_data_file_units %d\n", &input) != 1)
437     goto fmterr;
438     switch (input) { /* translate units */
439     case 1:
440     input = 1; /* radiometric quantities */
441     dir_is_horiz = 0; /* direct is perpendicular meas. */
442     break;
443     case 2:
444     input = 1; /* radiometric quantities */
445     dir_is_horiz = 1; /* solar measured horizontally */
446     break;
447     case 3:
448     input = 2; /* photometric quantities */
449     dir_is_horiz = 0; /* direct is perpendicular meas. */
450     break;
451     default:
452     goto fmterr;
453     }
454     rh_init(); /* initialize sky patches */
455     if (verbose) {
456     fprintf(stderr, "%s: location '%s'\n", progname, buf);
457     fprintf(stderr, "%s: (lat,long)=(%.1f,%.1f) degrees north, west\n",
458     progname, s_latitude, s_longitude);
459     fprintf(stderr, "%s: %d sky patches per time step\n",
460     progname, nskypatch);
461 greg 2.8 if (rotation != 0)
462     fprintf(stderr, "%s: rotating output %.0f degrees\n",
463     progname, rotation);
464 greg 2.1 }
465 greg 2.2 /* convert quantities to radians */
466     s_latitude = DegToRad(s_latitude);
467     s_longitude = DegToRad(s_longitude);
468     s_meridian = DegToRad(s_meridian);
469 greg 2.27 /* initial allocation */
470 greg 2.28 mtx_data = resize_dmatrix(mtx_data, tstorage=2, nskypatch);
471 greg 2.1 /* process each time step in tape */
472     while (scanf("%d %d %lf %lf %lf\n", &mo, &da, &hr, &dir, &dif) == 5) {
473     double sda, sta;
474 greg 2.27
475 greg 2.28 mtx_offset = 3*nskypatch*nstored;
476     nstored += !avgSky | !nstored;
477 greg 2.27 /* make space for next row */
478 greg 2.28 if (nstored > tstorage) {
479     tstorage += (tstorage>>1) + nstored + 7;
480     mtx_data = resize_dmatrix(mtx_data, tstorage, nskypatch);
481 greg 2.15 }
482 greg 2.27 ntsteps++; /* keep count of time steps */
483 greg 2.1 if (dif <= 1e-4) {
484 greg 2.28 if (!avgSky | !mtx_offset)
485 greg 2.27 memset(mtx_data+mtx_offset, 0, sizeof(float)*3*nskypatch);
486 greg 2.1 continue;
487     }
488     if (verbose && mo != last_monthly)
489     fprintf(stderr, "%s: stepping through month %d...\n",
490     progname, last_monthly=mo);
491     /* compute solar position */
492     julian_date = jdate(mo, da);
493     sda = sdec(julian_date);
494     sta = stadj(julian_date);
495     altitude = salt(sda, hr+sta);
496 greg 2.8 azimuth = sazi(sda, hr+sta) + PI - DegToRad(rotation);
497 greg 2.1 /* convert measured values */
498     if (dir_is_horiz && altitude > 0.)
499     dir /= sin(altitude);
500     if (input == 1) {
501     dir_irrad = dir;
502     diff_irrad = dif;
503     } else /* input == 2 */ {
504     dir_illum = dir;
505     diff_illum = dif;
506     }
507     /* compute sky patch values */
508     ComputeSky(mtx_data+mtx_offset);
509 greg 2.4 AddDirect(mtx_data+mtx_offset);
510 greg 2.27 /* update cumulative sky? */
511     for (i = 3*nskypatch*(avgSky&(ntsteps>1)); i--; )
512     mtx_data[i] += mtx_data[mtx_offset+i];
513 greg 2.1 }
514     /* check for junk at end */
515     while ((i = fgetc(stdin)) != EOF)
516     if (!isspace(i)) {
517     fprintf(stderr, "%s: warning - unexpected data past EOT: ",
518     progname);
519     buf[0] = i; buf[1] = '\0';
520     fgets(buf+1, sizeof(buf)-1, stdin);
521     fputs(buf, stderr); fputc('\n', stderr);
522     break;
523     }
524 greg 2.27 if (!ntsteps) {
525     fprintf(stderr, "%s: no valid time steps on input\n", progname);
526     exit(1);
527     }
528     dif = 1./(double)ntsteps; /* average sky? */
529     for (i = 3*nskypatch*(avgSky&(ntsteps>1)); i--; )
530     mtx_data[i] *= dif;
531 greg 2.1 /* write out matrix */
532 greg 2.14 if (outfmt != 'a')
533     SET_FILE_BINARY(stdout);
534 greg 2.1 #ifdef getc_unlocked
535     flockfile(stdout);
536     #endif
537     if (verbose)
538     fprintf(stderr, "%s: writing %smatrix with %d time steps...\n",
539 greg 2.28 progname, outfmt=='a' ? "" : "binary ", nstored);
540 greg 2.14 if (doheader) {
541     newheader("RADIANCE", stdout);
542     printargs(argc, argv, stdout);
543     printf("LATLONG= %.8f %.8f\n", RadToDeg(s_latitude),
544     -RadToDeg(s_longitude));
545     printf("NROWS=%d\n", nskypatch);
546 greg 2.28 printf("NCOLS=%d\n", nstored);
547 greg 2.14 printf("NCOMP=3\n");
548 greg 2.29 if ((outfmt == 'f') | (outfmt == 'd'))
549     fputendian(stdout);
550 greg 2.18 fputformat((char *)getfmtname(outfmt), stdout);
551 greg 2.14 putchar('\n');
552     }
553 greg 2.1 /* patches are rows (outer sort) */
554     for (i = 0; i < nskypatch; i++) {
555     mtx_offset = 3*i;
556     switch (outfmt) {
557     case 'a':
558 greg 2.28 for (j = 0; j < nstored; j++) {
559 greg 2.3 printf("%.3g %.3g %.3g\n", mtx_data[mtx_offset],
560 greg 2.1 mtx_data[mtx_offset+1],
561     mtx_data[mtx_offset+2]);
562     mtx_offset += 3*nskypatch;
563     }
564 greg 2.28 if (nstored > 1)
565 greg 2.2 fputc('\n', stdout);
566 greg 2.1 break;
567     case 'f':
568 greg 2.28 for (j = 0; j < nstored; j++) {
569 greg 2.23 putbinary(mtx_data+mtx_offset, sizeof(float), 3,
570 greg 2.1 stdout);
571     mtx_offset += 3*nskypatch;
572     }
573     break;
574     case 'd':
575 greg 2.28 for (j = 0; j < nstored; j++) {
576 greg 2.1 double ment[3];
577     ment[0] = mtx_data[mtx_offset];
578     ment[1] = mtx_data[mtx_offset+1];
579     ment[2] = mtx_data[mtx_offset+2];
580 greg 2.23 putbinary(ment, sizeof(double), 3, stdout);
581 greg 2.1 mtx_offset += 3*nskypatch;
582     }
583     break;
584     }
585     if (ferror(stdout))
586     goto writerr;
587     }
588     if (fflush(stdout) == EOF)
589     goto writerr;
590     if (verbose)
591     fprintf(stderr, "%s: done.\n", progname);
592     exit(0);
593     userr:
594 greg 2.27 fprintf(stderr, "Usage: %s [-v][-h][-A][-d|-s][-r deg][-m N][-g r g b][-c r g b][-o{f|d}][-O{0|1}] [tape.wea]\n",
595 greg 2.1 progname);
596     exit(1);
597     fmterr:
598     fprintf(stderr, "%s: input weather tape format error\n", progname);
599     exit(1);
600     writerr:
601     fprintf(stderr, "%s: write error on output\n", progname);
602     exit(1);
603     }
604    
605     /* Return maximum of two doubles */
606     double dmax( double a, double b )
607     { return (a > b) ? a : b; }
608    
609     /* Compute sky patch radiance values (modified by GW) */
610     void
611     ComputeSky(float *parr)
612     {
613     int index; /* Category index */
614     double norm_diff_illum; /* Normalized diffuse illuimnance */
615     int i;
616    
617     /* Calculate atmospheric precipitable water content */
618     apwc = CalcPrecipWater(dew_point);
619    
620 greg 2.6 /* Calculate sun zenith angle (don't let it dip below horizon) */
621     /* Also limit minimum angle to keep circumsolar off zenith */
622     if (altitude <= 0.0)
623     sun_zenith = DegToRad(90.0);
624     else if (altitude >= DegToRad(87.0))
625     sun_zenith = DegToRad(3.0);
626     else
627     sun_zenith = DegToRad(90.0) - altitude;
628 greg 2.1
629     /* Compute the inputs for the calculation of the sky distribution */
630    
631     if (input == 0) /* XXX never used */
632     {
633     /* Calculate irradiance */
634     diff_irrad = CalcDiffuseIrradiance();
635     dir_irrad = CalcDirectIrradiance();
636    
637     /* Calculate illuminance */
638     index = GetCategoryIndex();
639     diff_illum = diff_irrad * CalcDiffuseIllumRatio(index);
640     dir_illum = dir_irrad * CalcDirectIllumRatio(index);
641     }
642     else if (input == 1)
643     {
644     sky_brightness = CalcSkyBrightness();
645     sky_clearness = CalcSkyClearness();
646    
647 greg 2.9 /* Limit sky clearness */
648     if (sky_clearness > 11.9)
649     sky_clearness = 11.9;
650    
651     /* Limit sky brightness */
652     if (sky_brightness < 0.01)
653 greg 2.11 sky_brightness = 0.01;
654 greg 2.9
655 greg 2.1 /* Calculate illuminance */
656     index = GetCategoryIndex();
657     diff_illum = diff_irrad * CalcDiffuseIllumRatio(index);
658     dir_illum = dir_irrad * CalcDirectIllumRatio(index);
659     }
660     else if (input == 2)
661     {
662     /* Calculate sky brightness and clearness from illuminance values */
663     index = CalcSkyParamFromIllum();
664     }
665    
666 greg 2.11 if (output == 1) { /* hack for solar radiance */
667     diff_illum = diff_irrad * WHTEFFICACY;
668     dir_illum = dir_irrad * WHTEFFICACY;
669     }
670    
671 greg 2.2 if (bright(skycolor) <= 1e-4) { /* 0 sky component? */
672     memset(parr, 0, sizeof(float)*3*nskypatch);
673     return;
674     }
675 greg 2.1 /* Compute ground radiance (include solar contribution if any) */
676 greg 2.3 parr[0] = diff_illum;
677 greg 2.1 if (altitude > 0)
678 greg 2.3 parr[0] += dir_illum * sin(altitude);
679 greg 2.4 parr[2] = parr[1] = parr[0] *= (1./PI/WHTEFFICACY);
680     multcolor(parr, grefl);
681 greg 2.1
682     /* Calculate Perez sky model parameters */
683     CalcPerezParam(sun_zenith, sky_clearness, sky_brightness, index);
684    
685     /* Calculate sky patch luminance values */
686     CalcSkyPatchLumin(parr);
687    
688     /* Calculate relative horizontal illuminance */
689     norm_diff_illum = CalcRelHorzIllum(parr);
690    
691 greg 2.13 /* Check for zero sky -- make uniform in that case */
692     if (norm_diff_illum <= FTINY) {
693     for (i = 1; i < nskypatch; i++)
694     setcolor(parr+3*i, 1., 1., 1.);
695     norm_diff_illum = PI;
696     }
697 greg 2.1 /* Normalization coefficient */
698     norm_diff_illum = diff_illum / norm_diff_illum;
699    
700     /* Apply to sky patches to get absolute radiance values */
701     for (i = 1; i < nskypatch; i++) {
702 greg 2.7 scalecolor(parr+3*i, norm_diff_illum*(1./WHTEFFICACY));
703 greg 2.1 multcolor(parr+3*i, skycolor);
704     }
705     }
706    
707     /* Add in solar direct to nearest sky patches (GW) */
708     void
709     AddDirect(float *parr)
710     {
711     FVECT svec;
712 greg 2.3 double near_dprod[NSUNPATCH];
713     int near_patch[NSUNPATCH];
714     double wta[NSUNPATCH], wtot;
715 greg 2.1 int i, j, p;
716    
717 greg 2.4 if (dir_illum <= 1e-4 || bright(suncolor) <= 1e-4)
718 greg 2.1 return;
719 greg 2.10 /* identify nsuns closest patches */
720     if (nsuns > NSUNPATCH)
721     nsuns = NSUNPATCH;
722     else if (nsuns <= 0)
723     nsuns = 1;
724     for (i = nsuns; i--; )
725 greg 2.1 near_dprod[i] = -1.;
726     vector(svec, altitude, azimuth);
727     for (p = 1; p < nskypatch; p++) {
728     FVECT pvec;
729     double dprod;
730     rh_vector(pvec, p);
731     dprod = DOT(pvec, svec);
732 greg 2.10 for (i = 0; i < nsuns; i++)
733 greg 2.1 if (dprod > near_dprod[i]) {
734 greg 2.10 for (j = nsuns; --j > i; ) {
735 greg 2.1 near_dprod[j] = near_dprod[j-1];
736     near_patch[j] = near_patch[j-1];
737     }
738     near_dprod[i] = dprod;
739     near_patch[i] = p;
740     break;
741     }
742     }
743     wtot = 0; /* weight by proximity */
744 greg 2.10 for (i = nsuns; i--; )
745 greg 2.1 wtot += wta[i] = 1./(1.002 - near_dprod[i]);
746     /* add to nearest patch radiances */
747 greg 2.10 for (i = nsuns; i--; ) {
748 greg 2.2 float *pdest = parr + 3*near_patch[i];
749 greg 2.10 float val_add = wta[i] * dir_illum / (WHTEFFICACY * wtot);
750    
751     val_add /= (fixed_sun_sa > 0) ? fixed_sun_sa
752     : rh_dom[near_patch[i]] ;
753 greg 2.4 *pdest++ += val_add*suncolor[0];
754     *pdest++ += val_add*suncolor[1];
755     *pdest++ += val_add*suncolor[2];
756 greg 2.2 }
757 greg 2.1 }
758    
759     /* Initialize Reinhart sky patch positions (GW) */
760     int
761     rh_init(void)
762     {
763     #define NROW 7
764     static const int tnaz[NROW] = {30, 30, 24, 24, 18, 12, 6};
765     const double alpha = (PI/2.)/(NROW*rhsubdiv + .5);
766     int p, i, j;
767     /* allocate patch angle arrays */
768     nskypatch = 0;
769     for (p = 0; p < NROW; p++)
770     nskypatch += tnaz[p];
771     nskypatch *= rhsubdiv*rhsubdiv;
772     nskypatch += 2;
773     rh_palt = (float *)malloc(sizeof(float)*nskypatch);
774     rh_pazi = (float *)malloc(sizeof(float)*nskypatch);
775     rh_dom = (float *)malloc(sizeof(float)*nskypatch);
776     if ((rh_palt == NULL) | (rh_pazi == NULL) | (rh_dom == NULL)) {
777     fprintf(stderr, "%s: out of memory in rh_init()\n", progname);
778     exit(1);
779     }
780     rh_palt[0] = -PI/2.; /* ground & zenith patches */
781     rh_pazi[0] = 0.;
782     rh_dom[0] = 2.*PI;
783     rh_palt[nskypatch-1] = PI/2.;
784     rh_pazi[nskypatch-1] = 0.;
785     rh_dom[nskypatch-1] = 2.*PI*(1. - cos(alpha*.5));
786     p = 1; /* "normal" patches */
787     for (i = 0; i < NROW*rhsubdiv; i++) {
788     const float ralt = alpha*(i + .5);
789     const int ninrow = tnaz[i/rhsubdiv]*rhsubdiv;
790 greg 2.3 const float dom = 2.*PI*(sin(alpha*(i+1)) - sin(alpha*i)) /
791     (double)ninrow;
792 greg 2.1 for (j = 0; j < ninrow; j++) {
793     rh_palt[p] = ralt;
794     rh_pazi[p] = 2.*PI * j / (double)ninrow;
795     rh_dom[p++] = dom;
796     }
797     }
798     return nskypatch;
799     #undef NROW
800     }
801    
802     /* Resize daylight matrix (GW) */
803     float *
804     resize_dmatrix(float *mtx_data, int nsteps, int npatch)
805     {
806     if (mtx_data == NULL)
807     mtx_data = (float *)malloc(sizeof(float)*3*nsteps*npatch);
808     else
809     mtx_data = (float *)realloc(mtx_data,
810     sizeof(float)*3*nsteps*npatch);
811     if (mtx_data == NULL) {
812     fprintf(stderr, "%s: out of memory in resize_dmatrix(%d,%d)\n",
813     progname, nsteps, npatch);
814     exit(1);
815     }
816     return(mtx_data);
817     }
818    
819     /* Determine category index */
820     int GetCategoryIndex()
821     {
822     int index; /* Loop index */
823    
824     for (index = 0; index < 8; index++)
825     if ((sky_clearness >= SkyClearCat[index].lower) &&
826     (sky_clearness < SkyClearCat[index].upper))
827     break;
828    
829     return index;
830     }
831    
832     /* Calculate diffuse illuminance to diffuse irradiance ratio */
833    
834     /* Reference: Perez, R., P. Ineichen, R. Seals, J. Michalsky, and R. */
835     /* Stewart. 1990. ìModeling Daylight Availability and */
836     /* Irradiance Components from Direct and Global */
837     /* Irradiance,î Solar Energy 44(5):271-289, Eqn. 7. */
838    
839     double CalcDiffuseIllumRatio( int index )
840     {
841     ModelCoeff const *pnle; /* Category coefficient pointer */
842    
843     /* Get category coefficient pointer */
844     pnle = &(DiffuseLumEff[index]);
845    
846     return pnle->a + pnle->b * apwc + pnle->c * cos(sun_zenith) +
847     pnle->d * log(sky_brightness);
848     }
849    
850     /* Calculate direct illuminance to direct irradiance ratio */
851    
852     /* Reference: Perez, R., P. Ineichen, R. Seals, J. Michalsky, and R. */
853     /* Stewart. 1990. ìModeling Daylight Availability and */
854     /* Irradiance Components from Direct and Global */
855     /* Irradiance,î Solar Energy 44(5):271-289, Eqn. 8. */
856    
857     double CalcDirectIllumRatio( int index )
858     {
859     ModelCoeff const *pnle; /* Category coefficient pointer */
860    
861     /* Get category coefficient pointer */
862     pnle = &(DirectLumEff[index]);
863    
864     /* Calculate direct illuminance from direct irradiance */
865    
866     return dmax((pnle->a + pnle->b * apwc + pnle->c * exp(5.73 *
867     sun_zenith - 5.0) + pnle->d * sky_brightness),
868     0.0);
869     }
870    
871     /* Calculate sky brightness */
872    
873     /* Reference: Perez, R., P. Ineichen, R. Seals, J. Michalsky, and R. */
874     /* Stewart. 1990. ìModeling Daylight Availability and */
875     /* Irradiance Components from Direct and Global */
876     /* Irradiance,î Solar Energy 44(5):271-289, Eqn. 2. */
877    
878     double CalcSkyBrightness()
879     {
880     return diff_irrad * CalcAirMass() / (DC_SolarConstantE *
881     CalcEccentricity());
882     }
883    
884     /* Calculate sky clearness */
885    
886     /* Reference: Perez, R., P. Ineichen, R. Seals, J. Michalsky, and R. */
887     /* Stewart. 1990. ìModeling Daylight Availability and */
888     /* Irradiance Components from Direct and Global */
889     /* Irradiance,î Solar Energy 44(5):271-289, Eqn. 1. */
890    
891     double CalcSkyClearness()
892     {
893     double sz_cubed; /* Sun zenith angle cubed */
894    
895     /* Calculate sun zenith angle cubed */
896 greg 2.11 sz_cubed = sun_zenith*sun_zenith*sun_zenith;
897 greg 2.1
898     return ((diff_irrad + dir_irrad) / diff_irrad + 1.041 *
899     sz_cubed) / (1.0 + 1.041 * sz_cubed);
900     }
901    
902     /* Calculate diffuse horizontal irradiance from Perez sky brightness */
903    
904     /* Reference: Perez, R., P. Ineichen, R. Seals, J. Michalsky, and R. */
905     /* Stewart. 1990. ìModeling Daylight Availability and */
906     /* Irradiance Components from Direct and Global */
907     /* Irradiance,î Solar Energy 44(5):271-289, Eqn. 2 */
908     /* (inverse). */
909    
910     double CalcDiffuseIrradiance()
911     {
912     return sky_brightness * DC_SolarConstantE * CalcEccentricity() /
913     CalcAirMass();
914     }
915    
916     /* Calculate direct normal irradiance from Perez sky clearness */
917    
918     /* Reference: Perez, R., P. Ineichen, R. Seals, J. Michalsky, and R. */
919     /* Stewart. 1990. ìModeling Daylight Availability and */
920     /* Irradiance Components from Direct and Global */
921     /* Irradiance,î Solar Energy 44(5):271-289, Eqn. 1 */
922     /* (inverse). */
923    
924     double CalcDirectIrradiance()
925     {
926     return CalcDiffuseIrradiance() * ((sky_clearness - 1.0) * (1 + 1.041
927 greg 2.11 * sun_zenith*sun_zenith*sun_zenith));
928 greg 2.1 }
929    
930     /* Calculate sky brightness and clearness from illuminance values */
931     int CalcSkyParamFromIllum()
932     {
933     double test1 = 0.1;
934     double test2 = 0.1;
935     int counter = 0;
936     int index = 0; /* Category index */
937    
938     /* Convert illuminance to irradiance */
939     diff_irrad = diff_illum * DC_SolarConstantE /
940     (DC_SolarConstantL * 1000.0);
941     dir_irrad = dir_illum * DC_SolarConstantE /
942     (DC_SolarConstantL * 1000.0);
943    
944     /* Calculate sky brightness and clearness */
945     sky_brightness = CalcSkyBrightness();
946     sky_clearness = CalcSkyClearness();
947    
948     /* Limit sky clearness */
949     if (sky_clearness > 12.0)
950     sky_clearness = 12.0;
951    
952     /* Limit sky brightness */
953 greg 2.9 if (sky_brightness < 0.01)
954 greg 2.1 sky_brightness = 0.01;
955    
956     while (((fabs(diff_irrad - test1) > 10.0) ||
957     (fabs(dir_irrad - test2) > 10.0)) && !(counter == 5))
958     {
959     test1 = diff_irrad;
960     test2 = dir_irrad;
961     counter++;
962    
963     /* Convert illuminance to irradiance */
964     index = GetCategoryIndex();
965     diff_irrad = diff_illum / CalcDiffuseIllumRatio(index);
966 greg 2.26 dir_irrad = CalcDirectIllumRatio(index);
967     if (dir_irrad > 0.1)
968     dir_irrad = dir_illum / dir_irrad;
969 greg 2.1
970     /* Calculate sky brightness and clearness */
971     sky_brightness = CalcSkyBrightness();
972     sky_clearness = CalcSkyClearness();
973    
974     /* Limit sky clearness */
975     if (sky_clearness > 12.0)
976     sky_clearness = 12.0;
977    
978     /* Limit sky brightness */
979 greg 2.9 if (sky_brightness < 0.01)
980 greg 2.1 sky_brightness = 0.01;
981     }
982    
983     return GetCategoryIndex();
984     }
985    
986     /* Calculate relative luminance */
987    
988     /* Reference: Perez, R., R. Seals, and J. Michalsky. 1993. */
989     /* ìAll-Weather Model for Sky Luminance Distribution - */
990     /* Preliminary Configuration and Validation,î Solar Energy */
991     /* 50(3):235-245, Eqn. 1. */
992    
993     double CalcRelLuminance( double gamma, double zeta )
994     {
995     return (1.0 + perez_param[0] * exp(perez_param[1] / cos(zeta))) *
996     (1.0 + perez_param[2] * exp(perez_param[3] * gamma) +
997     perez_param[4] * cos(gamma) * cos(gamma));
998     }
999    
1000     /* Calculate Perez sky model parameters */
1001    
1002     /* Reference: Perez, R., R. Seals, and J. Michalsky. 1993. */
1003     /* ìAll-Weather Model for Sky Luminance Distribution - */
1004     /* Preliminary Configuration and Validation,î Solar Energy */
1005     /* 50(3):235-245, Eqns. 6 - 8. */
1006    
1007     void CalcPerezParam( double sz, double epsilon, double delta,
1008     int index )
1009     {
1010     double x[5][4]; /* Coefficents a, b, c, d, e */
1011     int i, j; /* Loop indices */
1012    
1013     /* Limit sky brightness */
1014     if (epsilon > 1.065 && epsilon < 2.8)
1015     {
1016     if (delta < 0.2)
1017     delta = 0.2;
1018     }
1019    
1020     /* Get Perez coefficients */
1021     for (i = 0; i < 5; i++)
1022     for (j = 0; j < 4; j++)
1023     x[i][j] = PerezCoeff[index][4 * i + j];
1024    
1025     if (index != 0)
1026     {
1027     /* Calculate parameter a, b, c, d and e (Eqn. 6) */
1028     for (i = 0; i < 5; i++)
1029     perez_param[i] = x[i][0] + x[i][1] * sz + delta * (x[i][2] +
1030     x[i][3] * sz);
1031     }
1032     else
1033     {
1034     /* Parameters a, b and e (Eqn. 6) */
1035     perez_param[0] = x[0][0] + x[0][1] * sz + delta * (x[0][2] +
1036     x[0][3] * sz);
1037     perez_param[1] = x[1][0] + x[1][1] * sz + delta * (x[1][2] +
1038     x[1][3] * sz);
1039     perez_param[4] = x[4][0] + x[4][1] * sz + delta * (x[4][2] +
1040     x[4][3] * sz);
1041    
1042     /* Parameter c (Eqn. 7) */
1043     perez_param[2] = exp(pow(delta * (x[2][0] + x[2][1] * sz),
1044     x[2][2])) - x[2][3];
1045    
1046     /* Parameter d (Eqn. 8) */
1047     perez_param[3] = -exp(delta * (x[3][0] + x[3][1] * sz)) +
1048     x[3][2] + delta * x[3][3];
1049     }
1050     }
1051    
1052     /* Calculate relative horizontal illuminance (modified by GW) */
1053    
1054     /* Reference: Perez, R., R. Seals, and J. Michalsky. 1993. */
1055     /* ìAll-Weather Model for Sky Luminance Distribution - */
1056     /* Preliminary Configuration and Validation,î Solar Energy */
1057     /* 50(3):235-245, Eqn. 3. */
1058    
1059     double CalcRelHorzIllum( float *parr )
1060     {
1061     int i;
1062     double rh_illum = 0.0; /* Relative horizontal illuminance */
1063    
1064     for (i = 1; i < nskypatch; i++)
1065 greg 2.7 rh_illum += parr[3*i+1] * rh_cos(i) * rh_dom[i];
1066 greg 2.1
1067 greg 2.7 return rh_illum;
1068 greg 2.1 }
1069    
1070     /* Calculate earth orbit eccentricity correction factor */
1071    
1072     /* Reference: Sen, Z. 2008. Solar Energy Fundamental and Modeling */
1073     /* Techniques. Springer, p. 72. */
1074    
1075     double CalcEccentricity()
1076     {
1077     double day_angle; /* Day angle (radians) */
1078     double E0; /* Eccentricity */
1079    
1080     /* Calculate day angle */
1081     day_angle = (julian_date - 1.0) * (2.0 * PI / 365.0);
1082    
1083     /* Calculate eccentricity */
1084     E0 = 1.00011 + 0.034221 * cos(day_angle) + 0.00128 * sin(day_angle)
1085     + 0.000719 * cos(2.0 * day_angle) + 0.000077 * sin(2.0 *
1086     day_angle);
1087    
1088     return E0;
1089     }
1090    
1091     /* Calculate atmospheric precipitable water content */
1092    
1093     /* Reference: Perez, R., P. Ineichen, R. Seals, J. Michalsky, and R. */
1094     /* Stewart. 1990. ìModeling Daylight Availability and */
1095     /* Irradiance Components from Direct and Global */
1096     /* Irradiance,î Solar Energy 44(5):271-289, Eqn. 3. */
1097    
1098     /* Note: The default surface dew point temperature is 11 deg. C */
1099     /* (52 deg. F). Typical values are: */
1100    
1101     /* Celsius Fahrenheit Human Perception */
1102     /* > 24 > 75 Extremely uncomfortable */
1103     /* 21 - 24 70 - 74 Very humid */
1104     /* 18 - 21 65 - 69 Somewhat uncomfortable */
1105     /* 16 - 18 60 - 64 OK for most people */
1106     /* 13 - 16 55 - 59 Comfortable */
1107     /* 10 - 12 50 - 54 Very comfortable */
1108     /* < 10 < 49 A bit dry for some */
1109    
1110     double CalcPrecipWater( double dpt )
1111     { return exp(0.07 * dpt - 0.075); }
1112    
1113     /* Calculate relative air mass */
1114    
1115     /* Reference: Kasten, F. 1966. "A New Table and Approximation Formula */
1116     /* for the Relative Optical Air Mass," Arch. Meteorol. */
1117     /* Geophys. Bioklimataol. Ser. B14, pp. 206-233. */
1118    
1119     /* Note: More sophisticated relative air mass models are */
1120     /* available, but they differ significantly only for */
1121     /* sun zenith angles greater than 80 degrees. */
1122    
1123     double CalcAirMass()
1124     {
1125     return (1.0 / (cos(sun_zenith) + 0.15 * pow(93.885 -
1126     RadToDeg(sun_zenith), -1.253)));
1127     }
1128    
1129     /* Calculate Perez All-Weather sky patch luminances (modified by GW) */
1130    
1131     /* NOTE: The sky patches centers are determined in accordance with the */
1132     /* BRE-IDMP sky luminance measurement procedures. (See for example */
1133     /* Mardaljevic, J. 2001. "The BRE-IDMP Dataset: A New Benchmark */
1134     /* for the Validation of Illuminance Prediction Techniques," */
1135     /* Lighting Research & Technology 33(2):117-136.) */
1136    
1137     void CalcSkyPatchLumin( float *parr )
1138     {
1139     int i;
1140     double aas; /* Sun-sky point azimuthal angle */
1141     double sspa; /* Sun-sky point angle */
1142     double zsa; /* Zenithal sun angle */
1143    
1144     for (i = 1; i < nskypatch; i++)
1145     {
1146     /* Calculate sun-sky point azimuthal angle */
1147     aas = fabs(rh_pazi[i] - azimuth);
1148    
1149     /* Calculate zenithal sun angle */
1150     zsa = PI * 0.5 - rh_palt[i];
1151    
1152     /* Calculate sun-sky point angle (Equation 8-20) */
1153     sspa = acos(cos(sun_zenith) * cos(zsa) + sin(sun_zenith) *
1154     sin(zsa) * cos(aas));
1155    
1156     /* Calculate patch luminance */
1157     parr[3*i] = CalcRelLuminance(sspa, zsa);
1158     if (parr[3*i] < 0) parr[3*i] = 0;
1159     parr[3*i+2] = parr[3*i+1] = parr[3*i];
1160     }
1161     }