updated top-level README and version_decl for V4.5 (#1847)
[WRF.git] / var / da / da_radiance / da_oi_stats_rad.inc
blob1cd0e2fe96334cf18990b603d3d08c2807f6d1be
1 subroutine da_oi_stats_rad ( stats_unit, iv )
3    !---------------------------------------------------------------------------
4    ! Purpose: Calculate statistics of obs minus background for radiance data.
5    !
6    ! METHOD:  average, rms, minimum, maximum of iv
7    !---------------------------------------------------------------------------
9    implicit none
11    integer,        intent (in)      :: stats_unit    ! Output unit for stats.
12    type (iv_type), intent (in)      :: iv            ! Innovation
14    type (stats_rad_type),    pointer  :: rad(:)
15    integer                          :: n, k, i
17    if (trace_use) call da_trace_entry("da_oi_stats_rad")
19    allocate ( rad(1:iv%num_inst) )
21    do i = 1, iv%num_inst                          !! loop for sensors
22       allocate ( rad(i)%ichan(1:iv%instid(i)%nchan) )
23       rad(i)%ichan(1:iv%instid(i)%nchan)%num  = 0
24       rad(i)%ichan(1:iv%instid(i)%nchan)%ave  = 0.0
25       rad(i)%ichan(1:iv%instid(i)%nchan)%rms  = 0.0
26       rad(i)%ichan(1:iv%instid(i)%nchan)%minimum%value  = -missing_r
27       rad(i)%ichan(1:iv%instid(i)%nchan)%maximum%value  =  missing_r
28       rad(i)%ichan(1:iv%instid(i)%nchan)%minimum%n      = 1
29       rad(i)%ichan(1:iv%instid(i)%nchan)%maximum%n      = 1
30       do k=1,iv%instid(i)%nchan
31          rad(i)%ichan(k)%minimum%l      = k
32          rad(i)%ichan(k)%maximum%l      = k
33       end do
35       if (iv%instid(i)%num_rad < 1) cycle
36       do k=1, iv%instid(i)%nchan               !! loop for channels
37          do n=1, iv%instid(i)%num_rad              !! loop for pixels
38             if (iv%instid(i)%info%proc_domain(1,n)) then
39                call da_stats_calculate( n,k,iv%instid(i)%tb_qc(k,n), &
40                   iv%instid(i)%tb_inv(k,n), rad(i)%ichan(k)%num, &
41                   rad(i)%ichan(k)%minimum, rad(i)%ichan(k)%maximum, &
42                   rad(i)%ichan(k)%ave, rad(i)%ichan(k)%rms)
44             end if          ! end if( iv%sound(n)%loc%proc_domain )
45          end do                                 !! end loop for pixels
46       end do                        !  end loop for channels
47    end do                         !  end loop for sensor
49    do i = 1, iv%num_inst                          !! loop for sensors
50       do k=1, iv%instid(i)%nchan               !! loop for channels
51          ! Do inter-processor communication to gather statistics.
52          call da_proc_sum_int ( rad(i)%ichan(k)%num )
53          call da_proc_stats_combine(rad(i)%ichan(k)%ave, rad(i)%ichan(k)%rms, &
54             rad(i)%ichan(k)%minimum%value, rad(i)%ichan(k)%maximum%value, &
55             rad(i)%ichan(k)%minimum%n, rad(i)%ichan(k)%maximum%n, &
56             rad(i)%ichan(k)%minimum%l, rad(i)%ichan(k)%maximum%l )
58       end do                        !  end loop for channels
60       if (rootproc) then
61          if ( any( rad(i)%ichan(:)%num /= 0 ) ) then
62             write(unit=stats_unit, fmt='(/a/)') &
63                ' Diagnostics of OI for radiance         '//iv%instid(i)%rttovid_string 
64             call da_print_stats_rad( stats_unit, iv%instid(i)%nchan, rad(i) )
65          end if
66       end if
67    end do                         !  end loop for sensor
69    do i = 1, iv%num_inst           ! loop for sensors
70       deallocate (rad(i)%ichan)
71    end do
73    deallocate (rad)
75    if (trace_use) call da_trace_exit("da_oi_stats_rad")
77 end subroutine da_oi_stats_rad