ViewVC Help
View File | Revision Log | Show Annotations | Download File | Root Listing
root/radiance/ray/src/gen/gendaymtx.c
Revision: 2.22
Committed: Fri Oct 30 17:06:34 2015 UTC (8 years, 6 months ago) by greg
Content type: text/plain
Branch: MAIN
Changes since 2.21: +12 -3 lines
Log Message:
Added warning about inconsistent option settings

File Contents

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