ViewVC Help
View File | Revision Log | Show Annotations | Download File | Root Listing
root/radiance/ray/src/gen/gendaymtx.c
Revision: 2.36
Committed: Mon Apr 13 17:12:19 2020 UTC (4 years ago) by greg
Content type: text/plain
Branch: MAIN
Changes since 2.35: +55 -21 lines
Log Message:
Updates to gendaymtx requested by Ladybug Tools

File Contents

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