ViewVC Help
View File | Revision Log | Show Annotations | Download File | Root Listing
root/radiance/ray/src/gen/gendaymtx.c
Revision: 2.7
Committed: Sat Jan 26 00:59:08 2013 UTC (11 years, 2 months ago) by greg
Content type: text/plain
Branch: MAIN
Changes since 2.6: +4 -11 lines
Log Message:
Fixed silly coding error that caused everything to be wrong

File Contents

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