updated top-level README and version_decl for V4.5 (#1847)
[WRF.git] / var / da / da_radiance / da_ao_stats_rad.inc
blob503c67c88f3df6b862f05ced5bc59e0ee545597e
1 subroutine da_ao_stats_rad ( stats_unit, iv, re )
3    !---------------------------------------------------------------------------
4    ! Purpose: Calculate statistics of obs minus analysis for radiance data.
5    !
6    ! METHOD:  average, rms, minimum, maximum of re
7    !---------------------------------------------------------------------------
9    implicit none
11    integer,        intent (in)    :: stats_unit    ! output unit for stats.
12    type (iv_type), intent (inout) :: iv            ! innovation
13    type (y_type),  intent (in)    :: re            ! o-a
15    type (stats_rad_type), pointer  :: rad(:)
16    integer                         :: n, k, i
18    iv%nstats(radiance) = 0
20    if ( iv%num_inst < 1 ) return
22    if (trace_use) call da_trace_entry("da_ao_stats_rad")
24    allocate ( rad(1:iv%num_inst) )
26    do i = 1, iv%num_inst                          !! loop for sensors
28       allocate (rad(i)%ichan(1:iv%instid(i)%nchan))
29       rad(i)%ichan(1:iv%instid(i)%nchan)%num  = 0
30       rad(i)%ichan(1:iv%instid(i)%nchan)%ave  = 0.0
31       rad(i)%ichan(1:iv%instid(i)%nchan)%rms  = 0.0
32       rad(i)%ichan(1:iv%instid(i)%nchan)%minimum%value  = -missing_r
33       rad(i)%ichan(1:iv%instid(i)%nchan)%maximum%value  =  missing_r
34       rad(i)%ichan(1:iv%instid(i)%nchan)%minimum%n      = 1
35       rad(i)%ichan(1:iv%instid(i)%nchan)%maximum%n      = 1
36       do k=1,iv%instid(i)%nchan
37          rad(i)%ichan(k)%minimum%l      = k
38          rad(i)%ichan(k)%maximum%l      = k
39       end do
41       if (iv%instid(i)%num_rad < 1) cycle
43       do k=1, iv%instid(i)%nchan               !! loop for channels
44          do n=1, iv%instid(i)%num_rad              !! loop for pixels
45             if (iv%instid(i)%info%proc_domain(1,n)) then
46                call da_stats_calculate( n,k,iv%instid(i)%tb_qc(k,n), &
47                                  re%instid(i)%tb(k,n), rad(i)%ichan(k)%num, &
48                                  rad(i)%ichan(k)%minimum, rad(i)%ichan(k)%maximum, &
49                                  rad(i)%ichan(k)%ave, rad(i)%ichan(k)%rms)
51             end if          ! end if( oi%sound(n)%loc%proc_domain )
52          end do                                 !! end loop for pixels
53       end do                        !  end loop for channels
54    end do                         !  end loop for sensor
56    do i = 1, iv%num_inst                          !! loop for sensors
57       do k=1, iv%instid(i)%nchan               !! loop for channels
58          ! FIX? generate 1D array to allow a da_proc_sum_ints call here
59          ! Do inter-processor communication to gather statistics.
60          call da_proc_sum_int ( rad(i)%ichan(k)%num )
61          call da_proc_stats_combine(rad(i)%ichan(k)%ave, rad(i)%ichan(k)%rms, &
62                            rad(i)%ichan(k)%minimum%value, rad(i)%ichan(k)%maximum%value, &
63                            rad(i)%ichan(k)%minimum%n, rad(i)%ichan(k)%maximum%n, &
64                            rad(i)%ichan(k)%minimum%l, rad(i)%ichan(k)%maximum%l )
66          iv%nstats(radiance) = iv%nstats(radiance) + rad(i)%ichan(k)%num
67       end do                        !  end loop for channels
69       if (rootproc) then
70          if (any(rad(i)%ichan(:)%num /= 0)) then
71             write(unit=stats_unit, fmt='(/a/)') ' Diagnostics of AO for radiance         '//iv%instid(i)%rttovid_string
72             call da_print_stats_rad( stats_unit, iv%instid(i)%nchan, rad(i) )
73          end if
74       end if
75    end do                         !  end loop for sensor
77    do i = 1, iv%num_inst           ! loop for sensors
78       deallocate (rad(i)%ichan)
79    end do
81    deallocate (rad)
83    if (trace_use) call da_trace_exit("da_ao_stats_rad")
85 end subroutine da_ao_stats_rad