ViewVC Help
View File | Revision Log | Show Annotations | Download File | Root Listing
root/radiance/ray/src/gen/gendaymtx.c
Revision: 2.33
Committed: Mon Jan 6 21:22:46 2020 UTC (4 years, 3 months ago) by greg
Content type: text/plain
Branch: MAIN
Changes since 2.32: +2 -14 lines
Log Message:
Made sure -d option works as before

File Contents

# Content
1 #ifndef lint
2 static const char RCSid[] = "$Id: gendaymtx.c,v 2.32 2020/01/06 21:02:57 greg Exp $";
3 #endif
4 /*
5 * gendaymtx.c
6 *
7 * Generate a daylight matrix based on Perez Sky Model.
8 *
9 * Most of this code is borrowed (see copyright below) from Ian Ashdown's
10 * excellent re-implementation of Jean-Jacques Delaunay's gendaylit.c
11 *
12 * Created by Greg Ward on 1/16/13.
13 */
14
15 /*********************************************************************
16 *
17 * H32_gendaylit.CPP - Perez Sky Model Calculation
18 *
19 * Version: 1.00A
20 *
21 * History: 09/10/01 - Created.
22 * 11/10/08 - Modified for Unix compilation.
23 * 11/10/12 - Fixed conditional __max directive.
24 * 1/11/13 - Tweaks and optimizations (G.Ward)
25 *
26 * Compilers: Microsoft Visual C/C++ Professional V10.0
27 *
28 * Author: Ian Ashdown, P.Eng.
29 * byHeart Consultants Limited
30 * 620 Ballantree Road
31 * West Vancouver, B.C.
32 * Canada V7S 1W3
33 * e-mail: [email protected]
34 *
35 * References: Perez, R., P. Ineichen, R. Seals, J. Michalsky, and R.
36 * Stewart. 1990. ìModeling Daylight Availability and
37 * Irradiance Components from Direct and Global
38 * Irradiance,î Solar Energy 44(5):271-289.
39 *
40 * Perez, R., R. Seals, and J. Michalsky. 1993.
41 * ìAll-Weather Model for Sky Luminance Distribution -
42 * Preliminary Configuration and Validation,î Solar Energy
43 * 50(3):235-245.
44 *
45 * Perez, R., R. Seals, and J. Michalsky. 1993. "ERRATUM to
46 * All-Weather Model for Sky Luminance Distribution -
47 * Preliminary Configuration and Validation,î Solar Energy
48 * 51(5):423.
49 *
50 * NOTE: This program is a completely rewritten version of
51 * gendaylit.c written by Jean-Jacques Delaunay (1994).
52 *
53 * Copyright 2009-2012 byHeart Consultants Limited. All rights
54 * reserved.
55 *
56 * Redistribution and use in source and binary forms, with or without
57 * modification, are permitted for personal and commercial purposes
58 * provided that redistribution of source code must retain the above
59 * copyright notice, this list of conditions and the following
60 * disclaimer:
61 *
62 * THIS SOFTWARE IS PROVIDED "AS IS" AND ANY EXPRESSED OR IMPLIED
63 * WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED WARRANTIES
64 * OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE
65 * DISCLAIMED. IN NO EVENT SHALL byHeart Consultants Limited OR
66 * ITS EMPLOYEES BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL,
67 * SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT
68 * LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF
69 * USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED
70 * AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
71 * LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
72 * ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
73 * POSSIBILITY OF SUCH DAMAGE.
74 *
75 *********************************************************************/
76
77 /* Zenith is along the Z-axis */
78 /* X-axis points east */
79 /* Y-axis points north */
80 /* azimuth is measured as degrees or radians east of North */
81
82 /* Include files */
83 #define _USE_MATH_DEFINES
84 #include <stdlib.h>
85 #include <ctype.h>
86 #include "platform.h"
87 #include "rtmath.h"
88 #include "rtio.h"
89 #include "color.h"
90 #include "sun.h"
91
92 char *progname; /* Program name */
93 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 int nsuns = NSUNPATCH; /* number of sun patches to use */
256 double fixed_sun_sa = -1; /* fixed solid angle per sun? */
257
258 int verbose = 0; /* progress reports to stderr? */
259
260 int outfmt = 'a'; /* output format */
261
262 int rhsubdiv = 1; /* Reinhart sky subdivisions */
263
264 COLOR skycolor = {.96, 1.004, 1.118}; /* sky coloration */
265 COLOR suncolor = {1., 1., 1.}; /* sun color */
266 COLOR grefl = {.2, .2, .2}; /* ground reflectance */
267
268 int nskypatch; /* number of Reinhart patches */
269 float *rh_palt; /* sky patch altitudes (radians) */
270 float *rh_pazi; /* sky patch azimuths (radians) */
271 float *rh_dom; /* sky patch solid angle (sr) */
272
273 #define vector(v,alt,azi) ( (v)[1] = tcos(alt), \
274 (v)[0] = (v)[1]*tsin(azi), \
275 (v)[1] *= tcos(azi), \
276 (v)[2] = tsin(alt) )
277
278 #define rh_vector(v,i) vector(v,rh_palt[i],rh_pazi[i])
279
280 #define rh_cos(i) tsin(rh_palt[i])
281
282 extern int rh_init(void);
283 extern float * resize_dmatrix(float *mtx_data, int nsteps, int npatch);
284 extern void AddDirect(float *parr);
285
286
287 static const char *
288 getfmtname(int fmt)
289 {
290 switch (fmt) {
291 case 'a':
292 return("ascii");
293 case 'f':
294 return("float");
295 case 'd':
296 return("double");
297 }
298 return("unknown");
299 }
300
301
302 int
303 main(int argc, char *argv[])
304 {
305 char buf[256];
306 int doheader = 1; /* output header? */
307 double rotation = 0; /* site rotation (degrees) */
308 double elevation; /* site elevation (meters) */
309 int dir_is_horiz; /* direct is meas. on horizontal? */
310 float *mtx_data = NULL; /* our matrix data */
311 int avgSky = 0; /* compute average sky r.t. matrix? */
312 int ntsteps = 0; /* number of time steps */
313 int tstorage = 0; /* number of allocated time steps */
314 int nstored = 0; /* number of time steps in matrix */
315 int last_monthly = 0; /* month of last report */
316 int mo, da; /* month (1-12) and day (1-31) */
317 double hr; /* hour (local standard time) */
318 double dir, dif; /* direct and diffuse values */
319 int mtx_offset;
320 int i, j;
321
322 progname = argv[0];
323 /* get options */
324 for (i = 1; i < argc && argv[i][0] == '-'; i++)
325 switch (argv[i][1]) {
326 case 'g': /* ground reflectance */
327 grefl[0] = atof(argv[++i]);
328 grefl[1] = atof(argv[++i]);
329 grefl[2] = atof(argv[++i]);
330 break;
331 case 'v': /* verbose progress reports */
332 verbose++;
333 break;
334 case 'h': /* turn off header */
335 doheader = 0;
336 break;
337 case 'o': /* output format */
338 switch (argv[i][2]) {
339 case 'f':
340 case 'd':
341 case 'a':
342 outfmt = argv[i][2];
343 break;
344 default:
345 goto userr;
346 }
347 break;
348 case 'O': /* output type */
349 switch (argv[i][2]) {
350 case '0':
351 output = 0;
352 break;
353 case '1':
354 output = 1;
355 break;
356 default:
357 goto userr;
358 }
359 if (argv[i][3])
360 goto userr;
361 break;
362 case 'm': /* Reinhart subdivisions */
363 rhsubdiv = atoi(argv[++i]);
364 break;
365 case 'c': /* sky color */
366 skycolor[0] = atof(argv[++i]);
367 skycolor[1] = atof(argv[++i]);
368 skycolor[2] = atof(argv[++i]);
369 break;
370 case 'd': /* solar (direct) only */
371 skycolor[0] = skycolor[1] = skycolor[2] = 0;
372 grefl[0] = grefl[1] = grefl[2] = 0;
373 break;
374 case 's': /* sky only (no direct) */
375 suncolor[0] = suncolor[1] = suncolor[2] = 0;
376 break;
377 case 'r': /* rotate distribution */
378 if (argv[i][2] && argv[i][2] != 'z')
379 goto userr;
380 rotation = atof(argv[++i]);
381 break;
382 case '5': /* 5-phase calculation */
383 nsuns = 1;
384 fixed_sun_sa = PI/360.*atof(argv[++i]);
385 if (fixed_sun_sa <= 0) {
386 fprintf(stderr, "%s: missing solar disk size argument for '-5' option\n",
387 argv[0]);
388 exit(1);
389 }
390 fixed_sun_sa *= fixed_sun_sa*PI;
391 break;
392 case 'A': /* compute average sky */
393 avgSky = 1;
394 break;
395 default:
396 goto userr;
397 }
398 if (i < argc-1)
399 goto userr;
400 if (i == argc-1 && freopen(argv[i], "r", stdin) == NULL) {
401 fprintf(stderr, "%s: cannot open '%s' for input\n",
402 progname, argv[i]);
403 exit(1);
404 }
405 if (verbose) {
406 if (i == argc-1)
407 fprintf(stderr, "%s: reading weather tape '%s'\n",
408 progname, argv[i]);
409 else
410 fprintf(stderr, "%s: reading weather tape from <stdin>\n",
411 progname);
412 }
413 /* read weather tape header */
414 if (scanf("place %[^\r\n] ", buf) != 1)
415 goto fmterr;
416 if (scanf("latitude %lf\n", &s_latitude) != 1)
417 goto fmterr;
418 if (scanf("longitude %lf\n", &s_longitude) != 1)
419 goto fmterr;
420 if (scanf("time_zone %lf\n", &s_meridian) != 1)
421 goto fmterr;
422 if (scanf("site_elevation %lf\n", &elevation) != 1)
423 goto fmterr;
424 if (scanf("weather_data_file_units %d\n", &input) != 1)
425 goto fmterr;
426 switch (input) { /* translate units */
427 case 1:
428 input = 1; /* radiometric quantities */
429 dir_is_horiz = 0; /* direct is perpendicular meas. */
430 break;
431 case 2:
432 input = 1; /* radiometric quantities */
433 dir_is_horiz = 1; /* solar measured horizontally */
434 break;
435 case 3:
436 input = 2; /* photometric quantities */
437 dir_is_horiz = 0; /* direct is perpendicular meas. */
438 break;
439 default:
440 goto fmterr;
441 }
442 rh_init(); /* initialize sky patches */
443 if (verbose) {
444 fprintf(stderr, "%s: location '%s'\n", progname, buf);
445 fprintf(stderr, "%s: (lat,long)=(%.1f,%.1f) degrees north, west\n",
446 progname, s_latitude, s_longitude);
447 fprintf(stderr, "%s: %d sky patches per time step\n",
448 progname, nskypatch);
449 if (rotation != 0)
450 fprintf(stderr, "%s: rotating output %.0f degrees\n",
451 progname, rotation);
452 }
453 /* convert quantities to radians */
454 s_latitude = DegToRad(s_latitude);
455 s_longitude = DegToRad(s_longitude);
456 s_meridian = DegToRad(s_meridian);
457 /* initial allocation */
458 mtx_data = resize_dmatrix(mtx_data, tstorage=2, nskypatch);
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
463 mtx_offset = 3*nskypatch*nstored;
464 nstored += !avgSky | !nstored;
465 /* make space for next row */
466 if (nstored > tstorage) {
467 tstorage += (tstorage>>1) + nstored + 7;
468 mtx_data = resize_dmatrix(mtx_data, tstorage, nskypatch);
469 }
470 ntsteps++; /* keep count of time steps */
471 if (dif <= 1e-4) {
472 if (!avgSky | !mtx_offset)
473 memset(mtx_data+mtx_offset, 0, sizeof(float)*3*nskypatch);
474 continue;
475 }
476 if (verbose && mo != last_monthly)
477 fprintf(stderr, "%s: stepping through month %d...\n",
478 progname, last_monthly=mo);
479 /* compute solar position */
480 julian_date = jdate(mo, da);
481 sda = sdec(julian_date);
482 sta = stadj(julian_date);
483 altitude = salt(sda, hr+sta);
484 azimuth = sazi(sda, hr+sta) + PI - DegToRad(rotation);
485 /* convert measured values */
486 if (dir_is_horiz && altitude > 0.)
487 dir /= sin(altitude);
488 if (input == 1) {
489 dir_irrad = dir;
490 diff_irrad = dif;
491 } else /* input == 2 */ {
492 dir_illum = dir;
493 diff_illum = dif;
494 }
495 /* compute sky patch values */
496 ComputeSky(mtx_data+mtx_offset);
497 AddDirect(mtx_data+mtx_offset);
498 /* update cumulative sky? */
499 for (i = 3*nskypatch*(avgSky&(ntsteps>1)); i--; )
500 mtx_data[i] += mtx_data[mtx_offset+i];
501 }
502 /* check for junk at end */
503 while ((i = fgetc(stdin)) != EOF)
504 if (!isspace(i)) {
505 fprintf(stderr, "%s: warning - unexpected data past EOT: ",
506 progname);
507 buf[0] = i; buf[1] = '\0';
508 fgets(buf+1, sizeof(buf)-1, stdin);
509 fputs(buf, stderr); fputc('\n', stderr);
510 break;
511 }
512 if (!ntsteps) {
513 fprintf(stderr, "%s: no valid time steps on input\n", progname);
514 exit(1);
515 }
516 dif = 1./(double)ntsteps; /* average sky? */
517 for (i = 3*nskypatch*(avgSky&(ntsteps>1)); i--; )
518 mtx_data[i] *= dif;
519 /* write out matrix */
520 if (outfmt != 'a')
521 SET_FILE_BINARY(stdout);
522 #ifdef getc_unlocked
523 flockfile(stdout);
524 #endif
525 if (verbose)
526 fprintf(stderr, "%s: writing %smatrix with %d time steps...\n",
527 progname, outfmt=='a' ? "" : "binary ", nstored);
528 if (doheader) {
529 newheader("RADIANCE", stdout);
530 printargs(argc, argv, stdout);
531 printf("LATLONG= %.8f %.8f\n", RadToDeg(s_latitude),
532 -RadToDeg(s_longitude));
533 printf("NROWS=%d\n", nskypatch);
534 printf("NCOLS=%d\n", nstored);
535 printf("NCOMP=3\n");
536 if ((outfmt == 'f') | (outfmt == 'd'))
537 fputendian(stdout);
538 fputformat((char *)getfmtname(outfmt), stdout);
539 putchar('\n');
540 }
541 /* patches are rows (outer sort) */
542 for (i = 0; i < nskypatch; i++) {
543 mtx_offset = 3*i;
544 switch (outfmt) {
545 case 'a':
546 for (j = 0; j < nstored; j++) {
547 printf("%.3g %.3g %.3g\n", mtx_data[mtx_offset],
548 mtx_data[mtx_offset+1],
549 mtx_data[mtx_offset+2]);
550 mtx_offset += 3*nskypatch;
551 }
552 if (nstored > 1)
553 fputc('\n', stdout);
554 break;
555 case 'f':
556 for (j = 0; j < nstored; j++) {
557 putbinary(mtx_data+mtx_offset, sizeof(float), 3,
558 stdout);
559 mtx_offset += 3*nskypatch;
560 }
561 break;
562 case 'd':
563 for (j = 0; j < nstored; j++) {
564 double ment[3];
565 ment[0] = mtx_data[mtx_offset];
566 ment[1] = mtx_data[mtx_offset+1];
567 ment[2] = mtx_data[mtx_offset+2];
568 putbinary(ment, sizeof(double), 3, stdout);
569 mtx_offset += 3*nskypatch;
570 }
571 break;
572 }
573 if (ferror(stdout))
574 goto writerr;
575 }
576 if (fflush(stdout) == EOF)
577 goto writerr;
578 if (verbose)
579 fprintf(stderr, "%s: done.\n", progname);
580 exit(0);
581 userr:
582 fprintf(stderr, "Usage: %s [-v][-h][-A][-d|-s][-r deg][-m N][-g r g b][-c r g b][-o{f|d}][-O{0|1}] [tape.wea]\n",
583 progname);
584 exit(1);
585 fmterr:
586 fprintf(stderr, "%s: input weather tape format error\n", progname);
587 exit(1);
588 writerr:
589 fprintf(stderr, "%s: write error on output\n", progname);
590 exit(1);
591 }
592
593 /* Return maximum of two doubles */
594 double dmax( double a, double b )
595 { return (a > b) ? a : b; }
596
597 /* Compute sky patch radiance values (modified by GW) */
598 void
599 ComputeSky(float *parr)
600 {
601 int index; /* Category index */
602 double norm_diff_illum; /* Normalized diffuse illuimnance */
603 int i;
604
605 /* Calculate atmospheric precipitable water content */
606 apwc = CalcPrecipWater(dew_point);
607
608 /* Calculate sun zenith angle (don't let it dip below horizon) */
609 /* Also limit minimum angle to keep circumsolar off zenith */
610 if (altitude <= 0.0)
611 sun_zenith = DegToRad(90.0);
612 else if (altitude >= DegToRad(87.0))
613 sun_zenith = DegToRad(3.0);
614 else
615 sun_zenith = DegToRad(90.0) - altitude;
616
617 /* Compute the inputs for the calculation of the sky distribution */
618
619 if (input == 0) /* XXX never used */
620 {
621 /* Calculate irradiance */
622 diff_irrad = CalcDiffuseIrradiance();
623 dir_irrad = CalcDirectIrradiance();
624
625 /* Calculate illuminance */
626 index = GetCategoryIndex();
627 diff_illum = diff_irrad * CalcDiffuseIllumRatio(index);
628 dir_illum = dir_irrad * CalcDirectIllumRatio(index);
629 }
630 else if (input == 1)
631 {
632 sky_brightness = CalcSkyBrightness();
633 sky_clearness = CalcSkyClearness();
634
635 /* Limit sky clearness */
636 if (sky_clearness > 11.9)
637 sky_clearness = 11.9;
638
639 /* Limit sky brightness */
640 if (sky_brightness < 0.01)
641 sky_brightness = 0.01;
642
643 /* Calculate illuminance */
644 index = GetCategoryIndex();
645 diff_illum = diff_irrad * CalcDiffuseIllumRatio(index);
646 dir_illum = dir_irrad * CalcDirectIllumRatio(index);
647 }
648 else if (input == 2)
649 {
650 /* Calculate sky brightness and clearness from illuminance values */
651 index = CalcSkyParamFromIllum();
652 }
653
654 if (output == 1) { /* hack for solar radiance */
655 diff_illum = diff_irrad * WHTEFFICACY;
656 dir_illum = dir_irrad * WHTEFFICACY;
657 }
658 /* Compute ground radiance (include solar contribution if any) */
659 parr[0] = diff_illum;
660 if (altitude > 0)
661 parr[0] += dir_illum * sin(altitude);
662 parr[2] = parr[1] = parr[0] *= (1./PI/WHTEFFICACY);
663 multcolor(parr, grefl);
664
665 if (bright(skycolor) <= 1e-4) { /* 0 sky component? */
666 memset(parr+3, 0, sizeof(float)*3*(nskypatch-1));
667 return;
668 }
669 /* Calculate Perez sky model parameters */
670 CalcPerezParam(sun_zenith, sky_clearness, sky_brightness, index);
671
672 /* Calculate sky patch luminance values */
673 CalcSkyPatchLumin(parr);
674
675 /* Calculate relative horizontal illuminance */
676 norm_diff_illum = CalcRelHorzIllum(parr);
677
678 /* Check for zero sky -- make uniform in that case */
679 if (norm_diff_illum <= FTINY) {
680 for (i = 1; i < nskypatch; i++)
681 setcolor(parr+3*i, 1., 1., 1.);
682 norm_diff_illum = PI;
683 }
684 /* Normalization coefficient */
685 norm_diff_illum = diff_illum / norm_diff_illum;
686
687 /* Apply to sky patches to get absolute radiance values */
688 for (i = 1; i < nskypatch; i++) {
689 scalecolor(parr+3*i, norm_diff_illum*(1./WHTEFFICACY));
690 multcolor(parr+3*i, skycolor);
691 }
692 }
693
694 /* Add in solar direct to nearest sky patches (GW) */
695 void
696 AddDirect(float *parr)
697 {
698 FVECT svec;
699 double near_dprod[NSUNPATCH];
700 int near_patch[NSUNPATCH];
701 double wta[NSUNPATCH], wtot;
702 int i, j, p;
703
704 if (dir_illum <= 1e-4 || bright(suncolor) <= 1e-4)
705 return;
706 /* identify nsuns closest patches */
707 if (nsuns > NSUNPATCH)
708 nsuns = NSUNPATCH;
709 else if (nsuns <= 0)
710 nsuns = 1;
711 for (i = nsuns; i--; )
712 near_dprod[i] = -1.;
713 vector(svec, altitude, azimuth);
714 for (p = 1; p < nskypatch; p++) {
715 FVECT pvec;
716 double dprod;
717 rh_vector(pvec, p);
718 dprod = DOT(pvec, svec);
719 for (i = 0; i < nsuns; i++)
720 if (dprod > near_dprod[i]) {
721 for (j = nsuns; --j > i; ) {
722 near_dprod[j] = near_dprod[j-1];
723 near_patch[j] = near_patch[j-1];
724 }
725 near_dprod[i] = dprod;
726 near_patch[i] = p;
727 break;
728 }
729 }
730 wtot = 0; /* weight by proximity */
731 for (i = nsuns; i--; )
732 wtot += wta[i] = 1./(1.002 - near_dprod[i]);
733 /* add to nearest patch radiances */
734 for (i = nsuns; i--; ) {
735 float *pdest = parr + 3*near_patch[i];
736 float val_add = wta[i] * dir_illum / (WHTEFFICACY * wtot);
737
738 val_add /= (fixed_sun_sa > 0) ? fixed_sun_sa
739 : rh_dom[near_patch[i]] ;
740 *pdest++ += val_add*suncolor[0];
741 *pdest++ += val_add*suncolor[1];
742 *pdest++ += val_add*suncolor[2];
743 }
744 }
745
746 /* Initialize Reinhart sky patch positions (GW) */
747 int
748 rh_init(void)
749 {
750 #define NROW 7
751 static const int tnaz[NROW] = {30, 30, 24, 24, 18, 12, 6};
752 const double alpha = (PI/2.)/(NROW*rhsubdiv + .5);
753 int p, i, j;
754 /* allocate patch angle arrays */
755 nskypatch = 0;
756 for (p = 0; p < NROW; p++)
757 nskypatch += tnaz[p];
758 nskypatch *= rhsubdiv*rhsubdiv;
759 nskypatch += 2;
760 rh_palt = (float *)malloc(sizeof(float)*nskypatch);
761 rh_pazi = (float *)malloc(sizeof(float)*nskypatch);
762 rh_dom = (float *)malloc(sizeof(float)*nskypatch);
763 if ((rh_palt == NULL) | (rh_pazi == NULL) | (rh_dom == NULL)) {
764 fprintf(stderr, "%s: out of memory in rh_init()\n", progname);
765 exit(1);
766 }
767 rh_palt[0] = -PI/2.; /* ground & zenith patches */
768 rh_pazi[0] = 0.;
769 rh_dom[0] = 2.*PI;
770 rh_palt[nskypatch-1] = PI/2.;
771 rh_pazi[nskypatch-1] = 0.;
772 rh_dom[nskypatch-1] = 2.*PI*(1. - cos(alpha*.5));
773 p = 1; /* "normal" patches */
774 for (i = 0; i < NROW*rhsubdiv; i++) {
775 const float ralt = alpha*(i + .5);
776 const int ninrow = tnaz[i/rhsubdiv]*rhsubdiv;
777 const float dom = 2.*PI*(sin(alpha*(i+1)) - sin(alpha*i)) /
778 (double)ninrow;
779 for (j = 0; j < ninrow; j++) {
780 rh_palt[p] = ralt;
781 rh_pazi[p] = 2.*PI * j / (double)ninrow;
782 rh_dom[p++] = dom;
783 }
784 }
785 return nskypatch;
786 #undef NROW
787 }
788
789 /* Resize daylight matrix (GW) */
790 float *
791 resize_dmatrix(float *mtx_data, int nsteps, int npatch)
792 {
793 if (mtx_data == NULL)
794 mtx_data = (float *)malloc(sizeof(float)*3*nsteps*npatch);
795 else
796 mtx_data = (float *)realloc(mtx_data,
797 sizeof(float)*3*nsteps*npatch);
798 if (mtx_data == NULL) {
799 fprintf(stderr, "%s: out of memory in resize_dmatrix(%d,%d)\n",
800 progname, nsteps, npatch);
801 exit(1);
802 }
803 return(mtx_data);
804 }
805
806 /* Determine category index */
807 int GetCategoryIndex()
808 {
809 int index; /* Loop index */
810
811 for (index = 0; index < 8; index++)
812 if ((sky_clearness >= SkyClearCat[index].lower) &&
813 (sky_clearness < SkyClearCat[index].upper))
814 break;
815
816 return index;
817 }
818
819 /* Calculate diffuse illuminance to diffuse irradiance ratio */
820
821 /* Reference: Perez, R., P. Ineichen, R. Seals, J. Michalsky, and R. */
822 /* Stewart. 1990. ìModeling Daylight Availability and */
823 /* Irradiance Components from Direct and Global */
824 /* Irradiance,î Solar Energy 44(5):271-289, Eqn. 7. */
825
826 double CalcDiffuseIllumRatio( int index )
827 {
828 ModelCoeff const *pnle; /* Category coefficient pointer */
829
830 /* Get category coefficient pointer */
831 pnle = &(DiffuseLumEff[index]);
832
833 return pnle->a + pnle->b * apwc + pnle->c * cos(sun_zenith) +
834 pnle->d * log(sky_brightness);
835 }
836
837 /* Calculate direct illuminance to direct irradiance ratio */
838
839 /* Reference: Perez, R., P. Ineichen, R. Seals, J. Michalsky, and R. */
840 /* Stewart. 1990. ìModeling Daylight Availability and */
841 /* Irradiance Components from Direct and Global */
842 /* Irradiance,î Solar Energy 44(5):271-289, Eqn. 8. */
843
844 double CalcDirectIllumRatio( int index )
845 {
846 ModelCoeff const *pnle; /* Category coefficient pointer */
847
848 /* Get category coefficient pointer */
849 pnle = &(DirectLumEff[index]);
850
851 /* Calculate direct illuminance from direct irradiance */
852
853 return dmax((pnle->a + pnle->b * apwc + pnle->c * exp(5.73 *
854 sun_zenith - 5.0) + pnle->d * sky_brightness),
855 0.0);
856 }
857
858 /* Calculate sky brightness */
859
860 /* Reference: Perez, R., P. Ineichen, R. Seals, J. Michalsky, and R. */
861 /* Stewart. 1990. ìModeling Daylight Availability and */
862 /* Irradiance Components from Direct and Global */
863 /* Irradiance,î Solar Energy 44(5):271-289, Eqn. 2. */
864
865 double CalcSkyBrightness()
866 {
867 return diff_irrad * CalcAirMass() / (DC_SolarConstantE *
868 CalcEccentricity());
869 }
870
871 /* Calculate sky clearness */
872
873 /* Reference: Perez, R., P. Ineichen, R. Seals, J. Michalsky, and R. */
874 /* Stewart. 1990. ìModeling Daylight Availability and */
875 /* Irradiance Components from Direct and Global */
876 /* Irradiance,î Solar Energy 44(5):271-289, Eqn. 1. */
877
878 double CalcSkyClearness()
879 {
880 double sz_cubed; /* Sun zenith angle cubed */
881
882 /* Calculate sun zenith angle cubed */
883 sz_cubed = sun_zenith*sun_zenith*sun_zenith;
884
885 return ((diff_irrad + dir_irrad) / diff_irrad + 1.041 *
886 sz_cubed) / (1.0 + 1.041 * sz_cubed);
887 }
888
889 /* Calculate diffuse horizontal irradiance from Perez sky brightness */
890
891 /* Reference: Perez, R., P. Ineichen, R. Seals, J. Michalsky, and R. */
892 /* Stewart. 1990. ìModeling Daylight Availability and */
893 /* Irradiance Components from Direct and Global */
894 /* Irradiance,î Solar Energy 44(5):271-289, Eqn. 2 */
895 /* (inverse). */
896
897 double CalcDiffuseIrradiance()
898 {
899 return sky_brightness * DC_SolarConstantE * CalcEccentricity() /
900 CalcAirMass();
901 }
902
903 /* Calculate direct normal irradiance from Perez sky clearness */
904
905 /* Reference: Perez, R., P. Ineichen, R. Seals, J. Michalsky, and R. */
906 /* Stewart. 1990. ìModeling Daylight Availability and */
907 /* Irradiance Components from Direct and Global */
908 /* Irradiance,î Solar Energy 44(5):271-289, Eqn. 1 */
909 /* (inverse). */
910
911 double CalcDirectIrradiance()
912 {
913 return CalcDiffuseIrradiance() * ((sky_clearness - 1.0) * (1 + 1.041
914 * sun_zenith*sun_zenith*sun_zenith));
915 }
916
917 /* Calculate sky brightness and clearness from illuminance values */
918 int CalcSkyParamFromIllum()
919 {
920 double test1 = 0.1;
921 double test2 = 0.1;
922 int counter = 0;
923 int index = 0; /* Category index */
924
925 /* Convert illuminance to irradiance */
926 diff_irrad = diff_illum * DC_SolarConstantE /
927 (DC_SolarConstantL * 1000.0);
928 dir_irrad = dir_illum * DC_SolarConstantE /
929 (DC_SolarConstantL * 1000.0);
930
931 /* Calculate sky brightness and clearness */
932 sky_brightness = CalcSkyBrightness();
933 sky_clearness = CalcSkyClearness();
934
935 /* Limit sky clearness */
936 if (sky_clearness > 12.0)
937 sky_clearness = 12.0;
938
939 /* Limit sky brightness */
940 if (sky_brightness < 0.01)
941 sky_brightness = 0.01;
942
943 while (((fabs(diff_irrad - test1) > 10.0) ||
944 (fabs(dir_irrad - test2) > 10.0)) && !(counter == 5))
945 {
946 test1 = diff_irrad;
947 test2 = dir_irrad;
948 counter++;
949
950 /* Convert illuminance to irradiance */
951 index = GetCategoryIndex();
952 diff_irrad = diff_illum / CalcDiffuseIllumRatio(index);
953 dir_irrad = CalcDirectIllumRatio(index);
954 if (dir_irrad > 0.1)
955 dir_irrad = dir_illum / dir_irrad;
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 }