updated top-level README and version_decl for V4.5 (#1847)
[WRF.git] / var / da / da_radiance / da_rttov_init.inc
blob6728b7adf3261455c0176a82f742ddabf8f4193a
1 subroutine da_rttov_init(iv,ob,nsensor,nchan)
3    !---------------------------------------------------------------------------
4    !  Purpose: interface to the initialization subroutine of RTTOV.
5    !
6    !  METHOD:  1. read RTTOV coefs files; 2. allocate radiance structure
7    !
8    !  HISTORY: 07/28/2005 - Creation                               Zhiquan Liu
9    !           03/22/2006   add error tuning factor read in        Zhiquan Liu
10    !           10/24/2007   limit this routine to RTTOV init       Tom Auligne
11    !---------------------------------------------------------------------------
13    implicit none 
15    type (iv_type), intent (inout) :: iv
16    type (y_type) , intent (inout) :: ob
17    integer ,       intent (in)    :: nsensor
18    integer ,       intent (in)    :: nchan(nsensor)
20    !  local arguments
21    !------------------- 
22    integer   :: i, j, n
24    !  input parameters of RTTOV_SETUP
25    !----------------------------------
26    integer :: err_unit        ! Logical error unit (<0 for default)
27    integer :: verbosity_level ! (<0 for default)
28    integer, allocatable  :: sensor(:,:) ! instrument id
29    integer, allocatable  :: coefs_channels (:,:)
30    character(len=256)    :: fmt_string
31    character(len=256)    :: file_path, file_suffix
32    character(len=256)    :: coef_prefix, scaer_prefix, sccld_prefix
33    character(len=256)    :: coef_filename, scaer_filename, sccld_filename
35    !  output parameters of RTTOV_SETUP
36    !-------------------------------------
37    integer :: errorstatus  ! return code
39    ! local variables
40    !----------------
41    integer               :: mxchn
42    integer(jpim)         :: id_sensor
44    if (trace_use) call da_trace_entry("da_rttov_init")
46    !--------------------------------------------------------------
47    !  1.0 setup RTTOV instrument triplets from namelist parameter
48    !--------------------------------------------------------------
50    mxchn           =  MAXVAL(nchan)
51    err_unit        =  stderr
52    verbosity_level =  rtminit_print
54    allocate (coefs(nsensor))
55    allocate (opts(nsensor))
56    allocate (opts_rt_ir(nsensor))      
57    allocate (sensor(3,nsensor))
58    allocate (coefs_channels(mxchn,nsensor))
59    if ( rttov_emis_atlas_ir > 0 .or. rttov_emis_atlas_mw > 0 ) then
60       allocate (atlas_type(nsensor))
61       allocate (atlas_id(nsensor))
62    end if
64    sensor (1,1:nsensor)  = rtminit_platform (1:nsensor) 
65    sensor (2,1:nsensor)  = rtminit_satid    (1:nsensor)
66    sensor (3,1:nsensor)  = rtminit_sensor   (1:nsensor)
68    coefs_channels(:,:) = 0
69    do n = 1, nsensor
70       coefs_channels(1:nchan(n),n) = iv%instid(n)%ichan(1:nchan(n))
71    end do
73    if (print_detail_rad) then
74       write(unit=message(1),fmt='(A,I5)') 'err_unit             = ', err_unit
75       write(unit=message(2),fmt='(A,I5)') 'verbosity_level      = ', verbosity_level
76       write(unit=message(3),fmt='(A,I5)') 'nsensor              = ', nsensor
77       write(unit=message(4),fmt='(A,10I5)') 'sensor (1,1:nsensor) = ', sensor (1,1:nsensor)
78       write(unit=message(5),fmt='(A,10I5)') 'sensor (2,1:nsensor) = ', sensor (2,1:nsensor)
79       write(unit=message(6),fmt='(A,10I5)') 'sensor (3,1:nsensor) = ', sensor (3,1:nsensor)
80       call da_message(message(1:6))
81    end if
83    !-----------------------------------------------------------
84    ! 2.0 read and initialize coefficients
85    !     rttov_read_coefs and rttov_init_coefs are called for each instrument
86    !-----------------------------------------------------------
88    file_path    = 'rttov_coeffs/'
89    file_suffix  = '.dat'
90    coef_prefix  = 'rtcoef_'
91    scaer_prefix = 'scaercoef_'  ! scattering coefficients - aerosols
92    sccld_prefix = 'sccldcoef_'  ! scattering coefficients - clouds
94    write(unit=message(1),fmt='(a)') 'Read in the RTTOV coef files for the following sensors'
95    call da_message(message(1:1))
97    do n = 1, nsensor
99       opts(n) % rt_all % switchrad  = .true.   ! brightness temperature radiance%bt is the input perturbation
100       opts(n) % rt_all % addrefrac  = .false.  ! refraction in path calc
101       opts(n) % interpolation % addinterp  = .false.  ! interpolation of input profile
102       opts(n) % rt_mw % clw_data = .false.
103       opts_rt_ir(n) % addsolar   = .false.  ! reflected solar
104       opts_rt_ir(n) % addclouds  = .false.  ! cloud effect
105       opts_rt_ir(n) % addaerosl  = .false.  ! aerosol effect
106       opts_rt_ir(n) % ozone_data = .false.
107       opts_rt_ir(n) % co2_data   = .false.
108       opts_rt_ir(n) % n2o_data   = .false.
109       opts_rt_ir(n) % ch4_data   = .false.
110       opts_rt_ir(n) % co_data    = .false.
112       ! construct the full path name to the coef file
113       coef_filename = trim(file_path)//trim(coef_prefix)//trim(iv%instid(n)%rttovid_string_coef)//trim(file_suffix)
114       scaer_filename = trim(file_path)//trim(scaer_prefix)//trim(iv%instid(n)%rttovid_string_coef)//trim(file_suffix)
115       sccld_filename = trim(file_path)//trim(sccld_prefix)//trim(iv%instid(n)%rttovid_string_coef)//trim(file_suffix)
117       write(unit=fmt_string, fmt='(a,i3,a)') '(2a,2x,a,', nchan(n), 'i5)'
118       write(unit=message(2),fmt=trim(fmt_string)) "   ", &
119       !write(unit=message(2),fmt='(2a,2x,a,(30i5))') "   ", &
120          trim(iv%instid(n)%rttovid_string), 'nchan = ', coefs_channels(1:nchan(n),n)
121       call da_message(message(2:2))
123       call rttov_read_coefs( &
124         & errorstatus,       &   ! out
125         & coefs(n),          &   ! out
126         & opts(n),           &   ! in
127         & channels     = coefs_channels(1:nchan(n),n),  &   ! list of channels to extract coefs
128         & file_coef    = trim(coef_filename),  &
129         & file_scaer   = trim(scaer_filename), &
130         & file_sccld   = trim(sccld_filename) )
131       if ( errorstatus /= errorstatus_success ) then
132          call da_error(__FILE__,__LINE__,(/"rttov_read_coefs fatal error"/))
133       end if
135       iv%instid(n)%nlevels = coefs(n)%coef%nlevels
137       if ( rttov_emis_atlas_ir > 0 .or. rttov_emis_atlas_mw > 0 ) then
138          id_sensor = coefs(n)%coef%id_sensor
139          atlas_type(n) = 0
140          if( id_sensor == sensor_id_ir .OR. id_sensor == sensor_id_hi ) then
141             atlas_type(n) = atlas_type_ir
142 !            atlas_id(n) = uwiremis_atlas_id    !(Previous WRFDA default)
143             atlas_id(n) = rttov_emis_atlas_ir !(namelist variable, can either be 1=uwiremis or 2=camel)
144          end if
145          if( id_sensor == sensor_id_mw .OR. id_sensor == sensor_id_po ) then
146             atlas_type(n) = atlas_type_mw
147             atlas_id(n) = rttov_emis_atlas_mw !(namelist variable, can either be 1=TELSEM2 or 2=CNRW)
148          end if
149       end if
150    end do
152    deallocate (sensor)
153    deallocate (coefs_channels)
155    if (trace_use) call da_trace_exit("da_rttov_init")
157 end subroutine da_rttov_init