Merge remote-tracking branch 'origin/release-v4.6.1'
[WRF.git] / var / da / da_obs_io / da_final_write_obs.inc
blob9c7f1453fae4026028a50a6a423220862149c594
1 subroutine da_final_write_obs(it,iv)
3    !-------------------------------------------------------------------------
4    ! Purpose: Writes full diagnostics for O, (O-B) & OMA together
5    !-------------------------------------------------------------------------   
7    implicit none
8  
9    integer,        intent(in)    :: it
10    type (iv_type), intent(in)    :: iv      ! O-B structure.
11    integer                       :: n, k, iunit
12    integer                       :: ios  ! Error code from MPI routines.
13    integer                       :: num_obs
14    logical                       :: if_wind_sd
15    character(len=filename_len), allocatable     :: filename(:)
16    character(len=filename_len)                  :: file
19    if (trace_use) call da_trace_entry("da_final_write_obs")
21 #ifdef DM_PARALLEL
22    ! Wait to ensure all temporary files have been written
23    call mpi_barrier(comm, ierr)
24 #endif
26    if (rootproc) then
27       call da_get_unit(iunit)
28       allocate (filename(0:num_procs-1))
29       do k = 0,num_procs-1
30          write(unit=filename(k),fmt ='(a,i2.2,a,i4.4)')'gts_omb_oma_',it,'.',k
31       end do 
32       call da_get_unit(omb_unit)
33        write(unit=file,fmt ='(a,i2.2)')'gts_omb_oma_',it
34       open(unit=omb_unit,file=trim(file),form='formatted', status='replace', iostat=ios) 
35       if (ios /= 0) call da_error(__FILE__,__LINE__, &
36          (/"Cannot open file "//file/))
37    end if
39    num_obs = 0
40    if (iv%info(synop)%nlocal > 0) then
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    end if
45    call da_proc_sum_int(num_obs)
46    if_wind_sd = .false.
47    if (wind_sd_synop) if_wind_sd = .true.
48    if (num_obs > 0 .and. rootproc) then
49       write(omb_unit,'(a20,i8)')'synop', num_obs  
50       num_obs = 0
51       do k = 0,num_procs-1
52          call da_read_omb_tmp(filename(k),iunit,num_obs,'synop',5,if_wind_sd)
53       end do
54    end if
56    !------------------------------------------------------------------
57    ! [2] writing Metar
58    !------------------------------------------------------------------
60    num_obs = 0
61    if (iv%info(metar)%nlocal > 0) then
62       do n = 1, iv%info(metar)%nlocal
63          if (iv%info(metar)%proc_domain(1,n)) num_obs = num_obs + 1
64       end do
65    end if
66    call da_proc_sum_int(num_obs)
67    if_wind_sd = .false.
68    if (wind_sd_metar) if_wind_sd = .true.
69    if (num_obs > 0 .and. rootproc) then
70       write(omb_unit,'(a20,20i8)')'metar', num_obs  
71       num_obs = 0
72       do k = 0,num_procs-1
73          call da_read_omb_tmp(filename(k),iunit,num_obs,'metar',5,if_wind_sd)    
74       end do
75    end if
77    !------------------------------------------------------------------
78    ! [3] writing Ships
79    !------------------------------------------------------------------
81    num_obs = 0
82    if (iv%info(ships)%nlocal > 0) then
83       do n = 1, iv%info(ships)%nlocal
84          if(iv%info(ships)%proc_domain(1,n)) num_obs = num_obs + 1
85       end do
86    end if
87    call da_proc_sum_int(num_obs)
88    if_wind_sd = .false.
89    if (wind_sd_ships) if_wind_sd = .true.
90    if (num_obs > 0 .and. rootproc) then
91       write(omb_unit,'(a20,i8)')'ships', num_obs  
92       num_obs = 0
93       do k = 0,num_procs-1
94          call da_read_omb_tmp(filename(k),iunit,num_obs,'ships',5,if_wind_sd)
95       end do
96    end if
98    !------------------------------------------------------------------
99    ! [4] writing GeoAMV
100    !------------------------------------------------------------------
102    num_obs = 0
103    if (iv%info(geoamv)%nlocal > 0) then
104        do n = 1, iv%info(geoamv)%nlocal
105           if (iv%info(geoamv)%proc_domain(1,n)) num_obs = num_obs + 1
106       end do
107    end if
108    call da_proc_sum_int(num_obs)
109    if_wind_sd = .false.
110    if (wind_sd_geoamv) if_wind_sd = .true.
111    if (num_obs > 0 .and. rootproc) then
112       write(omb_unit,'(a20,i8)')'geoamv', num_obs  
113       num_obs = 0
114       do k = 0,num_procs-1
115          call da_read_omb_tmp(filename(k),iunit,num_obs,'geoamv',6,if_wind_sd)
116       end do
117    end if
119    !------------------------------------------------------------------
120    ! [5] writing PolarAMV
121    !------------------------------------------------------------------
123    num_obs = 0
124    if (iv%info(polaramv)%nlocal > 0) then
125       do n = 1, iv%info(polaramv)%nlocal
126          if (iv%info(polaramv)%proc_domain(1,n)) num_obs = num_obs + 1
127       end do
128    end if
129    call da_proc_sum_int(num_obs)
130    if_wind_sd = .false.
131    if (wind_sd_polaramv) if_wind_sd = .true.
132    if (num_obs > 0 .and. rootproc) then
133       write(omb_unit,'(a20,i8)')'polaramv', num_obs  
134       num_obs = 0
135       do k = 0,num_procs-1
136          call da_read_omb_tmp(filename(k),iunit,num_obs,'polaramv',8,if_wind_sd) 
137       end do
138    end if
140    !------------------------------------------------------------------
141    ! [5] writing GPSPW  
142    !------------------------------------------------------------------
144    num_obs = 0
145    if (iv%info(gpspw)%nlocal > 0) then
146       do n = 1, iv%info(gpspw)%nlocal
147          if(iv%info(gpspw)%proc_domain(1,n)) num_obs = num_obs + 1
148       end do
149    end if
150    call da_proc_sum_int(num_obs)
151    if_wind_sd = .false.
152    if (num_obs > 0 .and. rootproc) then
153       write(omb_unit,'(a20,i8)')'gpspw', num_obs  
154       num_obs = 0
155       do k = 0,num_procs-1
156          call da_read_omb_tmp(filename(k),iunit,num_obs,'gpspw',5,if_wind_sd)     
157       end do
158    end if
160    !------------------------------------------------------------------
161    ! [6] writing Sonde  
162    !------------------------------------------------------------------
164    num_obs = 0
165    if (iv%info(sound)%nlocal > 0) then
166       do n = 1, iv%info(sound)%nlocal
167          if (iv%info(sound)%proc_domain(1,n)) num_obs = num_obs + 1
168       end do
169    end if
170    call da_proc_sum_int(num_obs)
171    if_wind_sd = .false.
172    if (wind_sd_sound) if_wind_sd = .true.
173    if (num_obs > 0 .and. rootproc) then
174       write(omb_unit,'(a20,i8)')'sound', num_obs  
175       num_obs = 0
176       do k = 0,num_procs-1
177          call da_read_omb_tmp(filename(k),iunit,num_obs,'sound',5,if_wind_sd) 
178       end do
179    end if
181 !  Now sonde_sfc
182    num_obs = 0
183    if (iv%info(sonde_sfc)%nlocal > 0) then
184       do n = 1, iv%info(sonde_sfc)%nlocal
185          if(iv%info(sonde_sfc)%proc_domain(1,n)) num_obs = num_obs + 1
186       end do
187    end if
188    call da_proc_sum_int(num_obs)
189    if (num_obs > 0 .and. rootproc) then
190       write(omb_unit,'(a20,i8)')'sonde_sfc', num_obs
191       num_obs = 0
192       do k = 0,num_procs-1
193          call da_read_omb_tmp(filename(k),iunit,num_obs,'sonde_sfc',9,if_wind_sd)
194       end do
195    end if
197    !------------------------------------------------------------------
198    ! [7] writing Airep  
199    !------------------------------------------------------------------
201    num_obs = 0
202    if (iv%info(airep)%nlocal > 0) then
203       do n = 1, iv%info(airep)%nlocal
204          if(iv%info(airep)%proc_domain(1,n)) num_obs = num_obs + 1
205       end do
206    end if
207    call da_proc_sum_int(num_obs)
208    if_wind_sd = .false.
209    if (wind_sd_airep) if_wind_sd = .true.
210    if (num_obs > 0 .and. rootproc) then
211       write(omb_unit,'(a20,i8)')'airep', num_obs  
212       num_obs = 0
213       do k = 0,num_procs-1
214          call da_read_omb_tmp(filename(k),iunit,num_obs,'airep',5,if_wind_sd)
215       end do
216    end if
218    !------------------------------------------------------------------
219    ! [8] writing Pilot  
220    !------------------------------------------------------------------
221    
222    num_obs = 0
223    if (iv%info(pilot)%nlocal > 0) then
224       do n = 1, iv%info(pilot)%nlocal
225          if(iv%info(pilot)%proc_domain(1,n)) num_obs = num_obs + 1
226       end do
227    end if
228    call da_proc_sum_int(num_obs)
229    if_wind_sd = .false.
230    if (wind_sd_pilot) if_wind_sd = .true.
231    if (num_obs > 0 .and. rootproc) then
232       write(omb_unit,'(a20,i8)')'pilot', num_obs  
233       num_obs = 0
234       do k = 0,num_procs-1
235          call da_read_omb_tmp(filename(k),iunit,num_obs,'pilot',5,if_wind_sd)    
236       end do
237    end if
239    !------------------------------------------------------------------
240    ! [9] writing ssmi_rv
241    !------------------------------------------------------------------
243    num_obs = 0
244    if (iv%info(ssmi_rv)%nlocal > 0) then
245       do n = 1, iv%info(ssmi_rv)%nlocal
246          if(iv%info(ssmi_rv)%proc_domain(1,n)) num_obs = num_obs + 1
247       end do
248    end if
249    call da_proc_sum_int(num_obs)
250    if_wind_sd = .false.
251    if (num_obs > 0 .and. rootproc) then
252       write(omb_unit,'(a20,i8)')'ssmir', num_obs  
253       num_obs = 0
254       do k = 0,num_procs-1
255          call da_read_omb_tmp(filename(k),iunit,num_obs,'ssmir',5,if_wind_sd)    
256       end do
257    end if
259    !------------------------------------------------------------------
260    ! [10] writing SSMITB
261    !------------------------------------------------------------------
262    
263    num_obs = 0
264    if (iv%info(ssmi_tb)%nlocal > 0) then
265       do n = 1, iv%info(ssmi_tb)%nlocal
266          if (iv%info(ssmi_tb)%proc_domain(1,n)) num_obs = num_obs + 1
267       end do
268    end if
269    call da_proc_sum_int(num_obs)
270    if_wind_sd = .false.
271    if (num_obs > 0 .and. rootproc) then
272       write(omb_unit,'(a20,i8)')'ssmiT', num_obs  
273       num_obs = 0
274       do k = 0,num_procs-1
275          call da_read_omb_tmp(filename(k),iunit,num_obs,'ssmiT',5,if_wind_sd)    
276       end do
277    end if
279    !------------------------------------------------------------------
280    ! [11] writing SATEM  
281    !------------------------------------------------------------------
283    num_obs = 0
284    if (iv%info(satem)%nlocal > 0) then
285       do n = 1, iv%info(satem)%nlocal
286          if(iv%info(satem)%proc_domain(1,n)) num_obs = num_obs + 1
287       end do
288    end if
289    call da_proc_sum_int(num_obs)
290    if_wind_sd = .false.
291    if (num_obs > 0 .and. rootproc) then
292       write(omb_unit,'(a20,i8)')'satem', num_obs  
293       num_obs = 0
294       do k = 0,num_procs-1
295          call da_read_omb_tmp(filename(k),iunit,num_obs,'satem',5,if_wind_sd)    
296       end do
297    end if
299    !------------------------------------------------------------------
300    ! [12] writing SSMT1  
301    !------------------------------------------------------------------
303    num_obs = 0
304    if (iv%info(ssmt1)%nlocal > 0) then
305       do n = 1, iv%info(ssmt1)%nlocal
306          if(iv%info(ssmt1)%proc_domain(1,n)) num_obs = num_obs + 1
307       end do
308    end if
309    call da_proc_sum_int(num_obs)
310    if_wind_sd = .false.
311    if (num_obs > 0 .and. rootproc) then
312       write(omb_unit,'(a20,i8)')'ssmt1', num_obs  
313       num_obs = 0
314       do k = 0,num_procs-1
315          call da_read_omb_tmp(filename(k),iunit,num_obs,'ssmt1',5,if_wind_sd)    
316       end do
317    end if
319    !------------------------------------------------------------------
320    ! [13] writing SSMT2  
321    !------------------------------------------------------------------
323    num_obs = 0
324    if (iv%info(ssmt2)%nlocal > 0) then
325       do n = 1, iv%info(ssmt2)%nlocal
326          if(iv%info(ssmt2)%proc_domain(1,n)) num_obs = num_obs + 1
327       end do
328    end if
329    call da_proc_sum_int(num_obs)
330    if_wind_sd = .false.
331    if (num_obs > 0 .and. rootproc) then
332       write(omb_unit,'(a20,i8)')'ssmt2', num_obs  
333       num_obs = 0
334       do k = 0,num_procs-1
335          call da_read_omb_tmp(filename(k),iunit,num_obs,'ssmt2',5,if_wind_sd)    
336       end do
337    end if
339    !------------------------------------------------------------------
340    ! [14] writing QSCAT  
341    !------------------------------------------------------------------
342    
343    num_obs = 0
344    if (iv%info(qscat)%nlocal > 0) then
345       do n = 1, iv%info(qscat)%nlocal
346          if(iv%info(qscat)%proc_domain(1,n)) num_obs = num_obs + 1
347       end do
348    end if
349    call da_proc_sum_int(num_obs)
350    if_wind_sd = .false.
351    if (wind_sd_qscat) if_wind_sd = .true.
352    if (num_obs > 0 .and. rootproc) then
353       write(omb_unit,'(a20,i8)')'qscat', num_obs  
354       num_obs = 0
355       do k = 0,num_procs-1
356          call da_read_omb_tmp(filename(k),iunit,num_obs,'qscat',5,if_wind_sd)    
357       end do
358    end if
360    !------------------------------------------------------------------
361    ! [15] writing Profiler
362    !------------------------------------------------------------------
364    num_obs = 0
365    if (iv%info(profiler)%nlocal > 0) then
366       do n = 1, iv%info(profiler)%nlocal
367          if(iv%info(profiler)%proc_domain(1,n)) num_obs = num_obs + 1
368       end do
369    end if
370    call da_proc_sum_int(num_obs)
371    if_wind_sd = .false.
372    if (wind_sd_profiler) if_wind_sd = .true.
373    if (num_obs > 0 .and. rootproc) then
374       write(omb_unit,'(a20,i8)')'profiler', num_obs  
375       num_obs = 0
376       do k = 0,num_procs-1
377          call da_read_omb_tmp(filename(k),iunit,num_obs,'profiler',8,if_wind_sd)    
378       end do
379    end if
381    !------------------------------------------------------------------
382    ! [16] writing Buoy 
383    !------------------------------------------------------------------
384    
385    num_obs = 0
386    if (iv%info(buoy)%nlocal > 0) then
387       do n = 1, iv%info(buoy)%nlocal
388          if(iv%info(buoy)%proc_domain(1,n)) num_obs = num_obs + 1
389       end do
390    end if
391    call da_proc_sum_int(num_obs)
392    if_wind_sd = .false.
393    if (wind_sd_buoy) if_wind_sd = .true.
394    if (num_obs > 0 .and. rootproc) then
395       write(omb_unit,'(a20,i8)')'buoy', num_obs  
396       num_obs = 0
397       do k = 0,num_procs-1
398          call da_read_omb_tmp(filename(k),iunit,num_obs,'buoy',4,if_wind_sd)    
399       end do
400    end if
402    !------------------------------------------------------------------
403    ! [17] writing Bogus 
404    !------------------------------------------------------------------
405    
406    num_obs = 0
407    if (iv%info(bogus)%nlocal > 0) then
408       do n = 1, iv%info(bogus)%nlocal
409          if(iv%info(bogus)%proc_domain(1,n)) num_obs = num_obs + 1
410       end do
411    end if
412    call da_proc_sum_int(num_obs)
413    if_wind_sd = .false.
414    if (num_obs > 0 .and. rootproc) then
415       write(omb_unit,'(a20,i8)')'bogus', num_obs  
416       num_obs = 0
417       do k = 0,num_procs-1
418          call da_read_omb_tmp(filename(k),iunit,num_obs,'bogus',5,if_wind_sd)    
419       end do
420    end if
422    !------------------------------------------------------------------
423    ! [18] writing Tamdar
424    !------------------------------------------------------------------
426    num_obs = 0
427    if (iv%info(tamdar)%nlocal > 0) then
428       do n = 1, iv%info(tamdar)%nlocal
429          if (iv%info(tamdar)%proc_domain(1,n)) num_obs = num_obs + 1
430       end do
431    end if
432    call da_proc_sum_int(num_obs)
433    if_wind_sd = .false.
434    if (wind_sd_tamdar) if_wind_sd = .true.
435    if (num_obs > 0 .and. rootproc) then
436       write(omb_unit,'(a20,i8)')'tamdar', num_obs
437       num_obs = 0
438       do k = 0,num_procs-1
439          call da_read_omb_tmp(filename(k),iunit,num_obs,'tamdar',6,if_wind_sd)
440       end do
442    end if
444 !  Now tamdar_sfc
445    num_obs = 0
446    if (iv%info(tamdar_sfc)%nlocal > 0) then
447       do n = 1, iv%info(tamdar_sfc)%nlocal
448          if(iv%info(tamdar_sfc)%proc_domain(1,n)) num_obs = num_obs + 1
449       end do
450    end if
451    call da_proc_sum_int(num_obs)
452    if (num_obs > 0 .and. rootproc) then
453       write(omb_unit,'(a20,i8)')'tamdar_sfc', num_obs  
454       num_obs = 0
455       do k = 0,num_procs-1
456          call da_read_omb_tmp(filename(k),iunit,num_obs,'tamdar_sfc',10,if_wind_sd)
457       end do
458    end if
460    !------------------------------------------------------------------
461    ! [19] writing AIRS retrievals:
462    !------------------------------------------------------------------
464    num_obs = 0
465    if (iv%info(airsr)%nlocal > 0) then
466       do n = 1, iv%info(airsr)%nlocal
467          if(iv%info(airsr)%proc_domain(1,n)) num_obs = num_obs + 1
468       end do
469    end if
470    call da_proc_sum_int(num_obs)
471    if_wind_sd = .false.
472    if (num_obs > 0 .and. rootproc) then
473       write(omb_unit,'(a20,i8)')'airsr', num_obs  
474       num_obs = 0
475       do k = 0,num_procs-1
476          call da_read_omb_tmp(filename(k),iunit,num_obs,'airsr',5,if_wind_sd)
477       end do
478    end if
480    !------------------------------------------------------------------
481    ! [20] writing GPS refractivity
482    !------------------------------------------------------------------
484    num_obs = 0
485    if (iv%info(gpsref)%nlocal > 0) then
486       do n = 1, iv%info(gpsref)%nlocal
487          if(iv%info(gpsref)%proc_domain(1,n)) num_obs = num_obs + 1
488       end do
489    end if
490    call da_proc_sum_int(num_obs)
491    if_wind_sd = .false.
492    if (num_obs > 0 .and. rootproc) then
493       write(omb_unit,'(a20,i8)')'gpsref', num_obs  
494       num_obs = 0
495       do k = 0,num_procs-1
496          call da_read_omb_tmp(filename(k),iunit,num_obs,'gpsref',6,if_wind_sd)    
497       end do
498    end if
500    !------------------------------------------------------------------
501    ! [20.1] writing GPS Excess Phase
502    !------------------------------------------------------------------
504    num_obs = 0
505    if (iv%info(gpseph)%nlocal > 0) then
506       do n = 1, iv%info(gpseph)%nlocal
507          if(iv%info(gpseph)%proc_domain(1,n)) num_obs = num_obs + 1
508       end do
509    end if
510    call da_proc_sum_int(num_obs)
511    if (num_obs > 0 .and. rootproc) then
512       write(omb_unit,'(a20,i8)')'gpseph', num_obs
513       num_obs = 0
514       do k = 0,num_procs-1
515          call da_read_omb_tmp(filename(k),iunit,num_obs,'gpseph',6,if_wind_sd)
516       end do
517    end if
519    !------------------------------------------------------------------
520    ! [21] writing rain 
521    !------------------------------------------------------------------
522    
523    num_obs = 0
524    if (iv%info(rain)%nlocal > 0) then
525       do n = 1, iv%info(rain)%nlocal
526          if(iv%info(rain)%proc_domain(1,n)) num_obs = num_obs + 1
527       end do
528    end if
529    call da_proc_sum_int(num_obs)
530    if_wind_sd = .false.
531    if (num_obs > 0 .and. rootproc) then
532       write(omb_unit,'(a20,i8)')'rain', num_obs  
533       num_obs = 0
534       do k = 0,num_procs-1
535          call da_read_omb_tmp(filename(k),iunit,num_obs,'rain',4,if_wind_sd)    
536       end do
537    end if
539    !------------------------------------------------------------------
540    ! [22] writing lightning
541    !------------------------------------------------------------------
543    num_obs = 0
544    if (iv%info(lightning)%nlocal > 0) then
545       do n = 1, iv%info(lightning)%nlocal
546          if(iv%info(lightning)%proc_domain(1,n)) num_obs = num_obs + 1
547       end do
548    end if
549    call da_proc_sum_int(num_obs)
550    if_wind_sd = .false.
551    if (num_obs > 0 .and. rootproc) then
552       write(omb_unit,'(a20,i8)')'lightning', num_obs
553       num_obs = 0
554       do k = 0,num_procs-1
555          call da_read_omb_tmp(filename(k),iunit,num_obs,'lightning',5,if_wind_sd)
556       end do
557    end if
560    if (rootproc) then
561       close(iunit)
562       close(omb_unit)
563       call da_free_unit(iunit)
564       call da_free_unit(omb_unit)
565       deallocate (filename)
566    end if
568    if (trace_use) call da_trace_exit("da_final_write_obs")
569    
570 end subroutine da_final_write_obs