NASA Logo
Ocean Color Science Software

ocssw V2022
l1_oci.c
Go to the documentation of this file.
1 /* ============================================================================ */
2 /* module l1b_oci.c - functions to read OCI L1B for l2gen */
3 /* Written By: Don Shea */
4 /* */
5 /* ============================================================================ */
6 
7 /* Issues:
8  -Assume num_pixels = ccd_pixels = SWIR_pixels
9 */
10 
11 /*
12  further var population by Tincho, 6/16->
13 */
14 
15 #include <netcdf.h>
16 #include "l1_oci.h"
17 
18 #include "l1.h"
19 #include <nc4utils.h>
20 #include "libnav.h"
21 #include <stdio.h>
22 #include <math.h>
23 #include <allocate2d.h>
24 
25 #define SATURATION_BIT 1
26 
27 //static float *Fobar; // reflectance to radiance conversion factors
28 static int extract_pixel_start = 0;
29 static int use_rhot = 0;
30 
31 static int bad_num_bands = 0;
32 
33 static short *tmpShort;
34 
35 // whole file stuff
36 static size_t expected_num_blue_bands = 119;
37 static size_t expected_num_red_bands = 163;
38 static size_t expected_num_SWIR_bands = 9;
39 
40 static size_t num_scans, num_pixels;
41 static size_t num_blue_bands, num_red_bands, num_SWIR_bands;
42 static size_t tot_num_bands = 286;
43 static int ncid_L1B;
44 
45 // scan line attributes
46 static double *scan_time; // seconds of day
47 static double file_start_day; // unix time of start day 00:00:00
48 static float *tilt;
49 static unsigned char *scanQual;
50 static uint8_t *hamside;
51 
52 // geolocation data
53 static int geolocationGrp; // netCDF groupid
54 static int lonId, latId, heightId, senzId, senaId, solzId, solaId; // netCDF varids
55 static float latFillValue = BAD_FLT;
56 static float lonFillValue = BAD_FLT;
57 static short heightFillValue = BAD_FLT;
58 static short senzFillValue = BAD_FLT;
59 static float senzScale = 0.01;
60 static float senzOffset = 0.0;
61 static short senaFillValue = BAD_FLT;
62 static float senaScale = 0.01;
63 static float senaOffset = 0.0;
64 static short solzFillValue = BAD_FLT;
65 static float solzScale = 0.01;
66 static float solzOffset = 0.0;
67 static short solaFillValue = BAD_FLT;
68 static float solaScale = 0.01;
69 static float solaOffset = 0.0;
70 
71 // Observation data
72 static int observationGrp,navigationGrp;
73 static int Lt_blueId, Lt_redId, Lt_SWIRId, tiltId, qual_SWIRId;
74 static float **Lt_blue; //[num_blue_bands][num_pixels], This scan
75 static float **Lt_red; //[num_red_bands][num_pixels], This scan
76 static float **Lt_SWIR; //[num_SWIR_bands][num_pixels], This scan
77 static float Lt_blueFillValue = BAD_FLT;
78 static float Lt_redFillValue = BAD_FLT;
79 static float Lt_SWIRFillValue = BAD_FLT;
80 static float tiltFillValue = BAD_FLT;
81 static uint8_t **qual_SWIR;
82 static float *blue_solar_irradiance; // [num_blue_bands]
83 static float *red_solar_irradiance; // [num_red_bands]
84 static float *SWIR_solar_irradiance; // [num_SWIR_bands]
85 static double earth_sun_distance_correction;
86 
99 int openl1_oci(filehandle * file) {
100  int dimid, status;
101 
102  // Open the netcdf4 input file
103  printf("Opening OCI L1B file\n");
104  status = nc_open(file->name, NC_NOWRITE, &ncid_L1B);
105  if (status != NC_NOERR) {
106  fprintf(stderr, "-E- %s line %d: nc_open(%s) failed.\n",
107  __FILE__, __LINE__, file->name);
108  return (1);
109  }
110 
111  // num_scans
112  status = nc_inq_dimid(ncid_L1B, "scans", &dimid);
113  if (status != NC_NOERR) {
114  // Try to look for the old name for this dimension
115  if (nc_inq_dimid(ncid_L1B, "number_of_scans", &dimid) != NC_NOERR) {
116  fprintf(stderr, "-E- Error reading scan dimension.\n");
117  exit(EXIT_FAILURE);
118  }
119  }
120  nc_inq_dimlen(ncid_L1B, dimid, &num_scans);
121 
122  // num_pixels
123  status = nc_inq_dimid(ncid_L1B, "pixels", &dimid);
124  if (status != NC_NOERR) {
125  // If "pixels" doesn't exist, it's the old "ccd_pixels"
126  if (nc_inq_dimid(ncid_L1B, "ccd_pixels", &dimid)) {
127  fprintf(stderr, "-E- Error reading num_pixels.\n");
128  exit(EXIT_FAILURE);
129  }
130  }
131  nc_inq_dimlen(ncid_L1B, dimid, &num_pixels);
132 
133  // num_blue_bands
134  status = nc_inq_dimid(ncid_L1B, "blue_bands", &dimid);
135  if (status != NC_NOERR) {
136  fprintf(stderr, "-E- Error reading num_blue_bands.\n");
137  exit(EXIT_FAILURE);
138  }
139  nc_inq_dimlen(ncid_L1B, dimid, &num_blue_bands);
140  if(num_blue_bands < expected_num_blue_bands) {
141  fprintf(stderr, "-W- Not enough blue bands, expecting %d, found %d.\n",
142  (int)expected_num_blue_bands, (int)num_blue_bands);
143  bad_num_bands = 1;
144  }
145 
146  // num_red_bands
147  status = nc_inq_dimid(ncid_L1B, "red_bands", &dimid);
148  if (status != NC_NOERR) {
149  fprintf(stderr, "-E- Error reading num_red_bands.\n");
150  exit(EXIT_FAILURE);
151  };
152  nc_inq_dimlen(ncid_L1B, dimid, &num_red_bands);
153  if(num_red_bands < expected_num_red_bands) {
154  fprintf(stderr, "-W- Not enough red bands, expecting %d, found %d.\n",
155  (int)expected_num_red_bands, (int)num_red_bands);
156  bad_num_bands = 1;
157  }
158 
159  // num_SWIR_bands
160  status = nc_inq_dimid(ncid_L1B, "SWIR_bands", &dimid);
161  if (status != NC_NOERR) {
162  fprintf(stderr, "-E- Error reading num_SWIR_bands.\n");
163  exit(EXIT_FAILURE);
164  };
165  nc_inq_dimlen(ncid_L1B, dimid, &num_SWIR_bands);
166  if(num_SWIR_bands < expected_num_SWIR_bands) {
167  fprintf(stderr, "-E- Not enough SWIR bands, expecting %d, found %d.\n",
168  (int)expected_num_SWIR_bands, (int)num_SWIR_bands);
169  exit(EXIT_FAILURE);
170  }
171 
172  if (want_verbose) {
173  printf("OCI L1B Npix :%d Nlines:%d\n", (int)num_pixels, (int)num_scans);
174  } // want_verbose
175 
176  // allocate all of the data
177  tmpShort = (short*) malloc(num_pixels * sizeof(short));
178  scan_time = (double*) malloc(num_scans * sizeof(double));
179 
180  Lt_blue = allocate2d_float(num_blue_bands, num_pixels);
181  Lt_red = allocate2d_float(num_red_bands, num_pixels);
182  Lt_SWIR = allocate2d_float(num_SWIR_bands, num_pixels);
183  tilt = (float*) malloc(num_scans * sizeof(float));
184  hamside = (uint8_t *) malloc(num_scans * sizeof(uint8_t));
185  scanQual = (unsigned char *) malloc(num_scans * sizeof(unsigned char));
186  qual_SWIR = allocate2d_uchar(num_SWIR_bands, num_pixels);
187 
188  // Get group id from L1B file for GROUP scan_line_attributes.
189  int scanLineGrp;
190  if ((nc_inq_grp_ncid(ncid_L1B, "scan_line_attributes", &scanLineGrp)) == NC_NOERR) {
191  } else {
192  fprintf(stderr, "-E- Error finding scan_line_attributes.\n");
193  exit(EXIT_FAILURE);
194  }
195  int varId;
196  double scan_timeFillValue = BAD_FLT;
197  status = nc_inq_varid(scanLineGrp, "time", &varId);
198  if(status == NC_NOERR) {
199  status = nc_inq_var_fill(scanLineGrp, varId, NULL, &scan_timeFillValue);
200  check_err(status, __LINE__, __FILE__);
201  status = nc_get_var_double(scanLineGrp, varId, scan_time);
202  check_err(status, __LINE__, __FILE__);
203  } else {
204  status = nc_inq_varid(scanLineGrp, "ev_mid_time", &varId);
205  check_err(status, __LINE__, __FILE__);
206  status = nc_inq_var_fill(scanLineGrp, varId, NULL, &scan_timeFillValue);
207  check_err(status, __LINE__, __FILE__);
208  status = nc_get_var_double(scanLineGrp, varId, scan_time);
209  check_err(status, __LINE__, __FILE__);
210  }
211 
212  // get HAM side info
213  status = nc_inq_varid(scanLineGrp, "HAM_side", &varId);
214  if( ( status = nc_inq_varid(scanLineGrp, "HAM_side", &varId) ) == NC_NOERR) {
215  status = nc_get_var_ubyte(scanLineGrp, varId, hamside);
216  check_err(status, __LINE__, __FILE__);
217  } else {
218  fprintf(stderr, "-E- Error reading the HAM side data.\n");
219  exit(EXIT_FAILURE);
220  }
221 
222  // get scan quality info
223  status = nc_inq_varid(scanLineGrp, "scan_quality_flags", &varId);
224  if( ( status = nc_inq_varid(scanLineGrp, "scan_quality_flags", &varId) ) == NC_NOERR) {
225  status = nc_get_var_uchar(scanLineGrp, varId, scanQual);
226  check_err(status, __LINE__, __FILE__);
227  } else {
228  fprintf(stderr, "-E- Error reading the scan quality flags.\n");
229  exit(EXIT_FAILURE);
230  }
231 
232  // get start time
233  size_t att_len;
234  status = nc_inq_attlen(ncid_L1B, NC_GLOBAL, "time_coverage_start", &att_len);
235  check_err(status, __LINE__, __FILE__);
236 
237  // allocate required space before retrieving values
238  char* time_str = (char *) malloc(att_len + 1); // + 1 for trailing null
239 
240  // get attribute values
241  status = nc_get_att_text(ncid_L1B, NC_GLOBAL, "time_coverage_start", time_str);
242  check_err(status, __LINE__, __FILE__);
243  time_str[att_len] = '\0';
244 
245  double start_time = isodate2unix(time_str);
246  int16_t syear, smon, sday;
247  double secs;
248  unix2ymds(start_time, &syear, &smon, &sday, &secs);
249  file_start_day = ymds2unix(syear, smon, sday, 0.0);
250 
251  free(time_str);
252 
253  for(int i=0; i<num_scans; i++) {
254  if(scan_time[i] == scan_timeFillValue)
255  scan_time[i] = BAD_FLT;
256  }
257 
258  // Setup geofile pointers
259  status = nc_inq_grp_ncid(ncid_L1B, "geolocation_data", &geolocationGrp);
260  check_err(status, __LINE__, __FILE__);
261  status = nc_inq_varid(geolocationGrp, "longitude", &lonId);
262  check_err(status, __LINE__, __FILE__);
263  status = nc_inq_var_fill(geolocationGrp, lonId, NULL, &lonFillValue);
264  check_err(status, __LINE__, __FILE__);
265  status = nc_inq_varid(geolocationGrp, "latitude", &latId);
266  check_err(status, __LINE__, __FILE__);
267  status = nc_inq_var_fill(geolocationGrp, latId, NULL, &latFillValue);
268  check_err(status, __LINE__, __FILE__);
269  status = nc_inq_varid(geolocationGrp, "height", &heightId);
270  check_err(status, __LINE__, __FILE__);
271  status = nc_inq_var_fill(geolocationGrp, heightId, NULL, &heightFillValue);
272  check_err(status, __LINE__, __FILE__);
273 
274  status = nc_inq_varid(geolocationGrp, "sensor_zenith", &senzId);
275  check_err(status, __LINE__, __FILE__);
276  status = nc_inq_var_fill(geolocationGrp, senzId, NULL, &senzFillValue);
277  check_err(status, __LINE__, __FILE__);
278  status = nc_get_att_float(geolocationGrp, senzId, "scale_factor", &senzScale);
279  // ignore missing
280  //check_err(status, __LINE__, __FILE__);
281  status = nc_get_att_float(geolocationGrp, senzId, "add_offset", &senzOffset);
282  // ignore missing
283  //check_err(status, __LINE__, __FILE__);
284 
285  status = nc_inq_varid(geolocationGrp, "sensor_azimuth", &senaId);
286  check_err(status, __LINE__, __FILE__);
287  status = nc_inq_var_fill(geolocationGrp, senaId, NULL, &senaFillValue);
288  check_err(status, __LINE__, __FILE__);
289  status = nc_get_att_float(geolocationGrp, senaId, "scale_factor", &senaScale);
290  // ignore missing
291  //check_err(status, __LINE__, __FILE__);
292  status = nc_get_att_float(geolocationGrp, senaId, "add_offset", &senaOffset);
293  // ignore missing
294  //check_err(status, __LINE__, __FILE__);
295 
296  status = nc_inq_varid(geolocationGrp, "solar_zenith", &solzId);
297  check_err(status, __LINE__, __FILE__);
298  status = nc_inq_var_fill(geolocationGrp, solzId, NULL, &solzFillValue);
299  check_err(status, __LINE__, __FILE__);
300  status = nc_get_att_float(geolocationGrp, solzId, "scale_factor", &solzScale);
301  // ignore missing
302  //check_err(status, __LINE__, __FILE__);
303  status = nc_get_att_float(geolocationGrp, solzId, "add_offset", &solzOffset);
304  // ignore missing
305  //check_err(status, __LINE__, __FILE__);
306 
307  status = nc_inq_varid(geolocationGrp, "solar_azimuth", &solaId);
308  check_err(status, __LINE__, __FILE__);
309  status = nc_inq_var_fill(geolocationGrp, solaId, NULL, &solaFillValue);
310  check_err(status, __LINE__, __FILE__);
311  status = nc_get_att_float(geolocationGrp, solaId, "scale_factor", &solaScale);
312  // ignore missing
313  //check_err(status, __LINE__, __FILE__);
314  status = nc_get_att_float(geolocationGrp, solaId, "add_offset", &solaOffset);
315  // ignore missing
316  //check_err(status, __LINE__, __FILE__);
317 
318 
319  // get IDs for the observations
320  status = nc_inq_grp_ncid(ncid_L1B, "observation_data", &observationGrp);
321  check_err(status, __LINE__, __FILE__);
322 
323  // Get varids for each of the Lt_*
324  status = nc_inq_varid(observationGrp, "Lt_blue", &Lt_blueId);
325  if(status == NC_NOERR) {
326  use_rhot = 0;
327 
328  status = nc_inq_var_fill(observationGrp, Lt_blueId, NULL, &Lt_blueFillValue);
329  check_err(status, __LINE__, __FILE__);
330 
331  status = nc_inq_varid(observationGrp, "Lt_red", &Lt_redId);
332  check_err(status, __LINE__, __FILE__);
333  status = nc_inq_var_fill(observationGrp, Lt_redId, NULL, &Lt_redFillValue);
334  check_err(status, __LINE__, __FILE__);
335 
336  status = nc_inq_varid(observationGrp, "Lt_SWIR", &Lt_SWIRId);
337  check_err(status, __LINE__, __FILE__);
338  status = nc_inq_var_fill(observationGrp, Lt_SWIRId, NULL, &Lt_SWIRFillValue);
339  check_err(status, __LINE__, __FILE__);
340  } else { // no Lt, so check for rhot
341  use_rhot = 1;
342 
343  status = nc_inq_varid(observationGrp, "rhot_blue", &Lt_blueId);
344  check_err(status, __LINE__, __FILE__);
345  status = nc_inq_var_fill(observationGrp, Lt_blueId, NULL, &Lt_blueFillValue);
346  check_err(status, __LINE__, __FILE__);
347 
348  status = nc_inq_varid(observationGrp, "rhot_red", &Lt_redId);
349  check_err(status, __LINE__, __FILE__);
350  status = nc_inq_var_fill(observationGrp, Lt_redId, NULL, &Lt_redFillValue);
351  check_err(status, __LINE__, __FILE__);
352 
353  status = nc_inq_varid(observationGrp, "rhot_SWIR", &Lt_SWIRId);
354  check_err(status, __LINE__, __FILE__);
355  status = nc_inq_var_fill(observationGrp, Lt_SWIRId, NULL, &Lt_SWIRFillValue);
356  check_err(status, __LINE__, __FILE__);
357  }
358 
359  status = nc_inq_grp_ncid(ncid_L1B, "navigation_data", &navigationGrp);
360  check_err(status, __LINE__, __FILE__);
361  status = nc_inq_varid(navigationGrp, "tilt_angle", &tiltId);
362  if(status){
363  status = nc_inq_varid(navigationGrp, "tilt", &tiltId);
364  if(status)
365  check_err(status, __LINE__, __FILE__);
366  }
367  status = nc_inq_var_fill(navigationGrp, tiltId, NULL, &tiltFillValue);
368  check_err(status, __LINE__, __FILE__);
369  status = nc_get_var_float(navigationGrp,tiltId,tilt);
370  check_err(status, __LINE__, __FILE__);
371 
372  status = nc_inq_varid(observationGrp, "qual_SWIR", &qual_SWIRId);
373  check_err(status, __LINE__, __FILE__);
374 
375  if(use_rhot) {
376  //rdsensorinfo(file->sensorID, l1_input->evalmask, "Fobar", (void **) &Fobar);
377 
378  status = nc_get_att_double(ncid_L1B, NC_GLOBAL, "earth_sun_distance_correction", &earth_sun_distance_correction);
379  check_err(status, __LINE__, __FILE__);
380 
381  int sensorBandGrp;
382  int tmpId;
383  blue_solar_irradiance = malloc(num_blue_bands * sizeof(float));
384  red_solar_irradiance = malloc(num_red_bands * sizeof(float));
385  SWIR_solar_irradiance = malloc(num_SWIR_bands * sizeof(float));
386  status = nc_inq_grp_ncid(ncid_L1B, "sensor_band_parameters", &sensorBandGrp);
387  check_err(status, __LINE__, __FILE__);
388 
389  status = nc_inq_varid(sensorBandGrp, "blue_solar_irradiance", &tmpId);
390  check_err(status, __LINE__, __FILE__);
391  status = nc_get_var_float(sensorBandGrp, tmpId, blue_solar_irradiance);
392  check_err(status, __LINE__, __FILE__);
393 
394  status = nc_inq_varid(sensorBandGrp, "red_solar_irradiance", &tmpId);
395  check_err(status, __LINE__, __FILE__);
396  status = nc_get_var_float(sensorBandGrp, tmpId, red_solar_irradiance);
397  check_err(status, __LINE__, __FILE__);
398 
399  status = nc_inq_varid(sensorBandGrp, "SWIR_solar_irradiance", &tmpId);
400  check_err(status, __LINE__, __FILE__);
401  status = nc_get_var_float(sensorBandGrp, tmpId, SWIR_solar_irradiance);
402  check_err(status, __LINE__, __FILE__);
403  }
404 
405  file->sd_id = ncid_L1B;
406  file->nbands = tot_num_bands;
407  file->npix = num_pixels;
408  file->nscan = num_scans;
409  file->ndets = 1;
410  file->terrain_corrected = 1; // presumed.
411  strcpy(file->spatialResolution, "1000 m");
412 
413  if (want_verbose)
414  printf("file->nbands = %d\n", (int) file->nbands);
415 
416  return (LIFE_IS_GOOD);
417 }
418 
419 
434 int readl1_oci(filehandle *file, int32_t line, l1str *l1rec, int lonlat) {
435 
436  int status;
437  size_t start[] = { 0, 0, 0 };
438  size_t count[] = { 1, 1, 1 };
439 
440  for (int ip = 0; ip < num_pixels; ip++) {
441  l1rec->pixnum[ip] = ip + extract_pixel_start;
442  }
443 
444  l1rec->npix = file->npix;
445  l1rec->scantime = BAD_FLT;
446  l1rec->mside = hamside[line];
447  if (scan_time[line] != BAD_FLT) {
448  l1rec->scantime = file_start_day + scan_time[line];
449  }
450 
451  l1rec->tilt=tilt[line];
452 
453  int16_t syear, sday;
454  double secs;
455  unix2yds(l1rec->scantime, &syear, &sday, &secs);
456 
457  int32_t yr = syear;
458  int32_t dy = sday;
459  int32_t msec = (int32_t) (secs * 1000.0);
460  double esdist = esdist_(&yr, &dy, &msec);
461 
462  l1rec->fsol = pow(1.0 / esdist, 2);
463 
464  start[0] = line;
465  start[1] = 0;
466  start[2] = 0;
467  count[0] = 1;
468  count[1] = num_pixels; // 1 line at a time
469  count[2] = 1;
470  status = nc_get_vara_float(geolocationGrp, latId, start, count, l1rec->lat);
471  check_err(status, __LINE__, __FILE__);
472  status = nc_get_vara_float(geolocationGrp, lonId, start, count, l1rec->lon);
473  check_err(status, __LINE__, __FILE__);
474 
475  // dont need height for lonlat, but do when it is not lonlat
476  if (!lonlat) {
477  status = nc_get_vara_float(geolocationGrp, heightId, start, count, l1rec->height);
478  check_err(status, __LINE__, __FILE__);
479  }
480 
481  for(int i=0; i<num_pixels; i++) {
482  if(l1rec->lat[i] == latFillValue)
483  l1rec->lat[i] = BAD_FLT;
484  if(l1rec->lon[i] == lonFillValue)
485  l1rec->lon[i] = BAD_FLT;
486  }
487 
488  if(scanQual[line] & 1) {
489  for(int i=0; i<num_pixels; i++)
490  l1rec->navwarn[i] = 1;
491  }
492 
493  status = nc_get_vara_short(geolocationGrp, solzId, start, count, tmpShort);
494  check_err(status, __LINE__, __FILE__);
495  for(int i=0; i<num_pixels; i++) {
496  if(tmpShort[i] == solzFillValue)
497  l1rec->solz[i] = BAD_FLT;
498  else
499  l1rec->solz[i] = tmpShort[i] * solzScale + solzOffset;
500  }
501 
502  // lon lat, only need solarz and nothing else
503  // at this point, it should have read in: time, lon, lat, solz
504  if (lonlat)
505  return (LIFE_IS_GOOD);
506 
507  if (bad_num_bands) {
508  fprintf(stderr, "-E- Bad number of bands\n");
509  exit(EXIT_FAILURE);
510  }
511 
512  status = nc_get_vara_short(geolocationGrp, senaId, start, count, tmpShort);
513  check_err(status, __LINE__, __FILE__);
514  for(int i=0; i<num_pixels; i++) {
515  if(tmpShort[i] == senaFillValue)
516  l1rec->sena[i] = BAD_FLT;
517  else
518  l1rec->sena[i] = tmpShort[i] * senaScale + senaOffset;
519  }
520 
521  status = nc_get_vara_short(geolocationGrp, senzId, start, count, tmpShort);
522  check_err(status, __LINE__, __FILE__);
523  for(int i=0; i<num_pixels; i++) {
524  if(tmpShort[i] == senzFillValue)
525  l1rec->senz[i] = BAD_FLT;
526  else
527  l1rec->senz[i] = tmpShort[i] * senzScale + senzOffset;
528  }
529 
530  status = nc_get_vara_short(geolocationGrp, solaId, start, count, tmpShort);
531  check_err(status, __LINE__, __FILE__);
532  for(int i=0; i<num_pixels; i++) {
533  if(tmpShort[i] == solaFillValue)
534  l1rec->sola[i] = BAD_FLT;
535  else
536  l1rec->sola[i] = tmpShort[i] * solaScale + solaOffset;
537  }
538 
539 
540  start[0] = 0;
541  start[1] = line;
542  start[2] = 0;
543  count[0] = expected_num_blue_bands;
544  count[1] = 1; // 1 line at a time
545  count[2] = num_pixels;
546  status = nc_get_vara_float(observationGrp, Lt_blueId, start, count, Lt_blue[0]);
547  check_err(status, __LINE__, __FILE__);
548 
549  count[0] = expected_num_red_bands;
550  status = nc_get_vara_float(observationGrp, Lt_redId, start, count, Lt_red[0]);
551  check_err(status, __LINE__, __FILE__);
552 
553  count[0] = expected_num_SWIR_bands;
554  status = nc_get_vara_float(observationGrp, Lt_SWIRId, start, count, Lt_SWIR[0]);
555  check_err(status, __LINE__, __FILE__);
556 
557  status = nc_get_vara_ubyte(observationGrp, qual_SWIRId, start, count, qual_SWIR[0]);
558  check_err(status, __LINE__, __FILE__);
559 
560  for(int ip=0; ip<num_pixels; ip++) {
561  int band;
562  int ib = 0;
563  int ipb = ip * tot_num_bands;
564 
565  // load up the blue bands skip last 3
566  for(band=0; band<expected_num_blue_bands-3; band++) {
567  if(Lt_blue[band][ip] == Lt_blueFillValue) {
568  l1rec->Lt[ipb] = BAD_FLT;
569  } else {
570  l1rec->Lt[ipb] = Lt_blue[band][ip];
571  if (use_rhot) {
572  //l1rec->Lt[ipb] *= Fobar[ib] * l1rec->fsol * cos(l1rec->solz[ip]/RADEG) / M_PI;
573  l1rec->Lt[ipb] *= blue_solar_irradiance[band] * cos(l1rec->solz[ip]/RADEG) / earth_sun_distance_correction / M_PI / 10.0;
574  } else {
575  l1rec->Lt[ipb] /= 10.; // the input is in W/m2 ...
576  }
577  }
578  ib++;
579  ipb++;
580  }
581 
582  // load up the red bands use all of them
583  for(band=0; band<expected_num_red_bands; band++) {
584  if(Lt_red[band][ip] == Lt_redFillValue) {
585  l1rec->Lt[ipb] = BAD_FLT;
586  } else {
587  l1rec->Lt[ipb] = Lt_red[band][ip];
588  if (use_rhot) {
589  //l1rec->Lt[ipb] *= Fobar[ib] * l1rec->fsol * cos(l1rec->solz[ip]/RADEG) / M_PI;
590  l1rec->Lt[ipb] *= red_solar_irradiance[band] * cos(l1rec->solz[ip]/RADEG) / earth_sun_distance_correction / M_PI / 10.0;
591  } else {
592  l1rec->Lt[ipb] /= 10.; // the input is in W/m2 ...
593  }
594  }
595  ib++;
596  ipb++;
597  }
598 
599  // load up the SWIR bands choose low gain if either high gain band saturates
600  // band 2 low gain
601  // band 3 high gain
602  // band 5 low gain
603  // band 6 high gain
604 
605  int saturated = 0;
606  if((qual_SWIR[3][ip] & SATURATION_BIT) || (qual_SWIR[6][ip] & SATURATION_BIT)) {
607  saturated = 1;
608  l1rec->hilt[ip] = 1;
609  }
610 
611  for(band=0; band<expected_num_SWIR_bands; band++) {
612 
613  if((band==2 || band==5) && !saturated)
614  continue;
615  else if((band==3 || band==6) && saturated)
616  continue;
617 
618  if(Lt_SWIR[band][ip] == Lt_SWIRFillValue) {
619  l1rec->Lt[ipb] = BAD_FLT;
620  } else {
621  l1rec->Lt[ipb] = Lt_SWIR[band][ip];
622  if (use_rhot) {
623  //l1rec->Lt[ipb] *= Fobar[ib] * l1rec->fsol * cos(l1rec->solz[ip]/RADEG) / M_PI;
624  l1rec->Lt[ipb] *= SWIR_solar_irradiance[band] * cos(l1rec->solz[ip]/RADEG) / earth_sun_distance_correction / M_PI / 10.0;
625  } else {
626  l1rec->Lt[ipb] /= 10.; // the input is in W/m2 ...
627  }
628  }
629  ib++;
630  ipb++;
631  }
632  }
633 
634  return (LIFE_IS_GOOD);
635 }
636 
637 
643 int closel1_oci(filehandle *file) {
644  int status;
645 
646  printf("Closing oci l1b file\n");
647  status = nc_close(file->sd_id);
648  check_err(status, __LINE__, __FILE__);
649 
650  // Free memory
651  // From openl1_oci
652  if (tmpShort) free(tmpShort);
653  if (scan_time) free(scan_time);
654  if (hamside) free(hamside);
655  if (tilt) free(tilt);
656  if (scanQual) free(scanQual);
657  if (Lt_blue) free2d_float(Lt_blue);
658  if (Lt_red) free2d_float(Lt_red);
659  if (Lt_SWIR) free2d_float(Lt_SWIR);
660  if (qual_SWIR) free2d_uchar(qual_SWIR);
661  if (blue_solar_irradiance) free(blue_solar_irradiance);
662  if (red_solar_irradiance) free(red_solar_irradiance);
663  if (SWIR_solar_irradiance) free(SWIR_solar_irradiance);
664 
665  return (LIFE_IS_GOOD);
666 }
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 NULL
Definition: decode_rs.h:63
read l1rec
PARAM_TYPE_NONE Default value No parameter is buried in the product name name_prefix is case insensitive string compared to the product name PARAM_TYPE_VIS_WAVE The visible wavelength bands from the sensor are buried in the product name The product name is compared by appending and name_suffix ie aph_412_giop where prod_ix will be set to PARAM_TYPE_IR_WAVE same search method as PARAM_TYPE_VIS_WAVE except only wavelength above are looped through but prod_ix is still based ie aph_2_giop for the second band
int openl1_oci(filehandle *file)
Definition: l1_oci.c:99
void unix2ymds(double usec, int16_t *year, int16_t *mon, int16_t *day, double *secs)
Definition: unix2ymds.c:8
int32 * msec
Definition: l1_czcs_hdf.c:31
int readl1_oci(filehandle *file, int32_t line, l1str *l1rec, int lonlat)
Definition: l1_oci.c:434
double esdist_(int32_t *year, int32_t *day, int32_t *msec)
int syear
Definition: l1_czcs_hdf.c:15
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 M_PI
Definition: dtranbrdf.cpp:19
void free2d_uchar(unsigned char **p)
Free a two-dimensional array created by allocate2d_uchar.
Definition: allocate2d.c:52
#define LIFE_IS_GOOD
Definition: passthebuck.h:4
int sday
Definition: l1_czcs_hdf.c:15
int time_str(short, short, int, char *)
Definition: time_str.c:3
int closel1_oci(filehandle *file)
Definition: l1_oci.c:643
subroutine lonlat(alon, alat, xlon, ylat)
Definition: lonlat.f:2
void free2d_float(float **p)
Free a two-dimensional array created by allocate2d_float.
Definition: allocate2d.c:140
int want_verbose
#define SATURATION_BIT
Definition: l1_oci.c:25
void unix2yds(double usec, short *year, short *day, double *secs)
unsigned char ** allocate2d_uchar(size_t h, size_t w)
Allocate a two-dimensional array of type unsigned char of a given size.
Definition: allocate2d.c:35
#define RADEG
Definition: czcs_ctl_pt.c:5
Utility functions for allocating and freeing two-dimensional arrays of various types.
#define BAD_FLT
Definition: jplaeriallib.h:19
real *8 function esdist(iyr, iday, msec)
Definition: esdist.f:3
float ** allocate2d_float(size_t h, size_t w)
Allocate a two-dimensional array of type float of a given size.
Definition: allocate2d.c:123
double ymds2unix(short year, short month, short day, double secs)
int16_t * tilt
Definition: l2bin.cpp:80
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")
double isodate2unix(const char *isodate)
Definition: unix2isodate.c:61
int count
Definition: decode_rs.h:79