OB.DAAC Logo
NASA Logo
Ocean Color Science Software

ocssw V2022
mng_ms.cpp
Go to the documentation of this file.
1 /*
2  * mng_ms.cpp - take over the management of the large multi-scattering
3  * reflectance arrays for the cloud processing
4  * W Robinson, SAIC, 17 Dec 2021
5  */
6 #include <string.h>
7 #include "mfhdf.h"
8 #include "dim_mgr.hpp"
9 
11  {
12  int fid_ms; /* fid of MS table file */
13  int fid_ms_std; /* fid of MS std file (for files over water */
14  int sds_id; /* the SDS id refl table */
15  int sds_id_std; /* the SDS id refl SD table */
16  int32_t *dim_siz; /* array of dim sizes */
17  double **dim_vals; /* arrays of each dim value set */
18  };
19 typedef struct ms_finfo_struc_def ms_finfo_struc;
20 
22 // for the water surface
23 ms_finfo_struc wtrcld_wtrsfc_info[3], icecld_wtrsfc_info[3];
24 
25 dim_mgr land_mgr( 3, 83 ); /* manage separately for over water, land */
26 dim_mgr wtr_mgr( 3, 83 );
27 
28 // storage for the size of each reflectance / sd table size
29 // the waater has 18 radii and the ice has 12
30 int32_t wtrcld_nelts, icecld_nelts; // size of water, ice cloud arrays
31 
32 /*
33  * mng_ms_init - start the management of the large tables
34  */
35 /*******************************************************************
36  mng_ms_init_ initialize the management of the large cloud multi-scatering
37  reflectance arrays
38 
39  Returns:
40  int32_t - nothing at the moment
41 
42  Parameters: (in calling order)
43  Type Name I/O Description
44  ---- ---- --- -----------
45  char * fil_wtrcld_landsfc I file name for water cloud
46  land sfc arrays
47  char * fil_icecld_landsfc I file name for ice cloud
48  land sfc arrays
49  int len1 I length of 1st char *
50  int len2 I length of 2nd char *
51  The above is for 2 land surface files and their lengths
52  The rest are the files for the water surface and the lengths after
53  them all (a fortran convention). below, I'll just list the file names
54  and their lengths same types as above
55  file name length Description
56  ----------------------- ---------- -----------------------------
57  fil_wtrcld_landsfc len1 water cloud, land sfc
58  fil_icecld_landsfc len2 ice cloud, land sfc
59  fil_wtrcld_wtrsfc1 len_ww1 water cloud, water sfc wind1
60  fil_wtrcld_wtrsfc1_sd len_ww1sd water cloud, water sfc wind1 sd
61  fil_icecld_wtrsfc1 len_iw1 ice cloud, water sfc wind1
62  fil_icecld_wtrsfc1_sd len_iw1sd ice cloud, water sfc wind1 sd
63  fil_wtrcld_wtrsfc2 len_ww2 water cloud, water sfc wind2
64  fil_wtrcld_wtrsfc2_sd len_ww2sd water cloud, water sfc wind2 sd
65  fil_icecld_wtrsfc2 len_iw2 ice cloud, water sfc wind2
66  fil_icecld_wtrsfc2_sd len_iw2sd ice cloud, water sfc wind sd
67  fil_wtrcld_wtrsfc3 len_ww3 water cloud, water sfc wind3
68  fil_wtrcld_wtrsfc3_sd len_ww3sd water cloud, water sfc wind3 sd
69  fil_icecld_wtrsfc3 len_iw3 ice cloud, water sfc wind3
70  fil_icecld_wtrsfc3_sd len_iw3sd ice cloud, water sfc wind3 sd
71 
72 ********************************************************************/
73 extern"C" int32_t mng_ms_init_( char *fil_wtrcld_landsfc,
74  char *fil_icecld_landsfc,
75  char *fil_wtrcld_wtrsfc1, char *fil_wtrcld_wtrsfc1_sd,
76  char *fil_icecld_wtrsfc1, char *fil_icecld_wtrsfc1_sd,
77  char *fil_wtrcld_wtrsfc2, char *fil_wtrcld_wtrsfc2_sd,
78  char *fil_icecld_wtrsfc2, char *fil_icecld_wtrsfc2_sd,
79  char *fil_wtrcld_wtrsfc3, char *fil_wtrcld_wtrsfc3_sd,
80  char *fil_icecld_wtrsfc3, char *fil_icecld_wtrsfc3_sd,
81  int len1, int len2, int len_ww1, int len_ww1sd, int len_iw1, int len_iw1sd,
82  int len_ww2, int len_ww2sd, int len_iw2, int len_iw2sd,
83  int len_ww3, int len_ww3sd, int len_iw3, int len_iw3sd )
84  {
85  char *fn, *fn2;
86  char *ftn_str_cond( char *, int );
87  int32_t status;
88  int32_t set_ms_file( char *, char *, int32_t , ms_finfo_struc * );
89 
90  /* set up the water and ice cloud over land file table info */
91  fn = ftn_str_cond( fil_wtrcld_landsfc, len1 );
93  free( fn );
94 
95  fn = ftn_str_cond( fil_icecld_landsfc, len2 );
97  free( fn );
98 
99  /* make sure the first 5 dims, less radius, are equal
100  Also, set up the # elements for ice, water in the last 3 dims */
101  status = 0;
102  wtrcld_nelts = wtrcld_landsfc_info.dim_siz[3] + 1; // for tau
103  icecld_nelts = icecld_landsfc_info.dim_siz[3] + 1;
104  for( int idim = 0; idim < 6; idim++ )
105  {
106  if( ( idim != 5 ) && ( wtrcld_landsfc_info.dim_siz[idim] !=
107  icecld_landsfc_info.dim_siz[idim] ) )
108  {
109  status = 1;
110  printf( "%s, %d, E: Water and ice cloud over land dimension mismatch\n",
111  __FILE__, __LINE__ );
112  printf( " file_water: %s\n", fil_wtrcld_landsfc );
113  printf( " file_ice : %s\n", fil_icecld_landsfc );
114  }
115  // set up the size of the water and ice # elements,
116  // that will go in the data blob
117  if( idim > 3 ) // for band and radius
118  {
119  wtrcld_nelts *= wtrcld_landsfc_info.dim_siz[idim];
120  icecld_nelts *= icecld_landsfc_info.dim_siz[idim];
121  }
122  }
123  if( status != 0 ) exit(27);
124 
125 /*
126  * setup the water and ice cloud over water surface tables
127  * loop through all 3 water surface types (3, 7, 15 m/s), get the SDS IDs
128  * and check consostency of dims with over-land info
129  */
130  status = 0;
131  for( int32_t isfc = 0; isfc < 3; isfc++ )
132  {
133  if( isfc == 0 ) {
134  fn = ftn_str_cond( fil_wtrcld_wtrsfc1, len_ww1 );
135  fn2 = ftn_str_cond( fil_wtrcld_wtrsfc1_sd, len_ww1sd );
136  } else if( isfc == 1 ) {
137  fn = ftn_str_cond( fil_wtrcld_wtrsfc2, len_ww2 );
138  fn2 = ftn_str_cond( fil_wtrcld_wtrsfc2_sd, len_ww2sd );
139  } else if( isfc == 2 ) {
140  fn = ftn_str_cond( fil_wtrcld_wtrsfc3, len_ww3 );
141  fn2 = ftn_str_cond( fil_wtrcld_wtrsfc3_sd, len_ww3sd );
142  }
143  set_ms_file( fn, fn2, 1, &wtrcld_wtrsfc_info[isfc] );
144  free( fn );
145  free( fn2 );
146 
147  if( isfc == 0 ) {
148  fn = ftn_str_cond( fil_icecld_wtrsfc1, len_iw1 );
149  fn2 = ftn_str_cond( fil_icecld_wtrsfc1_sd, len_iw1sd );
150  } else if( isfc == 1 ) {
151  fn = ftn_str_cond( fil_icecld_wtrsfc2, len_iw2 );
152  fn2 = ftn_str_cond( fil_icecld_wtrsfc2_sd, len_iw2sd );
153  } else if( isfc == 2 ) {
154  fn = ftn_str_cond( fil_icecld_wtrsfc3, len_iw3 );
155  fn2 = ftn_str_cond( fil_icecld_wtrsfc3_sd, len_iw3sd );
156  }
157  set_ms_file( fn, fn2, 1, &icecld_wtrsfc_info[isfc] );
158  free( fn );
159  free( fn2 );
160 
161  for( int idim = 0; idim < 6; idim++ )
162  {
163  if( ( ( idim == 3 ) && ( wtrcld_wtrsfc_info[isfc].dim_siz[idim] !=
164  wtrcld_landsfc_info.dim_siz[idim] + 1 ) )
165  || ( ( idim != 3 ) && (wtrcld_wtrsfc_info[isfc].dim_siz[idim] !=
166  wtrcld_landsfc_info.dim_siz[idim] ) ) )
167  {
168  status = 1;
169  printf(
170  "%s, %d, E: water cloud over-water MS file dimension mismatch\n",
171  __FILE__, __LINE__ );
172  printf( " dim # %d, water sfc # %d\n", idim, isfc );
173  }
174  if( status != 0 ) exit(27);
175 
176  if( ( ( idim == 3 ) && ( icecld_wtrsfc_info[isfc].dim_siz[idim] !=
177  icecld_landsfc_info.dim_siz[idim] + 1 ) )
178  || ( ( idim != 3 ) && (icecld_wtrsfc_info[isfc].dim_siz[idim] !=
179  icecld_landsfc_info.dim_siz[idim] ) ) )
180  {
181  status = 1;
182  printf(
183  "%s, %d, E: ice cloud over-water MS file dimension mismatch\n",
184  __FILE__, __LINE__ );
185  printf( " dim # %d, water sfc # %d\n", idim, isfc );
186  }
187  if( status != 0 ) exit(27);
188  }
189  }
190  /* set up the dimension manager for points over land */
191  //dim_mgr land_mgr( 3, 80 );
192  for( int idim = 0; idim < 3; idim++ )
193  {
194  land_mgr.init_dim( idim, wtrcld_landsfc_info.dim_siz[idim],
195  wtrcld_landsfc_info.dim_vals[idim] );
196  wtr_mgr.init_dim( idim, wtrcld_wtrsfc_info[0].dim_siz[idim],
197  wtrcld_wtrsfc_info[0].dim_vals[idim] );
198  }
199 
200  return 0;
201  }
202 
203 /********************************************************************
204  ftn_str_cond - Just condition a string to be null terminated
205 
206  Returns:
207  char * the good null terminated astring
208  Parameters: (in calling order)
209  Type Name I/O Description
210  ---- ---- --- -----------
211  char * str I initial fortran string
212  int len I string length
213 ********************************************************************/
214 extern"C" char *ftn_str_cond( char *str, int len )
215  {
216  char *fn;
217 
218  fn = (char *) malloc( ( len + 1 ) * sizeof( char ) );
219  strncpy( fn, str, len );
220  fn[ len ] = (char) 0;
221  return fn;
222  }
223 
224 extern"C" int32_t set_ms_file( char *fil_ms, char *fil_ms_std, int32_t fil_typ,
225  ms_finfo_struc *finfo )
226 /*******************************************************************
227 
228  set_ms_file - set up all the info for the multi-scattering file table
229 
230  Returns type: int - size of dimension or -1 if problem
231 
232  Parameters: (in calling order)
233  Type Name I/O Description
234  ---- ---- --- -----------
235  char * fil_ms I name of multi-scatering file
236  char * fil_ms_std I name of multi-scatering std
237  file
238  int32_t fil_typ I 0 for a land surface file
239  (has the MS and MSstd tables, 1 for the
240  water surface files (each file has the MS or MSstd)
241  ms_finfo_struc * finfo O all information of the file(s)
242 
243  Modification history:
244  Programmer Date Description of change
245  ---------- ---- ---------------------
246  W. Robinson 20 Dec 2017 original development
247 
248 *******************************************************************/
249  {
250  char bl_nam[500];
251  int rank, dim_lens[7], nattr, ret1, ret2, vtyp;
252  int start[7], end[7], sds_id;
253  int32_t trg_rank = 6;
254  float *fstore;
255 
256  char *dim_sds_names[] = { "ReflectanceSensorZenith", "ReflectanceSolarZenith",
257  "ReflectanceRelativeAzimuth", "OpticalThickness", "Wavelengths",
258  "ParticleRadius" };
259 
260  /*
261  * open the file and get the refl array SDS ID
262  */
263  finfo->fid_ms = SDstart( fil_ms, DFACC_READ );
264  finfo->sds_id = SDselect( finfo->fid_ms, SDnametoindex( finfo->fid_ms,
265  "MultiScatBDReflectance" ) );
266  SDgetinfo( finfo->sds_id, bl_nam, &rank, dim_lens, &vtyp, &nattr );
267  if( ( rank != trg_rank ) || ( finfo->sds_id == FAIL ) )
268  {
269  printf( "E: %s, %d: MS refl in cloud table file has wrong size\n",
270  __FILE__, __LINE__ );
271  printf( " file: %s", fil_ms );
272  printf( " SDS: MultiScatBDReflectance\n" );
273  exit(27);
274  }
275  /*
276  * for the wind-dependent tables, the SD is in separate file
277  */
278  if( fil_typ == 0 )
279  finfo->fid_ms_std = finfo->fid_ms;
280  else
281  finfo->fid_ms_std = SDstart( fil_ms_std, DFACC_READ );
282  /*
283  * get the Std Dev SDS ID
284  */
285  finfo->sds_id_std = SDselect( finfo->fid_ms_std,
286  SDnametoindex( finfo->fid_ms_std, "StdDevMultiScatBDReflectance" ) );
287  SDgetinfo( finfo->sds_id_std, bl_nam, &rank, dim_lens, &vtyp, &nattr );
288  if( ( rank != trg_rank ) || ( finfo->fid_ms_std == FAIL ) )
289  {
290  printf( "E: %s, %d: MS refl std in cloud table file has wrong size\n",
291  __FILE__, __LINE__ );
292  printf( " file: %s", fil_ms );
293  printf( " SDS: StdMultiScatBDReflectance\n" );
294  exit(27);
295  }
296  /* Note that the data in the file is float and we are putting it into
297  double, so some translation will be needed */
298  finfo->dim_siz = (int32_t *) malloc( trg_rank * sizeof( int32_t ) );
299  finfo->dim_vals = ( double **) malloc( trg_rank * sizeof( double * ) );
300 
301  /*
302  * get information for all 6 dimensions
303  */
304  for( int32_t isds = 0; isds < 6; isds++ )
305  {
306  if( ( sds_id = SDselect( finfo->fid_ms, SDnametoindex( finfo->fid_ms,
307  dim_sds_names[isds] ) ) ) == FAIL )
308  {
309  printf( "E: %s, %d: Can't read table SDS %s\n", __FILE__, __LINE__,
310  dim_sds_names[isds] );
311  exit(26);
312  }
313  SDgetinfo( sds_id, bl_nam, &rank, dim_lens, &vtyp, &nattr );
314  finfo->dim_siz[isds] = dim_lens[0];
315  finfo->dim_vals[isds] = (double *) malloc( dim_lens[0] * sizeof(double) );
316  fstore = (float *) malloc( dim_lens[0] * sizeof(float) );
317  start[0] = 0;
318  end[0] = dim_lens[0];
319  ret1 = SDreaddata( sds_id, start, NULL, end, (float *)fstore );
320  for( int32_t ipx = 0; ipx < end[0]; ipx++ )
321  finfo->dim_vals[isds][ipx] = *(fstore + ipx );
322  free( fstore );
323  ret2 = SDendaccess( sds_id );
324  if( ( ret1 !=0 ) || ( ret2 !=0 ) )
325  {
326  printf( "E: %s, %d: Can't read table SDS %s\n", __FILE__, __LINE__,
327  dim_sds_names[isds] );
328  exit(26);
329  }
330  }
331  return 0;
332  }
333 
334 /* the next is the retrieval of data side. It will have a part to get the
335  structure defining the set of data grid points and the weights to apply
336  And, a part to make the particular interpolated array requested
337 */
338 
339 extern"C" int mng_ms_get_lib_( float *sol_ang, float *sen_ang,
340  float *relaz_ang, int *sfc_typ, int *meas_typ, int *scan, float *wtr_int,
341  float *ice_int, int *stat )
342 /*******************************************************************
343 
344  mng_ms_get_lib - get the portion of the MS reflectance library
345  interpolated to the solar, sensor and relative azimuth angles
346 
347  Returns type: int - 0 good, otherwise an error
348  Note that in th ecase where the point is outside the table, a 2 status
349  is returned ans will be handled by making the final table all 0,
350  hopefully that will just make for fill values in cloud quantities
351  dealing with the calling code is horrendous.
352 
353  Parameters: (in calling order)
354  Type Name I/O Description
355  ---- ---- --- -----------
356  float * sol_ang I solar angle
357  float * sen_ang I sensor angle
358  float * relaz_ang I relative azimuth angle
359  int * sfc_typ I underlying surface type:
360  0 - land, 1 - water, 3 m/s
361  wind, 2 - water, 7 m/s,
362  3 - water, 15 m/s
363  int * meas_typ I measurement type to return:
364  0 - reflectance, 1 - Std Dev
365  int * scan I current scan line, use as
366  access id
367  float * wtr_int O interpolated water array
368  float * ice_int O interpolated ice array
369  int * stat O status - 0 is good, 5 is out
370  of table bounds
371 
372  Modification history:
373  Programmer Date Description of change
374  ---------- ---- ---------------------
375  W. Robinson 23 Dec 2021 original development
376 
377 NOTE - organization of data blobs
378 position
379  1 - reflectance wtrcld size wtrcld_nelts
380  2 - SD wtrcld size wtrcld_nelts
381  3 - reflectance icecld size icecld_nelts
382  2 - SD icecld size icecld_nelts
383  This repeats 2 more for the water surface arrays
384 
385 *******************************************************************/
386  {
387  /* follow the test program for dim_mgr */
388  /* */
389  int32_t access_id;
390  //static pt_info_struc *pt_info;
391  pt_info_struc *pt_info;
392  double pt[3];
393  int32_t status, ndat;
394  int32_t nsfc; // nsfc will be for the # of surface types (1 land, 3 water)
395  int32_t narr = 4; // # of arrays to get out from the wtrcld, wtrcld_sd,
396  // icecld, icecld_sd SDSes
397  int32_t ndim = 3, off_wtr, off_ice;
398  int32_t ntau, ntau1, nwav, nrad_wtr, nrad_ice, nrad;
399  static int32_t prune_called = 0;
400 
401  // ntau1 represents ntau + 1 the size in tau of the final storage
402  // the over-water data has ntau+1 while the over-land has ntau
403  // its just the (odd) way the cloud code, files were set up
404 
405  // WDR temp to report the points managed when the access_id changes
406  /* leave out normally */
407 //printf( "%s, %d, SCAN: %d\n", __FILE__, __LINE__, *scan );
408  //if( ( *scan % 20 ) == 19 )
409 /*
410  if( *scan == 300 )
411  {
412  double way = 0.;
413  printf( "%s, %d, LAND MANAGER DUMP\n", __FILE__, __LINE__ );
414  land_mgr.dump_mgr( &way );
415  //printf( "%s, %d, WATER MANAGER DUMP\n", __FILE__, __LINE__ );
416  // no need wtr_mgr.dump_mgr( &way );
417  exit( 0 );
418  }
419  */
420  // WDR end of the probe of points
421 
422  // WDR add the prune, every 200 lines prune all 100 less than current line
423  if( ( *scan % 200 ) == 199 )
424  {
425  if( prune_called == 0 )
426  {
427  //double way = 0.;
428  //printf( "LAND MANAGER DUMP PRE, scan %d\n", *scan );
429  //land_mgr.dump_mgr( &way );
430  //printf( "WATER MANAGER DUMP PRE, scan %d\n", *scan );
431  //wtr_mgr.dump_mgr( &way );
432  land_mgr.prune( *scan - 100 );
433  wtr_mgr.prune( *scan - 100 );
434  prune_called = 1;
435  //printf( "LAND MANAGER DUMP POST, scan %d\n", *scan );
436  //land_mgr.dump_mgr( &way );
437  }
438  }
439  else
440  prune_called = 0;
441 
442  access_id = *scan;
443  ndat = pow( 2, ndim );
444  // Note that the organization is senz, solz, relaz
445  pt[0] = *sen_ang;
446  pt[1] = *sol_ang;
447  pt[2] = *relaz_ang;
448 
449  // do over-land and over-water separately - no need to read land table
450  // over water or visa-versa
451  if( *sfc_typ == 0 ) // a request for over-land table data
452  {
453  pt_info = land_mgr.mng_pt( pt, access_id, &status );
454  nsfc = 1;
455  }
456  else
457  {
458  pt_info = wtr_mgr.mng_pt( pt, access_id, &status );
459  nsfc = 3;
460  }
461 
462  if( status != 0 )
463  {
464  // trying to deal with exceptions in the calling code is a mess
465  // so the error 2 will end up returning an unfilled table
466  if( status == 2 ) // out of table bounds
467  {
468  *stat = 5;
469  return 5;
470  }
471  printf( "%s, %d: mng_ms_get_lib error found of %d\n", __FILE__,
472  __LINE__, status );
473  exit( 27 );
474  }
475  else
476  {
477  ntau = wtrcld_landsfc_info.dim_siz[3];
478  ntau1 = ntau + 1;
479  nwav = wtrcld_landsfc_info.dim_siz[4];
480  nrad_wtr = wtrcld_landsfc_info.dim_siz[5];
481  nrad_ice = icecld_landsfc_info.dim_siz[5];
482 
483  if( pt_info->interval_needs_data == 1 )
484  {
485  /* we need to read data for some of the grid points in the interval */
486  int32_t n_new_blobs = 0;
487  int32_t offset[3], arr_off, tau_off;
488  int tsds_id;
489  int32_t istart[6] = { 0, 0, 0, 0, 0, 0 };
490  int32_t icount[6] = { 1, 1, 1, ntau, nwav, nrad_wtr };
491  if( *sfc_typ != 0 ) icount[3] = ntau1;
492 
493  for( int32_t idat = 0; idat < ndat; idat++ ) // loop interval corners
494  {
495  if( pt_info->pt_status[idat] == -1 )
496  {
497  n_new_blobs++;
498  /* do the read, reorg, store in struc and deposit in
499  a data blob pointer, use calloc to avoid needing to zero
500  the array for a particular tau index that changes */
501  float *dat_blob = (float *)calloc( ( wtrcld_nelts * 2 +
502  icecld_nelts * 2 ) * nsfc, sizeof(float) );
503  linear_to_offset( ndim, idat, offset );
504  for( int32_t idim = 0; idim < 3; idim++ )
505  istart[idim] = pt_info->pt_base_loc[idim] + offset[idim];
506 
507  // loop thru the nsfc surface types (3 for water sfc, 1 for land
508  for( int32_t isfc = 0; isfc < nsfc; isfc++ )
509  {
510  // loop thru the arrays to fill the data blob, do - 4 of them
511  // per surface: wtrcld ref, sd, icecld ref, sd (etc for more sfc)
512  for( int32_t iarr = 0; iarr < narr; iarr++ )
513  {
514  nrad = nrad_wtr;
515  tau_off = 0; /* for land sfc: is 0 to transfer data to tau
516  dim of blob starting at 0, 1 to start at
517  1 for the SD (to get expected arrays)
518  for wtr sfc: is always 0 */
519  if( *sfc_typ == 0 ) // the land sfc tables are 1 way
520  {
521  switch (iarr)
522  {
523  case 0: // water cloud refl
524  tsds_id = wtrcld_landsfc_info.sds_id;
525  arr_off = 0;
526  break;
527  case 1: // water cloud SD
528  tsds_id = wtrcld_landsfc_info.sds_id_std;
529  arr_off = wtrcld_nelts;
530  tau_off = 1;
531  break;
532  case 2: // ice cloud refl
533  tsds_id = icecld_landsfc_info.sds_id;
534  arr_off = 2 * wtrcld_nelts;
535  nrad = nrad_ice;
536  break;
537  case 3: // ice cloud SD
538  tsds_id = icecld_landsfc_info.sds_id_std;
539  arr_off = 2 * wtrcld_nelts + icecld_nelts;
540  nrad = nrad_ice;
541  tau_off = 1;
542  break;
543  }
544  }
545  else // the water surface tables are another way
546  {
547  tau_off = 0;
548  switch (iarr)
549  {
550  case 0: // water cloud refl
551  tsds_id = wtrcld_wtrsfc_info[isfc].sds_id;
552  arr_off = isfc * 2 * ( wtrcld_nelts + icecld_nelts );
553  break;
554  case 1: // water cloud SD
555  tsds_id = wtrcld_wtrsfc_info[isfc].sds_id_std;
556  arr_off = wtrcld_nelts +
557  isfc * 2 * ( wtrcld_nelts + icecld_nelts );
558  break;
559  case 2: // ice cloud refl
560  tsds_id = icecld_wtrsfc_info[isfc].sds_id;
561  arr_off = 2 * wtrcld_nelts +
562  isfc * 2 * ( wtrcld_nelts + icecld_nelts );
563  nrad = nrad_ice;
564  break;
565  case 3: // ice cloud SD
566  tsds_id = icecld_wtrsfc_info[isfc].sds_id_std;
567  arr_off = 2 * wtrcld_nelts + icecld_nelts +
568  isfc * 2 * ( wtrcld_nelts + icecld_nelts );
569  nrad = nrad_ice;
570  break;
571  }
572  }
573  icount[5] = nrad;
574  /* read from the file */
575  int32_t ntau_mod = ( *sfc_typ == 0 ) ? ntau : ntau1;
576  float *tmp_dat = (float *)malloc( ntau_mod * nwav * nrad *
577  sizeof(float) );
578  if( status = ( SDreaddata( tsds_id, istart, NULL, icount,
579  (void *)tmp_dat ) ) != SUCCEED )
580  {
581  printf(
582  "E: %s, %d, Read of cloud table data failed with status: %d\n",
583  __FILE__, __LINE__, status );
584  exit(27);
585  }
586  /* transfer water ref to data blob */
587  for( int32_t irad = 0; irad < nrad; irad++ )
588  {
589  for( int32_t iwav = 0; iwav < nwav; iwav++ )
590  {
591  // for water sfc, they have 1 more tau, so the setting 0
592  // not needed
593  if( *sfc_typ == 0 ) // = land, do this stuff
594  {
595  for( int32_t itau = 0; itau < ntau; itau++ )
596  {
597  *( dat_blob + arr_off + ( itau + tau_off ) + ntau1 *
598  ( iwav + nwav * irad ) ) =
599  *( tmp_dat + irad + nrad * ( iwav + nwav * itau ) );
600  }
601  }
602  else // - water surface, it has a 0 tau in the sds
603  {
604  for( int32_t itau = 0; itau < ntau + 1; itau++ )
605  {
606  *( dat_blob + arr_off + itau + ntau1 *
607  ( iwav + nwav * irad ) ) =
608  *( tmp_dat + irad + nrad * (iwav + nwav * itau ) );
609  }
610  }
611  }
612  }
613  free(tmp_dat);
614  }
615  }
616  /* place data blob pointer into the pt_info */
617  pt_info->dat_ptrs[idat] = (void *)dat_blob;
618  /* set status of the point to say data is there */
619  pt_info->pt_status[idat] = 1;
620  }
621  }
622 
623  /* update knowlege of point status and data location */
624  if( *sfc_typ == 0 ) // over land
625  {
627  land_mgr.add_pts( n_new_blobs );
628  }
629  else
630  {
632  wtr_mgr.add_pts( n_new_blobs );
633  }
634  }
635 
636  /* at this point, go through the interval grid points and weight the
637  proper data and sum to the interpolated array */
638  // NOTE - this assumes the interpolated array is ZERO on input
639  int32_t sfc_typ_mult [] = { 0, 0, 1, 2 };
640  int32_t sfc_typ_off = 2 * ( wtrcld_nelts + icecld_nelts );
641 
642  if( *meas_typ == 0 ) // refl
643  {
644  off_wtr = 0 + sfc_typ_mult[*sfc_typ] * sfc_typ_off;
645  off_ice = 2 * wtrcld_nelts + sfc_typ_mult[*sfc_typ] * sfc_typ_off;
646  }
647  else // SD
648  {
649  off_wtr = wtrcld_nelts + sfc_typ_mult[*sfc_typ] * sfc_typ_off;
650  off_ice = 2 * wtrcld_nelts + icecld_nelts +
651  sfc_typ_mult[*sfc_typ] * sfc_typ_off;
652  }
653  for( int idat = 0; idat < ndat; idat++ )
654  {
655  // this is just for doc, may go soon
656  //printf( "%s, %d, for interval at point %d, weight: %f\n", __FILE__,
657  // __LINE__, idat, pt_info->wt_pt[idat] );
658  for( int32_t itau = 0; itau < ntau1; itau++ )
659  {
660  for( int32_t iwav = 0; iwav < nwav; iwav++ )
661  {
662  for( int32_t irad = 0; irad < nrad_wtr; irad++ )
663  {
664  *( wtr_int + itau + ntau1 * ( iwav + nwav * irad ) ) +=
665  pt_info->wt_pt[idat] *
666  *( (float *)pt_info->dat_ptrs[idat] + off_wtr + itau + ntau1 *
667  ( iwav + nwav * irad ) );
668  }
669 
670  for( int32_t irad = 0; irad < nrad_ice; irad++ )
671  {
672  *( ice_int + itau + ntau1 * ( iwav + nwav * irad ) ) +=
673  pt_info->wt_pt[idat] *
674  *( (float *)pt_info->dat_ptrs[idat] + off_ice + itau + ntau1 *
675  ( iwav + nwav * irad ) );
676  }
677  }
678  }
679  }
680  }
681  /* at the end... */
682  //printf( "I %s, %d, interval loc: %d, %d, %d\n", __FILE__, __LINE__, pt_info->pt_base_loc[0], pt_info->pt_base_loc[1], pt_info->pt_base_loc[2] );
683  return 0;
684  }
int32_t init_dim(int32_t, int32_t, double *)
Definition: dim_mgr.cpp:172
int32_t wtrcld_nelts
Definition: mng_ms.cpp:30
void ** dat_ptrs
Definition: dim_mgr.hpp:70
ms_finfo_struc icecld_landsfc_info
Definition: mng_ms.cpp:21
int status
Definition: l1_czcs_hdf.c:32
#define FAIL
Definition: ObpgReadGrid.h:18
int32_t mng_ms_init_(char *fil_wtrcld_landsfc, char *fil_icecld_landsfc, char *fil_wtrcld_wtrsfc1, char *fil_wtrcld_wtrsfc1_sd, char *fil_icecld_wtrsfc1, char *fil_icecld_wtrsfc1_sd, char *fil_wtrcld_wtrsfc2, char *fil_wtrcld_wtrsfc2_sd, char *fil_icecld_wtrsfc2, char *fil_icecld_wtrsfc2_sd, char *fil_wtrcld_wtrsfc3, char *fil_wtrcld_wtrsfc3_sd, char *fil_icecld_wtrsfc3, char *fil_icecld_wtrsfc3_sd, int len1, int len2, int len_ww1, int len_ww1sd, int len_iw1, int len_iw1sd, int len_ww2, int len_ww2sd, int len_iw2, int len_iw2sd, int len_ww3, int len_ww3sd, int len_iw3, int len_iw3sd)
Definition: mng_ms.cpp:73
#define NULL
Definition: decode_rs.h:63
ms_finfo_struc wtrcld_wtrsfc_info[3]
Definition: mng_ms.cpp:23
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
char * ftn_str_cond(char *str, int len)
Definition: mng_ms.cpp:214
int32_t prune(int32_t)
Definition: dim_mgr.cpp:816
double * wt_pt
Definition: dim_mgr.hpp:69
int32_t * pt_status
Definition: dim_mgr.hpp:63
void add_pts(int32_t)
Definition: dim_mgr.cpp:113
ms_finfo_struc wtrcld_landsfc_info
Definition: mng_ms.cpp:21
int32_t icecld_nelts
Definition: mng_ms.cpp:30
int32_t * dim_siz
Definition: mng_ms.cpp:16
ms_finfo_struc icecld_wtrsfc_info[3]
Definition: mng_ms.cpp:23
pt_info_struc * mng_pt(double *, int32_t, int32_t *)
Definition: dim_mgr.cpp:241
void update_new_data()
Definition: dim_mgr.cpp:202
int32_t interval_needs_data
Definition: dim_mgr.hpp:62
const char * str
Definition: l1c_msi.cpp:35
int32_t linear_to_offset(int32_t n_dim, int32_t dat_ix, int32_t *offset)
Definition: dim_mgr.cpp:1010
double ** dim_vals
Definition: mng_ms.cpp:17
dim_mgr land_mgr(3, 83)
int mng_ms_get_lib_(float *sol_ang, float *sen_ang, float *relaz_ang, int *sfc_typ, int *meas_typ, int *scan, float *wtr_int, float *ice_int, int *stat)
Definition: mng_ms.cpp:339
int32_t set_ms_file(char *fil_ms, char *fil_ms_std, int32_t fil_typ, ms_finfo_struc *finfo)
Definition: mng_ms.cpp:224
Extra metadata that will be written to the HDF4 file l2prod rank
dim_mgr wtr_mgr(3, 83)
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 and prod_ix set to PARAM_TYPE_INT name_prefix is compared with the beginning of the product name If name_suffix is not empty the it must match the end of the product name The characters right after the prefix are read as an integer and prod_ix is set to that number strncpy(l2prod->name_prefix, "myprod", UNITLEN)
l2prod offset
int32_t * pt_base_loc
Definition: dim_mgr.hpp:66