ViewVC Help
View File | Revision Log | Show Annotations | Download File | Root Listing
root/radiance/ray/src/gen/gendaymtx.c
Revision: 2.42
Committed: Thu Feb 27 19:00:00 2025 UTC (2 months ago) by greg
Content type: text/plain
Branch: MAIN
CVS Tags: HEAD
Changes since 2.41: +14 -12 lines
Log Message:
fix(gendaymtx): Issue with interpretation of minutes in EPW input

File Contents

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