1 subroutine da_oi_stats_sound (stats_unit, iv, ob)
3 !-----------------------------------------------------------------------
5 !-----------------------------------------------------------------------
9 integer, intent (in) :: stats_unit ! Output unit for stats.
10 type (iv_type), intent (in) :: iv ! OI
11 type(y_type), intent (in) :: ob ! Observation structure.
13 type (stats_sound_type) :: stats
14 integer :: nu, nv, nt, nq
16 real :: u_inv, v_inv, u_obs, v_obs
18 if (trace_use_dull) call da_trace_entry("da_oi_stats_sound")
25 stats%maximum%u = maxmin_type(missing_r, 0, 0)
26 stats%maximum%v = maxmin_type(missing_r, 0, 0)
27 stats%maximum%t = maxmin_type(missing_r, 0, 0)
28 stats%maximum%q = maxmin_type(missing_r, 0, 0)
29 stats%minimum%u = maxmin_type(-missing_r, 0, 0)
30 stats%minimum%v = maxmin_type(-missing_r, 0, 0)
31 stats%minimum%t = maxmin_type(-missing_r, 0, 0)
32 stats%minimum%q = maxmin_type(-missing_r, 0, 0)
33 stats%average = residual_sound1_type(0.0, 0.0, 0.0, 0.0)
34 stats%rms_err = stats%average
36 do n=1, iv%info(sound)%nlocal
37 do k=1, iv%info(sound)%levels(n)
38 ! if (iv%info(sound)%proc_domain(k,n)) then
39 if (iv%info(sound)%proc_domain(1,n)) then
41 u_inv = iv%sound(n)%u(k)%inv
42 v_inv = iv%sound(n)%v(k)%inv
43 u_obs = ob%sound(n)%u(k)
44 v_obs = ob%sound(n)%v(k)
46 if (.not. wind_sd_sound .and. wind_stats_sd) &
47 call da_ffdduv_diagnose(u_obs, u_inv, u_obs, v_obs, v_inv, v_obs, &
48 iv%sound(n)%u(k)%qc, iv%sound(n)%v(k)%qc, convert_uv2fd)
49 if (wind_sd_sound .and. .not. wind_stats_sd) &
50 call da_ffdduv_diagnose(u_obs, u_inv, u_obs, v_obs, v_inv, v_obs, &
51 iv%sound(n)%u(k)%qc, iv%sound(n)%v(k)%qc, convert_fd2uv)
53 call da_stats_calculate(iv%info(sound)%obs_global_index(n), &
54 k, iv%sound(n)%u(k)%qc, u_inv, nu, &
55 stats%minimum%u, stats%maximum%u, stats%average%u, stats%rms_err%u)
56 call da_stats_calculate(iv%info(sound)%obs_global_index(n), &
57 k, iv%sound(n)%v(k)%qc, v_inv, nv, &
58 stats%minimum%v, stats%maximum%v, stats%average%v, stats%rms_err%v)
59 call da_stats_calculate(iv%info(sound)%obs_global_index(n), &
60 k, iv%sound(n)%t(k)%qc, iv%sound(n)%t(k)%inv, nt, &
61 stats%minimum%t, stats%maximum%t, stats%average%t, stats%rms_err%t)
62 call da_stats_calculate(iv%info(sound)%obs_global_index(n), &
63 k, iv%sound(n)%q(k)%qc, iv%sound(n)%q(k)%inv, nq, &
64 stats%minimum%q, stats%maximum%q, stats%average%q, stats%rms_err%q)
69 ! Do inter-processor communication to gather statistics.
70 call da_proc_sum_int(nu)
71 call da_proc_sum_int(nv)
72 call da_proc_sum_int(nt)
73 call da_proc_sum_int(nq)
75 call da_proc_stats_combine(stats%average%u, stats%rms_err%u, &
76 stats%minimum%u%value, stats%maximum%u%value, &
77 stats%minimum%u%n, stats%maximum%u%n, &
78 stats%minimum%u%l, stats%maximum%u%l)
79 call da_proc_stats_combine(stats%average%v, stats%rms_err%v, &
80 stats%minimum%v%value, stats%maximum%v%value, &
81 stats%minimum%v%n, stats%maximum%v%n, &
82 stats%minimum%v%l, stats%maximum%v%l)
83 call da_proc_stats_combine(stats%average%t, stats%rms_err%t, &
84 stats%minimum%t%value, stats%maximum%t%value, &
85 stats%minimum%t%n, stats%maximum%t%n, &
86 stats%minimum%t%l, stats%maximum%t%l)
87 call da_proc_stats_combine(stats%average%q, stats%rms_err%q, &
88 stats%minimum%q%value, stats%maximum%q%value, &
89 stats%minimum%q%n, stats%maximum%q%n, &
90 stats%minimum%q%l, stats%maximum%q%l)
93 if (nu /= 0 .or. nv /= 0 .or. nt /= 0 .or. nq /= 0) then
94 write(unit=stats_unit, fmt='(/a/)') ' Diagnostics of OI for sound'
95 call da_print_stats_sound(stats_unit, nu, nv, nt, nq, stats)
99 if (trace_use_dull) call da_trace_exit("da_oi_stats_sound")
101 end subroutine da_oi_stats_sound