da_oi_stats_airsr.inc
References to this file elsewhere.
1 subroutine da_oi_stats_airsr (stats_unit, oi)
2
3 !----------------------------------------------------------------------
4 ! Purpose: Prints out innov diagnostics for AIRS retrievals
5 !----------------------------------------------------------------------
6
7 implicit none
8
9 integer, intent (in) :: stats_unit ! Output unit for stats.
10 type (ob_type), intent (in) :: oi ! OI
11
12 type (stats_airsr_type) :: airsr
13 integer :: nt, nq
14 integer :: n, k
15
16 nt = 0
17 nq = 0
18
19 airsr%maximum%t = maxmin_type(missing_r, 0, 0)
20 airsr%maximum%q = maxmin_type(missing_r, 0, 0)
21 airsr%minimum%t = maxmin_type(-missing_r, 0, 0)
22 airsr%minimum%q = maxmin_type(-missing_r, 0, 0)
23 airsr%average = residual_airsr1_type(0.0, 0.0)
24 airsr%rms_err = airsr%average
25
26 if (oi%num_airsr > 0) then
27 do n=1, oi%num_airsr
28 if (oi%airsr(n)%loc%proc_domain) then
29 do k=1, oi%airsr(n)%info%levels
30 call da_stats_calculate(oi%airsr(n)%loc%obs_global_index, &
31 k, oi%airsr(n)%t(k)%qc, &
32 oi%airsr(n)%t(k)%inv, nt, &
33 airsr%minimum%t, airsr%maximum%t, &
34 airsr%average%t, airsr%rms_err%t)
35 call da_stats_calculate(oi%airsr(n)%loc%obs_global_index, &
36 k, oi%airsr(n)%q(k)%qc, &
37 oi%airsr(n)%q(k)%inv, nq, &
38 airsr%minimum%q, airsr%maximum%q, &
39 airsr%average%q, airsr%rms_err%q)
40 end do
41 end if ! end if (oi%airsr(n)%loc%proc_domain)
42 end do
43 end if
44
45 ! Do inter-processor communication to gather statistics.
46 call da_proc_sum_int(nt)
47 call da_proc_sum_int(nq)
48
49 call da_proc_stats_combine(airsr%average%t, airsr%rms_err%t, &
50 airsr%minimum%t%value, airsr%maximum%t%value, &
51 airsr%minimum%t%n, airsr%maximum%t%n, &
52 airsr%minimum%t%l, airsr%maximum%t%l)
53 call da_proc_stats_combine(airsr%average%q, airsr%rms_err%q, &
54 airsr%minimum%q%value, airsr%maximum%q%value, &
55 airsr%minimum%q%n, airsr%maximum%q%n, &
56 airsr%minimum%q%l, airsr%maximum%q%l)
57
58 if (rootproc) then
59 if (nt /= 0 .or. nq /= 0) then
60 write(unit=stats_unit, fmt='(/a/)') ' Diagnostics of OI for airs retrievals'
61 call da_print_stats_airsr(stats_unit, nt, nq, airsr)
62 end if
63 end if
64
65 end subroutine da_oi_stats_airsr
66
67