ViewVC Help
View File | Revision Log | Show Annotations | Download File | Root Listing
root/radiance/ray/src/gen/gendaymtx.c
Revision: 2.34
Committed: Tue Jan 7 01:42:30 2020 UTC (4 years, 3 months ago) by greg
Content type: text/plain
Branch: MAIN
Changes since 2.33: +46 -11 lines
Log Message:
Added -n and -D options to gendaymtx

File Contents

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