1 subroutine da_proc_maxmin_combine(n, max, min)
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
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
11 !---------------------------------------------------------------------------
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.
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.
27 if (trace_use_frequent) call da_trace_entry("da_proc_maxmin_combine")
29 ! Get minimum value and associated processor index.
31 in(2*i-1) = min(i)%value
35 #if ( RWORDSIZE == 4 )
36 call mpi_reduce(in, out, n, mpi_2real, mpi_minloc, root, comm, ierr)
38 call mpi_reduce(in, out, n, mpi_2double_precision, mpi_minloc, root, comm, ierr)
43 min(i)%value = out(2*i-1)
44 proc_id(i) = int(out(2*i))
48 call wrf_dm_bcast_integer (proc_id, n)
50 ! Get i and j where minimum occurs.
52 if (proc_id(i) .ne. 0) then
54 call mpi_recv(min(i)%i, 1, mpi_integer, proc_id(i), 10, comm, STATUS, ierr)
55 call mpi_recv(min(i)%j, 1, mpi_integer, proc_id(i), 11, comm, STATUS, ierr)
56 else if (myproc == proc_id(i)) then
57 call mpi_send(min(i)%i, 1, mpi_integer, root, 10, comm, ierr)
58 call mpi_send(min(i)%j, 1, mpi_integer, root, 11, comm, ierr)
63 ! Get maximum value and associated processor index.
65 in(2*i-1) = max(i)%value
68 #if ( RWORDSIZE == 4 )
69 call mpi_reduce(in, out, n, mpi_2real, mpi_maxloc, root, comm, ierr)
71 call mpi_reduce(in, out, n, mpi_2double_precision, mpi_maxloc, root, comm, ierr)
76 max(i)%value = out(2*i-1)
77 proc_id(i) = int(out(2*i))
81 call wrf_dm_bcast_integer (proc_id, n)
83 ! Get i and j where maximum occurs.
85 if (proc_id(i) .ne. root) then
87 call mpi_recv(max(i)%i, 1, mpi_integer, proc_id(i), 10, comm, STATUS, ierr)
88 call mpi_recv(max(i)%j, 1, mpi_integer, proc_id(i), 11, comm, STATUS, ierr)
89 else if (myproc == proc_id(i)) then
90 call mpi_send(max(i)%i, 1, mpi_integer, root, 10, comm, ierr)
91 call mpi_send(max(i)%j, 1, mpi_integer, root, 11, comm, ierr)
96 if (trace_use_frequent) call da_trace_exit("da_proc_maxmin_combine")
99 end subroutine da_proc_maxmin_combine