Merge remote-tracking branch 'origin/release-v4.6.1'
[WRF.git] / var / da / da_obs_io / da_read_omb_tmp.inc
blobff83dd642a5dd2065bde9d96bb3ba7dcd8fb78b9
1 subroutine da_read_omb_tmp(filename,unit_in,num,obs_type_in,nc,if_wind_sd)
3    !-------------------------------------------------------------------------
4    ! read diagnostics written to temporary file by WRFVAR
5    !
6    ! 07 MAR 2014 -- Variables of OMB/OMA for wind obs. in diagnostic files 
7    !                are optional, i.e. SPD/DIR or U/V        -- Feng Gao
8    !-------------------------------------------------------------------------
10    implicit none
12    integer      ,intent (in)    :: unit_in
13    integer      ,intent (inout) :: num
14    character*(*),intent (in)    :: obs_type_in, filename 
15    integer      ,intent (in)    :: nc
17    integer      :: num_obs, ios 
18    character*20 :: iv_type
19    logical      :: if_write, if_wind_sd
20    
21 !  data format from da_write_obs_chem_sfc.inc
22    character(len=120)  :: fmt_chem = '(i8,3x,a6,2f9.2,5(2f17.7,i8,2f17.7))'
23    character(len=120)  :: fmt_chem1 = '(i8,2x,a6,2f11.6,2f11.2,i8,2f11.2)'
24    character(len=120)  :: fmt_chem2 = '(i8,2x,a6,2f11.6,2(2f11.2,i8,2f11.2))'
25    character(len=120)  :: fmt_chem4 = '(i8,2x,a6,2f11.6,4(2f12.3,i8,2f12.3))'
26    character(len=120)  :: wrt_chem1 = '(a6,2f11.6,2f11.2,i8,2f11.2)'
27    character(len=120)  :: wrt_chem2 = '(a6,2f11.6,2(2f11.2,i8,2f11.2))'
28    character(len=120)  :: wrt_chem4 = '(a6,2f11.6,4(2f12.3,i8,2f12.3))'
30    character*6  :: stn_id
31    integer      :: n, k, kk, l, levels, dummy_i
32    real         :: lat, lon, press, height, dummy
33    real         :: tpw_obs, tpw_inv, tpw_err, tpw_inc
34    real         :: u_obs, u_inv, u_error, u_inc, & 
35                    v_obs, v_inv, v_error, v_inc, &
36                    t_obs, t_inv, t_error, t_inc, &
37                    p_obs, p_inv, p_error, p_inc, &
38                    q_obs, q_inv, q_error, q_inc, &
39                    spd_obs, spd_inv, spd_err, spd_inc,   &
40                    ref_obs, ref_inv, ref_error, ref_inc, &
41                    eph_obs, eph_inv, eph_error, eph_inc, &
42                    rain_obs, rain_inv, rain_error, rain_inc, zk, &
43                    w_obs, w_inv, w_error, w_inc, &         ! lightning
44                    div_obs, div_inv, div_error, div_inc    ! lightning
45    integer     :: u_qc, v_qc, t_qc, p_qc, q_qc, tpw_qc, spd_qc, ref_qc, rain_qc, w_qc, div_qc
46 #if (WRF_CHEM == 1)
47    real         :: chem_obs, chem_inv, chem_err, chem_inc, &
48                    chem_obs2, chem_inv2, chem_err2, chem_inc2, &
49                    chem_obs3, chem_inv3, chem_err3, chem_inc3, &
50                    chem_obs4, chem_inv4, chem_err4, chem_inc4, &
51                    chem_obs5, chem_inv5, chem_err5, chem_inc5, &
52                    chem_obs6, chem_inv6, chem_err6, chem_inc6
53    integer     :: chem_qc, chem_qc2, chem_qc3, chem_qc4, chem_qc5, chem_qc6
54 #endif
55    integer     :: eph_qc
56    integer     :: ifgat
58    if (trace_use_dull) call da_trace_entry("da_read_omb_tmp")
60    open(unit=unit_in,file=trim(filename),form='formatted',status='old',iostat=ios)
61    if (ios /= 0) then
62       call da_error(__FILE__,__LINE__, (/"Cannot open file"//filename/))
63    end if
65    reports: do
67       read(unit_in,'(a20,i8)', end = 999, err = 1000) iv_type,num_obs
68       if_write = .false.
69       if (index(iv_type,OBS_type_in(1:nc)) > 0) if_write = .true.
71       select case (trim(adjustl(iv_type)))
73       case ('synop', 'ships', 'buoy', 'metar', 'sonde_sfc', 'tamdar_sfc')
74          if (num_obs > 0) then
75             do n = 1, num_obs    
76                read(unit_in,'(2i8)')levels, ifgat
77                if (if_write) then
78                   write(omb_unit,'(2i8)')levels, ifgat
79                num = num + 1
80                end if
81                do k = 1, levels
82                   read(unit_in,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))', err= 1000)&
83                      kk,l, stn_id, &          ! Station
84                      lat, lon, press, &       ! Lat/lon, pressure
85                      u_obs, u_inv, u_qc, u_error, u_inc, & 
86                      v_obs, v_inv, v_qc, v_error, v_inc, &
87                      t_obs, t_inv, t_qc, t_error, t_inc, &
88                      p_obs, p_inv, p_qc, p_error, p_inc, &
89                      q_obs, q_inv, q_qc, q_error, q_inc
91                   if (.not. if_wind_sd .and. wind_stats_sd) & 
92                      call da_ffdduv_diagnose(u_obs, u_inv, u_inc, v_obs, v_inv, v_inc, u_qc, v_qc, convert_uv2fd)
93                   if (if_wind_sd .and. .not. wind_stats_sd) &
94                      call da_ffdduv_diagnose(u_obs, u_inv, u_inc, v_obs, v_inv, v_inc, u_qc, v_qc, convert_fd2uv)
96                   if (if_write) &
97                      write(omb_unit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))', err= 1000)&
98                         num, k, stn_id, &          ! Station
99                         lat, lon, press, &       ! Lat/lon, pressure
100                         u_obs, u_inv, u_qc, u_error, u_inc, & 
101                         v_obs, v_inv, v_qc, v_error, v_inc, &
102                         t_obs, t_inv, t_qc, t_error, t_inc, &
103                         p_obs, p_inv, p_qc, p_error, p_inc, &
104                         q_obs, q_inv, q_qc, q_error, q_inc
105                end do
106             end do
107          end if
109          if (if_write) exit reports
110          cycle reports
112       case ('pilot', 'profiler', 'geoamv', 'qscat', 'polaramv')
113          if (num_obs > 0) then
114             do n = 1, num_obs    
115                read(unit_in,'(2i8)')levels, ifgat
116                if (if_write) then
117                   write(omb_unit,'(2i8)')levels, ifgat
118                num = num + 1
119                end if
120                do k = 1, levels
121                   read(unit_in,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))', err= 1000)&
122                       kk, l, stn_id, &          ! Station
123                       lat, lon, press, &        ! Lat/lon, pressure
124                       u_obs, u_inv, u_qc, u_error, u_inc, & 
125                       v_obs, v_inv, v_qc, v_error, v_inc
127                   if (.not. if_wind_sd .and. wind_stats_sd) &
128                      call da_ffdduv_diagnose(u_obs, u_inv, u_inc, v_obs, v_inv, v_inc, u_qc, v_qc, convert_uv2fd)
129                   if (if_wind_sd .and. .not. wind_stats_sd) &
130                      call da_ffdduv_diagnose(u_obs, u_inv, u_inc, v_obs, v_inv, v_inc, u_qc, v_qc, convert_fd2uv)
132                   if (if_write) &
133                      write(omb_unit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))', err= 1000)&
134                         num, k, stn_id, &          ! Station
135                         lat, lon, press, &         ! Lat/lon, pressure
136                         u_obs, u_inv, u_qc, u_error, u_inc, & 
137                         v_obs, v_inv, v_qc, v_error, v_inc
139                end do 
140             end do
141          end if
142          if (if_write) exit reports
143          cycle reports
145       case ('gpspw' )
146          if (num_obs > 0) then
147             do n = 1, num_obs    
148                read(unit_in,'(2i8)')levels, ifgat
149                if (if_write) then
150                   write(omb_unit,'(2i8)')levels, ifgat
151                num = num + 1
152                end if
153                do k = 1, levels
154                   read(unit_in,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))', err= 1000)&
155                      kk,l, stn_id, &          ! Station
156                      lat, lon, dummy, &       ! Lat/lon, dummy    
157                      tpw_obs, tpw_inv, tpw_qc, tpw_err, tpw_inc
158                   if (if_write) &
159                      write(omb_unit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))', err= 1000)&
160                         num, k, stn_id,  &       ! Station
161                         lat, lon, dummy, &       ! Lat/lon, dummy    
162                         tpw_obs, tpw_inv, tpw_qc, tpw_err, tpw_inc
163                end do
164             end do
165          end if
166          if (if_write) exit reports
167          cycle reports
169       case ('sound', 'tamdar', 'airep')
170          if (num_obs > 0) then
171             do n = 1, num_obs    
172                read(unit_in,'(2i8)')levels, ifgat
173                if (if_write) then
174                    write(omb_unit,'(2i8)')levels, ifgat
175                    num = num + 1 
176                end if
177                do k = 1, levels
178                   read(unit_in,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))', err= 1000)&
179                      kk,l, stn_id, &          ! Station
180                      lat, lon, press, &       ! Lat/lon, dummy    
181                      u_obs, u_inv, u_qc, u_error, u_inc, & 
182                      v_obs, v_inv, v_qc, v_error, v_inc, &
183                      t_obs, t_inv, t_qc, t_error, t_inc, &
184                      q_obs, q_inv, q_qc, q_error, q_inc
186                   if (.not. if_wind_sd .and. wind_stats_sd) &
187                      call da_ffdduv_diagnose(u_obs, u_inv, u_inc, v_obs, v_inv, v_inc, u_qc, v_qc, convert_uv2fd)  
188                   if (if_wind_sd .and. .not. wind_stats_sd) &
189                      call da_ffdduv_diagnose(u_obs, u_inv, u_inc, v_obs, v_inv, v_inc, u_qc, v_qc, convert_fd2uv) 
191                   if (if_write) &
192                      write(omb_unit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))', err= 1000)&
193                         num, k, stn_id,  &       ! Station
194                         lat, lon, press, &       ! Lat/lon, dummy    
195                         u_obs, u_inv, u_qc, u_error, u_inc, & 
196                         v_obs, v_inv, v_qc, v_error, v_inc, &
197                         t_obs, t_inv, t_qc, t_error, t_inc, &
198                         q_obs, q_inv, q_qc, q_error, q_inc
199                end do 
200             end do
201          end if
202      if (if_write) exit reports
203      cycle reports
205       case ('ssmir' )
206          if (num_obs > 0) then
207             do n = 1, num_obs    
208                read(unit_in,'(i8)')levels
209                if (if_write) then
210                   write(omb_unit,'(i8)')levels
211                   num = num + 1 
212                end if
213                do k = 1, levels
214                   read(unit_in,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))', err= 1000)&
215                      kk,l, stn_id, &          ! Station
216                      lat, lon, dummy, &       ! Lat/lon, dummy    
217                      spd_obs, spd_inv, spd_qc, spd_err, spd_inc, &
218                      tpw_obs, tpw_inv, tpw_qc, tpw_err, tpw_inc
219                   if (if_write) &
220                      write(omb_unit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))', err= 1000)&
221                         num, k, stn_id,  &       ! Station
222                         lat, lon, dummy, &       ! Lat/lon, dummy    
223                         spd_obs, spd_inv, spd_qc, spd_err, spd_inc, &
224                         tpw_obs, tpw_inv, tpw_qc, tpw_err, tpw_inc
225                end do
226             end do
227          end if
228          if (if_write) exit reports
229          cycle reports
230    
231       case ('ssmit' )
232          if (num_obs > 0) then
233             do n = 1, num_obs    
234                read(unit_in,'(i8)')levels
235                if (if_write) then
236                   write(omb_unit,'(i8)')levels
237                   num = num + 1 
238                end if
239                do k = 1, levels
240                   read(unit_in,'(2i8,a5,2f9.2,f17.7,7(2f17.7,i8,2f17.7))', err= 1000)&
241                      kk,l, stn_id, &          ! Station
242                      lat, lon, dummy, &       ! Lat/lon, dummy    
243                      dummy, dummy, dummy_i, dummy, dummy, &    
244                      dummy, dummy, dummy_i, dummy, dummy, &    
245                      dummy, dummy, dummy_i, dummy, dummy, &    
246                      dummy, dummy, dummy_i, dummy, dummy, &    
247                      dummy, dummy, dummy_i, dummy, dummy, &    
248                      dummy, dummy, dummy_i, dummy, dummy, &    
249                      dummy, dummy, dummy_i, dummy, dummy
250                   if (if_write) &
251                      write(omb_unit,'(2i8,a5,2f9.2,f17.7,7(2f17.7,i8,2f17.7))', err= 1000)&
252                         num,k,stn_id, &          ! Station
253                         lat, lon, dummy, &       ! Lat/lon, dummy    
254                         dummy, dummy, dummy_i, dummy, dummy, &    
255                         dummy, dummy, dummy_i, dummy, dummy, &    
256                         dummy, dummy, dummy_i, dummy, dummy, &    
257                         dummy, dummy, dummy_i, dummy, dummy, &    
258                         dummy, dummy, dummy_i, dummy, dummy, &    
259                         dummy, dummy, dummy_i, dummy, dummy, &    
260                         dummy, dummy, dummy_i, dummy, dummy
261                end do
262             end do
263          end if
264          if (if_write) exit reports
265          cycle reports
267       case ('satem' )
268          if (num_obs > 0) then
269             do n = 1, num_obs    
270                read(unit_in,'(i8)') levels
271                if (if_write) then
272                   write(omb_unit,'(i8)')levels
273                   num = num + 1 
274                end if
275                do k = 1, levels
276                   read(unit_in,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))', err= 1000)&
277                      kk,l, stn_id, &          ! Station
278                      lat, lon, press, &       ! Lat/lon, dummy    
279                      tpw_obs, tpw_inv, tpw_qc, tpw_err, tpw_inc
280                   if (if_write) &
281                      write(omb_unit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))', err= 1000)&
282                         num,k,stn_id, &          ! Station
283                         lat, lon, press, &       ! Lat/lon, dummy    
284                         tpw_obs, tpw_inv, tpw_qc, tpw_err, tpw_inc
285                end do  
286             end do
287          end if
288          if (if_write) exit reports
289          cycle reports
291       case ('ssmt1' , 'ssmt2' )
292          if (num_obs > 0) then
293             do n = 1, num_obs    
294                read(unit_in,'(i8)') levels
295                if (if_write) then
296                   write(omb_unit,'(i8)')levels
297                   num = num + 1 
298                end if
299                do k = 1, levels
300                   read(unit_in,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))', err= 1000)&
301                      kk,l, stn_id, &          ! Station
302                      lat, lon, dummy, &       ! Lat/lon, dummy    
303                      dummy,dummy, dummy_i, dummy, dummy
304                   if (if_write) &
305                      write(omb_unit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))', err= 1000)&
306                         num,k,stn_id, &          ! Station
307                         lat, lon, dummy, &       ! Lat/lon, dummy    
308                         dummy,dummy, dummy_i, dummy, dummy
309                end do 
310             end do
311          end if
312          if (if_write) exit reports
313          cycle reports
315       case ('bogus' )          
316          ! TC Bogus data is written in two records
317          ! 1st record holds info about surface level
318          ! 2nd is for upper air
320          if (num_obs > 0) then
321             do n = 1, num_obs    
322                read(unit_in,'(i8)') levels
323                if (if_write) then
324                   write(omb_unit,'(i8)')levels
325                   num = num + 1 
326                end if
327                do k = 1, levels
328                   read(unit_in,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))', err= 1000)&
329                       kk,l, stn_id, &          ! Station
330                       lat, lon, press, &       ! Lat/lon, dummy    
331                       u_obs, u_inv, u_qc, u_error, u_inc, & 
332                       v_obs, v_inv, v_qc, v_error, v_inc
333                   if (if_write) &
334                      write(omb_unit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))', err= 1000)&
335                          num,l,stn_id, &          ! Station
336                          lat, lon, press, &       ! Lat/lon, dummy    
337                          u_obs, u_inv, u_qc, u_error, u_inc, & 
338                          v_obs, v_inv, v_qc, v_error, v_inc
339                end do
340                read(unit_in,'(i8)') levels
341                if (if_write) then
342                   write(omb_unit,'(i8)')levels
343                end if
344                do k = 1, levels
345                   read(unit_in,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))', err= 1000)&
346                      kk,l, stn_id, &          ! Station
347                      lat, lon, press, &       ! Lat/lon, dummy    
348                      u_obs, u_inv, u_qc, u_error, u_inc, & 
349                      v_obs, v_inv, v_qc, v_error, v_inc
350                   if (if_write) &
351                      write(omb_unit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))', err= 1000)&
352                          num,l,stn_id, &          ! Station
353                          lat, lon, press, &       ! Lat/lon, dummy    
354                          u_obs, u_inv, u_qc, u_error, u_inc, & 
355                          v_obs, v_inv, v_qc, v_error, v_inc
356                end do
357             end do
358          end if
359          if (if_write) exit reports
360          cycle reports
362       case ('airsr' )          
363          if (num_obs > 0) then
364             do n = 1, num_obs    
365                read(unit_in,'(2i8)') levels, ifgat
366                if (if_write) then
367                   write(omb_unit,'(2i8)')levels, ifgat
368                   num = num + 1
369                end if
370                do k = 1, levels
371                   read(unit_in,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))', err= 1000)&
372                      kk,l, stn_id, &          ! Station
373                      lat, lon, press, &       ! Lat/lon, dummy    
374                      t_obs, t_inv, t_qc, t_error, t_inc, & 
375                      q_obs, q_inv, q_qc, q_error, q_inc
376                   if (if_write) &
377                      write(omb_unit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))', err= 1000)&
378                          num,k,stn_id, &          ! Station
379                          lat, lon, press, &       ! Lat/lon, dummy    
380                          t_obs, t_inv, t_qc, t_error, t_inc, & 
381                          q_obs, q_inv, q_qc, q_error, q_inc
382                end do
383             end do
384          end if
385          if (if_write) exit reports
386          cycle reports
388       case ('gpsref' )          
389          if (num_obs > 0) then
390             do n = 1, num_obs    
391                read(unit_in,'(2i8)') levels, ifgat
392                if (if_write) then
393                   write(omb_unit,'(2i8)')levels, ifgat
394                   num = num + 1
395                end if
396                do k = 1, levels
397                   read(unit_in,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))', err= 1000)&
398                      kk,l, stn_id, &          ! Station
399                      lat, lon, height, &       ! Lat/lon, height   
400                      ref_obs, ref_inv, ref_qc, ref_error, ref_inc 
401                   if (if_write) &
402                      write(omb_unit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))', err= 1000)&
403                         num,k,stn_id, &          ! Station
404                         lat, lon, height, &       ! Lat/lon, height   
405                         ref_obs, ref_inv, ref_qc, ref_error, ref_inc 
406                end do
407             end do
408          end if
409          if (if_write) exit reports
410          cycle reports
412       case ('gpseph' )
413          if (num_obs > 0) then
414             do n = 1, num_obs
415                read(unit_in,'(2i8)') levels, ifgat
416                if (if_write) write(omb_unit,'(2i8)')levels, ifgat
417                num = num + 1
418                do k = 1, levels
419                   read(unit_in,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))', err= 1000)&
420                      kk, l, stn_id, &          ! Station
421                      lat, lon, height, &       ! Lat/lon, height
422                      eph_obs, eph_inv, eph_qc, eph_error, eph_inc
423                   if (if_write) &
424                      write(omb_unit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))', err= 1000)&
425                         num, l, stn_id, &      ! Station
426                         lat, lon, height, &    ! Lat/lon, height
427                         eph_obs, eph_inv, eph_qc, eph_error, eph_inc
428                end do
429             end do
430          end if
431          if (if_write) exit reports
432          cycle reports
434       case ('rain' )          
435          if (num_obs > 0) then
436             do n = 1, num_obs    
437                read(unit_in,'(2i8)') levels, ifgat
438                if (if_write) then
439                   write(omb_unit,'(2i8)')levels, ifgat
440                   num = num + 1
441                end if
442                do k = 1, levels
443                   read(unit_in,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))', err= 1000)&
444                      kk,l, stn_id, &          ! Station
445                      lat, lon, height, &       ! Lat/lon, height   
446                      rain_obs, rain_inv, rain_qc, rain_error, rain_inc 
447                   if (if_write) &
448                      write(omb_unit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))', err= 1000)&
449                         num,k,stn_id, &          ! Station
450                         lat, lon, height, &       ! Lat/lon, height   
451                         rain_obs, rain_inv, rain_qc, rain_error, rain_inc 
452                end do
453             end do
454          end if
455          if (if_write) exit reports
456          cycle reports
458      case ('lightning' )
459          if (num_obs > 0) then
460             do n = 1, num_obs
461                read(unit_in,'(2i8)') levels, ifgat
462                if (if_write) then
463                   write(omb_unit,'(2i8)')levels, ifgat
464                   num = num + 1
465                end if
466                do k = 1, levels
467                   read(unit_in,'(2i8,a5,2f9.2,f17.7,3(2f17.7,i8,2f17.7))', err= 1000)&
468                      kk,l, stn_id, &          ! Station
469                      lat, lon, height, &       ! Lat/lon, height
470                      w_obs, w_inv, w_qc, w_error, w_inc, & ! vertical velocity
471                      div_obs, div_inv, div_qc, div_error, div_inc, & ! divergence
472                      q_obs, q_inv, q_qc, q_error, q_inc ! water vapor
473                   if (if_write) &
474                      write(omb_unit,'(2i8,a5,2f9.2,f17.7,3(2f17.7,i8,2f17.7))', err= 1000)&
475                         num,k,stn_id, &          ! Station
476                         lat, lon, height, &       ! Lat/lon, height
477                         w_obs, w_inv, w_qc, w_error, w_inc, & ! vertical velocity
478                         div_obs, div_inv, div_qc, div_error, div_inc, & ! divergence
479                         q_obs, q_inv, q_qc, q_error, q_inc ! water vapor
480                end do
481             end do
482          end if
483          if (if_write) exit reports
484          cycle reports
486 #if (WRF_CHEM == 1)
487       case ('chem' )
488          if (num_obs > 0) then
489             do n = 1, num_obs
490                read(unit_in,'(2i8)')levels, ifgat
491                if (if_write) then
492                   write(omb_unit,'(2i8)')levels, ifgat
493                num = num + 1
494                end if
495                do k = 1, levels
496                  if (chemicda_opt == 1 .or. chemicda_opt ==2 ) then
497                   read(unit_in,fmt=fmt_chem1, err= 1000)&
498                      kk, stn_id, &       ! Station
499                      lat, lon,   &       ! Lat/lon
500                      chem_obs, chem_inv, chem_qc, chem_err, chem_inc
501                   if (if_write) &
502                      write(omb_unit,fmt=wrt_chem1, err= 1000)&
503                         stn_id,  &       ! Station
504                         lat, lon,  &          ! Lat/lon
505                         chem_obs, chem_inv, chem_qc, chem_err, chem_inc
506                  else if (chemicda_opt == 3) then
507                   read(unit_in,fmt=fmt_chem2, err= 1000)&
508                      kk, stn_id, &          ! Station
509                      lat, lon,  &       ! Lat/lon
510                      chem_obs, chem_inv, chem_qc, chem_err, chem_inc, &
511                      chem_obs2, chem_inv2, chem_qc2, chem_err2, chem_inc2
512                   if (if_write) &
513                      write(omb_unit,fmt=wrt_chem2, err= 1000)&
514                         stn_id,  &       ! Station
515                         lat, lon, &       ! Lat/lon
516                         chem_obs, chem_inv, chem_qc, chem_err, chem_inc, &
517                         chem_obs2, chem_inv2, chem_qc2, chem_err2, chem_inc2
518                  else if (chemicda_opt == 4) then
519                   read(unit_in,fmt=fmt_chem4, err= 1000)&
520                      kk, stn_id, &          ! Station
521                      lat, lon,  &       ! Lat/lon
522                      chem_obs, chem_inv, chem_qc, chem_err, chem_inc, &
523                      chem_obs2, chem_inv2, chem_qc2, chem_err2, chem_inc2, &
524                      chem_obs3, chem_inv3, chem_qc3, chem_err3, chem_inc3, &
525                      chem_obs4, chem_inv4, chem_qc4, chem_err4, chem_inc4
526                   if (if_write) &
527                      write(omb_unit,fmt=wrt_chem4, err= 1000)&
528                         stn_id,  &       ! Station
529                         lat, lon,  &       ! Lat/lon, dummy
530                         chem_obs, chem_inv, chem_qc, chem_err, chem_inc, &
531                         chem_obs2, chem_inv2, chem_qc2, chem_err2, chem_inc2, &
532                         chem_obs3, chem_inv3, chem_qc3, chem_err3, chem_inc3, &
533                         chem_obs4, chem_inv4, chem_qc4, chem_err4, chem_inc4
534                  else if (chemicda_opt == 5) then
535                   read(unit_in,fmt=fmt_chem2, err= 1000)&
536                      kk, stn_id, &          ! Station
537                      lat, lon,  &       ! Lat/lon
538                      chem_obs, chem_inv, chem_qc, chem_err, chem_inc, &
539                      chem_obs2, chem_inv2, chem_qc2, chem_err2, chem_inc2
540                   if (if_write) &
541                      write(omb_unit,fmt=wrt_chem2, err= 1000)&
542                         stn_id,  &       ! Station
543                         lat, lon, &       ! Lat/lon
544                         chem_obs, chem_inv, chem_qc, chem_err, chem_inc, &
545                         chem_obs2, chem_inv2, chem_qc2, chem_err2, chem_inc2
547                  end if
548                end do
549             end do
551          end if !  if (num_obs > 0) then
553          if (if_write) exit reports
554          cycle reports
556       case ('gas' )
557          if (num_obs > 0) then
558            do n = 1, num_obs
559               read(unit_in,'(2i8)')levels, ifgat
560               if (if_write) then
561                  write(omb_unit,'(2i8)')levels, ifgat
562                  num = num + 1
563               end if
564               if (chemicda_opt == 4 .or. chemicda_opt == 5) then
565                   read(unit_in, fmt = fmt_chem4, err= 1000)&
566                        kk,  stn_id, &       ! Station
567                        lat, lon,    &       ! Lat/lon
568                        chem_obs, chem_inv, chem_qc, chem_err, chem_inc, &
569                        chem_obs2, chem_inv2, chem_qc2, chem_err2, chem_inc2, &
570                        chem_obs3, chem_inv3, chem_qc3, chem_err3, chem_inc3, &
571                        chem_obs4, chem_inv4, chem_qc4, chem_err4, chem_inc4
572                   if (if_write) &
573                   write(omb_unit, fmt = wrt_chem4, err= 1000)&
574                        stn_id,  lat, lon,     &       ! station, lat/lon
575                        chem_obs, chem_inv, chem_qc, chem_err, chem_inc, &
576                        chem_obs2, chem_inv2, chem_qc2, chem_err2, chem_inc2, &
577                        chem_obs3, chem_inv3, chem_qc3, chem_err3, chem_inc3, &
578                        chem_obs4, chem_inv4, chem_qc4, chem_err4, chem_inc4
579               end if
580            end do   !  do n = 1, num_obs
581          end if     !  if (num_obs > 0) then
583          if (if_write) exit reports
584          cycle reports
585 #endif
587       case default;
589          write(unit=message(1), fmt='(a,a20,a,i3)') &
590             'Got unknown obs_type string:', trim(iv_type),' on unit ',unit_in
591          call da_error(__FILE__,__LINE__,message(1:1))
592       end select
593    end do reports 
595 999 continue
596    close (unit_in)
598    if (trace_use_dull) call da_trace_exit("da_read_omb_tmp")
599    return
601 1000 continue
602    write(unit=message(1), fmt='(a,i3)') &
603       'read error on unit: ',unit_in
604    call da_warning(__FILE__,__LINE__,message(1:1))
606 end subroutine da_read_omb_tmp