ViewVC Help
View File | Revision Log | Show Annotations | Download File | Root Listing
root/radiance/ray/src/gen/gendaymtx.c
Revision: 2.12
Committed: Thu Jun 6 20:21:10 2013 UTC (10 years, 9 months ago) by greg
Content type: text/plain
Branch: MAIN
Changes since 2.11: +2 -2 lines
Log Message:
Fixed potential bug pointed out by Ian Ashdown

File Contents

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