OB.DAAC Logo
NASA Logo
Ocean Color Science Software

ocssw V2022
main_l1det2det.c
Go to the documentation of this file.
1 /* =====================================================================*/
2 /* */
3 /* Program: l1det2det - multi-sensor detector-to-detector calibration */
4 /* */
5 /* */
6 /* Written By: */
7 /* */
8 /* Ewa J. Kwiatkowska */
9 /* SAIC */
10 /* NASA/Ocean Color Project */
11 /* June 2007 */
12 /* */
13 /* Programmer Organization Date Description of change*/
14 /* -------------- ------------ -------- ---------------------*/
15 /* Joel Gales Futuretech 08/01/13 Add support for */
16 /* NETCDF4 */
17 /* B. Franz NASA 09/10/13 Remove references to */
18 /* resolution and 1640 */
19 /* and 500m MODIS code */
20 /* Sean Bailey Futuretech 04/15/13 Entirely reworked to */
21 /* handle any mission, */
22 /* work with existing */
23 /* product definitions */
24 /* and require complete */
25 /* scans */
26 /* =====================================================================*/
27 
28 #include <stdio.h>
29 #include <sys/stat.h>
30 #include <stdlib.h>
31 #include <fcntl.h>
32 #include <string.h>
33 #include <errno.h>
34 #include <unistd.h>
35 #include <libgen.h>
36 #include <setupflags.h>
37 #include "l12_proto.h"
38 #include "calfile_utils.h"
39 #include <hdf5.h>
40 
41 int read9km_mask(char *file, char *mask);
42 int is_masked(int32_t bin, char *mask, int32_t nrows);
43 int16_t fileID = 0;
44 
45 float ***aots, ***nlws; //arrays to contain the values needed to avg over for inversion
46 
47 /* -------------------------------------------------------------------- */
48 /* main */
49 
50 /* -------------------------------------------------------------------- */
51 int main(int argc, char* argv[]) {
52 
53  static int32_t nrows = 2160;
54  long iscan = 0; /* input scan number */
55  long npix = 0; /* input number pixels per scan */
56  long spix = 0; /* start pixel for subscene process */
57  long epix = -1; /* end pixel for subscene process */
58  long sscan = 0; /* start scan for subscene process */
59  long escan = -1; /* end scan for subscene process */
60  long nscan;
61  long nbands;
62  int32_t ndets;
63  long nruns; /* number of output detector runs */
64  long ipix;
65  long jscan; /* loop variables */
66  int band;
67  long scanidx;
68  long runidx;
69  int initial_mirror_side;
70  int goodscan;
71  idDS ds_id;
72 
73  l1str l1rec; /* generic level-1b scan structure */
74  l2str l2rec; /* generic level-2 scan structure */
75  tgstr tgrec; /* structure to store target values */
76  aestr aerec; /* structure to store aerosol values */
77  calstr** detrecArray; /* detrec for each pixel in a scan */
78 
79  filehandle l1file; /* input l1 file handle */
80 
81  double start_time;
82  long i, j, ib, runCounter, nir_l, nprods;
83  int32 *numbin, *basebin;
84  int32_t bin, ibin, bin_row, bin_col;
85  long **valid_pixel_array; //Valid pixel array 0=invalid, 1=pixel valid
86  int numContiguousDetectorRuns; // number of scans for detector runs
87  int status;
88 
89  char buf[5000], *mask, *maskfile, maskname[1000];
90  float lat, lon, latbin;
91  int32_t flagusemask;
92  uint32 required;
93 
94  if (argc == 1) {
95  l2gen_usage("l1det2det");
96  return 0;
97  }
98 
99  want_verbose = 1;
100 
101  setvbuf(stdout, NULL, _IOLBF, 0);
102  setvbuf(stderr, NULL, _IOLBF, 0);
103 
104  /* Initialize file handles */
105  cdata_();
108 
109  /* Parse input parameters */
110  if (msl12_input(argc, argv, "l1det2det", &l1file) != 0) {
111  printf("-E- %s: Error parsing input parameters.\n", argv[0]);
112  exit(FATAL_ERROR);
113  }
114 
115  if (access(input->ifile[0], F_OK) || access(input->ifile[0], R_OK)) {
116  printf("-E- %s: Input file '%s' does not exist or cannot open.\n",
117  argv[0], input->ifile[0]);
118  exit(FATAL_ERROR);
119  }
120 
121  /* */
122  /* Open input file and get sensor */
123  /* and scan information from handle */
124  /* */
125  strcpy(l1file.name, input->ifile[0]);
126  l1file.geofile = input->geofile;
127 
128  if (openl1(&l1file) != 0) {
129  printf("-E- %s: Error opening %s for reading.\n",
130  argv[0], l1file.name);
131  exit(FATAL_ERROR);
132  }
133 
134  /*
135  * Check output file for matching sensorID
136  */
137  int existingSensorID;
138 
139  /* Save old HDF5 error handler */
140  H5E_auto_t old_func;
141  void *old_client_data;
142  H5Eget_auto(H5E_DEFAULT, &old_func, &old_client_data);
143 
144  /* Turn off error handling */
145  H5Eset_auto(H5E_DEFAULT, NULL, NULL);
146 
147  if (Hishdf(input->ofile[0]) == TRUE || H5Fis_hdf5(input->ofile[0]) == TRUE) {
148  ds_id = openDS(input->ofile[0]);
149  status = readAttr(ds_id, "sensorID", &existingSensorID);
150  if (status) {
151  printf("-E- %s: Problem reading output file: %s\n", argv[0], input->ofile[0]);
152  exit(EXIT_FAILURE);
153  }
154  endDS(ds_id);
155  if (l1file.sensorID != existingSensorID) {
156  printf("-E- %s: Mixing L2 data from different sensors, %s and %s.\n",
157  argv[0], sensorId2SensorName(l1file.sensorID), sensorId2SensorName(existingSensorID));
158  exit(EXIT_FAILURE);
159  }
160 
161  }
162  /* Restore previous HDF5 error handler */
163  H5Eset_auto(H5E_DEFAULT, old_func, old_client_data);
164 
165  npix = l1file.npix;
166  numContiguousDetectorRuns = input->ybox;
167  if (numContiguousDetectorRuns < 1)
168  numContiguousDetectorRuns = 1;
169 
170  if (numContiguousDetectorRuns > 2) {
171  printf("ybox > 2, sorry, only letting you use 2 - talk to management if you want more...\n");
172  numContiguousDetectorRuns = 2;
173  }
174  /* */
175  /* Allocate memory for L1 and L2 scan data */
176  /* */
177  if (alloc_l1(&l1file, &l1rec) == 0) {
178  printf("-E- %s: Unable to allocate L1 record.\n", argv[0]);
179  exit(FATAL_ERROR);
180  }
181  if (alloc_l2(&l1rec, &l2rec) == 0) {
182  printf("-E- %s: Unable to allocate L2 record.\n", argv[0]);
183  exit(FATAL_ERROR);
184  }
185 
186  /* Set the end pixel if it was not set by command argument */
187  if (l1_input->epixl == -1 || l1_input->epixl > l1file.npix)
188  l1_input->epixl = l1file.npix;
189  if (l1_input->eline == -1 || l1_input->eline > l1file.nscan)
190  l1_input->eline = l1file.nscan;
191  if (l1_input->spixl < 1)
192  l1_input->spixl = 1;
193  if (l1_input->sline < 1)
194  l1_input->sline = 1;
195 
196  spix = MAX(l1_input->spixl - 1, 0);
197  epix = MIN(l1_input->epixl - 1, l1file.npix - 1);
198  sscan = MAX(l1_input->sline - 1, 0);
199  escan = MIN(l1_input->eline - 1, l1file.nscan - 1);
200 
201  ndets = l1file.ndets;
202  // ensure there are full scans
203  if (sscan > 0) {
204  while (sscan % ndets)
205  sscan++;
206  }
207  int numContiguousDets = ndets * numContiguousDetectorRuns;
208  escan = ((escan + 1) / numContiguousDets) * numContiguousDets - 1;
209  //escan = ((escan + 1) / ndets) * ndets - 1;
210 
211  if (sscan > escan || spix > epix) {
212  printf("-E- %s: scan and pixel limits make no sense.\n", argv[0]);
213  printf(" start scan = %ld\n", sscan + 1);
214  printf(" end scan = %ld\n", escan + 1);
215  printf(" start pixel = %ld\n", spix + 1);
216  printf(" end pixel = %ld\n", epix + 1);
217  exit(FATAL_ERROR);
218  }
219 
220  /* Note: for the L1 file, npix is still the native scan pixel count */
221  l1file.spix = spix;
222  l1file.epix = epix;
223 
224  nscan = (escan - sscan) + 1;
225  npix = (epix - spix) + 1;
226  nbands = l1file.nbands;
227 
228  /* */
229  /* Transfer any additional header info to the L2 record headers */
230  /* */
231  l2rec.bindx = l1file.bindx;
232 
233  /* */
234  /* Set up parameters for strict ocean color forward processing */
235  /* */
236  input->proc_ocean = 1;
237  input->proc_land = 0;
238  input->proc_sst = 0;
239  input->atmocor = 1;
240  input->aer_opt = AERWANG;
241  l1_input->outband_opt = 0;
242  input->mode = FORWARD;
243  char outprod[L1_MAXPROD][32];
244  int nvars1d = 2;
245  char vars1Dnames[nvars1d][32]; // floating point one dimensional data to include in output - lon, lat
246  nprods = prodlist(l1file.sensorID, 0, (const char *) &(input->l2prod), (const char *) &(input->def_l2prod), outprod);
247 
248  /*
249  * Allocating memory for the multidimentional arrays of pixel index, nLw and aot
250  */
251  if ((valid_pixel_array = (long **) malloc(nscan * sizeof (long *))) == NULL) {
252  printf("-E- %s: Error allocating memory to the pixel index.\n", argv[0]);
253  exit(FATAL_ERROR);
254  }
255 
256  if ((aots = (float ***) malloc(nscan * sizeof (float *))) == NULL) {
257  printf("-E- : Error allocating memory to in situ inputs\n");
258  exit(FATAL_ERROR);
259  }
260  if ((nlws = (float ***) malloc(nscan * sizeof (float *))) == NULL) {
261  printf("-E- : Error allocating memory to in situ inputs\n");
262  exit(FATAL_ERROR);
263  }
264 
265  for (i = 0; i < nscan; i++) {
266  if ((valid_pixel_array[i] = (long *) calloc(npix, sizeof (long))) == NULL) {
267  printf("-E- %s: Error allocating memory to the pixel index.\n", argv[0]);
268  exit(FATAL_ERROR);
269  }
270 
271  if ((aots[i] = (float **) malloc(npix * sizeof (float *))) == NULL) {
272  printf("-E- %s: Error allocating memory to the AOT data.\n", argv[0]);
273  exit(FATAL_ERROR);
274  }
275  if ((nlws[i] = (float **) malloc(npix * sizeof (float *))) == NULL) {
276  printf("-E- %s: Error allocating memory to the nLw data.\n", argv[0]);
277  exit(FATAL_ERROR);
278  }
279  for (j = 0; j < npix; j++) {
280  if ((aots[i][j] = (float *) calloc(nbands, sizeof (float))) == NULL) {
281  printf("-E- %s: Error allocating memory to the aot data.\n", argv[0]);
282  exit(FATAL_ERROR);
283  }
284  if ((nlws[i][j] = (float *) calloc(nbands, sizeof (float))) == NULL) {
285  printf("-E- %s: Error allocating memory to the nLw data.\n", argv[0]);
286  exit(FATAL_ERROR);
287  }
288  }
289  }
290 
291  // Allocate memory for and set up deepwater mask
292  if ((mask = (char *) malloc(5940424 * sizeof (char))) == NULL) {
293  printf("-E- %s: Error allocating memory to the mask data\n", argv[0]);
294  exit(FATAL_ERROR);
295  }
296  if ((maskfile = getenv("OCDATAROOT")) == NULL) {
297  printf("-E- %s: OCDATAROOT directory undefined.\n", argv[0]);
298  exit(FATAL_ERROR);
299  }
300  strcpy(maskname, maskfile);
301  strcat(maskname, "/common/deep_water_mask_9km.dat\x0");
302 
303  if (read9km_mask(maskname, mask) > 0) {
304  printf("-E- %s: Failed reading the deep water mask file\n", argv[0]);
305  exit(FATAL_ERROR);
306  }
307 
308  /* ----------------- */
309  /* Set up bin arrays */
310  /* ----------------- */
311  if ((numbin = (int32 *) calloc(nrows, sizeof (int32))) == NULL) {
312  printf("-E- %s: Error allocating memory to numbin\n", argv[0]);
313  exit(FATAL_ERROR);
314  }
315  if ((basebin = (int32 *) calloc(nrows + 1, sizeof (int32))) == NULL) {
316  printf("-E- %s: Error allocating memory to basebin\n", argv[0]);
317  exit(FATAL_ERROR);
318  }
319 
320  bin = 0;
321  for (i = 0; i < nrows; i++) {
322  latbin = (i + 0.5) * (180.0 / nrows) - 90.0;
323  numbin[i] = (int32) (cos(latbin * PI / 180.0) * (2.0 * nrows) + 0.5);
324  bin += numbin[i];
325  }
326  basebin[0] = 1;
327  for (i = 1; i <= nrows; i++) {
328  basebin[i] = basebin[i - 1] + numbin[i - 1];
329  }
330 
331 
332  strcpy(buf, "\x0");
333  for (i = 0; i < L1_NFLAGS; i++) {
334  strcat(buf, l2_flag_lname[i]);
335  if (i < L1_NFLAGS - 1) strcat(buf, ",\x0");
336  }
337  setupflags(buf, input->flaguse, (uint32 *) & flagusemask, &required, &status);
338  if (status < 0) {
339  printf("-E- %s: Error reading L2 flags %s.\n", argv[0], input->ifile[0]);
340  exit(FATAL_ERROR);
341  }
342 
343  strcpy(vars1Dnames[0], "lon");
344  strcpy(vars1Dnames[1], "lat");
345 
346  nir_l = bindex_get(input->aer_wave_long);
347 
348  printf("\n\nBegin %s Processing\n", "l1det2det\x0");
349  printf("\nSensor is %s\n", sensorId2SensorName(l1file.sensorID));
350  printf("Sensor ID is %d\n", l1file.sensorID);
351  printf("Sensor has %d reflective bands\n", l1file.nbands);
352  printf("Sensor has %d emissive bands\n", l1file.nbandsir);
353  printf("Sensor resolution %d m\n", l1_input->resolution);
354  printf("Number of along-track detectors per band is %d\n", ndets);
355  printf("Number of input pixels per scan is %d\n", l1file.npix);
356  printf("Processing pixels %ld to %ld\n", spix + 1, epix + 1);
357  printf("Processing scans %ld to %ld\n", sscan + 1, escan + 1);
358 
359  start_time = now();
360  printf("\nBegin l1det2det processing at %s\n\n", ydhmsf(start_time, 'L'));
361 
362  /*
363  * Read file scan by scan, convert to L2
364  * determine valid detector "runs" and
365  * populate the nLw and tau arrays needed
366  * for averaging to insert into inverse
367  * processing.
368  */
369  for (iscan = sscan; iscan <= escan; iscan++) {
370 
371  scanidx = iscan - sscan;
372 
373  if (getl1rec(iscan, 1, &l1rec) != 0) {
374  exit(FATAL_ERROR);
375  }
376  // First Line!
377  if (iscan == sscan) {
378  initial_mirror_side = l1rec.mside;
379  }
380 
381  if ((iscan % 50) == 0)
382  printf("Processing scan #%6ld (%ld of %ld) after %6.0f seconds\n",
383  iscan, iscan - sscan + 1, escan - sscan + 1, now() - start_time);
384 
385  /* */
386  /* Convert the L1B radiances to L2 */
387  /* */
388  convl12(&l1rec, &l2rec, 0, l1rec.npix - 1, NULL);
389 
390  /* */
391  /* Loop through each pixel and check exclusion criteria */
392  /* on pixel flags and values */
393  /* */
394  for (ipix = 0; ipix < l1rec.npix; ipix++) {
395  ibin = 1;
396 
397  lat = l1rec.lat[ipix];
398  lon = l1rec.lon[ipix];
399  bin_row = (int32_t) ((90.0 + lat) * ((float64) nrows / (float64) 180.0));
400  bin_col = (int32_t) ((float64) numbin[bin_row] * (lon + 180.0) / (float64) 360.0);
401  bin = basebin[bin_row] + bin_col;
402  ibin = (int32_t) is_masked(bin, mask, nrows);
403 
404  if (((l1rec.flags[ipix] & flagusemask) == 0) &&
405  l2rec.chl[ipix] > 0.0 &&
406  l2rec.chl[ipix] <= input->chlthreshold &&
407  l2rec.taua[ipix * nbands + nir_l] > 0.0 &&
408  l2rec.taua[ipix * nbands + nir_l] <= input->aotthreshold && ibin > 0) {
409 
410  valid_pixel_array[scanidx][ipix] = 1;
411  for (ib = 0; ib < nbands; ib++) {
412  nlws[scanidx][ipix][ib] = l2rec.nLw [ipix * nbands + ib];
413  aots[scanidx][ipix][ib] = l2rec.taua[ipix * nbands + ib];
414  }
415  }
416  /*
417  * Once a detector run is complete,
418  * check the all lines in the scan to ensure all detectors are valid
419  */
420  if (l1rec.detnum == ndets - 1) {
421  if (numContiguousDetectorRuns == 1 || l1rec.mside != initial_mirror_side) {
422  for (i = 0; i < ndets * numContiguousDetectorRuns; i++) {
423  if (!valid_pixel_array[scanidx - i][ipix]) {
424  for (j = 0; j < ndets * numContiguousDetectorRuns; j++) {
425  valid_pixel_array[scanidx - j][ipix] = 0;
426  }
427  break;
428  }
429  }
430 
431  }
432  }
433  }
434  }
435 
436  printf("\nEnd of line-by-line processing at %s\n", ydhmsf(now(), 'L'));
437  printf("Processing Rate = %f scans/sec\n", nscan / (now() - start_time));
438 
439  // Free memory associated with deepwater mask
440  free(mask);
441  free(numbin);
442  free(basebin);
443 
444  /*
445  * Count the number of valid detector runs (i.e. scans with all detectors valid)
446  */
447  nruns = 0;
448  for (iscan = sscan; iscan <= escan; iscan += ndets) {
449  scanidx = iscan - sscan;
450  for (ipix = 0; ipix < l1rec.npix; ipix++) {
451  if (valid_pixel_array[scanidx][ipix])
452  nruns++;
453  }
454  }
455 
456  if (!nruns) {
457  printf("-W- %s: No analysis pixel data found in the L2 file %s\n", argv[0], input->ifile[0]);
458  } else {
459 
460  printf("\n\nNow extracting valid pixel runs\n");
461  printf("Along-track pixel runs of %d found: %ld\n", numContiguousDets, nruns);
462 
463  start_time = now();
464 
465  /*
466  * Allocate detrec, aer and tgt structures
467  */
468  if ((detrecArray = (calstr **) malloc(l1rec.npix * sizeof (calstr *))) == NULL) {
469  printf("-E- : Error allocating memory for detector record structures\n");
470  exit(EXIT_FAILURE);
471  }
472  for (i = 0; i < l1rec.npix; i++) {
473  detrecArray[i] = alloc_calrec(ndets, nbands, nprods, nvars1d);
474  }
475 
476 
477  alloc_aer(l1rec.npix, nbands, &aerec);
478  alloc_target(l1rec.npix, nbands, &tgrec);
479  /*
480  * set up the bits for inversion
481  */
482  input->aer_opt = FIXAOT;
483  input->mode = INVERSE_NLW;
484  aerec.mode = ON;
485  tgrec.mode = ON;
486  l2rec.tgrec = &tgrec;
487 
488  /*
489  * Open the output file. It will be created if necessary
490  */
491 
492  runCounter = 0;
493  ds_id = calfile_open(input->ofile[0], l1file.sensorID, nruns, ndets, nprods, nvars1d, outprod, vars1Dnames, &runCounter, DET2DET);
494 
495  /*
496  * Fill up detrec with detector info and grab avg aot/nlws
497  * for use in inversion to fill in vLt array in detrec
498  * along with the other
499  * detector values...
500  */
501 
502  runCounter--; //nruns index
503  int runNum = 0; //just for nice printing on output...
504 
505  for (iscan = sscan; iscan <= escan; iscan += ndets * numContiguousDetectorRuns) {
506  goodscan = 0;
507  runidx = iscan - sscan;
508  for (ipix = 0; ipix < l1rec.npix; ipix++) {
509  if (valid_pixel_array[runidx][ipix]) {
510  goodscan = 1;
511  break;
512  }
513  }
514  if (goodscan) {
515  // For "good scans", loop through detectors for the valid runs and get the pixel info for detrec
516  for (jscan = iscan; jscan < iscan + ndets * numContiguousDetectorRuns; jscan++) {
517  scanidx = jscan - sscan;
518  int detidx = jscan - iscan;
519  if (detidx >= ndets)
520  detidx -= ndets;
521 
522  if (getl1rec(jscan, 1, &l1rec) != 0) {
523  exit(EXIT_FAILURE);
524  }
525  // For the first line in the scan, fill up the per run values in detrec
526  if (detidx == 0) {
527  int16_t year, day;
528  double sec;
529  unix2yds(l1rec.scantime, &year, &day, &sec);
530 
531  for (ipix = 0; ipix < l1rec.npix; ipix++) {
532  if (valid_pixel_array[scanidx][ipix]) {
533  detrecArray[ipix]->sensorID = l1rec.l1file->sensorID;
534  detrecArray[ipix]->year = (int32_t) year;
535  detrecArray[ipix]->day = (int32_t) day;
536  detrecArray[ipix]->msec = (int32_t) (sec * 1.e3);
537  detrecArray[ipix]->mside = l1rec.mside;
538  detrecArray[ipix]->iscan = (int16_t) jscan;
539  detrecArray[ipix]->pixnum = l1rec.pixnum[ipix];
540  detrecArray[ipix]->vars1D[0] = l1rec.lon[ipix];
541  detrecArray[ipix]->vars1D[1] = l1rec.lat[ipix];
542  // average the nLws and aots over the pixels in ndets*numContiguousDetectorRuns
543  // and fills in the aerec and tgrec arrays needed for the inversion step
544  inversion_init(ndets*numContiguousDetectorRuns, runidx, nbands, ipix, &aerec, &tgrec);
545  }
546  }
547  }
548  // The Lt array is filled and the avg nLw and taua needed for the inversion
549  // are set up in this loop
550  for (ipix = 0; ipix < l1rec.npix; ipix++) {
551  if (valid_pixel_array[scanidx][ipix]) {
552  // average the nLws and aots over the pixels in ndets*numContiguousDetectorRuns
553  // and fills in the aerec and tgrec arrays needed for the inversion step
554  // inversion_init(ndets*numContiguousDetectorRuns, runidx, l1rec.nbands, ipix, &aerec, &tgrec);
555  for (band = 0; band < nbands; band++) {
556  ib = ipix * nbands + band;
557  detrecArray[ipix]->Lt[band][detidx] = (float) l1rec.Lt[ib];
558  }
559  }
560  }
561  /* */
562  /* Convert the L1B radiances to L2 */
563  /* */
564  convl12(&l1rec, &l2rec, 0, l1rec.npix - 1, &aerec);
565  /* */
566  /* Convert back to L1B to get vLts */
567  /* */
568  convl21(&l2rec, &tgrec, 0, l1rec.npix - 1, l1rec.Lt, NULL);
569 
570  // set detrec data and vLts obtained via inversion step
571  for (ipix = 0; ipix < l1rec.npix; ipix++) {
572  if (valid_pixel_array[scanidx][ipix]) {
573  // The vLts are filled here
574  for (band = 0; band < nbands; band++) {
575  ib = ipix * nbands + band;
576  detrecArray[ipix]->vLt[band][detidx] = (float) l1rec.Lt[ib];
577  }
578  // The detrec entries for the other request products are populated here
579  for (i = 0; i < nprods; i++) {
580  if (strcmp(outprod[i], "pixnum") != 0 &&
581  strcmp(outprod[i], "mside") != 0 &&
582  strncmp(outprod[i], "Lt", 2) != 0 &&
583  strncmp(outprod[i], "vLt", 3) != 0 &&
584  strcmp(outprod[i], "detnum") != 0 &&
585  strcmp(outprod[i], "l2_flags") != 0) {
586  // get the product index record
587  l2prodstr *p;
588 
589  if ((p = get_l2prod_index(outprod[i], l1rec.l1file->sensorID,
590  nbands + l1file.nbandsir, l1rec.npix, l1rec.l1file->nscan, l1rec.l1file->iwave)) == NULL) {
591  fprintf(stderr,
592  "-E- %s line %d: product index failure.\n",
593  __FILE__, __LINE__);
594  return (1);
595  };
596 
597  // extract the product
598  float *pbuf;
599  pbuf = prodgen(p, &l2rec);
600 
601  // fill up detrec data
602  detrecArray[ipix]->data[i][detidx] = pbuf[ipix];
603  }
604  }
605  }
606  }
607 
608  // Write out the valid runs for this scan
609  if (detidx == ndets - 1) {
610  for (ipix = 0; ipix < l1rec.npix; ipix++) {
611  if (valid_pixel_array[scanidx][ipix]) {
612  runCounter++;
613  runNum++;
614  if ((runNum % 500) == 0)
615  printf("Writing run %d of %ld\n", runNum, nruns);
616  calfile_write(ds_id, detrecArray[ipix], runCounter, nruns, ndets, nprods, nbands, nvars1d, outprod, vars1Dnames, DET2DET);
617  }
618  }
619  }
620  }
621  }
622  }
623  printf("\nEnd pixel-run extraction at %s, pixel-run number found: %ld\n\n\n", ydhmsf(now(), 'L'), nruns);
624  if (calfile_close(ds_id) != 0) {
625  printf("Trouble closing file %s\n", input->ofile[0]);
626  exit(EXIT_FAILURE);
627  }
628  // Free up some memory
629  for (i = 0; i < l1rec.npix; i++) {
630  free_calrec(detrecArray[i], nbands, nprods);
631  }
632 
633  }
634 
635  for (ib = 0; ib < nbands; ib++) {
636  free(aots[ib]);
637  free(nlws[ib]);
638  }
639  free(aots);
640  free(nlws);
641 
642  /* */
643  /* Close all files */
644  /* */
645  closel1(&l1file);
646 
647  if (nruns > 0) return (EXIT_SUCCESS);
648  else return (NOMATCH_ERROR);
649 }
char * ydhmsf(double dtime, char zone)
Definition: ydhmsf.c:12
#define L1_MAXPROD
Definition: filehandle.h:20
int32 l1file(int32 sdfid, int32 *nsamp, int32 *nscans, int16 *dtynum)
Definition: l1stat_chk.c:586
#define MAX(A, B)
Definition: swl0_utils.h:26
#define MIN(x, y)
Definition: rice.h:169
#define EXIT_SUCCESS
Definition: GEO_basic.h:72
int j
Definition: decode_rs.h:73
#define INVERSE_NLW
Definition: filehandle.h:13
int32_t day
int msl12_input(int argc, char *argv[], const char *progName, filehandle *l1file)
Definition: msl12_input.c:4201
int status
Definition: l1_czcs_hdf.c:32
float *** nlws
idDS openDS(const char *filename)
Definition: wrapper.c:616
int calfile_write(idDS ds_id, calstr *calrec, int recnum, int ydim, int xdim, int nprods, int nbands, int nvars1d, char l2prods[L1_MAXPROD][32], char vars1Dnames[L1_MAXPROD][32], caltype ctype)
int alloc_aer(int32_t npix, int32_t nbands, aestr *rec)
Definition: alloc_aer.c:12
#define NOMATCH_ERROR
Definition: l12_parms.h:82
#define NULL
Definition: decode_rs.h:63
void filehandle_init(filehandle *file)
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
#define TRUE
Definition: rice.h:165
float * lat
#define ON
Definition: l1.h:43
void msl12_input_init()
Definition: msl12_input.c:485
int32 nscan
Definition: l1_czcs_hdf.c:19
int l2gen_usage(const char *prog)
Definition: msl12_input.c:4254
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
instr * input
#define FORWARD
Definition: regen_attr.h:12
int readAttr(idDS ds_id, const char *nam, void *data)
#define PI
Definition: l3_get_org.c:6
int alloc_target(int32_t npix, int32_t nbands, tgstr *rec)
Definition: alloc_target.c:12
int32_t prodlist(int32_t sensorID, int32_t evalmask, const char *inprod, const char *defprod, char outprod[L1_MAXPROD][32])
int alloc_l2(l1str *l1rec, l2str *l2rec)
Definition: alloc_l2.c:17
int bindex_get(int32_t wave)
Definition: windex.c:45
#define AERWANG
Definition: l12_parms.h:23
void setupflags(char *flagdef, char *flaguse, uint32_t *flagusemask, uint32_t *required, int *status)
Definition: setupflags.c:5
int convl12(l1str *l1rec, l2str *l2rec, int32_t spix, int32_t epix, aestr *aerec)
Definition: convl12.c:16
void free_calrec(calstr *calrec, int nbands, int nprods)
void cdata_()
l1_input_t * l1_input
Definition: l1_options.c:9
#define FATAL_ERROR
Definition: swl0_parms.h:5
int32 nrows
INT32 getl1rec(INT16 sceneFrameNum, swl0scene *scene, swl0ctl *l0ctl, input_sType navinp[], navblk_sType navblk[], tilt_states_sType *tiltblk, swl1rec l1rec[])
Definition: getl1rec.c:45
int16_t fileID
int want_verbose
a context in which it is NOT documented to do so subscript which cannot be easily calculated when extracting TONS attitude data from the Terra L0 files Corrected several defects in extraction of entrained ephemeris and and as HDF file for both the L1A and Geolocation enabling retrieval of South Polar DEM data Resolved Bug by changing to opent the geolocation file only after a successful read of the L1A and also by checking for fatal errors from not restoring C5 and to report how many of those high resolution values were water in the new WaterPresent SDS Added valid_range attribute to Land SeaMask Changed to bilinearly interpolate the geoid_height to remove artifacts at one degree lines Made corrections to const qualification of pointers allowed by new version of M API library Removed casts that are no longer for same not the geoid Corrected off by one error in calculation of high resolution offsets Corrected parsing of maneuver list configuration parameter Corrected to set Height SDS to fill values when geolocation when for elevation and land water mask
Definition: HISTORY.txt:114
void unix2yds(double usec, short *year, short *day, double *secs)
@ DET2DET
Definition: calfile_utils.h:21
void closel1(filehandle *l1file)
Definition: l1_io.c:68
int main(int argc, char *argv[])
calstr * alloc_calrec(int ydim, int nbands, int nprods, int nvar1d)
void inversion_init(long ndets, long iscan, int nbands, long ipix, aestr *aerec, tgstr *tgrec)
int32_t nbands
const char * sensorId2SensorName(int sensorId)
Definition: sensorInfo.c:198
int32 spix
Definition: l1_czcs_hdf.c:21
int32_t iscan
void * prodgen(l2prodstr *p, l2str *l2rec)
Definition: prodgen.c:89
u5 which has been done in the LOCALGRANULEID metadata should have an extension NRT It is requested to identify the NRT production Changes from v6 which may affect scientific the sector rotation may actually occur during one of the scans earlier than the one where it is first reported As a the b1 values are about the LOCALGRANULEID metadata should have an extension NRT It is requested to identify the NRT to fill pixels affected by dead subframes with a special value Output the metadata of noisy and dead subframe Dead Subframe EV and Detector Quality Flag2 Removed the function call of Fill_Dead_Detector_SI to stop interpolating SI values for dead but also for all downstream products for science test only Changes from v5 which will affect scientific to conform to MODIS requirements Removed the Mixed option from the ScanType in the code because the L1A Scan Type is never Mixed Changed for ANSI C compliance and comments to better document the fact that when the HDF_EOS metadata is stricly the and products are off by and in the track respectively Corrected some misspelling of RCS swir_oob_sending_detector to the Reflective LUTs to enable the SWIR OOB correction detector so that if any of the sending detectors becomes noisy or non near by good detectors from the same sending band can be specified as the substitute in the new look up table Code change for adding an additional dimension of mirror side to the Band_21_b1 LUT to separate the coefficient of the two mirror sides for just like other thermal emissive so that the L1B code can calibrate Band scan to scan with mirror side dependency which leads better calibration result Changes which do not affect scientific when the EV data are not provided in this Crosstalk Correction will not be performed to the Band calibration data Changes which do not affect scientific and BB_500m in L1A Logic was added to turn off the or to spatial aggregation processes and the EV_250m_Aggr1km_RefSB and EV_500m_Aggr1km_RefSB fields were set to fill values when SDSs EV_250m and EV_500m are absent in L1A file Logic was added to skip the processing and turn off the output of the L1B QKM and HKM EV data when EV_250m and EV_500m are absent from L1A In this the new process avoids accessing and reading the and L1A EV skips and writing to the L1B and EV omits reading and subsampling SDSs from geolocation file and writing them to the L1B and omits writing metadata to L1B and EV and skips closing the L1A and L1B EV and SDSs Logic was added to turn off the L1B OBC output when the high resolution OBC SDSs are absent from L1A This is accomplished by skipping the openning the writing of metadata and the closing of the L1B OBC hdf which is Bit in the scan by scan bit QA has been changed Until now
Definition: HISTORY.txt:361
This should be set to the NetCDF standard name if exists for this product Create a function that computes your product edit get_myprod c add prototype to l12_proto h add get_myprod c to add_executable for l2gen and l3gen in CMakeLists txt Add an entry to the output routine to call your function edit prodgen c edit function prodgen() case CAT_myprod pbuf
Definition: dfutils.h:28
float * lon
int read9km_mask(char *file, char *mask)
Definition: read9km_mask.c:25
int convl21(l2str *l2rec, tgstr *tgrec, int32_t spix, int32_t epix, float *vLt, vcstr *vrec)
Definition: convl21.c:15
int32 epix
Definition: l1_czcs_hdf.c:23
int32_t alloc_l1(filehandle *l1file, l1str *l1rec)
Definition: alloc_l1.c:15
int calfile_close(idDS ds_id)
int endDS(idDS ds_id)
Definition: wrapper.c:634
idDS calfile_open(char *ofile, int sensorID, int ydim, int xdim, int nprods, int nvars1d, char l2prods[L1_MAXPROD][32], char vars1Dnames[L1_MAXPROD][32], long *numExistingRecords, caltype ctype)
Definition: calfile_utils.c:48
float *** aots
int i
Definition: decode_rs.h:71
int openl1(filehandle *l1file)
Definition: l1_io.c:207
How many dimensions is the output array Default is Not sure if anything above will work correctly strcpy(l2prod->title, "no title yet")
#define FIXAOT
Definition: l12_parms.h:30
#define L1_NFLAGS
Definition: filehandle.h:21
int npix
Definition: get_cmp.c:27
int is_masked(int32_t bin, char *mask, int32_t nrows)
Definition: read9km_mask.c:65
float p[MODELMAX]
Definition: atrem_corl1.h:131
int32_t get_l2prod_index(const l2_prod &l2, const char *prodname)
Definition: l2bin.cpp:345