Merge remote-tracking branch 'origin/release-v4.6.1'
[WRF.git] / var / da / da_obs_io / da_write_obs.inc
blob7148b480d27ad08705a9cbdd1d71c76481992b4d
1 subroutine da_write_obs(it,ob, iv, re)
3    !-------------------------------------------------------------------------
4    ! Purpose: Writes out components of iv=O-B structure.
5    !-------------------------------------------------------------------------   
7    implicit none
9    integer,        intent(in)    :: it
10    type (y_type),  intent(in)    :: ob      ! Observation structure.
11    type (iv_type), intent(in)    :: iv      ! O-B structure.
12    type (y_type),  intent(inout) :: re      ! residual vector.
13       
14    integer                     :: n, k, num_obs, ios
15    integer                     :: ounit     ! Output unit           
16    character(len=filename_len) :: filename
17    integer :: itime, ifgat
19    if (trace_use) call da_trace_entry("da_write_obs")
21    !-------------------------------------------------------------------------
22    ! Fix output unit
23    !-------------------------------------------------------------------------
25    call da_get_unit(ounit)
27 #ifdef DM_PARALLEL
28     write(unit=filename, fmt='(a,i2.2,a,i4.4)') 'gts_omb_oma_',it,'.', myproc
29 #else
30     write(unit=filename, fmt='(a,i2.2,a)') 'gts_omb_oma_',it,'.0000'
31 #endif
33    open (unit=ounit,file=trim(filename),form='formatted',status='replace', &
34       iostat=ios)
35    if (ios /= 0) then
36       call da_error(__FILE__,__LINE__, &
37          (/"Cannot open conventional observation omb and oma file"//filename/))
38    end if
40    num_obs = 0
41    do n = 1, iv%info(synop)%nlocal
42       if (iv%info(synop)%proc_domain(1,n)) num_obs = num_obs + 1
43    end do
44    if (num_obs > 0) then
45       write(ounit,'(a20,i8)')'synop', num_obs  
46       num_obs = 0
47       do n = 1, iv%info(synop)%nlocal  
48          do itime = 1, num_fgat_time
49             if ( n >= iv%info(synop)%plocal(itime-1)+1 .and. &
50                  n <= iv%info(synop)%plocal(itime) ) then
51                ifgat = itime
52                exit
53             end if
54          end do
55          if (iv%info(synop)%proc_domain(1,n)) then
56             num_obs = num_obs + 1
57             write(ounit,'(2i8)') 1, ifgat
58             write(ounit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))')&
59                num_obs , 1, iv%info(synop)%id(n), &  ! Station
60                iv%info(synop)%lat(1,n), &       ! Latitude
61                iv%info(synop)%lon(1,n), &       ! Longitude
62                ob%synop(n)%p,           &       ! Obs Pressure
63                ob%synop(n)%u,           & 
64                iv%synop(n)%u%inv, iv%synop(n)%u%qc, iv%synop(n)%u%error, &
65                re%synop(n)%u,           &
66                ob%synop(n)%v,           &
67                iv%synop(n)%v%inv, iv%synop(n)%v%qc, iv%synop(n)%v%error, &
68                re%synop(n)%v,           &
69                ob%synop(n)%t,           &
70                iv%synop(n)%t%inv, iv%synop(n)%t%qc, iv%synop(n)%t%error, &
71                re%synop(n)%t,           &
72                ob%synop(n)%p,           &
73                iv%synop(n)%p%inv, iv%synop(n)%p%qc, iv%synop(n)%p%error, &
74                re%synop(n)%p,           &
75                ob%synop(n)%q,           &
76                iv%synop(n)%q%inv, iv%synop(n)%q%qc, iv%synop(n)%q%error, & 
77                re%synop(n)%q
78          end if
79       end do
80    end if
82   num_obs = 0
83   do n = 1, iv%info(metar)%nlocal
84      if (iv%info(metar)%proc_domain(1,n)) num_obs = num_obs + 1
85   end do
86   if (num_obs > 0) then
87      write(ounit,'(a20,i8)')'metar', num_obs  
88      num_obs = 0
89      do n = 1, iv%info(metar)%nlocal  
90          do itime = 1, num_fgat_time
91             if ( n >= iv%info(metar)%plocal(itime-1)+1 .and. &
92                  n <= iv%info(metar)%plocal(itime) ) then
93                ifgat = itime
94                exit
95             end if
96          end do
97         if (iv%info(metar)%proc_domain(1,n)) then
98            num_obs = num_obs + 1
99            write(ounit,'(2i8)') 1, ifgat
100            write(ounit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))')&
101               num_obs  , 1, iv%info(metar)%id(n), &  ! Station
102               iv%info(metar)%lat(1,n), &       ! Latitude
103               iv%info(metar)%lon(1,n), &       ! Longitude
104               ob%metar(n)%p,           &       ! Obs Pressure
105               ob%metar(n)%u,           &
106               iv%metar(n)%u%inv, iv%metar(n)%u%qc, iv%metar(n)%u%error, &
107               re%metar(n)%u,           &
108               ob%metar(n)%v,           &
109               iv%metar(n)%v%inv, iv%metar(n)%v%qc, iv%metar(n)%v%error, &
110               re%metar(n)%v,           &
111               ob%metar(n)%t,           &
112               iv%metar(n)%t%inv, iv%metar(n)%t%qc, iv%metar(n)%t%error, &
113               re%metar(n)%t,           &
114               ob%metar(n)%p,           &
115               iv%metar(n)%p%inv, iv%metar(n)%p%qc, iv%metar(n)%p%error, &
116               re%metar(n)%p,           &
117               ob%metar(n)%q,           &
118               iv%metar(n)%q%inv, iv%metar(n)%q%qc, iv%metar(n)%q%error, &
119               re%metar(n)%q
120         end if
121      end do
122   end if
124   num_obs = 0
125   do n = 1, iv%info(ships)%nlocal
126      if (iv%info(ships)%proc_domain(1,n)) num_obs = num_obs + 1
127   end do
128   if (num_obs > 0) then
129      write(ounit,'(a20,i8)')'ships', num_obs    
130      num_obs = 0
131      do n = 1, iv%info(ships)%nlocal  
132          do itime = 1, num_fgat_time
133             if ( n >= iv%info(ships)%plocal(itime-1)+1 .and. &
134                  n <= iv%info(ships)%plocal(itime) ) then
135                ifgat = itime
136                exit
137             end if
138          end do
139         if (iv%info(ships)%proc_domain(1,n)) then
140            write(ounit,'(2i8)') 1, ifgat
141            num_obs = num_obs + 1
142            write(ounit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))')&
143               num_obs,1, iv%info(ships)%id(n), &  ! Station
144               iv%info(ships)%lat(1,n), &       ! Latitude
145               iv%info(ships)%lon(1,n), &       ! Longitude
146               ob%ships(n)%p,           &       ! Obs Pressure
147               ob%ships(n)%u,           &
148               iv%ships(n)%u%inv, iv%ships(n)%u%qc, iv%ships(n)%u%error, &
149               re%ships(n)%u,           &
150               ob%ships(n)%v,           &
151               iv%ships(n)%v%inv, iv%ships(n)%v%qc, iv%ships(n)%v%error, &
152               re%ships(n)%v,           &
153               ob%ships(n)%t,           &
154               iv%ships(n)%t%inv, iv%ships(n)%t%qc, iv%ships(n)%t%error, &
155               re%ships(n)%t,           &
156               ob%ships(n)%p,           &
157               iv%ships(n)%p%inv, iv%ships(n)%p%qc, iv%ships(n)%p%error, &
158               re%ships(n)%p,           &
159               ob%ships(n)%q,           &
160               iv%ships(n)%q%inv, iv%ships(n)%q%qc, iv%ships(n)%q%error, &
161               re%ships(n)%q
162         end if
163      end do
164   end if
166   num_obs = 0
167   do n = 1, iv%info(geoamv)%nlocal
168      if (iv%info(geoamv)%proc_domain(1,n)) num_obs = num_obs + 1
169   end do
170   if (num_obs > 0) then
171      write(ounit,'(a20,i8)')'geoamv', num_obs    
172      num_obs = 0
173      do n = 1, iv%info(geoamv)%nlocal
174          do itime = 1, num_fgat_time
175             if ( n >= iv%info(geoamv)%plocal(itime-1)+1 .and. &
176                  n <= iv%info(geoamv)%plocal(itime) ) then
177                ifgat = itime
178                exit
179             end if
180          end do
181         if (iv%info(geoamv)%proc_domain(1,n)) then                  
182            num_obs = num_obs + 1
183            write(ounit,'(2i8)')iv%info(geoamv)%levels(n), ifgat
184            do k = 1, iv%info(geoamv)%levels(n)
185                write(ounit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))')&
186                   num_obs, 1, iv%info(geoamv)%id(n), &  ! Station
187                   iv%info(geoamv)%lat(k,n), &       ! Latitude
188                   iv%info(geoamv)%lon(k,n), &       ! Longitude
189                   iv%geoamv(n)%p(k),        &       ! Obs Pressure
190                   ob%geoamv(n)%u(k),        &
191                   iv%geoamv(n)%u(k)%inv, iv%geoamv(n)%u(k)%qc, iv%geoamv(n)%u(k)%error, &
192                   re%geoamv(n)%u(k),        &
193                   ob%geoamv(n)%v(k),        &
194                   iv%geoamv(n)%v(k)%inv, iv%geoamv(n)%v(k)%qc, iv%geoamv(n)%v(k)%error, &
195                   re%geoamv(n)%v(k)
196            end do
197         end if
198      end do
199   end if
201    num_obs = 0
202    do n = 1, iv%info(polaramv)%nlocal
203       if (iv%info(polaramv)%proc_domain(1,n)) num_obs = num_obs + 1
204    end do
205    if (num_obs > 0) then
206       write(ounit,'(a20,i8)')'polaramv', num_obs      
207       num_obs = 0
208       do n = 1, iv%info(polaramv)%nlocal
209          do itime = 1, num_fgat_time
210             if ( n >= iv%info(polaramv)%plocal(itime-1)+1 .and. &
211                  n <= iv%info(polaramv)%plocal(itime) ) then
212                ifgat = itime
213                exit
214             end if
215          end do
216          if (iv%info(polaramv)%proc_domain(1,n)) then                    
217             num_obs = num_obs + 1
218             write(ounit,'(2i8)')iv%info(polaramv)%levels(n), ifgat
219             do k = 1, iv%info(polaramv)%levels(n)
220                 write(ounit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))')&
221                    num_obs, 1, iv%info(polaramv)%id(n), &  ! Station
222                    iv%info(polaramv)%lat(k,n), &       ! Latitude
223                    iv%info(polaramv)%lon(k,n), &       ! Longitude
224                    iv%polaramv(n)%p(k),        &       ! Obs Pressure
225                    ob%polaramv(n)%u(k),        &
226                    iv%polaramv(n)%u(k)%inv, iv%polaramv(n)%u(k)%qc, iv%polaramv(n)%u(k)%error, &
227                    re%polaramv(n)%u(k),        &
228                    ob%polaramv(n)%v(k),        &
229                    iv%polaramv(n)%v(k)%inv, iv%polaramv(n)%v(k)%qc, iv%polaramv(n)%v(k)%error, &
230                    re%polaramv(n)%v(k)
231             end do
232          end if
233       end do
234    end if
236    num_obs = 0
237    do n = 1, iv%info(gpspw)%nlocal   
238       if (iv%info(gpspw)%proc_domain(1,n)) num_obs = num_obs + 1
239    end do
240    if (num_obs > 0) then
241       write(ounit,'(a20,i8)')'gpspw', num_obs    
242       num_obs = 0
243       do n = 1, iv%info(gpspw)%nlocal
244          do itime = 1, num_fgat_time
245             if ( n >= iv%info(gpspw)%plocal(itime-1)+1 .and. &
246                  n <= iv%info(gpspw)%plocal(itime) ) then
247                ifgat = itime
248                exit
249             end if
250          end do
251          if (iv%info(gpspw)%proc_domain(1,n)) then
252             num_obs = num_obs + 1
253             write(ounit,'(2i8)') 1, ifgat
254             write(ounit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))')&
255                num_obs, 1, iv%info(gpspw)%id(n), &  ! Station
256                iv%info(gpspw)%lat(1,n), &       ! Latitude
257                iv%info(gpspw)%lon(1,n), &       ! Longitude
258                iv%info(gpspw)%elv(n)  , &
259                ob%gpspw(n)%tpw,         &
260                iv%gpspw(n)%tpw%inv, iv%gpspw(n)%tpw%qc, iv%gpspw(n)%tpw%error, &
261                re%gpspw(n)%tpw
262          end if
263       end do
264    end if
266    num_obs = 0
267    do n = 1, iv%info(sound)%nlocal
268      if (iv%info(sound)%proc_domain(1,n)) num_obs = num_obs + 1
269    end do
270    if (num_obs > 0) then
271       write(ounit,'(a20,i8)')'sound', num_obs    
272       num_obs = 0
273       do n = 1, iv%info(sound)%nlocal
274          do itime = 1, num_fgat_time
275             if ( n >= iv%info(sound)%plocal(itime-1)+1 .and. &
276                  n <= iv%info(sound)%plocal(itime) ) then
277                ifgat = itime
278                exit
279             end if
280          end do
281          if (iv%info(sound)%proc_domain(1,n)) then
282             num_obs = num_obs + 1
283             write(ounit,'(2i8)')iv%info(sound)%levels(n), ifgat
284             do k = 1, iv%info(sound)%levels(n)
285                write(ounit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))')&
286                   num_obs,k, iv%info(sound)%id(n), &  ! Station
287                   iv%info(sound)%lat(k,n), &       ! Latitude
288                   iv%info(sound)%lon(k,n), &       ! Longitude
289                   iv%sound(n)%p(k),        &       ! Obs Pressure
290                   ob%sound(n)%u(k),        &
291                   iv%sound(n)%u(k)%inv, iv%sound(n)%u(k)%qc, iv%sound(n)%u(k)%error, &
292                   re%sound(n)%u(k),        &
293                   ob%sound(n)%v(k),        &
294                   iv%sound(n)%v(k)%inv, iv%sound(n)%v(k)%qc, iv%sound(n)%v(k)%error, &
295                   re%sound(n)%v(k),        &
296                   ob%sound(n)%t(k),        &
297                   iv%sound(n)%t(k)%inv, iv%sound(n)%t(k)%qc, iv%sound(n)%t(k)%error, &
298                   re%sound(n)%t(k),        &
299                   ob%sound(n)%q(k),        &
300                   iv%sound(n)%q(k)%inv, iv%sound(n)%q(k)%qc, iv%sound(n)%q(k)%error, &
301                   re%sound(n)%q(k)
302             end do
303          end if
304       end do
305    end if
307    if (num_obs > 0) then
308       write(ounit,'(a20,i8)')'sonde_sfc', num_obs    
309       num_obs = 0
310       do n = 1, iv%info(sonde_sfc)%nlocal
311          do itime = 1, num_fgat_time
312             if ( n >= iv%info(sonde_sfc)%plocal(itime-1)+1 .and. &
313                  n <= iv%info(sonde_sfc)%plocal(itime) ) then
314                ifgat = itime
315                exit
316             end if
317          end do
318          if (iv%info(sound)%proc_domain(1,n)) then 
319             num_obs = num_obs + 1
320             write(ounit,'(2i8)') 1, ifgat
321             write(ounit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))')&
322                num_obs , 1, iv%info(sonde_sfc)%id(n), &  ! Station
323                iv%info(sonde_sfc)%lat(1,n), &       ! Latitude
324                iv%info(sonde_sfc)%lon(1,n), &       ! Longitude
325                ob%sonde_sfc(n)%p,           &       ! Obs Pressure
326                ob%sonde_sfc(n)%u,           &
327                iv%sonde_sfc(n)%u%inv, iv%sonde_sfc(n)%u%qc, iv%sonde_sfc(n)%u%error, &
328                re%sonde_sfc(n)%u,           &
329                ob%sonde_sfc(n)%v,           &
330                iv%sonde_sfc(n)%v%inv, iv%sonde_sfc(n)%v%qc, iv%sonde_sfc(n)%v%error, & 
331                re%sonde_sfc(n)%v,           &
332                ob%sonde_sfc(n)%t,           &
333                iv%sonde_sfc(n)%t%inv, iv%sonde_sfc(n)%t%qc, iv%sonde_sfc(n)%t%error, &
334                re%sonde_sfc(n)%t,           &
335                ob%sonde_sfc(n)%p,           &
336                iv%sonde_sfc(n)%p%inv, iv%sonde_sfc(n)%p%qc, iv%sonde_sfc(n)%p%error, & 
337                re%sonde_sfc(n)%p,           &
338                ob%sonde_sfc(n)%q,           &
339                iv%sonde_sfc(n)%q%inv, iv%sonde_sfc(n)%q%qc, iv%sonde_sfc(n)%q%error, &
340                re%sonde_sfc(n)%q
341          end if
342       end do
343    end if
345    num_obs = 0
346    do n = 1, iv%info(airep)%nlocal
347       if (iv%info(airep)%proc_domain(1,n)) num_obs = num_obs + 1
348    end do
349    if (num_obs > 0) then
350       write(ounit,'(a20,i8)')'airep', num_obs  
351       num_obs = 0
352       do n = 1, iv%info(airep)%nlocal
353          do itime = 1, num_fgat_time
354             if ( n >= iv%info(airep)%plocal(itime-1)+1 .and. &
355                  n <= iv%info(airep)%plocal(itime) ) then
356                ifgat = itime
357                exit
358             end if
359          end do
360          if (iv%info(airep)%proc_domain(1,n)) then                  
361             num_obs = num_obs + 1
362             write(ounit,'(2i8)')iv%info(airep)%levels(n), ifgat
363             do k = 1, iv%info(airep)%levels(n)
364                write(ounit,'(2i8,a5,2f9.2,f17.7,4(2f17.7,i8,2f17.7),f12.2)')&
365                   num_obs, k, iv%info(airep)%id(n), &  ! Station
366                   iv%info(airep)%lat(k,n), &       ! Latitude
367                   iv%info(airep)%lon(k,n), &       ! Longitude
368                   iv%airep(n)%p(k),        &       ! Obs pressure
369                   ob%airep(n)%u(k),        &
370                   iv%airep(n)%u(k)%inv, iv%airep(n)%u(k)%qc, iv%airep(n)%u(k)%error, & 
371                   re%airep(n)%u(k),        &
372                   ob%airep(n)%v(k),        &
373                   iv%airep(n)%v(k)%inv, iv%airep(n)%v(k)%qc, iv%airep(n)%v(k)%error, & 
374                   re%airep(n)%v(k),        &
375                   ob%airep(n)%t(k),        &
376                   iv%airep(n)%t(k)%inv, iv%airep(n)%t(k)%qc, iv%airep(n)%t(k)%error, & 
377                   re%airep(n)%t(k),        &
378                   ob%airep(n)%q(k),        &
379                   iv%airep(n)%q(k)%inv, iv%airep(n)%q(k)%qc, iv%airep(n)%q(k)%error, &
380                   re%airep(n)%q(k), iv%info(airep)%zk(k,n)
381             end do
382          end if
383       end do
384    end if
386    num_obs = 0
387    do n = 1, iv%info(pilot)%nlocal
388       if (iv%info(pilot)%proc_domain(1,n)) num_obs = num_obs + 1
389    end do
390    if (num_obs > 0) then
391       write(ounit,'(a20,i8)')'pilot', num_obs   
392       num_obs = 0
393       do n = 1, iv%info(pilot)%nlocal
394          do itime = 1, num_fgat_time
395             if ( n >= iv%info(pilot)%plocal(itime-1)+1 .and. &
396                  n <= iv%info(pilot)%plocal(itime) ) then
397                ifgat = itime
398                exit
399             end if
400          end do
401          if (iv%info(pilot)%proc_domain(1,n)) then
402             num_obs = num_obs + 1
403             write(ounit,'(2i8)')iv%info(pilot)%levels(n), ifgat
404             do k = 1, iv%info(pilot)%levels(n)
405                write(ounit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))')&
406                   num_obs, k, iv%info(pilot)%id(n), &  ! Station
407                   iv%info(pilot)%lat(1,n), &       ! Latitude
408                   iv%info(pilot)%lon(1,n), &       ! Longitude
409                   iv%pilot(n)%h(k),        &       ! Obs Height
410                   ob%pilot(n)%u(k),        &
411                   iv%pilot(n)%u(k)%inv, iv%pilot(n)%u(k)%qc, iv%pilot(n)%u(k)%error, & 
412                   re%pilot(n)%u(k),        &
413                   ob%pilot(n)%v(k),        &
414                   iv%pilot(n)%v(k)%inv, iv%pilot(n)%v(k)%qc, iv%pilot(n)%v(k)%error, & 
415                   re%pilot(n)%v(k)
416             end do
417          end if
418       end do
419    end if
421    num_obs = 0
422    do n = 1, iv%info(ssmi_rv)%nlocal
423       if (iv%info(ssmi_rv)%proc_domain(1,n)) num_obs = num_obs + 1
424    end do
425    if (num_obs > 0) then
426       write(ounit,'(a20,i8)')'ssmir',  num_obs
427       num_obs = 0
428       do n = 1, iv%info(ssmi_rv)%nlocal
429          if (iv%info(ssmi_rv)%proc_domain(1,n)) then
430             num_obs = num_obs + 1
431             write(ounit,'(i8)') 1
432             write(ounit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))')&
433                num_obs, 1, 'SSMIR',              &       ! Station
434                iv%info(ssmi_rv)%lat(1,n), &! Latitude
435                iv%info(ssmi_rv)%lon(1,n), &! Longitude
436                missing_r,                 &       ! Obs height
437                ob%ssmi_rv(n)%speed,       &
438                iv%ssmi_rv(n)%speed%inv, iv%ssmi_rv(n)%speed%qc, iv%ssmi_rv(n)%speed%error, & 
439                re%ssmi_rv(n)%speed,       &
440                ob%ssmi_rv(n)%tpw,         &
441                iv%ssmi_rv(n)%tpw%inv, iv%ssmi_rv(n)%tpw%qc, iv%ssmi_rv(n)%tpw%error, & 
442                re%ssmi_rv(n)%tpw
443          end if
444       end do
445    end if
447    num_obs = 0
448    do n = 1, iv%info(ssmi_tb)%nlocal
449       if (iv%info(ssmi_tb)%proc_domain(1,n)) num_obs = num_obs + 1
450    end do
451    if (num_obs > 0) then
452       write(ounit,'(a20,i8)')'ssmiT', num_obs     
453       num_obs = 0
454       do n = 1, iv%info(ssmi_tb)%nlocal
455          ! write(ounit,*)' SSMI radiance output not yet coded.'
456          if (iv%info(ssmi_tb)%proc_domain(1,n)) then
457             num_obs = num_obs + 1
458             write(ounit,'(i8)') 1
459             write(ounit,'(2i8,a5,2f9.2,f17.7,7(2f17.7,i8,2f17.7))')&
460                num_obs, 1, 'SSMIT',              &        ! Station
461                iv%info(ssmi_tb)%lat(1,n), &! Latitude
462                iv%info(ssmi_tb)%lon(1,n), &! Longitude
463                missing_r,                 &       ! Obs height
464                ob%ssmi_tb(n)%tb19h,       &
465                iv%ssmi_tb(n)%tb19h%inv, iv%ssmi_tb(n)%tb19h%qc, iv%ssmi_tb(n)%tb19h%error, & 
466                re%ssmi_tb(n)%tb19h,       &
467                ob%ssmi_tb(n)%tb19v,       &
468                iv%ssmi_tb(n)%tb19v%inv, iv%ssmi_tb(n)%tb19v%qc, iv%ssmi_tb(n)%tb19v%error, & 
469                re%ssmi_tb(n)%tb19v,       &
470                ob%ssmi_tb(n)%tb22v,       &
471                iv%ssmi_tb(n)%tb22v%inv, iv%ssmi_tb(n)%tb22v%qc, iv%ssmi_tb(n)%tb22v%error, &
472                re%ssmi_tb(n)%tb22v,       &
473                ob%ssmi_tb(n)%tb37h,       &
474                iv%ssmi_tb(n)%tb37h%inv, iv%ssmi_tb(n)%tb37h%qc, iv%ssmi_tb(n)%tb37h%error, &
475                re%ssmi_tb(n)%tb37h,       &
476                ob%ssmi_tb(n)%tb37v,       &
477                iv%ssmi_tb(n)%tb37v%inv, iv%ssmi_tb(n)%tb37v%qc, iv%ssmi_tb(n)%tb37v%error, & 
478                re%ssmi_tb(n)%tb37v,       &
479                ob%ssmi_tb(n)%tb85h,       &
480                iv%ssmi_tb(n)%tb85h%inv, iv%ssmi_tb(n)%tb85h%qc, iv%ssmi_tb(n)%tb85h%error, & 
481                re%ssmi_tb(n)%tb85h,       &
482                ob%ssmi_tb(n)%tb85v,       &
483                iv%ssmi_tb(n)%tb85v%inv, iv%ssmi_tb(n)%tb85v%qc, iv%ssmi_tb(n)%tb85v%error, & 
484                re%ssmi_tb(n)%tb85v
485          end if
486       end do
487    end if
489    num_obs = 0
490    do n = 1, iv%info(satem)%nlocal
491       if (iv%info(satem)%proc_domain(1,n)) num_obs = num_obs + 1
492    end do
493    if (num_obs > 0) then
494       write(ounit,'(a20,i8)')'satem', num_obs    
495       num_obs = 0
496       do n = 1, iv%info(satem)%nlocal
497          if (iv%info(satem)%proc_domain(1,n)) then
498             num_obs = num_obs + 1
499             write(ounit,'(i8)')iv%info(satem)%levels(n)
500             do k = 1, iv%info(satem)%levels(n)
501                write(ounit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))')&
502                   num_obs   , k, iv%info(satem)%id(n), &  ! Station
503                   iv%info(satem)%lat(1,n), &       ! Latitude
504                   iv%info(satem)%lon(1,n), &       ! Longitude
505                   iv%satem(n)%p(k),        &       ! Obs Pressure
506                   ob%satem(n)%thickness(k),       &
507                   iv%satem(n)%thickness(k)%inv,   &
508                   iv%satem(n)%thickness(k)%qc,    &
509                   iv%satem(n)%thickness(k)%error, &
510                   re%satem(n)%thickness(k)
511             end do
512          end if
513       end do
514    end if
516    num_obs = 0
517    do n = 1, iv%info(ssmt1)%nlocal
518       if (iv%info(ssmt1)%proc_domain(1,n)) num_obs = num_obs + 1
519    end do
520    if (num_obs > 0) then
521       write(ounit,'(a20,i8)')'ssmt1', num_obs
522       num_obs = 0
523       do n = 1, iv%info(ssmt1)%nlocal
524          if (iv%info(ssmt1)%proc_domain(1,n)) then
525             num_obs = num_obs + 1
526             write(ounit,'(i8)')iv%info(ssmt1)%levels(n)
527             do k = 1, iv%info(ssmt1)%levels(n)
528                write(ounit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))')&
529                   num_obs , k, iv%info(ssmt1)%id(n), &  ! Station
530                   iv%info(ssmt1)%lat(1,n), &       ! Latitude
531                   iv%info(ssmt1)%lon(1,n), &       ! Longitude
532                   iv%ssmt1(n)%h(k),       &        ! Obs height
533                   ob%ssmt1(n)%t(k),       &
534                   iv%ssmt1(n)%t(k)%inv,   &
535                   iv%ssmt1(n)%t(k)%qc,    &
536                   iv%ssmt1(n)%t(k)%error, &
537                   re%ssmt1(n)%t(k)
538             end do
539          end if
540       end do
541    end if
543    num_obs = 0
544    do n = 1, iv%info(ssmt2)%nlocal
545       if (iv%info(ssmt2)%proc_domain(1,n)) num_obs = num_obs + 1
546    end do
547    if (num_obs > 0) then
548       write(ounit,'(a20,i8)')'ssmt2', num_obs    
549       num_obs = 0
550       do n = 1, iv%info(ssmt2)%nlocal
551          if (iv%info(ssmt2)%proc_domain(1,n)) then                   
552             num_obs = num_obs + 1
553             write(ounit,'(i8)')iv%info(ssmt2)%levels(n)
554             do k = 1, iv%info(ssmt2)%levels(n)
555                write(ounit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))')&
556                   num_obs , k, iv%info(ssmt2)%id(n), &  ! Station
557                   iv%info(ssmt2)%lat(1,n), &       ! Latitude
558                   iv%info(ssmt2)%lon(1,n), &       ! Longitude
559                   iv%ssmt2(n)%h(k),        &       ! Obs height
560                   ob%ssmt2(n)%rh(k),       &
561                   iv%ssmt2(n)%rh(k)%inv,   &
562                   iv%ssmt2(n)%rh(k)%qc,    &
563                   iv%ssmt2(n)%rh(k)%error, &
564                   re%ssmt2(n)%rh(k)
565             end do
566          end if
567       end do
568    end if
570    num_obs = 0
571    do n = 1, iv%info(qscat)%nlocal  
572       if (iv%info(qscat)%proc_domain(1,n)) num_obs = num_obs + 1
573    end do
574    if (num_obs > 0) then
575       write(ounit,'(a20,i8)')'qscat', num_obs   
576       num_obs = 0
577       do n = 1, iv%info(qscat)%nlocal  
578          do itime = 1, num_fgat_time
579             if ( n >= iv%info(qscat)%plocal(itime-1)+1 .and. &
580                  n <= iv%info(qscat)%plocal(itime) ) then
581                ifgat = itime
582                exit
583             end if
584          end do
585          if (iv%info(qscat)%proc_domain(1,n)) then
586             num_obs = num_obs + 1
587             write(ounit,'(2i8)') 1, ifgat
588             write(ounit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))')&
589                 num_obs, 1, iv%info(qscat)%id(n), &  ! Station
590                 iv%info(qscat)%lat(1,n), &       ! Latitude
591                 iv%info(qscat)%lon(1,n), &       ! Longitude
592                 iv%qscat(n)%h,           &       ! Obs height
593                 ob%qscat(n)%u,           &
594                 iv%qscat(n)%u%inv, iv%qscat(n)%u%qc, iv%qscat(n)%u%error, &
595                 re%qscat(n)%u,           &
596                 ob%qscat(n)%v,           &
597                 iv%qscat(n)%v%inv, iv%qscat(n)%v%qc, iv%qscat(n)%v%error, & 
598                 re%qscat(n)%v
599          end if
600       end do
601    end if
603    num_obs = 0
604    do n = 1, iv%info(profiler)%nlocal
605       if (iv%info(profiler)%proc_domain(1,n)) num_obs = num_obs + 1
606    end do
607    if (num_obs > 0) then
608       write(ounit,'(a20,i8)')'profiler',  num_obs
609       num_obs = 0
610       do n = 1, iv%info(profiler)%nlocal
611          do itime = 1, num_fgat_time
612             if ( n >= iv%info(profiler)%plocal(itime-1)+1 .and. &
613                  n <= iv%info(profiler)%plocal(itime) ) then
614                ifgat = itime
615                exit
616             end if
617          end do
618          if (iv%info(profiler)%proc_domain(1,n)) then
619             num_obs = num_obs + 1
620             write(ounit,'(2i8)')iv%info(profiler)%levels(n), ifgat
621             do k = 1, iv%info(profiler)%levels(n)
622                write(ounit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))')&
623                   num_obs, k, iv%info(profiler)%id(n), &  ! Station
624                   iv%info(profiler)%lat(1,n), &       ! Latitude
625                   iv%info(profiler)%lon(1,n), &       ! Longitude
626                   iv%profiler(n)%h(k),        &       ! Obs Height 
627                   ob%profiler(n)%u(k),        &
628                   iv%profiler(n)%u(k)%inv, iv%profiler(n)%u(k)%qc, iv%profiler(n)%u(k)%error, &
629                   re%profiler(n)%u(k),        &
630                   ob%profiler(n)%v(k),        &
631                   iv%profiler(n)%v(k)%inv, iv%profiler(n)%v(k)%qc, iv%profiler(n)%v(k)%error, &
632                   re%profiler(n)%v(k) 
633             end do
634          end if 
635       end do
636    end if
638    num_obs = 0
639    do n = 1, iv%info(buoy)%nlocal  
640       if (iv%info(buoy)%proc_domain(1,n)) num_obs = num_obs + 1
641    end do
642    if (num_obs > 0) then
643       write(ounit,'(a20,i8)')'buoy', num_obs
644       num_obs = 0
645       do n = 1, iv%info(buoy)%nlocal  
646          do itime = 1, num_fgat_time
647             if ( n >= iv%info(buoy)%plocal(itime-1)+1 .and. &
648                  n <= iv%info(buoy)%plocal(itime) ) then
649                ifgat = itime
650                exit
651             end if
652          end do
653          if (iv%info(buoy)%proc_domain(1,n)) then
654             num_obs = num_obs + 1
655             write(ounit,'(2i8)') 1, ifgat
656             write(ounit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))')&
657                num_obs,1, iv%info(buoy)%id(n), &  ! Station
658                iv%info(buoy)%lat(1,n), &       ! Latitude
659                iv%info(buoy)%lon(1,n), &       ! Longitude
660                ob%buoy(n)%p,           &       ! Obs Pressure
661                ob%buoy(n)%u,           &
662                iv%buoy(n)%u%inv, iv%buoy(n)%u%qc, iv%buoy(n)%u%error, &
663                re%buoy(n)%u,           &
664                ob%buoy(n)%v,           &
665                iv%buoy(n)%v%inv, iv%buoy(n)%v%qc, iv%buoy(n)%v%error, &
666                re%buoy(n)%v,           &
667                ob%buoy(n)%t,           &
668                iv%buoy(n)%t%inv, iv%buoy(n)%t%qc, iv%buoy(n)%t%error, &
669                re%buoy(n)%t,           &
670                ob%buoy(n)%p,           &
671                iv%buoy(n)%p%inv, iv%buoy(n)%p%qc, iv%buoy(n)%p%error, &
672                re%buoy(n)%p,           &
673                ob%buoy(n)%q,           &
674                iv%buoy(n)%q%inv, iv%buoy(n)%q%qc, iv%buoy(n)%q%error, &
675                re%buoy(n)%q
676          end if
677       end do
678    end if
680    num_obs = 0
681    do n = 1, iv%info(bogus)%nlocal
682       if (iv%info(bogus)%proc_domain(1,n)) num_obs = num_obs + 1
683    end do
684    if (num_obs > 0) then
685       write(ounit,'(a20,i8)')'bogus', num_obs
686       num_obs = 0
687       do n = 1, iv%info(bogus)%nlocal
688          if (iv%info(bogus)%proc_domain(1,n)) then
689             num_obs = num_obs + 1
690             write(ounit,'(i8)') 1
691             write(ounit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))')&
692                 num_obs, 1, iv%info(bogus)%id(n), &  ! Station
693                 iv%info(bogus)%lat(1,n), &       ! Latitude
694                 iv%info(bogus)%lon(1,n), &       ! Longitude
695                 missing_r,               &
696                 ob%bogus(n)%slp,         &
697                 iv%bogus(n)%slp%inv, iv%bogus(n)%slp%qc, iv%bogus(n)%slp%error, &
698                 re%bogus(n)%slp    ! O, O-B, O-A p
699             write(ounit,'(i8)')iv%info(bogus)%levels(n)
700             do k = 1, iv%info(bogus)%levels(n)
701                write(ounit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))')&
702                   num_obs , k, iv%info(bogus)%id(n), &  ! Station
703                   iv%info(bogus)%lat(1,n), &       ! Latitude
704                   iv%info(bogus)%lon(1,n), &       ! Longitude
705                   iv%bogus(n)%p(k),        &       ! Obs Pressure
706                   ob%bogus(n)%u(k),        &
707                   iv%bogus(n)%u(k)%inv, iv%bogus(n)%u(k)%qc, iv%bogus(n)%u(k)%error, &
708                   re%bogus(n)%u(k),        &
709                   ob%bogus(n)%v(k),        &
710                   iv%bogus(n)%v(k)%inv, iv%bogus(n)%v(k)%qc, iv%bogus(n)%v(k)%error, &
711                   re%bogus(n)%v(k),        &
712                   ob%bogus(n)%t(k),        &
713                   iv%bogus(n)%t(k)%inv, iv%bogus(n)%t(k)%qc, iv%bogus(n)%t(k)%error, &
714                   re%bogus(n)%t(k),        &
715                   ob%bogus(n)%q(k),        &
716                   iv%bogus(n)%q(k)%inv, iv%bogus(n)%q(k)%qc, iv%bogus(n)%q(k)%error, &
717                   re%bogus(n)%q(k)
718             end do
719          end if
720       end do
721    end if
723    num_obs = 0
724    do n = 1, iv%info(airsr)%nlocal
725       if (iv%info(airsr)%proc_domain(1,n)) num_obs = num_obs + 1
726    end do
727    if (num_obs > 0) then
728       write(ounit,'(a20,i8)')'airsr', num_obs    
729       num_obs = 0
730       do n = 1, iv%info(airsr)%nlocal
731          do itime = 1, num_fgat_time
732             if ( n >= iv%info(airsr)%plocal(itime-1)+1 .and. &
733                  n <= iv%info(airsr)%plocal(itime) ) then
734                ifgat = itime
735                exit
736             end if
737          end do
738          if (iv%info(airsr)%proc_domain(1,n)) then
739             num_obs = num_obs + 1
740             write(ounit,'(2i8)')iv%info(airsr)%levels(n), ifgat
741             do k = 1, iv%info(airsr)%levels(n)
742                write(ounit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))')&
743                   num_obs,k, iv%info(airsr)%id(n), &  ! Station
744                   iv%info(airsr)%lat(1,n), &       ! Latitude
745                   iv%info(airsr)%lon(1,n), &       ! Longitude
746                   iv%airsr(n)%p(k),        &       ! Obs Pressure
747                   ob%airsr(n)%t(k),        &
748                   iv%airsr(n)%t(k)%inv, iv%airsr(n)%t(k)%qc, iv%airsr(n)%t(k)%error, &
749                   re%airsr(n)%t(k),        &
750                   ob%airsr(n)%q(k),        &
751                   iv%airsr(n)%q(k)%inv, iv%airsr(n)%q(k)%qc, iv%airsr(n)%q(k)%error, &
752                   re%airsr(n)%q(k)
753             end do
754          end if
755       end do
756    end if
758    num_obs = 0
759    do n = 1, iv%info(gpsref)%nlocal
760       if (iv%info(gpsref)%proc_domain(1,n)) num_obs = num_obs + 1
761    end do
762    if (num_obs > 0) then
763       write(ounit,'(a20,i8)')'gpsref', num_obs
764        num_obs = 0
765       do n = 1, iv%info(gpsref)%nlocal
766          do itime = 1, num_fgat_time
767             if ( n >= iv%info(gpsref)%plocal(itime-1)+1 .and. &
768                  n <= iv%info(gpsref)%plocal(itime) ) then
769                ifgat = itime
770                exit
771             end if
772          end do
773          if (iv%info(gpsref)%proc_domain(1,n)) then
774             num_obs = num_obs + 1
775             write(ounit,'(2i8)')iv%info(gpsref)%levels(n), ifgat
776             do k = 1, iv%info(gpsref)%levels(n)
777                write(ounit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))')&
778                   num_obs,k, iv%info(gpsref)%id(n), &  ! Station
779                   iv%info(gpsref)%lat(1,n), &       ! Latitude
780                   iv%info(gpsref)%lon(1,n), &       ! Longitude
781                   iv%gpsref(n)%h(k),        &       ! Obs Height    
782                   ob%gpsref(n)%ref(k),      &
783                   iv%gpsref(n)%ref(k)%inv, iv%gpsref(n)%ref(k)%qc, iv%gpsref(n)%ref(k)%error, &
784                   re%gpsref(n)%ref(k)
785             end do
786          end if
787       end do
788    end if
790    num_obs = 0
791    do n = 1, iv%info(gpseph)%nlocal
792       if (iv%info(gpseph)%proc_domain(1,n)) num_obs = num_obs + 1
793    end do
794    if (num_obs > 0) then
795       write(ounit,'(a20,i8)')'gpseph', num_obs
796       num_obs = 0
797       do n = 1, iv%info(gpseph)%nlocal
798          do itime = 1, num_fgat_time
799             if ( n >= iv%info(gpseph)%plocal(itime-1)+1 .and. &
800                  n <= iv%info(gpseph)%plocal(itime) ) then
801                ifgat = itime
802                exit
803             end if
804          end do
805          if (iv%info(gpseph)%proc_domain(1,n)) then
806             num_obs = num_obs + 1
807             !write(ounit,'(i8)')iv%info(gpseph)%levels(n)
808             write(ounit,'(2i8)') iv%gpseph(n)%level2 - iv%gpseph(n)%level1 + 1, ifgat
809             !do k = 1, iv%info(gpseph)%levels(n)
810             do k = iv%gpseph(n)%level1, iv%gpseph(n)%level2
811                write(ounit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))')&
812                   num_obs, k, iv%info(gpseph)%id(n), &  ! Station
813                   iv%gpseph(n)%lat(k),      &           ! Latitude
814                   iv%gpseph(n)%lon(k),      &           ! Longitude
815                   iv%gpseph(n)%h(k),        &           ! Obs Height
816                   ob%gpseph(n)%eph(k),      &
817                   iv%gpseph(n)%eph(k)%inv, iv%gpseph(n)%eph(k)%qc, iv%gpseph(n)%eph(k)%error, &
818                   re%gpseph(n)%eph(k)
819             end do
820          end if
821       end do
822    end if
824    num_obs = 0
825    do n = 1, iv%info(tamdar)%nlocal
826       if (iv%info(tamdar)%proc_domain(1,n)) num_obs = num_obs + 1
827    end do
828    if (num_obs > 0) then
829       write(ounit,'(a20,i8)')'tamdar', num_obs
830       num_obs = 0
831       do n = 1, iv%info(tamdar)%nlocal
832          do itime = 1, num_fgat_time
833             if ( n >= iv%info(tamdar)%plocal(itime-1)+1 .and. &
834                  n <= iv%info(tamdar)%plocal(itime) ) then
835                ifgat = itime
836                exit
837             end if
838          end do
839          if (iv%info(tamdar)%proc_domain(1,n)) then
840             num_obs = num_obs + 1
841             write(ounit,'(2i8)')iv%info(tamdar)%levels(n), ifgat
842             do k = 1, iv%info(tamdar)%levels(n)
843                write(ounit,'(2i8,a5,2f9.2,f17.7,4(2f17.7,i8,2f17.7),f12.2)')&
844                   num_obs,k, iv%info(tamdar)%id(n), &  ! Station
845                   iv%info(tamdar)%lat(k,n), &       ! Latitude
846                   iv%info(tamdar)%lon(k,n), &       ! Longitude
847                   iv%tamdar(n)%p(k),        &       ! Obs Pressure
848                   ob%tamdar(n)%u(k),        &
849                   iv%tamdar(n)%u(k)%inv, iv%tamdar(n)%u(k)%qc, iv%tamdar(n)%u(k)%error, &
850                   re%tamdar(n)%u(k),        &
851                   ob%tamdar(n)%v(k),        &
852                   iv%tamdar(n)%v(k)%inv, iv%tamdar(n)%v(k)%qc, iv%tamdar(n)%v(k)%error, &
853                   re%tamdar(n)%v(k),        &
854                   ob%tamdar(n)%t(k),        &
855                   iv%tamdar(n)%t(k)%inv, iv%tamdar(n)%t(k)%qc, iv%tamdar(n)%t(k)%error, &
856                   re%tamdar(n)%t(k),        &
857                   ob%tamdar(n)%q(k),        &
858                   iv%tamdar(n)%q(k)%inv, iv%tamdar(n)%q(k)%qc, iv%tamdar(n)%q(k)%error, &
859                   re%tamdar(n)%q(k), iv%info(tamdar)%zk(k,n)
860             end do
861          end if
862       end do
863    end if
865 ! tamdar_sfc
867    num_obs = 0
868    do n = 1, iv%info(tamdar_sfc)%nlocal
869       if (iv%info(tamdar_sfc)%proc_domain(1,n)) num_obs = num_obs + 1
870    end do
871    if (num_obs > 0) then
872       write(ounit,'(a20,i8)')'tamdar_sfc', num_obs  
873       num_obs = 0
874       do n = 1, iv%info(tamdar_sfc)%nlocal  
875          do itime = 1, num_fgat_time
876             if ( n >= iv%info(tamdar_sfc)%plocal(itime-1)+1 .and. &
877                  n <= iv%info(tamdar_sfc)%plocal(itime) ) then
878                ifgat = itime
879                exit
880             end if
881          end do
882          if (iv%info(tamdar_sfc)%proc_domain(1,n)) then
883             num_obs = num_obs + 1
884             write(ounit,'(2i8)') 1, ifgat
885             write(ounit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7),f12.2)')&
886                num_obs , 1, iv%info(tamdar_sfc)%id(n), &  ! Station
887                iv%info(tamdar_sfc)%lat(1,n), &       ! Latitude
888                iv%info(tamdar_sfc)%lon(1,n), &       ! Longitude
889                ob%tamdar_sfc(n)%p,           &       ! Obs Pressure
890                ob%tamdar_sfc(n)%u,           & 
891                iv%tamdar_sfc(n)%u%inv, iv%tamdar_sfc(n)%u%qc, iv%tamdar_sfc(n)%u%error, &
892                re%tamdar_sfc(n)%u,           &
893                ob%tamdar_sfc(n)%v,           &
894                iv%tamdar_sfc(n)%v%inv, iv%tamdar_sfc(n)%v%qc, iv%tamdar_sfc(n)%v%error, &
895                re%tamdar_sfc(n)%v,           &
896                ob%tamdar_sfc(n)%t,           &
897                iv%tamdar_sfc(n)%t%inv, iv%tamdar_sfc(n)%t%qc, iv%tamdar_sfc(n)%t%error, &
898                re%tamdar_sfc(n)%t,           &
899                ob%tamdar_sfc(n)%p,           &
900                iv%tamdar_sfc(n)%p%inv, iv%tamdar_sfc(n)%p%qc, iv%tamdar_sfc(n)%p%error, &
901                re%tamdar_sfc(n)%p,           &
902                ob%tamdar_sfc(n)%q,           &
903                iv%tamdar_sfc(n)%q%inv, iv%tamdar_sfc(n)%q%qc, iv%tamdar_sfc(n)%q%error, & 
904                re%tamdar_sfc(n)%q, iv%info(tamdar_sfc)%zk(1,n)
905          end if
906       end do
907    end if
909    num_obs = 0
910    do n = 1, iv%info(rain)%nlocal  
911       if (iv%info(rain)%proc_domain(1,n)) num_obs = num_obs + 1
912    end do
913    if (num_obs > 0) then
914       write(ounit,'(a20,i8)')'rain', num_obs
915       num_obs = 0
916       do n = 1, iv%info(rain)%nlocal  
917          do itime = 1, num_fgat_time
918             if ( n >= iv%info(rain)%plocal(itime-1)+1 .and. &
919                  n <= iv%info(rain)%plocal(itime) ) then
920                ifgat = itime
921                exit
922             end if
923          end do
924          if (iv%info(rain)%proc_domain(1,n)) then
925             num_obs = num_obs + 1
926             write(ounit,'(2i8)') 1, ifgat
927             write(ounit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))')&
928                num_obs,1, iv%info(rain)%id(n), &  ! Station  
929                iv%info(rain)%lat(1,n), &          ! Latitude
930                iv%info(rain)%lon(1,n), &          ! Longitude
931                iv%rain(n)%height,      &          ! Obs Pressure
932                ob%rain(n)%rain,        &
933                iv%rain(n)%rain%inv, iv%rain(n)%rain%qc, iv%rain(n)%rain%error, &
934                re%rain(n)%rain
935          end if
936       end do
937    end if
939    !! lightning
940    num_obs = 0
941    do n = 1, iv%info(lightning)%nlocal
942       if (iv%info(lightning)%proc_domain(1,n)) num_obs = num_obs + 1
943    end do
944    if (num_obs > 0) then
945       write(ounit,'(a20,i8)')'lightning', num_obs    
946       num_obs = 0
947       do n = 1, iv%info(lightning)%nlocal
948          do itime = 1, num_fgat_time
949             if ( n >= iv%info(lightning)%plocal(itime-1)+1 .and. &
950                  n <= iv%info(lightning)%plocal(itime) ) then
951                ifgat = itime
952                exit
953             end if
954          end do
955          if (iv%info(lightning)%proc_domain(1,n)) then
956             num_obs = num_obs + 1
957             write(ounit,'(2i8)')iv%info(lightning)%levels(n), ifgat
958             do k = 1, iv%info(lightning)%levels(n)
959                write(ounit,'(2i8,a5,2f9.2,f17.7,3(2f17.7,i8,2f17.7))')&
960                   num_obs   , k, iv%info(lightning)%id(n), &  ! Station
961                   iv%info(lightning)%lat(1,n), &       ! Latitude
962                   iv%info(lightning)%lon(1,n), &       ! Longitude
963                   iv%lightning(n)%height(k),   &       ! Obs Height
964                   ob%lightning(n)%w(k),        &
965                   iv%lightning(n)%w(k)%inv,    &
966                   iv%lightning(n)%w(k)%qc,     &
967                   iv%lightning(n)%w(k)%error,  &
968                   re%lightning(n)%w(k),        &
969                   ob%lightning(n)%div(k),      &
970                   iv%lightning(n)%div(k)%inv,  &
971                   iv%lightning(n)%div(k)%qc,   &
972                   iv%lightning(n)%div(k)%error,&
973                   re%lightning(n)%div(k),      &
974                   ob%lightning(n)%qv(k),       &
975                   iv%lightning(n)%qv(k)%inv,   &
976                   iv%lightning(n)%qv(k)%qc,    &
977                   iv%lightning(n)%qv(k)%error, &
978                   re%lightning(n)%qv(k)                           
979             end do
980          end if
981       end do
982    end if   
983    
984    close (ounit)
985    call da_free_unit(ounit)
987    if (trace_use) call da_trace_exit("da_write_obs")
989 end subroutine da_write_obs