da_scan_radar.inc

References to this file elsewhere.
1 subroutine da_scan_radar (ob, xp, filename)
2 
3    !---------------------------------------------------------------------------
4    ! Purpose: Scan the Radar observation file
5    !---------------------------------------------------------------------------
6 
7    implicit none
8 
9    type (xpose_type), intent(in)  :: xp       ! Domain decomposition vars.
10    type (ob_type), intent (inout) :: ob
11 
12    character(len=*), intent(in)   :: filename
13 
14    integer                       :: i, j, n, iost, nlevels, num_Radar, fm
15    integer                       :: total_Radar
16    integer                       :: iunit
17 
18    type (Radar_multi_level_type) :: platform
19 
20    character (LEN = 120)         :: char_total_Radar
21    character (LEN = 120)         :: char_ned
22 
23    logical                       :: outside
24    integer                        ::n_dup, ndup
25 
26    if (trace_use) call da_trace_entry("da_scan_radar")
27 
28    ! default value.
29    ob%ob_numb(ob%current_ob_time)%Radar = ob%num_Radar
30 
31    ! 1. open file
32    ! ============
33 
34    call da_get_unit(iunit)
35    open(unit   = iunit,     &
36         FILE   = trim(filename), &
37         FORM   = 'FORMATTED',  &
38         ACCESS = 'SEQUENTIAL', &
39         iostat =  iost,     &
40         STATUS = 'OLD')
41 
42    if (iost /= 0) then
43       ! Does not matter of radar file missing
44       call da_warning(__FILE__,__LINE__, &
45          (/"Cannot open radar file "//trim(filename)/))
46       call da_free_unit(iunit) 
47       if (trace_use) call da_trace_exit("da_scan_radar")
48       return
49    end if
50 
51    num_Radar = 0
52 
53    ! 2. read total radar
54    ! ===================
55 
56    ! 2.1 read first line
57    !     ---------------
58 
59    read (unit=iunit, fmt = '(A)', iostat = iost) char_total_Radar
60    if (iost /= 0) then
61       ! Does matter if present and unreadable
62       call da_error(__FILE__,__LINE__, &
63          (/"Cannot read radar file"/))
64    end if
65 
66    ! 2.3 total radar number
67 
68    read (unit=char_total_Radar (15:17),fmt='(I3)', iostat = iost) total_Radar
69 
70    ! 2.4 skip one lines
71 
72    read (unit=iunit, fmt = '(A)', iostat = iost)
73 
74    ! 3. read radar data
75 
76    do n = 1, total_Radar
77 
78       ! 3.1 skip one blank line
79 
80       read (unit=iunit, fmt = '(A)', iostat = iost)
81 
82       ! 3.2 read header
83 
84       read (unit=iunit, fmt = '(A)', iostat = iost) char_ned
85 
86       ! 3.3 read header information
87 
88       read (unit=char_ned (69:74), fmt='(I6)', iostat = iost) platform % stn % numObs
89 
90       ! 3.4 skip two lines
91 
92       read (unit=iunit, fmt = '(A)', iostat = iost)
93       read (unit=iunit, fmt = '(A)', iostat = iost)
94 
95       ! 3.5 loop over records
96 
97       reports: do j = 1, platform % stn % numObs
98 
99          ! 3.5.1 read station general info
100 
101          read (unit = iunit, iostat = iost, &
102                       fmt = '(A12,3X,A19,2X,2(F12.3,2X),F8.1,2X,I6)') &
103                       platform % info % platform,  &
104                       platform % info % date_char, &
105                       platform % info % lat,       &
106                       platform % info % lon,       &
107                       platform % info % elv,       &
108                       platform % info % levels
109 
110          read(unit=platform % info % platform (4:6), fmt='(I3)') fm
111 
112          ob % total_obs = ob % total_obs + 1
113 
114          !     3.5.2 read each level
115 
116          do i = 1, platform % info % levels
117             ! height
118             platform%each (i) = Radar_each_level_type(missing_r, missing, -1.0,&
119                field_type(missing_r, missing, missing_r), & ! rv
120                field_type(missing_r, missing, missing_r))   ! rf
121 
122             read (unit = iunit, fmt = '(3X, F12.1, 2(F12.3,I4,F12.3,2X))') &
123                              platform % each (i) % height,           &
124                              platform % each (i) % rv % inv,         &
125                              platform % each (i) % rv % qc,          &
126                              platform % each (i) % rv % error,       &
127                              platform % each (i) % rf % inv,         &
128                              platform % each (i) % rf % qc,          &
129                              platform % each (i) % rf % error
130          end do
131 
132          call da_ll_to_xy (platform%info, platform%loc, xp, outside)
133 
134          if (outside) then
135             cycle reports
136          end if
137 
138          nlevels = platform%info%levels
139 
140          if (nlevels > max_ob_levels) then
141              write(unit=message(1),fmt='(A,2I8)') &
142                 ' Radar=> nlevels > max_ob_levels:',nlevels, max_ob_levels
143              call da_warning(__FILE__,__LINE__,message(1:1))
144 
145              nlevels = max_ob_levels
146              platform%info%levels = nlevels
147          else if (nlevels < 1) then
148             cycle reports
149          end if
150 
151          ! Loop over duplicating obs for global
152          n_dup = 1
153          if (global .and. &
154             (platform%loc%i == xp%ids .or. platform%loc%i == xp%ide)) n_dup= 2
155    
156          do ndup = 1, n_dup
157             select case (fm)
158 
159             case (128)
160                num_Radar = num_Radar + 1
161 
162                if (num_Radar > max_Radar) then
163                   write(unit=message(1),fmt='(A,I6,A,I6)') &
164                      ' Radar #= ',ob % num_Radar, ' > max_Radar = ', max_Radar
165                   call da_error(__FILE__,__LINE__,message(1:1))
166                end if
167 
168             case default;
169                write(unit=stdout, fmt='(a)') 'Warining: unsaved obs found:'
170 
171                write(unit=stdout, fmt='(2a)') &
172                   'platform % info % platform=', platform % info % platform
173 
174                write(unit=stdout, fmt='(a, i3)') &
175                   'platform % info % levels=', platform % info % levels
176             end select
177          end do        !  loop over duplicate
178       end do reports
179    end do
180 
181    close (iunit)
182    call da_free_unit(iunit)
183  
184    if (print_detail_radar) write(unit=stdout, fmt='(/5x,a,i6,a)') &
185       ' Read:  ', num_Radar, ' Radar reports '
186 
187    ob%num_Radar = ob%num_Radar + num_Radar
188 
189    ob%ob_numb(ob%current_ob_time)%Radar = ob%num_Radar
190 
191    if (trace_use) call da_trace_exit("da_scan_radar")
192 
193 
194 end subroutine da_scan_radar
195 
196