ViewVC Help
View File | Revision Log | Show Annotations | Download File | Root Listing
root/radiance/ray/src/gen/gendaymtx.c
Revision: 2.3
Committed: Sat Jan 19 20:38:36 2013 UTC (11 years, 3 months ago) by greg
Content type: text/plain
Branch: MAIN
Changes since 2.2: +20 -15 lines
Log Message:
More bug fixes -- seems to be working OK for the most part

File Contents

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