1 subroutine da_write_noise_to_ob (iv)
3 !-------------------------------------------------------------------------
4 ! Purpose: Write consolidated obs noise created in da_add_noise_to_ob
5 !-------------------------------------------------------------------------
9 type (iv_type), intent(in) :: iv ! Obs and header structure.
10 integer :: n, k,kk,i, iunit
12 character(len=filename_len), allocatable :: filename(:)
13 character*20 :: ob_name
15 call da_trace_entry("da_write_noise_to_ob")
18 ! Ensure other processors have written their temporary files
19 call mpi_barrier(comm, ierr)
21 call da_get_unit (iunit)
22 allocate (filename(0:num_procs-1))
24 write(unit=filename(k), fmt='(a,i4.4)')'rand_obs_error.',k
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')
32 call da_error(__FILE__,__LINE__, (/"Cannot open file rand_obs_error"/))
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
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
46 call da_read_rand_unit(filename(k),iunit,num_obs,'synop',5)
50 !------------------------------------------------------------------
52 !------------------------------------------------------------------
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
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
65 call da_read_rand_unit(filename(k),iunit,num_obs,'metar',5)
69 !------------------------------------------------------------------
71 !------------------------------------------------------------------
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
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
84 call da_read_rand_unit(filename(k),iunit,num_obs,'ships',5)
88 !------------------------------------------------------------------
90 !------------------------------------------------------------------
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
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
103 call da_read_rand_unit(filename(k),iunit,num_obs,'geoamv',6)
107 !------------------------------------------------------------------
108 ! [5] writing PolarAMV
109 !------------------------------------------------------------------
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
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
122 call da_read_rand_unit(filename(k),iunit,num_obs,'polaramv',8)
126 !------------------------------------------------------------------
128 !------------------------------------------------------------------
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
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
141 call da_read_rand_unit(filename(k),iunit,num_obs,'gpspw',5)
145 !------------------------------------------------------------------
147 !------------------------------------------------------------------
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
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
160 call da_read_rand_unit(filename(k),iunit,num_obs,'sound',5)
163 if (rootproc) write(rand_unit,'(a20,i8)')'sonde_sfc', num_obs
166 call da_read_rand_unit(filename(k),iunit,num_obs,'sonde_sfc',9)
170 !------------------------------------------------------------------
172 !------------------------------------------------------------------
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
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
185 call da_read_rand_unit(filename(k),iunit,num_obs,'airep',5)
189 !------------------------------------------------------------------
191 !------------------------------------------------------------------
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
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
204 call da_read_rand_unit(filename(k),iunit,num_obs,'pilot',5)
208 !------------------------------------------------------------------
209 ! [9] writing ssmi_rv
210 !------------------------------------------------------------------
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
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
223 call da_read_rand_unit(filename(k),iunit,num_obs,'ssmir',5)
227 !------------------------------------------------------------------
228 ! [10] writing SSMITB
229 !------------------------------------------------------------------
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
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
242 call da_read_rand_unit(filename(k),iunit,num_obs,'ssmiT',5)
246 !------------------------------------------------------------------
248 !------------------------------------------------------------------
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
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
261 call da_read_rand_unit(filename(k),iunit,num_obs,'satem',5)
265 !------------------------------------------------------------------
267 !------------------------------------------------------------------
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
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
280 call da_read_rand_unit(filename(k),iunit,num_obs,'ssmt1',5)
284 !------------------------------------------------------------------
286 !------------------------------------------------------------------
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
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
299 call da_read_rand_unit(filename(k),iunit,num_obs,'ssmt2',5)
303 !------------------------------------------------------------------
305 !------------------------------------------------------------------
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
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
318 call da_read_rand_unit(filename(k),iunit,num_obs,'qscat',5)
322 !------------------------------------------------------------------
323 ! [15] writing Profiler
324 !------------------------------------------------------------------
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
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
337 call da_read_rand_unit(filename(k),iunit,num_obs,'profiler',8)
341 !------------------------------------------------------------------
343 !------------------------------------------------------------------
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
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
356 call da_read_rand_unit(filename(k),iunit,num_obs,'buoy',4)
360 !------------------------------------------------------------------
362 !------------------------------------------------------------------
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
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
375 call da_read_rand_unit(filename(k),iunit,num_obs,'bogus',5)
379 !------------------------------------------------------------------
380 ! [18] writing AIRS retrievals
381 !------------------------------------------------------------------
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
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
394 call da_read_rand_unit(filename(k),iunit,num_obs,'airsr',5)
398 !------------------------------------------------------------------
399 ! [19] writing gpsref
400 !------------------------------------------------------------------
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
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
413 call da_read_rand_unit(filename(k),iunit,num_obs,'gpsref',6)
417 !------------------------------------------------------------------
419 !------------------------------------------------------------------
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
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
432 call da_read_rand_unit(filename(k),iunit,num_obs,'gpseph',6)
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
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
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
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)))
463 end do ! end loop for channel
464 end do ! end loop for sensor
467 !------------------------------------------------------------------
469 !------------------------------------------------------------------
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
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
482 call da_read_rand_unit(filename(k),iunit,num_obs,'rain',4)
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")
493 end subroutine da_write_noise_to_ob