NASA Logo
Ocean Color Science Software

ocssw V2022
ch_cld_sci.f90
Go to the documentation of this file.
1 subroutine ch_cld_sci( tdat, nv, ubdat, nubyte, i32dat, ni32, &
2  sensor_id, cloud_hgt_file )
3 !
4 ! ch_cld_sci is the entry to start the CHIMAERA processing for the 3-
5 ! line radiance (and other) input data for 1 line of data
6 !
7 ! tdat, nv Float array of all real data and the number of values
8 ! ubdat, nubyte Unsigned byte values for processing and the number of values
9 ! i32dat, ni32 integer values for processing and the number of values
10 ! sensor_id integer sensor identification number from l2gen
11 ! cloud_hgt_file string name of the cloud height file
12 !
13 ! (see below for the unpacking of these values into all inputs)
14 !
15 ! W. Robinson, SAIC, Dec 2018
16 !
17 use ch_xfr
18 implicit none
19 
20 integer :: npix, nlin, nbnd, nbnd_albedo, nlvl_model, scan, st_samp, &
21  g_year, g_day
22 common /dim_ctl/ npix, nlin, nbnd, nbnd_albedo, nlvl_model, scan, st_samp, &
23  g_year, g_day
24 integer :: nv, i,j,k, na, off, nubyte, ni32, sensor_id, nbnd_sfc_albedo
25 real, dimension(nv) :: tdat
26 integer, dimension(3) :: dim_3
27 integer*1, dimension(nubyte) :: ubdat
28 integer, dimension(ni32) :: i32dat
29 character(len=500) :: cloud_hgt_file
30 real, parameter :: bad_float = -32767.
31 
32 c2_sensor_id = sensor_id ! a global sensor name
33 if( ( c2_sensor_id == oci_id ) .OR. ( c2_sensor_id == ocis_id ) ) then
34  nbnd_sfc_albedo = 6
35 else
36  nbnd_sfc_albedo = 5
37 endif
38 
39 ! set up the l2gen values that will go into Chimaera
40 if( .not. allocated( c2_refl ) ) then
41  allocate(c2_refl(npix, nlin, nbnd))
42  allocate(c2_bnd_unc(npix, nbnd_albedo, nlin))
43  allocate(c2_spec_unc(nbnd_albedo))
44  allocate(c2_unc_sf(nbnd_albedo))
45  allocate(c2_lat(npix, nlin))
46  allocate(c2_lon(npix, nlin))
47  allocate(c2_senz(npix, nlin))
48  allocate(c2_sena(npix, nlin))
49  allocate(c2_sola(npix, nlin))
50  allocate(c2_solz(npix, nlin))
51  allocate(c2_relaz(npix, nlin))
52  allocate(c2_prof_mixr(npix,nlin,nlvl_model))
53  allocate(c2_prof_t(npix,nlin,nlvl_model))
54  allocate(c2_prof_p(npix,nlin,nlvl_model))
55  allocate(c2_prof_hgt(npix,nlin,nlvl_model))
56  allocate(c2_tsfc(npix,nlin))
57  allocate(c2_psfc(npix,nlin))
58  allocate(c2_wind(npix,nlin))
59  allocate(c2_tot_o3(npix,nlin))
60  allocate(c2_ice_frac(npix,nlin))
61  allocate(c2_snow_frac(npix,nlin))
62  allocate(c2_alt_o3(npix,nlin))
63  allocate(c2_alt_icec(npix,nlin))
64  allocate(c2_sfc_albedo(npix,nlin,nbnd_sfc_albedo))
65  allocate(c2_sfc_lvl(npix,nlin))
66  allocate(c2_trop_lvl(npix,nlin))
67  allocate(c2_cld_det(npix,nlin))
68  allocate(c2_conf_cld(npix,nlin))
69  allocate(c2_clr_66(npix,nlin))
70  allocate(c2_clr_95(npix,nlin))
71  allocate(c2_clr_99(npix,nlin))
72  allocate(c2_sno_sfc(npix,nlin))
73  allocate(c2_wtr_sfc(npix,nlin))
74  allocate(c2_coast_sfc(npix,nlin))
75  allocate(c2_desert_sfc(npix,nlin))
76  allocate(c2_lnd_sfc(npix,nlin))
77  allocate(c2_night(npix,nlin))
78  allocate(c2_glint(npix,nlin))
79  allocate(c2_ocean_no_snow(npix,nlin))
80  allocate(c2_ocean_snow(npix,nlin))
81  allocate(c2_lnd_no_snow(npix,nlin))
82  allocate(c2_lnd_snow(npix,nlin))
83  allocate(c2_tst_h_138(npix,nlin))
84  allocate(c2_tst_vis_refl(npix,nlin))
85  allocate(c2_tst_vis_ratio(npix,nlin))
86  allocate(c2_vis_cld_250(npix,nlin))
87  allocate(c2_appl_hcld_138(npix,nlin))
88  allocate(c2_appl_vis_refl(npix,nlin))
89  allocate(c2_appl_vis_nir_ratio(npix,nlin))
90  allocate(c2_alt_snowice(npix,nlin))
91  allocate(c2_cmp_there(nbnd))
92  ! for out point diagnostics
93  allocate(prd_out_refl_loc_2100(npix,nlin,10))
94  allocate(prd_out_refl_loc_1600(npix,nlin,10))
95  allocate(prd_out_refl_loc_2200(npix,nlin,10))
96  allocate(prd_out_refl_loc_1621(npix,nlin,10))
97  !
98  allocate(c2_cth(npix,nlin))
99  allocate(c2_ctp(npix,nlin))
100  allocate(c2_ctt(npix,nlin))
101  endif
102 c2_npix = npix
103 c2_nbnd = nbnd
104 c2_nlin = nlin
105 c2_nbnd_albedo = nbnd_albedo
106 c2_nlvl_model = nlvl_model
107 c2_scan = scan
108 c2_st_samp = st_samp
109 c2_g_year = g_year
110 c2_g_day = g_day
111 !
112 ! Start with the L1b information, rads, geoloc, view angles
113 na = npix * nbnd * nlin
114 dim_3(1) = npix
115 dim_3(2) = nbnd
116 dim_3(3) = nlin
117 c2_refl = reshape( tdat(1:na), (/ npix, nbnd, nlin /) )
118 !
119 off = na
120 na = npix * nbnd_albedo * nlin
121 c2_bnd_unc = reshape( tdat(1+off:na+off), (/ npix, nbnd_albedo, nlin /) )
122 !
123 off = off + na
124 na = nbnd_albedo
125 c2_spec_unc = tdat(1+off:na+off)
126 !
127 off = off + na
128 c2_unc_sf = tdat(1+off:na+off)
129 !
130 off = off + na
131 na = npix * nlin
132 c2_lat = reshape( tdat(1+off:na+off), (/ npix, nlin /) )
133 !
134 off = off + na
135 c2_lon = reshape( tdat(1+off:na+off), (/ npix, nlin /) )
136 !
137 off = off + na
138 c2_senz = reshape( tdat(1+off:na+off), (/ npix, nlin /) )
139 where( c2_senz >= 90. ) c2_senz = bad_float
140 !
141 off = off + na
142 c2_sena = reshape( tdat(1+off:na+off), (/ npix, nlin /) )
143 !
144 off = off + na
145 c2_sola = reshape( tdat(1+off:na+off), (/ npix, nlin /) )
146 !
147 off = off + na
148 c2_solz = reshape( tdat(1+off:na+off), (/ npix, nlin /) )
149 where( c2_solz >= 90. ) c2_solz = bad_float
150 !
151 off = off + na
152 c2_relaz = reshape( tdat(1+off:na+off), (/ npix, nlin /) )
153 !
154 ! The met information, profiles, 2D fields
155 off = off + na
156 na = npix * nlin * nlvl_model
157 c2_prof_mixr = reshape( tdat(1+off:na+off), (/npix, nlin,nlvl_model/) )
158 !
159 off = off + na
160 c2_prof_t = reshape( tdat(1+off:na+off), (/npix, nlin,nlvl_model/) )
161 !
162 off = off + na
163 c2_prof_p = reshape( tdat(1+off:na+off), (/npix, nlin,nlvl_model/) )
164 !
165 off = off + na
166 c2_prof_hgt = reshape( tdat(1+off:na+off), (/npix, nlin,nlvl_model/) )
167 !
168 off = off + na
169 na = npix * nlin
170 c2_tsfc = reshape( tdat(1+off:na+off), (/npix, nlin/) )
171 !
172 off = off + na
173 c2_psfc = reshape( tdat(1+off:na+off), (/npix, nlin/) )
174 !
175 off = off + na
176 c2_wind = reshape( tdat(1+off:na+off), (/npix, nlin/) )
177 !
178 off = off + na
179 c2_tot_o3 = reshape( tdat(1+off:na+off), (/npix, nlin/) )
180 !
181 off = off + na
182 c2_ice_frac = reshape( tdat(1+off:na+off), (/npix, nlin/) )
183 !
184 off = off + na
185 c2_snow_frac = reshape( tdat(1+off:na+off), (/npix, nlin/) )
186 !
187 off = off + na
188 c2_alt_o3 = reshape( tdat(1+off:na+off), (/npix, nlin/) )
189 !
190 off = off + na
191 c2_alt_icec = reshape( tdat(1+off:na+off), (/npix, nlin/) )
192 !
193 off = off + na
194 c2_sfc_albedo = reshape( tdat(1+off: na*nbnd_sfc_albedo + off), &
195  (/npix, nlin, nbnd_sfc_albedo/))
196 !
197 ! the new cloud top height, press, temp
198 off = off + na * nbnd_sfc_albedo
199 c2_cth = reshape( tdat(1+off:na+off), (/npix, nlin/) )
200 off = off + na
201 c2_ctp = reshape( tdat(1+off:na+off), (/npix, nlin/) )
202 off = off + na
203 c2_ctt = reshape( tdat(1+off:na+off), (/npix, nlin/) )
204 !
205 ! byte met data
206 na = npix * nlin
207 off = 0
208 c2_sfc_lvl = reshape( ubdat(1+off:na+off), (/npix, nlin/) )
209 !
210 off = off + na
211 c2_trop_lvl = reshape( ubdat(1+off:na+off), (/npix, nlin/) )
212 !
213 ! the rest are the cloud mask byte data
214 off = off + na
215 c2_cld_det = reshape( ubdat(1+off:na+off), (/npix, nlin/) )
216 !
217 off = off + na
218 c2_conf_cld = reshape( ubdat(1+off:na+off), (/npix, nlin/) )
219 !
220 off = off + na
221 c2_clr_66 = reshape( ubdat(1+off:na+off), (/npix, nlin/) )
222 !
223 off = off + na
224 c2_clr_95 = reshape( ubdat(1+off:na+off), (/npix, nlin/) )
225 !
226 off = off + na
227 c2_clr_99 = reshape( ubdat(1+off:na+off), (/npix, nlin/) )
228 !
229 off = off + na
230 c2_sno_sfc = reshape( ubdat(1+off:na+off), (/npix, nlin/) )
231 !
232 off = off + na
233 c2_wtr_sfc = reshape( ubdat(1+off:na+off), (/npix, nlin/) )
234 !
235 off = off + na
236 c2_coast_sfc = reshape( ubdat(1+off:na+off), (/npix, nlin/) )
237 !
238 off = off + na
239 c2_desert_sfc = reshape( ubdat(1+off:na+off), (/npix, nlin/) )
240 !
241 off = off + na
242 c2_lnd_sfc = reshape( ubdat(1+off:na+off), (/npix, nlin/) )
243 !
244 off = off + na
245 c2_night = reshape( ubdat(1+off:na+off), (/npix, nlin/) )
246 !
247 off = off + na
248 c2_glint = reshape( ubdat(1+off:na+off), (/npix, nlin/) )
249 !
250 off = off + na
251 c2_ocean_no_snow = reshape( ubdat(1+off:na+off), (/npix, nlin/) )
252 !
253 off = off + na
254 c2_ocean_snow = reshape( ubdat(1+off:na+off), (/npix, nlin/) )
255 !
256 off = off + na
257 c2_lnd_no_snow = reshape( ubdat(1+off:na+off), (/npix, nlin/) )
258 !
259 off = off + na
260 c2_lnd_snow = reshape( ubdat(1+off:na+off), (/npix, nlin/) )
261 !
262 off = off + na
263 c2_tst_h_138 = reshape( ubdat(1+off:na+off), (/npix, nlin/) )
264 !
265 off = off + na
266 c2_tst_vis_refl = reshape( ubdat(1+off:na+off), (/npix, nlin/) )
267 !
268 off = off + na
269 c2_tst_vis_ratio = reshape( ubdat(1+off:na+off), (/npix, nlin/) )
270 !
271 off = off + na
272 c2_vis_cld_250 = reshape( ubdat(1+off:na+off), (/npix, nlin/) )
273 !
274 off = off + na
275 c2_appl_hcld_138 = reshape( ubdat(1+off:na+off), (/npix, nlin/) )
276 !
277 off = off + na
278 c2_appl_vis_refl = reshape( ubdat(1+off:na+off), (/npix, nlin/) )
279 !
280 off = off + na
281 c2_appl_vis_nir_ratio = reshape( ubdat(1+off:na+off), (/npix, nlin/) )
282 !
283 off = off + na
284 c2_cmp_there = ubdat( 1+off:nbnd+off )
285 !
286 !
287 ! do the integer data
288 na = npix * nlin
289 c2_alt_snowice = reshape( i32dat(1:na), (/npix, nlin/) )
290 !
291 ! set a global height file name
292 c2_cloud_hgt_file = cloud_hgt_file
293 !
294 ! Here's what needs to be done:
295 ! - Make local values of the dim sizes from the common values
296 ! - Extract real arrays to correctly allocated science-familliar arrays
297 ! - call modified scienceinterface with the right stuff
298 !i = 1
299 !k = 2
300 ! init the output diagnostics
301 ! WDR this done in general_science_module already prd_out_refl_loc_2100 = -999.
302 ! WDR rewind prd_out_refl_loc_1600 = -999.
303 !prd_out_refl_loc_1600 = -999.
304 !
305 ! we zero out the L1B data and met less the center line for test purposes
306 ! and to see what needs surrounding lines
307 ! WDR out now test done call zero_surround( )
308 !
309 ! OK, segway into the chimaera code
310 call driver_mod_pr06od( )
311  return
312 end subroutine ch_cld_sci
313 
314 !
315 ! so zero_surround will just zero the values in surrounding L1b lines
316 ! (for now) mainly to test that no area data is needed in doing the
317 ! cloud work except for the cloud stuff
318 !
319 subroutine zero_surround( )
321 integer :: clin
322 
323 clin = c2_nlin / 2 + 1
324 ! reflectance
325 c2_refl(:,:,1:clin - 1) = 0.
326 c2_refl(:,:,clin + 1:c2_nlin ) = 0.
327 ! uncertainty
328 c2_bnd_unc(:,:,1:clin - 1) = 0.
329 c2_bnd_unc(:,:,1:clin + 1:c2_nlin ) = 0.
330 ! lat ** leave out - quite possibly, used for the anc (which is still
331 ! processed in the CHIMAERA code)
332 !c2_lat(:,1:clin - 1) = 0.
333 !c2_lat(:,clin + 1:c2_nlin ) = 0.
334 ! lon
335 !c2_lon(:,1:clin - 1) = 0.
336 !c2_lon(:,clin + 1:c2_nlin ) = 0.
337 ! senz
338 !c2_senz(:,1:clin - 1) = 0.
339 !c2_senz(:,clin + 1:c2_nlin ) = 0.
340 ! sena
341 !c2_sena(:,1:clin - 1) = 0.
342 !c2_sena(:,clin + 1:c2_nlin ) = 0.
343 ! solz
344 !c2_solz(:,1:clin - 1) = 0.
345 !c2_solz(:,clin + 1:c2_nlin ) = 0.
346 ! sola
347 !c2_sola(:,1:clin - 1) = 0.
348 !c2_sola(:,clin + 1:c2_nlin ) = 0.
349 ! relaz
350 !c2_relaz(:,1:clin - 1) = 0.
351 !c2_relaz(:,clin + 1:c2_nlin ) = 0.
352 end subroutine zero_surround
353 
354 subroutine get_cmp_prod( prod_num, prod_array, n_prd )
355 !
356 ! get_cmp_prod will grab the requested product from the CHIMAERA final
357 ! data arrays
358 !
359 ! prod_num integer I The l2gen-defined product ID number
360 ! prod_array real O array to carry out the line of data
361 ! n_prd integer I # of sub-products in this product, ie the product
362 ! has c2_npix pixels and n_prd sub-products
363 !
364 ! W. Robinson, SAIC, 8 Feb 2019
380  integer :: prod_num, n_prd, iprd
381  !real, dimension(c2_npix, n_prd ) :: prod_array
382  real, dimension(n_prd,c2_npix) :: prod_array
383 ! the same catalog #s defined in l2prod.h
384 integer, parameter :: CAT_Cld_p = 469, cat_cld_t = 470, cat_cld_h = 493
385 integer, parameter :: CAT_CER_2100 = 445
386 integer, parameter :: CAT_CER_1600 = 446, cat_cot_2100 = 447
387 integer, parameter :: CAT_COT_1600 = 448, cat_cer_1621 = 449
388 integer, parameter :: CAT_COT_1621 = 450, cat_cwp_2100 = 451
389 integer, parameter :: CAT_CWP_1621 = 452, cat_cwp_1600 = 453
390 ! these are used in the byte call get_cmp_byte()
391 !integer, parameter :: CAT_Cld_Sfc_Type = 454, CAT_Cld_Phase_2100 = 455
392 !integer, parameter :: CAT_Cld_Non_Abs_Band = 456, CAT_Cld_Phase_1600 = 457
393 !integer, parameter :: CAT_Cld_Phase_1621 = 458
394 integer, parameter :: CAT_Cld_Top_Refl_650 = 459
395 integer, parameter :: CAT_Cld_Top_Refl_860 = 460, cat_cld_top_refl_1200 = 461
396 integer, parameter :: CAT_Cld_Top_Refl_1600 = 462, cat_cld_top_refl_2100 = 463
397 integer, parameter :: CAT_Surface_Albedo_650 = 464, cat_surface_albedo_860 = 465
398 integer, parameter :: CAT_Surface_Albedo_1200 = 466
399 integer, parameter :: CAT_Surface_Albedo_1600 = 467
400 integer, parameter :: CAT_Surface_Albedo_2100 = 468
401 
402 integer, parameter :: CAT_COT_fail_2100 = 471
403 integer, parameter :: CAT_COT_fail_1600 = 472
404 integer, parameter :: CAT_COT_fail_1621 = 473
405 integer, parameter :: CAT_CER_fail_2100 = 474
406 integer, parameter :: CAT_CER_fail_1600 = 475
407 integer, parameter :: CAT_CER_fail_1621 = 476
408 integer, parameter :: CAT_CMP_fail_pct_2100 = 477
409 integer, parameter :: CAT_CMP_fail_pct_1600 = 478
410 integer, parameter :: CAT_CMP_fail_pct_1621 = 479
411 integer, parameter :: CAT_refl_loc_1600 = 480
412 integer, parameter :: CAT_refl_loc_2100 = 481
413 integer, parameter :: CAT_refl_loc_1621 = 482
414 ! recent additions for the 2.2 um band
415 integer, parameter :: CAT_CER_2200 = 483
416 integer, parameter :: CAT_COT_2200 = 484
417 integer, parameter :: CAT_CWP_2200 = 485
418 ! byte integer, parameter :: CAT_Cld_Phase_2200 = 486
419 integer, parameter :: CAT_Cld_Top_Refl_2200 = 487
420 integer, parameter :: CAT_Surface_Albedo_2200 = 488
421 integer, parameter :: CAT_COT_fail_2200 = 489
422 integer, parameter :: CAT_CER_fail_2200 = 490
423 integer, parameter :: CAT_CMP_fail_pct_2200 = 491
424 integer, parameter :: CAT_refl_loc_2200 = 492
425 
426 real, parameter :: bad_float = -32767.
427 !
428 ! for the product type, transfer the center line
429 !
430 select case( prod_num )
431  case( cat_cld_p )
432  prod_array( 1, : ) = cloud_top_pressure( :, c2_nlin / 2 + 1 )
433  case( cat_cld_t )
434  prod_array( 1, : ) = cloud_top_temperature( :, c2_nlin / 2 + 1 )
435  case( cat_cld_h )
436  prod_array( 1, : ) = cloud_top_height( :, c2_nlin / 2 + 1 )
437  case( cat_cer_2100 )
438  prod_array( 1, : ) = effective_radius_21_final( :, c2_nlin / 2 + 1 )
439  case( cat_cer_1600 )
440  prod_array( 1, : ) = effective_radius_16_final( :, c2_nlin / 2 + 1 )
441  case( cat_cer_2200 )
442  prod_array( 1, : ) = effective_radius_22_final( :, c2_nlin / 2 + 1 )
443  case( cat_cot_2100 )
444  prod_array( 1, : ) = optical_thickness_final( :, c2_nlin / 2 + 1 )
445  case( cat_cot_1600 )
446  prod_array( 1, : ) = optical_thickness_16_final( :, c2_nlin / 2 + 1 )
447  case( cat_cot_2200 )
448  prod_array( 1, : ) = optical_thickness_22_final( :, c2_nlin / 2 + 1 )
449  case( cat_cer_1621 )
450  prod_array( 1, : ) = effective_radius_1621_final( :, c2_nlin / 2 + 1 )
451  case( cat_cot_1621 )
452  prod_array( 1, : ) = optical_thickness_1621_final( :, c2_nlin / 2 + 1 )
453  case( cat_cwp_2100 )
454  prod_array( 1, : ) = liquid_water_path( :, c2_nlin / 2 + 1 )
455  case( cat_cwp_1621 )
456  prod_array( 1, : ) = liquid_water_path_1621( :, c2_nlin / 2 + 1 )
457  case( cat_cwp_1600 )
458  prod_array( 1, : ) = liquid_water_path_16( :, c2_nlin / 2 + 1 )
459  case( cat_cwp_2200 )
460  prod_array( 1, : ) = liquid_water_path_22( :, c2_nlin / 2 + 1 )
461  case( cat_cld_top_refl_650 )
462  prod_array( 1, : ) = atm_corr_refl( 1, :, c2_nlin / 2 + 1 )
463  case( cat_cld_top_refl_860 )
464  prod_array( 1, : ) = atm_corr_refl( 2, :, c2_nlin / 2 + 1 )
465  case( cat_cld_top_refl_1200 )
466  prod_array( 1, : ) = atm_corr_refl( 3, :, c2_nlin / 2 + 1 )
467  case( cat_cld_top_refl_1600 )
468  prod_array( 1, : ) = atm_corr_refl( 4, :, c2_nlin / 2 + 1 )
469  case( cat_cld_top_refl_2100 )
470  prod_array( 1, : ) = atm_corr_refl( 5, :, c2_nlin / 2 + 1 )
471  case( cat_cld_top_refl_2200 )
472  prod_array( 1, : ) = atm_corr_refl( 6, :, c2_nlin / 2 + 1 )
473  case( cat_surface_albedo_650 )
474  prod_array( 1, : ) = 0.001 * surface_albedo( :, c2_nlin / 2 + 1, 1 )
475  case( cat_surface_albedo_860 )
476  prod_array( 1, : ) = 0.001 * surface_albedo( :, c2_nlin / 2 + 1, 2 )
477  case( cat_surface_albedo_1200 )
478  prod_array( 1, : ) = 0.001 * surface_albedo( :, c2_nlin / 2 + 1, 3 )
479  case( cat_surface_albedo_1600 )
480  prod_array( 1, : ) = 0.001 * surface_albedo( :, c2_nlin / 2 + 1, 4 )
481  case( cat_surface_albedo_2100 )
482  prod_array( 1, : ) = 0.001 * surface_albedo( :, c2_nlin / 2 + 1, 5 )
483  case( cat_surface_albedo_2200 )
484  prod_array( 1, : ) = 0.001 * surface_albedo( :, c2_nlin / 2 + 1, 6 )
485 
486  case( cat_cot_fail_2100 )
487  prod_array( 1, : ) = failure_metric( :, c2_nlin / 2 + 1 )%tau / 100.
488  where ( prod_array < -90. ) prod_array = bad_float
489  case( cat_cot_fail_1600 )
490  prod_array( 1, : ) = failure_metric_16( :, c2_nlin / 2 + 1 )%tau / 100.
491  where ( prod_array < -90. ) prod_array = bad_float
492  case( cat_cot_fail_1621 )
493  prod_array( 1, : ) = failure_metric_1621( :, c2_nlin / 2 + 1 )%tau / 100.
494  where ( prod_array < -90. ) prod_array = bad_float
495  case( cat_cot_fail_2200 )
496  prod_array( 1, : ) = failure_metric_22( :, c2_nlin / 2 + 1 )%tau / 100.
497  where ( prod_array < -90. ) prod_array = bad_float
498  case( cat_cer_fail_2100 )
499  prod_array( 1, : ) = failure_metric( :, c2_nlin / 2 + 1 )%re / 100.
500  where ( prod_array < -90. ) prod_array = bad_float
501  case( cat_cer_fail_1600 )
502  prod_array( 1, : ) = failure_metric_16( :, c2_nlin / 2 + 1 )%re / 100.
503  where ( prod_array < -90. ) prod_array = bad_float
504  case( cat_cer_fail_1621 )
505  prod_array( 1, : ) = failure_metric_1621( :, c2_nlin / 2 + 1 )%re / 100.
506  where ( prod_array < -90. ) prod_array = bad_float
507  case( cat_cer_fail_2200 )
508  prod_array( 1, : ) = failure_metric_22( :, c2_nlin / 2 + 1 )%re / 100.
509  where ( prod_array < -90. ) prod_array = bad_float
510  case( cat_cmp_fail_pct_2100 )
511  prod_array( 1, : ) = failure_metric( :, c2_nlin / 2 + 1 )%RSS / 100.
512  where ( prod_array < -90. ) prod_array = bad_float
513  case( cat_cmp_fail_pct_1600 )
514  prod_array( 1, : ) = failure_metric_16( :, c2_nlin / 2 + 1 )%RSS / 100.
515  where ( prod_array < -90. ) prod_array = bad_float
516  case( cat_cmp_fail_pct_1621 )
517  prod_array( 1, : ) = failure_metric_1621( :, c2_nlin / 2 + 1 )%RSS / 100.
518  where ( prod_array < -90. ) prod_array = bad_float
519  case( cat_cmp_fail_pct_2200 )
520  prod_array( 1, : ) = failure_metric_22( :, c2_nlin / 2 + 1 )%RSS / 100.
521  where ( prod_array < -90. ) prod_array = bad_float
522  case( cat_refl_loc_1600 )
523  do iprd=1,n_prd
524  prod_array( iprd, :) = prd_out_refl_loc_1600( :, c2_nlin / 2 + 1, iprd )
525  end do
526  !
527  case( cat_refl_loc_2100 )
528  do iprd=1,n_prd
529  prod_array( iprd, :) = prd_out_refl_loc_2100( :, c2_nlin / 2 + 1, iprd )
530  end do
531  case( cat_refl_loc_1621 )
532  do iprd=1,n_prd
533  prod_array( iprd, :) = prd_out_refl_loc_1621( :, c2_nlin / 2 + 1, iprd )
534  end do
535  case( cat_refl_loc_2200 )
536  do iprd=1,n_prd
537  prod_array( iprd, :) = prd_out_refl_loc_2200( :, c2_nlin / 2 + 1, iprd )
538  end do
539  case default
540  print*, "Improper product ID, case encountered, ID:", prod_num
541  prod_array = bad_float
542  end select
543 !
544 ! set bad value to obpg standard
545 where(prod_array < -900. )
546  prod_array = bad_float
547 end where
548 return
549 end subroutine get_cmp_prod
550 
551 subroutine get_cmp_byt( prod_num, bprod )
552 !
553 ! get_cmp_byt is for all the flags / indexes of byte size
554 !
555 ! prod_num integer I The l2gen-defined product ID number
556 ! bprod integer*1 O array to carry out the line of data
557 !
558 ! W. Robinson, SAIC, 5 Aug 2019
559 use ch_xfr, only: c2_nlin, c2_npix
561 integer :: prod_num, lin_mid
562 integer*1, dimension(c2_npix) :: bprod
563 integer, parameter :: CAT_Cld_Sfc_Type = 454, cat_cld_phase_2100 = 455
564 integer, parameter :: CAT_Cld_Non_Abs_Band = 456, cat_cld_phase_1600 = 457
565 integer, parameter :: CAT_Cld_Phase_1621 = 458
566 integer, parameter :: CAT_Cld_Phase_2200 = 486
567 integer, parameter :: CAT_Cld_water_cloud = 440, cat_cld_ice_cloud = 441
568 integer*1 :: iand_comp = 7
569 
570 !
571 ! for the product type, transfer the center line
572 
573 lin_mid = c2_nlin / 2 + 1
574 prod_sel: select case( prod_num )
575  case( cat_cld_sfc_type )
576  bprod = 0
577  where(cloudmask(:,lin_mid)%ocean_snow == 1 ) &
578  bprod = 1
579  where(cloudmask(:,lin_mid)%land_no_snow == 1 ) &
580  bprod = 2
581  where(cloudmask(:,lin_mid)%land_snow == 1 ) &
582  bprod = 3
583  case( cat_cld_phase_2100 )
584  bprod = processing_information(:,lin_mid)%path_and_outcome
585  bprod = ibclr( bprod, 3 )
586  case( cat_cld_non_abs_band )
587  bprod = processing_information(:,lin_mid)%band_used_for_optical_thickness
588  case( cat_cld_phase_1600 )
589  bprod = processing_information(:,lin_mid)%path_and_outcome_16
590  bprod = ibclr( bprod, 3 )
591  case( cat_cld_phase_1621 )
592  bprod = processing_information(:,lin_mid)%path_and_outcome_1621
593  bprod = ibclr( bprod, 3 )
594  case( cat_cld_phase_2200 )
595  bprod = processing_information(:,lin_mid)%path_and_outcome_22
596  bprod = ibclr( bprod, 3 )
597  case( cat_cld_water_cloud )
598  bprod = 0
599  where(iand(processing_information(:,lin_mid)%path_and_outcome,iand_comp) == 2 )
600  bprod = 1
601  elsewhere(iand(processing_information(:,lin_mid)%path_and_outcome,iand_comp) == 4 )
602  bprod = 1
603  endwhere
604  case( cat_cld_ice_cloud )
605  bprod = 0
606  where(iand(processing_information(:,lin_mid)%path_and_outcome,iand_comp) == 3 ) &
607  bprod = 1
608  case default
609  print*, "Improper product ID, case encountered, ID:", prod_num
610  return
611 
612  end select prod_sel
613  return
614 end subroutine get_cmp_byt
Definition: ch_xfr.f90:1
type(failed_type), dimension(:,:), allocatable failure_metric_22
real, dimension(:,:), allocatable c2_alt_icec
Definition: ch_xfr.f90:20
integer *1, dimension(:,:), allocatable c2_trop_lvl
Definition: ch_xfr.f90:26
real, dimension(:,:), allocatable c2_psfc
Definition: ch_xfr.f90:17
integer, dimension(:,:), allocatable c2_alt_snowice
Definition: ch_xfr.f90:27
integer *1, dimension(:,:), allocatable c2_wtr_sfc
Definition: ch_xfr.f90:31
real(single), dimension(:,:), allocatable liquid_water_path
Definition: core_arrays.f90:65
integer c2_g_year
Definition: ch_xfr.f90:46
integer *1, dimension(:,:), allocatable c2_clr_99
Definition: ch_xfr.f90:30
integer c2_nbnd_albedo
Definition: ch_xfr.f90:46
integer refl_abs
Definition: ch_xfr.f90:66
real, dimension(:,:,:), allocatable c2_prof_hgt
Definition: ch_xfr.f90:16
real, dimension(:,:,:), allocatable c2_bnd_unc
Definition: ch_xfr.f90:10
integer *1, dimension(:,:), allocatable c2_sno_sfc
Definition: ch_xfr.f90:31
integer ocis_id
Definition: ch_xfr.f90:51
real, dimension(:,:), allocatable c2_lon
Definition: ch_xfr.f90:12
type(cloudmask_type), dimension(:,:), allocatable cloudmask
integer *1, dimension(:,:), allocatable c2_appl_hcld_138
Definition: ch_xfr.f90:38
real(single), dimension(:,:), allocatable cloud_top_pressure
integer *1, dimension(:,:), allocatable c2_vis_cld_250
Definition: ch_xfr.f90:37
real, dimension(:,:), allocatable c2_sena
Definition: ch_xfr.f90:13
integer *1, dimension(:,:), allocatable c2_tst_h_138
Definition: ch_xfr.f90:36
integer *1, dimension(:,:), allocatable c2_ocean_snow
Definition: ch_xfr.f90:34
integer *1, dimension(:,:), allocatable c2_desert_sfc
Definition: ch_xfr.f90:32
integer tbl_lo_nabs_ice
Definition: ch_xfr.f90:66
real, dimension(:), allocatable c2_unc_sf
Definition: ch_xfr.f90:11
unsigned char * get_cmp_byt(l2str *l2rec, int prodnum)
Definition: get_cmp.c:128
real, dimension(:,:,:), allocatable c2_prof_mixr
Definition: ch_xfr.f90:15
real, dimension(:,:,:), allocatable prd_out_refl_loc_1600
Definition: ch_xfr.f90:63
real, dimension(:,:), allocatable c2_solz
Definition: ch_xfr.f90:13
real, dimension(:,:,:), allocatable c2_sfc_albedo
Definition: ch_xfr.f90:22
real, dimension(:,:,:), allocatable atm_corr_refl
type(failed_type), dimension(:,:), allocatable failure_metric_1621
real, dimension(:,:,:), allocatable surface_albedo
integer *1, dimension(:,:), allocatable c2_night
Definition: ch_xfr.f90:33
real, dimension(:,:), allocatable c2_lat
Definition: ch_xfr.f90:12
real, dimension(:,:), allocatable c2_ctp
Definition: ch_xfr.f90:24
integer c2_nlin
Definition: ch_xfr.f90:46
real(single), dimension(:,:), allocatable cloud_top_height
real, dimension(:,:), allocatable c2_tsfc
Definition: ch_xfr.f90:17
type(failed_type), dimension(:,:), allocatable failure_metric_16
real, dimension(:,:,:), allocatable prd_out_refl_loc_2200
Definition: ch_xfr.f90:65
integer tbl_hi_nabs_ice
Definition: ch_xfr.f90:66
integer refl_nabs
Definition: ch_xfr.f90:66
integer *1, dimension(:,:), allocatable c2_cld_det
Definition: ch_xfr.f90:29
integer *1, dimension(:,:), allocatable c2_appl_vis_refl
Definition: ch_xfr.f90:39
real, dimension(:,:), allocatable c2_relaz
Definition: ch_xfr.f90:13
real(single), dimension(:,:), allocatable optical_thickness_final
Definition: core_arrays.f90:38
void print(std::ostream &stream, const char *format)
Definition: PrintDebug.hpp:38
subroutine get_cmp_prod(prod_num, prod_array, n_prd)
Definition: ch_cld_sci.f90:355
integer c2_npix
Definition: ch_xfr.f90:46
real, dimension(:,:), allocatable c2_ice_frac
Definition: ch_xfr.f90:19
type(failed_type), dimension(:,:), allocatable failure_metric
real(single), dimension(:,:), allocatable optical_thickness_22_final
Definition: core_arrays.f90:39
integer *1, dimension(:,:), allocatable c2_tst_vis_refl
Definition: ch_xfr.f90:36
real(single), dimension(:,:), allocatable liquid_water_path_1621
Definition: core_arrays.f90:69
character(len=500) c2_cloud_hgt_file
Definition: ch_xfr.f90:57
subroutine ch_cld_sci(tdat, nv, ubdat, nubyte, i32dat, ni32, sensor_id, cloud_hgt_file)
Definition: ch_cld_sci.f90:3
integer c2_g_day
Definition: ch_xfr.f90:46
integer *1, dimension(:,:), allocatable c2_coast_sfc
Definition: ch_xfr.f90:32
integer c2_sensor_id
Definition: ch_xfr.f90:50
integer *1, dimension(:,:), allocatable c2_glint
Definition: ch_xfr.f90:33
integer *1, dimension(:,:), allocatable c2_lnd_snow
Definition: ch_xfr.f90:35
subroutine zero_surround()
Definition: ch_cld_sci.f90:320
real, dimension(:,:,:), allocatable c2_prof_p
Definition: ch_xfr.f90:16
integer tbl_hi_abs_ice
Definition: ch_xfr.f90:66
real(single), dimension(:,:), allocatable cloud_top_temperature
real, dimension(:,:,:), allocatable prd_out_refl_loc_1621
Definition: ch_xfr.f90:64
integer *1, dimension(:,:), allocatable c2_lnd_sfc
Definition: ch_xfr.f90:33
integer *1, dimension(:,:), allocatable c2_tst_vis_ratio
Definition: ch_xfr.f90:37
integer tbl_lo_abs_ice
Definition: ch_xfr.f90:66
real(single), dimension(:,:), allocatable optical_thickness_1621_final
Definition: core_arrays.f90:42
integer *1, dimension(:,:), allocatable c2_appl_vis_nir_ratio
Definition: ch_xfr.f90:40
integer *1, dimension(:,:), allocatable c2_clr_95
Definition: ch_xfr.f90:30
integer c2_scan
Definition: ch_xfr.f90:46
type(qualityanalysis), dimension(:,:), allocatable processing_information
integer tbl_lo_nabs
Definition: ch_xfr.f90:66
integer *1, dimension(:,:), allocatable c2_lnd_no_snow
Definition: ch_xfr.f90:35
real(single), dimension(:,:), allocatable liquid_water_path_16
Definition: core_arrays.f90:67
integer c2_nlvl_model
Definition: ch_xfr.f90:46
integer *1, dimension(:,:), allocatable c2_ocean_no_snow
Definition: ch_xfr.f90:34
integer *1, dimension(:,:), allocatable c2_sfc_lvl
Definition: ch_xfr.f90:26
real(single), dimension(:,:), allocatable effective_radius_22_final
Definition: core_arrays.f90:45
real(single), dimension(:,:), allocatable liquid_water_path_22
Definition: core_arrays.f90:66
integer oci_id
Definition: ch_xfr.f90:52
real, dimension(:,:), allocatable c2_alt_o3
Definition: ch_xfr.f90:20
real(single), dimension(:,:), allocatable optical_thickness_16_final
Definition: core_arrays.f90:40
real(single), dimension(:,:), allocatable effective_radius_16_final
Definition: core_arrays.f90:43
integer c2_nbnd
Definition: ch_xfr.f90:46
integer tbl_hi_nabs
Definition: ch_xfr.f90:66
integer tbl_lo_abs
Definition: ch_xfr.f90:66
integer *1, dimension(:,:), allocatable c2_clr_66
Definition: ch_xfr.f90:30
real, dimension(:,:), allocatable c2_ctt
Definition: ch_xfr.f90:24
integer c2_st_samp
Definition: ch_xfr.f90:46
integer tbl_hi_abs
Definition: ch_xfr.f90:66
real, dimension(:,:), allocatable c2_snow_frac
Definition: ch_xfr.f90:19
real, dimension(:,:), allocatable c2_wind
Definition: ch_xfr.f90:18
real, dimension(:,:), allocatable c2_senz
Definition: ch_xfr.f90:12
real, dimension(:,:,:), allocatable prd_out_refl_loc_2100
Definition: ch_xfr.f90:62
real, dimension(:,:,:), allocatable c2_prof_t
Definition: ch_xfr.f90:15
real(single), dimension(:,:), allocatable effective_radius_1621_final
Definition: core_arrays.f90:47
real, dimension(:,:), allocatable c2_sola
Definition: ch_xfr.f90:13
subroutine driver_mod_pr06od()
real, dimension(:), allocatable c2_spec_unc
Definition: ch_xfr.f90:11
real, dimension(:,:), allocatable c2_tot_o3
Definition: ch_xfr.f90:17
integer *1, dimension(:), allocatable c2_cmp_there
Definition: ch_xfr.f90:43
real, dimension(:,:), allocatable c2_cth
Definition: ch_xfr.f90:24
integer *1, dimension(:,:), allocatable c2_conf_cld
Definition: ch_xfr.f90:29
real, dimension(:,:,:), allocatable c2_refl
Definition: ch_xfr.f90:10
real(single), dimension(:,:), allocatable effective_radius_21_final
Definition: core_arrays.f90:44