ViewVC Help
View File | Revision Log | Show Annotations | Download File | Root Listing
root/radiance/ray/src/gen/gendaymtx.c
Revision: 2.9
Committed: Sat Feb 9 00:20:24 2013 UTC (11 years, 2 months ago) by greg
Content type: text/plain
Branch: MAIN
Changes since 2.8: +11 -3 lines
Log Message:
Fixed problem with sky clearness >= 12

File Contents

# User Rev Content
1 greg 2.1 #ifndef lint
2 greg 2.9 static const char RCSid[] = "$Id: gendaymtx.c,v 2.8 2013/02/05 06:00:19 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     int verbose = 0; /* progress reports to stderr? */
264    
265     int outfmt = 'a'; /* output format */
266    
267     int rhsubdiv = 1; /* Reinhart sky subdivisions */
268    
269 greg 2.4 COLOR skycolor = {.96, 1.004, 1.118}; /* sky coloration */
270     COLOR suncolor = {1., 1., 1.}; /* sun color */
271     COLOR grefl = {.2, .2, .2}; /* ground reflectance */
272 greg 2.1
273     int nskypatch; /* number of Reinhart patches */
274     float *rh_palt; /* sky patch altitudes (radians) */
275     float *rh_pazi; /* sky patch azimuths (radians) */
276     float *rh_dom; /* sky patch solid angle (sr) */
277    
278     #define vector(v,alt,azi) ( (v)[1] = tcos(alt), \
279     (v)[0] = (v)[1]*tsin(azi), \
280     (v)[1] *= tcos(azi), \
281     (v)[2] = tsin(alt) )
282    
283     #define rh_vector(v,i) vector(v,rh_palt[i],rh_pazi[i])
284    
285     #define rh_cos(i) tsin(rh_palt[i])
286    
287     extern int rh_init(void);
288     extern float * resize_dmatrix(float *mtx_data, int nsteps, int npatch);
289     extern void AddDirect(float *parr);
290    
291     int
292     main(int argc, char *argv[])
293     {
294     char buf[256];
295 greg 2.8 double rotation = 0; /* site rotation (degrees) */
296 greg 2.1 double elevation; /* site elevation (meters) */
297     int dir_is_horiz; /* direct is meas. on horizontal? */
298     float *mtx_data = NULL; /* our matrix data */
299     int ntsteps = 0; /* number of rows in matrix */
300     int last_monthly = 0; /* month of last report */
301     int mo, da; /* month (1-12) and day (1-31) */
302     double hr; /* hour (local standard time) */
303     double dir, dif; /* direct and diffuse values */
304     int mtx_offset;
305     int i, j;
306    
307     progname = argv[0];
308     /* get options */
309     for (i = 1; i < argc && argv[i][0] == '-'; i++)
310     switch (argv[i][1]) {
311 greg 2.4 case 'g': /* ground reflectance */
312     grefl[0] = atof(argv[++i]);
313     grefl[1] = atof(argv[++i]);
314     grefl[2] = atof(argv[++i]);
315 greg 2.1 break;
316 greg 2.4 case 'v': /* verbose progress reports */
317 greg 2.1 verbose++;
318     break;
319 greg 2.4 case 'o': /* output format */
320 greg 2.1 switch (argv[i][2]) {
321     case 'f':
322     case 'd':
323     case 'a':
324     outfmt = argv[i][2];
325     break;
326     default:
327     goto userr;
328     }
329     break;
330 greg 2.4 case 'm': /* Reinhart subdivisions */
331 greg 2.1 rhsubdiv = atoi(argv[++i]);
332     break;
333 greg 2.4 case 'c': /* sky color */
334 greg 2.1 skycolor[0] = atof(argv[++i]);
335     skycolor[1] = atof(argv[++i]);
336     skycolor[2] = atof(argv[++i]);
337     break;
338 greg 2.4 case 'd': /* solar (direct) only */
339 greg 2.1 skycolor[0] = skycolor[1] = skycolor[2] = 0;
340 greg 2.4 if (suncolor[1] <= 1e-4)
341     suncolor[0] = suncolor[1] = suncolor[2] = 1;
342 greg 2.1 break;
343 greg 2.4 case 's': /* sky only (no direct) */
344     suncolor[0] = suncolor[1] = suncolor[2] = 0;
345 greg 2.1 if (skycolor[1] <= 1e-4)
346     skycolor[0] = skycolor[1] = skycolor[2] = 1;
347     break;
348 greg 2.8 case 'r': /* rotate distribution */
349     if (argv[i][2] && argv[i][2] != 'z')
350     goto userr;
351     rotation = atof(argv[++i]);
352     break;
353 greg 2.1 default:
354     goto userr;
355     }
356     if (i < argc-1)
357     goto userr;
358     if (i == argc-1 && freopen(argv[i], "r", stdin) == NULL) {
359     fprintf(stderr, "%s: cannot open '%s' for input\n",
360     progname, argv[i]);
361     exit(1);
362     }
363     if (verbose) {
364     if (i == argc-1)
365     fprintf(stderr, "%s: reading weather tape '%s'\n",
366     progname, argv[i]);
367     else
368     fprintf(stderr, "%s: reading weather tape from <stdin>\n",
369     progname);
370     }
371     /* read weather tape header */
372 greg 2.2 if (scanf("place %[^\r\n] ", buf) != 1)
373 greg 2.1 goto fmterr;
374     if (scanf("latitude %lf\n", &s_latitude) != 1)
375     goto fmterr;
376     if (scanf("longitude %lf\n", &s_longitude) != 1)
377     goto fmterr;
378     if (scanf("time_zone %lf\n", &s_meridian) != 1)
379     goto fmterr;
380     if (scanf("site_elevation %lf\n", &elevation) != 1)
381     goto fmterr;
382     if (scanf("weather_data_file_units %d\n", &input) != 1)
383     goto fmterr;
384     switch (input) { /* translate units */
385     case 1:
386     input = 1; /* radiometric quantities */
387     dir_is_horiz = 0; /* direct is perpendicular meas. */
388     break;
389     case 2:
390     input = 1; /* radiometric quantities */
391     dir_is_horiz = 1; /* solar measured horizontally */
392     break;
393     case 3:
394     input = 2; /* photometric quantities */
395     dir_is_horiz = 0; /* direct is perpendicular meas. */
396     break;
397     default:
398     goto fmterr;
399     }
400     rh_init(); /* initialize sky patches */
401     if (verbose) {
402     fprintf(stderr, "%s: location '%s'\n", progname, buf);
403     fprintf(stderr, "%s: (lat,long)=(%.1f,%.1f) degrees north, west\n",
404     progname, s_latitude, s_longitude);
405     fprintf(stderr, "%s: %d sky patches per time step\n",
406     progname, nskypatch);
407 greg 2.8 if (rotation != 0)
408     fprintf(stderr, "%s: rotating output %.0f degrees\n",
409     progname, rotation);
410 greg 2.1 }
411 greg 2.2 /* convert quantities to radians */
412     s_latitude = DegToRad(s_latitude);
413     s_longitude = DegToRad(s_longitude);
414     s_meridian = DegToRad(s_meridian);
415 greg 2.1 /* process each time step in tape */
416     while (scanf("%d %d %lf %lf %lf\n", &mo, &da, &hr, &dir, &dif) == 5) {
417     double sda, sta;
418     /* make space for next time step */
419     mtx_offset = 3*nskypatch*ntsteps++;
420     mtx_data = resize_dmatrix(mtx_data, ntsteps, nskypatch);
421     if (dif <= 1e-4) {
422     memset(mtx_data+mtx_offset, 0, sizeof(float)*3*nskypatch);
423     continue;
424     }
425     if (verbose && mo != last_monthly)
426     fprintf(stderr, "%s: stepping through month %d...\n",
427     progname, last_monthly=mo);
428     /* compute solar position */
429     julian_date = jdate(mo, da);
430     sda = sdec(julian_date);
431     sta = stadj(julian_date);
432     altitude = salt(sda, hr+sta);
433 greg 2.8 azimuth = sazi(sda, hr+sta) + PI - DegToRad(rotation);
434 greg 2.1 /* convert measured values */
435     if (dir_is_horiz && altitude > 0.)
436     dir /= sin(altitude);
437     if (input == 1) {
438     dir_irrad = dir;
439     diff_irrad = dif;
440     } else /* input == 2 */ {
441     dir_illum = dir;
442     diff_illum = dif;
443     }
444     /* compute sky patch values */
445     ComputeSky(mtx_data+mtx_offset);
446 greg 2.4 AddDirect(mtx_data+mtx_offset);
447 greg 2.1 }
448     /* check for junk at end */
449     while ((i = fgetc(stdin)) != EOF)
450     if (!isspace(i)) {
451     fprintf(stderr, "%s: warning - unexpected data past EOT: ",
452     progname);
453     buf[0] = i; buf[1] = '\0';
454     fgets(buf+1, sizeof(buf)-1, stdin);
455     fputs(buf, stderr); fputc('\n', stderr);
456     break;
457     }
458     /* write out matrix */
459     #ifdef getc_unlocked
460     flockfile(stdout);
461     #endif
462     if (verbose)
463     fprintf(stderr, "%s: writing %smatrix with %d time steps...\n",
464     progname, outfmt=='a' ? "" : "binary ", ntsteps);
465     /* patches are rows (outer sort) */
466     for (i = 0; i < nskypatch; i++) {
467     mtx_offset = 3*i;
468     switch (outfmt) {
469     case 'a':
470     for (j = 0; j < ntsteps; j++) {
471 greg 2.3 printf("%.3g %.3g %.3g\n", mtx_data[mtx_offset],
472 greg 2.1 mtx_data[mtx_offset+1],
473     mtx_data[mtx_offset+2]);
474     mtx_offset += 3*nskypatch;
475     }
476 greg 2.2 if (ntsteps > 1)
477     fputc('\n', stdout);
478 greg 2.1 break;
479     case 'f':
480     for (j = 0; j < ntsteps; j++) {
481     fwrite(mtx_data+mtx_offset, sizeof(float), 3,
482     stdout);
483     mtx_offset += 3*nskypatch;
484     }
485     break;
486     case 'd':
487     for (j = 0; j < ntsteps; j++) {
488     double ment[3];
489     ment[0] = mtx_data[mtx_offset];
490     ment[1] = mtx_data[mtx_offset+1];
491     ment[2] = mtx_data[mtx_offset+2];
492     fwrite(ment, sizeof(double), 3, stdout);
493     mtx_offset += 3*nskypatch;
494     }
495     break;
496     }
497     if (ferror(stdout))
498     goto writerr;
499     }
500     if (fflush(stdout) == EOF)
501     goto writerr;
502     if (verbose)
503     fprintf(stderr, "%s: done.\n", progname);
504     exit(0);
505     userr:
506 greg 2.8 fprintf(stderr, "Usage: %s [-v][-d|-s][-r deg][-m N][-g r g b][-c r g b][-o{f|d}] [tape.wea]\n",
507 greg 2.1 progname);
508     exit(1);
509     fmterr:
510     fprintf(stderr, "%s: input weather tape format error\n", progname);
511     exit(1);
512     writerr:
513     fprintf(stderr, "%s: write error on output\n", progname);
514     exit(1);
515     }
516    
517     /* Return maximum of two doubles */
518     double dmax( double a, double b )
519     { return (a > b) ? a : b; }
520    
521     /* Compute sky patch radiance values (modified by GW) */
522     void
523     ComputeSky(float *parr)
524     {
525     int index; /* Category index */
526     double norm_diff_illum; /* Normalized diffuse illuimnance */
527     int i;
528    
529     /* Calculate atmospheric precipitable water content */
530     apwc = CalcPrecipWater(dew_point);
531    
532 greg 2.6 /* Calculate sun zenith angle (don't let it dip below horizon) */
533     /* Also limit minimum angle to keep circumsolar off zenith */
534     if (altitude <= 0.0)
535     sun_zenith = DegToRad(90.0);
536     else if (altitude >= DegToRad(87.0))
537     sun_zenith = DegToRad(3.0);
538     else
539     sun_zenith = DegToRad(90.0) - altitude;
540 greg 2.1
541     /* Compute the inputs for the calculation of the sky distribution */
542    
543     if (input == 0) /* XXX never used */
544     {
545     /* Calculate irradiance */
546     diff_irrad = CalcDiffuseIrradiance();
547     dir_irrad = CalcDirectIrradiance();
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 == 1)
555     {
556     sky_brightness = CalcSkyBrightness();
557     sky_clearness = CalcSkyClearness();
558    
559 greg 2.9 /* Limit sky clearness */
560     if (sky_clearness > 11.9)
561     sky_clearness = 11.9;
562    
563     /* Limit sky brightness */
564     if (sky_brightness < 0.01)
565     sky_brightness = 0.01;
566    
567 greg 2.1 /* Calculate illuminance */
568     index = GetCategoryIndex();
569     diff_illum = diff_irrad * CalcDiffuseIllumRatio(index);
570     dir_illum = dir_irrad * CalcDirectIllumRatio(index);
571     }
572     else if (input == 2)
573     {
574     /* Calculate sky brightness and clearness from illuminance values */
575     index = CalcSkyParamFromIllum();
576     }
577    
578 greg 2.2 if (bright(skycolor) <= 1e-4) { /* 0 sky component? */
579     memset(parr, 0, sizeof(float)*3*nskypatch);
580     return;
581     }
582 greg 2.1 /* Compute ground radiance (include solar contribution if any) */
583 greg 2.3 parr[0] = diff_illum;
584 greg 2.1 if (altitude > 0)
585 greg 2.3 parr[0] += dir_illum * sin(altitude);
586 greg 2.4 parr[2] = parr[1] = parr[0] *= (1./PI/WHTEFFICACY);
587     multcolor(parr, grefl);
588 greg 2.1
589     /* Calculate Perez sky model parameters */
590     CalcPerezParam(sun_zenith, sky_clearness, sky_brightness, index);
591    
592     /* Calculate sky patch luminance values */
593     CalcSkyPatchLumin(parr);
594    
595     /* Calculate relative horizontal illuminance */
596     norm_diff_illum = CalcRelHorzIllum(parr);
597    
598     /* Normalization coefficient */
599     norm_diff_illum = diff_illum / norm_diff_illum;
600    
601     /* Apply to sky patches to get absolute radiance values */
602     for (i = 1; i < nskypatch; i++) {
603 greg 2.7 scalecolor(parr+3*i, norm_diff_illum*(1./WHTEFFICACY));
604 greg 2.1 multcolor(parr+3*i, skycolor);
605     }
606     }
607    
608     /* Add in solar direct to nearest sky patches (GW) */
609     void
610     AddDirect(float *parr)
611     {
612     FVECT svec;
613 greg 2.3 double near_dprod[NSUNPATCH];
614     int near_patch[NSUNPATCH];
615     double wta[NSUNPATCH], wtot;
616 greg 2.1 int i, j, p;
617    
618 greg 2.4 if (dir_illum <= 1e-4 || bright(suncolor) <= 1e-4)
619 greg 2.1 return;
620 greg 2.3 /* identify NSUNPATCH closest patches */
621     for (i = NSUNPATCH; i--; )
622 greg 2.1 near_dprod[i] = -1.;
623     vector(svec, altitude, azimuth);
624     for (p = 1; p < nskypatch; p++) {
625     FVECT pvec;
626     double dprod;
627     rh_vector(pvec, p);
628     dprod = DOT(pvec, svec);
629 greg 2.3 for (i = 0; i < NSUNPATCH; i++)
630 greg 2.1 if (dprod > near_dprod[i]) {
631 greg 2.3 for (j = NSUNPATCH; --j > i; ) {
632 greg 2.1 near_dprod[j] = near_dprod[j-1];
633     near_patch[j] = near_patch[j-1];
634     }
635     near_dprod[i] = dprod;
636     near_patch[i] = p;
637     break;
638     }
639     }
640     wtot = 0; /* weight by proximity */
641 greg 2.3 for (i = NSUNPATCH; i--; )
642 greg 2.1 wtot += wta[i] = 1./(1.002 - near_dprod[i]);
643     /* add to nearest patch radiances */
644 greg 2.3 for (i = NSUNPATCH; i--; ) {
645 greg 2.2 float *pdest = parr + 3*near_patch[i];
646     float val_add = wta[i] * dir_illum /
647     (WHTEFFICACY * wtot * rh_dom[near_patch[i]]);
648 greg 2.4 *pdest++ += val_add*suncolor[0];
649     *pdest++ += val_add*suncolor[1];
650     *pdest++ += val_add*suncolor[2];
651 greg 2.2 }
652 greg 2.1 }
653    
654     /* Initialize Reinhart sky patch positions (GW) */
655     int
656     rh_init(void)
657     {
658     #define NROW 7
659     static const int tnaz[NROW] = {30, 30, 24, 24, 18, 12, 6};
660     const double alpha = (PI/2.)/(NROW*rhsubdiv + .5);
661     int p, i, j;
662     /* allocate patch angle arrays */
663     nskypatch = 0;
664     for (p = 0; p < NROW; p++)
665     nskypatch += tnaz[p];
666     nskypatch *= rhsubdiv*rhsubdiv;
667     nskypatch += 2;
668     rh_palt = (float *)malloc(sizeof(float)*nskypatch);
669     rh_pazi = (float *)malloc(sizeof(float)*nskypatch);
670     rh_dom = (float *)malloc(sizeof(float)*nskypatch);
671     if ((rh_palt == NULL) | (rh_pazi == NULL) | (rh_dom == NULL)) {
672     fprintf(stderr, "%s: out of memory in rh_init()\n", progname);
673     exit(1);
674     }
675     rh_palt[0] = -PI/2.; /* ground & zenith patches */
676     rh_pazi[0] = 0.;
677     rh_dom[0] = 2.*PI;
678     rh_palt[nskypatch-1] = PI/2.;
679     rh_pazi[nskypatch-1] = 0.;
680     rh_dom[nskypatch-1] = 2.*PI*(1. - cos(alpha*.5));
681     p = 1; /* "normal" patches */
682     for (i = 0; i < NROW*rhsubdiv; i++) {
683     const float ralt = alpha*(i + .5);
684     const int ninrow = tnaz[i/rhsubdiv]*rhsubdiv;
685 greg 2.3 const float dom = 2.*PI*(sin(alpha*(i+1)) - sin(alpha*i)) /
686     (double)ninrow;
687 greg 2.1 for (j = 0; j < ninrow; j++) {
688     rh_palt[p] = ralt;
689     rh_pazi[p] = 2.*PI * j / (double)ninrow;
690     rh_dom[p++] = dom;
691     }
692     }
693     return nskypatch;
694     #undef NROW
695     }
696    
697     /* Resize daylight matrix (GW) */
698     float *
699     resize_dmatrix(float *mtx_data, int nsteps, int npatch)
700     {
701     if (mtx_data == NULL)
702     mtx_data = (float *)malloc(sizeof(float)*3*nsteps*npatch);
703     else
704     mtx_data = (float *)realloc(mtx_data,
705     sizeof(float)*3*nsteps*npatch);
706     if (mtx_data == NULL) {
707     fprintf(stderr, "%s: out of memory in resize_dmatrix(%d,%d)\n",
708     progname, nsteps, npatch);
709     exit(1);
710     }
711     return(mtx_data);
712     }
713    
714     /* Determine category index */
715     int GetCategoryIndex()
716     {
717     int index; /* Loop index */
718    
719     for (index = 0; index < 8; index++)
720     if ((sky_clearness >= SkyClearCat[index].lower) &&
721     (sky_clearness < SkyClearCat[index].upper))
722     break;
723    
724     return index;
725     }
726    
727     /* Calculate diffuse illuminance to diffuse irradiance ratio */
728    
729     /* Reference: Perez, R., P. Ineichen, R. Seals, J. Michalsky, and R. */
730     /* Stewart. 1990. ìModeling Daylight Availability and */
731     /* Irradiance Components from Direct and Global */
732     /* Irradiance,î Solar Energy 44(5):271-289, Eqn. 7. */
733    
734     double CalcDiffuseIllumRatio( int index )
735     {
736     ModelCoeff const *pnle; /* Category coefficient pointer */
737    
738     /* Get category coefficient pointer */
739     pnle = &(DiffuseLumEff[index]);
740    
741     return pnle->a + pnle->b * apwc + pnle->c * cos(sun_zenith) +
742     pnle->d * log(sky_brightness);
743     }
744    
745     /* Calculate direct illuminance to direct irradiance ratio */
746    
747     /* Reference: Perez, R., P. Ineichen, R. Seals, J. Michalsky, and R. */
748     /* Stewart. 1990. ìModeling Daylight Availability and */
749     /* Irradiance Components from Direct and Global */
750     /* Irradiance,î Solar Energy 44(5):271-289, Eqn. 8. */
751    
752     double CalcDirectIllumRatio( int index )
753     {
754     ModelCoeff const *pnle; /* Category coefficient pointer */
755    
756     /* Get category coefficient pointer */
757     pnle = &(DirectLumEff[index]);
758    
759     /* Calculate direct illuminance from direct irradiance */
760    
761     return dmax((pnle->a + pnle->b * apwc + pnle->c * exp(5.73 *
762     sun_zenith - 5.0) + pnle->d * sky_brightness),
763     0.0);
764     }
765    
766     /* Calculate sky brightness */
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. 2. */
772    
773     double CalcSkyBrightness()
774     {
775     return diff_irrad * CalcAirMass() / (DC_SolarConstantE *
776     CalcEccentricity());
777     }
778    
779     /* Calculate sky clearness */
780    
781     /* Reference: Perez, R., P. Ineichen, R. Seals, J. Michalsky, and R. */
782     /* Stewart. 1990. ìModeling Daylight Availability and */
783     /* Irradiance Components from Direct and Global */
784     /* Irradiance,î Solar Energy 44(5):271-289, Eqn. 1. */
785    
786     double CalcSkyClearness()
787     {
788     double sz_cubed; /* Sun zenith angle cubed */
789    
790     /* Calculate sun zenith angle cubed */
791     sz_cubed = pow(sun_zenith, 3.0);
792    
793     return ((diff_irrad + dir_irrad) / diff_irrad + 1.041 *
794     sz_cubed) / (1.0 + 1.041 * sz_cubed);
795     }
796    
797     /* Calculate diffuse horizontal irradiance from Perez sky brightness */
798    
799     /* Reference: Perez, R., P. Ineichen, R. Seals, J. Michalsky, and R. */
800     /* Stewart. 1990. ìModeling Daylight Availability and */
801     /* Irradiance Components from Direct and Global */
802     /* Irradiance,î Solar Energy 44(5):271-289, Eqn. 2 */
803     /* (inverse). */
804    
805     double CalcDiffuseIrradiance()
806     {
807     return sky_brightness * DC_SolarConstantE * CalcEccentricity() /
808     CalcAirMass();
809     }
810    
811     /* Calculate direct normal irradiance from Perez sky clearness */
812    
813     /* Reference: Perez, R., P. Ineichen, R. Seals, J. Michalsky, and R. */
814     /* Stewart. 1990. ìModeling Daylight Availability and */
815     /* Irradiance Components from Direct and Global */
816     /* Irradiance,î Solar Energy 44(5):271-289, Eqn. 1 */
817     /* (inverse). */
818    
819     double CalcDirectIrradiance()
820     {
821     return CalcDiffuseIrradiance() * ((sky_clearness - 1.0) * (1 + 1.041
822     * pow(sun_zenith, 3.0)));
823     }
824    
825     /* Calculate sky brightness and clearness from illuminance values */
826     int CalcSkyParamFromIllum()
827     {
828     double test1 = 0.1;
829     double test2 = 0.1;
830     int counter = 0;
831     int index = 0; /* Category index */
832    
833     /* Convert illuminance to irradiance */
834     diff_irrad = diff_illum * DC_SolarConstantE /
835     (DC_SolarConstantL * 1000.0);
836     dir_irrad = dir_illum * DC_SolarConstantE /
837     (DC_SolarConstantL * 1000.0);
838    
839     /* Calculate sky brightness and clearness */
840     sky_brightness = CalcSkyBrightness();
841     sky_clearness = CalcSkyClearness();
842    
843     /* Limit sky clearness */
844     if (sky_clearness > 12.0)
845     sky_clearness = 12.0;
846    
847     /* Limit sky brightness */
848 greg 2.9 if (sky_brightness < 0.01)
849 greg 2.1 sky_brightness = 0.01;
850    
851     while (((fabs(diff_irrad - test1) > 10.0) ||
852     (fabs(dir_irrad - test2) > 10.0)) && !(counter == 5))
853     {
854     test1 = diff_irrad;
855     test2 = dir_irrad;
856     counter++;
857    
858     /* Convert illuminance to irradiance */
859     index = GetCategoryIndex();
860     diff_irrad = diff_illum / CalcDiffuseIllumRatio(index);
861     dir_irrad = dir_illum / CalcDirectIllumRatio(index);
862    
863     /* Calculate sky brightness and clearness */
864     sky_brightness = CalcSkyBrightness();
865     sky_clearness = CalcSkyClearness();
866    
867     /* Limit sky clearness */
868     if (sky_clearness > 12.0)
869     sky_clearness = 12.0;
870    
871     /* Limit sky brightness */
872 greg 2.9 if (sky_brightness < 0.01)
873 greg 2.1 sky_brightness = 0.01;
874     }
875    
876     return GetCategoryIndex();
877     }
878    
879     /* Calculate relative luminance */
880    
881     /* Reference: Perez, R., R. Seals, and J. Michalsky. 1993. */
882     /* ìAll-Weather Model for Sky Luminance Distribution - */
883     /* Preliminary Configuration and Validation,î Solar Energy */
884     /* 50(3):235-245, Eqn. 1. */
885    
886     double CalcRelLuminance( double gamma, double zeta )
887     {
888     return (1.0 + perez_param[0] * exp(perez_param[1] / cos(zeta))) *
889     (1.0 + perez_param[2] * exp(perez_param[3] * gamma) +
890     perez_param[4] * cos(gamma) * cos(gamma));
891     }
892    
893     /* Calculate Perez sky model parameters */
894    
895     /* Reference: Perez, R., R. Seals, and J. Michalsky. 1993. */
896     /* ìAll-Weather Model for Sky Luminance Distribution - */
897     /* Preliminary Configuration and Validation,î Solar Energy */
898     /* 50(3):235-245, Eqns. 6 - 8. */
899    
900     void CalcPerezParam( double sz, double epsilon, double delta,
901     int index )
902     {
903     double x[5][4]; /* Coefficents a, b, c, d, e */
904     int i, j; /* Loop indices */
905    
906     /* Limit sky brightness */
907     if (epsilon > 1.065 && epsilon < 2.8)
908     {
909     if (delta < 0.2)
910     delta = 0.2;
911     }
912    
913     /* Get Perez coefficients */
914     for (i = 0; i < 5; i++)
915     for (j = 0; j < 4; j++)
916     x[i][j] = PerezCoeff[index][4 * i + j];
917    
918     if (index != 0)
919     {
920     /* Calculate parameter a, b, c, d and e (Eqn. 6) */
921     for (i = 0; i < 5; i++)
922     perez_param[i] = x[i][0] + x[i][1] * sz + delta * (x[i][2] +
923     x[i][3] * sz);
924     }
925     else
926     {
927     /* Parameters a, b and e (Eqn. 6) */
928     perez_param[0] = x[0][0] + x[0][1] * sz + delta * (x[0][2] +
929     x[0][3] * sz);
930     perez_param[1] = x[1][0] + x[1][1] * sz + delta * (x[1][2] +
931     x[1][3] * sz);
932     perez_param[4] = x[4][0] + x[4][1] * sz + delta * (x[4][2] +
933     x[4][3] * sz);
934    
935     /* Parameter c (Eqn. 7) */
936     perez_param[2] = exp(pow(delta * (x[2][0] + x[2][1] * sz),
937     x[2][2])) - x[2][3];
938    
939     /* Parameter d (Eqn. 8) */
940     perez_param[3] = -exp(delta * (x[3][0] + x[3][1] * sz)) +
941     x[3][2] + delta * x[3][3];
942     }
943     }
944    
945     /* Calculate relative horizontal illuminance (modified by GW) */
946    
947     /* Reference: Perez, R., R. Seals, and J. Michalsky. 1993. */
948     /* ìAll-Weather Model for Sky Luminance Distribution - */
949     /* Preliminary Configuration and Validation,î Solar Energy */
950     /* 50(3):235-245, Eqn. 3. */
951    
952     double CalcRelHorzIllum( float *parr )
953     {
954     int i;
955     double rh_illum = 0.0; /* Relative horizontal illuminance */
956    
957     for (i = 1; i < nskypatch; i++)
958 greg 2.7 rh_illum += parr[3*i+1] * rh_cos(i) * rh_dom[i];
959 greg 2.1
960 greg 2.7 return rh_illum;
961 greg 2.1 }
962    
963     /* Calculate earth orbit eccentricity correction factor */
964    
965     /* Reference: Sen, Z. 2008. Solar Energy Fundamental and Modeling */
966     /* Techniques. Springer, p. 72. */
967    
968     double CalcEccentricity()
969     {
970     double day_angle; /* Day angle (radians) */
971     double E0; /* Eccentricity */
972    
973     /* Calculate day angle */
974     day_angle = (julian_date - 1.0) * (2.0 * PI / 365.0);
975    
976     /* Calculate eccentricity */
977     E0 = 1.00011 + 0.034221 * cos(day_angle) + 0.00128 * sin(day_angle)
978     + 0.000719 * cos(2.0 * day_angle) + 0.000077 * sin(2.0 *
979     day_angle);
980    
981     return E0;
982     }
983    
984     /* Calculate atmospheric precipitable water content */
985    
986     /* Reference: Perez, R., P. Ineichen, R. Seals, J. Michalsky, and R. */
987     /* Stewart. 1990. ìModeling Daylight Availability and */
988     /* Irradiance Components from Direct and Global */
989     /* Irradiance,î Solar Energy 44(5):271-289, Eqn. 3. */
990    
991     /* Note: The default surface dew point temperature is 11 deg. C */
992     /* (52 deg. F). Typical values are: */
993    
994     /* Celsius Fahrenheit Human Perception */
995     /* > 24 > 75 Extremely uncomfortable */
996     /* 21 - 24 70 - 74 Very humid */
997     /* 18 - 21 65 - 69 Somewhat uncomfortable */
998     /* 16 - 18 60 - 64 OK for most people */
999     /* 13 - 16 55 - 59 Comfortable */
1000     /* 10 - 12 50 - 54 Very comfortable */
1001     /* < 10 < 49 A bit dry for some */
1002    
1003     double CalcPrecipWater( double dpt )
1004     { return exp(0.07 * dpt - 0.075); }
1005    
1006     /* Calculate relative air mass */
1007    
1008     /* Reference: Kasten, F. 1966. "A New Table and Approximation Formula */
1009     /* for the Relative Optical Air Mass," Arch. Meteorol. */
1010     /* Geophys. Bioklimataol. Ser. B14, pp. 206-233. */
1011    
1012     /* Note: More sophisticated relative air mass models are */
1013     /* available, but they differ significantly only for */
1014     /* sun zenith angles greater than 80 degrees. */
1015    
1016     double CalcAirMass()
1017     {
1018     return (1.0 / (cos(sun_zenith) + 0.15 * pow(93.885 -
1019     RadToDeg(sun_zenith), -1.253)));
1020     }
1021    
1022     /* Calculate Perez All-Weather sky patch luminances (modified by GW) */
1023    
1024     /* NOTE: The sky patches centers are determined in accordance with the */
1025     /* BRE-IDMP sky luminance measurement procedures. (See for example */
1026     /* Mardaljevic, J. 2001. "The BRE-IDMP Dataset: A New Benchmark */
1027     /* for the Validation of Illuminance Prediction Techniques," */
1028     /* Lighting Research & Technology 33(2):117-136.) */
1029    
1030     void CalcSkyPatchLumin( float *parr )
1031     {
1032     int i;
1033     double aas; /* Sun-sky point azimuthal angle */
1034     double sspa; /* Sun-sky point angle */
1035     double zsa; /* Zenithal sun angle */
1036    
1037     for (i = 1; i < nskypatch; i++)
1038     {
1039     /* Calculate sun-sky point azimuthal angle */
1040     aas = fabs(rh_pazi[i] - azimuth);
1041    
1042     /* Calculate zenithal sun angle */
1043     zsa = PI * 0.5 - rh_palt[i];
1044    
1045     /* Calculate sun-sky point angle (Equation 8-20) */
1046     sspa = acos(cos(sun_zenith) * cos(zsa) + sin(sun_zenith) *
1047     sin(zsa) * cos(aas));
1048    
1049     /* Calculate patch luminance */
1050     parr[3*i] = CalcRelLuminance(sspa, zsa);
1051     if (parr[3*i] < 0) parr[3*i] = 0;
1052     parr[3*i+2] = parr[3*i+1] = parr[3*i];
1053     }
1054     }