NASA Logo
Ocean Color Science Software

ocssw V2022
rtm_support.f90
Go to the documentation of this file.
1 module rtm_support
2 !
3 ! Disclaimer: All this stuff is courtesy of UW-Madison. Author unknown.
4 ! There are no comments in the original code.
5 !
6 ! WDR Note that for the transition to l2gen, we changed the
7 ! model_info(...)%<val> to c2_model_info%<val> for l2gen anc values
8 !
9 !
10 ! have to use MASTER channel numbers in modis_planck, modis_bright, but MAS channels in PFAAST
11 
13  use names
14  use science_parameters, only : model_levels, d2r
16  use ch_xfr, only : c2_sensor_id, oci_id, ocis_id
17 #if !CT_CODE
19 #else
21 #endif
22 
23 ! use FASCODE_routines
24  use pfaast
25 #ifdef VIIRS_INST
26  use pfaast_modis
27 #endif
28 
29  implicit none
30 
31  private
32 
33 
34  integer :: start_level
35 
36 #ifdef CT_CODE
37 #if !MAS_INST & !EMAS_INST
38  integer, parameter :: nchan = set_number_of_bands
39 #else
40  integer, parameter :: nchan = set_number_of_bands-3
41 #endif
42 #else
43  integer, parameter :: nchan = 2
44 #endif
45 
46  real :: rtm_trans_2way(nchan, model_levels)
48  real :: rtm_rad_atm_clr(nchan, model_levels)
49  real :: rtm_cloud_prof(nchan, model_levels)
50 
51  real :: b_profile(nchan, model_levels)
52 
53 #ifndef CT_CODE
54 
56 
61 
66 
67 
68 #endif
69 
70 
71  integer :: last_i, last_j
72  real :: last_2way_angle
73  real :: last_zenith
74 
75 
76 
78  public :: init_rtm_vars
79 #ifndef CT_CODE
83 #endif
84 
85 contains
86 
87 
88  subroutine init_rtm_vars()
89 
90  last_zenith = -999.
91  last_2way_angle = -999.
92 
93  last_i = -1
94  last_j = -1
95 
96  start_level = 1
97 
98 #if MAS_INST || EMAS_INST
99  start_level = 36
100 #endif
101 
102  end subroutine init_rtm_vars
103 
104  subroutine get_rtm_parameters (platform, surface_emissivity, view_zenith, sun_zenith, i, j, x, y)
105 
106 ! WDR - this looks to be entirely for 3.7 um and up. Strictly speaking,
107 ! OCI should not need this, so try removing call to MODIS_fascode
108 ! True - OCI is unchanged, but Aqua does change - we'll have to do sensor
109 ! specific calling
111 
112  real, intent(in) :: view_zenith, surface_emissivity(:), sun_zenith
113  character(*), intent(in) :: platform
114  integer, intent(in) :: i, j, x, y
115 
116  real :: ozone(model_levels)
117  integer :: platform_id
118  integer :: idet, iok, sfc_level
119  real :: rtm_rad_clr, rtm_bt_clr
120  integer :: channels(nchan)
121  integer :: k, jj
122  integer :: emis_idx, path_len, channel_use
123  real :: angle_two_way, miu, miu0
124  logical :: do_2way, newatm, newang, new_2way
125 
126  real :: temp_mixr(model_levels)
127 
128 ! channels(1) = bands(band_0370)
129 #ifndef CT_CODE
130  channels(1) = set_bands(band_1100)
131  channels(2) = set_bands(band_0370)
132 #ifdef VIIRS_INST
133  if (modis_mode) then
134  channels(1) = set_bands_modis(band_1100)
135  channels(2) = set_bands_modis(band_0370)
136  endif
137 #endif
138 
139 
140  ozone = 0.
141 
142 #else
143  ! WDR, if I understand right, for modis, this is not used
144  ozone = c2_model_info%o3_profile
145 
146 #if !MAS_INST & !EMAS_INST
147  channels = set_bands
148 #else
149  channels = set_bands(1:set_number_of_bands-3)
150 #endif
151 
152 #endif
153 
154  idet = 0
155 
156 
157  ! re-execute the transmittance calculations only if it makes sense to do so.
158  ! There is more than one MODIS pixel that has the same vza and fits within same 1-degree box
159  ! so we will recycle the transmittance for points like that.
160 
161 
162  path_len = len(trim(act_lib_path))
163 
164  platform_id = -1
165  if (platform(1:5) == "Terra" .or. platform(1:5) == "terra") then
166  platform_id = 0
167  endif
168  if (platform(1:4) == "Aqua" .or. platform(1:4) == "aqua" ) then
169  platform_id = 1
170  endif
171 
172  miu0 = cos(d2r*sun_zenith)
173  miu = cos(d2r*view_zenith)
174 
175 
176  angle_two_way = acos( miu*miu0 / (miu+miu0)) / d2r
177 
178  newatm = .false.
179  newang = .false.
180  new_2way = .false.
181 
182  if (abs(view_zenith - last_zenith) > 0.005) newang = .true.
183  if (angle_two_way /= last_2way_angle) new_2way = .true.
184  if (last_i /= i .or. last_j /= j) newatm = .true.
185 
186  do_2way = .true.
187 
188 #ifdef CT_CODE
189  do_2way = .false.
190  new_2way = .false.
191 #endif
192  if (newatm) then
193  do k=1, nchan
194  do jj=start_level, model_levels
195  b_profile(k,jj) = modis_planck(platform, c2_model_info%temp_profile(jj), channels(k), 1)
196  end do
197  end do
198 
199  end if
200 
201  do k=1, nchan
202 
203 
204 #ifdef MODIS_INST
205 
206  if( ( c2_sensor_id .ne. oci_id ) .and. &
207  ( c2_sensor_id .ne. ocis_id ) ) then
208  call modis_fascode(trim(act_lib_path) // char(0), myyear, &
209  myday, c2_model_info%temp_profile, c2_model_info%mixr_profile, &
210  ozone, view_zenith, angle_two_way, platform_id, channels(k), &
211  idet, rtm_trans_atm_clr(k,:), rtm_trans_2way(k,:), &
212  newang, newatm, new_2way, do_2way, iok, x, y)
213  endif
214 
215 #endif
216 
217 #if MAS_INST || EMAS_INST
218  if (platform == "MASTER") then
219 #ifdef CT_CODE
220 ! channel_use = set_bands_mas(k)
221  channel_use = set_bands_master(k)
222 #else
223  if (k==1) then
224 ! channel_use = set_bands_mas(band_1100)
225  channel_use = set_bands_master(k)
226  else
227  channel_use = channels(k)
228  endif
229 #endif
230 
231  else
232  channel_use = channels(k)
233  endif
234 
235  call mas_fascode(trim(act_lib_path), path_len, myyear, myday, c2_model_info%temp_profile(36:model_levels), &
236  c2_model_info%mixr_profile(36:model_levels), ozone(36:model_levels), view_zenith, &
237  angle_two_way, mymission, channel_use, &
238  rtm_trans_atm_clr(k,36:model_levels), rtm_trans_2way(k, 36:model_levels), newang, newatm, new_2way, &
239  do_2way, iok)
240 
241 
242 #endif
243 
244 ! need to get (and keep) information about which MSG satellite we have: 8 or 9
245 ! this information must be passed to the code somehow by the production crew
246 #ifdef SEVIRI_INST
247 
248 
249  call seviri_fascode(trim(act_lib_path), path_len, myyear, myday, c2_model_info%temp_profile, &
250  c2_model_info%mixr_profile, ozone, view_zenith, angle_two_way, mymsg, channels(k), &
251  rtm_trans_atm_clr(k,:), rtm_trans_2way(k,:), newang, newatm, new_2way, do_2way, iok, x, y)
252 
253 #endif
254 
255 #ifdef AHI_INST
256 
257  call ahi_fascode(trim(act_lib_path), path_len, myyear, myday, c2_model_info%temp_profile, &
258  c2_model_info%mixr_profile, ozone, view_zenith, angle_two_way, channels(k), &
259  rtm_trans_atm_clr(k,:), rtm_trans_2way(k,:), newang, newatm, new_2way, do_2way, iok, x, y)
260 
261 #endif
262 
263 
264 
265 #ifdef VIIRS_INST
266  if (modis_mode) then
267 
268  call modis_fascode(trim(act_lib_path) // char(0), myyear, myday, c2_model_info%temp_profile, &
269  c2_model_info%mixr_profile, ozone, view_zenith, angle_two_way, &
270  platform_id, channels(k), idet, rtm_trans_atm_clr(k,:), rtm_trans_2way(k,:), &
271  newang, newatm, new_2way, do_2way, iok, x, y)
272 
273  else
274 
275 
276  call viirs_fascode(trim(act_lib_path), myyear, myday, c2_model_info%temp_profile, &
277  c2_model_info%mixr_profile, ozone, view_zenith, angle_two_way, channels(k), &
278  rtm_trans_atm_clr(k,:), rtm_trans_2way(k,:), newang, newatm, new_2way, do_2way, iok, x, y)
279  endif
280 #endif
281 
282 
283  end do
284 
285 #ifndef CT_CODE
286 
287 #ifdef MODIS_INST
288 ! we only need all this stuff for 3.7um channel, ever, so only run the 3.7um channel, don't waste time.
289  newatm = .true.
290  temp_mixr = c2_model_info%mixr_profile*0.8
291  do k=1, nchan
292 ! k = nchan
293  if( ( c2_sensor_id .ne. oci_id ) .and. ( c2_sensor_id .ne. ocis_id ) ) then
294  call modis_fascode(trim(act_lib_path) // char(0), myyear, myday, &
295  c2_model_info%temp_profile, temp_mixr, ozone, view_zenith, &
296  angle_two_way, platform_id, channels(k), idet, &
298  newang, newatm, new_2way, do_2way, iok, x, y)
299  endif
300  end do
301  newatm = .true.
302  temp_mixr = c2_model_info%mixr_profile*1.2
303  do k=1, nchan
304 ! k = nchan
305  if( ( c2_sensor_id .ne. oci_id ) .and. &
306  ( c2_sensor_id .ne. ocis_id ) ) then
307  call modis_fascode(trim(act_lib_path) // char(0), myyear, myday, &
308  c2_model_info%temp_profile, temp_mixr, ozone, view_zenith, &
309  angle_two_way, platform_id, channels(k), idet, &
311  newang, newatm, new_2way, do_2way, iok, x, y)
312  endif
313  end do
314 
316  abs(rtm_trans_2way_high(2,:) - rtm_trans_2way(2,:)) ) / 2.
317 
318 #endif
319 
320 #if MAS_INST || EMAS_INST
321 ! we only need all this stuff for 3.7um channel, ever, so only run the 3.7um channel, don't waste time.
322 ! the 3.7um channel is the same for MAS and MASTER
323  newatm = .true.
324  temp_mixr = c2_model_info%mixr_profile*0.8
325  do k=1, nchan
326 ! k = nchan
327 
328  call mas_fascode(trim(act_lib_path), path_len, myyear, myday, c2_model_info%temp_profile(36:model_levels), &
329  temp_mixr(36:model_levels), ozone(36:model_levels), view_zenith, &
330  angle_two_way, mymission, channel_use, &
331  rtm_trans_atm_clr_low(k,36:model_levels), rtm_trans_2way_low(k, 36:model_levels), newang, newatm, new_2way, &
332  do_2way, iok)
333  end do
334 
335  newatm = .true.
336  temp_mixr = c2_model_info%mixr_profile*1.2
337  do k=1, nchan
338 ! k = nchan
339  call mas_fascode(trim(act_lib_path), path_len, myyear, myday, c2_model_info%temp_profile(36:model_levels), &
340  temp_mixr(36:model_levels), ozone(36:model_levels), view_zenith, &
341  angle_two_way, mymission, channel_use, &
342  rtm_trans_atm_clr_high(k,36:model_levels), rtm_trans_2way_high(k, 36:model_levels), newang, newatm, new_2way, &
343  do_2way, iok)
344  end do
345 
347  abs(rtm_trans_2way_high(2,:) - rtm_trans_2way(2,:)) ) / 2.
348 
349 
350 #endif
351 
352 
353 
354 #ifdef SEVIRI_INST
355 ! we only need all this stuff for 3.7um channel, ever, so only run the 3.7um channel, don't waste time.
356  newatm = .true.
357  temp_mixr = c2_model_info%mixr_profile*0.8
358  do k=1, nchan
359 ! k = nchan
360  call seviri_fascode(trim(act_lib_path), path_len, myyear, myday, c2_model_info%temp_profile, &
361  temp_mixr, ozone, view_zenith, angle_two_way, mymsg, channels(k), &
362  rtm_trans_atm_clr_low(k,:), rtm_trans_2way_low(k,:), newang, newatm, new_2way, do_2way, iok, x, y)
363  end do
364 
365  newatm = .true.
366  temp_mixr = c2_model_info%mixr_profile*1.2
367  do k=1, nchan
368 ! k = nchan
369  call seviri_fascode(trim(act_lib_path), path_len, myyear, myday, c2_model_info%temp_profile, &
370  temp_mixr, ozone, view_zenith, angle_two_way, mymsg, channels(k), &
371  rtm_trans_atm_clr_high(k,:), rtm_trans_2way_high(k,:), newang, newatm, new_2way, do_2way, iok, x, y)
372  end do
373 
375  abs(rtm_trans_2way_high(2,:) - rtm_trans_2way(2,:)) ) / 2.
376 
377 
378 #endif
379 
380 #ifdef AHI_INST
381 ! we only need all this stuff for 3.7um channel, ever, so only run the 3.7um channel, don't waste time.
382  newatm = .true.
383  temp_mixr = c2_model_info%mixr_profile*0.8
384  do k=1, nchan
385 ! k = nchan
386  call ahi_fascode(trim(act_lib_path), path_len, myyear, myday, c2_model_info%temp_profile, &
387  temp_mixr, ozone, view_zenith, angle_two_way, channels(k), &
388  rtm_trans_atm_clr_low(k,:), rtm_trans_2way_low(k,:), newang, newatm, new_2way, do_2way, iok, x, y)
389 
390  end do
391 
392  newatm = .true.
393  temp_mixr = c2_model_info%mixr_profile*1.2
394  do k=1, nchan
395 ! k = nchan
396 
397  call ahi_fascode(trim(act_lib_path), path_len, myyear, myday, c2_model_info%temp_profile, &
398  temp_mixr, ozone, view_zenith, angle_two_way, channels(k), &
399  rtm_trans_atm_clr_high(k,:), rtm_trans_2way_high(k,:), newang, newatm, new_2way, do_2way, iok, x, y)
400 
401  end do
402 
404  abs(rtm_trans_2way_high(2,:) - rtm_trans_2way(2,:)) ) / 2.
405 
406 
407 #endif
408 
409 
410 #ifdef VIIRS_INST
411 ! we only need all this stuff for 3.7um channel, ever, so only run the 3.7um channel, don't waste time.
412 
413  newatm = .true.
414  temp_mixr = c2_model_info%mixr_profile*0.8
415  do k=1, nchan
416 ! k = nchan
417 
418  if (modis_mode) then
419  call modis_fascode(trim(act_lib_path) // char(0), myyear, myday, c2_model_info%temp_profile, &
420  temp_mixr, ozone, view_zenith, angle_two_way, &
421  platform_id, channels(k), idet, rtm_trans_atm_clr_low(k,:), rtm_trans_2way_low(k,:), &
422  newang, newatm, new_2way, do_2way, iok, x, y)
423  else
424  call viirs_fascode(trim(act_lib_path), myyear, myday, c2_model_info%temp_profile, &
425  temp_mixr, ozone, view_zenith, angle_two_way, channels(k), &
426  rtm_trans_atm_clr_low(k,:), rtm_trans_2way_low(k,:), newang, newatm, new_2way, do_2way, iok, x, y)
427  endif
428  end do
429  newatm = .true.
430  temp_mixr = c2_model_info%mixr_profile*1.2
431  do k=1, nchan
432 ! k = nchan
433 
434  if (modis_mode) then
435  call modis_fascode(trim(act_lib_path) // char(0), myyear, myday, c2_model_info%temp_profile, &
436  temp_mixr, ozone, view_zenith, angle_two_way, &
437  platform_id, channels(k), idet, rtm_trans_atm_clr_high(k,:), rtm_trans_2way_high(k,:), &
438  newang, newatm, new_2way, do_2way, iok, x, y)
439  else
440  call viirs_fascode(trim(act_lib_path), myyear, myday, c2_model_info%temp_profile, &
441  temp_mixr, ozone, view_zenith, angle_two_way, channels(k), &
442  rtm_trans_atm_clr_high(k,:), rtm_trans_2way_high(k,:), newang, newatm, new_2way, do_2way, iok, x, y)
443  endif
444  end do
445 
447  abs(rtm_trans_2way_high(2,:) - rtm_trans_2way(2,:)) ) / 2.
448 #endif
449 
450 
451 
452 #endif
453 
454 
455  sfc_level = c2_model_info%surface_level
456 
457  last_zenith = view_zenith
458  last_2way_angle = angle_two_way
459  last_i = i
460  last_j = j
461 
462  do k=1, nchan
463 
464 #ifndef CT_CODE
465  emis_idx = nchan - k + 1
466 #else
467  emis_idx = k
468 #endif
469 
470  call clear_atm_rad(platform, b_profile(k,:), c2_model_info%Ts, sfc_level, &
471  surface_emissivity(emis_idx), &
473  channels(k), rtm_rad_clr, rtm_bt_clr)
474 
475 #ifndef CT_CODE
476 
477  call clear_atm_rad(platform, b_profile(k,:), c2_model_info%Ts, sfc_level, &
478  surface_emissivity(emis_idx), &
480  channels(k), rtm_rad_clr, rtm_bt_clr)
481 
482  call clear_atm_rad(platform, b_profile(k,:), c2_model_info%Ts, sfc_level, &
483  surface_emissivity(emis_idx), &
485  channels(k), rtm_rad_clr, rtm_bt_clr)
486 
487 
488 #endif
489 
490 
491  end do
492 
493  end subroutine get_rtm_parameters
494 
495 
496  subroutine clear_atm_rad(platform, B_prof, Tsfc, sfc_level, esfc, &
497  rt_trans_atm, rt_rad_atm_clr, rt_cloud_prof, channel, rt_rad_clr, rt_bt_clr)
498 
499  character(*), intent(in) :: platform
500  real, intent(in) :: B_prof(:), rt_trans_atm(:)
501  real, intent(in) :: Tsfc, esfc
502  integer, intent(in) :: sfc_level, channel
503  real, intent(inout) :: rt_rad_clr, rt_bt_clr, rt_rad_atm_clr(:), rt_cloud_prof(:)
504 
505  integer :: i
506  real :: B1, B2, dtrn
507 
508  rt_rad_atm_clr(1:start_level) = 0.
509 
510 ! B1 = modis_planck(platform, temp_profile(start_level), channel, 1)
511  b1 = b_prof(start_level)
512  rt_rad_atm_clr(start_level) = 0.
513  rt_cloud_prof(start_level) = b1*rt_trans_atm(start_level)
514 
515  do i=start_level+1, sfc_level
516 
517 ! B2 = modis_planck(platform, temp_profile(i), channel, 1)
518  b2 = b_prof(i)
519  dtrn = -(rt_trans_atm(i) - rt_trans_atm(i-1))
520  rt_rad_atm_clr(i) = rt_rad_atm_clr(i-1) + (b1 + b2) / 2.0 * dtrn
521  rt_cloud_prof(i) = rt_rad_atm_clr(i) + b2*rt_trans_atm(i)
522  b1 = b2
523 
524  end do
525 
526 ! call clear_toa_rad(platform, rt_rad_atm_clr(sfc_level), rt_trans_atm(sfc_level), Tsfc, &
527 ! esfc, channel, rt_rad_clr, rt_bt_clr, .false.)
528 
529 
530  end subroutine clear_atm_rad
531 
532  subroutine get_clear_toa_rad(platform, Tsfc, esfc, sfc_level, rad_clr, bt_clr, clear_rad_table, clear_trans_table, PRN)
533 
535 
536  character(*), intent(in) :: platform
537  real, intent(in) :: esfc(:)
538  integer, intent(in) :: sfc_level
539  real, intent(inout) :: rad_clr(:), bt_clr(:)
540  real, intent(in) :: clear_rad_table(:,:), clear_trans_table(:,:)
541  integer :: channels(nchan)
542  integer :: k
543  real, intent(in) :: tsfc
544  logical, intent(in) :: prn
545 
546 ! channels(1) = bands(band_0370)
547 #ifndef CT_CODE
548  channels(1) = set_bands(band_1100)
549  channels(2) = set_bands(band_0370)
550 #ifdef VIIRS_INST
551  if (modis_mode) then
552  channels(1) = set_bands_modis(band_1100)
553  channels(2) = set_bands_modis(band_0370)
554  endif
555 #endif
556 
557 
558 #else
559 
560 #if !MAS_INST & !EMAS_INST
561  channels = set_bands
562 #else
563  channels = set_bands(1:set_number_of_bands-3)
564 #endif
565 
566 #endif
567 
568 ! rtm_rad_atm_clr, rtm_trans_atm_clr
569  do k=1, nchan
570  call clear_toa_rad(platform, clear_rad_table(k, sfc_level), clear_trans_table(k, sfc_level), tsfc, &
571  esfc(k), channels(k), rad_clr(k), bt_clr(k), prn )
572  end do
573 
574 
575  end subroutine get_clear_toa_rad
576 
577 
578  subroutine clear_toa_rad(platform, rad_atm, tau_atm, tsfc, &
579  esfc, channel, rad_clr, bt_clr, PRN)
580 
581  character(*), intent(in) :: platform
582  real, intent(in) :: rad_atm, tau_atm, tsfc, esfc
583  real, intent(inout) :: rad_clr, bt_clr
584  integer, intent(in) :: channel
585  logical, intent(in) :: PRN
586 
587  real :: rad_temp
588 
589 
590  rad_temp = modis_planck(platform, tsfc, channel, 1)
591 
592  rad_clr = rad_atm + esfc*rad_temp*tau_atm
593 
594  bt_clr = modis_bright(platform, rad_clr, channel, 1)
595 
596 
597  end subroutine clear_toa_rad
598 
599 end module rtm_support
Definition: ch_xfr.f90:1
Definition: pfaast.f90:1
real, dimension(nchan, model_levels), public rtm_trans_atm_clr
Definition: rtm_support.f90:47
real, dimension(nchan, model_levels), public rtm_rad_atm_clr_low
Definition: rtm_support.f90:59
character(len=2) mymission
Definition: names.f90:11
real, dimension(nchan, model_levels), public rtm_trans_atm_clr_low
Definition: rtm_support.f90:58
integer ocis_id
Definition: ch_xfr.f90:51
subroutine clear_atm_rad(platform, B_prof, Tsfc, sfc_level, esfc, rt_trans_atm, rt_rad_atm_clr, rt_cloud_prof, channel, rt_rad_clr, rt_bt_clr)
real, dimension(nchan, model_levels), public rtm_cloud_prof_high
Definition: rtm_support.f90:65
character(len=500) act_lib_path
Definition: names.f90:44
character *15 platform_name
integer, parameter band_0370
subroutine, public init_rtm_vars()
Definition: rtm_support.f90:89
integer, parameter band_1100
real, dimension(nchan, model_levels), public rtm_trans_atm_clr_high
Definition: rtm_support.f90:63
real, dimension(nchan, model_levels), public rtm_rad_atm_clr
Definition: rtm_support.f90:48
subroutine clear_toa_rad(platform, rad_atm, tau_atm, tsfc, esfc, channel, rad_clr, bt_clr, PRN)
real, dimension(nchan, model_levels), public rtm_trans_2way_high
Definition: rtm_support.f90:62
real function, public modis_bright(platform_name, RAD, BAND, UNITS, cwn_array, tcs_array, tci_array)
real, dimension(nchan, model_levels), public rtm_cloud_prof_low
Definition: rtm_support.f90:60
real, dimension(nchan, model_levels), public rtm_trans_2way_low
Definition: rtm_support.f90:57
real, dimension(nchan, model_levels), public rtm_cloud_prof
Definition: rtm_support.f90:49
integer mymsg
Definition: names.f90:9
real function, public modis_planck(platform_name, TEMP, BAND, UNITS)
subroutine, public get_rtm_parameters(platform, surface_emissivity, view_zenith, sun_zenith, i, j, x, y)
integer c2_sensor_id
Definition: ch_xfr.f90:50
subroutine, public get_clear_toa_rad(platform, Tsfc, esfc, sfc_level, rad_clr, bt_clr, clear_rad_table, clear_trans_table, PRN)
integer myyear
Definition: names.f90:9
integer myday
Definition: names.f90:9
string & trim(string &s, const string &delimiters)
Definition: EnvsatUtil.cpp:29
integer, dimension(set_number_of_bands), parameter set_bands
type(ancillary_type) c2_model_info
integer, parameter model_levels
Definition: names.f90:1
integer, parameter set_number_of_bands
real, dimension(nchan, model_levels), public rtm_trans_2way
Definition: rtm_support.f90:46
integer oci_id
Definition: ch_xfr.f90:52
real, dimension(nchan, model_levels), public rtm_rad_atm_clr_high
Definition: rtm_support.f90:64
real, dimension(model_levels), public rtm_trans_2way_mean
Definition: rtm_support.f90:55
#define abs(a)
Definition: misc.h:90
subroutine, public modis_fascode(coeff_dir_path, year, jday, temp, wvmr, ozmr, theta, ang_2way, platform, kban, jdet, taut, taut_2way, newang, newatm, new_2way, do_2way, iok, xxx, yyy)
Definition: pfaast.f90:23