ViewVC Help
View File | Revision Log | Show Annotations | Download File | Root Listing
root/radiance/ray/src/gen/gendaymtx.c
Revision: 2.14
Committed: Fri May 30 00:00:54 2014 UTC (9 years, 10 months ago) by greg
Content type: text/plain
Branch: MAIN
Changes since 2.13: +36 -2 lines
Log Message:
Added NROWS, NCOLS and NCOMP to matrix headers

File Contents

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