OB.DAAC Logo
NASA Logo
Ocean Color Science Software

ocssw V2022
l1b_oci.c
Go to the documentation of this file.
1 /* ============================================================================ */
2 /* module l1b_oci.c - functions to read OCI L1B for MSL12 */
3 /* Written By: Steve Lockhart SAIC August 2018 */
4 /* -Started with l1a_hawkeye.c and made changes to support OCI. */
5 /* */
6 /* ============================================================================ */
7 
8 /* Issues:
9  -Assume num_pixels = ccd_pixels = SWIR_pixels
10  -Temporarily set tot_num_bands = 62 (since temporarily using copy of ocia share stuff)
11  -No scan_delta_time_ms, so use secs of day (with wrap-around?)
12  -Get Lt
13  -Skip over first two red bands (600, 605), as they overlap with last two blue bands
14  i.e. start at the third red band. So, this affects the total number of bands.
15  -Combine the high-res and low-res SWIR bands?
16  -Get sensor_band_parameters
17  -Review Ltir: what I carried over from VIIRS reader may not be appropriate for OCI.
18  -Make sure default (FILL) values are correct wrt what is in the input file.
19 */
20 
21 // Includes from l1_viirs_nc.c
22 #include <nc4utils.h>
23 #include "l1.h"
24 #include "libnav.h"
25 #include <productInfo.h>
26 #include <allocate2d.h>
27 
28 // Includes from demo_cal.c
29 #include <stdio.h>
30 #include <sys/types.h>
31 #include <sys/stat.h>
32 #include <stdlib.h>
33 #include <fcntl.h>
34 #include <string.h>
35 #include <errno.h>
36 #include <unistd.h>
37 #include <libgen.h>
38 #include <math.h>
39 #include <timeutils.h>
40 #include <gsl/gsl_math.h>
41 #include <gsl/gsl_spline.h>
42 #include <gsl/gsl_sort.h>
43 #include <gsl/gsl_statistics.h>
44 #include <stdbool.h>
45 
46 // New includes
47 #include "l1b_oci.h"
48 
49 
50 // MBAND_NUM_DETECTORS changed to 1
51 #define MBAND_NUM_DETECTORS 1
52 
53 
54 // Declarations from l1_viirs_nc.c
55 static int32_t prevScan = -1;
56 
57 //static int geoFileId;
58 static int geoGeolocationGrp;
59 static int l1bScanLineGrp;
60 static int l1bObservationGrp;
61 static int lonId, latId, senzId, senaId, solzId, solaId, HAMSideId, scanQualityId;
62 
63 static float *Fobar; // reflectance to radiance conversion factors
64 static int extract_pixel_start = 0;
65 
66 static short *tmpShort;
67 static unsigned char *tmpByte;
68 static int nline;
69 static size_t num_scans, num_pixels;
70 static size_t num_blue_bands, num_red_bands, num_SWIR_bands, tot_num_bands;
71 
72 static int firstCall = 1;
73 static double starttime;
74 static double lastvalidtime;
75 static int lastvalidscan = 0;
76 static double time_interval;
77 
78 static float latGeoFillValue = -999.9;
79 static float lonGeoFillValue = -999.9;
80 static short senzGeoFillValue = -32768;
81 static short senaGeoFillValue = -32768;
82 static short solzGeoFillValue = -32768;
83 static short solaGeoFillValue = -32768;
84 
85 static float latL2FillValue = -999.0;
86 static float lonL2FillValue = -999.0;
87 static float senzL2FillValue = -32767;
88 static float senaL2FillValue = -32767;
89 static float solzL2FillValue = -32767;
90 static float solaL2FillValue = -32767;
91 
92 // Declarations added for OCI L1B
93 static int Lt_blue_varid, Lt_red_varid, Lt_SWIR_varid;
94 static double **Lt_blue; //[num_blue_bands][num_pixels], This scan
95 static double **Lt_red; //[num_red_bands][num_pixels], This scan
96 static double **Lt_SWIR; //[num_SWIR_bands][num_pixels], This scan
97 static double **Lt; //[tot_num_bands][num_pixels], This scan
98 
99 
113 int openl1b_oci(filehandle * file) {
114  char *fltime;
115 
116  int ncid_L1B, dimid, status;
117  size_t att_len;
118  int orbit_number;
119 
120  // Open the netcdf4 input file
121  printf("Opening oci l1b file\n");
122  status = nc_open(file->name, NC_NOWRITE, &ncid_L1B);
123  if (status != NC_NOERR) {
124  fprintf(stderr, "-E- %s line %d: nc_open(%s) failed.\n",
125  __FILE__, __LINE__, file->name);
126  return (1);
127  }
128 
129  // Get dims from L1B file: num_scans, num_bands, num_pixels
130  // num_scans
131  status = nc_inq_dimid(ncid_L1B, "number_of_scans", &dimid);
132  if (status != NC_NOERR) {
133  fprintf(stderr, "-E- Error reading num_scans.\n");
134  exit(EXIT_FAILURE);
135  };
136  nc_inq_dimlen(ncid_L1B, dimid, &num_scans);
137  // num_blue_bands
138  status = nc_inq_dimid(ncid_L1B, "blue_bands", &dimid);
139  if (status != NC_NOERR) {
140  fprintf(stderr, "-E- Error reading num_blue_bands.\n");
141  exit(EXIT_FAILURE);
142  };
143  nc_inq_dimlen(ncid_L1B, dimid, &num_blue_bands);
144  // num_red_bands
145  status = nc_inq_dimid(ncid_L1B, "red_bands", &dimid);
146  if (status != NC_NOERR) {
147  fprintf(stderr, "-E- Error reading num_red_bands.\n");
148  exit(EXIT_FAILURE);
149  };
150  nc_inq_dimlen(ncid_L1B, dimid, &num_red_bands);
151  // num_SWIR_bands
152  status = nc_inq_dimid(ncid_L1B, "SWIR_bands", &dimid);
153  if (status != NC_NOERR) {
154  fprintf(stderr, "-E- Error reading num_SWIR_bands.\n");
155  exit(EXIT_FAILURE);
156  };
157  nc_inq_dimlen(ncid_L1B, dimid, &num_SWIR_bands);
158  // num_pixels
159  status = nc_inq_dimid(ncid_L1B, "ccd_pixels", &dimid);
160  if (status != NC_NOERR) {
161  fprintf(stderr, "-E- Error reading num_pixels.\n");
162  exit(EXIT_FAILURE);
163  };
164  nc_inq_dimlen(ncid_L1B, dimid, &num_pixels);
165 
166  // Derive nline (See issue above.) This line is from l1_viirs_nc.c.
167  nline = num_scans * MBAND_NUM_DETECTORS;
168 
169  // Now that we know dims, prep for additional one-time reads by allocating memory
170  // for arrays (declared static above). (Memory is freed in closel1b_oci.)
171  //tot_num_bands = num_blue_bands + num_red_bands + num_SWIR_bands;
172  tot_num_bands = 62;
173 
174  // Get group id from L1B file for GROUP scan_line_attributes.
175  if ((nc_inq_grp_ncid(ncid_L1B, "scan_line_attributes", &l1bScanLineGrp)) == NC_NOERR) {
176  } else {
177  fprintf(stderr, "-E- Error finding scan_line_attributes.\n");
178  exit(EXIT_FAILURE);
179  }
180 
181  // get mirror side from this GROUP
182  status = nc_inq_varid(l1bScanLineGrp, "HAM_side", &HAMSideId);
183  check_err(status, __LINE__, __FILE__);
184  // get scan quality from this group
185  status = nc_inq_varid(l1bScanLineGrp, "scan_quality_flags", &scanQualityId);
186  check_err(status, __LINE__, __FILE__);
187 
188  // Get attribute values (from l1_viirs_nc.c) e.g. "time_coverage_start" and "time_coverage_end" to derive
189  // time_interval etc. Note that time_coverage_end is the START of the last scan.
190  nc_type vr_type; // attribute type
191  size_t vr_len; // attribute length
192  // Comment out extract_pixel_start/stop stuff
193  /*
194  if ((nc_inq_att(fileID, NC_GLOBAL, "extract_pixel_start", &vr_type, &vr_len) == 0)) {
195  status = nc_get_att_int(fileID, NC_GLOBAL, "extract_pixel_start", &extract_pixel_start);
196  check_err(status, __LINE__, __FILE__);
197  extract_pixel_start--; // Attribute is one-based
198  status = nc_get_att_int(fileID, NC_GLOBAL, "extract_pixel_stop", &extract_pixel_stop);
199  check_err(status, __LINE__, __FILE__);
200  extract_pixel_stop--; // Attribute is one-based
201  if (npix != (extract_pixel_stop - extract_pixel_start + 1)) {
202  fprintf(stderr, "-E- Problem with the extracted L1B file pixel dimension.\n");
203  printf(" npix(%d), extract_pixel_stop(%d), extract_pixel_start(%d) do not work together.\n",
204  npix, extract_pixel_stop, extract_pixel_start);
205  exit(EXIT_FAILURE);
206  }
207  }
208  */
209 
210  if (want_verbose) {
211  printf("OCI L1B Npix :%d Nlines:%d\n", (int)num_pixels, nline);
212  } // want_verbose
213 
214  // get start and end time
215  status = nc_inq_attlen(ncid_L1B, NC_GLOBAL, "time_coverage_start", &att_len);
216  check_err(status, __LINE__, __FILE__);
217 
218  // allocate required space before retrieving values
219  fltime = (char *) malloc(att_len + 1); // + 1 for trailing null
220 
221  // get attribute values
222  status = nc_get_att_text(ncid_L1B, NC_GLOBAL, "time_coverage_start", fltime);
223  check_err(status, __LINE__, __FILE__);
224  fltime[att_len] = '\0';
225  // isodate2ydmsec(fltime, (int32_t *) &year,(int32_t *) &day, (int32_t *) &msec);
226 
227  // Convert "time_coverage_start" ISO string to unix (seconds since 1/1/1970)
228  starttime = lastvalidtime = isodate2unix(fltime);
229  free(fltime);
230 
231  status = nc_inq_attlen(ncid_L1B, NC_GLOBAL, "time_coverage_end", &att_len);
232  check_err(status, __LINE__, __FILE__);
233 
234  // allocate required space before retrieving values
235  fltime = (char *) malloc(att_len + 1); // + 1 for trailing null
236 
237  // get attribute values
238  status = nc_get_att_text(ncid_L1B, NC_GLOBAL, "time_coverage_end", fltime);
239  check_err(status, __LINE__, __FILE__);
240  fltime[att_len] = '\0';
241 
242  // Convert "time_coverage_stop" ISO string to unix (seconds since 1/1/1970)
243  double stoptime = isodate2unix(fltime);
244  free(fltime);
245 
246  // time_interval may be used in readl1b_oci if there is not a good scan time (per scan)
247  time_interval = (stoptime - starttime)/(num_scans-1); // secs per scan
248 
249  if ((nc_inq_att(ncid_L1B, NC_GLOBAL, "orbit_number", &vr_type, &vr_len) == NC_NOERR)) {
250  status = nc_get_att_int(ncid_L1B, NC_GLOBAL, "orbit_number", &orbit_number);
251  check_err(status, __LINE__, __FILE__);
252  } else {
253  orbit_number = 0;
254  }
255 
256 
257  // Identify the "observation_data" GROUP and its "Lt_" vars, to be used later by readl1b_oci
258  // Store the ids in static variables so we don't have to do nc_inq_grp_ncid per scan.
259  if ((nc_inq_grp_ncid(ncid_L1B, "observation_data", &l1bObservationGrp)) == NC_NOERR) {
260  // Get varids for each of the Lt_*
261  status = nc_inq_varid(l1bObservationGrp, "Lt_blue", &Lt_blue_varid);
262  check_err(status, __LINE__, __FILE__);
263  status = nc_inq_varid(l1bObservationGrp, "Lt_red", &Lt_red_varid);
264  check_err(status, __LINE__, __FILE__);
265  status = nc_inq_varid(l1bObservationGrp, "Lt_SWIR", &Lt_SWIR_varid);
266  check_err(status, __LINE__, __FILE__);
267  } else {
268  fprintf(stderr, "-E- Error finding observation_data.\n");
269  exit(EXIT_FAILURE);
270  }
271 
272  file->sd_id = ncid_L1B;
273  file->nbands = tot_num_bands;
274  file->npix = num_pixels;
275  file->nscan = nline;
276  file->ndets = MBAND_NUM_DETECTORS;
277  file->terrain_corrected = 1; // presumed.
278  file->orbit_number = orbit_number;
279  strcpy(file->spatialResolution, "??? m");
280 
281  rdsensorinfo(file->sensorID, l1_input->evalmask,
282  "Fobar", (void **) &Fobar);
283 
284  if (want_verbose)
285  printf("file->nbands = %d\n", (int) file->nbands);
286 
287  // Setup geofile pointers
288  status = nc_inq_grp_ncid(ncid_L1B, "geolocation_data", &geoGeolocationGrp);
289  check_err(status, __LINE__, __FILE__);
290  status = nc_inq_varid(geoGeolocationGrp, "longitude", &lonId);
291  check_err(status, __LINE__, __FILE__);
292  status = nc_inq_var_fill(geoGeolocationGrp, lonId, NULL, &lonGeoFillValue);
293  check_err(status, __LINE__, __FILE__);
294  status = nc_inq_varid(geoGeolocationGrp, "latitude", &latId);
295  check_err(status, __LINE__, __FILE__);
296  status = nc_inq_var_fill(geoGeolocationGrp, latId, NULL, &latGeoFillValue);
297  check_err(status, __LINE__, __FILE__);
298  status = nc_inq_varid(geoGeolocationGrp, "sensor_zenith", &senzId);
299  check_err(status, __LINE__, __FILE__);
300  status = nc_inq_var_fill(geoGeolocationGrp, senzId, NULL, &senzGeoFillValue);
301  check_err(status, __LINE__, __FILE__);
302  status = nc_inq_varid(geoGeolocationGrp, "sensor_azimuth", &senaId);
303  check_err(status, __LINE__, __FILE__);
304  status = nc_inq_var_fill(geoGeolocationGrp, senaId, NULL, &senaGeoFillValue);
305  check_err(status, __LINE__, __FILE__);
306  status = nc_inq_varid(geoGeolocationGrp, "solar_zenith", &solzId);
307  check_err(status, __LINE__, __FILE__);
308  status = nc_inq_var_fill(geoGeolocationGrp, solzId, NULL, &solzGeoFillValue);
309  check_err(status, __LINE__, __FILE__);
310  status = nc_inq_varid(geoGeolocationGrp, "solar_azimuth", &solaId);
311  check_err(status, __LINE__, __FILE__);
312  status = nc_inq_var_fill(geoGeolocationGrp, solaId, NULL, &solaGeoFillValue);
313  check_err(status, __LINE__, __FILE__);
314  //status = nc_inq_varid(geoGeolocationGrp, "quality_flag", &pixelQualityId);
315  //check_err(status, __LINE__, __FILE__);
316 
317  //status = nc_inq_grp_ncid(ncid_L1B, "navigation_data", &geoNavigationGrp);
318  //check_err(status, __LINE__, __FILE__);
319  //angId not yet in GEO file
320  //status = nc_inq_varid(geoNavigationGrp, "att_ang_mid", &angId);
321  //check_err(status, __LINE__, __FILE__);
322  //status = nc_inq_varid(geoNavigationGrp, "orb_pos", &posId);
323  //check_err(status, __LINE__, __FILE__);
324  //status = nc_inq_varid(geoNavigationGrp, "orb_vel", &velId);
325  //check_err(status, __LINE__, __FILE__);
326 
327  //status = nc_inq_grp_ncid(ncid_L1B, "scan_line_attributes", &geoScanLineGrp);
328  //check_err(status, __LINE__, __FILE__);
329  // The following field is not yet in the GEO file
330  //status = nc_inq_varid(geoScanLineGrp, "scan_quality", &scanQualityId);
331  //check_err(status, __LINE__, __FILE__);
332 
333 
334  // Setup the fill values for the geo products
335  // Replaced HAWKEYE with OCI.
336  productInfo_t* info = allocateProductInfo();
337  status = findProductInfo("lat", OCI, info);
338  if (status)
339  latL2FillValue = info->fillValue;
340  status = findProductInfo("lon", OCI, info);
341  if (status)
342  lonL2FillValue = info->fillValue;
343  status = findProductInfo("sena", OCI, info);
344  if (status)
345  senaL2FillValue = info->fillValue;
346  status = findProductInfo("senz", OCI, info);
347  if (status)
348  senzL2FillValue = info->fillValue;
349  status = findProductInfo("sola", OCI, info);
350  if (status)
351  solaL2FillValue = info->fillValue;
352  status = findProductInfo("solz", OCI, info);
353  if (status)
354  solzL2FillValue = info->fillValue;
355  freeProductInfo(info);
356 
357  return (LIFE_IS_GOOD);
358 }
359 
360 
373 int readl1b_oci(filehandle *file, int32_t line, l1str *l1rec) {
374 
375  // Declarations from l1_viirs_nc.c
376  int32_t ip;
377  int i;
378  double scan_sec;
379  int16_t scan_year, scan_day;
380 
381  int status;
382  size_t start[] = { 0, 0, 0 };
383  size_t count[] = { 1, 1, 1 };
384 
385  int32_t scan = line / MBAND_NUM_DETECTORS;
386 
387  // Additional declarations
388  int band_num, pixel_num;
389  int blue_band_num, red_band_num, SWIR_band_num;
390  int16_t scan_month, scan_dom;
391  double scan_delta_time_ms = 0; // ms since start of image per scan
392 
393  //printf("reading oci l1b file\n");
394  for (ip = 0; ip < num_pixels; ip++) {
395  l1rec->pixnum[ip] = ip + extract_pixel_start;
396  }
397 
398  // If first call,
399  if (firstCall) {
400  firstCall = 0;
401 
402  // One-time memory allocations
403  tmpShort = (short *) malloc(num_pixels * sizeof(short));
404  tmpByte = (unsigned char *) malloc(num_pixels);
405  Lt = allocate2d_double(tot_num_bands, num_pixels);
406  Lt_blue = allocate2d_double(num_blue_bands, num_pixels);
407  Lt_red = allocate2d_double(num_red_bands, num_pixels);
408  Lt_SWIR = allocate2d_double(num_SWIR_bands, num_pixels);
409  }
410 
411  // l1rec->sensorID = file->sensorID;
412  l1rec->npix = file->npix;
413 
414  // Time
415  // Get delta_time
416  //start[0] = line;
417  //count[0] = 1; // 1 scan at a time
418  //count[1] = 0;
419  //count[2] = 0;
420  //status = nc_inq_varid(l1bScanLineGrp, "delta_time", &varid);
421  //if (status != NC_NOERR) {
422  // fprintf(stderr, "-E- Error finding scan_delta_time_ms.\n");
423  // exit(EXIT_FAILURE);
424  //}
425  //status = nc_get_vara_double(l1bScanLineGrp, varid, start, count, &scan_delta_time_ms);
426  //if (status != NC_NOERR) {
427  // fprintf(stderr, "-E- Error reading scan_delta_time_ms.\n");
428  // exit(EXIT_FAILURE);
429  //};
430  // Set lastvalidtime, the start of this scan (in secs since 1/1/1970)
431  if (scan_delta_time_ms > 0) {
432  // starttime is the start of this granule (L1B "time_coverage_start"), converted to secs since 1/1/1970
433  // scan_delta_time_ms is when this scan starts (number of milliseconds since the start of this granule)
434  lastvalidtime = starttime + scan_delta_time_ms/1000;
435  } else {
436  lastvalidtime = lastvalidtime + (time_interval * (line - lastvalidscan));
437  }
438  // Set scan_year, scan_day, scan_sec
439  unix2yds(lastvalidtime, &scan_year, &scan_day, &scan_sec);
440  // Store lastvalidtime in l1rec->scantime
441  l1rec->scantime = lastvalidtime;
442  // Convert lastvalidtime to julian date
443  yd2md(scan_year, scan_day, &scan_month, &scan_dom);
444 
445  // GROUP earth_view_data
446  // Get Lt_blue[num_blue_bands][num_pixels] from this GROUP just for THIS scan
447  start[1] = line;
448  count[1] = 1; // 1 line at a time
449  start[2] = 0;
450  count[2] = num_pixels;
451  for (blue_band_num=0; blue_band_num<num_blue_bands; blue_band_num++) {
452  start[0] = blue_band_num;
453  count[0] = 1;
454  status = nc_get_vara_double(l1bObservationGrp, Lt_blue_varid, start, count, &Lt_blue[blue_band_num][0]);
455  if (status != NC_NOERR) {
456  fprintf(stderr, "-E- Error reading Lt_blue.\n");
457  exit(EXIT_FAILURE);
458  }
459  }
460 
461 
462  // Get Lt_red[num_red_bands][num_pixels] from this GROUP just for THIS scan
463  start[1] = line;
464  count[1] = 1; // 1 line at a time
465  start[2] = 0;
466  count[2] = num_pixels;
467  for (red_band_num=0; red_band_num<num_red_bands; red_band_num++) {
468  start[0] = red_band_num;
469  count[0] = 1;
470  status = nc_get_vara_double(l1bObservationGrp, Lt_red_varid, start, count, &Lt_red[red_band_num][0]);
471  if (status != NC_NOERR) {
472  fprintf(stderr, "-E- Error reading Lt_red.\n");
473  exit(EXIT_FAILURE);
474  }
475  }
476 
477  // Get Lt_SWIR[num_SWIR_bands][num_pixels] from this GROUP just for THIS scan
478  start[1] = line;
479  count[1] = 1; // 1 line at a time
480  start[2] = 0;
481  count[2] = num_pixels;
482  for (SWIR_band_num=0; SWIR_band_num<num_SWIR_bands; SWIR_band_num++) {
483  start[0] = SWIR_band_num;
484  count[0] = 1;
485  status = nc_get_vara_double(l1bObservationGrp, Lt_SWIR_varid, start, count, &Lt_SWIR[SWIR_band_num][0]);
486  if (status != NC_NOERR) {
487  fprintf(stderr, "-E- Error reading Lt_SWIR.\n");
488  exit(EXIT_FAILURE);
489  }
490  }
491 
492  // Set l1rec->Lt (fudged for now)
493  for (pixel_num = 0; pixel_num < num_pixels; pixel_num++) {
494  blue_band_num = 0;
495  red_band_num = 0;
496  SWIR_band_num = 0;
497  for (band_num = 0; band_num < tot_num_bands; band_num++) {
498  //ipb = ip * nbands + ib
499  // Fudge for now the map from source bands to target bands e.g. take every other source band
500  if (band_num < 30) {
501  l1rec->Lt[pixel_num*tot_num_bands + band_num] = (float)Lt_blue[blue_band_num][pixel_num];
502  blue_band_num = blue_band_num + 2;
503  }
504  if ((band_num >= 30) && (band_num < 60)) {
505  l1rec->Lt[pixel_num*tot_num_bands + band_num] = (float)Lt_red[red_band_num][pixel_num];
506  red_band_num = red_band_num + 2;
507  }
508  if (band_num >= 60) {
509  l1rec->Lt[pixel_num*tot_num_bands + band_num] = (float)Lt_SWIR[SWIR_band_num][pixel_num];
510  SWIR_band_num = SWIR_band_num + 2;
511  }
512  }
513  }
514 
515 
516 
517 
518 
519  // first check the scan quality flag
520  // 1 SCE_side_A_B
521  // 2 SCE_side_invalid
522  // 4 Sector_rotation
523  // 8 Encoder_degraded
524  // 16 SAA
525  // 32 Solar_eclipse
526  // 64 Lunar_eclipse
527  // 128 HAM_side
528  //
529  // Sector_rotation
530  short scanQualityWarnMask = 2 | 8 | 128;
531  short scanQualityFailMask = 4;
532  static short scanQualityFlag = 0;
533 
534  if (scan != prevScan) {
535  start[0] = scan;
536  status = nc_get_var1_short(l1bScanLineGrp, scanQualityId, start, &scanQualityFlag);
537  check_err(status, __LINE__, __FILE__);
538  }
539  if (scanQualityFlag & scanQualityFailMask) {
540  for (ip = 0; ip < num_pixels; ip++)
541  l1rec->flags[ip] |= NAVFAIL;
542  return 0;
543  }
544  if (scanQualityFlag & scanQualityWarnMask) {
545  for (ip = 0; ip < num_pixels; ip++)
546  l1rec->flags[ip] |= NAVWARN;
547  }
548 
549 
550  static unsigned char HAMSideVal = 0;
551  if (scan != prevScan) {
552  start[0] = scan;
553  status = nc_get_var1_uchar(l1bScanLineGrp, HAMSideId, start, &HAMSideVal);
554  check_err(status, __LINE__, __FILE__);
555  }
556  l1rec->mside = HAMSideVal;
557 
558  // set up to read all pixels of the line.
559  start[0] = line;
560  start[1] = 0;
561  count[0] = 1;
562  count[1] = num_pixels; // read all pixels
563 
564  status = nc_get_vara_float(geoGeolocationGrp, latId, start, count, l1rec->lat);
565  check_err(status, __LINE__, __FILE__);
566  for (i = 0; i < num_pixels; i++)
567  if (l1rec->lat[i] == latGeoFillValue)
568  l1rec->lat[i] = latL2FillValue;
569 
570  status = nc_get_vara_float(geoGeolocationGrp, lonId, start, count, l1rec->lon);
571  check_err(status, __LINE__, __FILE__);
572  for (i = 0; i < num_pixels; i++)
573  if (l1rec->lon[i] == lonGeoFillValue)
574  l1rec->lon[i] = lonL2FillValue;
575 
576  status = nc_get_vara_short(geoGeolocationGrp, solzId, start, count, tmpShort);
577  check_err(status, __LINE__, __FILE__);
578  for (i = 0; i < num_pixels; i++)
579  if (tmpShort[i] == solzGeoFillValue)
580  l1rec->solz[i] = solzL2FillValue;
581  else
582  l1rec->solz[i] = tmpShort[i] * 0.01;
583 
584  status = nc_get_vara_short(geoGeolocationGrp, solaId, start, count, tmpShort);
585  check_err(status, __LINE__, __FILE__);
586  for (i = 0; i < num_pixels; i++)
587  if (tmpShort[i] == solaGeoFillValue)
588  l1rec->sola[i] = solaL2FillValue;
589  else
590  l1rec->sola[i] = tmpShort[i] * 0.01;
591 
592  status = nc_get_vara_short(geoGeolocationGrp, senzId, start, count, tmpShort);
593  check_err(status, __LINE__, __FILE__);
594  for (i = 0; i < num_pixels; i++)
595  if (tmpShort[i] == senzGeoFillValue)
596  l1rec->senz[i] = senzL2FillValue;
597  else
598  l1rec->senz[i] = tmpShort[i] * 0.01;
599 
600  status = nc_get_vara_short(geoGeolocationGrp, senaId, start, count, tmpShort);
601  check_err(status, __LINE__, __FILE__);
602  for (i = 0; i < num_pixels; i++)
603  if (tmpShort[i] == senaGeoFillValue)
604  l1rec->sena[i] = senaL2FillValue;
605  else
606  l1rec->sena[i] = tmpShort[i] * 0.01;
607 
608  //status = nc_get_vara_float(geoNavigationGrp, angId, s3, c3, ang);
609  //check_err(status, __LINE__, __FILE__);
610  //status = nc_get_vara_float(geoNavigationGrp, posId, s3, c3, pos);
611  //check_err(status, __LINE__, __FILE__);
612  //status = nc_get_vara_float(geoNavigationGrp, velId, s3, c3, vel);
613  //check_err(status, __LINE__, __FILE__);
614  //for (i = 0; i < 3; i++) {
615  // pos[i] /= 1000.; // m -> km
616  // vel[i] /= 1000.; // m/s -> km/s
617  //}
618 
619  // Compute polarization rotation angles
620  // Skipping for now since no ang [SBL]
621  /*
622  float sen_mat[3][3], coeff[10];
623  double mnorm[3];
624  ocorient_(pos, vel, ang, sen_mat, coeff);
625  for (i = 0; i < 3; i++)
626  mnorm[i] = sen_mat[i][0];
627  compute_alpha(l1rec->lon, l1rec->lat,
628  l1rec->senz, l1rec->sena,
629  mnorm, l1rec->npix, l1rec->alpha);
630  */
631 
632  // Check pixel values
633  //status = nc_get_vara_uchar(geoGeolocationGrp, pixelQualityId, start, count, tmpByte);
634  //check_err(status, __LINE__, __FILE__);
635  // 1 Input_invalid
636  // 2 Pointing_bad
637  // 4 Terrain_bad
638  unsigned char qualityFailMask = 1 | 2;
639  unsigned char qualityWarnMask = 4;
640  for (i = 0; i < num_pixels; i++) {
641  if (tmpByte[i] & qualityFailMask)
642  l1rec->flags[i] |= NAVFAIL;
643  if (tmpByte[i] & qualityWarnMask)
644  l1rec->flags[i] |= NAVWARN;
645  }
646 
647  // Earth-sun distance correction for this scan
648  static double esdist = -999.9;
649  if (scan != prevScan) {
650  int16_t year, day;
651  double dsec;
652  unix2yds(l1rec->scantime, &year, &day, &dsec);
653  int32_t yr = (int32_t) year;
654  int32_t dy = (int32_t) day;
655  int32_t msec = (int32_t) (dsec * 1000.0);
656  esdist = esdist_(&yr, &dy, &msec);
657  }
658  l1rec->fsol = pow(1.0 / esdist, 2);
659 
660  // Skipping this VIIRS-related stuff
661  /*
662  int nbands = 16; //, nRSBbands = 10, nCIRbands = 1, nTEBbands = 5;
663 
664  // read in calibrated L1B data
665  static float *l1bptrs[16] = {0x0, 0x0, 0x0, 0x0, 0x0, 0x0, 0x0, 0x0,
666  0x0, 0x0, 0x0, 0x0, 0x0, 0x0, 0x0, 0x0};
667 
668  static char *bandType[16] = {"RSB", "RSB", "RSB", "RSB", "RSB", "RSB",
669  "RSB", "RSB", "CIR", "RSB", "RSB", "TEB",
670  "TEB", "TEB", "TEB", "TEB"};
671 
672  // Note: l1bptrs arrays are 3200 pixels wide
673  int oldVerbose = want_verbose;
674  want_verbose = 0;
675  VcstViirsCal_calibrateMOD(line, nbands, l1bptrs);
676  want_verbose = oldVerbose;
677  */
678 
679  l1rec->detnum = line % file->ndets;
680 
681  // Skipping this VIIRS-related stuff
682  /*
683  int irsb = 0, iteb = 0;
684  int l1bptrs_scan_idx;
685  for (ib = 0; ib < nbands; ib++) {
686 
687  // get specific f table cal correction
688  f_corr = (f_cal_corr == NULL) ? 1.0
689  : f_cal_corr[ib][l1rec->detnum][l1rec->mside];
690 
691  if (strcmp(bandType[ib], "TEB") == 0) {
692 
693  for (ip = 0; ip < num_pixels; ip++) {
694  ipb = ip * NBANDSIR + iteb;
695  l1rec->Ltir[ipb] = 0;
696 
697  l1bptrs_scan_idx = l1rec->detnum * 3200 + ip + extract_pixel_start;
698 
699  if (l1bptrs[ib][l1bptrs_scan_idx] != -32767) {
700  l1rec->Ltir[ipb] = l1bptrs[ib][l1bptrs_scan_idx] / 10.0;
701 
702  // Apply F-factor
703  l1rec->Ltir[ipb] *= f_corr;
704  }
705 
706  }
707  iteb++;
708 
709  } else if (strcmp(bandType[ib], "CIR") == 0) {
710 
711  for (ip = 0; ip < num_pixels; ip++) {
712 
713  l1bptrs_scan_idx = l1rec->detnum * 3200 + ip + extract_pixel_start;
714 
715  if (l1bptrs[ib][l1bptrs_scan_idx] != -32767) {
716  l1rec->rho_cirrus[ip] = l1bptrs[ib][l1bptrs_scan_idx];
717 
718  // Normalize reflectance by solar zenith angle
719  l1rec->rho_cirrus[ip] /= cos(l1rec->solz[ip] / RADEG);
720 
721  // Apply F-factor
722  l1rec->rho_cirrus[ip] *= f_corr;
723  }
724 
725  }
726 
727  } else if (strcmp(bandType[ib], "RSB") == 0) {
728 
729  l1rec->Fo[irsb] = Fobar[irsb] * l1rec->fsol;
730 
731  // copy to Lt record.
732  for (ip = 0; ip < num_pixels; ip++) {
733  ipb = ip * l1rec->l1file->nbands + irsb;
734 
735  l1bptrs_scan_idx = l1rec->detnum * 3200 + ip + extract_pixel_start;
736 
737  if (l1bptrs[ib][l1bptrs_scan_idx] != -32767) {
738  l1rec->Lt[ipb] = l1bptrs[ib][l1bptrs_scan_idx];
739 
740  // convert from reflectance to radiance
741  l1rec->Lt[ipb] *= l1rec->Fo[irsb] / PI;
742 
743  // Apply F-factor
744  l1rec->Lt[ipb] *= f_corr;
745  }
746 
747  }
748  irsb++;
749  } // if RSB
750 
751  } // for ib
752  */
753 
754  // Skipping for now to avoid "No brightness temperature conversion provided for this sensor" [SBL]
755  //radiance2bt(l1rec, -1); // calculate brightness temperature
756 
757  // Bowtie stuff is specific to VIIRS, so skip it.
758  //for (ip = 0; ip < num_pixels; ip++) {
759  // flag_bowtie_deleted(l1rec, ip, extract_pixel_start);
760  //}
761 
762  prevScan = scan;
763 
764 
765  return (LIFE_IS_GOOD);
766 }
767 
768 
774 int closel1b_oci(filehandle *file) {
775  int status;
776 
777  printf("Closing oci l1b file\n");
778  status = nc_close(file->sd_id);
779  check_err(status, __LINE__, __FILE__);
780 
781  // Free memory
782  // From readl1b_oci
783  if (tmpShort) free(tmpShort);
784  if (tmpByte) free(tmpByte);
785  if (Lt) free2d_double(Lt);
786  if (Lt_blue) free2d_double(Lt_blue);
787  if (Lt_red) free2d_double(Lt_red);
788  if (Lt_SWIR) free2d_double(Lt_SWIR);
789  // From openl1b_oci
790 
791 
792  return (LIFE_IS_GOOD);
793 }
794 
795 
796 
797 
798 
void free2d_double(double **p)
Free a two-dimensional array created by allocate2d_double.
Definition: allocate2d.c:170
void freeProductInfo(productInfo_t *info)
int32_t day
int status
Definition: l1_czcs_hdf.c:32
void check_err(const int stat, const int line, const char *file)
Definition: nc4utils.c:35
#define OCI
Definition: sensorDefs.h:42
#define NULL
Definition: decode_rs.h:63
int openl1b_oci(filehandle *file)
Definition: l1b_oci.c:113
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
int32 * msec
Definition: l1_czcs_hdf.c:31
double esdist_(int32_t *year, int32_t *day, int32_t *msec)
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
#define LIFE_IS_GOOD
Definition: passthebuck.h:4
int readl1b_oci(filehandle *file, int32_t line, l1str *l1rec)
Definition: l1b_oci.c:373
void yd2md(int16_t year, int16_t doy, int16_t *month, int16_t *dom)
Definition: yd2md.c:6
productInfo_t * allocateProductInfo()
int closel1b_oci(filehandle *file)
Definition: l1b_oci.c:774
l1_input_t * l1_input
Definition: l1_options.c:9
int want_verbose
void unix2yds(double usec, short *year, short *day, double *secs)
Utility functions for allocating and freeing two-dimensional arrays of various types.
int findProductInfo(const char *productName, int sensorId, productInfo_t *info)
double ** allocate2d_double(size_t h, size_t w)
Allocate a two-dimensional array of type double of a given size.
Definition: allocate2d.c:153
real *8 function esdist(iyr, iday, msec)
Definition: esdist.f:3
int32_t rdsensorinfo(int32_t, int32_t, const char *, void **)
Definition: rdsensorinfo.c:69
#define MBAND_NUM_DETECTORS
Definition: l1b_oci.c:51
#define NAVWARN
Definition: l2_flags.h:27
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")
#define NAVFAIL
Definition: l2_flags.h:36
double isodate2unix(const char *isodate)
Definition: unix2isodate.c:61
int count
Definition: decode_rs.h:79