<HTML> <BODY BGCOLOR=#ccccdd LINK=#0000aa VLINK=#0000ff ALINK=#ff0000 ><BASE TARGET="bottom_target"><PRE>
<A NAME='DA_AO_STATS_AIREP'><A href='../../html_code/airep/da_ao_stats_airep.inc.html#DA_AO_STATS_AIREP' TARGET='top_target'><IMG SRC="../../gif/bar_red.gif" border=0></A>

subroutine da_ao_stats_airep (stats_unit, iv, re, ob) 1,17

   !-----------------------------------------------------------------------
   ! Purpose: TBD
   !-----------------------------------------------------------------------

   implicit none

   integer,        intent(in)    :: stats_unit    ! Output unit for stats.
   type (iv_type), intent(inout) :: iv            ! iv
   type  (y_type), intent(in)    :: re            ! A - O
   type(y_type),   intent (in)   :: ob            ! Observation structure.

   type (stats_airep_type) :: stats
   integer                 :: nu, nv, nt, nq
   integer                 :: n, k
   real                    :: u_inc, v_inc, u_obs, v_obs
   
   if (trace_use_dull) call da_trace_entry("da_ao_stats_airep")

   nu = 0
   nv = 0
   nt = 0
   nq = 0

   stats%maximum%u = maxmin_type (missing_r, 0, 0)
   stats%maximum%v = maxmin_type (missing_r, 0, 0)
   stats%maximum%t = maxmin_type (missing_r, 0, 0)
   stats%maximum%q = maxmin_type (missing_r, 0, 0)
   stats%minimum%u = maxmin_type(-missing_r, 0, 0)
   stats%minimum%v = maxmin_type(-missing_r, 0, 0)
   stats%minimum%t = maxmin_type(-missing_r, 0, 0)
   stats%minimum%q = maxmin_type(-missing_r, 0, 0)
   stats%average = residual_airep1_type(0.0, 0.0, 0.0, 0.0)
   stats%rms_err = stats%average

   do n=1, iv%info(airep)%nlocal
      do k=1, iv%info(airep)%levels(n)
         if (iv%info(airep)%proc_domain(k,n)) then

            u_inc = re%airep(n)%u(k)
            v_inc = re%airep(n)%v(k)
            u_obs = ob%airep(n)%u(k)
            v_obs = ob%airep(n)%v(k)

            if (.not. wind_sd_airep .and. wind_stats_sd) &amp;
               call da_ffdduv_diagnose(u_obs,u_obs,u_inc,v_obs,v_obs,v_inc, &amp;
                                       iv%airep(n)%u(k)%qc,iv%airep(n)%v(k)%qc, convert_uv2fd)
            if (wind_sd_airep .and. .not. wind_stats_sd) &amp;
               call da_ffdduv_diagnose(u_obs,u_obs,u_inc,v_obs,v_obs,v_inc, &amp;
                                       iv%airep(n)%u(k)%qc,iv%airep(n)%v(k)%qc, convert_fd2uv)

            call da_stats_calculate (n, k, iv%airep(n)%u(k)%qc,  &amp; 
               u_inc, nu, &amp;
               stats%minimum%u, stats%maximum%u, &amp;
               stats%average%u, stats%rms_err%u)
            call da_stats_calculate (n, k, iv%airep(n)%v(k)%qc,  &amp; 
               v_inc, nv, &amp;
               stats%minimum%v, stats%maximum%v, &amp;
               stats%average%v, stats%rms_err%v)
            call da_stats_calculate (n, k, iv%airep(n)%t(k)%qc,  &amp; 
               re%airep(n)%t(k), nt, &amp;
               stats%minimum%t, stats%maximum%t, &amp;
               stats%average%t, stats%rms_err%t)
            call da_stats_calculate (n, k, iv%airep(n)%q(k)%qc,  &amp;
               re%airep(n)%q(k), nq, &amp;
               stats%minimum%q, stats%maximum%q, &amp;
               stats%average%q, stats%rms_err%q)
        end if
      end do
   end do

   ! Do inter-processor communication to gather statistics.
   call da_proc_sum_int (nu)
   call da_proc_sum_int (nv)
   call da_proc_sum_int (nt)
   call da_proc_sum_int (nq)
   iv%nstats(airep) = nu + nv + nt + nq

   call da_proc_stats_combine(stats%average%u, stats%rms_err%u, &amp;
      stats%minimum%u%value, stats%maximum%u%value, &amp;
      stats%minimum%u%n, stats%maximum%u%n, &amp;
      stats%minimum%u%l, stats%maximum%u%l)
   call da_proc_stats_combine(stats%average%v, stats%rms_err%v, &amp;
      stats%minimum%v%value, stats%maximum%v%value, &amp;
      stats%minimum%v%n, stats%maximum%v%n, &amp;
      stats%minimum%v%l, stats%maximum%v%l)
   call da_proc_stats_combine(stats%average%t, stats%rms_err%t, &amp;
      stats%minimum%t%value, stats%maximum%t%value, &amp;
      stats%minimum%t%n, stats%maximum%t%n, &amp;
      stats%minimum%t%l, stats%maximum%t%l)
   call da_proc_stats_combine(stats%average%q, stats%rms_err%q, &amp;
      stats%minimum%q%value, stats%maximum%q%value, &amp;
      stats%minimum%q%n, stats%maximum%q%n, &amp;
      stats%minimum%q%l, stats%maximum%q%l)

   if (rootproc) then
      if (nu /= 0 .or. nv /= 0 .or. nt /= 0 .or. nq /= 0) then
         write(unit=stats_unit, fmt='(/a/)') ' Diagnostics of AO for airep'
         call da_print_stats_airep(stats_unit, nu, nv, nt, nq, stats)
      end if
   end if
   
   if (trace_use_dull) call da_trace_exit("da_ao_stats_airep")

 end subroutine da_ao_stats_airep