ViewVC Help
View File | Revision Log | Show Annotations | Download File | Root Listing
root/radiance/ray/src/gen/gensdaymtx.c
Revision: 1.1
Committed: Fri Aug 2 18:59:26 2024 UTC (8 months, 4 weeks ago) by greg
Content type: text/plain
Branch: MAIN
Log Message:
feat(gensdaymtx,epw2wea,genssky): Taoning added new gensdaymtx and updated others

File Contents

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