da_proc_stats_combine.inc

References to this file elsewhere.
1 subroutine da_proc_stats_combine(proc_ave, proc_err, proc_min, proc_max, &
2    nobs_min, nobs_max, klev_min, klev_max)
3 
4    !---------------------------------------------------------------------------
5    !  Purpose: Do MPI reduction operations across processors to get the average, 
6    !           rms error, minimum, and maximum values for an observation field.
7    !           These are stored only on the root processor, i.e., processor 0.
8    !           (In this way, we do not have to do all-to-all communication.)
9    !---------------------------------------------------------------------------
10 
11    implicit none
12 
13    real,      intent(inout)      :: proc_ave       ! Processor average.
14    real,      intent(inout)      :: proc_err       ! Processor rms error.
15    real,      intent(inout)      :: proc_min       ! Processor minumum.
16    real,      intent(inout)      :: proc_max       ! Processor maximum.
17    integer,   intent(inout)      :: nobs_min       ! Obs number of minimum.
18    integer,   intent(inout)      :: nobs_max       ! Obs number of maximum.
19    integer,   intent(inout)      :: klev_min       ! Level of minimum.
20    integer,   intent(inout)      :: klev_max       ! Level of maximum.
21 
22    real    :: average            ! Global average.
23    real    :: rms_err            ! Global rms_error.
24    real    :: in(2)              ! mpi_reduce input value with processor rank.
25    real    :: out(2)             ! mpi_reduce output min/max with associated processor.
26    integer :: proc_id            ! Id of processor with max or min value.
27    integer :: status(mpi_status_size) ! MPI status.
28    real    :: buf(1)
29 
30 #ifdef DM_PARALLEL
31 
32    if (trace_use_frequent) call da_trace_entry("da_proc_stats_combine")
33 
34    ! Sum average and rms error over all processors and store on monitor processor.
35    call mpi_reduce(proc_ave, average, 1, true_mpi_real, mpi_sum, root, comm, ierr)
36    call mpi_reduce(proc_err, rms_err, 1, true_mpi_real, mpi_sum, root, comm, ierr)
37 
38    if (rootproc) then
39       proc_ave = average
40       proc_err = rms_err
41    end if
42 
43    ! Get minimum value and associated processor index.
44    in(1) = proc_min
45    in(2) = myproc
46    call mpi_reduce(in, out, 1, mpi_2double_precision, mpi_minloc, root, comm, ierr)
47 
48    if (myproc == root) then
49       proc_min = out(1)
50       proc_id = inT(out(2))
51    end if
52 
53    buf(1) = real(proc_id)
54    call wrf_dm_bcast_real (buf, 1)
55    proc_id=int(buf(1))
56 
57    ! Get obs number and k-level where minimum occurs.
58    if (proc_id .ne. root) then
59       if (rootproc) then
60          call mpi_recv(nobs_min, 1, mpi_integer, proc_id, 10, comm, STATUS, ierr)
61          call mpi_recv(klev_min, 1, mpi_integer, proc_id, 11, comm, STATUS, ierr)
62       else if (myproc == proc_id) then
63          call mpi_send(nobs_min, 1, mpi_integer, root, 10, comm, ierr)
64          call mpi_send(klev_min, 1, mpi_integer, root, 11, comm, ierr)
65       end if
66    end if
67 
68    ! Get maximum value and associated processor index.
69    in(1) = proc_max
70    in(2) = myproc
71 
72    call mpi_reduce(in, out, 1, mpi_2double_precision, mpi_maxloc, root, comm, ierr)
73 
74    if (rootproc) then
75       proc_max = out(1)
76       proc_id = int(out(2))
77    end if
78 
79    buf(1) = real(proc_id)
80    call wrf_dm_bcast_real (buf, 1)
81    proc_id=int(buf(1))
82 
83    ! Get obs number and k-level where maximum occurs.
84    if (proc_id .ne. root) then
85       if (rootproc) then
86          call mpi_recv(nobs_max, 1, mpi_integer, proc_id, 10, comm, STATUS, ierr)
87          call mpi_recv(klev_max, 1, mpi_integer, proc_id, 11, comm, STATUS, ierr)
88       else if (myproc == proc_id) then
89          call mpi_send(nobs_max, 1, mpi_integer, root, 10, comm, ierr)
90          call mpi_send(klev_max, 1, mpi_integer, root, 11, comm, ierr)
91       end if
92    end if
93 
94    if (trace_use_frequent) call da_trace_exit("da_proc_stats_combine")
95 #endif
96 
97 end subroutine da_proc_stats_combine
98 
99