ViewVC Help
View File | Revision Log | Show Annotations | Download File | Root Listing
root/radiance/ray/src/gen/gendaymtx.c
Revision: 2.40
Committed: Fri Apr 26 23:10:59 2024 UTC (6 days, 23 hours ago) by greg
Content type: text/plain
Branch: MAIN
CVS Tags: HEAD
Changes since 2.39: +63 -3 lines
Error occurred while calculating annotation data.
Log Message:
fix(gendaymtx): Added -i option and improved consistency with gendaylit, thanks to Yongqing

File Contents

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