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    ! Sum average and rms error over all processors and store on monitor processor.
33    call mpi_reduce(proc_ave, average, 1, true_mpi_real, mpi_sum, root, comm, ierr)
34    call mpi_reduce(proc_err, rms_err, 1, true_mpi_real, mpi_sum, root, comm, ierr)
35 
36    if (rootproc) then
37       proc_ave = average
38       proc_err = rms_err
39    end if
40 
41    ! Get minimum value and associated processor index.
42    in(1) = proc_min
43    in(2) = myproc
44    call mpi_reduce(in, out, 1, mpi_2double_precision, mpi_minloc, root, comm, ierr)
45 
46    if (myproc == root) then
47       proc_min = out(1)
48       proc_id = inT(out(2))
49    end if
50 
51    buf(1) = real(proc_id)
52    call wrf_dm_bcast_real (buf, 1)
53    proc_id=int(buf(1))
54 
55    ! Get obs number and k-level where minimum occurs.
56    if (proc_id .ne. root) then
57       if (rootproc) then
58          call mpi_recv(nobs_min, 1, mpi_integer, proc_id, 10, comm, STATUS, ierr)
59          call mpi_recv(klev_min, 1, mpi_integer, proc_id, 11, comm, STATUS, ierr)
60       else if (myproc == proc_id) then
61          call mpi_send(nobs_min, 1, mpi_integer, root, 10, comm, ierr)
62          call mpi_send(klev_min, 1, mpi_integer, root, 11, comm, ierr)
63       end if
64    end if
65 
66    ! Get maximum value and associated processor index.
67    in(1) = proc_max
68    in(2) = myproc
69 
70    call mpi_reduce(in, out, 1, mpi_2double_precision, mpi_maxloc, root, comm, ierr)
71 
72    if (rootproc) then
73       proc_max = out(1)
74       proc_id = int(out(2))
75    end if
76 
77    buf(1) = real(proc_id)
78    call wrf_dm_bcast_real (buf, 1)
79    proc_id=int(buf(1))
80 
81    ! Get obs number and k-level where maximum occurs.
82    if (proc_id .ne. root) then
83       if (rootproc) then
84          call mpi_recv(nobs_max, 1, mpi_integer, proc_id, 10, comm, STATUS, ierr)
85          call mpi_recv(klev_max, 1, mpi_integer, proc_id, 11, comm, STATUS, ierr)
86       else if (myproc == proc_id) then
87          call mpi_send(nobs_max, 1, mpi_integer, root, 10, comm, ierr)
88          call mpi_send(klev_max, 1, mpi_integer, root, 11, comm, ierr)
89       end if
90    end if
91 #endif
92 
93 end subroutine da_proc_stats_combine
94 
95