da_write_obs.inc

References to this file elsewhere.
1 subroutine da_write_obs(ob, iv, re)
2 
3    !-------------------------------------------------------------------------
4    ! Purpose: Writes out components of iv=O-B structure.
5    !-------------------------------------------------------------------------   
6 
7    implicit none
8 
9    type (y_type), intent(in)     :: ob      ! Observation structure.
10    type (ob_type), intent(in)    :: iv      ! O-B structure.
11    type (y_type), intent(inout)  :: re      ! residual vector.
12       
13    integer                       :: n, k, num_obs, ios
14    integer                       :: ounit     ! Output unit           
15    character(len=filename_len)   :: filename
16 
17    if (trace_use) call da_trace_entry("da_write_obs")
18 
19    !-------------------------------------------------------------------------
20    ! Fix output unit
21    !-------------------------------------------------------------------------
22 
23    call da_get_unit(ounit)
24 
25 #ifdef DM_PARALLEL
26     write(unit=filename, fmt='(a,i3.3)') 'gts_omb_oma.', myproc
27 #else
28     write(unit=filename, fmt='(a)') 'gts_omb_oma.000'
29 #endif
30 
31    open (unit=ounit,file=trim(filename),form='formatted',status='replace', &
32       iostat=ios)
33    if (ios /= 0) then
34       call da_error(__FILE__,__LINE__, &
35          (/"Cannot open conventional observation omb and oma file"//filename/))
36    end if
37 
38 
39    ! [1] Transfer surface obs:
40 
41    if (iv % num_synop > 0) then
42       num_obs = 0
43       do n = 1, iv % num_synop
44          if (iv%synop(n)%loc%proc_domain) num_obs = num_obs + 1
45       end do
46       if (num_obs > 0) then
47          write(ounit,'(a20,i8)')'synop', num_obs  
48          num_obs = 0
49          do n = 1, iv % num_synop  
50             if (iv%synop(n)%loc%proc_domain) then
51                num_obs = num_obs + 1
52                write(ounit,'(i8)') 1                                 
53                write(ounit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))')&
54                   num_obs , 1, iv % synop(n) % info % id, &  ! Station
55                   iv % synop(n) % info % lat, &       ! Latitude
56                   iv % synop(n) % info % lon, &       ! Longitude
57                   ob % synop(n) % p, &                ! Obs Pressure
58                   ob%synop(n)%u, iv%synop(n)%u, re%synop(n)%u, &! O, O-B, O-A u
59                   ob%synop(n)%v, iv%synop(n)%v, re%synop(n)%v, &! O, O-B, O-A v
60                   ob%synop(n)%t, iv%synop(n)%t, re%synop(n)%t, &! O, O-B, O-A t
61                   ob%synop(n)%p, iv%synop(n)%p, re%synop(n)%p, &! O, O-B, O-A p
62                   ob%synop(n)%q, iv%synop(n)%q, re%synop(n)%q   ! O, O-B, O-A q
63             end if
64          end do
65       end if
66    end if
67 
68    ! [2] Transfer metar obs:
69 
70    if (iv % num_metar > 0) then
71       num_obs = 0
72       do n = 1, iv % num_metar
73          if (iv%metar(n)%loc%proc_domain) num_obs = num_obs + 1
74       end do
75       if (num_obs > 0) then
76       write(ounit,'(a20,i8)')'metar', num_obs  
77       num_obs = 0
78       do n = 1, iv % num_metar  
79        if (iv%metar(n)%loc%proc_domain) then
80        num_obs = num_obs + 1
81        write(ounit,'(i8)') 1                                 
82          write(ounit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))')&
83                  num_obs  , 1, iv % metar(n) % info % id, &  ! Station
84                          iv % metar(n) % info % lat, &       ! Latitude
85                          iv % metar(n) % info % lon, &       ! Longitude
86                          ob % metar(n) % p, &                ! Obs Pressure
87                          ob%metar(n)%u, iv%metar(n)%u, re%metar(n)%u, &! O, O-B, O-A u
88                          ob%metar(n)%v, iv%metar(n)%v, re%metar(n)%v, &! O, O-B, O-A v
89                          ob%metar(n)%t, iv%metar(n)%t, re%metar(n)%t, &! O, O-B, O-A t
90                          ob%metar(n)%p, iv%metar(n)%p, re%metar(n)%p, &! O, O-B, O-A p
91                          ob%metar(n)%q, iv%metar(n)%q, re%metar(n)%q   ! O, O-B, O-A q
92        end if
93       end do
94      end if
95    end if
96 
97    ! [3] Transfer ships obs:
98 
99    if (iv % num_ships > 0) then
100       num_obs = 0
101       do n = 1, iv % num_ships
102        if (iv%ships(n)%loc%proc_domain) num_obs = num_obs + 1
103       end do
104       if (num_obs > 0) then
105       write(ounit,'(a20,i8)')'ships', num_obs    
106       num_obs = 0
107       do n = 1, iv % num_ships  
108        if (iv%ships(n)%loc%proc_domain) then
109        write(ounit,'(i8)') 1                                 
110        num_obs = num_obs + 1
111          write(ounit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))')&
112                     num_obs,1, iv % ships(n) % info % id, &  ! Station
113                          iv % ships(n) % info % lat, &       ! Latitude
114                          iv % ships(n) % info % lon, &       ! Longitude
115                          ob % ships(n) % p, &                ! Obs Pressure
116                          ob%ships(n)%u, iv%ships(n)%u, re%ships(n)%u, &! O, O-B, O-A u
117                          ob%ships(n)%v, iv%ships(n)%v, re%ships(n)%v, &! O, O-B, O-A v
118                          ob%ships(n)%t, iv%ships(n)%t, re%ships(n)%t, &! O, O-B, O-A t
119                          ob%ships(n)%p, iv%ships(n)%p, re%ships(n)%p, &! O, O-B, O-A p
120                          ob%ships(n)%q, iv%ships(n)%q, re%ships(n)%q   ! O, O-B, O-A q
121       end if
122       end do
123      end if
124    end if
125 
126    ! [4.1] Transfer Geo AMVs Obs:
127 
128    if (iv % num_geoamv > 0) then
129       num_obs = 0
130       do n = 1, iv % num_geoamv
131          if (iv%geoamv(n)%loc%proc_domain) num_obs = num_obs + 1
132       end do
133       if (num_obs > 0) then
134          write(ounit,'(a20,i8)')'geoamv', num_obs    
135          num_obs = 0
136          do n = 1, iv % num_geoamv
137             if (iv%geoamv(n)%loc%proc_domain) then                  
138                num_obs = num_obs + 1
139                write(ounit,'(i8)')iv % geoamv(n) % info % levels
140                do k = 1, iv % geoamv(n) % info % levels
141                    write(ounit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))')&
142                       num_obs, 1, iv % geoamv(n) % info % id, &  ! Station
143                       iv % geoamv(n) % info % lat, &       ! Latitude
144                       iv % geoamv(n) % info % lon, &       ! Longitude
145                       iv % geoamv(n) % p(k), &                ! Obs Pressure
146                       ob%geoamv(n)%u(k), iv%geoamv(n)%u(k), re%geoamv(n)%u(k), &! O, O-B, O-A u
147                       ob%geoamv(n)%v(k), iv%geoamv(n)%v(k), re%geoamv(n)%v(k)
148                end do
149             end if
150          end do
151       end if
152    end if
153 
154    ! [4.2] Transfer Polar AMVs Obs:
155 
156    if (iv % num_polaramv > 0) then
157       num_obs = 0
158       do n = 1, iv % num_polaramv
159          if (iv%polaramv(n)%loc%proc_domain) num_obs = num_obs + 1
160       end do
161       if (num_obs > 0) then
162          write(ounit,'(a20,i8)')'polaramv', num_obs      
163          num_obs = 0
164          do n = 1, iv % num_polaramv
165             if (iv%polaramv(n)%loc%proc_domain) then                    
166                num_obs = num_obs + 1
167                write(ounit,'(i8)')iv % polaramv(n) % info % levels
168                do k = 1, iv % polaramv(n) % info % levels
169                    write(ounit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))')&
170                       num_obs, 1, iv % polaramv(n) % info % id, &  ! Station
171                       iv % polaramv(n) % info % lat, &       ! Latitude
172                       iv % polaramv(n) % info % lon, &       ! Longitude
173                       iv % polaramv(n) % p(k), &                ! Obs Pressure
174                       ob%polaramv(n)%u(k), iv%polaramv(n)%u(k), re%polaramv(n)%u(k), &! O, O-B, O-A u
175                       ob%polaramv(n)%v(k), iv%polaramv(n)%v(k), re%polaramv(n)%v(k)
176                end do
177             end if
178          end do
179       end if
180    end if
181 
182    ! [5] Transfer gpspw obs:
183 
184    if (iv % num_gpspw > 0) then
185       num_obs = 0
186       do n = 1, iv % num_gpspw   
187          if (iv%gpspw(n)%loc%proc_domain) num_obs = num_obs + 1
188       end do
189       if (num_obs > 0) then
190          write(ounit,'(a20,i8)')'gpspw', num_obs    
191          num_obs = 0
192          do n = 1, iv % num_gpspw
193             if (iv%gpspw(n)%loc%proc_domain) then
194                num_obs = num_obs + 1
195                write(ounit,'(i8)') 1                                 
196                write(ounit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))')&
197                   num_obs, 1, iv % gpspw(n) % info % id, &  ! Station
198                   iv % gpspw(n) % info % lat, &       ! Latitude
199                   iv % gpspw(n) % info % lon, &       ! Longitude
200                   missing_r, &
201                   ob%gpspw(n)%tpw, iv%gpspw(n)%tpw, re%gpspw(n)%tpw ! O, O-B, O-A u
202             end if
203          end do
204       end if
205    end if
206 
207    ! [6] Transfer sonde obs:
208 
209    if (iv % num_sound > 0) then
210       num_obs = 0
211       do n = 1, iv % num_sound
212         if (iv%sound(n)%loc%proc_domain) num_obs = num_obs + 1
213       end do
214       if (num_obs > 0) then
215          write(ounit,'(a20,i8)')'sound', num_obs    
216          num_obs = 0
217          do n = 1, iv % num_sound
218             if (iv%sound(n)%loc%proc_domain) then
219                num_obs = num_obs + 1
220                write(ounit,'(i8)')iv % sound(n) % info % levels
221                do k = 1, iv % sound(n) % info % levels
222                   write(ounit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))')&
223                      num_obs,k, iv % sound(n) % info % id, &  ! Station
224                      iv % sound(n) % info % lat, &       ! Latitude
225                      iv % sound(n) % info % lon, &       ! Longitude
226                      iv % sound(n) % p(k), &             ! Obs Pressure
227                      ob%sound(n)%u(k), iv%sound(n)%u(k), re%sound(n)%u(k), &! O,O-B,O-A u
228                      ob%sound(n)%v(k), iv%sound(n)%v(k), re%sound(n)%v(k), &! O,O-B,O-A v
229                      ob%sound(n)%t(k), iv%sound(n)%t(k), re%sound(n)%t(k), &! O,O-B,O-A t
230                      ob%sound(n)%q(k), iv%sound(n)%q(k), re%sound(n)%q(k)   ! O,O-B,O-A q
231                end do
232             end if
233          end do
234       end if
235    end if
236 
237    if (iv % num_sound > 0) then
238       if (num_obs > 0) then
239          write(ounit,'(a20,i8)')'sonde_sfc', num_obs    
240          num_obs = 0
241          do n = 1, iv % num_sound
242             if (iv%sound(n)%loc%proc_domain) then 
243                num_obs = num_obs + 1
244                write(ounit,'(i8)') 1
245                write(ounit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))')&
246                   num_obs , 1, iv % sonde_sfc(n) % info % id, &  ! Station
247                   iv % sonde_sfc(n) % info % lat, &       ! Latitude
248                   iv % sonde_sfc(n) % info % lon, &       ! Longitude
249                   ob % sonde_sfc(n) % p, &                ! Obs Pressure
250                   ob%sonde_sfc(n)%u, iv%sonde_sfc(n)%u, re%sonde_sfc(n)%u, &! O, O-B, O-A u
251                   ob%sonde_sfc(n)%v, iv%sonde_sfc(n)%v, re%sonde_sfc(n)%v, &! O, O-B, O-A v
252                   ob%sonde_sfc(n)%t, iv%sonde_sfc(n)%t, re%sonde_sfc(n)%t, &! O, O-B, O-A t
253                   ob%sonde_sfc(n)%p, iv%sonde_sfc(n)%p, re%sonde_sfc(n)%p, &! O, O-B, O-A p
254                   ob%sonde_sfc(n)%q, iv%sonde_sfc(n)%q, re%sonde_sfc(n)%q   ! O, O-B, O-A q
255             end if
256          end do
257       end if
258    end if
259 
260    ! [7] Transfer airep obs:
261 
262    if (iv % num_airep > 0) then
263       num_obs = 0
264       do n = 1, iv % num_airep
265          if (iv%airep(n)%loc%proc_domain) num_obs = num_obs + 1
266       end do
267       if (num_obs > 0) then
268          write(ounit,'(a20,i8)')'airep', num_obs  
269          num_obs = 0
270          do n = 1, iv % num_airep
271             if (iv%airep(n)%loc%proc_domain) then                  
272                num_obs = num_obs + 1
273                write(ounit,'(i8)')iv % airep(n) % info % levels
274                do k = 1, iv % airep(n) % info % levels
275                   write(ounit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))')&
276                      num_obs, k, iv % airep(n) % info % id, &  ! Station
277                      iv % airep(n) % info % lat, &       ! Latitude
278                      iv % airep(n) % info % lon, &       ! Longitude
279                      iv % airep(n) % p(k), &             ! Obs pressure
280                      ob%airep(n)%u(k), iv%airep(n)%u(k), re%airep(n)%u(k), &! O,O-B,O-A u
281                      ob%airep(n)%v(k), iv%airep(n)%v(k), re%airep(n)%v(k), &! O,O-B,O-A v
282                      ob%airep(n)%t(k), iv%airep(n)%t(k), re%airep(n)%t(k)
283                end do
284             end if
285          end do
286       end if
287    end if
288 
289    ! [8] Transfer pilot obs:
290 
291    if (iv % num_pilot > 0) then
292       num_obs = 0
293       do n = 1, iv % num_pilot
294          if (iv%pilot(n)%loc%proc_domain) num_obs = num_obs + 1
295       end do
296       if (num_obs > 0) then
297          write(ounit,'(a20,i8)')'pilot', num_obs   
298          num_obs = 0
299          do n = 1, iv % num_pilot
300             if (iv%pilot(n)%loc%proc_domain) then
301                num_obs = num_obs + 1
302                write(ounit,'(i8)')iv % pilot(n) % info % levels
303                do k = 1, iv % pilot(n) % info % levels
304                   write(ounit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))')&
305                      num_obs, k, iv % pilot(n) % info % id, &  ! Station
306                      iv % pilot(n) % info % lat, &       ! Latitude
307                      iv % pilot(n) % info % lon, &       ! Longitude
308                      iv % pilot(n) % p(k), &             ! Obs Pressure
309                      ob%pilot(n)%u(k), iv%pilot(n)%u(k), re%pilot(n)%u(k), &! O,O-B,O-A u
310                      ob%pilot(n)%v(k), iv%pilot(n)%v(k), re%pilot(n)%v(k) ! O,O-B,O-A v
311                end do
312             end if
313          end do
314       end if
315    end if
316 
317    ! [9] Transfer SSM/I obs:SSMI:
318 
319    if (iv % num_ssmi_retrieval > 0) then
320       num_obs = 0
321       do n = 1, iv % num_ssmi_retrieval
322          if (iv%ssmi_retrieval(n)%loc%proc_domain) num_obs = num_obs + 1
323       end do
324       if (num_obs > 0) then
325          write(ounit,'(a20,i8)')'ssmir',  num_obs
326          num_obs = 0
327          do n = 1, iv % num_ssmi_retrieval
328             if (iv%ssmi_retrieval(n)%loc%proc_domain) then
329                num_obs = num_obs + 1
330                write(ounit,'(i8)') 1
331                write(ounit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))')&
332                   num_obs, 1, 'SSMIR',              &       ! Station
333                   iv % ssmi_retrieval(n) % info % lat, &! Latitude
334                   iv % ssmi_retrieval(n) % info % lon, &! Longitude
335                   missing_r,                  &       ! Obs height
336                   ob%ssmi_retrieval(n)%speed, iv%ssmi_retrieval(n)%speed, re%ssmi_retrieval(n)%speed, & ! O, O-B, O-A speed
337                   ob%ssmi_retrieval(n)%tpw, iv%ssmi_retrieval(n)%tpw, re%ssmi_retrieval(n)%tpw ! O, O-B, O-A tpw
338             end if
339          end do
340       end if
341    end if
342 
343    if (iv % num_ssmi_tb > 0) then
344       num_obs = 0
345       do n = 1, iv % num_ssmi_tb
346          if (iv%ssmi_tb(n)%loc%proc_domain) num_obs = num_obs + 1
347       end do
348       if (num_obs > 0) then
349          write(ounit,'(a20,i8)')'ssmiT', num_obs     
350          num_obs = 0
351          do n = 1, iv % num_ssmi_tb
352             write(ounit,*)' SSMI radiance output not yet coded.'
353             if (iv%ssmi_tb(n)%loc%proc_domain) then
354                num_obs = num_obs + 1
355                write(ounit,'(i8)') 1
356                write(ounit,'(2i8,a5,2f9.2,f17.7,7(2f17.7,i8,2f17.7))')&
357                   num_obs, 1, 'SSMIT',              &        ! Station
358                   iv % ssmi_tb(n) % info % lat, &! Latitude
359                   iv % ssmi_tb(n) % info % lon, &! Longitude
360                   missing_r,                  &       ! Obs height
361                   ob%ssmi_tb(n)%tb19h, iv%ssmi_tb(n)%tb19h, re%ssmi_tb(n)%tb19h, & ! O, O-B, O-A Tb19h
362                   ob%ssmi_tb(n)%tb19v, iv%ssmi_tb(n)%tb19v, re%ssmi_tb(n)%tb19v, & ! O, O-B, O-A Tb19v
363                   ob%ssmi_tb(n)%tb22v, iv%ssmi_tb(n)%tb22v, re%ssmi_tb(n)%tb22v, & ! O, O-B, O-A Tb22v
364                   ob%ssmi_tb(n)%tb37h, iv%ssmi_tb(n)%tb37h, re%ssmi_tb(n)%tb37h, & ! O, O-B, O-A Tb37h
365                   ob%ssmi_tb(n)%tb37v, iv%ssmi_tb(n)%tb37v, re%ssmi_tb(n)%tb37v, & ! O, O-B, O-A Tb37v
366                   ob%ssmi_tb(n)%tb85h, iv%ssmi_tb(n)%tb85h, re%ssmi_tb(n)%tb85h, & ! O, O-B, O-A Tb85h
367                   ob%ssmi_tb(n)%tb85v, iv%ssmi_tb(n)%tb85v, re%ssmi_tb(n)%tb85v    ! O, O-B, O-A Tb85v
368             end if
369          end do
370       end if
371    end if
372 
373    ! [10] Transfer satem obs:
374 
375    if (iv % num_satem > 0) then
376       num_obs = 0
377       do n = 1, iv % num_satem
378          if (iv%satem(n)%loc%proc_domain) num_obs = num_obs + 1
379       end do
380       if (num_obs > 0) then
381          write(ounit,'(a20,i8)')'satem', num_obs    
382          num_obs = 0
383          do n = 1, iv % num_satem
384             if (iv%satem(n)%loc%proc_domain) then
385                num_obs = num_obs + 1
386                write(ounit,'(i8)')iv % satem(n) % info % levels
387                do k = 1, iv % satem(n) % info % levels
388                   write(ounit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))')&
389                      num_obs   , k, iv % satem(n) % info % id, &  ! Station
390                      iv % satem(n) % info % lat, &       ! Latitude
391                      iv % satem(n) % info % lon, &       ! Longitude
392                      iv % satem(n) % p(k), &             ! Obs Pressure
393                      ob%satem(n)%thickness(k), &
394                      iv%satem(n)%thickness(k), &
395                      re%satem(n)%thickness(k)            ! O,O-B,O-A
396                end do
397             end if
398          end do
399       end if
400    end if
401    
402    ! [11] Transfer ssmt1 obs:
403 
404    if (iv % num_ssmt1 > 0) then
405       num_obs = 0
406       do n = 1, iv % num_ssmt1
407          if (iv%ssmt1(n)%loc%proc_domain) num_obs = num_obs + 1
408       end do
409       if (num_obs > 0) then
410          write(ounit,'(a20,i8)')'ssmt1', num_obs
411          num_obs = 0
412          do n = 1, iv % num_ssmt1
413             if (iv%ssmt1(n)%loc%proc_domain) then
414                num_obs = num_obs + 1
415                write(ounit,'(i8)')iv % ssmt1(n) % info % levels
416                do k = 1, iv % ssmt1(n) % info % levels
417                   write(ounit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))')&
418                      num_obs , k, iv % ssmt1(n) % info % id, &  ! Station
419                      iv % ssmt1(n) % info % lat, &       ! Latitude
420                      iv % ssmt1(n) % info % lon, &       ! Longitude
421                      iv % ssmt1(n) % h(k), &             ! Obs height
422                      ob%ssmt1(n)%t(k), &
423                      iv%ssmt1(n)%t(k), &
424                      re%ssmt1(n)%t(k)                    ! O,O-B,O-A u
425                end do
426             end if
427          end do
428       end if
429    end if
430 
431    ! [12] Transfer ssmt2 obs:
432 
433    if (iv % num_ssmt2 > 0) then
434       num_obs = 0
435       do n = 1, iv % num_ssmt2
436          if (iv%ssmt2(n)%loc%proc_domain) num_obs = num_obs + 1
437       end do
438       if (num_obs > 0) then
439          write(ounit,'(a20,i8)')'ssmt2', num_obs    
440          num_obs = 0
441          do n = 1, iv % num_ssmt2
442             if (iv%ssmt2(n)%loc%proc_domain) then                   
443                num_obs = num_obs + 1
444                write(ounit,'(i8)')iv % ssmt2(n) % info % levels
445                do k = 1, iv % ssmt2(n) % info % levels
446                   write(ounit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))')&
447                      num_obs , k, iv % ssmt2(n) % info % id, &  ! Station
448                      iv % ssmt2(n) % info % lat, &       ! Latitude
449                      iv % ssmt2(n) % info % lon, &       ! Longitude
450                      iv % ssmt2(n) % h(k), &             ! Obs height
451                      ob%ssmt2(n)%rh(k), &
452                      iv%ssmt2(n)%rh(k), &
453                      re%ssmt2(n)%rh(k)                   ! O,O-B,O-A u
454                end do
455             end if
456          end do
457       end if
458    end if
459    
460    ! [3] Transfer scatterometer obs:
461 
462    if (iv % num_qscat > 0) then
463       num_obs = 0
464       do n = 1, iv % num_qscat  
465          if (iv%qscat(n)%loc%proc_domain) num_obs = num_obs + 1
466       end do
467       if (num_obs > 0) then
468          write(ounit,'(a20,i8)')'qscat', num_obs   
469          num_obs = 0
470          do n = 1, iv % num_qscat  
471             if (iv%qscat(n)%loc%proc_domain) then
472                num_obs = num_obs + 1
473                write(ounit,'(i8)') 1
474                write(ounit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))')&
475                    num_obs, 1, iv % qscat(n) % info % id, &  ! Station
476                    iv % qscat(n) % info % lat, &       ! Latitude
477                    iv % qscat(n) % info % lon, &       ! Longitude
478                    iv % qscat(n) % h, &                ! Obs height
479                    ob%qscat(n)%u, iv%qscat(n)%u, re%qscat(n)%u, &! O, O-B, O-A u
480                    ob%qscat(n)%v, iv%qscat(n)%v, re%qscat(n)%v   ! O, O-B, O-A v
481             end if
482          end do
483       end if
484    end if
485 
486    ! Transfer profiler obs:
487 
488    if (iv % num_profiler > 0) then
489       num_obs = 0
490       do n = 1, iv % num_profiler
491          if (iv%profiler(n)%loc%proc_domain) num_obs = num_obs + 1
492       end do
493       if (num_obs > 0) then
494          write(ounit,'(a20,i8)')'profiler',  num_obs
495          num_obs = 0
496          do n = 1, iv % num_profiler
497             if (iv%profiler(n)%loc%proc_domain) then
498                num_obs = num_obs + 1
499                write(ounit,'(i8)')iv % profiler(n) % info % levels
500                do k = 1, iv % profiler(n) % info % levels
501                   write(ounit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))')&
502                      num_obs, k, iv % profiler(n) % info % id, &  ! Station
503                      iv % profiler(n) % info % lat, &       ! Latitude
504                      iv % profiler(n) % info % lon, &       ! Longitude
505                      iv % profiler(n) % p(k), &             ! Obs Pressure
506                      ob%profiler(n)%u(k), iv%profiler(n)%u(k), re%profiler(n)%u(k), &! O,O-B,O-A u
507                      ob%profiler(n)%v(k), iv%profiler(n)%v(k), re%profiler(n)%v(k) ! O,O-B,O-A v
508                end do
509             end if 
510          end do
511       end if
512    end if
513 
514    ! Transfer Buoy obs:
515 
516    if (iv % num_buoy > 0) then
517       num_obs = 0
518       do n = 1, iv % num_buoy  
519          if (iv%buoy(n)%loc%proc_domain) num_obs = num_obs + 1
520       end do
521       if (num_obs > 0) then
522          write(ounit,'(a20,i8)')'buoy', num_obs
523          num_obs = 0
524          do n = 1, iv % num_buoy  
525             if (iv%buoy(n)%loc%proc_domain) then
526                num_obs = num_obs + 1
527                write(ounit,'(i8)') 1
528                write(ounit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))')&
529                   num_obs,1, iv % buoy(n) % info % id, &  ! Station
530                   iv % buoy(n) % info % lat, &       ! Latitude
531                   iv % buoy(n) % info % lon, &       ! Longitude
532                   ob % buoy(n) % p, &                ! Obs Pressure
533                   ob%buoy(n)%u, iv%buoy(n)%u, re%buoy(n)%u, &! O, O-B, O-A u
534                   ob%buoy(n)%v, iv%buoy(n)%v, re%buoy(n)%v, &! O, O-B, O-A v
535                   ob%buoy(n)%t, iv%buoy(n)%t, re%buoy(n)%t, &! O, O-B, O-A t
536                   ob%buoy(n)%p, iv%buoy(n)%p, re%buoy(n)%p, &! O, O-B, O-A p
537                   ob%buoy(n)%q, iv%buoy(n)%q, re%buoy(n)%q   ! O, O-B, O-A q
538             end if
539          end do
540       end if
541    end if
542 
543    ! Transfer TC bogus obs:
544 
545    if (iv % num_bogus > 0) then
546       num_obs = 0
547       do n = 1, iv % num_bogus
548          if (iv%bogus(n)%loc%proc_domain) num_obs = num_obs + 1
549       end do
550       if (num_obs > 0) then
551          write(ounit,'(a20,i8)')'bogus', num_obs
552          num_obs = 0
553          do n = 1, iv % num_bogus
554             if (iv%bogus(n)%loc%proc_domain) then
555                num_obs = num_obs + 1
556                write(ounit,'(i8)') 1
557                write(ounit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))')&
558                    num_obs, 1, iv % bogus(n) % info % id, &  ! Station
559                    iv % bogus(n) % info % lat, &       ! Latitude
560                    iv % bogus(n) % info % lon, &       ! Longitude
561                    missing_r,                  &
562                    ob%bogus(n)%slp, iv%bogus(n)%slp, re%bogus(n)%slp    ! O, O-B, O-A p
563                write(ounit,'(i8)')iv % bogus(n) % info % levels
564                do k = 1, iv % bogus(n) % info % levels
565                   write(ounit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))')&
566                      num_obs , k, iv % bogus(n) % info % id, &  ! Station
567                      iv % bogus(n) % info % lat, &       ! Latitude
568                      iv % bogus(n) % info % lon, &       ! Longitude
569                      iv % bogus(n) % p(k), &             ! Obs Pressure
570                      ob%bogus(n)%u(k), iv%bogus(n)%u(k), re%bogus(n)%u(k), &! O,O-B,O-A u
571                      ob%bogus(n)%v(k), iv%bogus(n)%v(k), re%bogus(n)%v(k), &! O,O-B,O-A v
572                      ob%bogus(n)%t(k), iv%bogus(n)%t(k), re%bogus(n)%t(k), &! O,O-B,O-A t
573                      ob%bogus(n)%q(k), iv%bogus(n)%q(k), re%bogus(n)%q(k)   ! O,O-B,O-A q 
574                end do
575             end if
576          end do
577       end if
578    end if
579 
580    ! Transfer AIRS retrievals:
581 
582    if (iv % num_airsr > 0) then
583       num_obs = 0
584       do n = 1, iv % num_airsr
585          if (iv%airsr(n)%loc%proc_domain) num_obs = num_obs + 1
586       end do
587       if (num_obs > 0) then
588          write(ounit,'(a20,i8)')'airsr', num_obs    
589          num_obs = 0
590          do n = 1, iv % num_airsr
591             if (iv%airsr(n)%loc%proc_domain) then
592                num_obs = num_obs + 1
593                write(ounit,'(i8)')iv % airsr(n) % info % levels
594                do k = 1, iv % airsr(n) % info % levels
595                   write(ounit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))')&
596                      num_obs,k, iv % airsr(n) % info % id, &  ! Station
597                      iv % airsr(n) % info % lat, &       ! Latitude
598                      iv % airsr(n) % info % lon, &       ! Longitude
599                      iv % airsr(n) % p(k), &             ! Obs Pressure
600                      ob%airsr(n)%t(k), iv%airsr(n)%t(k), re%airsr(n)%t(k), &! O,O-B,O-A t
601                      ob%airsr(n)%q(k), iv%airsr(n)%q(k), re%airsr(n)%q(k)   ! O,O-B,O-A q
602                end do
603             end if
604          end do
605       end if
606    end if
607 
608    ! Transfer gpsref obs:
609 
610    if (iv % num_gpsref > 0) then
611       num_obs = 0
612       do n = 1, iv % num_gpsref
613          if (iv%gpsref(n)%loc%proc_domain) num_obs = num_obs + 1
614       end do
615       if (num_obs > 0) then
616          write(ounit,'(a20,i8)')'gpsref', num_obs
617          num_obs = 0
618          do n = 1, iv % num_gpsref
619             if (iv%gpsref(n)%loc%proc_domain) then
620                num_obs = num_obs + 1
621                write(ounit,'(i8)')iv % gpsref(n) % info % levels
622                do k = 1, iv % gpsref(n) % info % levels
623                   write(ounit,'(2i8,a5,2f9.2,f17.7,5(2f17.7,i8,2f17.7))')&
624                      num_obs,k, iv % gpsref(n) % info % id, &  ! Station
625                      iv % gpsref(n) % info % lat, &       ! Latitude
626                      iv % gpsref(n) % info % lon, &       ! Longitude
627                      iv % gpsref(n) % h(k), &             ! Obs Height    
628                      ob%gpsref(n)%ref(k), iv%gpsref(n)%ref(k), re%gpsref(n)%ref(k) ! O, O-B, O-A ref
629                end do
630             end if
631          end do
632       end if
633    end if
634 
635    close (ounit)
636    call da_free_unit(ounit)
637 
638    if (trace_use) call da_trace_exit("da_write_obs")
639 
640 end subroutine da_write_obs
641 
642