updated top-level README and version_decl for V4.5 (#1847)
[WRF.git] / var / da / da_radiance / da_sort_rad.inc
bloba8c0da807abd05e95636c32b62f322240b6dcf12
1 subroutine da_sort_rad (iv)
3    !---------------------------------------------------------------------------
4    ! Purpose: sorting radiance innovation to FGAT time bins
5    !---------------------------------------------------------------------------
7    implicit none
9    type (iv_type), intent (inout) :: iv
11    integer                           :: i,j, n,t, error
12    integer, allocatable              :: ifgat(:),landsea_mask(:)
13    integer, allocatable              :: loc_i(:,:), loc_j(:,:),loc_k(:,:)
14    real, allocatable                 :: loc_dx(:,:),loc_dy(:,:),loc_dz(:,:)
15    real, allocatable                 :: loc_dxm(:,:),loc_dym(:,:),loc_dzm(:,:)
16    character (len = 40), allocatable :: name(:)       
17    character (len = 12), allocatable :: platform(:)   
18    character (len =  5), allocatable :: id(:)          
19    character (len = 19), allocatable :: date_char(:)   
20    integer, allocatable              :: levels(:)      
21    real, allocatable                 :: lat(:,:)         
22    real, allocatable                 :: lon(:,:)         
23    real, allocatable                 :: elv(:)        
24    real, allocatable                 :: pstar(:)      
25    real, allocatable                 :: scanpos(:), satzen(:), satazi(:), solzen(:)
26    real, allocatable                 :: tb_inv(:,:)
28    if (trace_use) call da_trace_entry("da_sort_rad")
30    iv%info(radiance)%plocal(:) = 0
31    if (num_fgat_time == 1) then
32       do i=1,rtminit_nsensor
33          iv%instid(i)%info%plocal(:) = 0
34          iv%instid(i)%info%plocal(num_fgat_time) = iv%instid(i)%num_rad
35          iv%info(radiance)%plocal(num_fgat_time) = iv%info(radiance)%plocal(num_fgat_time) + iv%instid(i)%num_rad
36       end do
37       if (trace_use) call da_trace_exit("da_sort_rad")
38       return
39    end if
41    do i=1,rtminit_nsensor
42       iv%instid(i)%info%plocal(:) = 0
43       if (iv%instid(i)%num_rad < 1) cycle
45       allocate (ifgat       (iv%instid(i)%num_rad)) 
46       allocate (landsea_mask(iv%instid(i)%num_rad))
47       allocate (loc_i       (kts:kte, iv%instid(i)%num_rad))
48       allocate (loc_j       (kts:kte, iv%instid(i)%num_rad))
49       allocate (loc_k       (iv%instid(i)%nlevels,iv%instid(i)%num_rad))
50       allocate (loc_dx      (kts:kte, iv%instid(i)%num_rad))
51       allocate (loc_dy      (kts:kte, iv%instid(i)%num_rad))
52       allocate (loc_dz      (iv%instid(i)%nlevels,iv%instid(i)%num_rad))
53       allocate (loc_dxm     (kts:kte, iv%instid(i)%num_rad))
54       allocate (loc_dym     (kts:kte, iv%instid(i)%num_rad))
55       allocate (loc_dzm     (iv%instid(i)%nlevels,iv%instid(i)%num_rad))
56       allocate (name        (iv%instid(i)%num_rad))
57       allocate (platform    (iv%instid(i)%num_rad))
58       allocate (id          (iv%instid(i)%num_rad))
59       allocate (date_char   (iv%instid(i)%num_rad))
60       allocate (levels      (iv%instid(i)%num_rad))
61       allocate (lat         (kts:kte,iv%instid(i)%num_rad))
62       allocate (lon         (kts:kte,iv%instid(i)%num_rad))
63       allocate (elv         (iv%instid(i)%num_rad))
64       allocate (pstar       (iv%instid(i)%num_rad))
65       allocate (scanpos     (iv%instid(i)%num_rad))
66       allocate (satzen      (iv%instid(i)%num_rad))
67       allocate (satazi      (iv%instid(i)%num_rad))
68       allocate (solzen      (iv%instid(i)%num_rad))
69       allocate (tb_inv      (iv%instid(i)%nchan,iv%instid(i)%num_rad))
71       j = 0
72       do t = 1,num_fgat_time
73          do n = 1,iv%instid(i)%num_rad
74             if (iv%instid(i)%ifgat(n) /= t) cycle
75             j = j + 1
76             ifgat(j)        = iv%instid(i)%ifgat(n)
77             name(j)         = iv%instid(i)%info%name(n)
78             platform(j)     = iv%instid(i)%info%platform(n)
79             id(j)           = iv%instid(i)%info%id(n)
80             date_char(j)    = iv%instid(i)%info%date_char(n)
81             levels(j)       = iv%instid(i)%info%levels(n)
82             elv(j)          = iv%instid(i)%info%elv(n)
83             pstar(j)        = iv%instid(i)%info%pstar(n)
85             lat    (:,j) = iv%instid(i)%info%lat(:,n)
86             lon    (:,j) = iv%instid(i)%info%lon(:,n)
87             loc_i  (:,j) = iv%instid(i)%info%i(:,n)
88             loc_j  (:,j) = iv%instid(i)%info%j(:,n)
89             loc_k  (:,j) = iv%instid(i)%info%k(:,n)
90             loc_dx (:,j) = iv%instid(i)%info%dx(:,n)
91             loc_dy (:,j) = iv%instid(i)%info%dy(:,n)
92             loc_dxm(:,j) = iv%instid(i)%info%dxm(:,n)
93             loc_dym(:,j) = iv%instid(i)%info%dym(:,n)
94             loc_dz (:,j) = iv%instid(i)%info%dz(:,n)
95             loc_dzm(:,j) = iv%instid(i)%info%dzm(:,n)
97             landsea_mask(j) = iv%instid(i)%landsea_mask(n)
98             scanpos(j)      = iv%instid(i)%scanpos(n)
99             satzen(j)       = iv%instid(i)%satzen(n)
100             satazi(j)       = iv%instid(i)%satazi(n)
101             solzen(j)       = iv%instid(i)%solzen(n)
103             tb_inv(1:iv%instid(i)%nchan,j) = iv%instid(i)%tb_inv(1:iv%instid(i)%nchan,n) 
104          end do
105          iv%instid(i)%info%plocal(t) = j
106          !write (0,*) __FILE__,__LINE__,"i,t,iv%instid(i)%info%plocal(t)",i,t,iv%instid(i)%info%plocal(t)
107       end do
109       iv%info(radiance)%plocal = iv%info(radiance)%plocal + iv%instid(i)%info%plocal
111       write(unit=stdout,fmt='(a,2x,a,2x,240i7)') &
112          ' FGAT: ',iv%instid(i)%rttovid_string, iv%instid(i)%info%plocal(1:num_fgat_time)
114       do n = 1,iv%instid(i)%num_rad
115          iv%instid(i)%ifgat(n)        = ifgat(n)
117          iv%instid(i)%info%name(n)          = name(n)
118          iv%instid(i)%info%platform(n)  = platform(n)
119          iv%instid(i)%info%id(n)            = id(n)
120          iv%instid(i)%info%date_char(n) = date_char(n)
121          iv%instid(i)%info%levels(n)    = levels(n)
122          iv%instid(i)%info%lat(:,n)     = lat(:,n)
123          iv%instid(i)%info%lon(:,n)     = lon(:,n)
124          iv%instid(i)%info%elv(n)           = elv(n)
125          iv%instid(i)%info%pstar(n)     = pstar(n)
127          iv%instid(i)%info%i(:,n)    = loc_i(:,n)
128          iv%instid(i)%info%j(:,n)    = loc_j(:,n)
129          iv%instid(i)%info%k(:,n)    = loc_k(:,n)
130          iv%instid(i)%info%dx(:,n)   = loc_dx(:,n)
131          iv%instid(i)%info%dy(:,n)   = loc_dy(:,n)
132          iv%instid(i)%info%dz(:,n)   = loc_dz(:,n)
133          iv%instid(i)%info%dxm(:,n)  = loc_dxm(:,n)
134          iv%instid(i)%info%dym(:,n)  = loc_dym(:,n)
135          iv%instid(i)%info%dzm(:,n)  = loc_dzm(:,n)
136          iv%instid(i)%landsea_mask(n) = landsea_mask(n)
137          iv%instid(i)%scanpos(n)      = scanpos(n)
138          iv%instid(i)%satzen(n)       = satzen(n)
139          iv%instid(i)%satazi(n)       = satazi(n)
140          iv%instid(i)%solzen(n)       = solzen(n)
142          iv%instid(i)%tb_inv(1:iv%instid(i)%nchan,n) = tb_inv(1:iv%instid(i)%nchan,n)
143       end do
145       deallocate (ifgat) 
146       deallocate (landsea_mask)
147       deallocate (name)
148       deallocate (platform)
149       deallocate (id)
150       deallocate (date_char)
151       deallocate (levels)
152       deallocate (lat)
153       deallocate (lon)
154       deallocate (elv)
155       deallocate (pstar)
156       deallocate (loc_i)
157       deallocate (loc_j)
158       deallocate (loc_k)
159       deallocate (loc_dx)
160       deallocate (loc_dy)
161       deallocate (loc_dz)
162       deallocate (loc_dxm)
163       deallocate (loc_dym)
164       deallocate (loc_dzm)
165       deallocate (scanpos)
166       deallocate (satzen)
167       deallocate (satazi)
168       deallocate (solzen)
169       deallocate (tb_inv)
170    end do
172    if (trace_use) call da_trace_exit("da_sort_rad")
174 end subroutine da_sort_rad