ViewVC Help
View File | Revision Log | Show Annotations | Download File | Root Listing
root/radiance/ray/src/gen/gendaymtx.c
Revision: 2.44
Committed: Sat Jun 7 05:09:45 2025 UTC (2 weeks, 4 days ago) by greg
Content type: text/plain
Branch: MAIN
CVS Tags: HEAD
Changes since 2.43: +1 -2 lines
Log Message:
refactor: Put some declarations into "paths.h" and included in "platform.h"

File Contents

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