ViewVC Help
View File | Revision Log | Show Annotations | Download File | Root Listing
root/radiance/ray/src/gen/gendaymtx.c
Revision: 2.37
Committed: Sat Aug 15 03:28:56 2020 UTC (3 years, 8 months ago) by greg
Content type: text/plain
Branch: MAIN
CVS Tags: rad5R3
Changes since 2.36: +23 -19 lines
Log Message:
feat(gendaymtx): added -u option for daylight hours only

File Contents

# Content
1 #ifndef lint
2 static const char RCSid[] = "$Id: gendaymtx.c,v 2.36 2020/04/13 17:12:19 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 /* Perez sky model coefficients */
139
140 /* Reference: Perez, R., R. Seals, and J. Michalsky, 1993. "All- */
141 /* Weather Model for Sky Luminance Distribution - */
142 /* Preliminary Configuration and Validation," Solar */
143 /* Energy 50(3):235-245, Table 1. */
144
145 static const double PerezCoeff[8][20] =
146 {
147 /* Sky clearness (epsilon): 1.000 to 1.065 */
148 { 1.3525, -0.2576, -0.2690, -1.4366, -0.7670,
149 0.0007, 1.2734, -0.1233, 2.8000, 0.6004,
150 1.2375, 1.0000, 1.8734, 0.6297, 0.9738,
151 0.2809, 0.0356, -0.1246, -0.5718, 0.9938 },
152 /* Sky clearness (epsilon): 1.065 to 1.230 */
153 { -1.2219, -0.7730, 1.4148, 1.1016, -0.2054,
154 0.0367, -3.9128, 0.9156, 6.9750, 0.1774,
155 6.4477, -0.1239, -1.5798, -0.5081, -1.7812,
156 0.1080, 0.2624, 0.0672, -0.2190, -0.4285 },
157 /* Sky clearness (epsilon): 1.230 to 1.500 */
158 { -1.1000, -0.2515, 0.8952, 0.0156, 0.2782,
159 -0.1812, - 4.5000, 1.1766, 24.7219, -13.0812,
160 -37.7000, 34.8438, -5.0000, 1.5218, 3.9229,
161 -2.6204, -0.0156, 0.1597, 0.4199, -0.5562 },
162 /* Sky clearness (epsilon): 1.500 to 1.950 */
163 { -0.5484, -0.6654, -0.2672, 0.7117, 0.7234,
164 -0.6219, -5.6812, 2.6297, 33.3389, -18.3000,
165 -62.2500, 52.0781, -3.5000, 0.0016, 1.1477,
166 0.1062, 0.4659, -0.3296, -0.0876, -0.0329 },
167 /* Sky clearness (epsilon): 1.950 to 2.800 */
168 { -0.6000, -0.3566, -2.5000, 2.3250, 0.2937,
169 0.0496, -5.6812, 1.8415, 21.0000, -4.7656 ,
170 -21.5906, 7.2492, -3.5000, -0.1554, 1.4062,
171 0.3988, 0.0032, 0.0766, -0.0656, -0.1294 },
172 /* Sky clearness (epsilon): 2.800 to 4.500 */
173 { -1.0156, -0.3670, 1.0078, 1.4051, 0.2875,
174 -0.5328, -3.8500, 3.3750, 14.0000, -0.9999,
175 -7.1406, 7.5469, -3.4000, -0.1078, -1.0750,
176 1.5702, -0.0672, 0.4016, 0.3017, -0.4844 },
177 /* Sky clearness (epsilon): 4.500 to 6.200 */
178 { -1.0000, 0.0211, 0.5025, -0.5119, -0.3000,
179 0.1922, 0.7023, -1.6317, 19.0000, -5.0000,
180 1.2438, -1.9094, -4.0000, 0.0250, 0.3844,
181 0.2656, 1.0468, -0.3788, -2.4517, 1.4656 },
182 /* Sky clearness (epsilon): 6.200 to ... */
183 { -1.0500, 0.0289, 0.4260, 0.3590, -0.3250,
184 0.1156, 0.7781, 0.0025, 31.0625, -14.5000,
185 -46.1148, 55.3750, -7.2312, 0.4050, 13.3500,
186 0.6234, 1.5000, -0.6426, 1.8564, 0.5636 }
187 };
188
189 /* Perez irradiance component model coefficients */
190
191 /* Reference: Perez, R., P. Ineichen, R. Seals, J. Michalsky, and R. */
192 /* Stewart. 1990. ìModeling Daylight Availability and */
193 /* Irradiance Components from Direct and Global */
194 /* Irradiance,î Solar Energy 44(5):271-289. */
195
196 typedef struct
197 {
198 double lower; /* Lower bound */
199 double upper; /* Upper bound */
200 } CategoryBounds;
201
202 /* Perez sky clearness (epsilon) categories (Table 1) */
203 static const CategoryBounds SkyClearCat[8] =
204 {
205 { 1.000, 1.065 }, /* Overcast */
206 { 1.065, 1.230 },
207 { 1.230, 1.500 },
208 { 1.500, 1.950 },
209 { 1.950, 2.800 },
210 { 2.800, 4.500 },
211 { 4.500, 6.200 },
212 { 6.200, 12.01 } /* Clear */
213 };
214
215 /* Luminous efficacy model coefficients */
216 typedef struct
217 {
218 double a;
219 double b;
220 double c;
221 double d;
222 } ModelCoeff;
223
224 /* Diffuse luminous efficacy model coefficients (Table 4, Eqn. 7) */
225 static const ModelCoeff DiffuseLumEff[8] =
226 {
227 { 97.24, -0.46, 12.00, -8.91 },
228 { 107.22, 1.15, 0.59, -3.95 },
229 { 104.97, 2.96, -5.53, -8.77 },
230 { 102.39, 5.59, -13.95, -13.90 },
231 { 100.71, 5.94, -22.75, -23.74 },
232 { 106.42, 3.83, -36.15, -28.83 },
233 { 141.88, 1.90, -53.24, -14.03 },
234 { 152.23, 0.35, -45.27, -7.98 }
235 };
236
237 /* Direct luminous efficacy model coefficients (Table 4, Eqn. 8) */
238 static const ModelCoeff DirectLumEff[8] =
239 {
240 { 57.20, -4.55, -2.98, 117.12 },
241 { 98.99, -3.46, -1.21, 12.38 },
242 { 109.83, -4.90, -1.71, -8.81 },
243 { 110.34, -5.84, -1.99, -4.56 },
244 { 106.36, -3.97, -1.75, -6.16 },
245 { 107.19, -1.25, -1.51, -26.73 },
246 { 105.75, 0.77, -1.26, -34.44 },
247 { 101.18, 1.58, -1.10, -8.29 }
248 };
249
250 #ifndef NSUNPATCH
251 #define NSUNPATCH 4 /* max. # patches to spread sun into */
252 #endif
253
254 #define SUN_ANG_DEG 0.533 /* sun full-angle in degrees */
255
256 int nsuns = NSUNPATCH; /* number of sun patches to use */
257 double fixed_sun_sa = -1; /* fixed solid angle per sun? */
258
259 int verbose = 0; /* progress reports to stderr? */
260
261 int outfmt = 'a'; /* output format */
262
263 int rhsubdiv = 1; /* Reinhart sky subdivisions */
264
265 COLOR skycolor = {.96, 1.004, 1.118}; /* sky coloration */
266 COLOR suncolor = {1., 1., 1.}; /* sun color */
267 COLOR grefl = {.2, .2, .2}; /* ground reflectance */
268
269 int nskypatch; /* number of Reinhart patches */
270 float *rh_palt; /* sky patch altitudes (radians) */
271 float *rh_pazi; /* sky patch azimuths (radians) */
272 float *rh_dom; /* sky patch solid angle (sr) */
273
274 #define vector(v,alt,azi) ( (v)[1] = cos(alt), \
275 (v)[0] = (v)[1]*sin(azi), \
276 (v)[1] *= cos(azi), \
277 (v)[2] = sin(alt) )
278
279 #define rh_vector(v,i) vector(v,rh_palt[i],rh_pazi[i])
280
281 #define rh_cos(i) tsin(rh_palt[i])
282
283 #define solar_minute(jd,hr) ((24*60)*((jd)-1)+(int)((hr)*60.+.5))
284
285 extern int rh_init(void);
286 extern float * resize_dmatrix(float *mtx_data, int nsteps, int npatch);
287 extern void OutputSun(int id, int goodsun, FILE *fp, FILE *mfp);
288 extern void AddDirect(float *parr);
289
290
291 static const char *
292 getfmtname(int fmt)
293 {
294 switch (fmt) {
295 case 'a':
296 return("ascii");
297 case 'f':
298 return("float");
299 case 'd':
300 return("double");
301 }
302 return("unknown");
303 }
304
305
306 int
307 main(int argc, char *argv[])
308 {
309 char buf[256];
310 int doheader = 1; /* output header? */
311 double rotation = 0; /* site rotation (degrees) */
312 double elevation; /* site elevation (meters) */
313 int leap_day = 0; /* add leap day? */
314 int sun_hours_only = 0; /* only output sun hours? */
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 'u': /* solar hours only */
410 sun_hours_only = 1;
411 break;
412 case 'r': /* rotate distribution */
413 if (argv[i][2] && argv[i][2] != 'z')
414 goto userr;
415 rotation = atof(argv[++i]);
416 break;
417 case '5': /* 5-phase calculation */
418 nsuns = 1;
419 fixed_sun_sa = PI/360.*atof(argv[++i]);
420 if (fixed_sun_sa <= 0) {
421 fprintf(stderr, "%s: missing solar disk size argument for '-5' option\n",
422 progname);
423 exit(1);
424 }
425 fixed_sun_sa *= fixed_sun_sa*PI;
426 break;
427 case 'A': /* compute average sky */
428 avgSky = 1;
429 break;
430 default:
431 goto userr;
432 }
433 if (i < argc-1)
434 goto userr;
435 if (i == argc-1 && freopen(argv[i], "r", stdin) == NULL) {
436 fprintf(stderr, "%s: cannot open '%s' for input\n",
437 progname, argv[i]);
438 exit(1);
439 }
440 if ((modsfp != NULL) & (sunsfp == NULL))
441 fprintf(stderr, "%s: warning -M output will be empty without -D\n",
442 progname);
443 if (verbose) {
444 if (i == argc-1)
445 fprintf(stderr, "%s: reading weather tape '%s'\n",
446 progname, argv[i]);
447 else
448 fprintf(stderr, "%s: reading weather tape from <stdin>\n",
449 progname);
450 }
451 /* read weather tape header */
452 if (scanf("place %[^\r\n] ", buf) != 1)
453 goto fmterr;
454 if (scanf("latitude %lf\n", &s_latitude) != 1)
455 goto fmterr;
456 if (scanf("longitude %lf\n", &s_longitude) != 1)
457 goto fmterr;
458 if (scanf("time_zone %lf\n", &s_meridian) != 1)
459 goto fmterr;
460 if (scanf("site_elevation %lf\n", &elevation) != 1)
461 goto fmterr;
462 if (scanf("weather_data_file_units %d\n", &input) != 1)
463 goto fmterr;
464 switch (input) { /* translate units */
465 case 1:
466 input = 1; /* radiometric quantities */
467 dir_is_horiz = 0; /* direct is perpendicular meas. */
468 break;
469 case 2:
470 input = 1; /* radiometric quantities */
471 dir_is_horiz = 1; /* solar measured horizontally */
472 break;
473 case 3:
474 input = 2; /* photometric quantities */
475 dir_is_horiz = 0; /* direct is perpendicular meas. */
476 break;
477 default:
478 goto fmterr;
479 }
480 rh_init(); /* initialize sky patches */
481 if (verbose) {
482 fprintf(stderr, "%s: location '%s'\n", progname, buf);
483 fprintf(stderr, "%s: (lat,long)=(%.1f,%.1f) degrees north, west\n",
484 progname, s_latitude, s_longitude);
485 if (avgSky >= 0)
486 fprintf(stderr, "%s: %d sky patches\n",
487 progname, nskypatch);
488 if (sunsfp)
489 fprintf(stderr, "%s: outputting suns to file\n",
490 progname);
491 if (rotation != 0)
492 fprintf(stderr, "%s: rotating output %.0f degrees\n",
493 progname, rotation);
494 }
495 /* convert quantities to radians */
496 s_latitude = DegToRad(s_latitude);
497 s_longitude = DegToRad(s_longitude);
498 s_meridian = DegToRad(s_meridian);
499 /* initial allocation */
500 mtx_data = resize_dmatrix(mtx_data, tstorage=2, nskypatch);
501 /* process each time step in tape */
502 while (scanf("%d %d %lf %lf %lf\n", &mo, &da, &hr, &dir, &dif) == 5) {
503 double sda, sta;
504 int sun_in_sky;
505 /* compute solar position */
506 if ((mo == 2) & (da == 29)) {
507 julian_date = 60;
508 leap_day = 1;
509 } else
510 julian_date = jdate(mo, da) + leap_day;
511 sda = sdec(julian_date);
512 sta = stadj(julian_date);
513 altitude = salt(sda, hr+sta);
514 sun_in_sky = (altitude > -DegToRad(SUN_ANG_DEG/2.));
515 if (sun_hours_only && !sun_in_sky)
516 continue; /* skipping nighttime points */
517 azimuth = sazi(sda, hr+sta) + PI - DegToRad(rotation);
518
519 mtx_offset = 3*nskypatch*nstored;
520 nstored += !avgSky | !nstored;
521 /* make space for next row */
522 if (nstored > tstorage) {
523 tstorage += (tstorage>>1) + nstored + 7;
524 mtx_data = resize_dmatrix(mtx_data, tstorage, nskypatch);
525 }
526 ntsteps++; /* keep count of time steps */
527
528 if (dir+dif <= 1e-4) { /* effectively nighttime? */
529 if (!avgSky | !mtx_offset)
530 memset(mtx_data+mtx_offset, 0,
531 sizeof(float)*3*nskypatch);
532 /* output black sun? */
533 if (sunsfp && sun_in_sky)
534 OutputSun(solar_minute(julian_date,hr), 0,
535 sunsfp, modsfp);
536 continue;
537 }
538 /* convert measured values */
539 if (dir_is_horiz && altitude > 0.)
540 dir /= sin(altitude);
541 if (input == 1) {
542 dir_irrad = dir;
543 diff_irrad = dif;
544 } else /* input == 2 */ {
545 dir_illum = dir;
546 diff_illum = dif;
547 }
548 /* compute sky patch values */
549 ComputeSky(mtx_data+mtx_offset);
550 /* output sun if requested */
551 if (sunsfp && sun_in_sky)
552 OutputSun(solar_minute(julian_date,hr), 1,
553 sunsfp, modsfp);
554
555 if (avgSky < 0) /* no matrix? */
556 continue;
557
558 AddDirect(mtx_data+mtx_offset);
559 /* update cumulative sky? */
560 for (i = 3*nskypatch*(avgSky&(ntsteps>1)); i--; )
561 mtx_data[i] += mtx_data[mtx_offset+i];
562 /* monthly reporting */
563 if (verbose && mo != last_monthly)
564 fprintf(stderr, "%s: stepping through month %d...\n",
565 progname, last_monthly=mo);
566 /* note whether leap-day was given */
567 }
568 if (!ntsteps) {
569 fprintf(stderr, "%s: no valid time steps on input\n", progname);
570 exit(1);
571 }
572 /* check for junk at end */
573 while ((i = fgetc(stdin)) != EOF)
574 if (!isspace(i)) {
575 fprintf(stderr, "%s: warning - unexpected data past EOT: ",
576 progname);
577 buf[0] = i; buf[1] = '\0';
578 fgets(buf+1, sizeof(buf)-1, stdin);
579 fputs(buf, stderr); fputc('\n', stderr);
580 break;
581 }
582
583 if (avgSky < 0) /* no matrix output? */
584 goto alldone;
585
586 dif = 1./(double)ntsteps; /* average sky? */
587 for (i = 3*nskypatch*(avgSky&(ntsteps>1)); i--; )
588 mtx_data[i] *= dif;
589 /* write out matrix */
590 if (outfmt != 'a')
591 SET_FILE_BINARY(stdout);
592 #ifdef getc_unlocked
593 flockfile(stdout);
594 #endif
595 if (verbose)
596 fprintf(stderr, "%s: writing %smatrix with %d time steps...\n",
597 progname, outfmt=='a' ? "" : "binary ", nstored);
598 if (doheader) {
599 newheader("RADIANCE", stdout);
600 printargs(argc, argv, stdout);
601 printf("LATLONG= %.8f %.8f\n", RadToDeg(s_latitude),
602 -RadToDeg(s_longitude));
603 printf("NROWS=%d\n", nskypatch);
604 printf("NCOLS=%d\n", nstored);
605 printf("NCOMP=3\n");
606 if ((outfmt == 'f') | (outfmt == 'd'))
607 fputendian(stdout);
608 fputformat((char *)getfmtname(outfmt), stdout);
609 putchar('\n');
610 }
611 /* patches are rows (outer sort) */
612 for (i = 0; i < nskypatch; i++) {
613 mtx_offset = 3*i;
614 switch (outfmt) {
615 case 'a':
616 for (j = 0; j < nstored; j++) {
617 printf("%.3g %.3g %.3g\n", mtx_data[mtx_offset],
618 mtx_data[mtx_offset+1],
619 mtx_data[mtx_offset+2]);
620 mtx_offset += 3*nskypatch;
621 }
622 if (nstored > 1)
623 fputc('\n', stdout);
624 break;
625 case 'f':
626 for (j = 0; j < nstored; j++) {
627 putbinary(mtx_data+mtx_offset, sizeof(float), 3,
628 stdout);
629 mtx_offset += 3*nskypatch;
630 }
631 break;
632 case 'd':
633 for (j = 0; j < nstored; j++) {
634 double ment[3];
635 ment[0] = mtx_data[mtx_offset];
636 ment[1] = mtx_data[mtx_offset+1];
637 ment[2] = mtx_data[mtx_offset+2];
638 putbinary(ment, sizeof(double), 3, stdout);
639 mtx_offset += 3*nskypatch;
640 }
641 break;
642 }
643 if (ferror(stdout))
644 goto writerr;
645 }
646 alldone:
647 if (fflush(NULL) == EOF)
648 goto writerr;
649 if (verbose)
650 fprintf(stderr, "%s: done.\n", progname);
651 exit(0);
652 userr:
653 fprintf(stderr, "Usage: %s [-v][-h][-A][-d|-s|-n][-u][-D file [-M modfile]][-r deg][-m N][-g r g b][-c r g b][-o{f|d}][-O{0|1}] [tape.wea]\n",
654 progname);
655 exit(1);
656 fmterr:
657 fprintf(stderr, "%s: weather tape format error in header\n", progname);
658 exit(1);
659 writerr:
660 fprintf(stderr, "%s: write error on output\n", progname);
661 exit(1);
662 }
663
664 /* Return maximum of two doubles */
665 double dmax( double a, double b )
666 { return (a > b) ? a : b; }
667
668 /* Compute sky patch radiance values (modified by GW) */
669 void
670 ComputeSky(float *parr)
671 {
672 int index; /* Category index */
673 double norm_diff_illum; /* Normalized diffuse illuimnance */
674 int i;
675
676 /* Calculate atmospheric precipitable water content */
677 apwc = CalcPrecipWater(dew_point);
678
679 /* Calculate sun zenith angle (don't let it dip below horizon) */
680 /* Also limit minimum angle to keep circumsolar off zenith */
681 if (altitude <= 0.0)
682 sun_zenith = DegToRad(90.0);
683 else if (altitude >= DegToRad(87.0))
684 sun_zenith = DegToRad(3.0);
685 else
686 sun_zenith = DegToRad(90.0) - altitude;
687
688 /* Compute the inputs for the calculation of the sky distribution */
689
690 if (input == 0) /* XXX never used */
691 {
692 /* Calculate irradiance */
693 diff_irrad = CalcDiffuseIrradiance();
694 dir_irrad = CalcDirectIrradiance();
695
696 /* Calculate illuminance */
697 index = GetCategoryIndex();
698 diff_illum = diff_irrad * CalcDiffuseIllumRatio(index);
699 dir_illum = dir_irrad * CalcDirectIllumRatio(index);
700 }
701 else if (input == 1)
702 {
703 sky_brightness = CalcSkyBrightness();
704 sky_clearness = CalcSkyClearness();
705
706 /* Limit sky clearness */
707 if (sky_clearness > 11.9)
708 sky_clearness = 11.9;
709
710 /* Limit sky brightness */
711 if (sky_brightness < 0.01)
712 sky_brightness = 0.01;
713
714 /* Calculate illuminance */
715 index = GetCategoryIndex();
716 diff_illum = diff_irrad * CalcDiffuseIllumRatio(index);
717 dir_illum = dir_irrad * CalcDirectIllumRatio(index);
718 }
719 else if (input == 2)
720 {
721 /* Calculate sky brightness and clearness from illuminance values */
722 index = CalcSkyParamFromIllum();
723 }
724
725 if (output == 1) { /* hack for solar radiance */
726 diff_illum = diff_irrad * WHTEFFICACY;
727 dir_illum = dir_irrad * WHTEFFICACY;
728 }
729 /* Compute ground radiance (include solar contribution if any) */
730 parr[0] = diff_illum;
731 if (altitude > 0)
732 parr[0] += dir_illum * sin(altitude);
733 parr[2] = parr[1] = parr[0] *= (1./PI/WHTEFFICACY);
734 multcolor(parr, grefl);
735
736 if (bright(skycolor) <= 1e-4) { /* 0 sky component? */
737 memset(parr+3, 0, sizeof(float)*3*(nskypatch-1));
738 return;
739 }
740 /* Calculate Perez sky model parameters */
741 CalcPerezParam(sun_zenith, sky_clearness, sky_brightness, index);
742
743 /* Calculate sky patch luminance values */
744 CalcSkyPatchLumin(parr);
745
746 /* Calculate relative horizontal illuminance */
747 norm_diff_illum = CalcRelHorzIllum(parr);
748
749 /* Check for zero sky -- make uniform in that case */
750 if (norm_diff_illum <= FTINY) {
751 for (i = 1; i < nskypatch; i++)
752 setcolor(parr+3*i, 1., 1., 1.);
753 norm_diff_illum = PI;
754 }
755 /* Normalization coefficient */
756 norm_diff_illum = diff_illum / norm_diff_illum;
757
758 /* Apply to sky patches to get absolute radiance values */
759 for (i = 1; i < nskypatch; i++) {
760 scalecolor(parr+3*i, norm_diff_illum*(1./WHTEFFICACY));
761 multcolor(parr+3*i, skycolor);
762 }
763 }
764
765 /* Add in solar direct to nearest sky patches (GW) */
766 void
767 AddDirect(float *parr)
768 {
769 FVECT svec;
770 double near_dprod[NSUNPATCH];
771 int near_patch[NSUNPATCH];
772 double wta[NSUNPATCH], wtot;
773 int i, j, p;
774
775 if (dir_illum <= 1e-4 || bright(suncolor) <= 1e-4)
776 return;
777 /* identify nsuns closest patches */
778 if (nsuns > NSUNPATCH)
779 nsuns = NSUNPATCH;
780 else if (nsuns <= 0)
781 nsuns = 1;
782 for (i = nsuns; i--; )
783 near_dprod[i] = -1.;
784 vector(svec, altitude, azimuth);
785 for (p = 1; p < nskypatch; p++) {
786 FVECT pvec;
787 double dprod;
788 rh_vector(pvec, p);
789 dprod = DOT(pvec, svec);
790 for (i = 0; i < nsuns; i++)
791 if (dprod > near_dprod[i]) {
792 for (j = nsuns; --j > i; ) {
793 near_dprod[j] = near_dprod[j-1];
794 near_patch[j] = near_patch[j-1];
795 }
796 near_dprod[i] = dprod;
797 near_patch[i] = p;
798 break;
799 }
800 }
801 wtot = 0; /* weight by proximity */
802 for (i = nsuns; i--; )
803 wtot += wta[i] = 1./(1.002 - near_dprod[i]);
804 /* add to nearest patch radiances */
805 for (i = nsuns; i--; ) {
806 float *pdest = parr + 3*near_patch[i];
807 float val_add = wta[i] * dir_illum / (WHTEFFICACY * wtot);
808
809 val_add /= (fixed_sun_sa > 0) ? fixed_sun_sa
810 : rh_dom[near_patch[i]] ;
811 *pdest++ += val_add*suncolor[0];
812 *pdest++ += val_add*suncolor[1];
813 *pdest++ += val_add*suncolor[2];
814 }
815 }
816
817 /* Output a sun to indicated file if appropriate for this time step */
818 void
819 OutputSun(int id, int goodsun, FILE *fp, FILE *mfp)
820 {
821 double srad;
822 FVECT sv;
823
824 srad = DegToRad(SUN_ANG_DEG/2.);
825 srad = goodsun ? dir_illum/(WHTEFFICACY * PI*srad*srad) : 0;
826 vector(sv, altitude, azimuth);
827 fprintf(fp, "\nvoid light solar%d\n0\n0\n", id);
828 fprintf(fp, "3 %.3e %.3e %.3e\n", srad*suncolor[0],
829 srad*suncolor[1], srad*suncolor[2]);
830 fprintf(fp, "\nsolar%d source sun%d\n0\n0\n", id, id);
831 fprintf(fp, "4 %.6f %.6f %.6f %.4f\n", sv[0], sv[1], sv[2], SUN_ANG_DEG);
832
833 if (mfp != NULL) /* saving modifier IDs? */
834 fprintf(mfp, "solar%d\n", id);
835 }
836
837 /* Initialize Reinhart sky patch positions (GW) */
838 int
839 rh_init(void)
840 {
841 #define NROW 7
842 static const int tnaz[NROW] = {30, 30, 24, 24, 18, 12, 6};
843 const double alpha = (PI/2.)/(NROW*rhsubdiv + .5);
844 int p, i, j;
845 /* allocate patch angle arrays */
846 nskypatch = 0;
847 for (p = 0; p < NROW; p++)
848 nskypatch += tnaz[p];
849 nskypatch *= rhsubdiv*rhsubdiv;
850 nskypatch += 2;
851 rh_palt = (float *)malloc(sizeof(float)*nskypatch);
852 rh_pazi = (float *)malloc(sizeof(float)*nskypatch);
853 rh_dom = (float *)malloc(sizeof(float)*nskypatch);
854 if ((rh_palt == NULL) | (rh_pazi == NULL) | (rh_dom == NULL)) {
855 fprintf(stderr, "%s: out of memory in rh_init()\n", progname);
856 exit(1);
857 }
858 rh_palt[0] = -PI/2.; /* ground & zenith patches */
859 rh_pazi[0] = 0.;
860 rh_dom[0] = 2.*PI;
861 rh_palt[nskypatch-1] = PI/2.;
862 rh_pazi[nskypatch-1] = 0.;
863 rh_dom[nskypatch-1] = 2.*PI*(1. - cos(alpha*.5));
864 p = 1; /* "normal" patches */
865 for (i = 0; i < NROW*rhsubdiv; i++) {
866 const float ralt = alpha*(i + .5);
867 const int ninrow = tnaz[i/rhsubdiv]*rhsubdiv;
868 const float dom = 2.*PI*(sin(alpha*(i+1)) - sin(alpha*i)) /
869 (double)ninrow;
870 for (j = 0; j < ninrow; j++) {
871 rh_palt[p] = ralt;
872 rh_pazi[p] = 2.*PI * j / (double)ninrow;
873 rh_dom[p++] = dom;
874 }
875 }
876 return nskypatch;
877 #undef NROW
878 }
879
880 /* Resize daylight matrix (GW) */
881 float *
882 resize_dmatrix(float *mtx_data, int nsteps, int npatch)
883 {
884 if (mtx_data == NULL)
885 mtx_data = (float *)malloc(sizeof(float)*3*nsteps*npatch);
886 else
887 mtx_data = (float *)realloc(mtx_data,
888 sizeof(float)*3*nsteps*npatch);
889 if (mtx_data == NULL) {
890 fprintf(stderr, "%s: out of memory in resize_dmatrix(%d,%d)\n",
891 progname, nsteps, npatch);
892 exit(1);
893 }
894 return(mtx_data);
895 }
896
897 /* Determine category index */
898 int GetCategoryIndex()
899 {
900 int index; /* Loop index */
901
902 for (index = 0; index < 8; index++)
903 if ((sky_clearness >= SkyClearCat[index].lower) &&
904 (sky_clearness < SkyClearCat[index].upper))
905 break;
906
907 return index;
908 }
909
910 /* Calculate diffuse illuminance to diffuse irradiance ratio */
911
912 /* Reference: Perez, R., P. Ineichen, R. Seals, J. Michalsky, and R. */
913 /* Stewart. 1990. ìModeling Daylight Availability and */
914 /* Irradiance Components from Direct and Global */
915 /* Irradiance,î Solar Energy 44(5):271-289, Eqn. 7. */
916
917 double CalcDiffuseIllumRatio( int index )
918 {
919 ModelCoeff const *pnle; /* Category coefficient pointer */
920
921 /* Get category coefficient pointer */
922 pnle = &(DiffuseLumEff[index]);
923
924 return pnle->a + pnle->b * apwc + pnle->c * cos(sun_zenith) +
925 pnle->d * log(sky_brightness);
926 }
927
928 /* Calculate direct illuminance to direct irradiance ratio */
929
930 /* Reference: Perez, R., P. Ineichen, R. Seals, J. Michalsky, and R. */
931 /* Stewart. 1990. ìModeling Daylight Availability and */
932 /* Irradiance Components from Direct and Global */
933 /* Irradiance,î Solar Energy 44(5):271-289, Eqn. 8. */
934
935 double CalcDirectIllumRatio( int index )
936 {
937 ModelCoeff const *pnle; /* Category coefficient pointer */
938
939 /* Get category coefficient pointer */
940 pnle = &(DirectLumEff[index]);
941
942 /* Calculate direct illuminance from direct irradiance */
943
944 return dmax((pnle->a + pnle->b * apwc + pnle->c * exp(5.73 *
945 sun_zenith - 5.0) + pnle->d * sky_brightness),
946 0.0);
947 }
948
949 /* Calculate sky brightness */
950
951 /* Reference: Perez, R., P. Ineichen, R. Seals, J. Michalsky, and R. */
952 /* Stewart. 1990. ìModeling Daylight Availability and */
953 /* Irradiance Components from Direct and Global */
954 /* Irradiance,î Solar Energy 44(5):271-289, Eqn. 2. */
955
956 double CalcSkyBrightness()
957 {
958 return diff_irrad * CalcAirMass() / (DC_SolarConstantE *
959 CalcEccentricity());
960 }
961
962 /* Calculate sky clearness */
963
964 /* Reference: Perez, R., P. Ineichen, R. Seals, J. Michalsky, and R. */
965 /* Stewart. 1990. ìModeling Daylight Availability and */
966 /* Irradiance Components from Direct and Global */
967 /* Irradiance,î Solar Energy 44(5):271-289, Eqn. 1. */
968
969 double CalcSkyClearness()
970 {
971 double sz_cubed; /* Sun zenith angle cubed */
972
973 /* Calculate sun zenith angle cubed */
974 sz_cubed = sun_zenith*sun_zenith*sun_zenith;
975
976 return ((diff_irrad + dir_irrad) / diff_irrad + 1.041 *
977 sz_cubed) / (1.0 + 1.041 * sz_cubed);
978 }
979
980 /* Calculate diffuse horizontal irradiance from Perez sky brightness */
981
982 /* Reference: Perez, R., P. Ineichen, R. Seals, J. Michalsky, and R. */
983 /* Stewart. 1990. ìModeling Daylight Availability and */
984 /* Irradiance Components from Direct and Global */
985 /* Irradiance,î Solar Energy 44(5):271-289, Eqn. 2 */
986 /* (inverse). */
987
988 double CalcDiffuseIrradiance()
989 {
990 return sky_brightness * DC_SolarConstantE * CalcEccentricity() /
991 CalcAirMass();
992 }
993
994 /* Calculate direct normal irradiance from Perez sky clearness */
995
996 /* Reference: Perez, R., P. Ineichen, R. Seals, J. Michalsky, and R. */
997 /* Stewart. 1990. ìModeling Daylight Availability and */
998 /* Irradiance Components from Direct and Global */
999 /* Irradiance,î Solar Energy 44(5):271-289, Eqn. 1 */
1000 /* (inverse). */
1001
1002 double CalcDirectIrradiance()
1003 {
1004 return CalcDiffuseIrradiance() * ((sky_clearness - 1.0) * (1 + 1.041
1005 * sun_zenith*sun_zenith*sun_zenith));
1006 }
1007
1008 /* Calculate sky brightness and clearness from illuminance values */
1009 int CalcSkyParamFromIllum()
1010 {
1011 double test1 = 0.1;
1012 double test2 = 0.1;
1013 int counter = 0;
1014 int index = 0; /* Category index */
1015
1016 /* Convert illuminance to irradiance */
1017 diff_irrad = diff_illum * DC_SolarConstantE /
1018 (DC_SolarConstantL * 1000.0);
1019 dir_irrad = dir_illum * DC_SolarConstantE /
1020 (DC_SolarConstantL * 1000.0);
1021
1022 /* Calculate sky brightness and clearness */
1023 sky_brightness = CalcSkyBrightness();
1024 sky_clearness = CalcSkyClearness();
1025
1026 /* Limit sky clearness */
1027 if (sky_clearness > 12.0)
1028 sky_clearness = 12.0;
1029
1030 /* Limit sky brightness */
1031 if (sky_brightness < 0.01)
1032 sky_brightness = 0.01;
1033
1034 while (((fabs(diff_irrad - test1) > 10.0) ||
1035 (fabs(dir_irrad - test2) > 10.0)) && !(counter == 5))
1036 {
1037 test1 = diff_irrad;
1038 test2 = dir_irrad;
1039 counter++;
1040
1041 /* Convert illuminance to irradiance */
1042 index = GetCategoryIndex();
1043 diff_irrad = diff_illum / CalcDiffuseIllumRatio(index);
1044 dir_irrad = CalcDirectIllumRatio(index);
1045 if (dir_irrad > 0.1)
1046 dir_irrad = dir_illum / dir_irrad;
1047
1048 /* Calculate sky brightness and clearness */
1049 sky_brightness = CalcSkyBrightness();
1050 sky_clearness = CalcSkyClearness();
1051
1052 /* Limit sky clearness */
1053 if (sky_clearness > 12.0)
1054 sky_clearness = 12.0;
1055
1056 /* Limit sky brightness */
1057 if (sky_brightness < 0.01)
1058 sky_brightness = 0.01;
1059 }
1060
1061 return GetCategoryIndex();
1062 }
1063
1064 /* Calculate relative luminance */
1065
1066 /* Reference: Perez, R., R. Seals, and J. Michalsky. 1993. */
1067 /* ìAll-Weather Model for Sky Luminance Distribution - */
1068 /* Preliminary Configuration and Validation,î Solar Energy */
1069 /* 50(3):235-245, Eqn. 1. */
1070
1071 double CalcRelLuminance( double gamma, double zeta )
1072 {
1073 return (1.0 + perez_param[0] * exp(perez_param[1] / cos(zeta))) *
1074 (1.0 + perez_param[2] * exp(perez_param[3] * gamma) +
1075 perez_param[4] * cos(gamma) * cos(gamma));
1076 }
1077
1078 /* Calculate Perez sky model parameters */
1079
1080 /* Reference: Perez, R., R. Seals, and J. Michalsky. 1993. */
1081 /* ìAll-Weather Model for Sky Luminance Distribution - */
1082 /* Preliminary Configuration and Validation,î Solar Energy */
1083 /* 50(3):235-245, Eqns. 6 - 8. */
1084
1085 void CalcPerezParam( double sz, double epsilon, double delta,
1086 int index )
1087 {
1088 double x[5][4]; /* Coefficents a, b, c, d, e */
1089 int i, j; /* Loop indices */
1090
1091 /* Limit sky brightness */
1092 if (epsilon > 1.065 && epsilon < 2.8)
1093 {
1094 if (delta < 0.2)
1095 delta = 0.2;
1096 }
1097
1098 /* Get Perez coefficients */
1099 for (i = 0; i < 5; i++)
1100 for (j = 0; j < 4; j++)
1101 x[i][j] = PerezCoeff[index][4 * i + j];
1102
1103 if (index != 0)
1104 {
1105 /* Calculate parameter a, b, c, d and e (Eqn. 6) */
1106 for (i = 0; i < 5; i++)
1107 perez_param[i] = x[i][0] + x[i][1] * sz + delta * (x[i][2] +
1108 x[i][3] * sz);
1109 }
1110 else
1111 {
1112 /* Parameters a, b and e (Eqn. 6) */
1113 perez_param[0] = x[0][0] + x[0][1] * sz + delta * (x[0][2] +
1114 x[0][3] * sz);
1115 perez_param[1] = x[1][0] + x[1][1] * sz + delta * (x[1][2] +
1116 x[1][3] * sz);
1117 perez_param[4] = x[4][0] + x[4][1] * sz + delta * (x[4][2] +
1118 x[4][3] * sz);
1119
1120 /* Parameter c (Eqn. 7) */
1121 perez_param[2] = exp(pow(delta * (x[2][0] + x[2][1] * sz),
1122 x[2][2])) - x[2][3];
1123
1124 /* Parameter d (Eqn. 8) */
1125 perez_param[3] = -exp(delta * (x[3][0] + x[3][1] * sz)) +
1126 x[3][2] + delta * x[3][3];
1127 }
1128 }
1129
1130 /* Calculate relative horizontal illuminance (modified by GW) */
1131
1132 /* Reference: Perez, R., R. Seals, and J. Michalsky. 1993. */
1133 /* ìAll-Weather Model for Sky Luminance Distribution - */
1134 /* Preliminary Configuration and Validation,î Solar Energy */
1135 /* 50(3):235-245, Eqn. 3. */
1136
1137 double CalcRelHorzIllum( float *parr )
1138 {
1139 int i;
1140 double rh_illum = 0.0; /* Relative horizontal illuminance */
1141
1142 for (i = 1; i < nskypatch; i++)
1143 rh_illum += parr[3*i+1] * rh_cos(i) * rh_dom[i];
1144
1145 return rh_illum;
1146 }
1147
1148 /* Calculate earth orbit eccentricity correction factor */
1149
1150 /* Reference: Sen, Z. 2008. Solar Energy Fundamental and Modeling */
1151 /* Techniques. Springer, p. 72. */
1152
1153 double CalcEccentricity()
1154 {
1155 double day_angle; /* Day angle (radians) */
1156 double E0; /* Eccentricity */
1157
1158 /* Calculate day angle */
1159 day_angle = (julian_date - 1.0) * (2.0 * PI / 365.0);
1160
1161 /* Calculate eccentricity */
1162 E0 = 1.00011 + 0.034221 * cos(day_angle) + 0.00128 * sin(day_angle)
1163 + 0.000719 * cos(2.0 * day_angle) + 0.000077 * sin(2.0 *
1164 day_angle);
1165
1166 return E0;
1167 }
1168
1169 /* Calculate atmospheric precipitable water content */
1170
1171 /* Reference: Perez, R., P. Ineichen, R. Seals, J. Michalsky, and R. */
1172 /* Stewart. 1990. ìModeling Daylight Availability and */
1173 /* Irradiance Components from Direct and Global */
1174 /* Irradiance,î Solar Energy 44(5):271-289, Eqn. 3. */
1175
1176 /* Note: The default surface dew point temperature is 11 deg. C */
1177 /* (52 deg. F). Typical values are: */
1178
1179 /* Celsius Fahrenheit Human Perception */
1180 /* > 24 > 75 Extremely uncomfortable */
1181 /* 21 - 24 70 - 74 Very humid */
1182 /* 18 - 21 65 - 69 Somewhat uncomfortable */
1183 /* 16 - 18 60 - 64 OK for most people */
1184 /* 13 - 16 55 - 59 Comfortable */
1185 /* 10 - 12 50 - 54 Very comfortable */
1186 /* < 10 < 49 A bit dry for some */
1187
1188 double CalcPrecipWater( double dpt )
1189 { return exp(0.07 * dpt - 0.075); }
1190
1191 /* Calculate relative air mass */
1192
1193 /* Reference: Kasten, F. 1966. "A New Table and Approximation Formula */
1194 /* for the Relative Optical Air Mass," Arch. Meteorol. */
1195 /* Geophys. Bioklimataol. Ser. B14, pp. 206-233. */
1196
1197 /* Note: More sophisticated relative air mass models are */
1198 /* available, but they differ significantly only for */
1199 /* sun zenith angles greater than 80 degrees. */
1200
1201 double CalcAirMass()
1202 {
1203 return (1.0 / (cos(sun_zenith) + 0.15 * pow(93.885 -
1204 RadToDeg(sun_zenith), -1.253)));
1205 }
1206
1207 /* Calculate Perez All-Weather sky patch luminances (modified by GW) */
1208
1209 /* NOTE: The sky patches centers are determined in accordance with the */
1210 /* BRE-IDMP sky luminance measurement procedures. (See for example */
1211 /* Mardaljevic, J. 2001. "The BRE-IDMP Dataset: A New Benchmark */
1212 /* for the Validation of Illuminance Prediction Techniques," */
1213 /* Lighting Research & Technology 33(2):117-136.) */
1214
1215 void CalcSkyPatchLumin( float *parr )
1216 {
1217 int i;
1218 double aas; /* Sun-sky point azimuthal angle */
1219 double sspa; /* Sun-sky point angle */
1220 double zsa; /* Zenithal sun angle */
1221
1222 for (i = 1; i < nskypatch; i++)
1223 {
1224 /* Calculate sun-sky point azimuthal angle */
1225 aas = fabs(rh_pazi[i] - azimuth);
1226
1227 /* Calculate zenithal sun angle */
1228 zsa = PI * 0.5 - rh_palt[i];
1229
1230 /* Calculate sun-sky point angle (Equation 8-20) */
1231 sspa = acos(cos(sun_zenith) * cos(zsa) + sin(sun_zenith) *
1232 sin(zsa) * cos(aas));
1233
1234 /* Calculate patch luminance */
1235 parr[3*i] = CalcRelLuminance(sspa, zsa);
1236 if (parr[3*i] < 0) parr[3*i] = 0;
1237 parr[3*i+2] = parr[3*i+1] = parr[3*i];
1238 }
1239 }