da_proc_maxmin_combine.inc

References to this file elsewhere.
1 subroutine da_proc_maxmin_combine(n, max, min)
2 
3    !---------------------------------------------------------------------------
4    !  Purpose: Do MPI reduction operations across processors to get the minimum
5    !           and maximum values for an observation field of length n. The
6    !           i, j location of the minimum and maximum, for each n, is also
7    !           communicated.
8    !           The return values are stored only on the root processor, i.e., 
9    !           processor 0.  (In this way, we do not have to do all-to-all 
10    !           communication.)
11    !---------------------------------------------------------------------------
12    
13    implicit none
14 
15    integer,                  intent(in)     :: n       ! Length of input fields.
16    type (maxmin_field_type), intent(inout)  :: max(n)  ! Max values over proc.
17    type (maxmin_field_type), intent(inout)  :: min(n)  ! Min values over proc.
18 
19    real    :: in(2*n)            ! mpi_reduce input value with processor rank.
20    real    :: out(2*n)           ! mpi_reduce output min/max with associated processor.
21    integer :: i                  ! Loop counter.
22    integer :: proc_id(n)         ! Id of processor with max or min value.
23    integer :: status(mpi_status_size) ! MPI status.
24 
25 #ifdef DM_PARALLEL
26 
27    if (trace_use_frequent) call da_trace_entry("da_proc_maxmin_combine")
28 
29    ! Get minimum value and associated processor index.
30    do i = 1, n
31       in(2*i-1) = min(i)%value
32       in(2*i)   = myproc
33    end do
34 
35    call mpi_reduce(in, out, n, mpi_2double_precision, mpi_minloc, root, comm, ierr)
36 
37    if (rootproc) then
38       do i = 1, n
39          min(i)%value = out(2*i-1)
40          proc_id(i)   = int(out(2*i))
41       end do
42    end if
43 
44    call wrf_dm_bcast_integer (proc_id, n)
45 
46    ! Get i and j where minimum occurs.
47    do i = 1, n
48       if (proc_id(i) .ne. 0) then
49          if (rootproc) then
50             call mpi_recv(min(i)%i, 1, mpi_integer, proc_id(i), 10, comm, STATUS, ierr)
51             call mpi_recv(min(i)%j, 1, mpi_integer, proc_id(i), 11, comm, STATUS, ierr)
52          else if (myproc == proc_id(i)) then
53             call mpi_send(min(i)%i, 1, mpi_integer, root, 10, comm, ierr)
54             call mpi_send(min(i)%j, 1, mpi_integer, root, 11, comm, ierr)
55          end if
56       end if
57    end do
58 
59    ! Get maximum value and associated processor index.
60    do i = 1, n
61       in(2*i-1) = max(i)%value
62       in(2*i)   = myproc
63    end do
64    call mpi_reduce(in, out, n, mpi_2double_precision, mpi_maxloc, root, comm, ierr)
65 
66    if (rootproc) then
67       do i = 1, n
68          max(i)%value = out(2*i-1)
69          proc_id(i)   = int(out(2*i))
70       end do
71    end if
72 
73    call wrf_dm_bcast_integer (proc_id, n)
74 
75    ! Get i and j where maximum occurs.
76    do i = 1, n
77       if (proc_id(i) .ne. root) then
78          if (rootproc) then
79             call mpi_recv(max(i)%i, 1, mpi_integer, proc_id(i), 10, comm, STATUS, ierr)
80             call mpi_recv(max(i)%j, 1, mpi_integer, proc_id(i), 11, comm, STATUS, ierr)
81          else if (myproc == proc_id(i)) then
82             call mpi_send(max(i)%i, 1, mpi_integer, root, 10, comm, ierr)
83             call mpi_send(max(i)%j, 1, mpi_integer, root, 11, comm, ierr)
84          end if
85       end if
86    end do
87 
88    if (trace_use_frequent) call da_trace_exit("da_proc_maxmin_combine")
89 #endif
90 
91 end subroutine da_proc_maxmin_combine
92 
93