ViewVC Help
View File | Revision Log | Show Annotations | Download File | Root Listing
root/radiance/ray/src/gen/gendaymtx.c
Revision: 2.16
Committed: Tue Jun 17 21:01:21 2014 UTC (9 years, 10 months ago) by greg
Content type: text/plain
Branch: MAIN
CVS Tags: rad4R2, rad4R2P1
Changes since 2.15: +3 -3 lines
Log Message:
Fixed operator precedence error

File Contents

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