da_ao_stats_gpsref.inc

References to this file elsewhere.
1 subroutine da_ao_stats_gpsref (stats_unit, iv, re)
2 
3    !--------------------------------------------------------------------
4    ! Purpose: Called by da_minimisation/da_write_diagnostics.inc.
5    !--------------------------------------------------------------------
6 
7    implicit none
8 
9    integer,        intent(in)    :: stats_unit    ! Output unit for stats.
10    type (iv_type), intent(inout) :: iv            ! iv
11    type  (y_type), intent(in)    :: re            ! A - O
12 
13    type (stats_gpsref_type)         :: stats
14    integer                          :: ngpsref
15    integer                          :: n, k
16    real                             :: o_minus_b, o_minus_a, sigma_o, sigma_b
17    real                             :: o_minus_b_0, o_minus_a_0, sigma_o_0, sigma_b_0
18 
19    if (trace_use_dull) call da_trace_entry("da_ao_stats_gpsref")
20 
21    ngpsref = 0
22    o_minus_b_0 = 0.0; o_minus_a_0 = 0.0; sigma_o_0 = 0.0; sigma_b_0 = 0.0
23    
24    stats%maximum%ref = maxmin_type (missing_r, 0, 0)
25    stats%minimum%ref = maxmin_type(-missing_r, 0, 0)
26 
27    stats%average = residual_gpsref1_type(0.0,0.0,0.0,0.0)
28    stats%rms_err = stats%average
29 
30    do n=1, iv%info(gpsref)%nlocal
31       if (iv%info(gpsref)%proc_domain(1,n)) then
32          do k=1, iv%info(gpsref)%levels(n)
33             call da_stats_calculate (n, k, iv%gpsref(n)%ref(k)%qc, & 
34                re%gpsref(n)%ref(k), ngpsref, &
35                stats%minimum%ref, stats%maximum%ref, &
36                stats%average%ref, stats%rms_err%ref)
37 
38             if (pseudo_var(1:3) == 'ref' .and. num_pseudo > 0) then
39                o_minus_b = iv%GPSRef(n)%ref(k)%inv
40                o_minus_a = re%gpsref(n)%ref(k)
41                sigma_o   = iv%gpsref(n)%ref(k)%error
42             end if
43          end do
44 
45          if (pseudo_var(1:3) == 'ref' .and. num_pseudo > 0) then
46             ! Calculate equivalent sigma_b using
47             ! O-A=(O-B)*sigma_o/(sigma_o+sigma_b)
48 
49             sigma_b = sqrt ((o_minus_b - o_minus_a) &
50                / o_minus_a) * sigma_o
51             o_minus_b_0 = wrf_dm_sum_real (o_minus_b)
52             o_minus_a_0 = wrf_dm_sum_real (o_minus_a)
53             sigma_o_0 = wrf_dm_sum_real (sigma_o)
54             sigma_b_0 = wrf_dm_sum_real (sigma_b)
55             write (unit=stdout,fmt='(A,F10.2)') &
56                'TEST_COVERAGE_da_ao_stats_gpsref:  o_minus_b_0 = ', o_minus_b_0 
57             write (unit=stdout,fmt='(A,F10.2)') &
58                'TEST_COVERAGE_da_ao_stats_gpsref:  o_minus_a_0 = ', o_minus_a_0
59             write (unit=stdout,fmt='(A,F10.2)') & 
60                'TEST_COVERAGE_da_ao_stats_gpsref:  sigma_o_0 = ', sigma_o_0
61             write (unit=stdout,fmt='(A,F10.2)') &
62                'TEST_COVERAGE_da_ao_stats_gpsref:  sigma_b_0 = ', sigma_b_0
63             if (rootproc) then 
64                write(stats_unit,'(/A,A3,A,f12.3)')  & 
65                   ' Pseudo ', pseudo_var, ' O-B: ', o_minus_b_0 
66                write(stats_unit,' (A,A3,A,f12.3)')  & 
67                   ' Pseudo ', pseudo_var, ' O-A: ', o_minus_a_0 
68                write(stats_unit,' (A,A3,A,f12.3)')  & 
69                   ' Pseudo ', pseudo_var, ' sigma_o: ', sigma_o_0 
70                write(stats_unit,'(A,A3,A,f12.3)')  & 
71                   ' Pseudo ', pseudo_var, ' sigma_b: ', sigma_b_0
72             end if 
73          end if
74       end if    ! end if (iv%info(gpsref)%proc_domain(1,n))
75    end do
76 
77    ! Do inter-processor communication to gather statistics.
78 
79    call da_proc_sum_int (ngpsref)
80    iv%nstats(gpsref) = ngpsref
81     
82    call da_proc_stats_combine(stats%average%ref, stats%rms_err%ref, &
83       stats%minimum%ref%value, stats%maximum%ref%value, &
84       stats%minimum%ref%n, stats%maximum%ref%n, &
85       stats%minimum%ref%l, stats%maximum%ref%l)
86    
87    if (rootproc) then
88       if (ngpsref > 0) then
89          write(unit=stats_unit, fmt='(/a/)') ' Diagnostics of AO for gpsref'
90             call da_print_stats_gpsref(stats_unit, ngpsref, stats)
91       end if
92    end if
93 
94    if (trace_use_dull) call da_trace_exit("da_ao_stats_gpsref")
95 
96 end subroutine da_ao_stats_gpsref
97 
98