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 ! Get minimum value and associated processor index.
28 do i = 1, n
29 in(2*i-1) = min(i)%value
30 in(2*i) = myproc
31 end do
32
33 call mpi_reduce(in, out, n, mpi_2double_precision, mpi_minloc, &
34 root, comm, ierr)
35
36 if (rootproc) then
37 do i = 1, n
38 min(i)%value = out(2*i-1)
39 proc_id(i) = INT(out(2*i))
40 end do
41 end if
42
43 call wrf_dm_bcast_integer (proc_id, n)
44
45 ! Get i and j where minimum occurs.
46 do i = 1, n
47 if (proc_id(i) .ne. 0) then
48 if (rootproc) then
49 call mpi_recv(min(i)%i, 1, mpi_integer, proc_id(i), &
50 10, comm, STATUS, ierr)
51 call mpi_recv(min(i)%j, 1, mpi_integer, proc_id(i), &
52 11, comm, STATUS, ierr)
53 else if (myproc == proc_id(i)) then
54 call mpi_send(min(i)%i, 1, mpi_integer, root, &
55 10, comm, ierr)
56 call mpi_send(min(i)%j, 1, mpi_integer, root, &
57 11, comm, ierr)
58 end if
59 end if
60 end do
61
62 ! Get maximum value and associated processor index.
63 do i = 1, n
64 in(2*i-1) = max(i)%value
65 in(2*i) = myproc
66 end do
67 call mpi_reduce(in, out, n, mpi_2double_precision, mpi_maxloc, &
68 root, comm, ierr)
69
70 if (rootproc) then
71 do i = 1, n
72 max(i)%value = out(2*i-1)
73 proc_id(i) = INT(out(2*i))
74 end do
75 end if
76
77 call wrf_dm_bcast_integer (proc_id, n)
78
79 ! Get i and j where maximum occurs.
80 do i = 1, n
81 if (proc_id(i) .ne. root) then
82 if (rootproc) then
83 call mpi_recv(max(i)%i, 1, mpi_integer, proc_id(i), &
84 10, comm, STATUS, ierr)
85 call mpi_recv(max(i)%j, 1, mpi_integer, proc_id(i), &
86 11, comm, STATUS, ierr)
87 else if (myproc == proc_id(i)) then
88 call mpi_send(max(i)%i, 1, mpi_integer, root, &
89 10, comm, ierr)
90 call mpi_send(max(i)%j, 1, mpi_integer, root, &
91 11, comm, ierr)
92 end if
93 end if
94 end do
95 #endif
96
97 end subroutine da_proc_maxmin_combine
98
99