ViewVC Help
View File | Revision Log | Show Annotations | Download File | Root Listing
root/radiance/ray/src/gen/gendaymtx.c
Revision: 2.39
Committed: Fri Mar 11 18:04:39 2022 UTC (2 years, 1 month ago) by greg
Content type: text/plain
Branch: MAIN
CVS Tags: rad5R4
Changes since 2.38: +3 -3 lines
Log Message:
feat(gendaymtx): Minor improvement in verbose output

File Contents

# User Rev Content
1 greg 2.1 #ifndef lint
2 greg 2.39 static const char RCSid[] = "$Id: gendaymtx.c,v 2.38 2020/09/11 16:50:50 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 <stdlib.h>
85     #include <ctype.h>
86 greg 2.30 #include "platform.h"
87 greg 2.1 #include "rtmath.h"
88 greg 2.23 #include "rtio.h"
89 greg 2.1 #include "color.h"
90 greg 2.30 #include "sun.h"
91 greg 2.1
92 greg 2.38 char *progname; /* Program name */
93 greg 2.1 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     /* 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 greg 2.12 { 6.200, 12.01 } /* Clear */
212 greg 2.1 };
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 greg 2.10 #define NSUNPATCH 4 /* max. # patches to spread sun into */
251 greg 2.3 #endif
252    
253 greg 2.35 #define SUN_ANG_DEG 0.533 /* sun full-angle in degrees */
254    
255 greg 2.10 int nsuns = NSUNPATCH; /* number of sun patches to use */
256     double fixed_sun_sa = -1; /* fixed solid angle per sun? */
257    
258 greg 2.1 int verbose = 0; /* progress reports to stderr? */
259    
260     int outfmt = 'a'; /* output format */
261    
262     int rhsubdiv = 1; /* Reinhart sky subdivisions */
263    
264 greg 2.4 COLOR skycolor = {.96, 1.004, 1.118}; /* sky coloration */
265     COLOR suncolor = {1., 1., 1.}; /* sun color */
266     COLOR grefl = {.2, .2, .2}; /* ground reflectance */
267 greg 2.1
268     int nskypatch; /* number of Reinhart patches */
269     float *rh_palt; /* sky patch altitudes (radians) */
270     float *rh_pazi; /* sky patch azimuths (radians) */
271     float *rh_dom; /* sky patch solid angle (sr) */
272    
273 greg 2.36 #define vector(v,alt,azi) ( (v)[1] = cos(alt), \
274     (v)[0] = (v)[1]*sin(azi), \
275     (v)[1] *= cos(azi), \
276     (v)[2] = sin(alt) )
277 greg 2.1
278     #define rh_vector(v,i) vector(v,rh_palt[i],rh_pazi[i])
279    
280     #define rh_cos(i) tsin(rh_palt[i])
281    
282 greg 2.36 #define solar_minute(jd,hr) ((24*60)*((jd)-1)+(int)((hr)*60.+.5))
283    
284 greg 2.1 extern int rh_init(void);
285     extern float * resize_dmatrix(float *mtx_data, int nsteps, int npatch);
286 greg 2.36 extern void OutputSun(int id, int goodsun, FILE *fp, FILE *mfp);
287 greg 2.1 extern void AddDirect(float *parr);
288    
289 greg 2.14
290     static const char *
291     getfmtname(int fmt)
292     {
293     switch (fmt) {
294     case 'a':
295     return("ascii");
296     case 'f':
297     return("float");
298     case 'd':
299     return("double");
300     }
301     return("unknown");
302     }
303    
304    
305 greg 2.1 int
306     main(int argc, char *argv[])
307     {
308     char buf[256];
309 greg 2.14 int doheader = 1; /* output header? */
310 greg 2.8 double rotation = 0; /* site rotation (degrees) */
311 greg 2.1 double elevation; /* site elevation (meters) */
312 greg 2.36 int leap_day = 0; /* add leap day? */
313 greg 2.37 int sun_hours_only = 0; /* only output sun hours? */
314 greg 2.1 int dir_is_horiz; /* direct is meas. on horizontal? */
315 greg 2.34 FILE *sunsfp = NULL; /* output file for individual suns */
316 greg 2.36 FILE *modsfp = NULL; /* modifier output file */
317 greg 2.1 float *mtx_data = NULL; /* our matrix data */
318 greg 2.28 int avgSky = 0; /* compute average sky r.t. matrix? */
319 greg 2.27 int ntsteps = 0; /* number of time steps */
320 greg 2.28 int tstorage = 0; /* number of allocated time steps */
321     int nstored = 0; /* number of time steps in matrix */
322 greg 2.1 int last_monthly = 0; /* month of last report */
323     int mo, da; /* month (1-12) and day (1-31) */
324     double hr; /* hour (local standard time) */
325     double dir, dif; /* direct and diffuse values */
326     int mtx_offset;
327     int i, j;
328    
329     progname = argv[0];
330     /* get options */
331     for (i = 1; i < argc && argv[i][0] == '-'; i++)
332     switch (argv[i][1]) {
333 greg 2.4 case 'g': /* ground reflectance */
334     grefl[0] = atof(argv[++i]);
335     grefl[1] = atof(argv[++i]);
336     grefl[2] = atof(argv[++i]);
337 greg 2.1 break;
338 greg 2.4 case 'v': /* verbose progress reports */
339 greg 2.1 verbose++;
340     break;
341 greg 2.14 case 'h': /* turn off header */
342     doheader = 0;
343     break;
344 greg 2.4 case 'o': /* output format */
345 greg 2.1 switch (argv[i][2]) {
346     case 'f':
347     case 'd':
348     case 'a':
349     outfmt = argv[i][2];
350     break;
351     default:
352     goto userr;
353     }
354     break;
355 greg 2.11 case 'O': /* output type */
356     switch (argv[i][2]) {
357     case '0':
358     output = 0;
359     break;
360     case '1':
361     output = 1;
362     break;
363     default:
364     goto userr;
365     }
366     if (argv[i][3])
367     goto userr;
368     break;
369 greg 2.4 case 'm': /* Reinhart subdivisions */
370 greg 2.1 rhsubdiv = atoi(argv[++i]);
371     break;
372 greg 2.4 case 'c': /* sky color */
373 greg 2.1 skycolor[0] = atof(argv[++i]);
374     skycolor[1] = atof(argv[++i]);
375     skycolor[2] = atof(argv[++i]);
376     break;
377 greg 2.36 case 'D': /* output suns to file */
378     if (strcmp(argv[++i], "-")) {
379     sunsfp = fopen(argv[i], "w");
380     if (sunsfp == NULL) {
381     fprintf(stderr,
382     "%s: cannot open '%s' for output\n",
383     progname, argv[i]);
384     exit(1);
385     }
386     break; /* still may output matrix */
387     }
388     sunsfp = stdout; /* sending to stdout, so... */
389     /* fall through */
390 greg 2.34 case 'n': /* no matrix output */
391     avgSky = -1;
392     rhsubdiv = 1;
393     /* fall through */
394 greg 2.4 case 'd': /* solar (direct) only */
395 greg 2.1 skycolor[0] = skycolor[1] = skycolor[2] = 0;
396 greg 2.33 grefl[0] = grefl[1] = grefl[2] = 0;
397 greg 2.1 break;
398 greg 2.36 case 'M': /* send sun modifiers to file */
399     if ((modsfp = fopen(argv[++i], "w")) == NULL) {
400 greg 2.34 fprintf(stderr, "%s: cannot open '%s' for output\n",
401 greg 2.35 progname, argv[i]);
402 greg 2.34 exit(1);
403     }
404     break;
405 greg 2.4 case 's': /* sky only (no direct) */
406     suncolor[0] = suncolor[1] = suncolor[2] = 0;
407 greg 2.1 break;
408 greg 2.37 case 'u': /* solar hours only */
409     sun_hours_only = 1;
410     break;
411 greg 2.8 case 'r': /* rotate distribution */
412     if (argv[i][2] && argv[i][2] != 'z')
413     goto userr;
414     rotation = atof(argv[++i]);
415     break;
416 greg 2.10 case '5': /* 5-phase calculation */
417     nsuns = 1;
418 greg 2.19 fixed_sun_sa = PI/360.*atof(argv[++i]);
419 greg 2.21 if (fixed_sun_sa <= 0) {
420     fprintf(stderr, "%s: missing solar disk size argument for '-5' option\n",
421 greg 2.35 progname);
422 greg 2.21 exit(1);
423     }
424 greg 2.19 fixed_sun_sa *= fixed_sun_sa*PI;
425 greg 2.10 break;
426 greg 2.27 case 'A': /* compute average sky */
427     avgSky = 1;
428     break;
429 greg 2.1 default:
430     goto userr;
431     }
432     if (i < argc-1)
433     goto userr;
434     if (i == argc-1 && freopen(argv[i], "r", stdin) == NULL) {
435     fprintf(stderr, "%s: cannot open '%s' for input\n",
436     progname, argv[i]);
437     exit(1);
438     }
439 greg 2.36 if ((modsfp != NULL) & (sunsfp == NULL))
440     fprintf(stderr, "%s: warning -M output will be empty without -D\n",
441     progname);
442 greg 2.1 if (verbose) {
443     if (i == argc-1)
444     fprintf(stderr, "%s: reading weather tape '%s'\n",
445     progname, argv[i]);
446     else
447     fprintf(stderr, "%s: reading weather tape from <stdin>\n",
448     progname);
449     }
450     /* read weather tape header */
451 greg 2.2 if (scanf("place %[^\r\n] ", buf) != 1)
452 greg 2.1 goto fmterr;
453     if (scanf("latitude %lf\n", &s_latitude) != 1)
454     goto fmterr;
455     if (scanf("longitude %lf\n", &s_longitude) != 1)
456     goto fmterr;
457     if (scanf("time_zone %lf\n", &s_meridian) != 1)
458     goto fmterr;
459     if (scanf("site_elevation %lf\n", &elevation) != 1)
460     goto fmterr;
461     if (scanf("weather_data_file_units %d\n", &input) != 1)
462     goto fmterr;
463     switch (input) { /* translate units */
464     case 1:
465     input = 1; /* radiometric quantities */
466     dir_is_horiz = 0; /* direct is perpendicular meas. */
467     break;
468     case 2:
469     input = 1; /* radiometric quantities */
470     dir_is_horiz = 1; /* solar measured horizontally */
471     break;
472     case 3:
473     input = 2; /* photometric quantities */
474     dir_is_horiz = 0; /* direct is perpendicular meas. */
475     break;
476     default:
477     goto fmterr;
478     }
479     rh_init(); /* initialize sky patches */
480     if (verbose) {
481     fprintf(stderr, "%s: location '%s'\n", progname, buf);
482     fprintf(stderr, "%s: (lat,long)=(%.1f,%.1f) degrees north, west\n",
483     progname, s_latitude, s_longitude);
484 greg 2.35 if (avgSky >= 0)
485     fprintf(stderr, "%s: %d sky patches\n",
486     progname, nskypatch);
487     if (sunsfp)
488 greg 2.39 fprintf(stderr, "%s: outputting suns to %s\n",
489     progname, sunsfp==stdout ? "stdout" : "file");
490 greg 2.8 if (rotation != 0)
491     fprintf(stderr, "%s: rotating output %.0f degrees\n",
492     progname, rotation);
493 greg 2.1 }
494 greg 2.2 /* convert quantities to radians */
495     s_latitude = DegToRad(s_latitude);
496     s_longitude = DegToRad(s_longitude);
497     s_meridian = DegToRad(s_meridian);
498 greg 2.27 /* initial allocation */
499 greg 2.28 mtx_data = resize_dmatrix(mtx_data, tstorage=2, nskypatch);
500 greg 2.1 /* process each time step in tape */
501     while (scanf("%d %d %lf %lf %lf\n", &mo, &da, &hr, &dir, &dif) == 5) {
502     double sda, sta;
503 greg 2.37 int sun_in_sky;
504 greg 2.1 /* compute solar position */
505 greg 2.36 if ((mo == 2) & (da == 29)) {
506     julian_date = 60;
507     leap_day = 1;
508     } else
509     julian_date = jdate(mo, da) + leap_day;
510 greg 2.1 sda = sdec(julian_date);
511     sta = stadj(julian_date);
512     altitude = salt(sda, hr+sta);
513 greg 2.37 sun_in_sky = (altitude > -DegToRad(SUN_ANG_DEG/2.));
514     if (sun_hours_only && !sun_in_sky)
515     continue; /* skipping nighttime points */
516 greg 2.8 azimuth = sazi(sda, hr+sta) + PI - DegToRad(rotation);
517 greg 2.36
518 greg 2.37 mtx_offset = 3*nskypatch*nstored;
519     nstored += !avgSky | !nstored;
520     /* make space for next row */
521     if (nstored > tstorage) {
522     tstorage += (tstorage>>1) + nstored + 7;
523     mtx_data = resize_dmatrix(mtx_data, tstorage, nskypatch);
524     }
525     ntsteps++; /* keep count of time steps */
526    
527 greg 2.36 if (dir+dif <= 1e-4) { /* effectively nighttime? */
528     if (!avgSky | !mtx_offset)
529     memset(mtx_data+mtx_offset, 0,
530     sizeof(float)*3*nskypatch);
531 greg 2.37 /* output black sun? */
532     if (sunsfp && sun_in_sky)
533 greg 2.36 OutputSun(solar_minute(julian_date,hr), 0,
534     sunsfp, modsfp);
535     continue;
536     }
537 greg 2.38 if (!sun_in_sky && dir > (input==1 ? 20. : 20.*WHTEFFICACY))
538     fprintf(stderr,
539     "%s: warning - unusually bright at %.1f on %d-%d\n",
540     progname, hr, mo, da);
541 greg 2.1 /* convert measured values */
542 greg 2.38 if (dir_is_horiz && altitude > FTINY)
543 greg 2.1 dir /= sin(altitude);
544     if (input == 1) {
545     dir_irrad = dir;
546     diff_irrad = dif;
547     } else /* input == 2 */ {
548     dir_illum = dir;
549     diff_illum = dif;
550     }
551     /* compute sky patch values */
552     ComputeSky(mtx_data+mtx_offset);
553 greg 2.37 /* output sun if requested */
554     if (sunsfp && sun_in_sky)
555 greg 2.36 OutputSun(solar_minute(julian_date,hr), 1,
556     sunsfp, modsfp);
557 greg 2.35
558 greg 2.34 if (avgSky < 0) /* no matrix? */
559     continue;
560    
561 greg 2.4 AddDirect(mtx_data+mtx_offset);
562 greg 2.27 /* update cumulative sky? */
563     for (i = 3*nskypatch*(avgSky&(ntsteps>1)); i--; )
564     mtx_data[i] += mtx_data[mtx_offset+i];
565 greg 2.34 /* monthly reporting */
566     if (verbose && mo != last_monthly)
567     fprintf(stderr, "%s: stepping through month %d...\n",
568     progname, last_monthly=mo);
569 greg 2.36 /* note whether leap-day was given */
570 greg 2.34 }
571     if (!ntsteps) {
572     fprintf(stderr, "%s: no valid time steps on input\n", progname);
573     exit(1);
574 greg 2.1 }
575     /* check for junk at end */
576     while ((i = fgetc(stdin)) != EOF)
577     if (!isspace(i)) {
578     fprintf(stderr, "%s: warning - unexpected data past EOT: ",
579     progname);
580     buf[0] = i; buf[1] = '\0';
581     fgets(buf+1, sizeof(buf)-1, stdin);
582     fputs(buf, stderr); fputc('\n', stderr);
583     break;
584     }
585 greg 2.34
586     if (avgSky < 0) /* no matrix output? */
587     goto alldone;
588    
589 greg 2.27 dif = 1./(double)ntsteps; /* average sky? */
590     for (i = 3*nskypatch*(avgSky&(ntsteps>1)); i--; )
591     mtx_data[i] *= dif;
592 greg 2.1 /* write out matrix */
593 greg 2.14 if (outfmt != 'a')
594     SET_FILE_BINARY(stdout);
595 greg 2.1 #ifdef getc_unlocked
596     flockfile(stdout);
597     #endif
598     if (verbose)
599     fprintf(stderr, "%s: writing %smatrix with %d time steps...\n",
600 greg 2.28 progname, outfmt=='a' ? "" : "binary ", nstored);
601 greg 2.14 if (doheader) {
602     newheader("RADIANCE", stdout);
603     printargs(argc, argv, stdout);
604     printf("LATLONG= %.8f %.8f\n", RadToDeg(s_latitude),
605     -RadToDeg(s_longitude));
606     printf("NROWS=%d\n", nskypatch);
607 greg 2.28 printf("NCOLS=%d\n", nstored);
608 greg 2.14 printf("NCOMP=3\n");
609 greg 2.29 if ((outfmt == 'f') | (outfmt == 'd'))
610     fputendian(stdout);
611 greg 2.18 fputformat((char *)getfmtname(outfmt), stdout);
612 greg 2.14 putchar('\n');
613     }
614 greg 2.1 /* patches are rows (outer sort) */
615     for (i = 0; i < nskypatch; i++) {
616     mtx_offset = 3*i;
617     switch (outfmt) {
618     case 'a':
619 greg 2.28 for (j = 0; j < nstored; j++) {
620 greg 2.3 printf("%.3g %.3g %.3g\n", mtx_data[mtx_offset],
621 greg 2.1 mtx_data[mtx_offset+1],
622     mtx_data[mtx_offset+2]);
623     mtx_offset += 3*nskypatch;
624     }
625 greg 2.28 if (nstored > 1)
626 greg 2.2 fputc('\n', stdout);
627 greg 2.1 break;
628     case 'f':
629 greg 2.28 for (j = 0; j < nstored; j++) {
630 greg 2.23 putbinary(mtx_data+mtx_offset, sizeof(float), 3,
631 greg 2.1 stdout);
632     mtx_offset += 3*nskypatch;
633     }
634     break;
635     case 'd':
636 greg 2.28 for (j = 0; j < nstored; j++) {
637 greg 2.1 double ment[3];
638     ment[0] = mtx_data[mtx_offset];
639     ment[1] = mtx_data[mtx_offset+1];
640     ment[2] = mtx_data[mtx_offset+2];
641 greg 2.23 putbinary(ment, sizeof(double), 3, stdout);
642 greg 2.1 mtx_offset += 3*nskypatch;
643     }
644     break;
645     }
646     if (ferror(stdout))
647     goto writerr;
648     }
649 greg 2.34 alldone:
650     if (fflush(NULL) == EOF)
651 greg 2.1 goto writerr;
652     if (verbose)
653     fprintf(stderr, "%s: done.\n", progname);
654     exit(0);
655     userr:
656 greg 2.37 fprintf(stderr, "Usage: %s [-v][-h][-A][-d|-s|-n][-u][-D file [-M modfile]][-r deg][-m N][-g r g b][-c r g b][-o{f|d}][-O{0|1}] [tape.wea]\n",
657 greg 2.1 progname);
658     exit(1);
659     fmterr:
660 greg 2.34 fprintf(stderr, "%s: weather tape format error in header\n", progname);
661 greg 2.1 exit(1);
662     writerr:
663     fprintf(stderr, "%s: write error on output\n", progname);
664     exit(1);
665     }
666    
667     /* Return maximum of two doubles */
668     double dmax( double a, double b )
669     { return (a > b) ? a : b; }
670    
671     /* Compute sky patch radiance values (modified by GW) */
672     void
673     ComputeSky(float *parr)
674     {
675     int index; /* Category index */
676     double norm_diff_illum; /* Normalized diffuse illuimnance */
677     int i;
678    
679     /* Calculate atmospheric precipitable water content */
680     apwc = CalcPrecipWater(dew_point);
681    
682 greg 2.6 /* Calculate sun zenith angle (don't let it dip below horizon) */
683     /* Also limit minimum angle to keep circumsolar off zenith */
684     if (altitude <= 0.0)
685     sun_zenith = DegToRad(90.0);
686     else if (altitude >= DegToRad(87.0))
687     sun_zenith = DegToRad(3.0);
688     else
689     sun_zenith = DegToRad(90.0) - altitude;
690 greg 2.1
691     /* Compute the inputs for the calculation of the sky distribution */
692    
693     if (input == 0) /* XXX never used */
694     {
695     /* Calculate irradiance */
696     diff_irrad = CalcDiffuseIrradiance();
697     dir_irrad = CalcDirectIrradiance();
698    
699     /* Calculate illuminance */
700     index = GetCategoryIndex();
701     diff_illum = diff_irrad * CalcDiffuseIllumRatio(index);
702     dir_illum = dir_irrad * CalcDirectIllumRatio(index);
703     }
704     else if (input == 1)
705     {
706     sky_brightness = CalcSkyBrightness();
707     sky_clearness = CalcSkyClearness();
708    
709 greg 2.9 /* Limit sky clearness */
710     if (sky_clearness > 11.9)
711     sky_clearness = 11.9;
712    
713     /* Limit sky brightness */
714     if (sky_brightness < 0.01)
715 greg 2.11 sky_brightness = 0.01;
716 greg 2.9
717 greg 2.1 /* Calculate illuminance */
718     index = GetCategoryIndex();
719     diff_illum = diff_irrad * CalcDiffuseIllumRatio(index);
720     dir_illum = dir_irrad * CalcDirectIllumRatio(index);
721     }
722     else if (input == 2)
723     {
724     /* Calculate sky brightness and clearness from illuminance values */
725     index = CalcSkyParamFromIllum();
726     }
727    
728 greg 2.11 if (output == 1) { /* hack for solar radiance */
729     diff_illum = diff_irrad * WHTEFFICACY;
730     dir_illum = dir_irrad * WHTEFFICACY;
731     }
732 greg 2.1 /* Compute ground radiance (include solar contribution if any) */
733 greg 2.3 parr[0] = diff_illum;
734 greg 2.1 if (altitude > 0)
735 greg 2.3 parr[0] += dir_illum * sin(altitude);
736 greg 2.4 parr[2] = parr[1] = parr[0] *= (1./PI/WHTEFFICACY);
737     multcolor(parr, grefl);
738 greg 2.1
739 greg 2.32 if (bright(skycolor) <= 1e-4) { /* 0 sky component? */
740     memset(parr+3, 0, sizeof(float)*3*(nskypatch-1));
741     return;
742     }
743 greg 2.1 /* Calculate Perez sky model parameters */
744     CalcPerezParam(sun_zenith, sky_clearness, sky_brightness, index);
745    
746     /* Calculate sky patch luminance values */
747     CalcSkyPatchLumin(parr);
748    
749     /* Calculate relative horizontal illuminance */
750     norm_diff_illum = CalcRelHorzIllum(parr);
751    
752 greg 2.13 /* Check for zero sky -- make uniform in that case */
753     if (norm_diff_illum <= FTINY) {
754     for (i = 1; i < nskypatch; i++)
755     setcolor(parr+3*i, 1., 1., 1.);
756     norm_diff_illum = PI;
757     }
758 greg 2.1 /* Normalization coefficient */
759     norm_diff_illum = diff_illum / norm_diff_illum;
760    
761     /* Apply to sky patches to get absolute radiance values */
762     for (i = 1; i < nskypatch; i++) {
763 greg 2.7 scalecolor(parr+3*i, norm_diff_illum*(1./WHTEFFICACY));
764 greg 2.1 multcolor(parr+3*i, skycolor);
765     }
766     }
767    
768     /* Add in solar direct to nearest sky patches (GW) */
769     void
770     AddDirect(float *parr)
771     {
772     FVECT svec;
773 greg 2.3 double near_dprod[NSUNPATCH];
774     int near_patch[NSUNPATCH];
775     double wta[NSUNPATCH], wtot;
776 greg 2.1 int i, j, p;
777    
778 greg 2.4 if (dir_illum <= 1e-4 || bright(suncolor) <= 1e-4)
779 greg 2.1 return;
780 greg 2.10 /* identify nsuns closest patches */
781     if (nsuns > NSUNPATCH)
782     nsuns = NSUNPATCH;
783     else if (nsuns <= 0)
784     nsuns = 1;
785     for (i = nsuns; i--; )
786 greg 2.1 near_dprod[i] = -1.;
787     vector(svec, altitude, azimuth);
788     for (p = 1; p < nskypatch; p++) {
789     FVECT pvec;
790     double dprod;
791     rh_vector(pvec, p);
792     dprod = DOT(pvec, svec);
793 greg 2.10 for (i = 0; i < nsuns; i++)
794 greg 2.1 if (dprod > near_dprod[i]) {
795 greg 2.10 for (j = nsuns; --j > i; ) {
796 greg 2.1 near_dprod[j] = near_dprod[j-1];
797     near_patch[j] = near_patch[j-1];
798     }
799     near_dprod[i] = dprod;
800     near_patch[i] = p;
801     break;
802     }
803     }
804     wtot = 0; /* weight by proximity */
805 greg 2.10 for (i = nsuns; i--; )
806 greg 2.1 wtot += wta[i] = 1./(1.002 - near_dprod[i]);
807     /* add to nearest patch radiances */
808 greg 2.10 for (i = nsuns; i--; ) {
809 greg 2.2 float *pdest = parr + 3*near_patch[i];
810 greg 2.10 float val_add = wta[i] * dir_illum / (WHTEFFICACY * wtot);
811    
812     val_add /= (fixed_sun_sa > 0) ? fixed_sun_sa
813     : rh_dom[near_patch[i]] ;
814 greg 2.4 *pdest++ += val_add*suncolor[0];
815     *pdest++ += val_add*suncolor[1];
816     *pdest++ += val_add*suncolor[2];
817 greg 2.2 }
818 greg 2.1 }
819    
820 greg 2.35 /* Output a sun to indicated file if appropriate for this time step */
821     void
822 greg 2.36 OutputSun(int id, int goodsun, FILE *fp, FILE *mfp)
823 greg 2.35 {
824     double srad;
825     FVECT sv;
826    
827 greg 2.36 srad = DegToRad(SUN_ANG_DEG/2.);
828     srad = goodsun ? dir_illum/(WHTEFFICACY * PI*srad*srad) : 0;
829 greg 2.35 vector(sv, altitude, azimuth);
830     fprintf(fp, "\nvoid light solar%d\n0\n0\n", id);
831     fprintf(fp, "3 %.3e %.3e %.3e\n", srad*suncolor[0],
832     srad*suncolor[1], srad*suncolor[2]);
833     fprintf(fp, "\nsolar%d source sun%d\n0\n0\n", id, id);
834     fprintf(fp, "4 %.6f %.6f %.6f %.4f\n", sv[0], sv[1], sv[2], SUN_ANG_DEG);
835 greg 2.36
836     if (mfp != NULL) /* saving modifier IDs? */
837     fprintf(mfp, "solar%d\n", id);
838 greg 2.35 }
839    
840 greg 2.1 /* Initialize Reinhart sky patch positions (GW) */
841     int
842     rh_init(void)
843     {
844     #define NROW 7
845     static const int tnaz[NROW] = {30, 30, 24, 24, 18, 12, 6};
846     const double alpha = (PI/2.)/(NROW*rhsubdiv + .5);
847     int p, i, j;
848     /* allocate patch angle arrays */
849     nskypatch = 0;
850     for (p = 0; p < NROW; p++)
851     nskypatch += tnaz[p];
852     nskypatch *= rhsubdiv*rhsubdiv;
853     nskypatch += 2;
854     rh_palt = (float *)malloc(sizeof(float)*nskypatch);
855     rh_pazi = (float *)malloc(sizeof(float)*nskypatch);
856     rh_dom = (float *)malloc(sizeof(float)*nskypatch);
857     if ((rh_palt == NULL) | (rh_pazi == NULL) | (rh_dom == NULL)) {
858     fprintf(stderr, "%s: out of memory in rh_init()\n", progname);
859     exit(1);
860     }
861     rh_palt[0] = -PI/2.; /* ground & zenith patches */
862     rh_pazi[0] = 0.;
863     rh_dom[0] = 2.*PI;
864     rh_palt[nskypatch-1] = PI/2.;
865     rh_pazi[nskypatch-1] = 0.;
866     rh_dom[nskypatch-1] = 2.*PI*(1. - cos(alpha*.5));
867     p = 1; /* "normal" patches */
868     for (i = 0; i < NROW*rhsubdiv; i++) {
869     const float ralt = alpha*(i + .5);
870     const int ninrow = tnaz[i/rhsubdiv]*rhsubdiv;
871 greg 2.3 const float dom = 2.*PI*(sin(alpha*(i+1)) - sin(alpha*i)) /
872     (double)ninrow;
873 greg 2.1 for (j = 0; j < ninrow; j++) {
874     rh_palt[p] = ralt;
875     rh_pazi[p] = 2.*PI * j / (double)ninrow;
876     rh_dom[p++] = dom;
877     }
878     }
879     return nskypatch;
880     #undef NROW
881     }
882    
883     /* Resize daylight matrix (GW) */
884     float *
885     resize_dmatrix(float *mtx_data, int nsteps, int npatch)
886     {
887     if (mtx_data == NULL)
888     mtx_data = (float *)malloc(sizeof(float)*3*nsteps*npatch);
889     else
890     mtx_data = (float *)realloc(mtx_data,
891     sizeof(float)*3*nsteps*npatch);
892     if (mtx_data == NULL) {
893     fprintf(stderr, "%s: out of memory in resize_dmatrix(%d,%d)\n",
894     progname, nsteps, npatch);
895     exit(1);
896     }
897     return(mtx_data);
898     }
899    
900     /* Determine category index */
901     int GetCategoryIndex()
902     {
903     int index; /* Loop index */
904    
905     for (index = 0; index < 8; index++)
906     if ((sky_clearness >= SkyClearCat[index].lower) &&
907     (sky_clearness < SkyClearCat[index].upper))
908     break;
909    
910     return index;
911     }
912    
913     /* Calculate diffuse illuminance to diffuse irradiance ratio */
914    
915     /* Reference: Perez, R., P. Ineichen, R. Seals, J. Michalsky, and R. */
916     /* Stewart. 1990. ìModeling Daylight Availability and */
917     /* Irradiance Components from Direct and Global */
918     /* Irradiance,î Solar Energy 44(5):271-289, Eqn. 7. */
919    
920     double CalcDiffuseIllumRatio( int index )
921     {
922     ModelCoeff const *pnle; /* Category coefficient pointer */
923    
924     /* Get category coefficient pointer */
925     pnle = &(DiffuseLumEff[index]);
926    
927     return pnle->a + pnle->b * apwc + pnle->c * cos(sun_zenith) +
928     pnle->d * log(sky_brightness);
929     }
930    
931     /* Calculate direct illuminance to direct irradiance ratio */
932    
933     /* Reference: Perez, R., P. Ineichen, R. Seals, J. Michalsky, and R. */
934     /* Stewart. 1990. ìModeling Daylight Availability and */
935     /* Irradiance Components from Direct and Global */
936     /* Irradiance,î Solar Energy 44(5):271-289, Eqn. 8. */
937    
938     double CalcDirectIllumRatio( int index )
939     {
940     ModelCoeff const *pnle; /* Category coefficient pointer */
941    
942     /* Get category coefficient pointer */
943     pnle = &(DirectLumEff[index]);
944    
945     /* Calculate direct illuminance from direct irradiance */
946    
947     return dmax((pnle->a + pnle->b * apwc + pnle->c * exp(5.73 *
948     sun_zenith - 5.0) + pnle->d * sky_brightness),
949     0.0);
950     }
951    
952     /* Calculate sky brightness */
953    
954     /* Reference: Perez, R., P. Ineichen, R. Seals, J. Michalsky, and R. */
955     /* Stewart. 1990. ìModeling Daylight Availability and */
956     /* Irradiance Components from Direct and Global */
957     /* Irradiance,î Solar Energy 44(5):271-289, Eqn. 2. */
958    
959     double CalcSkyBrightness()
960     {
961     return diff_irrad * CalcAirMass() / (DC_SolarConstantE *
962     CalcEccentricity());
963     }
964    
965     /* Calculate sky clearness */
966    
967     /* Reference: Perez, R., P. Ineichen, R. Seals, J. Michalsky, and R. */
968     /* Stewart. 1990. ìModeling Daylight Availability and */
969     /* Irradiance Components from Direct and Global */
970     /* Irradiance,î Solar Energy 44(5):271-289, Eqn. 1. */
971    
972     double CalcSkyClearness()
973     {
974     double sz_cubed; /* Sun zenith angle cubed */
975    
976     /* Calculate sun zenith angle cubed */
977 greg 2.11 sz_cubed = sun_zenith*sun_zenith*sun_zenith;
978 greg 2.1
979     return ((diff_irrad + dir_irrad) / diff_irrad + 1.041 *
980     sz_cubed) / (1.0 + 1.041 * sz_cubed);
981     }
982    
983     /* Calculate diffuse horizontal irradiance from Perez sky brightness */
984    
985     /* Reference: Perez, R., P. Ineichen, R. Seals, J. Michalsky, and R. */
986     /* Stewart. 1990. ìModeling Daylight Availability and */
987     /* Irradiance Components from Direct and Global */
988     /* Irradiance,î Solar Energy 44(5):271-289, Eqn. 2 */
989     /* (inverse). */
990    
991     double CalcDiffuseIrradiance()
992     {
993     return sky_brightness * DC_SolarConstantE * CalcEccentricity() /
994     CalcAirMass();
995     }
996    
997     /* Calculate direct normal irradiance from Perez sky clearness */
998    
999     /* Reference: Perez, R., P. Ineichen, R. Seals, J. Michalsky, and R. */
1000     /* Stewart. 1990. ìModeling Daylight Availability and */
1001     /* Irradiance Components from Direct and Global */
1002     /* Irradiance,î Solar Energy 44(5):271-289, Eqn. 1 */
1003     /* (inverse). */
1004    
1005     double CalcDirectIrradiance()
1006     {
1007     return CalcDiffuseIrradiance() * ((sky_clearness - 1.0) * (1 + 1.041
1008 greg 2.11 * sun_zenith*sun_zenith*sun_zenith));
1009 greg 2.1 }
1010    
1011     /* Calculate sky brightness and clearness from illuminance values */
1012     int CalcSkyParamFromIllum()
1013     {
1014     double test1 = 0.1;
1015     double test2 = 0.1;
1016     int counter = 0;
1017     int index = 0; /* Category index */
1018    
1019     /* Convert illuminance to irradiance */
1020     diff_irrad = diff_illum * DC_SolarConstantE /
1021     (DC_SolarConstantL * 1000.0);
1022     dir_irrad = dir_illum * DC_SolarConstantE /
1023     (DC_SolarConstantL * 1000.0);
1024    
1025     /* Calculate sky brightness and clearness */
1026     sky_brightness = CalcSkyBrightness();
1027     sky_clearness = CalcSkyClearness();
1028    
1029     /* Limit sky clearness */
1030     if (sky_clearness > 12.0)
1031     sky_clearness = 12.0;
1032    
1033     /* Limit sky brightness */
1034 greg 2.9 if (sky_brightness < 0.01)
1035 greg 2.1 sky_brightness = 0.01;
1036    
1037     while (((fabs(diff_irrad - test1) > 10.0) ||
1038     (fabs(dir_irrad - test2) > 10.0)) && !(counter == 5))
1039     {
1040     test1 = diff_irrad;
1041     test2 = dir_irrad;
1042     counter++;
1043    
1044     /* Convert illuminance to irradiance */
1045     index = GetCategoryIndex();
1046     diff_irrad = diff_illum / CalcDiffuseIllumRatio(index);
1047 greg 2.26 dir_irrad = CalcDirectIllumRatio(index);
1048     if (dir_irrad > 0.1)
1049     dir_irrad = dir_illum / dir_irrad;
1050 greg 2.1
1051     /* Calculate sky brightness and clearness */
1052     sky_brightness = CalcSkyBrightness();
1053     sky_clearness = CalcSkyClearness();
1054    
1055     /* Limit sky clearness */
1056     if (sky_clearness > 12.0)
1057     sky_clearness = 12.0;
1058    
1059     /* Limit sky brightness */
1060 greg 2.9 if (sky_brightness < 0.01)
1061 greg 2.1 sky_brightness = 0.01;
1062     }
1063    
1064     return GetCategoryIndex();
1065     }
1066    
1067     /* Calculate relative luminance */
1068    
1069     /* Reference: Perez, R., R. Seals, and J. Michalsky. 1993. */
1070     /* ìAll-Weather Model for Sky Luminance Distribution - */
1071     /* Preliminary Configuration and Validation,î Solar Energy */
1072     /* 50(3):235-245, Eqn. 1. */
1073    
1074     double CalcRelLuminance( double gamma, double zeta )
1075     {
1076     return (1.0 + perez_param[0] * exp(perez_param[1] / cos(zeta))) *
1077     (1.0 + perez_param[2] * exp(perez_param[3] * gamma) +
1078     perez_param[4] * cos(gamma) * cos(gamma));
1079     }
1080    
1081     /* Calculate Perez sky model parameters */
1082    
1083     /* Reference: Perez, R., R. Seals, and J. Michalsky. 1993. */
1084     /* ìAll-Weather Model for Sky Luminance Distribution - */
1085     /* Preliminary Configuration and Validation,î Solar Energy */
1086     /* 50(3):235-245, Eqns. 6 - 8. */
1087    
1088     void CalcPerezParam( double sz, double epsilon, double delta,
1089     int index )
1090     {
1091     double x[5][4]; /* Coefficents a, b, c, d, e */
1092     int i, j; /* Loop indices */
1093    
1094     /* Limit sky brightness */
1095     if (epsilon > 1.065 && epsilon < 2.8)
1096     {
1097     if (delta < 0.2)
1098     delta = 0.2;
1099     }
1100    
1101     /* Get Perez coefficients */
1102     for (i = 0; i < 5; i++)
1103     for (j = 0; j < 4; j++)
1104     x[i][j] = PerezCoeff[index][4 * i + j];
1105    
1106     if (index != 0)
1107     {
1108     /* Calculate parameter a, b, c, d and e (Eqn. 6) */
1109     for (i = 0; i < 5; i++)
1110     perez_param[i] = x[i][0] + x[i][1] * sz + delta * (x[i][2] +
1111     x[i][3] * sz);
1112     }
1113     else
1114     {
1115     /* Parameters a, b and e (Eqn. 6) */
1116     perez_param[0] = x[0][0] + x[0][1] * sz + delta * (x[0][2] +
1117     x[0][3] * sz);
1118     perez_param[1] = x[1][0] + x[1][1] * sz + delta * (x[1][2] +
1119     x[1][3] * sz);
1120     perez_param[4] = x[4][0] + x[4][1] * sz + delta * (x[4][2] +
1121     x[4][3] * sz);
1122    
1123     /* Parameter c (Eqn. 7) */
1124     perez_param[2] = exp(pow(delta * (x[2][0] + x[2][1] * sz),
1125     x[2][2])) - x[2][3];
1126    
1127     /* Parameter d (Eqn. 8) */
1128     perez_param[3] = -exp(delta * (x[3][0] + x[3][1] * sz)) +
1129     x[3][2] + delta * x[3][3];
1130     }
1131     }
1132    
1133     /* Calculate relative horizontal illuminance (modified by GW) */
1134    
1135     /* Reference: Perez, R., R. Seals, and J. Michalsky. 1993. */
1136     /* ìAll-Weather Model for Sky Luminance Distribution - */
1137     /* Preliminary Configuration and Validation,î Solar Energy */
1138     /* 50(3):235-245, Eqn. 3. */
1139    
1140     double CalcRelHorzIllum( float *parr )
1141     {
1142     int i;
1143     double rh_illum = 0.0; /* Relative horizontal illuminance */
1144    
1145     for (i = 1; i < nskypatch; i++)
1146 greg 2.7 rh_illum += parr[3*i+1] * rh_cos(i) * rh_dom[i];
1147 greg 2.1
1148 greg 2.7 return rh_illum;
1149 greg 2.1 }
1150    
1151     /* Calculate earth orbit eccentricity correction factor */
1152    
1153     /* Reference: Sen, Z. 2008. Solar Energy Fundamental and Modeling */
1154     /* Techniques. Springer, p. 72. */
1155    
1156     double CalcEccentricity()
1157     {
1158     double day_angle; /* Day angle (radians) */
1159     double E0; /* Eccentricity */
1160    
1161     /* Calculate day angle */
1162     day_angle = (julian_date - 1.0) * (2.0 * PI / 365.0);
1163    
1164     /* Calculate eccentricity */
1165     E0 = 1.00011 + 0.034221 * cos(day_angle) + 0.00128 * sin(day_angle)
1166     + 0.000719 * cos(2.0 * day_angle) + 0.000077 * sin(2.0 *
1167     day_angle);
1168    
1169     return E0;
1170     }
1171    
1172     /* Calculate atmospheric precipitable water content */
1173    
1174     /* Reference: Perez, R., P. Ineichen, R. Seals, J. Michalsky, and R. */
1175     /* Stewart. 1990. ìModeling Daylight Availability and */
1176     /* Irradiance Components from Direct and Global */
1177     /* Irradiance,î Solar Energy 44(5):271-289, Eqn. 3. */
1178    
1179     /* Note: The default surface dew point temperature is 11 deg. C */
1180     /* (52 deg. F). Typical values are: */
1181    
1182     /* Celsius Fahrenheit Human Perception */
1183     /* > 24 > 75 Extremely uncomfortable */
1184     /* 21 - 24 70 - 74 Very humid */
1185     /* 18 - 21 65 - 69 Somewhat uncomfortable */
1186     /* 16 - 18 60 - 64 OK for most people */
1187     /* 13 - 16 55 - 59 Comfortable */
1188     /* 10 - 12 50 - 54 Very comfortable */
1189     /* < 10 < 49 A bit dry for some */
1190    
1191     double CalcPrecipWater( double dpt )
1192     { return exp(0.07 * dpt - 0.075); }
1193    
1194     /* Calculate relative air mass */
1195    
1196     /* Reference: Kasten, F. 1966. "A New Table and Approximation Formula */
1197     /* for the Relative Optical Air Mass," Arch. Meteorol. */
1198     /* Geophys. Bioklimataol. Ser. B14, pp. 206-233. */
1199    
1200     /* Note: More sophisticated relative air mass models are */
1201     /* available, but they differ significantly only for */
1202     /* sun zenith angles greater than 80 degrees. */
1203    
1204     double CalcAirMass()
1205     {
1206     return (1.0 / (cos(sun_zenith) + 0.15 * pow(93.885 -
1207     RadToDeg(sun_zenith), -1.253)));
1208     }
1209    
1210     /* Calculate Perez All-Weather sky patch luminances (modified by GW) */
1211    
1212     /* NOTE: The sky patches centers are determined in accordance with the */
1213     /* BRE-IDMP sky luminance measurement procedures. (See for example */
1214     /* Mardaljevic, J. 2001. "The BRE-IDMP Dataset: A New Benchmark */
1215     /* for the Validation of Illuminance Prediction Techniques," */
1216     /* Lighting Research & Technology 33(2):117-136.) */
1217    
1218     void CalcSkyPatchLumin( float *parr )
1219     {
1220     int i;
1221     double aas; /* Sun-sky point azimuthal angle */
1222     double sspa; /* Sun-sky point angle */
1223     double zsa; /* Zenithal sun angle */
1224    
1225     for (i = 1; i < nskypatch; i++)
1226     {
1227     /* Calculate sun-sky point azimuthal angle */
1228     aas = fabs(rh_pazi[i] - azimuth);
1229    
1230     /* Calculate zenithal sun angle */
1231     zsa = PI * 0.5 - rh_palt[i];
1232    
1233     /* Calculate sun-sky point angle (Equation 8-20) */
1234     sspa = acos(cos(sun_zenith) * cos(zsa) + sin(sun_zenith) *
1235     sin(zsa) * cos(aas));
1236    
1237     /* Calculate patch luminance */
1238     parr[3*i] = CalcRelLuminance(sspa, zsa);
1239     if (parr[3*i] < 0) parr[3*i] = 0;
1240     parr[3*i+2] = parr[3*i+1] = parr[3*i];
1241     }
1242     }