ViewVC Help
View File | Revision Log | Show Annotations | Download File | Root Listing
root/radiance/ray/src/gen/gendaymtx.c
Revision: 2.26
Committed: Fri Apr 28 16:07:34 2017 UTC (6 years, 11 months ago) by greg
Content type: text/plain
Branch: MAIN
CVS Tags: rad5R2, rad5R1
Changes since 2.25: +4 -2 lines
Log Message:
Eliminated NaN error pointed out by David Geisler-Moroder

File Contents

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