OB.DAAC Logo
NASA Logo
Ocean Color Science Software

ocssw V2022
GEO_landsea_mask.c
Go to the documentation of this file.
1 #include "GEO_DEM.h"
2 #include "GEO_output.h"
3 #include "GEO_earth.h"
4 #include "PGS_MODIS_35251.h"
5 
6 PGSt_SMF_status GEO_landsea_mask(
7  int num_samples,
8  int num_detectors,
9  double terrain_sample_position[][MAX_PADDED][3],
10  uint8 sample_flags[][MAX_PADDED],
11  uint8 *land_seamask_qaflag,
12  uint8 sample_landsea[][MAX_PADDED]
13 )
14 
15 /*
16 !C***************************************************************************
17 !Description:
18  subroutine in output group of the Level-1A geolocation software to
19  retrieve the EOS format land/sea mask data for the scan from the SDP
20  Toolkit.
21 
22 !Input Parameters:
23  num_samples The number of samples to be processed.
24  terrain_sample_position High resolution ground locations in geodetic
25  coordinates.
26  sample_flags High resolution pixel flags.
27 
28 !Output Parameters:
29  landseamask_qaflag QA flag for land/water mask issues.
30  sample_landsea High resolution Land/Water mask values
31 
32 Return Values:
33  MODIS_E_BAD_INPUT_ARG If any argument is invalid
34  MODIS_E_GEO If any subroutine fails
35  PGS_S_SUCCESS Otherwise
36 
37 Externally Defined:
38  DEM_resolutions "GEO_DEM.h"
39  DETECTORS_QKM "GEO_geo.h"
40  L_SMASK_FVALUE "GEO_geo.h"
41  LAT_IDX "GEO_earth.h"
42  LON_IDX "GEO_earth.h"
43  max_lon "GEO_DEM.h"
44  MAX_PADDED "GEO_geo.h"
45  min_lat "GEO_DEM.h"
46  MODIS_E_BAD_INPUT_ARG "PGS_MODIS_35251.h"
47  MODIS_E_GEO "PGS_MODIS_35251.h"
48  NO_ELLIPSE_INTERSECT "GEO_geo.h"
49  PGS_FALSE "PGS_SMF.h"
50  PGS_S_SUCCESS "PGS_SMF.h"
51 
52 Called by:
53  GEO_locate_one_scan "GEO_earth.h"
54 
55 Routines Called:
56  PGS_DEM_GetPoint "PGS_DEM.h"
57  PGS_SMF_TestErrorLevel "PGS_SMF.h"
58 
59 !Revision History:
60  $Log: GEO_landsea_mask.c,v $
61  Revision 6.6 2011/02/09 19:42:14 kuyper
62  Changed to work at high resolution, rather than 1km resolution.
63  Change to pass landsea information directly, rather than through scan_data.
64 
65  Revision 6.5 2010/12/16 22:08:37 kuyper
66  Switched from using *_index to *_IDX.
67 
68  Revision 6.4 2010/08/03 20:18:18 kuyper
69  Changed to use 15 arc second DEM files, falling back to 30 arc second if
70  necessary.
71  Changed to use GEO_DEM.h
72 
73  Revision 6.3 2009/09/11 18:28:54 kuyper
74  Corrected location of calculation of the pixel index.
75 
76  Revision 6.2 2009/05/31 01:30:11 ltan
77  Minor corrections.
78 
79  Revision 6.1 2009/05/29 15:03:25 ltan
80  Moved inputs to proper section of prolog.
81  Changed MAX_SCAN_SAMPLE to MAX_FRAMES, MAX_DETECTORS to DETECTORS_1KM.
82  Changed num_samples to num_frames.
83  Changed to use terrain_frame_position and frame_flags, rather than the corresponding members of scan_data.
84  Changed to return status codes.
85 
86  Revision 5.1 2004/11/03 17:53:01 vlin
87  unsigned integer conversion for scan_data->land_seamask.
88  vlin@saicmodis.com
89 
90  Revision 4.2 2003/08/15 14:23:53 vlin
91  Changed to process by detector within a frame, rather than by frame within
92  a detector, in hopes of improving tile caching by PGS_DEM_GetPoint().
93 
94  Revision 4.1 2003/02/21 22:07:48 kuyper
95  Corrected to use void* pointers with %p format code.
96 
97 Requirements: PR03-F-3.8-1
98 
99 !Team-unique Header:
100  This software is developed by the MODIS Science Data Support
101  Team for the National Aeronautics and Space Administration,
102  Goddard Space Flight Center, under contract NAS5-32373.
103 
104 !END***************************************************************************
105 */
106 
107 {
108  PGSt_double *latitude = (PGSt_double*) malloc(DETECTORS_QKM*MAX_PADDED*sizeof(PGSt_double));
109  PGSt_double *longitude = (PGSt_double*) malloc(DETECTORS_QKM*MAX_PADDED*sizeof(PGSt_double));
110  PGSt_integer numPoints=0; /* The number of point to retrieve.*/
111  int samp, det;
112  int8 *land_seamask = (int8*) malloc(DETECTORS_QKM*MAX_PADDED*sizeof(int8));
113  char msgbuf[256];
114  char filefunc[] = __FILE__ ", GEO_landsea_mask";
115 
116  if (terrain_sample_position == NULL || sample_landsea == NULL ||
117  sample_flags == NULL)
118  {
119  sprintf(msgbuf, "terrain_sample_position:%p sample_landsea:%p\n"
120  "sample_flags:%p\n", (void*)terrain_sample_position,
121  (void*)sample_landsea, (void*)sample_flags);
122  modsmf(MODIS_E_BAD_INPUT_ARG, msgbuf, filefunc);
123 
124  return MODIS_E_BAD_INPUT_ARG;
125  }
126  else if (num_samples <0 || num_samples>MAX_PADDED ||
127  num_detectors <=0 || num_detectors>DETECTORS_QKM)
128  {
129  sprintf(msgbuf, "samples:%d detectors:%d", num_samples, num_detectors);
130  modsmf(MODIS_E_BAD_INPUT_ARG, msgbuf, filefunc);
131 
132  return MODIS_E_BAD_INPUT_ARG;
133  }
134 
135  numPoints = 0;
136 
137 /* The LandWater SDS in the DEM files is chunked with a chunk size of 600x800,
138  * corresponding to 20 degrees by 26 degrees 40 arcmin. Therefore, since MODIS
139  * scans are long and skinny, processing by detector within a sample, rather
140  * than by sample within a detector, improves chunk caching by
141  * PGS_DEM_GetPoint().
142  */
143  for (samp = 0; samp < num_samples; samp++)
144  {
145  for (det = 0; det < num_detectors; det++)
146  {
147  if (sample_flags[det][samp] < NO_ELLIPSE_INTERSECT)
148  { /* Convert lat/lon to degrees and pack earth-viewing locations
149  for input to PGS_DEM_GetPoint: */
150  latitude[numPoints] =
151  terrain_sample_position[det][samp][LAT_IDX] * RAD2DEG;
152  longitude[numPoints] =
153  terrain_sample_position[det][samp][LON_IDX] * RAD2DEG;
154 
155  /* Adjust to avoid region that the SDP Toolkit handles incorrectly*/
156  if (latitude[numPoints] < min_lat)
157  latitude[numPoints] = min_lat;
158  if (longitude[numPoints] > max_lon)
159  longitude[numPoints] = max_lon;
160  numPoints++;
161  }
162  }
163  }
164 
165  if (numPoints == 0)
166  {
167  *land_seamask_qaflag = 1; /* indicating no valid data */
168 
169  /* Normal processing may include 0-sample scans, or scans with bad data*/
170  return PGS_S_SUCCESS;
171  }
172 
173  if (PGS_SMF_TestErrorLevel( PGS_DEM_GetPoint(DEM_resolutions, RESOLUTIONS,
174  PGSd_DEM_WATER_LAND, PGSd_DEM_DEGREE, latitude, longitude, numPoints,
175  PGSd_DEM_NEAREST_NEIGHBOR, land_seamask)) != PGS_FALSE)
176  {
177  *land_seamask_qaflag = 1; /* indicating invalid data. */
178  return MODIS_E_GEO;
179  }
180  else
181  { /* indicate (at least some) valid land/sea mask data */
182  *land_seamask_qaflag = 0;
183  numPoints = 0;
184  for (samp = 0; samp < num_samples; samp++)
185  for (det = 0; det < num_detectors; det++)
186  { /* unpack the land_seamask data to match the other scan data
187  * for output.
188  */
189  if (sample_flags[det][samp] >= NO_ELLIPSE_INTERSECT){
190  sample_landsea[det][samp] = (uint8)L_SMASK_FVALUE;
191  }
192  else {
193  if (land_seamask[numPoints] >= 0)
194  sample_landsea[det][samp] = (uint8)land_seamask[numPoints];
195  numPoints++;
196  }
197  }
198  }
199 
200  free(latitude);
201  free(longitude);
202  free(land_seamask);
203 
204  return PGS_S_SUCCESS;
205 }
real(single), dimension(:,:), allocatable longitude
#define NULL
Definition: decode_rs.h:63
#define MODIS_E_BAD_INPUT_ARG
const unsigned char NO_ELLIPSE_INTERSECT
real(single), dimension(:,:), allocatable latitude
#define MAX_PADDED
Definition: GEO_geo.h:84
@ L_SMASK_FVALUE
Definition: GEO_geo.h:130
#define MODIS_E_GEO
PGSt_SMF_status GEO_landsea_mask(int num_samples, int num_detectors, double terrain_sample_position[][MAX_PADDED][3], uint8 sample_flags[][MAX_PADDED], uint8 *land_seamask_qaflag, uint8 sample_landsea[][MAX_PADDED])
#define RAD2DEG
Definition: GEO_geo.h:173
PGSt_DEM_Tag DEM_resolutions[RESOLUTIONS]
Definition: GEO_DEM.c:62
@ LON_IDX
@ LAT_IDX
double max_lon
Definition: GEO_DEM.c:64
double min_lat
Definition: GEO_DEM.c:63
#define DETECTORS_QKM
Definition: GEO_geo.h:87
#define RESOLUTIONS
Definition: GEO_DEM.h:52