ViewVC Help
View File | Revision Log | Show Annotations | Download File | Root Listing
root/radiance/ray/src/gen/gensdaymtx.c
Revision: 1.4
Committed: Thu Aug 8 02:00:20 2024 UTC (8 months, 3 weeks ago) by greg
Content type: text/plain
Branch: MAIN
Changes since 1.3: +23 -13 lines
Log Message:
fix(gensdaymtx): Fixed sky value accumulation in first column/timestep

File Contents

# Content
1 #ifndef lint
2 static const char RCSid[] = "$Id$";
3 #endif
4 #include "atmos.h"
5 #include "copyright.h"
6 #include "data.h"
7 #include "platform.h"
8 #include "rtio.h"
9 #include <ctype.h>
10 #include <stdlib.h>
11 #ifdef _WIN32
12 #include <windows.h>
13 #else
14 #include <errno.h>
15 #include <sys/stat.h>
16 #include <sys/types.h>
17 #endif
18
19 char *progname;
20
21 double altitude; /* Solar altitude (radians) */
22 double azimuth; /* Solar azimuth (radians) */
23 int julian_date; /* Julian date */
24 double sun_zenith; /* Sun zenith angle (radians) */
25 int input = 0; /* Input type */
26 int output = 0; /* Output type */
27 FVECT sundir;
28
29 const double ARCTIC_LAT = 67.;
30 const double TROPIC_LAT = 23.;
31 const int SUMMER_START = 4;
32 const int SUMMER_END = 9;
33 const double GNORM = 0.777778;
34
35 const double D65EFF = 203.; /* standard illuminant D65 */
36
37 /* Mean normalized relative daylight spectra where CCT = 6415K for overcast */
38 const double D6415[NSSAMP] = {0.63231, 1.06171, 1.00779, 1.36423, 1.34133,
39 1.27258, 1.26276, 1.26352, 1.22201, 1.13246,
40 1.0434, 1.05547, 0.98212, 0.94445, 0.9722,
41 0.82387, 0.87853, 0.82559, 0.75111, 0.78925};
42 /* Degrees into radians */
43 #define DegToRad(deg) ((deg) * (PI / 180.))
44
45 /* Radiuans into degrees */
46 #define RadToDeg(rad) ((rad) * (180. / PI))
47
48 #ifndef NSUNPATCH
49 #define NSUNPATCH 4 /* max. # patches to spread sun into */
50 #endif
51
52 #define SUN_ANG_DEG 0.533 /* sun full-angle in degrees */
53
54 int nsuns = NSUNPATCH; /* number of sun patches to use */
55 double fixed_sun_sa = -1; /* fixed solid angle per sun? */
56
57 int verbose = 0; /* progress reports to stderr? */
58
59 int outfmt = 'a'; /* output format */
60
61 int rhsubdiv = 1; /* Reinhart sky subdivisions */
62
63 COLOR skycolor = {.96, 1.004, 1.118}; /* sky coloration */
64 COLOR suncolor = {1., 1., 1.}; /* sun color */
65 double grefl = .2; /* ground reflectance */
66
67 int nskypatch; /* number of Reinhart patches */
68 float *rh_palt; /* sky patch altitudes (radians) */
69 float *rh_pazi; /* sky patch azimuths (radians) */
70 float *rh_dom; /* sky patch solid angle (sr) */
71
72 double sun_ct;
73
74 #define vector(v, alt, azi) \
75 ((v)[1] = cos(alt), (v)[0] = (v)[1] * sin(azi), (v)[1] *= cos(azi), \
76 (v)[2] = sin(alt))
77
78 #define rh_vector(v, i) vector(v, rh_palt[i], rh_pazi[i])
79
80 #define rh_cos(i) tsin(rh_palt[i])
81
82 #define solar_minute(jd, hr) ((24 * 60) * ((jd) - 1) + (int)((hr) * 60. + .5))
83
84 inline void vectorize(double altitude, double azimuth, FVECT v) {
85 v[1] = cos(altitude);
86 v[0] = (v)[1] * sin(azimuth);
87 v[1] *= cos(azimuth);
88 v[2] = sin(altitude);
89 }
90
91 static int make_directory(const char *path) {
92 #ifdef _WIN32
93 if (CreateDirectory(path, NULL) || GetLastError() == ERROR_ALREADY_EXISTS) {
94 return 1;
95 }
96 return 0;
97 #else
98 if (mkdir(path, 0777) == 0 || errno == EEXIST) {
99 return 1;
100 }
101 return 0;
102 #endif
103 }
104
105 static const char *getfmtname(int fmt) {
106 switch (fmt) {
107 case 'a':
108 return ("ascii");
109 case 'f':
110 return ("float");
111 case 'd':
112 return ("double");
113 }
114 return ("unknown");
115 }
116
117 static inline double wmean2(const double a, const double b, const double x) {
118 return a * (1 - x) + b * x;
119 }
120
121 static inline double wmean(const double a, const double x, const double b,
122 const double y) {
123 return (a * x + b * y) / (a + b);
124 }
125
126 static double get_zenith_brightness(const double sundir[3]) {
127 double zenithbr;
128 if (sundir[2] < 0) {
129 zenithbr = 0;
130 } else {
131 zenithbr = (8.6 * sundir[2] + .123) * 1000.0 / D65EFF;
132 }
133 return zenithbr;
134 }
135
136 /* from gensky.c */
137 static double get_overcast_brightness(const double dz, const double zenithbr) {
138 double groundbr = zenithbr * GNORM;
139 return wmean(pow(dz + 1.01, 10), zenithbr * (1 + 2 * dz) / 3,
140 pow(dz + 1.01, -10), groundbr);
141 }
142
143 int rh_init(void) {
144 #define NROW 7
145 static const int tnaz[NROW] = {30, 30, 24, 24, 18, 12, 6};
146 const double alpha = (PI / 2.) / (NROW * rhsubdiv + .5);
147 int p, i, j;
148 /* allocate patch angle arrays */
149 nskypatch = 0;
150 for (p = 0; p < NROW; p++)
151 nskypatch += tnaz[p];
152 nskypatch *= rhsubdiv * rhsubdiv;
153 nskypatch += 2;
154 rh_palt = (float *)malloc(sizeof(float) * nskypatch);
155 rh_pazi = (float *)malloc(sizeof(float) * nskypatch);
156 rh_dom = (float *)malloc(sizeof(float) * nskypatch);
157 if ((rh_palt == NULL) | (rh_pazi == NULL) | (rh_dom == NULL)) {
158 fprintf(stderr, "%s: out of memory in rh_init()\n", progname);
159 exit(1);
160 }
161 rh_palt[0] = -PI / 2.; /* ground & zenith patches */
162 rh_pazi[0] = 0.;
163 rh_dom[0] = 2. * PI;
164 rh_palt[nskypatch - 1] = PI / 2.;
165 rh_pazi[nskypatch - 1] = 0.;
166 rh_dom[nskypatch - 1] = 2. * PI * (1. - cos(alpha * .5));
167 p = 1; /* "normal" patches */
168 for (i = 0; i < NROW * rhsubdiv; i++) {
169 const float ralt = alpha * (i + .5);
170 const int ninrow = tnaz[i / rhsubdiv] * rhsubdiv;
171 const float dom =
172 2. * PI * (sin(alpha * (i + 1)) - sin(alpha * i)) / (double)ninrow;
173 for (j = 0; j < ninrow; j++) {
174 rh_palt[p] = ralt;
175 rh_pazi[p] = 2. * PI * j / (double)ninrow;
176 rh_dom[p++] = dom;
177 }
178 }
179 return nskypatch;
180 #undef NROW
181 }
182
183 /* Resize daylight matrix (GW) */
184 float *resize_dmatrix(float *mtx_data, int nsteps, int npatch) {
185 if (mtx_data == NULL)
186 mtx_data = (float *)malloc(sizeof(float) * NSSAMP * nsteps * npatch);
187 else
188 mtx_data =
189 (float *)realloc(mtx_data, sizeof(float) * NSSAMP * nsteps * npatch);
190 if (mtx_data == NULL) {
191 fprintf(stderr, "%s: out of memory in resize_dmatrix(%d,%d)\n", progname,
192 nsteps, npatch);
193 exit(1);
194 }
195 return (mtx_data);
196 }
197
198 static Atmosphere init_atmos(const double aod, const double grefl) {
199 Atmosphere atmos = {.ozone_density = {.layers =
200 {
201 {.width = 25000.0,
202 .exp_term = 0.0,
203 .exp_scale = 0.0,
204 .linear_term = 1.0 / 15000.0,
205 .constant_term = -2.0 / 3.0},
206 {.width = AH,
207 .exp_term = 0.0,
208 .exp_scale = 0.0,
209 .linear_term = -1.0 / 15000.0,
210 .constant_term = 8.0 / 3.0},
211 }},
212 .rayleigh_density = {.layers =
213 {
214 {.width = AH,
215 .exp_term = 1.0,
216 .exp_scale = -1.0 / HR_MS,
217 .linear_term = 0.0,
218 .constant_term = 0.0},
219 }},
220 .beta_r0 = BR0_MS,
221 .beta_scale = aod / AOD0_CA,
222 .beta_m = NULL,
223 .grefl = grefl};
224 return atmos;
225 }
226
227 static DpPaths get_dppaths(const char *dir, const double aod, const char *mname,
228 const char *tag) {
229 DpPaths paths;
230
231 snprintf(paths.tau, PATH_MAX, "%s%ctau_%s_%s_%.2f.dat", dir, DIRSEP, tag,
232 mname, aod);
233 snprintf(paths.scat, PATH_MAX, "%s%cscat_%s_%s_%.2f.dat", dir, DIRSEP, tag,
234 mname, aod);
235 snprintf(paths.scat1m, PATH_MAX, "%s%cscat1m_%s_%s_%.2f.dat", dir, DIRSEP,
236 tag, mname, aod);
237 snprintf(paths.irrad, PATH_MAX, "%s%cirrad_%s_%s_%.2f.dat", dir, DIRSEP, tag,
238 mname, aod);
239
240 return paths;
241 }
242 static void set_rayleigh_density_profile(Atmosphere *atmos, char *tag,
243 const int is_summer,
244 const double s_latitude) {
245 /* Set rayleigh density profile */
246 if (fabs(s_latitude * 180.0 / PI) > ARCTIC_LAT) {
247 tag[0] = 's';
248 if (is_summer) {
249 tag[1] = 's';
250 atmos->rayleigh_density.layers[0].exp_scale = -1.0 / HR_SS;
251 atmos->beta_r0 = BR0_SS;
252 } else {
253 tag[1] = 'w';
254 atmos->rayleigh_density.layers[0].exp_scale = -1.0 / HR_SW;
255 atmos->beta_r0 = BR0_SW;
256 }
257 } else if (fabs(s_latitude * 180.0 / PI) > TROPIC_LAT) {
258 tag[0] = 'm';
259 if (is_summer) {
260 tag[1] = 's';
261 atmos->rayleigh_density.layers[0].exp_scale = -1.0 / HR_MS;
262 atmos->beta_r0 = BR0_MS;
263 } else {
264 tag[1] = 'w';
265 atmos->rayleigh_density.layers[0].exp_scale = -1.0 / HR_MW;
266 atmos->beta_r0 = BR0_MW;
267 }
268 } else {
269 tag[0] = 't';
270 tag[1] = 'r';
271 atmos->rayleigh_density.layers[0].exp_scale = -1.0 / HR_T;
272 atmos->beta_r0 = BR0_T;
273 }
274 tag[2] = '\0';
275 }
276 /* Add in solar direct to nearest sky patches (GW) */
277 void add_direct(DATARRAY *tau, DATARRAY *scat, DATARRAY *scat1m,
278 DATARRAY *irrad, double ccover, float *parr) {
279 FVECT svec;
280 double near_dprod[NSUNPATCH];
281 int near_patch[NSUNPATCH];
282 double wta[NSUNPATCH], wtot;
283 int i, j, p;
284
285 /* identify nsuns closest patches */
286 for (i = nsuns; i--;)
287 near_dprod[i] = -1.;
288 vectorize(altitude, azimuth, svec);
289 for (p = 1; p < nskypatch; p++) {
290 FVECT pvec;
291 double dprod;
292 vectorize(rh_palt[p], rh_pazi[p], pvec);
293 dprod = DOT(pvec, svec);
294 for (i = 0; i < nsuns; i++)
295 if (dprod > near_dprod[i]) {
296 for (j = nsuns; --j > i;) {
297 near_dprod[j] = near_dprod[j - 1];
298 near_patch[j] = near_patch[j - 1];
299 }
300 near_dprod[i] = dprod;
301 near_patch[i] = p;
302 break;
303 }
304 }
305 /* Get solar radiance */
306 double sun_radiance[NSSAMP] = {0};
307 get_solar_radiance(tau, scat, scat1m, sundir, ER, sun_ct, sun_radiance);
308 if (ccover > 0) {
309 double zenithbr = get_zenith_brightness(sundir);
310 double skybr = get_overcast_brightness(sundir[2], zenithbr);
311 int l;
312 for (l = 0; l < NSSAMP; ++l) {
313 sun_radiance[l] =
314 wmean2(sun_radiance[l], D6415[l] * skybr / WVLSPAN, ccover);
315 }
316 }
317 /* weight by proximity */
318 wtot = 0;
319 for (i = nsuns; i--;)
320 wtot += wta[i] = 1. / (1.002 - near_dprod[i]);
321 /* add to nearest patch radiances */
322 for (i = nsuns; i--;) {
323 float *pdest = parr + NSSAMP * near_patch[i];
324 int k;
325 for (k = 0; k < NSSAMP; k++) {
326 *pdest++ = sun_radiance[k] * wta[i] / wtot;
327 }
328 }
329 }
330
331 void calc_sky_patch_radiance(DATARRAY *scat, DATARRAY *scat1m, double ccover,
332 float *parr) {
333 int i;
334 double mu_sky; /* Sun-sky point azimuthal angle */
335 double sspa; /* Sun-sky point angle */
336 FVECT view_point = {0, 0, ER};
337 for (i = 1; i < nskypatch; i++) {
338 FVECT rdir_sky;
339 int k;
340 vectorize(rh_palt[i], rh_pazi[i], rdir_sky);
341 mu_sky = fdot(view_point, rdir_sky) / ER;
342 sspa = fdot(rdir_sky, sundir);
343 SCOLOR sky_radiance = {0};
344
345 get_sky_radiance(scat, scat1m, ER, mu_sky, sun_ct, sspa, sky_radiance);
346 for (k = 0; k < NSSAMP; ++k) {
347 sky_radiance[k] *= WVLSPAN;
348 }
349
350 if (ccover > 0) {
351 double zenithbr = get_zenith_brightness(sundir);
352 double grndbr = zenithbr * GNORM;
353 double skybr = get_overcast_brightness(rdir_sky[2], zenithbr);
354 int k;
355 for (k = 0; k < NSSAMP; ++k) {
356 sky_radiance[k] = wmean2(sky_radiance[k], skybr * D6415[k], ccover);
357 }
358 }
359
360 for (k = 0; k < NSSAMP; ++k) {
361 parr[NSSAMP * i + k] = sky_radiance[k];
362 }
363 }
364 }
365
366 /* Return maximum of two doubles */
367 static inline double dmax(double a, double b) { return (a > b) ? a : b; }
368
369 /* Compute sky patch radiance values (modified by GW) */
370 void compute_sky(DATARRAY *tau, DATARRAY *scat, DATARRAY *scat1m,
371 DATARRAY *irrad, double ccover, float *parr) {
372 int index; /* Category index */
373 int i;
374 float sun_zenith;
375 SCOLOR sky_radiance = {0};
376 SCOLOR ground_radiance = {0};
377 SCOLR sky_sclr = {0};
378 SCOLR ground_sclr = {0};
379 FVECT view_point = {0, 0, ER};
380 const double radius = VLEN(view_point);
381 const double sun_ct = fdot(view_point, sundir) / radius;
382 const FVECT rdir_grnd = {0, 0, -1};
383 const double mu_grnd = fdot(view_point, rdir_grnd) / radius;
384 const double nu_grnd = fdot(rdir_grnd, sundir);
385 int j;
386
387 /* Calculate sun zenith angle (don't let it dip below horizon) */
388 /* Also limit minimum angle to keep circumsolar off zenith */
389 if (altitude <= 0.0)
390 sun_zenith = DegToRad(90.0);
391 else if (altitude >= DegToRad(87.0))
392 sun_zenith = DegToRad(3.0);
393 else
394 sun_zenith = DegToRad(90.0) - altitude;
395
396 /* Compute ground radiance (include solar contribution if any) */
397 get_ground_radiance(tau, scat, scat1m, irrad, view_point, rdir_grnd, radius,
398 mu_grnd, sun_ct, nu_grnd, grefl, sundir, parr);
399 for (j = 0; j < NSSAMP; j++) {
400 parr[j] *= WVLSPAN;
401 }
402 calc_sky_patch_radiance(scat, scat1m, ccover, parr);
403 }
404
405 int main(int argc, char *argv[]) {
406
407 char buf[256];
408 int doheader = 1; /* output header? */
409 double rotation = 0.0;
410 double elevation = 0;
411 int leap_day = 0; /* add leap day? */
412 int sun_hours_only = 0; /* only output sun hours? */
413 float *mtx_data = NULL;
414 int ntsteps = 0; /* number of time steps */
415 int tstorage = 0; /* number of allocated time steps */
416 int nstored = 0; /* number of time steps in matrix */
417 int last_monthly = 0; /* month of last report */
418 int mo, da;
419 double hr, aod, cc;
420 double dni, dhi;
421 int mtx_offset = 0;
422 int i, j;
423 char lstag[3];
424 char *mie_path = getpath("mie_ca.dat", getrlibpath(), R_OK);
425 char *ddir = ".";
426 char mie_name[20] = "mie_ca";
427 int num_threads = 1;
428 int sorder = 4;
429 int solar_only = 0;
430 int sky_only = 0;
431 FVECT view_point = {0, 0, ER};
432
433 progname = argv[0];
434
435 for (i = 1; i < argc && argv[i][0] == '-'; i++) {
436 switch (argv[i][1]) {
437 case 'd': /* solar (direct) only */
438 solar_only = 1;
439 break;
440 case 's': /* sky only (no direct) */
441 sky_only = 1;
442 break;
443 case 'g':
444 grefl = atof(argv[++i]);
445 break;
446 case 'm':
447 rhsubdiv = atoi(argv[++i]);
448 break;
449 case 'n':
450 num_threads = atoi(argv[++i]);
451 break;
452 case 'r': /* rotate distribution */
453 if (argv[i][2] && argv[i][2] != 'z')
454 goto userr;
455 rotation = atof(argv[++i]);
456 break;
457 case 'u': /* solar hours only */
458 sun_hours_only = 1;
459 break;
460 case 'p':
461 ddir = argv[++i];
462 break;
463 case 'v': /* verbose progress reports */
464 verbose++;
465 break;
466 case 'h': /* turn off header */
467 doheader = 0;
468 break;
469 case '5': /* 5-phase calculation */
470 nsuns = 1;
471 fixed_sun_sa = PI / 360. * atof(argv[++i]);
472 if (fixed_sun_sa <= 0) {
473 fprintf(stderr,
474 "%s: missing solar disk size argument for '-5' option\n",
475 progname);
476 exit(1);
477 }
478 fixed_sun_sa *= fixed_sun_sa * PI;
479 break;
480 case 'o': /* output format */
481 switch (argv[i][2]) {
482 case 'f':
483 case 'd':
484 case 'a':
485 outfmt = argv[i][2];
486 break;
487 default:
488 goto userr;
489 }
490 break;
491 default:
492 goto userr;
493 }
494 }
495 if (i < argc - 1)
496 goto userr;
497 if (i == argc - 1 && freopen(argv[i], "r", stdin) == NULL) {
498 fprintf(stderr, "%s: cannot open '%s' for input\n", progname, argv[i]);
499 exit(1);
500 }
501 if (verbose) {
502 if (i == argc - 1)
503 fprintf(stderr, "%s: reading weather tape '%s'\n", progname, argv[i]);
504 else
505 fprintf(stderr, "%s: reading weather tape from <stdin>\n", progname);
506 }
507 /* read weather tape header */
508 if (scanf("place %[^\r\n] ", buf) != 1)
509 goto fmterr;
510 if (scanf("latitude %lf\n", &s_latitude) != 1)
511 goto fmterr;
512 if (scanf("longitude %lf\n", &s_longitude) != 1)
513 goto fmterr;
514 if (scanf("time_zone %lf\n", &s_meridian) != 1)
515 goto fmterr;
516 if (scanf("site_elevation %lf\n", &elevation) != 1)
517 goto fmterr;
518 if (scanf("weather_data_file_units %d\n", &input) != 1)
519 goto fmterr;
520
521 rh_init();
522 if (verbose) {
523 fprintf(stderr, "%s: location '%s'\n", progname, buf);
524 fprintf(stderr, "%s: (lat,long)=(%.1f,%.1f) degrees north, west\n",
525 progname, s_latitude, s_longitude);
526 if (rotation != 0)
527 fprintf(stderr, "%s: rotating output %.0f degrees\n", progname, rotation);
528 }
529
530 s_latitude = DegToRad(s_latitude);
531 s_longitude = DegToRad(s_longitude);
532 s_meridian = DegToRad(s_meridian);
533 /* initial allocation */
534 mtx_data = resize_dmatrix(mtx_data, tstorage = 2, nskypatch);
535
536 /* Load mie density data */
537 DATARRAY *mie_dp = getdata(mie_path);
538 if (mie_dp == NULL) {
539 fprintf(stderr, "Error reading mie data\n");
540 return 0;
541 }
542
543 while (scanf("%d %d %lf %lf %lf %lf %lf\n", &mo, &da, &hr, &dni, &dhi, &aod,
544 &cc) == 7) {
545 if (aod == 0.0) {
546 aod = AOD0_CA;
547 fprintf(stderr, "aod is zero, using default value %.3f\n", AOD0_CA);
548 }
549 double sda, sta;
550 int sun_in_sky;
551 /* compute solar position */
552 if ((mo == 2) & (da == 29)) {
553 julian_date = 60;
554 leap_day = 1;
555 } else
556 julian_date = jdate(mo, da) + leap_day;
557 sda = sdec(julian_date);
558 sta = stadj(julian_date);
559 altitude = salt(sda, hr + sta);
560 sun_in_sky = (altitude > -DegToRad(SUN_ANG_DEG / 2.));
561
562 azimuth = sazi(sda, hr + sta) + PI - DegToRad(rotation);
563
564 vectorize(altitude, azimuth, sundir);
565 if (sun_hours_only && sundir[2] <= 0.) {
566 continue; /* skipping nighttime points */
567 }
568 sun_ct = fdot(view_point, sundir) / ER;
569
570 mtx_offset = NSSAMP * nskypatch * nstored;
571 nstored += 1;
572 printf("mtx_offset = %d nstored = %d nskypatch = %d\n", mtx_offset, nstored,
573 nskypatch);
574 /* make space for next row */
575 if (nstored > tstorage) {
576 printf("make space for next row nstored = %d tstorage = %d\n", nstored,
577 tstorage);
578 tstorage += (tstorage >> 1) + nstored + 7;
579 mtx_data = resize_dmatrix(mtx_data, tstorage, nskypatch);
580 }
581 ntsteps++; /* keep count of time steps */
582 /* compute sky patch values */
583 Atmosphere clear_atmos = init_atmos(aod, grefl);
584 int is_summer = (mo >= SUMMER_START && mo <= SUMMER_END);
585 if (s_latitude < 0) {
586 is_summer = !is_summer;
587 }
588 set_rayleigh_density_profile(&clear_atmos, lstag, is_summer, s_latitude);
589
590 clear_atmos.beta_m = mie_dp;
591
592 char gsdir[PATH_MAX];
593 size_t siz = strlen(ddir);
594 if (ISDIRSEP(ddir[siz - 1]))
595 ddir[siz - 1] = '\0';
596 snprintf(gsdir, PATH_MAX, "%s%catmos_data", ddir, DIRSEP);
597 if (!make_directory(gsdir)) {
598 fprintf(stderr, "Failed creating atmos_data directory");
599 exit(1);
600 }
601 DpPaths clear_paths = get_dppaths(gsdir, aod, mie_name, lstag);
602
603 if (getpath(clear_paths.tau, ".", R_OK) == NULL ||
604 getpath(clear_paths.scat, ".", R_OK) == NULL ||
605 getpath(clear_paths.scat1m, ".", R_OK) == NULL ||
606 getpath(clear_paths.irrad, ".", R_OK) == NULL) {
607 printf("# Pre-computing...\n");
608 if (!precompute(sorder, clear_paths, &clear_atmos, num_threads)) {
609 fprintf(stderr, "Pre-compute failed\n");
610 return 0;
611 }
612 }
613
614 DATARRAY *tau_clear_dp = getdata(clear_paths.tau);
615 DATARRAY *irrad_clear_dp = getdata(clear_paths.irrad);
616 DATARRAY *scat_clear_dp = getdata(clear_paths.scat);
617 DATARRAY *scat1m_clear_dp = getdata(clear_paths.scat1m);
618
619 if (!solar_only)
620 compute_sky(tau_clear_dp, scat_clear_dp, scat1m_clear_dp, irrad_clear_dp,
621 cc, mtx_data + mtx_offset);
622 if (!sky_only)
623 add_direct(tau_clear_dp, scat_clear_dp, scat1m_clear_dp, irrad_clear_dp,
624 cc, mtx_data + mtx_offset);
625 /* monthly reporting */
626 if (verbose && mo != last_monthly)
627 fprintf(stderr, "%s: stepping through month %d...\n", progname,
628 last_monthly = mo);
629 }
630 freedata(mie_dp);
631 if (!ntsteps) {
632 fprintf(stderr, "%s: no valid time steps on input\n", progname);
633 exit(1);
634 }
635 /* check for junk at end */
636 while ((i = fgetc(stdin)) != EOF)
637 if (!isspace(i)) {
638 fprintf(stderr, "%s: warning - unexpected data past EOT: ", progname);
639 buf[0] = i;
640 buf[1] = '\0';
641 fgets(buf + 1, sizeof(buf) - 1, stdin);
642 fputs(buf, stderr);
643 fputc('\n', stderr);
644 break;
645 }
646 /* write out matrix */
647 if (outfmt != 'a')
648 SET_FILE_BINARY(stdout);
649 #ifdef getc_unlocked
650 flockfile(stdout);
651 #endif
652 if (verbose)
653 fprintf(stderr, "%s: writing %smatrix with %d time steps...\n", progname,
654 outfmt == 'a' ? "" : "binary ", nstored);
655 if (doheader) {
656 newheader("RADIANCE", stdout);
657 printargs(argc, argv, stdout);
658 printf("LATLONG= %.8f %.8f\n", RadToDeg(s_latitude),
659 -RadToDeg(s_longitude));
660 printf("NROWS=%d\n", nskypatch);
661 printf("NCOLS=%d\n", nstored);
662 printf("NCOMP=%d\n", NSSAMP);
663 if ((outfmt == 'f') | (outfmt == 'd'))
664 fputendian(stdout);
665 fputformat((char *)getfmtname(outfmt), stdout);
666 putchar('\n');
667 }
668 /* patches are rows (outer sort) */
669 for (i = 0; i < nskypatch; i++) {
670 mtx_offset = NSSAMP * i;
671 switch (outfmt) {
672 case 'a':
673 for (j = 0; j < nstored; j++) {
674 int k;
675 for (k = 0; k < NSSAMP; k++) {
676 printf("%.3g ", mtx_data[mtx_offset + k]);
677 }
678 printf("\n");
679 mtx_offset += NSSAMP * nskypatch;
680 }
681 if (nstored > 1)
682 fputc('\n', stdout);
683 break;
684 case 'f':
685 for (j = 0; j < nstored; j++) {
686 putbinary(mtx_data + mtx_offset, sizeof(float), NSSAMP, stdout);
687 mtx_offset += NSSAMP * nskypatch;
688 }
689 break;
690 case 'd':
691 for (j = 0; j < nstored; j++) {
692 double ment[NSSAMP];
693 for (j = 0; j < NSSAMP; j++)
694 ment[j] = mtx_data[mtx_offset + j];
695 putbinary(ment, sizeof(double), NSSAMP, stdout);
696 mtx_offset += NSSAMP * nskypatch;
697 }
698 break;
699 }
700 if (ferror(stdout))
701 goto writerr;
702 }
703 alldone:
704 if (fflush(NULL) == EOF)
705 goto writerr;
706 if (verbose)
707 fprintf(stderr, "%s: done.\n", progname);
708 exit(0);
709 userr:
710 fprintf(stderr,
711 "Usage: %s [-v][-h][-A][-d|-s|-n][-u][-D file [-M modfile]][-r "
712 "deg][-m N][-g r g b][-c r g b][-o{f|d}][-O{0|1}] [tape.wea]\n",
713 progname);
714 exit(1);
715 fmterr:
716 fprintf(stderr, "%s: weather tape format error in header\n", progname);
717 exit(1);
718 writerr:
719 fprintf(stderr, "%s: write error on output\n", progname);
720 exit(1);
721 }