Merge remote-tracking branch 'origin/release-v4.5.2'
[WRF.git] / var / da / da_par_util / da_proc_stats_combine.inc
blobd8fd895f57c615f7bc37caef75f562de5d31e31c
1 subroutine da_proc_stats_combine(proc_ave, proc_err, proc_min, proc_max, &
2    nobs_min, nobs_max, klev_min, klev_max)
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    !---------------------------------------------------------------------------
11    implicit none
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.
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)
30 #ifdef DM_PARALLEL
32    if (trace_use_frequent) call da_trace_entry("da_proc_stats_combine")
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)
38    if (rootproc) then
39       proc_ave = average
40       proc_err = rms_err
41    end if
43    ! Get minimum value and associated processor index.
44    in(1) = proc_min
45    in(2) = myproc
46 #if ( RWORDSIZE == 4 )
47    call mpi_reduce(in, out, 1, mpi_2real, mpi_minloc, root, comm, ierr)
48 #else
49    call mpi_reduce(in, out, 1, mpi_2double_precision, mpi_minloc, root, comm, ierr)
50 #endif
52    if (myproc == root) then
53       proc_min = out(1)
54       proc_id = inT(out(2))
55       buf(1) = real(proc_id)
56    end if
58    call wrf_dm_bcast_real (buf, 1)
59    proc_id=int(buf(1))
61    ! Get obs number and k-level where minimum occurs.
62    if (proc_id .ne. root) then
63       if (rootproc) then
64          call mpi_recv(nobs_min, 1, mpi_integer, proc_id, 10, comm, STATUS, ierr)
65          call mpi_recv(klev_min, 1, mpi_integer, proc_id, 11, comm, STATUS, ierr)
66       else if (myproc == proc_id) then
67          call mpi_send(nobs_min, 1, mpi_integer, root, 10, comm, ierr)
68          call mpi_send(klev_min, 1, mpi_integer, root, 11, comm, ierr)
69       end if
70    end if
72    ! Get maximum value and associated processor index.
73    in(1) = proc_max
74    in(2) = myproc
76 #if ( RWORDSIZE == 4 )
77    call mpi_reduce(in, out, 1, mpi_2real, mpi_maxloc, root, comm, ierr)
78 #else
79    call mpi_reduce(in, out, 1, mpi_2double_precision, mpi_maxloc, root, comm, ierr)
80 #endif
82    if (rootproc) then
83       proc_max = out(1)
84       proc_id = int(out(2))
85       buf(1) = real(proc_id)
86    end if
88    call wrf_dm_bcast_real (buf, 1)
89    proc_id=int(buf(1))
91    ! Get obs number and k-level where maximum occurs.
92    if (proc_id .ne. root) then
93       if (rootproc) then
94          call mpi_recv(nobs_max, 1, mpi_integer, proc_id, 10, comm, STATUS, ierr)
95          call mpi_recv(klev_max, 1, mpi_integer, proc_id, 11, comm, STATUS, ierr)
96       else if (myproc == proc_id) then
97          call mpi_send(nobs_max, 1, mpi_integer, root, 10, comm, ierr)
98          call mpi_send(klev_max, 1, mpi_integer, root, 11, comm, ierr)
99       end if
100    end if
102    if (trace_use_frequent) call da_trace_exit("da_proc_stats_combine")
103 #endif
105 end subroutine da_proc_stats_combine