OB.DAAC Logo
NASA Logo
Ocean Color Science Software

ocssw V2022
l1_hico_h5.c
Go to the documentation of this file.
1 /*
2  * l1_hico_hdf.c
3  *
4  */
5 
6 #include "l1_hico_h5.h"
7 
8 #include <string.h>
9 #include "h5io.h"
10 #include "l1.h"
11 
12 #define NUM_GEO_DATA 6
13 
14 static char *geo_name[] = {"latitude", "longitude", "sensor_azimuth",
15  "sensor_zenith", "solar_azimuth", "solar_zenith"};
16 
17 // HICO file private information stored in filehandle
18 
19 typedef struct hico_private_struct {
20  int numBands; // number of bands in the hico file
21  int* wave_ix; // array that maps the sensor bands into the hico file
22  int32_t syear;
23  int32_t sday;
24  int32_t smsec;
25  float lt_slope;
26  float lt_intercept;
27  uint16_t *rad_data;
28  h5io_str* fileID;
29  h5io_str* ds_id;
31  int orientation; // is HICO orientation flipped
32  float *fwhm;
33 } hico_private_t;
34 
35 static hico_private_t* allocatePrivateData() {
36  int i;
37  hico_private_t* pData = (hico_private_t*) calloc(1, sizeof (hico_private_t));
38  pData->lt_slope = 1.0;
39  pData->fileID = (h5io_str*) calloc(1, sizeof (h5io_str));
40  pData->ds_id = (h5io_str*) calloc(1, sizeof (h5io_str));
41  for (i = 0; i < NUM_GEO_DATA; i++)
42  pData->geo_dat_id[i] = (h5io_str*) calloc(1, sizeof (h5io_str));
43  pData->orientation = 0;
44  return pData;
45 }
46 
47 static void freePrivateData(hico_private_t* pData) {
48  int i;
49 
50  if (pData->rad_data)
51  free(pData->rad_data);
52  if (pData->wave_ix)
53  free(pData->wave_ix);
54  free(pData->fileID);
55  free(pData->ds_id);
56  for (i = 0; i < NUM_GEO_DATA; i++)
57  free(pData->geo_dat_id[i]);
58  free(pData);
59 }
60 
61 int openl1_hico_h5(filehandle * file) {
62  h5io_str gid;
63  int year, month, day, hour, minute;
64  int ids;
65  int ndim;
66  int dims[3];
67  int sto_len;
68  int sec;
69  char sdate[16];
70  char edate[16];
71  char stime[16];
72  char etime[16];
73  char g_path[100];
74  H5T_class_t f_class;
75  hid_t f_native_typ;
76  float* wvls;
77  hico_private_t* pData;
78  int32_t *Lambda;
79  int i;
80  char orientationStr[10];
81 
82  pData = allocatePrivateData();
83  file->private_data = pData;
84 
85  // get rid of the silly HDF5 error printing
86  H5Eset_auto(H5E_DEFAULT, NULL, NULL);
87 
88  if (h5io_openr(file->name, 0, pData->fileID)) {
89  fprintf(stderr, "-E- %s Line %d: Failure to open %s\n",
90  __FILE__, __LINE__, file->name);
91  freePrivateData(pData);
92  return 1;
93  }
94 
95  // set Lt dataset pointer
96  if (h5io_set_grp(pData->fileID, "products", &gid) != 0) {
97  fprintf(stderr, "-E- %s Line %d: Unable to open group products\n",
98  __FILE__, __LINE__);
99  h5io_close(pData->fileID);
100  freePrivateData(pData);
101  return 1;
102  }
103 
104  if (h5io_set_ds(&gid, "Lt", pData->ds_id) != 0) {
105  fprintf(stderr, "-E- %s Line %d: Unable to open data set: Lt\n",
106  __FILE__, __LINE__);
107  h5io_close(&gid);
108  h5io_close(pData->fileID);
109  freePrivateData(pData);
110  return 1;
111  }
112  h5io_close(&gid);
113 
114  // get number of wavelengths in file
115  if (h5io_info(pData->ds_id, NULL, &f_class, &f_native_typ, &ndim, dims, &sto_len)
116  != 0) {
117  fprintf(stderr, "-E- %s Line %d: Unable to get data set dimensions\n",
118  __FILE__, __LINE__);
119  h5io_close(pData->ds_id);
120  h5io_close(pData->fileID);
121  freePrivateData(pData);
122  return 1;
123  }
124  pData->numBands = dims[2];
125  // allocate space for wavelength array
126  pData->fwhm = (float*) malloc(sizeof (float)*pData->numBands);
127  if (pData->fwhm == NULL) {
128  fprintf(stderr,
129  "-E- %s Line %d: Unable to allocate data for wavelength bandwidth\n",
130  __FILE__, __LINE__);
131  h5io_close(pData->ds_id);
132  h5io_close(pData->fileID);
133  freePrivateData(pData);
134  return 1;
135  }
136  wvls = (float*) malloc(sizeof (float) * pData->numBands);
137  if (wvls == NULL) {
138  fprintf(stderr,
139  "-E- %s Line %d: Unable to allocate data for wavelengths\n",
140  __FILE__, __LINE__);
141  h5io_close(pData->ds_id);
142  h5io_close(pData->fileID);
143  freePrivateData(pData);
144  return 1;
145  }
146 
147  // Get the wavelength bandwidth list
148  if (h5io_rd_attr(pData->ds_id, "fwhm", (void *) pData->fwhm) != 0) {
149  fprintf(stderr,
150  "-E- %s Line %d: Unable to read the wavelength bandwidth (fwhm) attribute\n",
151  __FILE__, __LINE__);
152  free(pData->fwhm);
153  h5io_close(pData->ds_id);
154  h5io_close(pData->fileID);
155  freePrivateData(pData);
156  return 1;
157  }
158  // Get the wavelength list
159  if (h5io_rd_attr(pData->ds_id, "wavelengths", (void *) wvls) != 0) {
160  fprintf(stderr,
161  "-E- %s Line %d: Unable to read the wavelength attribute\n",
162  __FILE__, __LINE__);
163  free(wvls);
164  h5io_close(pData->ds_id);
165  h5io_close(pData->fileID);
166  freePrivateData(pData);
167  return 1;
168  }
169 
170  // allocate space for wavelength index array
171  pData->wave_ix = (int*) malloc(sizeof (int) * file->nbands);
172  if (pData->wave_ix == NULL) {
173  fprintf(stderr,
174  "-E- %s Line %d: Unable to allocate data for wavelength index array\n",
175  __FILE__, __LINE__);
176  free(wvls);
177  h5io_close(pData->ds_id);
178  h5io_close(pData->fileID);
179  freePrivateData(pData);
180  return 1;
181  }
182 
183  // fill up the wavelength index array
184  rdsensorinfo(file->sensorID, l1_input->evalmask, "Lambda", (void **) &Lambda);
185  for (i = 0; i < file->nbands; i++) {
186  pData->wave_ix[i] = windex(Lambda[i], wvls, pData->numBands);
187  }
188  free(wvls);
189 
190  // Get the start and end date,time
191  if (h5io_set_grp(pData->fileID,
192  "metadata/FGDC/Identification_Information/Time_Period_of_Content",
193  &gid) != 0) {
194  fprintf(stderr,
195  "-E- %s Line %d: Unable to open group Time_Period_of_Content\n",
196  __FILE__, __LINE__);
197  h5io_close(pData->ds_id);
198  h5io_close(pData->fileID);
199  freePrivateData(pData);
200  return 1;
201  }
202 
203  bzero(&sdate, sizeof (sdate));
204  if (h5io_rd_attr(&gid, "Beginning_Date", (void *) &sdate) != 0) {
205  fprintf(stderr,
206  "-E- %s Line %d: Unable to read the Beginning_Date attribute:\n",
207  __FILE__, __LINE__);
208  h5io_close(&gid);
209  h5io_close(pData->ds_id);
210  h5io_close(pData->fileID);
211  freePrivateData(pData);
212  return 1;
213  }
214 
215  bzero(&stime, sizeof (stime));
216  if (h5io_rd_attr(&gid, "Beginning_Time", (void *) &stime) != 0) {
217  fprintf(stderr,
218  "-E- %s Line %d: Unable to read the Beginning_Date attribute:\n",
219  __FILE__, __LINE__);
220  h5io_close(&gid);
221  h5io_close(pData->ds_id);
222  h5io_close(pData->fileID);
223  freePrivateData(pData);
224  return 1;
225  }
226 
227  bzero(&edate, sizeof (edate));
228  if (h5io_rd_attr(&gid, "Ending_Date", (void *) &edate) != 0) {
229  fprintf(stderr,
230  "-E- %s Line %d: Unable to read the Beginning_Date attribute:\n",
231  __FILE__, __LINE__);
232  h5io_close(&gid);
233  h5io_close(pData->ds_id);
234  h5io_close(pData->fileID);
235  freePrivateData(pData);
236  return 1;
237  }
238 
239  bzero(&etime, sizeof (etime));
240  if (h5io_rd_attr(&gid, "Ending_Time", (void *) &etime) != 0) {
241  fprintf(stderr,
242  "-E- %s Line %d: Unable to read the Beginning_Date attribute:\n",
243  __FILE__, __LINE__);
244  h5io_close(&gid);
245  h5io_close(pData->ds_id);
246  h5io_close(pData->fileID);
247  freePrivateData(pData);
248  return 1;
249  }
250  h5io_close(&gid);
251 
252  // Parse the date/time strings
253  sscanf(sdate, "%04d%02d%02d", &year, &month, &day);
254  sscanf(stime, "%02d%02d%02d", &hour, &minute, &sec);
255  ymdhms2ydmsec(year, month, day, hour, minute, sec, &pData->syear, &pData->sday, &pData->smsec);
256 
257  if (h5io_rd_attr(pData->ds_id, "scale_factor", (void *) &pData->lt_slope) != 0) {
258  if (h5io_rd_attr(pData->ds_id, "slope", (void *) &pData->lt_slope) != 0) {
259  fprintf(stderr,
260  "-E- %s Line %d: Unable to read the Lt slope attribute:\n",
261  __FILE__, __LINE__);
262  h5io_close(pData->ds_id);
263  h5io_close(pData->fileID);
264  freePrivateData(pData);
265  return 1;
266  }
267  }
268 
269  if (h5io_rd_attr(pData->ds_id, "add_offset", (void *) &pData->lt_intercept) != 0) {
270  if (h5io_rd_attr(pData->ds_id, "intercept", (void *) &pData->lt_intercept) != 0) {
271  fprintf(stderr,
272  "-E- %s Line %d: Unable to read the Lt intercept attribute:\n",
273  __FILE__, __LINE__);
274  h5io_close(pData->ds_id);
275  h5io_close(pData->fileID);
276  freePrivateData(pData);
277  return 1;
278  }
279  }
280 
281  // Set to the geolocation datasets for the lat, lon, and view angles
282  // TODO: Get real with the navigation code ;)
283 
284  for (ids = 0; ids < NUM_GEO_DATA; ids++) {
285  sprintf(g_path, "navigation/%s", geo_name[ids]);
286  if (h5io_set_ds(pData->fileID, g_path, pData->geo_dat_id[ids]) != 0) {
287 
288  // look for latitude/longitude and latitudes/longitudes
289  if(ids == 0 || ids == 1) {
290  sprintf(g_path, "navigation/%ss", geo_name[ids]);
291  if (h5io_set_ds(pData->fileID, g_path, pData->geo_dat_id[ids]) != 0) {
292  fprintf(stderr,
293  "-E- %s Line %d: Unable to set ds # %d in geolocation file\n",
294  __FILE__, __LINE__, ids);
295  h5io_close(pData->ds_id);
296  h5io_close(pData->fileID);
297  freePrivateData(pData);
298  return 1;
299  }
300  }
301  }
302  }
303 
304  // Get the HICO orientation
305  if (h5io_set_grp(pData->fileID, "metadata/HICO/Calibration", &gid) != 0) {
306  fprintf(stderr,
307  "-E- %s Line %d: Unable to open group metadata/HICO/Calibration\n",
308  __FILE__, __LINE__);
309  h5io_close(pData->ds_id);
310  h5io_close(pData->fileID);
311  freePrivateData(pData);
312  return 1;
313  }
314 
315  // #### note that there was a space in the attribute name
316  // #### how annoying.
317  bzero(&orientationStr, sizeof (orientationStr));
318  if (h5io_rd_attr(&gid, "hico_orientation_from_quaternion", (void *) &orientationStr) != 0) {
319  if (h5io_rd_attr(&gid, "hico_orientation_from_quaternion ", (void *) &orientationStr) != 0) {
320  fprintf(stderr,
321  "-E- %s Line %d: Unable to read the hico_orientation_from_quaternion attribute:\n",
322  __FILE__, __LINE__);
323  h5io_close(&gid);
324  h5io_close(pData->ds_id);
325  h5io_close(pData->fileID);
326  freePrivateData(pData);
327  return 1;
328  }
329  }
330  h5io_close(&gid);
331 
332  if (strstr(orientationStr, "-XVV")) {
333  pData->orientation = 1;
334  } else {
335  pData->orientation = 0;
336  }
337 
338  file->npix = dims[1];
339  file->nscan = dims[0];
340  file->sensorID = HICO;
341  strcpy(file->spatialResolution, "100 m");
342 
343  pData->rad_data = (uint16_t *) calloc(file->npix, sizeof (uint16_t));
344 
345  return (LIFE_IS_GOOD);
346 }
347 
348 int readl1_hico_h5(filehandle *file, int32_t scan, l1str *l1rec, int lonlat) {
349  int igeo;
350  int start[3], count[3];
351  int32_t ib, ip, ipb;
352  float *iptr;
353  int32_t dmsec = 13;
354  hico_private_t* pData = (hico_private_t*) file->private_data;
355 
356  l1rec->scantime = yds2unix(pData->syear, pData->sday, (double) ((pData->smsec + scan * dmsec) / 1.e3));
357 
358  // get the location and view information
359  for (igeo = 0; igeo < NUM_GEO_DATA; igeo++) {
360 
361  if (lonlat && igeo > 1)
362  return (LIFE_IS_GOOD);
363 
364  switch (igeo) {
365  case 0:
366  iptr = l1rec->lat;
367  break;
368  case 1:
369  iptr = l1rec->lon;
370  break;
371  case 2:
372  iptr = l1rec->sena;
373  break;
374  case 3:
375  iptr = l1rec->senz;
376  break;
377  case 4:
378  iptr = l1rec->sola;
379  break;
380  case 5:
381  iptr = l1rec->solz;
382  break;
383  }
384  if (pData->orientation) {
385  start[0] = file->nscan - scan - 1;
386  } else {
387  start[0] = scan;
388  }
389  start[1] = 0;
390  count[0] = 1;
391  count[1] = file->npix;
392  if (h5io_rd_ds_slice(pData->geo_dat_id[igeo], start, count, (void *) iptr)
393  != 0) {
394  fprintf(stderr,
395  "-E- %s Line %d: Failed to geonav scan %d of parm %d\n",
396  __FILE__, __LINE__, scan, igeo);
397  return 1;
398  }
399  }
400 
401  // read in radiance data
402  if (pData->orientation) {
403  start[0] = file->nscan - scan - 1;
404  } else {
405  start[0] = scan;
406  }
407  start[1] = 0;
408  count[0] = 1;
409  count[1] = file->npix;
410  count[2] = 1;
411 
412  for (ib = 0; ib < l1rec->l1file->nbands; ib++) {
413  start[2] = pData->wave_ix[ib];
414  if (h5io_rd_ds_slice(pData->ds_id, start, count, (void *) pData->rad_data) != 0) {
415  fprintf(stderr,
416  "-E- %s Line %d: Failed to read scan %d of band %d\n",
417  __FILE__, __LINE__, scan, ib);
418  return 1;
419  }
420 
421  for (ip = 0; ip < file->npix; ip++) {
422  ipb = ip * l1rec->l1file->nbands + ib;
423 
424  if (pData->rad_data[ip] > 0) {
425  l1rec->Lt[ipb] = ((float) (pData->rad_data[ip]) * pData->lt_slope
426  + pData->lt_intercept) / 10.0;
427  } else {
428  l1rec->Lt[ipb] = BAD_FLT;
429  // l1rec->navfail[ip] = 1;
430  }
431  }
432  }
433 
434  return (LIFE_IS_GOOD);
435 }
436 
437 int closel1_hico_h5(filehandle *file) {
438  int i;
439  hico_private_t* pData = (hico_private_t*) file->private_data;
440 
441  for (i = 0; i < NUM_GEO_DATA; i++) {
442  if (h5io_close(pData->geo_dat_id[i]) != 0) {
443  fprintf(stderr,
444  "-E- %s Line %d: Failed to close navigation data %s\n",
445  __FILE__, __LINE__, geo_name[i]);
446  return 1;
447  }
448  }
449 
450  if (h5io_close(pData->ds_id) != 0) {
451  fprintf(stderr, "-E- %s Line %d: Failed to close Lt data set\n",
452  __FILE__, __LINE__);
453  return 1;
454  }
455 
456  if (h5io_close(pData->fileID) != 0) {
457  fprintf(stderr, "-E- %s Line %d: Failed to close HDF5 file %s\n",
458  __FILE__, __LINE__, file->name);
459  return 1;
460  }
461 
462  freePrivateData(pData);
463  file->private_data = NULL;
464 
465  return (LIFE_IS_GOOD);
466 }
h5io_str * ds_id
Definition: l1_hico_h5.c:29
int32_t day
double yds2unix(int16_t year, int16_t day, double secs)
Definition: yds2unix.c:7
int h5io_openr(char *file, int opt, h5io_str *id)
Definition: h5io.c:4
#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 h5io_info(h5io_str *id, char *attr_name, H5T_class_t *class, hid_t *native_typ, int *ndim, int *dim_siz, int *sto_len)
Definition: h5io.c:173
int h5io_close(h5io_str *id)
Definition: h5io.c:115
int readl1_hico_h5(filehandle *file, int32_t scan, l1str *l1rec, int lonlat)
Definition: l1_hico_h5.c:348
void bzero()
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
int h5io_set_grp(h5io_str *id, char *path_name, h5io_str *grp_id)
Definition: h5io.c:369
#define LIFE_IS_GOOD
Definition: passthebuck.h:4
int h5io_set_ds(h5io_str *id, char *path_name, h5io_str *ds_id)
Definition: h5io.c:324
subroutine lonlat(alon, alat, xlon, ylat)
Definition: lonlat.f:2
int h5io_rd_ds_slice(h5io_str *ds_id, int *start, int *count, void *data)
Definition: h5io.c:602
int h5io_rd_attr(h5io_str *id, char *attr_name, void *data)
Definition: h5io.c:412
void freePrivateData(oli_t *data)
Definition: l1_oli.c:92
l1_input_t * l1_input
Definition: l1_options.c:9
int openl1_hico_h5(filehandle *file)
Definition: l1_hico_h5.c:61
#define NUM_GEO_DATA
Definition: l1_hico_h5.c:12
h5io_str * fileID
Definition: l1_hico_h5.c:28
#define BAD_FLT
Definition: jplaeriallib.h:19
void ymdhms2ydmsec(int yy, int mm, int dd, int hh, int mn, int sc, int32_t *year, int32_t *day, int32_t *msec)
Definition: ydhms2ydmsec.c:3
uint16_t * rad_data
Definition: l1_hico_h5.c:27
h5io_str * geo_dat_id[NUM_GEO_DATA]
Definition: l1_hico_h5.c:30
int windex(float wave, float twave[], int ntwave)
Definition: windex.c:73
int32_t rdsensorinfo(int32_t, int32_t, const char *, void **)
Definition: rdsensorinfo.c:69
for(i=0;i< NROOTS;i++) s[i]
Definition: decode_rs.h:85
#define HICO
Definition: sensorDefs.h:25
int closel1_hico_h5(filehandle *file)
Definition: l1_hico_h5.c:437
int i
Definition: decode_rs.h:71
How many dimensions is the output array Default is Not sure if anything above will work correctly strcpy(l2prod->title, "no title yet")
int count
Definition: decode_rs.h:79