Merge remote-tracking branch 'origin/release-v4.6.1'
[WRF.git] / var / da / da_radiance / da_read_pseudo_rad.inc
blob1308822ec27c5fd6f01f384cc7f0e8fa1c64758c
1 subroutine da_read_pseudo_rad (iv)
3    !---------------------------------------------------------------------------
4    !  Purpose: assign pseudo single radiance obs
5    !            one sensor, one pixel, one channel
6    !---------------------------------------------------------------------------
8    use da_control
10    implicit none
12    type (iv_type),  intent (inout) :: iv
14    ! Instrument triplet, follow the convension of RTTOV 
15    integer   :: platform_id, satellite_id, sensor_id
17    real, allocatable    :: tb_inv(:)  !  bright temperatures
18    type (datalink_type) :: p
20    logical :: outside, outside_all
21    integer :: i,k,n,nchan,inst, alloc_stat
22    type(info_type)       ::  info
23    type(model_loc_type)  ::  loc
25    call da_trace_entry("da_read_pseudo_rad")
27    ! Initialize variables
29    platform_id  = pseudo_rad_platid
30    satellite_id = pseudo_rad_satid
31    sensor_id    = pseudo_rad_senid
32    if (sensor_id == 0) then       ! hirs
33       nchan=19
34    else if (sensor_id == 1) then  ! msu
35       nchan=nchan_msu
36    else if (sensor_id == 3) then  ! amsua
37       nchan=nchan_amsua
38    else if (sensor_id == 4)  then ! amsub
39       nchan=nchan_amsub
40    else if (sensor_id == 15)  then ! mhs
41       nchan=nchan_mhs
42    else if (sensor_id == 10)  then ! ssmis
43       nchan=nchan_ssmis
44    else if (sensor_id == 11)  then ! airs
45       nchan=nchan_airs
46    else if (sensor_id == 16)  then ! iasi
47       nchan= 616
48    else if (sensor_id == 21)  then ! seviri
49       nchan = 8 
50    else if (sensor_id == 63)  then ! amsr2
51       nchan = 14
52    else if (sensor_id == 22)  then ! goes-imager
53       nchan = 4
54    else if (sensor_id == 23)  then ! mwhs2
55       nchan = 15
56    else
57       call da_error(__FILE__,__LINE__, &
58            (/"Error in setting up pseudo radiance, unknown sensor_id."/))
59    end if
61    inst = 1    ! single instrument
63    allocate (tb_inv(nchan))
65    info%lat  =  pseudo_rad_lat ! in degree
66    info%lon  =  pseudo_rad_lon
67    info%date_char = "0000-00-00_00:00:00"
68    call da_llxy (info, loc, outside, outside_all)
70    if (outside) then
71       iv%info(radiance)%nlocal          = 0
72       iv%info(radiance)%ntotal          = 0
73       iv%info(radiance)%ptotal(1)       = 0
74       iv%instid(inst)%num_rad           = 0
75       iv%instid(inst)%info%nlocal       = 0
76    else  
77       iv%info(radiance)%nlocal          = 1
78       iv%info(radiance)%ntotal          = 1
79       iv%info(radiance)%ptotal(1)       = 1
80       iv%instid(inst)%num_rad           = 1
81       iv%instid(inst)%info%nlocal       = 1
82    end if
84    do k = 1, nchan
85       tb_inv(k) = missing_r
86    end do
87    tb_inv(pseudo_rad_ichan)=pseudo_rad_inv
89    !  5.0 allocate innovation radiance structure
90    !----------------------------------------------------------------  
91    
92    i = 1
94    if ( iv%instid(inst)%num_rad > 0 ) then
95       write(UNIT=stdout,FMT='(a,i3,2x,a,3x,i10)') &
96         'Allocating space for radiance innov structure', &
97          i, iv%instid(inst)%rttovid_string, iv%instid(inst)%num_rad
99       call da_allocate_rad_iv(i,nchan,iv)
101       !  6.0 assign sequential structure to innovation structure
102       !-------------------------------------------------------------
104       n = 1     ! single obs
106       allocate (p%tb_inv(1:nchan), stat=alloc_stat)
107       if ( alloc_stat /= 0 ) CALL da_error(__FILE__,__LINE__,(/"error allocating"/))
109       p%info             = info
110       p%loc              = loc
111       p%landsea_mask     = 1
112       p%scanpos          = 1
113       p%satzen           = 0.0
114       p%satazi           = 0.0
115       p%solzen           = 0.0
116       p%tb_inv(1:nchan)  = tb_inv(1:nchan)
117       p%sensor_index     = inst          
118       p%ifgat            = 1
119      
120       call da_initialize_rad_iv (inst, n, iv, p)
122       iv%instid(inst)%tb_qc(:,n)                   = qc_bad
123       iv%instid(inst)%tb_error(pseudo_rad_ichan,n) = pseudo_rad_err
124       iv%instid(inst)%tb_qc(pseudo_rad_ichan,n)    = qc_good
126       deallocate(p%tb_inv)
127    end if
129    deallocate(tb_inv)
131    call da_trace_exit("da_read_pseudo_rad")
132   
134 end subroutine da_read_pseudo_rad