OB.DAAC Logo
NASA Logo
Ocean Color Science Software

ocssw V2022
l1_aci_hdf.c
Go to the documentation of this file.
1 /* ============================================================================ */
2 /* module l1_aci_hdf.c - functions to read Pathfinder Class data for MSL12 */
3 /* */
4 /* Modified By: S. Walsh, RSMAS, January 2008 to read Pathfinder */
5 /* */
6 /* ============================================================================ */
7 
8 #include "l1_aci_hdf.h"
9 #include <math.h>
10 #include "l1.h"
11 #include "rawcal.h"
12 #include "runcal.h"
13 #include <hdf4utils.h>
14 
15 #include <hdf.h>
16 #include <mfhdf.h>
17 
18 #define NAVBANDS 5
19 #define NVIRBANDS 2
20 #define NTIRBANDS 3
21 
22 #define TEMP_REF 273.15 /* convert Bt's from Deg K to Deg C */
23 
24 #define NRADTAB 100 /* Radiance to temperature table size */
25 
26 #define SDNAME "Channel_"
27 //#define SDNAME ""
28 
29 static int32_t sd_id;
30 static int32_t sd_id_g;
31 static int32_t spix = 0;
32 //static int32_t sscan = 0; /* start scan in original granule */
33 
34 static float T_LOW = 185.; /* Limit of 180K in etbpsub.f */
35 static float T_HIGH = 370.; /* Limit of 375K in etbpsub.f */
36 
37 static int pdbg1 = -1;
38 static int ldbg1 = -1;
39 //static int pdbg1 = 174; /* 2003131211434.N16 */
40 //static int ldbg1 = 1952;
41 //static int pdbg1 = 141; /* 2009200225057.N16 */
42 //static int ldbg1 = 122;
43 
44 /* ----------------------------------------------------------------------------------- */
45 /* openl1_aci_hdf() - opens a ACI L1B file for reading. */
46 /* */
47 /* B. Franz, SAIC, February 2003. */
48 
49 /* ----------------------------------------------------------------------------------- */
50 int openl1_aci_hdf(filehandle *file) {
51  int32_t npix;
52  int32_t nscan;
53  const char xsatida[10];
54 
55  int32_t sds_id;
56  int32_t rank;
57  int32_t dims[3];
58  int32_t type;
59  int32_t numattr;
60 
61  /* Open the HDF input file */
62  sd_id = SDstart(file->name, DFACC_RDONLY);
63  if (sd_id == FAIL) {
64  fprintf(stderr, "-E- %s line %d: SDstart(%s, %d) failed.\n",
65  __FILE__, __LINE__, file->name, DFACC_RDONLY);
66  return (HDF_FUNCTION_ERROR);
67  }
68 
69  /* Get pixel and scan dimensions */
70  sds_id = SDselect(sd_id, SDnametoindex(sd_id, (SDNAME"1")));
71  if (SDgetinfo(sds_id, NULL, &rank, dims, &type, &numattr) == -1) {
72  fprintf(stderr, "-E- %s line %d: error getting dimension info.\n",
73  __FILE__, __LINE__);
74  return (HDF_FUNCTION_ERROR);
75  }
76  npix = dims[1];
77  /* can't use dims for nscan because avhrr file is created before the */
78  /* actual number of scans is known - because missing scans need to */
79  /* be filled in */
80  // nscan = dims[0];
81  if (getHDFattr(sd_id, "Number of scans", "", (VOIDP) & nscan) != 0) {
82  printf("-E- %s line %d: Error reading Number of scans attribute.\n",
83  __FILE__, __LINE__);
84  return (1);
85  }
86 
87  if (getHDFattr(sd_id, "Orbit Number", "", (VOIDP) & file->orbit_number)
88  != 0) {
89  printf("-E- %s line %d: Error reading Orbit attribute.\n",
90  __FILE__, __LINE__);
91  return (1);
92  }
93 
94  /* Open the HDF geolocation input file.*/
95 
96  sd_id_g = SDstart(file->geofile, DFACC_RDONLY);
97  if (sd_id_g == FAIL) {
98  printf("Error opening geolocation file.\n");
99  fprintf(stderr, "-E- %s line %d: SDstart(%s, %d) failed.\n",
100  __FILE__, __LINE__, file->geofile, DFACC_RDONLY);
101  return (HDF_FUNCTION_ERROR);
102  }
103 
104  if (getHDFattr(sd_id, "Satellite", "", (VOIDP) xsatida) != 0) {
105  printf("-E- %s line %d: Error reading Satellite attribute.\n",
106  __FILE__, __LINE__);
107  return (1);
108  }
109 
110  if (getHDFattr(sd_id, "spatial_resolution", "", (VOIDP) & file->spatialResolution) != 0) {
111  printf("-W- %s line %d: Error reading spatial_resolution attribute.\n",
112  __FILE__, __LINE__);
113  }
114  file->subsensorID = satname2xsatid(xsatida);
115 
116  file->npix = npix;
117  file->ndets = 1;
118  file->nscan = nscan;
119  file->sd_id = sd_id;
120 
121  return (LIFE_IS_GOOD);
122 }
123 
124 /* ----------------------------------------------------------------------------------- */
125 /* readl1_aci_hdf() - reads 1 line (scan) from a ACI ingested CLASS file, loads l1rec. */
126 /* */
127 
128 /* ----------------------------------------------------------------------------------- */
129 int readl1_aci_hdf(filehandle *file, int32_t scan, l1str *l1rec) {
130  static int firstCall = 1;
131  static int firstScan = 1;
132 
133  static int16_t *ch1 = NULL;
134  static int16_t *ch2 = NULL;
135  static int16_t *ch3 = NULL;
136  static int16_t *ch4 = NULL;
137  static int16_t *ch5 = NULL;
138 
139  static float m[NAVBANDS];
140  static float b[NAVBANDS];
141 
142  static float rayly[2];
143  static float aersol[2];
144  static int year, jday, RelDay, LY, leap, month, day, hh, mm, ss;
145  static int32_t ms = 0;
146  static char starttime[17];
147 
148  static int startOfMonth[2][12] = {
149  { 0, 31, 59, 90, 120, 151, 181, 212, 243,
150  273, 304, 334},
151  { 0, 31, 60, 91, 121, 152, 182, 213, 244, 274, 305,
152  335}
153  };
154 
155  /* Jon Mittaz's calibration for ch4 and ch5 */
156  static double ch4_coeff[5] = {0.0351618, 1.24004, 0.979493, 8.6071114e-06,
157  8.6071114e-06};
158  static double ch5_coeff[5] = {0.0326272, 1.06256, 0.986256, 4.0256298e-06,
159  4.0256298e-06};
160 
161  static char ingester[5];
162 
163  int16_t *data = NULL;
164  int32_t npix = (int32_t) file->npix;
165  int32_t nbands = (int32_t) file->nbands;
166  int32_t nbandsir = (int32_t) file->nbandsir;
167  int32_t detnum = 1;
168  int32_t ich, iw, ip, ipb;
169  int32_t ii, jj;
170  int32_t dims[3];
171  /* lunin is used to open the noaa*.cal file. It's an almost random
172  * number chosen by searching for 'open' in the other fortran routines
173  * and selecting an unused number */
174  int32_t lunin = 23;
175 
176  static float radinc[3];
177  static float radbas[3];
178  static float radinv[3][NRADTAB];
179  static float teminc[3];
180  static float tembas[3];
181  static float teminv[3][NRADTAB];
182 
183  static int vbcksc[30];
184  static int vspace[50];
185  float rspc[3];
186  float radprt[3];
187  float avgbck[3];
188  float avgspc[3];
189  float64 sumbck[2][3];
190  float64 sumspc[2][5];
191 
192  char scaletype[128];
193  char startdate[15];
194  char satnam[12];
195 
196  static RAW_CAL *rawrec;
197  static RUN_CAL *runrec;
198 
199  static int32_t *scantime;
200 
201  char tmpstr[80];
202 
203  int32_t status;
204  int32_t irun, jrun, idif, iraw, jraw, ni;
205  float ai;
206  float ascan;
207  float percnt;
208  float trad, ttem, tnlc;
209  float tcor;
210  float prtemp;
211 
212  if (firstCall) {
213 
214  /* initialize index arrays */
215 
216  firstCall = 0;
217 
218  jj = 0;
219  for (ii = 0; ii < 50; ii++) {
220  vspace[ii] = jj;
221  jj += 5;
222  if (jj >= 50) {
223  jj = jj - 50 + 1;
224  }
225  }
226  jj = 0;
227  for (ii = 0; ii < 30; ii++) {
228  vbcksc[ii] = jj;
229  jj += 3;
230  if (jj >= 30) {
231  jj = jj - 30 + 1;
232  }
233  }
234 
235  if ((ch1 = (int16_t *) calloc(npix, sizeof (int16_t))) == NULL) {
236  printf("-E- %s line %d: Error allocating data space.\n",
237  __FILE__, __LINE__);
238  return (1);
239  }
240 
241  if (getHDFattr(sd_id, "Slope", SDNAME"1", (VOIDP) & m[0]) != 0) {
242  printf("-E- %s line %d: Error reading Channel_1 Slope attribute.\n",
243  __FILE__, __LINE__);
244  return (1);
245  }
246 
247  if (getHDFattr(sd_id, "Intercept", SDNAME"1", (VOIDP) & b[0]) != 0) {
248  printf(
249  "-E- %s line %d: Error reading Channel_1 Intercept attribute.\n",
250  __FILE__, __LINE__);
251  return (1);
252  }
253  if (getHDFattr(sd_id, "Scale_type", SDNAME"1", (VOIDP) & scaletype)
254  != 0) {
255  printf(
256  "-E- %s line %d: Error reading Channel_1 Scale_type attribute.\n",
257  __FILE__, __LINE__);
258  return (1);
259  }
260  if (strcmp(scaletype, "y= Slope * x + Intercept;") != 0) {
261  printf("-E- %s line %d: Channel_1 Scale_type must be linear.\n",
262  __FILE__, __LINE__);
263  return (1);
264  }
265 
266  if ((ch2 = (int16_t *) calloc(npix, sizeof (int16_t))) == NULL) {
267  printf("-E- %s line %d: Error allocating data space.\n",
268  __FILE__, __LINE__);
269  return (1);
270  }
271 
272  if (getHDFattr(sd_id, "Slope", SDNAME"2", (VOIDP) & m[1]) != 0) {
273  printf("-E- %s line %d: Error reading Channel_2 Slope attribute.\n",
274  __FILE__, __LINE__);
275  return (1);
276  }
277 
278  if (getHDFattr(sd_id, "Intercept", SDNAME"2", (VOIDP) & b[1]) != 0) {
279  printf(
280  "-E- %s line %d: Error reading Channel_2 Intercept attribute.\n",
281  __FILE__, __LINE__);
282  return (1);
283  }
284  if (getHDFattr(sd_id, "Scale_type", SDNAME"2", (VOIDP) & scaletype)
285  != 0) {
286  printf(
287  "-E- %s line %d: Error reading Channel_2 Scale_type attribute.\n",
288  __FILE__, __LINE__);
289  return (1);
290  }
291  if (strcmp(scaletype, "y= Slope * x + Intercept;") != 0) {
292  printf("-E- %s line %d: Channel_2 Scale_type must be linear.\n",
293  __FILE__, __LINE__);
294  return (1);
295  }
296 
297  if ((ch3 = (int16_t *) calloc(npix, sizeof (int16_t))) == NULL) {
298  printf("-E- %s line %d: Error allocating data space.\n",
299  __FILE__, __LINE__);
300  return (1);
301  }
302 
303  if (getHDFattr(sd_id, "Slope", SDNAME"3", (VOIDP) & m[2]) != 0) {
304  printf("-E- %s line %d: Error reading Channel_3 Slope attribute.\n",
305  __FILE__, __LINE__);
306  return (1);
307  }
308 
309  if (getHDFattr(sd_id, "Intercept", SDNAME"3", (VOIDP) & b[2]) != 0) {
310  printf(
311  "-E- %s line %d: Error reading Channel_3 Intercept attribute.\n",
312  __FILE__, __LINE__);
313  return (1);
314  }
315  if (getHDFattr(sd_id, "Scale_type", SDNAME"3", (VOIDP) & scaletype)
316  != 0) {
317  printf(
318  "-E- %s line %d: Error reading Channel_3 Scale_type attribute.\n",
319  __FILE__, __LINE__);
320  return (1);
321  }
322  if (strcmp(scaletype, "y= Slope * x + Intercept;") != 0) {
323  printf("-E- %s line %d: Channel_3 Scale_type must be linear.\n",
324  __FILE__, __LINE__);
325  return (1);
326  }
327 
328  if ((ch4 = (int16_t *) calloc(npix, sizeof (int16_t))) == NULL) {
329  printf("-E- %s line %d: Error allocating data space.\n",
330  __FILE__, __LINE__);
331  return (1);
332  }
333 
334  if (getHDFattr(sd_id, "Slope", SDNAME"4", (VOIDP) & m[3]) != 0) {
335  printf("-E- %s line %d: Error reading Channel_4 Slope attribute.\n",
336  __FILE__, __LINE__);
337  return (1);
338  }
339 
340  if (getHDFattr(sd_id, "Intercept", SDNAME"4", (VOIDP) & b[3]) != 0) {
341  printf(
342  "-E- %s line %d: Error reading Channel_4 Intercept attribute.\n",
343  __FILE__, __LINE__);
344  return (1);
345  }
346  if (getHDFattr(sd_id, "Scale_type", SDNAME"4", (VOIDP) & scaletype)
347  != 0) {
348  printf(
349  "-E- %s line %d: Error reading Channel_4 Scale_type attribute.\n",
350  __FILE__, __LINE__);
351  return (1);
352  }
353  if (strcmp(scaletype, "y= Slope * x + Intercept;") != 0) {
354  printf("-E- %s line %d: Channel_4 Scale_type must be linear.\n",
355  __FILE__, __LINE__);
356  return (1);
357  }
358 
359  if ((ch5 = (int16_t *) calloc(npix, sizeof (int16_t))) == NULL) {
360  printf("-E- %s line %d: Error allocating data space.\n",
361  __FILE__, __LINE__);
362  return (1);
363  }
364 
365  if (getHDFattr(sd_id, "Slope", SDNAME"5", (VOIDP) & m[4]) != 0) {
366  printf("-E- %s line %d: Error reading Channel_5 Slope attribute.\n",
367  __FILE__, __LINE__);
368  return (1);
369  }
370 
371  if (getHDFattr(sd_id, "Intercept", SDNAME"5", (VOIDP) & b[4]) != 0) {
372  printf(
373  "-E- %s line %d: Error reading Channel_5 Intercept attribute.\n",
374  __FILE__, __LINE__);
375  return (1);
376  }
377  if (getHDFattr(sd_id, "Scale_type", SDNAME"5", (VOIDP) & scaletype)
378  != 0) {
379  printf(
380  "-E- %s line %d: Error reading Channel_5 Scale_type attribute.\n",
381  __FILE__, __LINE__);
382  return (1);
383  }
384  if (strcmp(scaletype, "y= Slope * x + Intercept;") != 0) {
385  printf("-E- %s line %d: Channel_5 Scale_type must be linear.\n",
386  __FILE__, __LINE__);
387  return (1);
388  }
389 
390  if (getHDFattr(sd_id, "start date", "", (VOIDP) & startdate) != 0) {
391  printf("-E- %s line %d: Error reading start date attribute.\n",
392  __FILE__, __LINE__);
393  return (1);
394  }
395  /* convert yyyy-mm-dd to year, day */
396  sscanf(startdate, "%04d-%02d-%02d", &year, &month, &day);
397  if (((year % 400) == 0) || (((year % 4) == 0) && ((year % 100) != 0))) {
398  leap = 1;
399  LY = 366;
400  } else {
401  leap = 0;
402  LY = 365;
403  }
404 
405  jday = startOfMonth[leap][month - 1] + day;
406 
407  if (getHDFattr(sd_id, "start time", "", (VOIDP) & starttime) != 0) {
408  printf("-E- %s line %d: Error reading start time attribute.\n",
409  __FILE__, __LINE__);
410  return (1);
411  }
412  /* starttime is a float with a value in this form: hhmmss.ss */
413  // hh = starttime/10000.;
414  // mm = ((int)(starttime/100.) % 100);
415  // ss = ((int)(starttime-(100.0 * (mm+ (100.*hh)))) % 100);
416  // ms = (int32_t)(1000. * (starttime- (int)starttime));
417  /* starttime is "hh:mm:ss.sss UTC" */
418  sscanf(starttime, "%02d:%02d:%02d.%03d", &hh, &mm, &ss, &ms);
419 
420  if (getHDFattr(sd_id, "Satellite Name", "", (VOIDP) satnam) != 0) {
421  printf("-E- %s line %d: Error reading Satellite Name attribute.\n",
422  __FILE__, __LINE__);
423  return (1);
424  }
425 
426  if (getHDFattr(sd_id, "PRTEMP", "", (VOIDP) & prtemp) != 0) {
427  printf("-E- %s line %d: Error reading PRTEMP attribute.\n",
428  __FILE__, __LINE__);
429  return (1);
430  }
431  //printf(" I Mean baseplate temperature =%6.1fK\n",prtemp);
432  // Sue Walsh said to comment the ACI name out...
433  // strcpy(ingester, "ACI ");
434  // if (getDims(sd_id, "ACI Raw Running Calibration", dims) != 0) {
435  // printf(
436  // "-E- %s line %d: Error reading ACI Raw Running Calibration dim.\n",
437  // __FILE__, __LINE__);
438  // return (1);
439  // }
440  /* allocate number of bytes needed for data in this file */
441  if ((rawrec = (RAW_CAL *) malloc(dims[0])) == NULL) {
442  printf("-E- %s line %d: Error allocating rawcal data space.\n",
443  __FILE__, __LINE__);
444  return (1);
445  }
446  strcpy(tmpstr, ingester);
447  strcat(tmpstr, "Raw Running Calibration");
448  status = rdSDS(sd_id, tmpstr, 0, 0, 0, 0, (VOIDP) rawrec);
449  if (status != 0) {
450  printf("-E- %s line %d: Error reading %s.\n",
451  __FILE__, __LINE__, tmpstr);
452  return (1);
453  }
454 
455  strcpy(tmpstr, ingester);
456  strcat(tmpstr, "Running Calibration");
457  if (getDims(sd_id, tmpstr, dims) != 0) {
458  printf("-E- %s line %d: Error reading %s dim.\n",
459  __FILE__, __LINE__, tmpstr);
460  return (1);
461  }
462  /* allocate the number of bytes needed for the data in this file */
463  if ((runrec = (RUN_CAL *) malloc(dims[0])) == NULL) {
464  printf("-E- %s line %d: Error allocating runcal data space.\n",
465  __FILE__, __LINE__);
466  return (1);
467  }
468  status = rdSDS(sd_id, tmpstr, 0, 0, 0, 0, (VOIDP) runrec);
469 
470  /* Work-around for single-scan time errors */
471  // if (TAIsec < 0.0) {
472  // printf("-W- %s: bad time in geolocation file at frame %d, using previous.\n",
473  // __FILE__,(int)frame);
474  // TAIsec = lastTAIsec;
475  // } else
476  // lastTAIsec = TAIsec;
477  //
478  // usec = TAIsec + 725846400.0;
479  // unix2yds(usec,&year,&day,&dsec);
480  if (getDims(sd_id, "Scan start time", dims) != 0) {
481  printf("-E- %s line %d: Error reading Scan start time dim.\n",
482  __FILE__, __LINE__);
483  return (1);
484  }
485  if (dims[0] != l1rec->l1file->nscan) {
486  printf(
487  "-E- %s line %d: Error: Number of Scan start times (%d) must match Number of Scans (%d).\n",
488  __FILE__, __LINE__, dims[0], l1rec->l1file->nscan);
489  return (1);
490  }
491  if ((scantime = (int32_t *) malloc(dims[0] * sizeof (int32_t))) == NULL) {
492  printf("-E- %s line %d: Error allocating scantime data space.\n",
493  __FILE__, __LINE__);
494  return (1);
495  }
496  status = rdSDS(sd_id, "Scan start time", 0, 0, 0, 0, (VOIDP) scantime);
497  if (status != 0) {
498  printf("-E- %s line %d: Error reading scan start time\n",
499  __FILE__, __LINE__);
500  return (1);
501  }
502 
503  // /* degc option is now ONLY for sstref in nlsst* equations */
504  // if (file->input->degc == 1) {
505  // /* Convert to Deg C for calc */
506  // rabias = TEMP_REF;
507  // } else {
508  // /* Do calc in Deg K */
509  // rabias = 0.;
510  // }
511 
512  switch (file->subsensorID) {
513  case NO07:
514  printf("I Applying NOAA-7 visible degradation.\n");
515  m[0] = m[0] / (0.916 - 0.049 * ((float) year + ((float) jday / (float) LY) - 1981.5)
516  + 0.0050 * pow((float) year + ((float) jday / (float) LY) - 1981.5, 2));
517  b[0] = b[0] / (0.916 - 0.049 * ((float) year + ((float) jday / (float) LY) - 1981.5)
518  + 0.0050 * pow((float) year + ((float) jday / (float) LY) - 1981.5, 2));
519  m[1] = m[1] / (0.882 - 0.080 * ((float) year + ((float) jday / (float) LY) - 1981.5)
520  + 0.0110 * pow((float) year + ((float) jday / (float) LY) - 1981.5, 2));
521  b[1] = b[1] / (0.882 - 0.080 * ((float) year + ((float) jday / (float) LY) - 1981.5)
522  + 0.0110 * pow((float) year + ((float) jday / (float) LY) - 1981.5, 2));
523  break;
524  case NO09:
525  printf("I Applying NOAA-9 visible degradation.\n");
526  m[0] = m[0] / (0.953 - 0.051 * ((float) year + ((float) jday / (float) LY) - 1985.0));
527  b[0] = b[0] / (0.953 - 0.051 * ((float) year + ((float) jday / (float) LY) - 1985.0));
528  m[1] = m[1] / (0.866 - 0.026 * ((float) year + ((float) jday / (float) LY) - 1985.0));
529  b[1] = b[1] / (0.866 - 0.026 * ((float) year + ((float) jday / (float) LY) - 1985.0));
530  break;
531  case NO11:
532  printf("I Applying NOAA-11 visible degradation.\n");
533  m[0] = m[0] / (0.797 - 0.010 * ((float) year + ((float) jday / (float) LY) - 1989.0));
534  b[0] = b[0] / (0.797 - 0.010 * ((float) year + ((float) jday / (float) LY) - 1989.0));
535  m[1] = m[1] / (0.683 - 0.020 * ((float) year + ((float) jday / (float) LY) - 1989.0));
536  b[1] = b[1] / (0.683 - 0.020 * ((float) year + ((float) jday / (float) LY) - 1989.0));
537  break;
538  case NO14:
539  /* Dec 30, 1994 is base date */
540  RelDay = ((year - 1994) * 365 + (jday - 364));
541  printf("I Applying NOAA-14 visible degradation.\n");
542  /* Post-Launch Calibration for Noaa-14 */
543  /* C R Nagaraja Rao and Jianhua Chen */
544  /* NOAA/NESDIS */
545  /* (from web page, March 1997) */
546  m[0] = 0.0000232 * (float) RelDay + 0.109;
547  b[0] = -41.0 * m[0];
548  m[1] = 0.0000373 * (float) RelDay + 0.129;
549  b[1] = -41.0 * m[1];
550  break;
551  default:
552  printf("I No visible degradation available for NOAA-%2.0d.\n",
553  file->subsensorID);
554  break;
555 
556  }
557 
558  /* start running calibration */
559 
560  etloadresp_(&lunin, satnam);
561  /* for debug / info
562  printf(" t_low %10.4f\n",T_LOW);
563  printf(" t_high %10.4f\n",T_HIGH);
564  */
565  for (iw = 0; iw <= 2; iw++) {
566  ich = iw + 3;
567  radinc[iw] = (float) ((log10(etintegrate_(&ich, &T_HIGH))
568  - log10(etintegrate_(&ich, &T_LOW)))
569  / (float) (NRADTAB - 1));
570  radbas[iw] = (float) (log10(etintegrate_(&ich, &T_LOW)));
571  /* for debug / info
572  printf(" Ch %1d radbas %14.6E radinc %14.6E\n",ich,
573  radbas[iw],radinc[iw]);
574  printf(" Ch Idx radinv (LOG)\n");
575  */
576  radinv[iw][0] = T_LOW;
577  /* for debug / info
578  printf(" %2d %3d %14.6E\n",ich, 1, radinv[iw][0]);
579  */
580  for (ip = 1; ip < NRADTAB; ip++) {
581  ttem = pow(10., radbas[iw] + radinc[iw] * (float) (ip));
582  radinv[iw][ip] = etinvert_(&ich, &ttem);
583  /* for debug / info
584  printf(" %2d %3d %14.6E\n",ich, ip+1, radinv[iw][ip]);
585  */
586  }
587  teminc[iw] = (T_HIGH - T_LOW) / (float) (NRADTAB - 1);
588  tembas[iw] = T_LOW;
589  /* for debug / info
590  printf(" Ch %1d tembas %14.6E teminc %14.6E\n",ich,
591  tembas[iw],teminc[iw]);
592  printf(" Ch Idx teminv\n");
593  */
594  teminv[iw][0] = etintegrate_(&ich, &radinv[iw][0]);
595  /* for debug / info
596  printf(" %2d %3d %14.6E\n",ich, 1, teminv[iw][0]);
597  */
598  for (ip = 1; ip < NRADTAB; ip++) {
599  ttem = tembas[iw] + (teminc[iw] * (float) (ip));
600  teminv[iw][ip] = etintegrate_(&ich, &ttem);
601  /* for debug / info
602  printf(" %2d %3d %14.6E\n",ich, ip+1, teminv[iw][ip]);
603  */
604  }
605  }
606 
607  status = avconsh_(&lunin, &npix, &jday);
608 
609  /* end first */
610  }
611 
612  /* if (scan == ldbg1) { */
613  /* printf(" got to scan 1780\n"); */
614  /* } */
615  /* check to see if this scan line is a different day than the previous */
616  if (scan > 0 && (scantime[scan] < scantime[scan - 1])) {
617  /* crossed to the next day */
618  jday = jday + 1;
619  if (jday > LY) {
620  /* crossed to the next year */
621  jday -= LY;
622  year = year + 1;
623  }
624  }
625 
626  l1rec->scantime = yds2unix(year, jday, (double) (scantime[scan] / 1.e3));
627 
628  /* Get position and path geometry */
629  READ_SDS_ID(sd_id_g, "Longitude", l1rec->lon, scan, spix, 0, 0, 1, npix, 1,
630  1);
631  READ_SDS_ID(sd_id_g, "Latitude", l1rec->lat, scan, spix, 0, 0, 1, npix, 1,
632  1);
633  READ_SDS_ID(sd_id_g, "Solar Zenith Angle", l1rec->solz, scan, spix, 0, 0, 1,
634  npix, 1, 1);
635  READ_SDS_ID(sd_id_g, "Sensor Zenith Angle", l1rec->senz, scan, spix, 0, 0,
636  1, npix, 1, 1);
637  READ_SDS_ID(sd_id_g, "Relative Azimuth Angle", l1rec->delphi, scan, spix, 0,
638  0, 1, npix, 1, 1);
639 
640  for (ip = 0; ip < npix; ip++) {
641 
642  l1rec->pixnum[ip] = spix + ip;
643 
644  if (l1rec->lon[ip] < -181.0 || l1rec->lon[ip] > 181.0
645  || l1rec->lat[ip] < -91.0 || l1rec->lat[ip] > 91.0)
646  l1rec->navfail[ip] = 1;
647  }
648 
649  /* Read L1B data, scale to radiance, and copy relevant bands to L1 record */
650  READ_SDS_ID(sd_id, SDNAME"1", ch1, scan, spix, 0, 0, 1, npix, 1, 1);
651  READ_SDS_ID(sd_id, SDNAME"2", ch2, scan, spix, 0, 0, 1, npix, 1, 1);
652  READ_SDS_ID(sd_id, SDNAME"3", ch3, scan, spix, 0, 0, 1, npix, 1, 1);
653  READ_SDS_ID(sd_id, SDNAME"4", ch4, scan, spix, 0, 0, 1, npix, 1, 1);
654  READ_SDS_ID(sd_id, SDNAME"5", ch5, scan, spix, 0, 0, 1, npix, 1, 1);
655 
656  /* This check should somehow check the scan line number to see if it's
657  * the first or last boxsiz/2 lines at the beginning or end of the file.
658  *
659  * I'm assuming the boxsiz is 3 so only the first line doesn't get
660  * processed and that this routine is called for each scan line in order
661  * from first to last.
662  *
663  */
664 
665  if (firstScan) {
666  firstScan = 0;
667  /* what else happens on first line? */
668  return (LIFE_IS_GOOD);
669  }
670 
671  ascan = (float) scan + 1.0; /* ascan is one based */
672  irun = -1;
673  idif = 99999;
674 
675  /* numcal and numraw should always be the same */
676  for (ip = 0; ip < runrec->numcal; ip++) {
677  if (abs((int)(runrec->runcal[ip].cenlin - ascan)) < idif) {
678  idif = abs((int)(runrec->runcal[ip].cenlin - ascan));
679  irun = ip;
680  }
681  }
682 
683  if (irun == -1 || idif > runrec->intrvl) {
684  printf("-W- %s line %d: Hey, I can't find the RUNCAL entry!\n",
685  __FILE__, __LINE__);
686  printf("-W- scan: %f; first center scan: %f; last: %f", ascan,
687  runrec->runcal[0].cenlin,
688  runrec->runcal[runrec->numcal - 1].cenlin);
689  if (irun == -1) {
690  if (ascan < runrec->runcal[1].cenlin) {
691  irun = 0;
692  } else {
693  irun = runrec->numcal - 1;
694  }
695  }
696  }
697  if (ascan > runrec->runcal[irun].cenlin) {
698  if (irun < runrec->numcal - 1) {
699  jrun = irun + 1; /* next available entry */
700  } else {
701  jrun = irun; /* duplicate */
702  }
703  } else {
704  jrun = irun; /* next entry */
705  if (irun > 0) {
706  irun = irun - 1; /* previous entry */
707  }
708  }
709 
710  iraw = irun; /* should be the same */
711  jraw = jrun;
712 
713  if (runrec->runcal[irun].cenlin != runrec->runcal[jrun].cenlin) {
714  percnt = (ascan - runrec->runcal[irun].cenlin)
715  / (runrec->runcal[jrun].cenlin - runrec->runcal[irun].cenlin);
716  if (percnt < 0.0) {
717  percnt = 0.0; /* earlier than first set */
718  } else if (percnt > 1.0) {
719  percnt = 1.0; /* later than last set */
720  }
721  } else {
722  percnt = 0.0; /* use first set, duplicates */
723  }
724 
725  /* Decommutate raw telemetry data for sensor calibration */
726 
727  iw = 0; /* Force initialization for group */
728  for (ii = 1; ii <= 30; ii++) {
729  if (ii > iw * 10) {
730  iw++; /* Initialize for next group */
731  for (jj = 0; jj < 2; jj++) {
732  sumbck[jj][iw - 1] = 0.0;
733  }
734  }
735  sumbck[0][iw - 1] += rawrec->rawcal[iraw].bckscn[0][vbcksc[ii - 1]]
736  * (1.0 - percnt)
737  + rawrec->rawcal[jraw].bckscn[0][vbcksc[ii - 1]] * (percnt);
738  sumbck[1][iw - 1] += rawrec->rawcal[iraw].bckscn[1][vbcksc[ii - 1]]
739  * (1.0 - percnt)
740  + rawrec->rawcal[jraw].bckscn[1][vbcksc[ii - 1]] * (percnt);
741  }
742  for (iw = 0; iw < 3; iw++) {
743  if (sumbck[0][iw] > 0.0) {
744  avgbck[iw] = sumbck[1][iw] / sumbck[0][iw];
745  } else {
746  avgbck[iw] = 0.0;
747  }
748  }
749 
750  iw = 0; /* Force initialization for group 1 */
751  for (ii = 1; ii <= 50; ii++) {
752  if (ii > iw * 10) {
753  iw++; /* Initialize for next group */
754  for (jj = 0; jj < 2; jj++) {
755  sumspc[jj][iw - 1] = 0.0;
756  }
757  }
758  sumspc[0][iw - 1] += rawrec->rawcal[iraw].space[0][vspace[ii - 1]]
759  * (1.0 - percnt)
760  + rawrec->rawcal[jraw].space[0][vspace[ii - 1]] * (percnt);
761  sumspc[1][iw - 1] += rawrec->rawcal[iraw].space[1][vspace[ii - 1]]
762  * (1.0 - percnt)
763  + rawrec->rawcal[jraw].space[1][vspace[ii - 1]] * (percnt);
764  }
765  for (iw = 0; iw < 3; iw++) {
766  if (sumspc[0][iw + 2] > 0.0) {
767  avgspc[iw] = sumspc[1][iw + 2] / sumspc[0][iw + 2];
768 
769  } else {
770  avgspc[iw] = 0.0;
771  }
772  }
773 
774  prtemp = (1.0 - percnt) * runrec->runcal[irun].prtemp
775  + (percnt) * runrec->runcal[jrun].prtemp;
776 
777  for (iw = 0; iw < 3; iw++) {
778  /* Convert temperature to radiance */
779  ttem = prtemp;
780  if (ttem < tembas[iw])
781  ttem = tembas[iw];
782  ai = (ttem - tembas[iw]) / teminc[iw];
783  if (ai < 0)
784  ai = 0;
785  else if (ai > NRADTAB)
786  ai = NRADTAB - 2;
787  ni = ai; /* whole part */
788  ai = ai - ni; /* fractional part */
789  radprt[iw] = (1.0 - ai) * teminv[iw][ni] + ai * teminv[iw][ni + 1];
790  }
791 
792  rspc[0] = 0.0;
793  switch (file->subsensorID) {
794  case NO07:
795  rspc[1] = -5.86;
796  rspc[2] = -4.95;
797  break;
798  case NO09:
799  rspc[1] = -5.53; /* Based on Table 2A; edited data, Apr 93 */
800  rspc[2] = -3.06;
801  break;
802  case NO10:
803  rspc[1] = -7.29; /* Based on Rao memo, Jan 4, 94 */
804  rspc[2] = -7.29;
805  break;
806  case NO11:
807  rspc[1] = -8.05; /* Based on Rao memo, Jan 4, 94 */
808  rspc[2] = -3.51;
809  break;
810  case NO12:
811  rspc[1] = -5.51; /* Based on Rao memo, Jan 4, 94 */
812  rspc[2] = -2.51;
813  break;
814  case NO14:
815  rspc[1] = -4.05; /* Based on S. Brown memo, Nov 22, 94 */
816  rspc[2] = -2.29;
817  break;
818  case NO15:
819  rspc[1] = -4.50; /* Based on Draft KLM, app D.1-14 */
820  rspc[2] = -3.61;
821  break;
822  case NO16:
823  rspc[1] = -2.467; /* Based on Draft KLM, app D.2-15 */
824  rspc[2] = -2.009;
825  break;
826  case NO17:
827  rspc[1] = -8.55; /* Based on Draft KLM, app D.3-2 */
828  rspc[2] = -3.97;
829  break;
830  case NO18:
831  rspc[1] = -5.53; /* Based on Draft KLM, app D.4-2 */
832  rspc[2] = -2.22;
833  break;
834  case NO19:
835  rspc[1] = -5.49; /* Based on Draft KLM, app D.6-2 */
836  rspc[2] = -3.39;
837  break;
838  default:
839  rspc[1] = -99.;
840  rspc[2] = -99.;
841  }
842 
843  /* noaa-16,17,18,19 have channel 3a and 3b */
844  /* 3a is 1.6 um during the day */
845  /* 3b is 3.75 um at night */
846  /* earlier avhrr satellites just have the 3.75 um (day and night) */
847  /* so process channel 3 as visible during the day, and IR at night */
848  for (iw = 0; iw < nbands + 1; iw++) {
849 
850  switch (iw) {
851  case 0:
852  data = ch1;
853  break;
854  case 1:
855  data = ch2;
856  break;
857  case 2:
858  data = ch3;
859  break;
860  }
861 
862  for (ip = 0; ip < npix; ip++) {
863 
864  ipb = ip * nbands + iw;
865  l1rec->Lt[ipb] = 0.0;
866 
867  /* check for sentinel values and flag as appropriate */
868  if (data[ip] >= 32000) {
869 
870  switch (data[ip]) {
871 
872  default:
873  l1rec->hilt[ip] = 1;
874  l1rec->Lt[ipb] = 1000.0;
875  break;
876  }
877 
878  } else if (l1rec->solz[ip] < SOLZNIGHTA) {
879  /* only do channels 1 and 2 for daytime pixels */
880  /* only do channel 3 for daytime for noaa-16,17,18,19 */
881  if (iw == 2 && file->subsensorID != NO16 && file->subsensorID != NO17 &&
882  file->subsensorID != NO18 && file->subsensorID != NO19) {
883  /* this channel doesn't make the pixel bad */
884  // l1rec->hilt[ip] = 1;
885  l1rec->Lt[ipb] = 1000.0;
886  } else {
887  status = avlooph_(&l1rec->senz[ip], &l1rec->solz[ip],
888  &l1rec->delphi[ip], rayly, aersol,
889  &l1rec->glint_coef[ip]);
890 
891  trad = (data[ip] * m[iw]) + b[iw];
892  trad -= rayly[iw];
893 
894  if (aersol[iw] != 0.0) {
895  trad /= aersol[iw];
896  }
897  /* aermlt in pathnlc is 1.0 so we don't need it here */
898 
899  l1rec->Lt[ipb] = trad;
900  }
901  }
902  }
903  }
904 
905  for (iw = 0; iw < nbandsir; iw++) {
906 
907  switch (iw) {
908  case 0:
909  data = ch3;
910  break;
911  case 1:
912  data = ch4;
913  break;
914  case 2:
915  data = ch5;
916  break;
917  }
918 
919  for (ip = 0; ip < npix; ip++) {
920 
921  ipb = ip * NBANDSIR + iw;
922  l1rec->Ltir[ipb] = 0.0;
923  l1rec->Bt[ipb] = BT_LO;
924 
925  /* don't do channel 3 for daytime for noaa-16,17,18,19 */
926  if (iw == 0 && l1rec->solz[ip] < SOLZNIGHTA
927  && (file->subsensorID == NO16
928  || file->subsensorID == NO17
929  || file->subsensorID == NO18
930  || file->subsensorID == NO19)) {
931  continue;
932  } else {
933 
934  /* check for sentinel values and flag as appropriate */
935  if (data[ip] >= 32000) {
936 
937  switch (data[ip]) {
938 
939  default:
940  l1rec->hilt[ip] = 1;
941  l1rec->Lt[ipb] = 1000.0;
942  break;
943  }
944 
945  } else {
946  // if (ip == pdbg1 && scan == ldbg1) {
947  // printf(" got to 320,40\n");
948  // }
949  /* Miami calibration */
950  if (avgbck[iw] != avgspc[iw]) {
951  if (l1_input->newavhrrcal == 1
952  && file->subsensorID == NO16
953  && iw > 0) {
954  /* Jon Mittaz's calibration for ch4 and ch5 */
955  /* T_inst = prtemp, R_ICT = radprt[iw] */
956  /* C_S = avgspc[iw], C_BB = avgbck[iw] */
957  /* C = data[ip] */
958  if (iw == 1) {
959  /* avhrr channel 4 Bt11 */
960  if (ip == pdbg1 && scan == ldbg1) {
961  printf(" ch4_coeffs=%f %f %f %f %f\n",
962  ch4_coeff[0], ch4_coeff[1],
963  ch4_coeff[2], ch4_coeff[3],
964  ch4_coeff[4]);
965  printf(" T_inst, prtemp = %f \n", prtemp);
966  printf(" R_ICT, radprt = %f \n",
967  radprt[iw]);
968  printf(" C_S, avgspc = %f \n", avgspc[iw]);
969  printf(" C_BB, avgbck = %f \n", avgbck[iw]);
970  printf(" C, data = %d \n", data[ip]);
971  }
972  trad = ch4_coeff[0] * (prtemp - 288.) + ch4_coeff[1]
973  + ((((ch4_coeff[2] * radprt[iw]) - (ch4_coeff[3]
974  * (pow(avgspc[iw] - avgbck[iw], 2.0))))
975  / (avgspc[iw] - avgbck[iw])) * (avgspc[iw] - data[ip]))
976  + ch4_coeff[4] * (pow(avgspc[iw] - data[ip], 2.0));
977  } else if (iw == 2) {
978  /* avhrr channel 5 Bt12 */
979  if (ip == pdbg1 && scan == ldbg1) {
980  printf(" ch5_coeffs=%f %f %f %f %f\n",
981  ch5_coeff[0], ch5_coeff[1],
982  ch5_coeff[2], ch5_coeff[3],
983  ch5_coeff[4]);
984  printf(" T_inst, prtemp = %f \n", prtemp);
985  printf(" R_ICT, radprt = %f \n",
986  radprt[iw]);
987  printf(" C_S, avgspc = %f \n", avgspc[iw]);
988  printf(" C_BB, avgbck = %f \n", avgbck[iw]);
989  printf(" C, data = %d \n", data[ip]);
990  }
991  trad = ch5_coeff[0] * (prtemp - 288.) + ch5_coeff[1]
992  + ((((ch5_coeff[2] * radprt[iw]) - (ch5_coeff[3]
993  * (pow(avgspc[iw] - avgbck[iw], 2.0))))
994  / (avgspc[iw] - avgbck[iw])) * (avgspc[iw] - data[ip]))
995  + ch5_coeff[4] * (pow(avgspc[iw] - data[ip], 2.0));
996  }
997 
998  } else {
999  trad = radprt[iw] * (avgspc[iw] - data[ip])
1000  / (avgspc[iw] - avgbck[iw])
1001  + rspc[iw] * (avgbck[iw] - data[ip])
1002  / (avgbck[iw] - avgspc[iw]);
1003  if (iw == 1) {
1004  /* apply non-linearity correction for channel 4 */
1005  switch (file->subsensorID) {
1006  case NO07:
1007  tnlc = 5.7843 - 1.0754e-1 * trad + 4.8042e-4 * pow(trad, 2);
1008  break;
1009  case NO09: /* Based on Table 2A; edited data, Apr 93 */
1010  tnlc = 5.24 + (.88643 - 1.) * trad + 6.033e-4 * pow(trad, 2);
1011  break;
1012  case NO10: /* Based on Rao memo, Jan 4, 94 */
1013  tnlc = 5.76 + (.88428 - 1.) * trad + 5.882e-4 * pow(trad, 2);
1014  break;
1015  case NO11: /* Based on Rao memo, Jan 4, 94 */
1016  tnlc = 7.21 + (.8412 - 1.) * trad + 8.739e-4 * pow(trad, 2);
1017  break;
1018  case NO12: /* Based on Rao memo, Jan 4, 94 */
1019  tnlc = 5.11 + (.88929 - 1.) * trad + 5.968e-4 * pow(trad, 2);
1020  break;
1021  case NO14: /* Based on S. Brown memo, Nov 22, 94 */
1022  tnlc = 3.72 + (.92378 - 1.) * trad
1023  + 3.822e-4 * pow(trad, 2);
1024  break;
1025  case NO15: /* Based on Draft KLM, app D.1-14 */
1026  tnlc = 4.76 + (-0.0932) * trad + 4.524e-4 * pow(trad, 2);
1027  break;
1028  case NO16: /* Based on Draft KLM, app D.2-15 */
1029  tnlc = 2.96 + (-0.05411) * trad + 2.4532e-4 * pow(trad, 2);
1030  break;
1031  case NO17: /* Based on Draft KLM, app D.3-2 */
1032  tnlc = 8.22 + (-0.15795) * trad + 7.5579e-4 * pow(trad, 2);
1033  break;
1034  case NO18: /* Based on Draft KLM, app D.4-2 */
1035  tnlc = 5.82 + (-0.11069) * trad + 5.2337e-4 * pow(trad, 2);
1036  break;
1037  case NO19: /* Based on Draft KLM, app D.6-2 */
1038  tnlc = 5.70 + (-0.11187) * trad + 5.4668e-4 * pow(trad, 2);
1039  break;
1040  default:
1041  tnlc = -99.;
1042  }
1043  trad += tnlc;
1044  } else if (iw == 2) {
1045  /* apply non-linearity correction for channel 5 */
1046  switch (file->subsensorID) {
1047  case NO07:
1048  tnlc = 4.4035 - 6.9038e-2 * trad + 2.5741e-4 * pow(trad, 2);
1049  break;
1050  case NO09: /* Based on Table 2A; edited data, Apr 93 */
1051  tnlc = 2.42 + (.95311 - 1.) * trad + 2.198e-4 * pow(trad, 2);
1052  break;
1053  case NO10: /* Based on Rao memo, Jan 4, 94 */
1054  tnlc = 5.76 + (.88428 - 1.) * trad + 5.882e-4 * pow(trad, 2);
1055  break;
1056  case NO11: /* Based on Rao memo, Jan 4, 94 */
1057  tnlc = 2.92 + (.94598 - 1.) * trad + 2.504e-4 * pow(trad, 2);
1058  break;
1059  case NO12: /* Based on Rao memo, Jan 4, 94 */
1060  tnlc = 1.91 + (.96299 - 1.) * trad + 1.775e-4 * pow(trad, 2);
1061  break;
1062  case NO14: /* Based on S. Brown memo, Nov 22, 94 */
1063  tnlc = 2.00 + (.96194 - 1.) * trad + 1.742e-4 * pow(trad, 2);
1064  break;
1065  case NO15: /* Based on Draft KLM, app D.1-14 */
1066  tnlc = 3.83 + (-0.0659) * trad + 2.811e-4 * pow(trad, 2);
1067  break;
1068  case NO16: /* Based on Draft KLM, app D.2-15 */
1069  tnlc = 2.25 + (-0.03665) * trad + 1.4854e-4 * pow(trad, 2);
1070  break;
1071  case NO17: /* Based on Draft KLM, app D.3-2 */
1072  tnlc = 4.31 + (-0.07318) * trad + 3.0976e-4 * pow(trad, 2);
1073  break;
1074  case NO18: /* Based on Draft KLM, app D.4-2 */
1075  tnlc = 2.67 + (-0.04360) * trad + 1.7715e-4 * pow(trad, 2);
1076  break;
1077  case NO19: /* Based on Draft KLM, app D.6-2 */
1078  tnlc = 3.58 + (-0.05991) * trad + 2.4985e-4 * pow(trad, 2);
1079  break;
1080  default:
1081  tnlc = -99.;
1082  }
1083  trad += tnlc;
1084  }
1085  }
1086 
1087  if (ip == pdbg1 && scan == ldbg1) {
1088  printf(" ch %d trad = %f\n", iw + 3, trad);
1089  printf(" newavhrrcal = %d\n",
1090  l1_input->newavhrrcal);
1091  printf(" xsatid = %d NO16=%d\n",
1092  file->subsensorID, NO16);
1093  }
1094  /* convert radiance to temperature */
1095  if (trad < 1.e-10)
1096  trad = 1.e-10;
1097  if (radinc[iw] != 0.0)
1098  ai = (log10(trad) - radbas[iw]) / radinc[iw];
1099  else
1100  ai = NRADTAB - 1;
1101  if (ai < 0.0)
1102  ai = 0.0;
1103  else if (ai >= NRADTAB - 1)
1104  ai = NRADTAB - 2;
1105  ni = ai; /* whole part */
1106  ai = ai - ni; /* fractional part */
1107  ttem = (1.0 - ai) * radinv[iw][ni] + ai * radinv[iw][ni + 1];
1108  if (ip == pdbg1 && scan == ldbg1) {
1109  printf(" ch %d ttem = %f\n", iw + 3, ttem);
1110  }
1111  if (l1_input->newavhrrcal
1112  == 1 && file->subsensorID == NO16) {
1113  /* convert Mittaz Bt's to old so we can test sst with old coeffs */
1114  // if (iw == 1) {
1115  // /* convert Mittaz Bt11 to old so we can test sst with old coeffs */
1116  // tcor = -10.306526 + 0.068776865*ttem - 0.00010926801*pow(ttem,2);
1117  // /* tcor is Walton-Mittaz so add Mittaz to get Walton */
1118  // ttem = ttem + tcor;
1119  // } else if (iw == 2) {}
1120  if (iw == 2) {
1121  /* use ch5 (Bt12) correction for ch4 (Bt11) also */
1122  /* convert Mittaz Bt12 to old so we can test sst with old coeffs */
1123  tcor = -6.9176245 + 0.038554339 * ttem - 4.7661371e-05 * pow(ttem, 2);
1124  /* tcor is Walton-Mittaz so add Mittaz to get Walton */
1125  ttem = ttem + tcor;
1126  /* add ch5 correction to ch4 also */
1127  /* ipb is ip*NBANDSIR+iw, so ch4 should be ipb-1 */
1128  l1rec->Bt[ipb - 1] = l1rec->Bt[ipb - 1] + tcor;
1129  if (ip == pdbg1 && scan == ldbg1) {
1130  printf(" real ch = %d, ttem = %f, tcor = %f\n", iw + 3 - 1, ttem, tcor);
1131  }
1132  }
1133  if (ip == pdbg1 && scan == ldbg1) {
1134  printf(" ch = %d, ttem = %f, tcor = %f\n", iw + 3, ttem, tcor);
1135  }
1136  }
1137  ttem = ttem - TEMP_REF; /* convert from Deg K to Deg C */
1138  l1rec->Bt[ipb] = ttem;
1139  }
1140  /* put top of the atmosphere counts in Ltir */
1141  l1rec->Ltir[ipb] = data[ip];
1142  // l1rec->prtemp[ip] = prtemp - TEMP_REF;
1143 
1144  }
1145  }
1146  }
1147  }
1148 
1149  l1rec->npix = file->npix;
1150  l1rec->detnum = (int32_t) detnum;
1151  l1rec->mside = 0;
1152 
1153  /* Convert IR bands to brightness temperature */
1154  // radiance2bt(l1rec,-1);
1155  return (LIFE_IS_GOOD);
1156 }
1157 
1158 int closel1_aci_hdf(filehandle *file) {
1159  if ((sd_id_g != sd_id) && SDend(sd_id_g)) {
1160  fprintf(stderr, "-E- %s line %d: SDend(%d) failed for file, %s.\n",
1161  __FILE__, __LINE__, sd_id, file->geofile);
1162  return (HDF_FUNCTION_ERROR);
1163  }
1164  if (SDend(sd_id)) {
1165  fprintf(stderr, "-E- %s line %d: SDend(%d) failed for file, %s.\n",
1166  __FILE__, __LINE__, sd_id, file->name);
1167  return (HDF_FUNCTION_ERROR);
1168  }
1169 
1170  return (LIFE_IS_GOOD);
1171 }
1172 
1173 const char* xsatid2name(int xsatid) {
1174  switch (xsatid) {
1175  case NO06:
1176  return "NO06";
1177  case NO07:
1178  return "NO07";
1179  case NO08:
1180  return "NO08";
1181  case NO09:
1182  return "NO09";
1183  case NO10:
1184  return "NO10";
1185  case NO11:
1186  return "NO11";
1187  case NO12:
1188  return "NO12";
1189  case NO14:
1190  return "NO14";
1191  case NO15:
1192  return "NO15";
1193  case NO16:
1194  return "NO16";
1195  case NO17:
1196  return "NO17";
1197  case NO18:
1198  return "NO18";
1199  case NO19:
1200  return "NO19";
1201  default:
1202  return "UNKNOWN";
1203  }
1204 }
1205 
1206 int satname2xsatid(const char* satname) {
1207  if (strcmp(satname, "NOA6") == 0)
1208  return NO06;
1209  else if (strcmp(satname, "NOA7") == 0)
1210  return NO07;
1211  else if (strcmp(satname, "NOA8") == 0)
1212  return NO08;
1213  else if (strcmp(satname, "NOA9") == 0)
1214  return NO09;
1215  else if (strcmp(satname, "NO10") == 0)
1216  return NO10;
1217  else if (strcmp(satname, "NO11") == 0)
1218  return NO11;
1219  else if (strcmp(satname, "NO12") == 0)
1220  return NO12;
1221  else if (strcmp(satname, "NO14") == 0)
1222  return NO14;
1223  else if (strcmp(satname, "NO15") == 0)
1224  return NO15;
1225  else if (strcmp(satname, "NO16") == 0)
1226  return NO16;
1227  else if (strcmp(satname, "NO17") == 0)
1228  return NO17;
1229  else if (strcmp(satname, "NO18") == 0)
1230  return NO18;
1231  else if (strcmp(satname, "NO19") == 0)
1232  return NO19;
1233  else
1234  return -999;
1235 }
1236 
void etloadresp_(int32_t *lin, char *cal)
#define NO17
Definition: l1_aci_hdf.h:15
#define TEMP_REF
Definition: l1_aci_hdf.c:22
Definition: rawcal.h:84
int rdSDS(int32_t fileID, const char sdsname[], int32_t start1, int32_t start2, int32_t edges1, int32_t edges2, void *array_data)
int32_t day
int status
Definition: l1_czcs_hdf.c:32
double yds2unix(int16_t year, int16_t day, double secs)
Definition: yds2unix.c:7
#define NO11
Definition: l1_aci_hdf.h:10
#define NBANDSIR
Definition: filehandle.h:23
#define FAIL
Definition: ObpgReadGrid.h:18
#define NULL
Definition: decode_rs.h:63
MOD_PR01 Production producing one five minute granule of output data in each run It can be configured to produce as many as three five minute granules per run Each execution with one construction record and one date file for each dataset In normal these are created by which splits them out of the hour datasets For LANCE they are created by which merges all session MODIS L0 datasets overlapping the requested time and extracts from the merged data those packets which fall within that time period Each scan of data is stored in the L1A granule that covers the start time of that scan
Definition: MOD_PR01_pr.txt:19
read l1rec
int readl1_aci_hdf(filehandle *file, int32_t scan, l1str *l1rec)
Definition: l1_aci_hdf.c:129
#define NO08
Definition: l1_aci_hdf.h:7
#define NAVBANDS
Definition: l1_aci_hdf.c:18
float etintegrate_(int32_t *ich, float *etemp)
int32_t jday(int16_t i, int16_t j, int16_t k)
Definition: jday.c:4
int32 nscan
Definition: l1_czcs_hdf.c:19
#define READ_SDS_ID(sd_id, nam, ptr, s0, s1, s2, s3, e0, e1, e2, e3)
Definition: hdf4utils.h:109
no change in intended resolving MODur00064 Corrected handling of bad ephemeris attitude resolving resolving GSFcd00179 Corrected handling of fill values for[Sensor|Solar][Zenith|Azimuth] resolving MODxl01751 Changed to validate LUT version against a value retrieved from the resolving MODxl02056 Changed to calculate Solar Diffuser angles without adjustment for estimated post launch changes in the MODIS orientation relative to incidentally resolving defects MODxl01766 Also resolves MODxl01947 Changed to ignore fill values in SCI_ABNORM and SCI_STATE rather than treating them as resolving MODxl01780 Changed to use spacecraft ancillary data to recognise when the mirror encoder data is being set by side A or side B and to change calculations accordingly This removes the need for seperate LUTs for Side A and Side B data it makes the new LUTs incompatible with older versions of the and vice versa Also resolves MODxl01685 A more robust GRing algorithm is being which will create a non default GRing anytime there s even a single geolocated pixel in a granule Removed obsolete messages from seed file
Definition: HISTORY.txt:413
character(len=1000) if
Definition: names.f90:13
#define LIFE_IS_GOOD
Definition: passthebuck.h:4
#define NO06
Definition: l1_aci_hdf.h:5
#define NO10
Definition: l1_aci_hdf.h:9
int getDims(int32_t fileID, const char sdsname[], int32_t dims[])
#define NO12
Definition: l1_aci_hdf.h:11
#define NO09
Definition: l1_aci_hdf.h:8
Definition: runcal.h:81
l1_input_t * l1_input
Definition: l1_options.c:9
#define NO16
Definition: l1_aci_hdf.h:14
int getHDFattr(int32_t fileID, const char attrname[], const char sdsname[], void *data)
#define NO15
Definition: l1_aci_hdf.h:13
int32_t avlooph_(float *satz, float *solz, float *delphi, float *rayly, float *aersol, float *aglint)
const char * xsatid2name(int xsatid)
Definition: l1_aci_hdf.c:1173
no change in intended resolving MODur00064 Corrected handling of bad ephemeris attitude data
Definition: HISTORY.txt:356
float etinvert_(int32_t *ich, float *radi)
#define NO14
Definition: l1_aci_hdf.h:12
#define NO19
Definition: l1_aci_hdf.h:17
#define SDNAME
Definition: l1_aci_hdf.c:26
data_t b[NROOTS+1]
Definition: decode_rs.h:77
int closel1_aci_hdf(filehandle *file)
Definition: l1_aci_hdf.c:1158
int32_t nbands
int satname2xsatid(const char *satname)
Definition: l1_aci_hdf.c:1206
Extra metadata that will be written to the HDF4 file l2prod rank
#define BT_LO
Definition: l1.h:54
int32_t ch4
Definition: atrem_corl1.h:118
#define NRADTAB
Definition: l1_aci_hdf.c:24
#define SOLZNIGHTA
Definition: l1.h:58
#define NO07
Definition: l1_aci_hdf.h:6
logical function leap(YEAR)
Definition: leap.f:10
#define abs(a)
Definition: misc.h:90
How many dimensions is the output array Default is Not sure if anything above will work correctly strcpy(l2prod->title, "no title yet")
int npix
Definition: get_cmp.c:27
#define HDF_FUNCTION_ERROR
Definition: passthebuck.h:7
int openl1_aci_hdf(filehandle *file)
Definition: l1_aci_hdf.c:50
int32_t avconsh_(int32_t *lunin, int32_t *npix, int32_t *jday)
#define NO18
Definition: l1_aci_hdf.h:16