da_oi_stats_rad.inc

References to this file elsewhere.
1 subroutine da_oi_stats_rad ( stats_unit, iv )
2 
3    !---------------------------------------------------------------------------
4    ! Purpose: Calculate statistics of obs minus background for radiance data.
5    !
6    ! METHOD:  average, rms, minimum, maximum of iv
7    !---------------------------------------------------------------------------
8 
9    implicit none
10 
11    integer,        intent (in)      :: stats_unit    ! Output unit for stats.
12    type (ob_type), intent (in)      :: iv            ! Innovation
13 
14    type (stats_rad_type),    pointer  :: rad(:)
15    integer                          :: n, k, i
16 
17    if ( iv%num_inst < 1 ) RETURN
18 
19    if (trace_use) call da_trace_entry("da_oi_stats_rad")
20 
21    Allocate ( rad(1:iv%num_inst) )
22 
23    do i = 1, iv%num_inst                          !! loop for sensors
24 
25       Allocate ( rad(i)%ichan(1:iv%instid(i)%nchan) )
26       rad(i)%ichan(1:iv%instid(i)%nchan)%num  = 0
27       rad(i)%ichan(1:iv%instid(i)%nchan)%ave  = 0.0
28       rad(i)%ichan(1:iv%instid(i)%nchan)%rms  = 0.0
29       rad(i)%ichan(1:iv%instid(i)%nchan)%minimum%value  = -missing_r
30       rad(i)%ichan(1:iv%instid(i)%nchan)%maximum%value  =  missing_r
31       rad(i)%ichan(1:iv%instid(i)%nchan)%minimum%n      = 1
32       rad(i)%ichan(1:iv%instid(i)%nchan)%maximum%n      = 1
33       do k=1,iv%instid(i)%nchan
34          rad(i)%ichan(k)%minimum%l      = k
35          rad(i)%ichan(k)%maximum%l      = k
36       end do
37 
38       if (iv%instid(i)%num_rad < 1) cycle
39       do k=1, iv%instid(i)%nchan               !! loop for channels
40          do n=1, iv%instid(i)%num_rad              !! loop for pixels
41             if (iv%instid(i)%proc_domain(n)) then
42                call da_stats_calculate( n,k,iv%instid(i)%tb_qc(k,n), &
43                                   iv%instid(i)%tb_inv(k,n), rad(i)%ichan(k)%num, &
44                                   rad(i)%ichan(k)%minimum, rad(i)%ichan(k)%maximum, &
45                                   rad(i)%ichan(k)%ave, rad(i)%ichan(k)%rms)
46 
47             end if          ! end if( oi%sound(n)%loc%proc_domain )
48          end do                                 !! end loop for pixels
49       end do                        !  end loop for channels
50    end do                         !  end loop for sensor
51 
52    do i = 1, iv%num_inst                          !! loop for sensors
53       do k=1, iv%instid(i)%nchan               !! loop for channels
54          ! Do inter-processor communication to gather statistics.
55          call da_proc_sum_int ( rad(i)%ichan(k)%num )
56          call da_proc_stats_combine(rad(i)%ichan(k)%ave, rad(i)%ichan(k)%rms, &
57                            rad(i)%ichan(k)%minimum%value, rad(i)%ichan(k)%maximum%value, &
58                            rad(i)%ichan(k)%minimum%n, rad(i)%ichan(k)%maximum%n, &
59                            rad(i)%ichan(k)%minimum%l, rad(i)%ichan(k)%maximum%l )
60 
61       end do                        !  end loop for channels
62 
63       if (rootproc) then
64          if ( any( rad(i)%ichan(:)%num /= 0 ) ) then
65             write(unit=stats_unit, fmt='(/a/)') ' Diagnostics of OI for radiance         '//iv%instid(i)%rttovid_string 
66             call da_print_stats_rad( stats_unit, iv%instid(i)%nchan, rad(i) )
67          end if
68       end if
69    end do                         !  end loop for sensor
70 
71    do i = 1, iv%num_inst           ! loop for sensors
72       Deallocate ( rad(i)%ichan )
73    end do
74 
75    Deallocate ( rad )
76 
77    if (trace_use) call da_trace_exit("da_oi_stats_rad")
78 
79 end subroutine da_oi_stats_rad
80 
81