da_oi_stats_airsr.inc
References to this file elsewhere.
1 subroutine da_oi_stats_airsr (stats_unit, iv)
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 (iv_type), intent (in) :: iv ! OI
11
12 type (stats_airsr_type) :: stats
13 integer :: nt, nq
14 integer :: n, k
15
16 if (trace_use_dull) call da_trace_entry("da_oi_stats_airsr")
17
18 nt = 0
19 nq = 0
20
21 stats%maximum%t = maxmin_type(missing_r, 0, 0)
22 stats%maximum%q = maxmin_type(missing_r, 0, 0)
23 stats%minimum%t = maxmin_type(-missing_r, 0, 0)
24 stats%minimum%q = maxmin_type(-missing_r, 0, 0)
25
26 stats%average = residual_airsr1_type(0.0, 0.0)
27 stats%rms_err = stats%average
28
29 do n=1, iv%info(airsr)%nlocal
30 if (iv%info(airsr)%proc_domain(1,n)) then
31 do k=1, iv%info(airsr)%levels(n)
32 call da_stats_calculate(iv%info(airsr)%obs_global_index(n), &
33 k, iv%airsr(n)%t(k)%qc, &
34 iv%airsr(n)%t(k)%inv, nt, &
35 stats%minimum%t, stats%maximum%t, &
36 stats%average%t, stats%rms_err%t)
37 call da_stats_calculate(iv%info(airsr)%obs_global_index(n), &
38 k, iv%airsr(n)%q(k)%qc, &
39 iv%airsr(n)%q(k)%inv, nq, &
40 stats%minimum%q, stats%maximum%q, &
41 stats%average%q, stats%rms_err%q)
42 end do
43 end if ! end if (iv%info(airsr)%proc_domain(1,n))
44 end do
45
46 ! Do inter-processor communication to gather statistics.
47 call da_proc_sum_int(nt)
48 call da_proc_sum_int(nq)
49
50 call da_proc_stats_combine(stats%average%t, stats%rms_err%t, &
51 stats%minimum%t%value, stats%maximum%t%value, &
52 stats%minimum%t%n, stats%maximum%t%n, &
53 stats%minimum%t%l, stats%maximum%t%l)
54 call da_proc_stats_combine(stats%average%q, stats%rms_err%q, &
55 stats%minimum%q%value, stats%maximum%q%value, &
56 stats%minimum%q%n, stats%maximum%q%n, &
57 stats%minimum%q%l, stats%maximum%q%l)
58
59 if (rootproc) then
60 if (nt /= 0 .or. nq /= 0) then
61 write(unit=stats_unit, fmt='(/a/)') ' Diagnostics of OI for airs retrievals'
62 call da_print_stats_airsr(stats_unit, nt, nq, stats)
63 end if
64 end if
65
66 if (trace_use_dull) call da_trace_exit("da_oi_stats_airsr")
67
68 end subroutine da_oi_stats_airsr
69
70