Merge remote-tracking branch 'origin/release-v4.6.1'
[WRF.git] / var / da / da_obs_io / da_write_noise_to_ob.inc
blobcc800e35c4de0f92efc78c6d618406d1f33ec4fa
1 subroutine da_write_noise_to_ob (iv) 
3    !-------------------------------------------------------------------------
4    ! Purpose: Write consolidated obs noise created in da_add_noise_to_ob
5    !-------------------------------------------------------------------------
7    implicit none
9    type (iv_type), intent(in)    :: iv   ! Obs and header structure.
10    integer                       :: n, k,kk,i, iunit
11    integer                       :: num_obs
12    character(len=filename_len), allocatable     :: filename(:)   
13    character*20                  :: ob_name   
15    call da_trace_entry("da_write_noise_to_ob")
17 #ifdef DM_PARALLEL
18    ! Ensure other processors have written their temporary files
19    call mpi_barrier(comm, ierr)
20 #endif
21    call da_get_unit (iunit)
22    allocate (filename(0:num_procs-1))
23    do k = 0,num_procs-1
24       write(unit=filename(k), fmt='(a,i4.4)')'rand_obs_error.',k
25    end do
26    if (rootproc) then
27       call da_get_unit (rand_unit)
28       open(unit=rand_unit,file='rand_obs_error',form='formatted', &
29          iostat=ierr,status='unknown')
30 !         iostat=ierr,status='new')
31       if (ierr /= 0) &
32          call da_error(__FILE__,__LINE__, (/"Cannot open file rand_obs_error"/))
33    end if
35    num_obs = 0
36    if (iv%info(synop)%nlocal > 0) then
37       do n = 1, iv%info(synop)%nlocal
38          if (iv%info(synop)%proc_domain(1,n)) num_obs = num_obs + 1
39       end do
40    end if
41    call da_proc_sum_int(num_obs)
42    if (num_obs > 0 .and. rootproc) then
43       write(rand_unit,'(a20,i8)')'synop', num_obs  
44       num_obs = 0
45       do k = 0,num_procs-1
46          call da_read_rand_unit(filename(k),iunit,num_obs,'synop',5)
47       end do
48    end if
50    !------------------------------------------------------------------
51    ! [2] writing Metar
52    !------------------------------------------------------------------
54    num_obs = 0
55    if (iv%info(metar)%nlocal > 0) then
56       do n = 1, iv%info(metar)%nlocal
57          if (iv%info(metar)%proc_domain(1,n)) num_obs = num_obs + 1
58       end do
59    end if
60    call da_proc_sum_int(num_obs)
61    if (num_obs > 0 .and. rootproc) then
62       write(rand_unit,'(a20,20i8)')'metar', num_obs  
63       num_obs = 0
64       do k = 0,num_procs-1
65          call da_read_rand_unit(filename(k),iunit,num_obs,'metar',5)    
66       end do
67    end if
69    !------------------------------------------------------------------
70    ! [3] writing Ships
71    !------------------------------------------------------------------
73    num_obs = 0
74    if (iv%info(ships)%nlocal > 0) then
75       do n = 1, iv%info(ships)%nlocal
76          if (iv%info(ships)%proc_domain(1,n)) num_obs = num_obs + 1
77       end do
78    end if
79    call da_proc_sum_int(num_obs)
80    if (num_obs > 0 .and. rootproc) then
81       write(rand_unit,'(a20,i8)')'ships', num_obs  
82       num_obs = 0
83       do k = 0,num_procs-1
84          call da_read_rand_unit(filename(k),iunit,num_obs,'ships',5)
85       end do
86    end if
88    !------------------------------------------------------------------
89    ! [4] writing GeoAMV
90    !------------------------------------------------------------------
91    
92    num_obs = 0
93    if (iv%info(geoamv)%nlocal > 0) then
94       do n = 1, iv%info(geoamv)%nlocal
95          if (iv%info(geoamv)%proc_domain(1,n)) num_obs = num_obs + 1
96       end do
97    end if
98    call da_proc_sum_int(num_obs)
99    if (num_obs > 0 .and. rootproc) then
100       write(rand_unit,'(a20,i8)')'geoamv', num_obs  
101       num_obs = 0
102       do k = 0,num_procs-1
103          call da_read_rand_unit(filename(k),iunit,num_obs,'geoamv',6)
104       end do
105    end if
107    !------------------------------------------------------------------
108    ! [5] writing PolarAMV
109    !------------------------------------------------------------------
111    num_obs = 0
112    if (iv%info(polaramv)%nlocal > 0) then
113       do n = 1, iv%info(polaramv)%nlocal
114          if (iv%info(polaramv)%proc_domain(1,n)) num_obs = num_obs + 1
115       end do
116    end if
117    call da_proc_sum_int(num_obs)
118    if (num_obs > 0 .and. rootproc) then
119       if (rootproc) write(rand_unit,'(a20,i8)')'polaramv', num_obs  
120       num_obs = 0
121       do k = 0,num_procs-1
122          call da_read_rand_unit(filename(k),iunit,num_obs,'polaramv',8) 
123       end do
124    end if
126    !------------------------------------------------------------------
127    ! [5] writing GPSPW  
128    !------------------------------------------------------------------
129    
130    num_obs = 0
131    if (iv%info(gpspw)%nlocal > 0) then
132       do n = 1, iv%info(gpspw)%nlocal
133          if (iv%info(gpspw)%proc_domain(1,n)) num_obs = num_obs + 1
134       end do
135    end if
136    call da_proc_sum_int(num_obs)
137    if (num_obs > 0 .and. rootproc) then
138       write(rand_unit,'(a20,i8)')'gpspw', num_obs  
139       num_obs = 0
140       do k = 0,num_procs-1
141          call da_read_rand_unit(filename(k),iunit,num_obs,'gpspw',5)     
142       end do
143    end if
145    !------------------------------------------------------------------
146    ! [6] writing Sonde  
147    !------------------------------------------------------------------
148    
149    num_obs = 0
150    if (iv%info(sound)%nlocal > 0) then
151       do n = 1, iv%info(sound)%nlocal
152          if (iv%info(sound)%proc_domain(1,n)) num_obs = num_obs + 1
153       end do
154    end if
155    call da_proc_sum_int(num_obs)
156    if (num_obs > 0 .and. rootproc) then
157       write(rand_unit,'(a20,i8)')'sound', num_obs  
158       num_obs = 0
159       do k = 0,num_procs-1
160          call da_read_rand_unit(filename(k),iunit,num_obs,'sound',5)     
161       end do
162       !  writing Sonde_sfc  
163       if (rootproc) write(rand_unit,'(a20,i8)')'sonde_sfc', num_obs  
164       num_obs = 0
165       do k = 0,num_procs-1
166          call da_read_rand_unit(filename(k),iunit,num_obs,'sonde_sfc',9)
167       end do
168    end if
170    !------------------------------------------------------------------
171    ! [7] writing Airep  
172    !------------------------------------------------------------------
173    
174    num_obs = 0
175    if (iv%info(airep)%nlocal > 0) then
176       do n = 1, iv%info(airep)%nlocal
177          if (iv%info(airep)%proc_domain(1,n)) num_obs = num_obs + 1
178       end do
179    end if
180    call da_proc_sum_int(num_obs)
181    if (num_obs > 0 .and. rootproc) then
182       write(rand_unit,'(a20,i8)')'airep', num_obs  
183       num_obs = 0
184       do k = 0,num_procs-1
185          call da_read_rand_unit(filename(k),iunit,num_obs,'airep',5)   
186       end do
187    end if
189    !------------------------------------------------------------------
190    ! [8] writing   
191    !------------------------------------------------------------------
192    
193    num_obs = 0
194    if (iv%info(pilot)%nlocal > 0) then
195       do n = 1, iv%info(pilot)%nlocal
196          if (iv%info(pilot)%proc_domain(1,n)) num_obs = num_obs + 1
197       end do
198    end if
199    call da_proc_sum_int(num_obs)
200    if (num_obs > 0 .and. rootproc) then
201       write(rand_unit,'(a20,i8)')'pilot', num_obs  
202       num_obs = 0
203       do k = 0,num_procs-1
204          call da_read_rand_unit(filename(k),iunit,num_obs,'pilot',5)    
205       end do
206    end if
208    !------------------------------------------------------------------
209    ! [9] writing ssmi_rv
210    !------------------------------------------------------------------
211    
212    num_obs = 0
213    if (iv%info(ssmi_rv)%nlocal > 0) then
214       do n = 1, iv%info(ssmi_rv)%nlocal
215          if (iv%info(ssmi_rv)%proc_domain(1,n)) num_obs = num_obs + 1
216       end do
217    end if
218    call da_proc_sum_int(num_obs)
219    if (num_obs > 0 .and. rootproc) then
220       write(rand_unit,'(a20,i8)')'ssmir', num_obs  
221       num_obs = 0
222       do k = 0,num_procs-1
223          call da_read_rand_unit(filename(k),iunit,num_obs,'ssmir',5)    
224      end do
225    end if
227    !------------------------------------------------------------------
228    ! [10] writing SSMITB
229    !------------------------------------------------------------------
231    num_obs = 0
232    if (iv%info(ssmi_tb)%nlocal > 0) then
233       do n = 1, iv%info(ssmi_tb)%nlocal
234          if (iv%info(ssmi_tb)%proc_domain(1,n)) num_obs = num_obs + 1
235       end do
236    end if
237    call da_proc_sum_int(num_obs)
238    if (num_obs > 0 .and. rootproc) then
239       write(rand_unit,'(a20,i8)')'ssmiT', num_obs  
240       num_obs = 0
241       do k = 0,num_procs-1
242          call da_read_rand_unit(filename(k),iunit,num_obs,'ssmiT',5)    
243       end do
244    end if
246    !------------------------------------------------------------------
247    ! [11] writing SATEM  
248    !------------------------------------------------------------------
249    
250    num_obs = 0
251    if (iv%info(satem)%nlocal > 0) then
252       do n = 1, iv%info(satem)%nlocal
253          if (iv%info(satem)%proc_domain(1,n)) num_obs = num_obs + 1
254       end do
255    end if
256    call da_proc_sum_int(num_obs)
257    if (num_obs > 0 .and. rootproc) then
258       write(rand_unit,'(a20,i8)')'satem', num_obs  
259       num_obs = 0
260       do k = 0,num_procs-1
261          call da_read_rand_unit(filename(k),iunit,num_obs,'satem',5)    
262       end do
263    end if
265    !------------------------------------------------------------------
266    ! [12] writing SSMT1  
267    !------------------------------------------------------------------
269    num_obs = 0
270    if (iv%info(ssmt1)%nlocal > 0) then
271       do n = 1, iv%info(ssmt1)%nlocal
272          if (iv%info(ssmt1)%proc_domain(1,n)) num_obs = num_obs + 1
273       end do
274    end if
275    call da_proc_sum_int(num_obs)
276    if (num_obs > 0 .and. rootproc) then
277       write(rand_unit,'(a20,i8)')'ssmt1', num_obs  
278       num_obs = 0
279       do k = 0,num_procs-1
280          call da_read_rand_unit(filename(k),iunit,num_obs,'ssmt1',5)    
281       end do
282    end if
284    !------------------------------------------------------------------
285    ! [13] writing SSMT2  
286    !------------------------------------------------------------------
287    
288    num_obs = 0
289    if (iv%info(ssmt2)%nlocal > 0) then
290       do n = 1, iv%info(ssmt2)%nlocal
291          if (iv%info(ssmt2)%proc_domain(1,n)) num_obs = num_obs + 1
292       end do
293    end if
294    call da_proc_sum_int(num_obs)
295    if (num_obs > 0 .and. rootproc) then
296       write(rand_unit,'(a20,i8)')'ssmt2', num_obs  
297       num_obs = 0
298       do k = 0,num_procs-1
299          call da_read_rand_unit(filename(k),iunit,num_obs,'ssmt2',5)    
300       end do
301    end if
303    !------------------------------------------------------------------
304    ! [14] writing QSCAT  
305    !------------------------------------------------------------------
306     
307    num_obs = 0
308    if (iv%info(qscat)%nlocal > 0) then
309       do n = 1, iv%info(qscat)%nlocal
310          if (iv%info(qscat)%proc_domain(1,n)) num_obs = num_obs + 1
311       end do
312    end if
313    call da_proc_sum_int(num_obs)
314    if (num_obs > 0 .and. rootproc) then
315       write(rand_unit,'(a20,i8)')'qscat', num_obs  
316       num_obs = 0
317       do k = 0,num_procs-1
318          call da_read_rand_unit(filename(k),iunit,num_obs,'qscat',5)    
319       end do
320    end if
322    !------------------------------------------------------------------
323    ! [15] writing Profiler
324    !------------------------------------------------------------------
325    
326    num_obs = 0
327    if (iv%info(profiler)%nlocal > 0) then
328       do n = 1, iv%info(profiler)%nlocal
329          if (iv%info(profiler)%proc_domain(1,n)) num_obs = num_obs + 1
330       end do
331    end if
332    call da_proc_sum_int(num_obs)
333    if (num_obs > 0 .and. rootproc) then
334       write(rand_unit,'(a20,i8)')'profiler', num_obs  
335       num_obs = 0
336       do k = 0,num_procs-1
337          call da_read_rand_unit(filename(k),iunit,num_obs,'profiler',8)    
338       end do
339    end if
341    !------------------------------------------------------------------
342    ! [16] writing Buoy 
343    !------------------------------------------------------------------
344    
345    num_obs = 0
346    if (iv%info(buoy)%nlocal > 0) then
347       do n = 1, iv%info(buoy)%nlocal
348          if (iv%info(buoy)%proc_domain(1,n)) num_obs = num_obs + 1
349       end do
350    end if
351    call da_proc_sum_int(num_obs)
352    if (num_obs > 0 .and. rootproc) then
353       write(rand_unit,'(a20,i8)')'buoy', num_obs  
354       num_obs = 0
355       do k = 0,num_procs-1
356          call da_read_rand_unit(filename(k),iunit,num_obs,'buoy',4)    
357       end do
358    end if
360    !------------------------------------------------------------------
361    ! [17] writing  Bogus 
362    !------------------------------------------------------------------
363   
364    num_obs = 0
365    if (iv%info(bogus)%nlocal > 0) then
366       do n = 1, iv%info(bogus)%nlocal
367          if (iv%info(bogus)%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 (num_obs > 0 .and. rootproc) then
372       write(rand_unit,'(a20,i8)')'bogus', num_obs  
373       num_obs = 0
374       do k = 0,num_procs-1
375          call da_read_rand_unit(filename(k),iunit,num_obs,'bogus',5)    
376       end do
377    end if
379    !------------------------------------------------------------------
380    ! [18] writing AIRS retrievals
381    !------------------------------------------------------------------
382    
383    num_obs = 0
384    if (iv%info(airsr)%nlocal > 0) then
385       do n = 1, iv%info(airsr)%nlocal
386          if (iv%info(airsr)%proc_domain(1,n)) num_obs = num_obs + 1
387       end do
388    end if
389    call da_proc_sum_int(num_obs)
390    if (num_obs > 0 .and. rootproc) then
391        write(rand_unit,'(a20,i8)')'airsr', num_obs  
392        num_obs = 0
393        do k = 0,num_procs-1
394           call da_read_rand_unit(filename(k),iunit,num_obs,'airsr',5)    
395        end do
396     end if
398    !------------------------------------------------------------------
399    ! [19] writing gpsref
400    !------------------------------------------------------------------
402    num_obs = 0
403    if (iv%info(gpsref)%nlocal > 0) then
404       do n = 1, iv%info(gpsref)%nlocal
405          if (iv%info(gpsref)%proc_domain(1,n)) num_obs = num_obs + 1
406       end do
407    end if
408    call da_proc_sum_int(num_obs)
409    if (num_obs > 0 .and. rootproc) then
410       write(rand_unit,'(a20,20i8)')'gpsref', num_obs  
411       num_obs = 0
412       do k = 0,num_procs-1
413          call da_read_rand_unit(filename(k),iunit,num_obs,'gpsref',6)    
414       end do
415    end if
417    !------------------------------------------------------------------
418    ! writing gpseph
419    !------------------------------------------------------------------
421    num_obs = 0
422    if (iv%info(gpseph)%nlocal > 0) then
423       do n = 1, iv%info(gpseph)%nlocal
424          if (iv%info(gpseph)%proc_domain(1,n)) num_obs = num_obs + 1
425       end do
426    end if
427    call da_proc_sum_int(num_obs)
428    if (num_obs > 0 .and. rootproc) then
429       write(rand_unit,'(a20,20i8)')'gpseph', num_obs
430       num_obs = 0
431       do k = 0,num_procs-1
432          call da_read_rand_unit(filename(k),iunit,num_obs,'gpseph',6)
433       end do
434    end if
436    !------------------------------------------------------------------
437    ! [20] writing Radiance data:  
438    !------------------------------------------------------------------
440    if (iv%num_inst > 0) then
441       do i = 1, iv%num_inst                 ! loop for sensor
442          !if (iv%instid(i)%num_rad < 1) cycle ! may broken da_proc_sum_int(num_obs) if some PE num_rad=0
443          do k = 1,iv%instid(i)%nchan        ! loop for channel
444             ! Counting number of obs for channle k
445             num_obs = 0
446             do n = 1,iv%instid(i)%num_rad      ! loop for pixel
447                if (iv%instid(i)%info%proc_domain(1,n) .and. &
448                   (iv%instid(i)%tb_qc(k,n) >= obs_qc_pointer)) then
449                   num_obs = num_obs + 1
450                end if
451             end do                                ! end loop for pixel
452             call da_proc_sum_int(num_obs)
453             if (rootproc .and. num_obs > 0) then
454                write (ob_name,'(a,a,i4.4)') &
455                   trim(iv%instid(i)%rttovid_string),'-',iv%instid(i)%ichan(k)
456                write (rand_unit,'(a20,i8)')  ob_name,num_obs
457                num_obs = 0
458                do kk = 0,num_procs-1
459                   call da_read_rand_unit(filename(kk),iunit,num_obs, &
460                      trim(ob_name),len(trim(ob_name)))
461                end do
462             end if
463          end do                           ! end loop for channel
464       end do                            ! end loop for sensor
465    end if
467    !------------------------------------------------------------------
468    ! [21] writing RAIN 
469    !------------------------------------------------------------------
471    num_obs = 0
472    if (iv%info(rain)%nlocal > 0) then
473       do n = 1, iv%info(rain)%nlocal
474          if (iv%info(rain)%proc_domain(1,n)) num_obs = num_obs + 1
475       end do
476    end if
477    call da_proc_sum_int(num_obs)
478    if (num_obs > 0 .and. rootproc) then
479       write(rand_unit,'(a20,i8)')'rain', num_obs
480       num_obs = 0
481       do k = 0,num_procs-1
482          call da_read_rand_unit(filename(k),iunit,num_obs,'rain',4)
483       end do
484    end if
487    call da_free_unit (iunit)
488 !rizvi   call da_free_unit (rand_unit)
489    if (rootproc ) call da_free_unit (rand_unit)
490    deallocate (filename)
491    call da_trace_exit("da_write_noise_to_ob")
492    
493 end subroutine da_write_noise_to_ob