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(*), 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       return
47    end if
48 
49    num_Radar = 0
50 
51    ! 2. read total radar
52    ! ===================
53 
54    ! 2.1 read first line
55    !     ---------------
56 
57    read (unit=iunit, fmt = '(A)', iostat = iost) char_total_Radar
58    if (iost /= 0) then
59       ! Does matter if present and unreadable
60       call da_error(__FILE__,__LINE__, &
61          (/"Cannot read radar file"/))
62    end if
63 
64    ! 2.3 total radar number
65 
66    read (unit=char_total_Radar (15:17),fmt='(I3)', iostat = iost) total_Radar
67 
68    ! 2.4 skip one lines
69 
70    read (unit=iunit, fmt = '(A)', iostat = iost)
71 
72    ! 3. read radar data
73 
74    do n = 1, total_Radar
75 
76       ! 3.1 skip one blank line
77 
78       read (unit=iunit, fmt = '(A)', iostat = iost)
79 
80       ! 3.2 read header
81 
82       read (unit=iunit, fmt = '(A)', iostat = iost) char_ned
83 
84       ! 3.3 read header information
85 
86       read (unit=char_ned (69:74), fmt='(I6)', iostat = iost) platform % stn % numObs
87 
88       ! 3.4 skip two lines
89 
90       read (unit=iunit, fmt = '(A)', iostat = iost)
91       read (unit=iunit, fmt = '(A)', iostat = iost)
92 
93       ! 3.5 loop over records
94 
95       reports: do j = 1, platform % stn % numObs
96 
97          ! 3.5.1 read station general info
98 
99          read (unit = iunit, iostat = iost, &
100                       fmt = '(A12,3X,A19,2X,2(F12.3,2X),F8.1,2X,I6)') &
101                       platform % info % platform,  &
102                       platform % info % date_char, &
103                       platform % info % lat,       &
104                       platform % info % lon,       &
105                       platform % info % elv,       &
106                       platform % info % levels
107 
108          read(unit=platform % info % platform (4:6), fmt='(I3)') fm
109 
110          ob % total_obs = ob % total_obs + 1
111 
112          !     3.5.2 read each level
113 
114          do i = 1, platform % info % levels
115             ! height
116             platform%each (i) = Radar_each_level_type(missing_r, missing, -1.0,&
117                field_type(missing_r, missing, missing_r), & ! rv
118                field_type(missing_r, missing, missing_r))   ! rf
119 
120             read (unit = iunit, fmt = '(3X, F12.1, 2(F12.3,I4,F12.3,2X))') &
121                              platform % each (i) % height,           &
122                              platform % each (i) % rv % inv,         &
123                              platform % each (i) % rv % qc,          &
124                              platform % each (i) % rv % error,       &
125                              platform % each (i) % rf % inv,         &
126                              platform % each (i) % rf % qc,          &
127                              platform % each (i) % rf % error
128          end do
129 
130          call da_ll_to_xy (platform%info, platform%loc, xp, outside)
131 
132          if (outside) then
133             cycle reports
134          end if
135 
136          nlevels = platform%info%levels
137 
138          if (nlevels > max_ob_levels) then
139              write(unit=message(1),fmt='(A,2I8)') &
140                 ' Radar=> nlevels > max_ob_levels:',nlevels, max_ob_levels
141              call da_warning(__FILE__,__LINE__,message(1:1))
142 
143              nlevels = max_ob_levels
144              platform%info%levels = nlevels
145          else if (nlevels < 1) then
146             cycle reports
147          end if
148 
149          ! Loop over duplicating obs for global
150          n_dup = 1
151          if (global .and. &
152             (platform%loc%i == xp%ids .or. platform%loc%i == xp%ide)) n_dup= 2
153    
154          do ndup = 1, n_dup
155             select case (fm)
156 
157             case (128)
158                num_Radar = num_Radar + 1
159 
160                if (num_Radar > max_Radar) then
161                   write(unit=message(1),fmt='(A,I6,A,I6)') &
162                      ' Radar #= ',ob % num_Radar, ' > max_Radar = ', max_Radar
163                   call da_error(__FILE__,__LINE__,message(1:1))
164                end if
165 
166             case default;
167                write(unit=stdout, fmt='(a)') 'Warining: unsaved obs found:'
168 
169                write(unit=stdout, fmt='(2a)') &
170                   'platform % info % platform=', platform % info % platform
171 
172                write(unit=stdout, fmt='(a, i3)') &
173                   'platform % info % levels=', platform % info % levels
174             end select
175          end do        !  loop over duplicate
176       end do reports
177    end do
178 
179    close (iunit)
180    call da_free_unit(iunit)
181  
182    if (print_detail_radar) write(unit=stdout, fmt='(/5x,a,i6,a)') &
183       ' Read:  ', num_Radar, ' Radar reports '
184 
185    ob%num_Radar = ob%num_Radar + num_Radar
186 
187    ob%ob_numb(ob%current_ob_time)%Radar = ob%num_Radar
188 
189    if (trace_use) call da_trace_exit("da_scan_radar")
190 
191 
192 end subroutine da_scan_radar
193 
194