ViewVC Help
View File | Revision Log | Show Annotations | Download File | Root Listing
root/radiance/ray/src/gen/gendaymtx.c
Revision: 2.35
Committed: Tue Jan 7 18:26:55 2020 UTC (4 years, 3 months ago) by greg
Content type: text/plain
Branch: MAIN
Changes since 2.34: +36 -18 lines
Log Message:
Code clean-up from last change

File Contents

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