NASA Logo
Ocean Color Science Software

ocssw V2022
corescience_module.f90
Go to the documentation of this file.
2 
3  implicit none
4  private
5 
6  public :: corescience
7 
8  contains
9 
10 
11 subroutine corescience(xpoint,ypoint,process,measurements, Tc_liquid, &
12  Tc_ice, debug, na_band_used, nearest_liq, nearest_ice, RSS_liq, RSS_ice, &
13  alt_ray_liq, alt_ray_ice, do_retrievals_liq, do_retrievals_ice, status)
14 
15 ! Arguments
16 ! xpoint IN pixel location in radiance chunk array
17 ! ypoint IN line location in radiance chunk array
18 ! process IN structure of pixel's cloud status, type processflag
19 ! measurements IN vector of top of cloud reflectances
20 ! Tc_liquid IN/OUT cloud temperature for liquid cloud gotten from
21 ! cloud_top_temperature (IN) and may be re-calculated (OUT) - only
22 ! used for 3.7 um band
23 ! Tc_ice IN/OUT cloud temperature for ice cloud - same in/out as above
24 ! - only for 3.7 um band
25 ! debug IN debugging switch USED?
26 ! na_band_used OUT non-absorbing band used, 650 (land), 860 (water), or
27 ! 1240 (snow)
28 ! nearest_liq OUT size 4 array of 1 - nearest used, 0 - not used
29 ! for use of absorbing band: 1 - 1.6 um, 2 - 2.1 um, 3 - 3.7 | 2.2 um,
30 ! 4 - 1.6 / 2.1 um bands
31 ! nearest_ice OUT size 4 array of 1 - nearest used, 0 - not used
32 ! RSS_liq IN?OUT size 4 array of some root sum square uncertainty for
33 ! each of the abs bands - liquid --OR-- a distance away from
34 ! the closest solution in the alternate solution subplot/logic
35 ! RSS_ice IN?OUT size 4 array of some root sum square uncertainty for
36 ! each of the abs bands - ice
37 ! --- WDR lets try understanding std algorithm before alternate one (ASL)
38 ! alt_ray_liq
39 ! alt_ray_ice
40 ! do_retrievals_liq
41 ! do_retrievals_ice
42 ! status
43 
44  use generalauxtype
45  use core_arrays
46  use libraryarrays
56 
57  use specific_other
58  ! WDR for output of the some diagnostics
59  use ch_xfr, only: prd_out_refl_loc_2100, &
65 
66  implicit none
67 
68  logical, intent(in) :: debug
69  type(processflag), intent(in) :: process
70  real, dimension(:), intent(in) :: measurements
71  logical, dimension(:), intent(in) :: do_retrievals_liq, do_retrievals_ice
72  real, intent(inout) :: tc_liquid, tc_ice, alt_ray_liq, alt_ray_ice
73  integer, intent(in) :: xpoint,ypoint
74  integer, intent(out) :: status
75  integer, intent(out) :: na_band_used
76 
77  integer, dimension(:), intent(inout) :: nearest_liq, nearest_ice
78  real, dimension(:), intent(inout) :: rss_liq, rss_ice
79 
80  type(cloudphase) :: local_process
81  integer :: nonabsorbingband_index, absorbingband_index, maxradii, &
82  lib_vnir_index
83  real :: thermal_trans_1way, thermal_trans_2way
84  real, allocatable :: optical_thickness_vector(:), residual(:)
85  real :: retrievalopticalthickness, retrievalopticalthickness16, &
86  retrievalopticalthickness37, retrievalopticalthickness1621,&
87  retrievalradius21, retrievalradius16, retrievalradius37, &
88  retrievalradius1621, retrievalradius22, retrievalopticalthickness22
89 
90  integer :: idx16, idx21, idx37, channel_37, idx11, i, idx_alb37, &
91  idx_alb16, quality_in, idx22
92  real :: curtc, newtc, ray_temp
93  logical :: prn, use_nearest
94  integer :: iii, idx1621(2)
95 
96  status = 0
97 
98  thermal_trans_1way = thermal_correction_oneway(1) ! only the 3.7um is used
99  thermal_trans_2way = thermal_correction_twoway(1) ! only the 3.7um is used
100 
101  ! local band naming - var names are mostly obvious, more || less
102  call get_band_idx(idx16, idx21, idx37, channel_37, idx11, idx_alb37, &
103  idx_alb16)
104  ! for OCI, get the 2.2
105  if( ( c2_sensor_id == oci_id ) .OR. ( c2_sensor_id == ocis_id ) ) &
106  idx22 = band_0226 - 1
107 
108  nonabsorbingband_index = band_0065
109  processing_information(xpoint, ypoint)%band_used_for_optical_thickness = 0
110  !
111  ! set non-absorbing band based on the surface type
112  !
113  if (process%ocean_surface .or. process%coastal_surface) then!{
114  nonabsorbingband_index = band_0086
115  endif!}
116 
117  if (process%land_surface .or. process%desert_surface) then!{
118  nonabsorbingband_index = band_0065
119  endif!}
120  if (process%snowice_surface) then!{
121  ! 1.2 um
122  nonabsorbingband_index = band_0124
123  endif!}
124 
125  if (platform_name == "AVIRIS") then
126  if (nonabsorbingband_index == 4) &
127  processing_information(xpoint, ypoint)%band_used_for_optical_thickness = 3
128  if (nonabsorbingband_index == 3) &
129  processing_information(xpoint, ypoint)%band_used_for_optical_thickness = 2
130  if (nonabsorbingband_index == 1) &
131  processing_information(xpoint, ypoint)%band_used_for_optical_thickness = 1
132  else
133  if (nonabsorbingband_index > 3) then
134  processing_information(xpoint, ypoint)%band_used_for_optical_thickness &
135  = nonabsorbingband_index-1
136  else
137  processing_information(xpoint, ypoint)%band_used_for_optical_thickness &
138  = nonabsorbingband_index
139  endif
140  endif
141 
142  na_band_used = nonabsorbingband_index
143 
144  local_process%watercloud = 0
145  local_process%icecloud = 0
146 
147  ! initialize all retrievals to fillvalue
148  retrievalopticalthickness = fillvalue_real
149  retrievalopticalthickness1621 = fillvalue_real
150  retrievalopticalthickness16 = fillvalue_real
151  retrievalopticalthickness37 = fillvalue_real
152  retrievalradius16 = fillvalue_real
153  retrievalradius21 = fillvalue_real
154  retrievalradius1621 = fillvalue_real
155  retrievalradius37 = fillvalue_real
156  if( ( c2_sensor_id == oci_id ) .OR. ( c2_sensor_id == ocis_id ) )then
157  retrievalradius22 = fillvalue_real
158  retrievalopticalthickness22 = fillvalue_real
159  endif
160 
161  nearest_liq = 0
162  nearest_ice = 0
163 
164  local_process%watercloud = 1
165 
166  ! set surface QA
167 
168  ! if the measured surface reflectance in the non absorbing bands is
169  ! greater than the surface albedo , and we have a cloudy pixel then
170  ! we can try a retrieval
171 
172  if (process%cloudobserved ) then!{
173 
174  ! check for .86um saturation
175  ! if the .86um saturates (wisconsin reader sets band_meas to neg.
176  ! if saturation) or there can be bad noisy neg. reflectance on low end,
177  ! try the other .65um band (and change the QA flag)
178  if (nonabsorbingband_index == band_0086 .and. &
179  measurements(nonabsorbingband_index) < 0.) then!{
180  if (measurements(band_0065) > 0.) nonabsorbingband_index = band_0065
181  na_band_used = nonabsorbingband_index
182  processing_information(xpoint, ypoint)%band_used_for_optical_thickness &
183  = nonabsorbingband_index
184  endif!}
185 
186  nonabsorbingband_index = na_band_used
187 
188  ! WDR save the used reflectance - we still need the table ranges found
189  prd_out_refl_loc_2100(xpoint,ypoint,refl_nabs) = &
190  measurements(nonabsorbingband_index)
191  prd_out_refl_loc_1600(xpoint,ypoint,refl_nabs) = &
192  measurements(nonabsorbingband_index)
193  prd_out_refl_loc_1621(xpoint,ypoint,refl_nabs) = &
194  measurements(band_0163)
195  prd_out_refl_loc_2200(xpoint,ypoint,refl_nabs) = &
196  measurements(nonabsorbingband_index)
197 !
198  prd_out_refl_loc_2100(xpoint,ypoint,refl_abs) = &
199  measurements(band_0213)
200  prd_out_refl_loc_1600(xpoint,ypoint,refl_abs) = &
201  measurements(band_0163)
202  prd_out_refl_loc_1621(xpoint,ypoint,refl_abs) = &
203  measurements(band_0213)
204  prd_out_refl_loc_2200(xpoint,ypoint,refl_abs) = &
205  measurements(band_0226)
206 
207  ! this is for rotten snow albedoes
208  if (albedo_real4( nonabsorbingband_index) < 0. .or. &
209  albedo_real4( band_0163) < 0. .or. &
210  albedo_real4( band_0213) < 0.) return
211 
212  ! if the shortwave saturates, try the other shortwave band
213  allocate(optical_thickness_vector(number_waterradii))
214  allocate(residual(number_waterradii))
215  residual = 0 ! WDR UIV
216 
217  allocate(refliba(number_taus+1, number_waterradii), &
219  refliba = 0 ! WDR-UIV
220  reflibb = 0 ! WDR-UIV
222 
223  rad37lib = -999.0
224 
225  lib_vnir_index = nonabsorbingband_index
226  !
227  ! For water clouds, for the non-absorbing reflectance, find the
228  ! tau (optical thickness) values in the table for every eff. radius
229  !
230  call vis_nonabsorbing_science(measurements(nonabsorbingband_index), &
231  lib_vnir_index, albedo_real4(nonabsorbingband_index), &
234  int_reflectance_water, sensor_zenith_angle(xpoint,ypoint), &
235  solar_zenith_angle(xpoint,ypoint), &
236  relative_azimuth_angle(xpoint,ypoint), &
237  cloud_top_pressure(xpoint,ypoint), local_process, &
238  optical_thickness_vector)
239 
240  ! an option to knock out the SWIR channels completely so that we only
241  ! do optical thickness retrieval assuming an effective radius
242  ! This appears to handle the other COT band values too
243 #if NOSWIR
244  ! Retrieval using assumed effective radius
245  retrievalradius16 = 14.0
246  call solve_retrieval_noswir(optical_thickness_vector, water_radii, &
247  retrievalradius16, retrievalopticalthickness16)
248 
249  retrievalopticalthickness = retrievalopticalthickness16
250  retrievalradius21 = retrievalradius16
251 
252  ! Retrieval using assumed effective radius, minus 1-sigma CER
253  retrievalradius16 = 6.0
254  call solve_retrieval_noswir(optical_thickness_vector, water_radii, &
255  retrievalradius16, retrievalopticalthickness16)
256 
257  optical_thickness_37_final(xpoint,ypoint) = &
258  abs(retrievalopticalthickness-retrievalopticalthickness16)
259 
260  ! Retrieval using assumed effective radius, plus 1-sigma CER
261  retrievalradius16 = 22.0
262  call solve_retrieval_noswir(optical_thickness_vector, water_radii, &
263  retrievalradius16, retrievalopticalthickness16)
264 
265  optical_thickness_37_final(xpoint,ypoint) = &
266  optical_thickness_37_final(xpoint,ypoint) + &
267  abs(retrievalopticalthickness-retrievalopticalthickness16)
268 
269  ! Calculate mean COT error due to unknown CER
270  optical_thickness_37_final(xpoint,ypoint) = &
271  optical_thickness_37_final(xpoint,ypoint) / 2.0
272 
273  retrievalopticalthickness1621 = fillvalue_real
274  retrievalradius1621 = fillvalue_real
275 
276 #else
277  ! below is with SWIR
278  !
279  ! For 1.6 um, continue and find the re, tau using the absorbing
280  ! band reflectance for water clouds
281  !
282  if (do_retrievals_liq(1)) then
283  absorbingband_index = band_0163 ! 1.6um band
284  if(measurements(absorbingband_index) > 0. ) then!{
285 
286  call vis_absorbing_science(optical_thickness_vector, &
287  measurements(absorbingband_index), idx16, &
288  albedo_real4(idx_alb16), library_taus, water_radii, &
291  residual,maxradii, debug)
292 
293  if (maxradii > 2) then!{
294  call solveretrieval(residual(1:maxradii), &
295  optical_thickness_vector(1:maxradii), &
296  water_radii(1:maxradii), retrievalradius16, &
297  retrievalopticalthickness16, debug, use_nearest, quality_in)
298  if (use_nearest) then
299  nearest_liq(1) = 1
300  call solveretrieval_nearest(xpoint,ypoint, &
301  measurements(nonabsorbingband_index), &
302  measurements(absorbingband_index), &
303  (/nonabsorbingband_index, absorbingband_index/), &
304  water_radii, retrievalopticalthickness16, &
305  retrievalradius16, rss_liq(re16), &
306  .true., ray_temp, quality_in )
307  endif
308 
309  else!}{
310  retrievalradius16 = fillvalue_real
311  retrievalopticalthickness16 = fillvalue_real
312  endif!}
313 
314  else!}{
315  retrievalradius16 = fillvalue_real
316  retrievalopticalthickness16 = fillvalue_real
317  endif!}
318  else
319  retrievalradius16 = fillvalue_real
320  retrievalopticalthickness16 = fillvalue_real
321  endif
322  !
323  ! WDR get the table range
324  prd_out_refl_loc_1600(xpoint,ypoint,tbl_lo_nabs) = minval(refliba)
325  prd_out_refl_loc_1600(xpoint,ypoint,tbl_hi_nabs) = maxval(refliba)
326  prd_out_refl_loc_1600(xpoint,ypoint,tbl_lo_abs) = minval(reflibb)
327  prd_out_refl_loc_1600(xpoint,ypoint,tbl_hi_abs) = maxval(reflibb)
328  !
329  ! For 2.1 um, continue and find the re, tau using the absorbing
330  ! band reflectance for water clouds
331  !
332 
333  if (do_retrievals_liq(2) ) then
334  absorbingband_index = band_0213 ! 2.1um
335 
336  call vis_absorbing_science(optical_thickness_vector, &
337  measurements(absorbingband_index), idx21, &
338  albedo_real4(absorbingband_index), library_taus, water_radii, &
341  residual, maxradii, debug)
342 
343  if (maxradii > 2) then!{
344  call solveretrieval(residual(1:maxradii), &
345  optical_thickness_vector(1:maxradii), water_radii(1:maxradii), &
346  retrievalradius21, retrievalopticalthickness, &
347  debug, use_nearest, quality_in)
348 
349  if (use_nearest) then
350  nearest_liq(2) = 1
351  call solveretrieval_nearest(xpoint,ypoint, &
352  measurements(nonabsorbingband_index), &
353  measurements(absorbingband_index), &
354  (/nonabsorbingband_index, absorbingband_index/), &
355  water_radii, retrievalopticalthickness, retrievalradius21, &
356  rss_liq(re21), .true., alt_ray_liq, quality_in)
357  endif
358 
359  else!}{
360  retrievalradius21 = fillvalue_real
361  retrievalopticalthickness = fillvalue_real
362  endif!}
363  else
364  retrievalradius21 = fillvalue_real
365  retrievalopticalthickness = fillvalue_real
366  endif
367  !
368  ! WDR get the table range
369  prd_out_refl_loc_2100(xpoint,ypoint,tbl_lo_nabs) = minval(refliba)
370  prd_out_refl_loc_2100(xpoint,ypoint,tbl_hi_nabs) = maxval(refliba)
371  prd_out_refl_loc_2100(xpoint,ypoint,tbl_lo_abs) = minval(reflibb)
372  prd_out_refl_loc_2100(xpoint,ypoint,tbl_hi_abs) = maxval(reflibb)
373  !
374 
375 #endif
376  ! WDR best do 2.2 here - same as 1.6 and 2.1
377  if( ( c2_sensor_id == oci_id ) .OR. ( c2_sensor_id == ocis_id ) )then
378  if (do_retrievals_liq(3) ) then
379  absorbingband_index = band_0226 ! 2.2um
380 
381  call vis_absorbing_science(optical_thickness_vector, &
382  measurements(absorbingband_index), idx22, &
383  albedo_real4(absorbingband_index), library_taus, water_radii, &
386  residual, maxradii, debug)
387 
388  if (maxradii > 2) then!{
389  call solveretrieval(residual(1:maxradii), &
390  optical_thickness_vector(1:maxradii), water_radii(1:maxradii), &
391  retrievalradius22, retrievalopticalthickness22, &
392  debug, use_nearest, quality_in)
393 
394  if (use_nearest) then
395  nearest_liq(3) = 1
396  call solveretrieval_nearest(xpoint,ypoint, &
397  measurements(nonabsorbingband_index), &
398  measurements(absorbingband_index), &
399  (/nonabsorbingband_index, absorbingband_index/), &
400  water_radii, retrievalopticalthickness22, retrievalradius22, &
401  rss_liq(re22), .true., alt_ray_liq, quality_in)
402  endif
403 
404  else!}{
405  retrievalradius22 = fillvalue_real
406  retrievalopticalthickness22 = fillvalue_real
407  endif!}
408  else
409  retrievalradius22 = fillvalue_real
410  retrievalopticalthickness22 = fillvalue_real
411  endif
412  !
413  ! WDR get the table range
414  prd_out_refl_loc_2200(xpoint,ypoint,tbl_lo_nabs) = minval(refliba)
415  prd_out_refl_loc_2200(xpoint,ypoint,tbl_hi_nabs) = maxval(refliba)
416  prd_out_refl_loc_2200(xpoint,ypoint,tbl_lo_abs) = minval(reflibb)
417  prd_out_refl_loc_2200(xpoint,ypoint,tbl_hi_abs) = maxval(reflibb)
418  else ! END 2.2 um
419  !
420  !
421  ! For 3.7 um, continue and find the re, tau using the absorbing
422  ! band reflectance for water clouds
423  !
424  if (do_retrievals_liq(3)) then
425  absorbingband_index = band_0370 ! 3.7um
426 
427  ! sep, 4 May: absorbingband_index - 1 being passed in for all libarary
428  ! indices because
429  ! band_0370 = 7 and library index corresponding to this band
430  ! is 6. In the libraries, band_0935 is index 7.
431  ! wind, 7 Dec: absorbingband_index - 1 has to be passed
432  !in for albedo as well.
433 
434  call nir_absorbing_science(platform_name, optical_thickness_vector, &
435  measurements(absorbingband_index), idx37, albedo_real4( idx_alb37), &
436  xpoint, ypoint, cloud_top_temperature(xpoint, ypoint), &
437  thermal_trans_1way, thermal_trans_2way, library_taus, &
441  int_surface_emissivity_water, residual, maxradii, channel_37, &
443  sigma_r37_pw_liq,debug)
444 
445  if ( maxradii > 2 ) then!{
446  call solveretrieval(residual(1:maxradii), &
447  optical_thickness_vector(1:maxradii), &
448  water_radii(1:maxradii), retrievalradius37, &
449  retrievalopticalthickness37, debug, use_nearest, quality_in)
450 
451  if (use_nearest) then
452  nearest_liq(3) = 1
453  call solveretrieval_nearest(xpoint,ypoint, &
454  measurements(nonabsorbingband_index), &
455  measurements(absorbingband_index), &
456  (/nonabsorbingband_index, absorbingband_index/), &
457  water_radii, retrievalopticalthickness37, retrievalradius37, &
458  rss_liq(re37), .true., ray_temp, quality_in, &
459  ch37_idx = idx37, ctopt = cloud_top_temperature(xpoint, ypoint), &
460  ch37_num =channel_37 , platformname= platform_name)
461  endif
462 
463  else!}{
464  retrievalradius37 = fillvalue_real
465  retrievalopticalthickness37 = fillvalue_real
466  endif!}
467  !
468  ! now we iterate the retrieval. I believe the emission in 3.7 is
469  ! enough to require this iteration
470  !
471 #if !AMS_INST
472 
473  ! if the cloud too thick, don't bother, the iteration will not
474  ! accomplish absolutely anything
475  if ( retrievalopticalthickness37 > 0. &
476  .and. retrievalradius37 > 0. .and. &
477  irw_temperature(xpoint, ypoint) > 0. .and. nearest_liq(3) == 0 ) then
478 
479  ! if the cloud is too thick, then don't bother doing anything, but
480  ! still set the temperature
481  if (retrievalopticalthickness37 >= 10.) then
482  tc_liquid = cloud_top_temperature(xpoint, ypoint)
483  else
484  curtc = irw_temperature(xpoint, ypoint)
485 
486  do i=1, 5
487  if (retrievalopticalthickness37 < 0. .or. &
488  retrievalradius37 < 0. ) then
489  retrievalradius37 = fillvalue_real
490  retrievalopticalthickness37 = fillvalue_real
491  tc_liquid = curtc
492  exit
493  endif
494 
496  irw_temperature(xpoint, ypoint), &
497  surface_temperature(xpoint, ypoint), &
498  1.- surface_emissivity_land(xpoint, ypoint, 2), idx11, &
499  retrievalopticalthickness37, retrievalradius37, library_taus, &
502  int_surface_emissivity_water, newtc, .false.)
503 
504  if (newtc < 0.) then
505  tc_liquid = curtc
506  exit
507  endif
508 
510  optical_thickness_vector, &
511  measurements(absorbingband_index), idx37, &
512  albedo_real4( idx_alb37), &
513  xpoint, ypoint, newtc, thermal_trans_1way, thermal_trans_2way, &
518  residual, maxradii, channel_37, emission_uncertainty_pw_liq, &
520 
521  if ( maxradii > 2 ) then!{
522  call solveretrieval(residual(1:maxradii), &
523  optical_thickness_vector(1:maxradii), water_radii(1:maxradii), &
524  retrievalradius37, retrievalopticalthickness37, &
525  debug, use_nearest, quality_in)
526  if (use_nearest) then
527  nearest_liq(3) = 1
528  call solveretrieval_nearest(xpoint,ypoint, &
529  measurements(nonabsorbingband_index), &
530  measurements(absorbingband_index), &
531  (/nonabsorbingband_index, absorbingband_index/), &
532  water_radii, retrievalopticalthickness37, retrievalradius37, &
533  rss_liq(re37), .true., ray_temp, quality_in, &
534  ch37_idx = idx37, &
535  ctopt = cloud_top_temperature(xpoint, ypoint), &
536  ch37_num =channel_37 , platformname= platform_name)
537 
538  endif
539  else!}{
540  retrievalradius37 = fillvalue_real
541  retrievalopticalthickness37 = fillvalue_real
542  exit
543  endif!}
544 
545  if ( abs(curtc - newtc) < 0.01 .or. retrievalradius37 < 0.) then
546  tc_liquid = newtc
547  exit
548  endif
549 
550  curtc = newtc
551  end do
552  endif
553 
554  endif
555 #endif
556  else
557  retrievalradius37 = fillvalue_real
558  retrievalopticalthickness37 = fillvalue_real
559  endif ! END 3.7
560  endif
561 
562 #if NOSWIR
563  retrievalradius1621 = fillvalue_real
564  retrievalopticalthickness1621 = fillvalue_real
565 #else
566  !
567  ! For 1.6 and 2.1 um, continue and find the re, tau using the absorbing
568  ! band reflectance for water clouds
569  !
570  if (do_retrievals_liq(4)) then
571  retrievalradius1621 = fillvalue_real
572  retrievalopticalthickness1621 = fillvalue_real
573 
574  if( (process%ocean_surface .or. process%snowice_surface) .and. &
575  measurements(band_0163) > 0. ) then!{
576  nonabsorbingband_index = band_0163
577  absorbingband_index = band_0213
578 
579  ! get the COT for all radii that gives the non-abs refl at this point
580  call vis_nonabsorbing_science(measurements(nonabsorbingband_index), &
581  idx16, albedo_real4(nonabsorbingband_index), library_taus, &
584  sensor_zenith_angle(xpoint,ypoint), &
585  solar_zenith_angle(xpoint,ypoint), &
586  relative_azimuth_angle(xpoint,ypoint), &
587  cloud_top_pressure(xpoint,ypoint), local_process, &
588  optical_thickness_vector)
589 
590  call vis_absorbing_science(optical_thickness_vector, &
591  measurements(absorbingband_index), idx21, &
592  albedo_real4(absorbingband_index), library_taus, water_radii, &
595  residual,maxradii, debug)
596  if (maxradii > 2) then!{
597  call solveretrieval(residual(1:maxradii), &
598  optical_thickness_vector(1:maxradii), water_radii(1:maxradii), &
599  retrievalradius1621, retrievalopticalthickness1621, debug, &
600  use_nearest, quality_in)
601  if (use_nearest) then
602  nearest_liq(4) = 1
603 
604  call solveretrieval_nearest(xpoint,ypoint, &
605  measurements(nonabsorbingband_index), &
606  measurements(absorbingband_index), &
607  (/nonabsorbingband_index, absorbingband_index/), &
608  water_radii, retrievalopticalthickness1621, &
609  retrievalradius1621, rss_liq(re1621), &
610  .true., ray_temp, quality_in)
611  endif
612  endif!}
613 
614  endif!}
615  else
616  retrievalradius1621 = fillvalue_real
617  retrievalopticalthickness1621 = fillvalue_real
618  endif
619  !
620  ! WDR get the table range
621  prd_out_refl_loc_1621(xpoint,ypoint,tbl_lo_nabs) = minval(refliba)
622  prd_out_refl_loc_1621(xpoint,ypoint,tbl_hi_nabs) = maxval(refliba)
623  prd_out_refl_loc_1621(xpoint,ypoint,tbl_lo_abs) = minval(reflibb)
624  prd_out_refl_loc_1621(xpoint,ypoint,tbl_hi_abs) = maxval(reflibb)
625  !
626 #endif
627  !
628  ! End setting of water path for water phase clouds
629  ! Repeat entire process for ice clouds
630  !
631  deallocate(rad37lib)
632 
633  deallocate(optical_thickness_vector, residual, refliba, reflibb)
634 
635  optical_thickness_liquid = retrievalopticalthickness
636  optical_thickness_1621_liquid = retrievalopticalthickness1621
637  optical_thickness_16_liquid = retrievalopticalthickness16
638  optical_thickness_37_liquid = retrievalopticalthickness37
639  effective_radius_16_liquid = retrievalradius16
640  effective_radius_21_liquid = retrievalradius21
641  effective_radius_1621_liquid = retrievalradius1621
642  effective_radius_37_liquid = retrievalradius37
643  optical_thickness_22_liquid = retrievalopticalthickness22
644  effective_radius_22_liquid = retrievalradius22
645 
646  ! This was commented out in initiall chimaera code
647  ! if (nearest_liq(1) == 1) effective_radius_16_liquid(xpoint, ypoint) = &
648  ! fillvalue_real
649  ! if (nearest_liq(2) == 1) then
650  ! effective_radius_21_liquid(xpoint, ypoint) = fillvalue_real
651  ! optical_thickness_liquid(xpoint, ypoint) = fillvalue_real
652  ! endif
653  ! if (nearest_liq(3) == 1) effective_radius_37_liquid(xpoint, ypoint) = &
654  ! fillvalue_real
655  ! if (nearest_liq(4) == 1) then
656  ! effective_radius_1621_liquid(xpoint, ypoint) = fillvalue_real
657  ! optical_thickness_1621_liquid(xpoint, ypoint) = fillvalue_real
658  ! endif
659  !
660  !
661  ! now reset everything so the ice cloud can be processed
662  !
663  local_process%watercloud = 0
664  local_process%icecloud = 1
665 
666  nonabsorbingband_index = na_band_used
667  absorbingband_index = band_0163
668 
669  allocate(optical_thickness_vector(number_iceradii))
670  allocate(residual(number_iceradii))
671  residual = 0 ! WDR UIV
672  allocate(refliba(number_taus+1, number_iceradii), &
674  allocate(rad37lib(number_taus+1, number_iceradii))
675  ! WDR uiv sort of
676  refliba = 0
677  reflibb = 0
678 
679  rad37lib = -999.0
680 
681  ! initialize all retrievals to fillvalue
682  retrievalopticalthickness = fillvalue_real
683  retrievalopticalthickness1621 = fillvalue_real
684  retrievalopticalthickness16 = fillvalue_real
685  retrievalopticalthickness37 = fillvalue_real
686  retrievalradius16 = fillvalue_real
687  retrievalradius21 = fillvalue_real
688  retrievalradius1621 = fillvalue_real
689  retrievalradius37 = fillvalue_real
690  retrievalopticalthickness22 = fillvalue_real
691  retrievalradius22 = fillvalue_real
692 
693  lib_vnir_index = nonabsorbingband_index
694  !
695  ! For ice clouds, find the tau (optical thickness) values in the
696  ! table for the non-absorbing reflectance
697  !
698  call vis_nonabsorbing_science(measurements(nonabsorbingband_index), &
699  lib_vnir_index, albedo_real4(nonabsorbingband_index), &
702  sensor_zenith_angle(xpoint,ypoint), &
703  solar_zenith_angle(xpoint,ypoint), &
704  relative_azimuth_angle(xpoint,ypoint), &
705  cloud_top_pressure(xpoint,ypoint), &
706  local_process, optical_thickness_vector)
707 
708  ! code for work without SWIR bands
709 #if NOSWIR
710  ! Retrieval using assumed effective radius
711  retrievalradius16 = 30.0
712  call solve_retrieval_noswir(optical_thickness_vector, ice_radii, &
713  retrievalradius16, retrievalopticalthickness16)
714 
715  retrievalopticalthickness = retrievalopticalthickness16
716  retrievalradius21 = retrievalradius16
717 
718  ! Retrieval using assumed effective radius, minus 1-sigma CER
719  retrievalradius16 = 13.0
720  call solve_retrieval_noswir(optical_thickness_vector, ice_radii, &
721  retrievalradius16, retrievalopticalthickness16)
722 
723  optical_thickness_1621_final(xpoint,ypoint) = &
724  abs(retrievalopticalthickness-retrievalopticalthickness16)
725 
726  ! Retrieval using assumed effective radius, plus 1-sigma CER
727  retrievalradius16 = 37.0
728  call solve_retrieval_noswir(optical_thickness_vector, ice_radii, &
729  retrievalradius16, retrievalopticalthickness16)
730 
731  optical_thickness_1621_final(xpoint,ypoint) = &
732  optical_thickness_1621_final(xpoint,ypoint) + &
733  abs(retrievalopticalthickness-retrievalopticalthickness16)
734 
735  ! Calculate mean COT error due to unknown CER
736  optical_thickness_1621_final(xpoint,ypoint) = &
737  optical_thickness_1621_final(xpoint,ypoint) / 2.0
738 
739  retrievalopticalthickness1621 = fillvalue_real
740  retrievalradius1621 = fillvalue_real
741 #else
742  !
743  ! For 1.6 um, continue and find the re, tau using the absorbing
744  ! band reflectance for ice clouds
745  !
746  if (do_retrievals_ice(1)) then
747  if (measurements(absorbingband_index) > 0.) then!{
748  call vis_absorbing_science(optical_thickness_vector, &
749  measurements(absorbingband_index), idx16, albedo_real4(idx_alb16), &
752  int_reflectance_ice, residual,maxradii, debug)
753  if ( maxradii > 2) then!{
754  call solveretrieval(residual(1:maxradii), &
755  optical_thickness_vector(1:maxradii), &
756  ice_radii(1:maxradii), retrievalradius16, &
757  retrievalopticalthickness16, debug, use_nearest, quality_in)
758  if (use_nearest) then
759  nearest_ice(1) = 1
760 
761  call solveretrieval_nearest(xpoint,ypoint, &
762  measurements(nonabsorbingband_index), &
763  measurements(absorbingband_index), &
764  (/nonabsorbingband_index, absorbingband_index/), &
765  ice_radii, retrievalopticalthickness16, retrievalradius16, &
766  rss_ice(re16), .false., ray_temp, quality_in)
767  endif
768  else!}{
769  retrievalradius16 = fillvalue_real
770  retrievalopticalthickness16 = fillvalue_real
771  endif!}
772 
773  else!}{
774  retrievalradius16 = fillvalue_real
775  retrievalopticalthickness16 = fillvalue_real
776  endif!}
777  else
778  retrievalradius16 = fillvalue_real
779  retrievalopticalthickness16 = fillvalue_real
780  endif
781  !
782  ! WDR get the table range
783  prd_out_refl_loc_1600(xpoint,ypoint,tbl_lo_nabs_ice) = minval(refliba)
784  prd_out_refl_loc_1600(xpoint,ypoint,tbl_hi_nabs_ice) = maxval(refliba)
785  prd_out_refl_loc_1600(xpoint,ypoint,tbl_lo_abs_ice) = minval(reflibb)
786  prd_out_refl_loc_1600(xpoint,ypoint,tbl_hi_abs_ice) = maxval(reflibb)
787  !
788  !
789  ! For 2.1 um, continue and find the re, tau using the absorbing
790  ! band reflectance
791  !
792  if (do_retrievals_ice(2)) then
793  absorbingband_index = band_0213
794  call vis_absorbing_science(optical_thickness_vector, &
795  measurements(absorbingband_index), idx21, &
796  albedo_real4(absorbingband_index), library_taus, ice_radii, &
798  int_reflectance_ice, residual,maxradii, debug)
799 
800  if ( maxradii > 2) then!{
801  call solveretrieval(residual(1:maxradii), &
802  optical_thickness_vector(1:maxradii), ice_radii(1:maxradii), &
803  retrievalradius21, retrievalopticalthickness, &
804  debug, use_nearest, quality_in)
805 
806  if (use_nearest) then
807  nearest_ice(2) = 1
808  call solveretrieval_nearest(xpoint,ypoint, &
809  measurements(nonabsorbingband_index), &
810  measurements(absorbingband_index), &
811  (/nonabsorbingband_index, absorbingband_index/), ice_radii, &
812  retrievalopticalthickness, retrievalradius21, rss_ice(re21), &
813  .false., alt_ray_ice, quality_in)
814  endif
815  else!}{
816  retrievalradius21 = fillvalue_real
817  retrievalopticalthickness = fillvalue_real
818  endif!}
819  else
820  retrievalradius21 = fillvalue_real
821  retrievalopticalthickness = fillvalue_real
822  endif
823  !
824  ! WDR get the table range
825  prd_out_refl_loc_2100(xpoint,ypoint,tbl_lo_nabs_ice) = minval(refliba)
826  prd_out_refl_loc_2100(xpoint,ypoint,tbl_hi_nabs_ice) = maxval(refliba)
827  prd_out_refl_loc_2100(xpoint,ypoint,tbl_lo_abs_ice) = minval(reflibb)
828  prd_out_refl_loc_2100(xpoint,ypoint,tbl_hi_abs_ice) = maxval(reflibb)
829  !
830 #endif
831  ! the 2.2 um
832  if( ( c2_sensor_id == oci_id ) .OR. ( c2_sensor_id == ocis_id ) )then
833  if (do_retrievals_ice(3)) then
834  absorbingband_index = band_0226
835  call vis_absorbing_science(optical_thickness_vector, &
836  measurements(absorbingband_index), idx22, &
837  albedo_real4(absorbingband_index), library_taus, ice_radii, &
839  int_reflectance_ice, residual,maxradii, debug)
840 
841  if ( maxradii > 2) then!{
842  call solveretrieval(residual(1:maxradii), &
843  optical_thickness_vector(1:maxradii), ice_radii(1:maxradii), &
844  retrievalradius22, retrievalopticalthickness22, &
845  debug, use_nearest, quality_in)
846 
847  if (use_nearest) then
848  nearest_ice(3) = 1
849  call solveretrieval_nearest(xpoint,ypoint, &
850  measurements(nonabsorbingband_index), &
851  measurements(absorbingband_index), &
852  (/nonabsorbingband_index, absorbingband_index/), ice_radii, &
853  retrievalopticalthickness22, retrievalradius22, rss_ice(re22), &
854  .false., alt_ray_ice, quality_in)
855  endif
856  else!}{
857  retrievalradius22 = fillvalue_real
858  retrievalopticalthickness = fillvalue_real
859  endif!}
860  else
861  retrievalradius22 = fillvalue_real
862  retrievalopticalthickness = fillvalue_real
863  endif
864  !
865  ! WDR get the table range
866  prd_out_refl_loc_2200(xpoint,ypoint,tbl_lo_nabs_ice) = minval(refliba)
867  prd_out_refl_loc_2200(xpoint,ypoint,tbl_hi_nabs_ice) = maxval(refliba)
868  prd_out_refl_loc_2200(xpoint,ypoint,tbl_lo_abs_ice) = minval(reflibb)
869  prd_out_refl_loc_2200(xpoint,ypoint,tbl_hi_abs_ice) = maxval(reflibb)
870  else
871  !
872  ! For 3.7 um, continue and find the re, tau using the absorbing
873  ! band reflectance for ice clouds
874  !
875  if (do_retrievals_ice(3)) then
876  absorbingband_index = band_0370
877 
878  ! sep, 4 May: absorbingband_index - 1 being passed in for all
879  ! libarary indices because
880  ! band_0370 = 7 and library index corresponding to this band
881  ! is 6. In the
882  ! libraries, band_0935 is index 7.
883  !
884  ! wind, 7 Dec: absorbingband_index - 1 has to be passed in for
885  ! albedo as well.
886  !
887  call nir_absorbing_science(platform_name, optical_thickness_vector, &
888  measurements(absorbingband_index), idx37, albedo_real4( idx_alb37), &
889  xpoint, ypoint, cloud_top_temperature(xpoint, ypoint), &
890  thermal_trans_1way, thermal_trans_2way, library_taus, &
894  residual, maxradii, channel_37, emission_uncertainty_pw_ice, &
896 
897  if ( maxradii > 2 ) then!{
898  call solveretrieval(residual(1:maxradii), &
899  optical_thickness_vector(1:maxradii), ice_radii(1:maxradii), &
900  retrievalradius37, retrievalopticalthickness37, debug, use_nearest, &
901  quality_in)
902  if (use_nearest) then
903  nearest_ice(3) = 1
904 
905  call solveretrieval_nearest(xpoint,ypoint, &
906  measurements(nonabsorbingband_index), &
907  measurements(absorbingband_index), &
908  (/nonabsorbingband_index, absorbingband_index/), &
909  ice_radii, retrievalopticalthickness37, &
910  retrievalradius37, rss_ice(re37), .false., ray_temp, quality_in, &
911  ch37_idx = idx37, ctopt = cloud_top_temperature(xpoint, ypoint), &
912  ch37_num =channel_37 , platformname= platform_name)
913  endif
914 
915  else!}{
916  retrievalradius37 = fillvalue_real
917  retrievalopticalthickness37 = fillvalue_real
918  endif!}
919 
920  ! now we iterate the retrieval. I believe the emission in 3.7 is
921  ! enough to require this iteration
922 #if !AMS_INST
923  if ( retrievalopticalthickness37 > 0. &
924  .and. retrievalradius37 > 0. .and. &
925  irw_temperature(xpoint, ypoint) > 0. .and. nearest_ice(3) == 0 ) then
926 
927  ! if the cloud is too thick, then don't bother doing anything,
928  ! 1 but still set the temperature
929  if (retrievalopticalthickness37 > 10.) then
930  tc_ice = cloud_top_temperature(xpoint, ypoint)
931  else
932  curtc = irw_temperature(xpoint, ypoint)
933 
934  do i=1, 5
935 
936  if (retrievalopticalthickness37 < 0. .or. &
937  retrievalradius37 < 0. ) then
938  retrievalradius37 = fillvalue_real
939  retrievalopticalthickness37 = fillvalue_real
940  tc_ice = curtc
941  exit
942  endif
943 
945  irw_temperature(xpoint, ypoint), &
946  surface_temperature(xpoint, ypoint), &
947  1.- surface_emissivity_land(xpoint, ypoint, 2), &
948  idx11, retrievalopticalthickness37, retrievalradius37, &
952  newtc, .false.)
953 
954  if (newtc < 0.) then
955  tc_ice = curtc
956  exit
957  endif
958 
960  optical_thickness_vector, measurements(absorbingband_index), &
961  idx37, albedo_real4( idx_alb37), xpoint, ypoint, &
962  newtc, thermal_trans_1way, thermal_trans_2way, &
967  residual, maxradii, channel_37, emission_uncertainty_pw_ice, &
969 
970  if ( maxradii > 2 ) then!{
971  call solveretrieval(residual(1:maxradii), &
972  optical_thickness_vector(1:maxradii), &
973  ice_radii(1:maxradii), retrievalradius37, &
974  retrievalopticalthickness37, &
975  debug, use_nearest, quality_in)
976  if (use_nearest) then
977  nearest_ice(3) = 1
978 
979  call solveretrieval_nearest(xpoint,ypoint, &
980  measurements(nonabsorbingband_index), &
981  measurements(absorbingband_index), &
982  (/nonabsorbingband_index, absorbingband_index/), &
983  ice_radii, retrievalopticalthickness37, retrievalradius37, &
984  rss_ice(re37), .false., ray_temp, quality_in, &
985  ch37_idx = idx37, &
986  ctopt = cloud_top_temperature(xpoint, ypoint), &
987  ch37_num =channel_37 , platformname= platform_name)
988  endif
989 
990  else!}{
991  retrievalradius37 = fillvalue_real
992  retrievalopticalthickness37 = fillvalue_real
993  endif!}
994 
995  if ( abs(curtc - newtc) < 0.01 .or. retrievalradius37 < 0.) then
996  tc_ice = newtc
997  exit
998  endif
999  curtc = newtc
1000  end do
1001 
1002  endif
1003  endif
1004 #endif
1005 
1006  else
1007  retrievalradius37 = fillvalue_real
1008  retrievalopticalthickness37 = fillvalue_real
1009  endif ! end 3.7
1010  endif
1011 
1012 #if NOSWIR
1013  retrievalradius1621 = fillvalue_real
1014  retrievalopticalthickness1621 = fillvalue_real
1015 #else
1016  !
1017  ! For 1.6, 2.1 um, continue and find the re, tau using the absorbing
1018  ! band reflectance for ice clouds
1019  !
1020  if (do_retrievals_ice(4)) then
1021  retrievalradius1621 = fillvalue_real
1022  retrievalopticalthickness1621 = fillvalue_real
1023 
1024  if( (process%ocean_surface .or. process%snowice_surface) .and. &
1025  measurements(band_0163) > 0. ) then!{
1026 
1027  nonabsorbingband_index = band_0163
1028  absorbingband_index = band_0213
1029  idx1621(1) = idx16
1030  idx1621(2) = idx21
1031 
1032  ! the Nakajima-King space for VIIRS (and AHI) 1.6-2.2 um ice
1033  ! retrieval appears to be reversed
1034  ! where the 2.2um channel contains the tau information and the
1035  ! 1.6um channel contains the
1036  ! re information. So we play along and reverse the channels.
1037 #ifdef VIIRS_INST
1038  if ( .not. modis_mode ) then
1039  nonabsorbingband_index = band_0213
1040  absorbingband_index = band_0163
1041  idx1621(1) = idx21
1042  idx1621(2) = idx16
1043  endif
1044 #endif
1045 #if AHI_INST | AMS_INST
1046  nonabsorbingband_index = band_0213
1047  absorbingband_index = band_0163
1048  idx1621(1) = idx21
1049  idx1621(2) = idx16
1050 #endif
1051 #if AVIRIS_INST
1052  if (set_bands(absorbingband_index) == 14) then
1053  nonabsorbingband_index = band_0213
1054  absorbingband_index = band_0163
1055  idx1621(1) = idx21
1056  idx1621(2) = idx16
1057  endif
1058 #endif
1059  call vis_nonabsorbing_science(measurements(nonabsorbingband_index), &
1060  idx1621(1), albedo_real4(nonabsorbingband_index), library_taus, &
1063  sensor_zenith_angle(xpoint,ypoint), &
1064  solar_zenith_angle(xpoint,ypoint), &
1065  relative_azimuth_angle(xpoint,ypoint), &
1066  cloud_top_pressure(xpoint,ypoint), &
1067  local_process, optical_thickness_vector)
1068 
1069  call vis_absorbing_science(optical_thickness_vector, &
1070  measurements(absorbingband_index), idx1621(2), &
1071  albedo_real4(absorbingband_index), library_taus, ice_radii, &
1074  residual,maxradii,debug)
1075 
1076  if (maxradii > 2 ) then!{
1077  call solveretrieval(residual(1:maxradii), &
1078  optical_thickness_vector(1:maxradii), ice_radii(1:maxradii), &
1079  retrievalradius1621, retrievalopticalthickness1621, &
1080  debug, use_nearest, quality_in)
1081  if (use_nearest) then
1082  nearest_ice(4) = 1
1083  call solveretrieval_nearest(xpoint,ypoint, &
1084  measurements(nonabsorbingband_index), &
1085  measurements(absorbingband_index), &
1086  (/nonabsorbingband_index, absorbingband_index/), ice_radii, &
1087  retrievalopticalthickness1621, retrievalradius1621, &
1088  rss_ice(re1621), .false., ray_temp, quality_in)
1089 
1090  endif
1091  endif!}
1092  endif!}
1093  else
1094  retrievalradius1621 = fillvalue_real
1095  retrievalopticalthickness1621 = fillvalue_real
1096  endif
1097  !
1098  ! WDR get the table range
1099  prd_out_refl_loc_1621(xpoint,ypoint,tbl_lo_nabs_ice) = minval(refliba)
1100  prd_out_refl_loc_1621(xpoint,ypoint,tbl_hi_nabs_ice) = maxval(refliba)
1101  prd_out_refl_loc_1621(xpoint,ypoint,tbl_lo_abs_ice) = minval(reflibb)
1102  prd_out_refl_loc_1621(xpoint,ypoint,tbl_hi_abs_ice) = maxval(reflibb)
1103  !
1104 #endif
1105  !
1106  ! end of ice cloud processing
1107  !
1108  deallocate(rad37lib)
1109  deallocate(optical_thickness_vector, residual, refliba, reflibb)
1110 
1111  else!}{
1112  ! there is no cloud
1113  retrievalopticalthickness = fillvalue_real
1114  retrievalopticalthickness1621 = fillvalue_real
1115  retrievalopticalthickness16 = fillvalue_real
1116  retrievalopticalthickness37 = fillvalue_real
1117  retrievalradius16 = fillvalue_real
1118  retrievalradius1621 = fillvalue_real
1119  retrievalradius37 = fillvalue_real
1120  if( ( c2_sensor_id == oci_id ) .OR. ( c2_sensor_id == ocis_id ) ) then
1121  retrievalopticalthickness22 = fillvalue_real
1122  retrievalradius22 = fillvalue_real
1123  endif
1124  cloud_layer_flag(xpoint, ypoint) = 0
1125  status = 1
1126  endif!}
1127 
1128  ! so the statistics are computed properly. The retrieval is only
1129  ! successful for statistical purposes if the main
1130  ! re retrieval is successful.
1131  if (retrievalradius21 == fillvalue_real) then
1132  status = 1
1133  endif
1134 
1135  !assign the retrievals to the relevant arrays
1136  optical_thickness_ice = retrievalopticalthickness
1137  optical_thickness_1621_ice = retrievalopticalthickness1621
1138  optical_thickness_16_ice = retrievalopticalthickness16
1139  optical_thickness_37_ice = retrievalopticalthickness37
1140  effective_radius_16_ice = retrievalradius16
1141  effective_radius_21_ice = retrievalradius21
1142  effective_radius_1621_ice = retrievalradius1621
1143  effective_radius_37_ice = retrievalradius37
1144  if( ( c2_sensor_id == oci_id ) .OR. ( c2_sensor_id == ocis_id ) ) then
1145  optical_thickness_22_ice = retrievalopticalthickness22
1146  effective_radius_22_ice = retrievalradius22
1147  endif
1148 
1149  ! reset the extra retrievals for now.
1150  ! if (nearest_ice(1) == 1) effective_radius_16_ice(xpoint, ypoint) = &
1151  ! fillvalue_real
1152  ! if (nearest_ice(2) == 1) then
1153  ! effective_radius_21_ice(xpoint, ypoint) = fillvalue_real
1154  ! optical_thickness_ice(xpoint, ypoint) = fillvalue_real
1155  ! endif
1156  ! if (nearest_ice(3) == 1) effective_radius_37_ice(xpoint, ypoint) = &
1157  ! fillvalue_real
1158  ! if (nearest_ice(4) == 1) then
1159  ! effective_radius_1621_ice(xpoint, ypoint) = fillvalue_real
1160  ! optical_thickness_1621_ice(xpoint, ypoint) = fillvalue_real
1161  ! endif
1162 
1163 end subroutine corescience
1164 
1165 end module corescience_module
Definition: ch_xfr.f90:1
integer, parameter re21
real optical_thickness_16_liquid
Definition: core_arrays.f90:24
integer, parameter band_0124
integer refl_abs
Definition: ch_xfr.f90:66
subroutine vis_nonabsorbing_science(reflectance_nonabsorbing, nonabsorbing_index, nonabsorbing_albedo, library_taus, library_radii, sfr, fti1, fti0, rfi, theta, theta0, phi, cloudtop_pressure, process, optical_thickness_vector)
integer(integer_onebyte), dimension(:,:), allocatable cloud_layer_flag
Definition: core_arrays.f90:92
real, dimension(set_number_of_bands) thermal_correction_oneway
real, dimension(:,:), allocatable rad37lib
subroutine solveretrieval(residual, optical_thickness_vector, library_radii, effective_radius, optical_thickness, debug, use_nearest, quality_out)
integer ocis_id
Definition: ch_xfr.f90:51
real, dimension(20) sigma_r37_pw_liq
real optical_thickness_1621_ice
Definition: core_arrays.f90:26
real(single), dimension(:,:), allocatable optical_thickness_37_final
Definition: core_arrays.f90:41
subroutine nir_absorbing_science(platform_name, optical_thickness_vector, reflectance_absorbing, absorbing_index, absorbing_albedo, xpoint, ypoint, CTT, thermal_trans_1way, thermal_trans_2way, library_taus, library_radii, sfr, fti1, fti0, fri1, rfi, cl_emis, sf_emis, residual, maxradii, channel_number_37, emission_uncertainty_pw, emission_uncertainty_Tc, sigma_R37_PW, debug)
real, dimension(20) sigma_r37_pw_ice
real(single), dimension(:,:), allocatable cloud_top_pressure
real effective_radius_21_ice
Definition: core_arrays.f90:27
integer, parameter re37
real function, public linearinterpolation(X, Y, XX)
real(single), dimension(:,:,:), allocatable int_reflectance_ice
subroutine solve_retrieval_noswir(optical_thickness_vector, library_radii, effective_radius, optical_thickness)
real, dimension(:,:), allocatable irw_temperature
integer tbl_lo_nabs_ice
Definition: ch_xfr.f90:66
real optical_thickness_ice
Definition: core_arrays.f90:23
real effective_radius_21_liquid
Definition: core_arrays.f90:27
real, dimension(:,:,:), allocatable prd_out_refl_loc_1600
Definition: ch_xfr.f90:63
real optical_thickness_22_liquid
Definition: core_arrays.f90:29
real(single), dimension(:,:,:), allocatable int_fluxdownwater_solar
character *15 platform_name
integer, parameter band_0370
real optical_thickness_16_ice
Definition: core_arrays.f90:24
real, dimension(set_number_of_bands) thermal_correction_twoway
real effective_radius_1621_ice
Definition: core_arrays.f90:31
integer, parameter band_0226
integer, parameter band_0086
real(single), dimension(:,:,:), allocatable int_cloud_emissivity_ice
integer, parameter band_0213
real(single), dimension(:,:,:), allocatable int_cloud_emissivity_water
real effective_radius_22_ice
Definition: core_arrays.f90:28
real effective_radius_37_liquid
Definition: core_arrays.f90:32
real(single), dimension(:), allocatable water_radii
real(single), dimension(:), allocatable library_taus
real, dimension(:,:), allocatable refliba
real, dimension(20) emission_uncertainty_tc_liq
real optical_thickness_37_liquid
Definition: core_arrays.f90:25
subroutine vis_absorbing_science(optical_thickness_vector, reflectance_absorbing, absorbing_index, absorbing_albedo, library_taus, library_radii, sfr, fti1, fti0, rfi, residual, maxradii, debug)
real optical_thickness_1621_liquid
Definition: core_arrays.f90:26
real, dimension(:,:,:), allocatable prd_out_refl_loc_2200
Definition: ch_xfr.f90:65
real(single), dimension(:,:,:), allocatable spherical_albedo_water
integer tbl_hi_nabs_ice
Definition: ch_xfr.f90:66
real, dimension(:,:), allocatable solar_zenith_angle
Definition: core_arrays.f90:6
integer refl_nabs
Definition: ch_xfr.f90:66
real(single), dimension(:,:,:), allocatable int_surface_emissivity_ice
real(single), dimension(:,:,:), allocatable int_fluxdownice_solar
real(single), dimension(:,:,:), allocatable int_surface_emissivity_water
real(single), dimension(:,:,:), allocatable spherical_albedo_ice
real optical_thickness_37_ice
Definition: core_arrays.f90:25
endif() set(LIBS $
Definition: CMakeLists.txt:6
real optical_thickness_liquid
Definition: core_arrays.f90:23
integer, parameter re22
real(single), dimension(:), allocatable ice_radii
real effective_radius_22_liquid
Definition: core_arrays.f90:28
real, dimension(:,:), allocatable reflibb
real, dimension(20) emission_uncertainty_pw_liq
integer, parameter band_0163
real effective_radius_37_ice
Definition: core_arrays.f90:32
integer c2_sensor_id
Definition: ch_xfr.f90:50
integer tbl_hi_abs_ice
Definition: ch_xfr.f90:66
real(single), dimension(:,:), allocatable cloud_top_temperature
real(single), dimension(:,:,:), allocatable surface_emissivity_land
real, dimension(:,:,:), allocatable prd_out_refl_loc_1621
Definition: ch_xfr.f90:64
integer, dimension(set_number_of_bands), parameter set_bands
real effective_radius_16_liquid
Definition: core_arrays.f90:30
real, dimension(:,:), allocatable sensor_zenith_angle
real(single), dimension(:,:), allocatable relative_azimuth_angle
real(single), dimension(:,:,:), allocatable int_fluxupice_sensor
real(single), dimension(:,:), allocatable surface_temperature
integer tbl_lo_abs_ice
Definition: ch_xfr.f90:66
subroutine calculate_new_tc(platform_name, Tc, Tg, galbedo, wlen, tau, re, lib_taus, lib_res, sph_albedo, down_flux_sensor, up_flux_sensor, cloud_emiss, surface_emiss, newTc, PRN)
real(single), dimension(:,:,:), allocatable int_reflectance_water
integer, parameter band_0065
integer number_waterradii
real(single), dimension(:,:), allocatable optical_thickness_1621_final
Definition: core_arrays.f90:42
subroutine solveretrieval_nearest(xx_pt, yy_pt, Ram, Rbm, twobands, radii, tau, re, lib_dist, phase_liquid, Ram_corr, quality_in, CH37_IDX, CTopT, CH37_NUM, platFormName)
subroutine, public corescience(xpoint, ypoint, process, measurements, Tc_liquid, Tc_ice, debug, na_band_used, nearest_liq, nearest_ice, RSS_liq, RSS_ice, alt_ray_liq, alt_ray_ice, do_retrievals_liq, do_retrievals_ice, status)
real effective_radius_16_ice
Definition: core_arrays.f90:30
type(qualityanalysis), dimension(:,:), allocatable processing_information
integer tbl_lo_nabs
Definition: ch_xfr.f90:66
real(single), dimension(:,:,:), allocatable int_fluxdownwater_sensor
integer number_iceradii
real, dimension(20) emission_uncertainty_pw_ice
integer oci_id
Definition: ch_xfr.f90:52
subroutine get_band_idx(idx16, idx21, idx37, channel_37, idx_11, idx_alb37, idx_alb16)
integer tbl_hi_nabs
Definition: ch_xfr.f90:66
integer tbl_lo_abs
Definition: ch_xfr.f90:66
real(single), dimension(:,:,:), allocatable int_fluxdownice_sensor
integer tbl_hi_abs
Definition: ch_xfr.f90:66
real, dimension(:,:,:), allocatable prd_out_refl_loc_2100
Definition: ch_xfr.f90:62
real, dimension(set_albedo_bands) albedo_real4
#define abs(a)
Definition: misc.h:90
real optical_thickness_22_ice
Definition: core_arrays.f90:29
integer, parameter re1621
real, dimension(20) emission_uncertainty_tc_ice
real(single), dimension(:,:,:), allocatable int_fluxupwater_sensor
integer, parameter re16
real effective_radius_1621_liquid
Definition: core_arrays.f90:31
integer number_taus