1 | SUBROUTINE global_min_max( i1, i2, j1, j2, k1, k2, ar, mode, value, & |
---|
2 | value_ijk, value1, value1_ijk ) |
---|
3 | |
---|
4 | !-------------------------------------------------------------------------------! |
---|
5 | ! Actual revisions: |
---|
6 | ! ----------------- |
---|
7 | ! |
---|
8 | ! |
---|
9 | ! Former revisions: |
---|
10 | ! ----------------- |
---|
11 | ! $Log: global_min_max.f90,v $ |
---|
12 | ! Revision 1.11 2003/04/16 12:56:58 raasch |
---|
13 | ! Index values of the extrema are limited to the range 0..nx, 0..ny |
---|
14 | ! |
---|
15 | ! Revision 1.10 2003/03/16 09:39:21 raasch |
---|
16 | ! Two underscores (_) are placed in front of all define-strings |
---|
17 | ! |
---|
18 | ! Revision 1.9 2003/03/14 13:42:50 raasch |
---|
19 | ! IF clause extracted from loop which determines the absolut array maximum, |
---|
20 | ! because otherwise it prevents vectorization of that loop |
---|
21 | ! |
---|
22 | ! Revision 1.8 2002/12/19 14:51:03 raasch |
---|
23 | ! Speed optimization by removing MINVAL/MAXVAL calls and by handling |
---|
24 | ! the "abs" case in a different way than the min/max cases, |
---|
25 | ! translation of remaining variables |
---|
26 | ! |
---|
27 | ! Revision 1.7 2002/06/11 13:08:16 raasch |
---|
28 | ! Last change limited to IBM since it does not work on 64-bit machines |
---|
29 | ! |
---|
30 | ! Revision 1.6 2002/05/02 18:50:36 18:50:36 raasch (Siegfried Raasch) |
---|
31 | ! Type of fmax, fmax_l, fmin and fmin_l changed to REAL (KIND=4). |
---|
32 | ! |
---|
33 | ! Revision 1.5 2001/01/22 06:51:02 raasch |
---|
34 | ! Module test_variables removed |
---|
35 | ! |
---|
36 | ! Revision 1.4 2000/12/20 12:19:57 letzel |
---|
37 | ! All comments translated into English. |
---|
38 | ! |
---|
39 | ! Revision 1.3 1998/07/06 12:14:37 raasch |
---|
40 | ! + USE test_variables |
---|
41 | ! |
---|
42 | ! Revision 1.2 1997/08/11 06:16:54 raasch |
---|
43 | ! Indices des ermittelten Minumums/Maximums werden an PEs gesendet |
---|
44 | ! |
---|
45 | ! Revision 1.1 1997/07/24 11:14:03 raasch |
---|
46 | ! Initial revision |
---|
47 | ! |
---|
48 | ! |
---|
49 | ! Description: |
---|
50 | ! ------------ |
---|
51 | ! Determine the array minimum/maximum and the corresponding indices. |
---|
52 | !-------------------------------------------------------------------------------! |
---|
53 | |
---|
54 | USE indices |
---|
55 | USE pegrid |
---|
56 | |
---|
57 | IMPLICIT NONE |
---|
58 | |
---|
59 | CHARACTER (LEN=*) :: mode |
---|
60 | |
---|
61 | INTEGER :: i, i1, i2, id_fmax, id_fmin, j, j1, j2, k, k1, k2, & |
---|
62 | fmax_ijk(3), fmax_ijk_l(3), fmin_ijk(3), & |
---|
63 | fmin_ijk_l(3), value_ijk(3) |
---|
64 | INTEGER, OPTIONAL :: value1_ijk(3) |
---|
65 | REAL :: value, & |
---|
66 | ar(i1:i2,j1:j2,k1:k2) |
---|
67 | #if defined( __ibm ) |
---|
68 | REAL (KIND=4) :: fmax(2), fmax_l(2), fmin(2), fmin_l(2) ! on 32bit- |
---|
69 | ! machines MPI_2REAL must not be replaced by |
---|
70 | ! MPI_2DOUBLE_PRECISION |
---|
71 | #else |
---|
72 | REAL :: fmax(2), fmax_l(2), fmin(2), fmin_l(2) |
---|
73 | #endif |
---|
74 | REAL, OPTIONAL :: value1 |
---|
75 | |
---|
76 | |
---|
77 | ! |
---|
78 | !-- Determine array minimum |
---|
79 | IF ( mode == 'min' .OR. mode == 'minmax' ) THEN |
---|
80 | |
---|
81 | ! |
---|
82 | !-- Determine the local minimum |
---|
83 | fmin_ijk_l = MINLOC( ar ) |
---|
84 | fmin_ijk_l(1) = i1 + fmin_ijk_l(1) - 1 ! MINLOC assumes lowerbound = 1 |
---|
85 | fmin_ijk_l(2) = j1 + fmin_ijk_l(2) - 1 |
---|
86 | fmin_ijk_l(3) = k1 + fmin_ijk_l(3) - 1 |
---|
87 | fmin_l(1) = ar(fmin_ijk_l(1),fmin_ijk_l(2),fmin_ijk_l(3)) |
---|
88 | |
---|
89 | #if defined( __parallel ) |
---|
90 | fmin_l(2) = myid |
---|
91 | CALL MPI_ALLREDUCE( fmin_l, fmin, 1, MPI_2REAL, MPI_MINLOC, comm2d, ierr ) |
---|
92 | |
---|
93 | ! |
---|
94 | !-- Determine the global minimum. Result stored on PE0. |
---|
95 | id_fmin = fmin(2) |
---|
96 | IF ( id_fmin /= 0 ) THEN |
---|
97 | IF ( myid == 0 ) THEN |
---|
98 | CALL MPI_RECV( fmin_ijk, 3, MPI_INTEGER, id_fmin, 0, comm2d, & |
---|
99 | status, ierr ) |
---|
100 | ELSEIF ( myid == id_fmin ) THEN |
---|
101 | CALL MPI_SEND( fmin_ijk_l, 3, MPI_INTEGER, 0, 0, comm2d, ierr ) |
---|
102 | ENDIF |
---|
103 | ELSE |
---|
104 | fmin_ijk = fmin_ijk_l |
---|
105 | ENDIF |
---|
106 | ! |
---|
107 | !-- Send the indices of the just determined array minimum to other PEs |
---|
108 | CALL MPI_BCAST( fmin_ijk, 3, MPI_INTEGER, 0, comm2d, ierr ) |
---|
109 | #else |
---|
110 | fmin(1) = fmin_l(1) |
---|
111 | fmin_ijk = fmin_ijk_l |
---|
112 | #endif |
---|
113 | |
---|
114 | ENDIF |
---|
115 | |
---|
116 | ! |
---|
117 | !-- Determine array maximum |
---|
118 | IF ( mode == 'max' .OR. mode == 'minmax' ) THEN |
---|
119 | |
---|
120 | ! |
---|
121 | !-- Determine the local maximum |
---|
122 | fmax_ijk_l = MAXLOC( ar ) |
---|
123 | fmax_ijk_l(1) = i1 + fmax_ijk_l(1) - 1 ! MAXLOC assumes lowerbound = 1 |
---|
124 | fmax_ijk_l(2) = j1 + fmax_ijk_l(2) - 1 |
---|
125 | fmax_ijk_l(3) = k1 + fmax_ijk_l(3) - 1 |
---|
126 | fmax_l(1) = ar(fmax_ijk_l(1),fmax_ijk_l(2),fmax_ijk_l(3)) |
---|
127 | |
---|
128 | #if defined( __parallel ) |
---|
129 | fmax_l(2) = myid |
---|
130 | CALL MPI_ALLREDUCE( fmax_l, fmax, 1, MPI_2REAL, MPI_MAXLOC, comm2d, ierr ) |
---|
131 | |
---|
132 | ! |
---|
133 | !-- Determine the global maximum. Result stored on PE0. |
---|
134 | id_fmax = fmax(2) |
---|
135 | IF ( id_fmax /= 0 ) THEN |
---|
136 | IF ( myid == 0 ) THEN |
---|
137 | CALL MPI_RECV( fmax_ijk, 3, MPI_INTEGER, id_fmax, 0, comm2d, & |
---|
138 | status, ierr ) |
---|
139 | ELSEIF ( myid == id_fmax ) THEN |
---|
140 | CALL MPI_SEND( fmax_ijk_l, 3, MPI_INTEGER, 0, 0, comm2d, ierr ) |
---|
141 | ENDIF |
---|
142 | ELSE |
---|
143 | fmax_ijk = fmax_ijk_l |
---|
144 | ENDIF |
---|
145 | ! |
---|
146 | !-- send the indices of the just determined array maximum to other PEs |
---|
147 | CALL MPI_BCAST( fmax_ijk, 3, MPI_INTEGER, 0, comm2d, ierr ) |
---|
148 | #else |
---|
149 | fmax(1) = fmax_l(1) |
---|
150 | fmax_ijk = fmax_ijk_l |
---|
151 | #endif |
---|
152 | |
---|
153 | ENDIF |
---|
154 | |
---|
155 | ! |
---|
156 | !-- Determine absolute array maximum |
---|
157 | IF ( mode == 'abs' ) THEN |
---|
158 | |
---|
159 | ! |
---|
160 | !-- Determine the local absolut maximum |
---|
161 | fmax_l(1) = 0.0 |
---|
162 | fmax_ijk_l(1) = i1 |
---|
163 | fmax_ijk_l(2) = j1 |
---|
164 | fmax_ijk_l(3) = k1 |
---|
165 | DO k = k1, k2 |
---|
166 | DO j = j1, j2 |
---|
167 | DO i = i1, i2 |
---|
168 | IF ( ABS( ar(i,j,k) ) > fmax_l(1) ) THEN |
---|
169 | fmax_l(1) = ABS( ar(i,j,k) ) |
---|
170 | fmax_ijk_l(1) = i |
---|
171 | fmax_ijk_l(2) = j |
---|
172 | fmax_ijk_l(3) = k |
---|
173 | ENDIF |
---|
174 | ENDDO |
---|
175 | ENDDO |
---|
176 | ENDDO |
---|
177 | |
---|
178 | ! |
---|
179 | !-- Set a flag in case that the determined value is negative. |
---|
180 | !-- A constant offset has to be subtracted in order to handle the special |
---|
181 | !-- case i=0 correctly |
---|
182 | IF ( ar(fmax_ijk_l(1),fmax_ijk_l(2),fmax_ijk_l(3)) < 0.0 ) THEN |
---|
183 | fmax_ijk_l(1) = -fmax_ijk_l(1) - 10 |
---|
184 | ENDIF |
---|
185 | |
---|
186 | #if defined( __parallel ) |
---|
187 | fmax_l(2) = myid |
---|
188 | CALL MPI_ALLREDUCE( fmax_l, fmax, 1, MPI_2REAL, MPI_MAXLOC, comm2d, & |
---|
189 | ierr ) |
---|
190 | |
---|
191 | ! |
---|
192 | !-- Determine the global absolut maximum. Result stored on PE0. |
---|
193 | id_fmax = fmax(2) |
---|
194 | IF ( id_fmax /= 0 ) THEN |
---|
195 | IF ( myid == 0 ) THEN |
---|
196 | CALL MPI_RECV( fmax_ijk, 3, MPI_INTEGER, id_fmax, 0, comm2d, & |
---|
197 | status, ierr ) |
---|
198 | ELSEIF ( myid == id_fmax ) THEN |
---|
199 | CALL MPI_SEND( fmax_ijk_l, 3, MPI_INTEGER, 0, 0, comm2d, ierr ) |
---|
200 | ENDIF |
---|
201 | ELSE |
---|
202 | fmax_ijk = fmax_ijk_l |
---|
203 | ENDIF |
---|
204 | ! |
---|
205 | !-- Send the indices of the just determined absolut maximum to other PEs |
---|
206 | CALL MPI_BCAST( fmax_ijk, 3, MPI_INTEGER, 0, comm2d, ierr ) |
---|
207 | #else |
---|
208 | fmax(1) = fmax_l(1) |
---|
209 | fmax_ijk = fmax_ijk_l |
---|
210 | #endif |
---|
211 | |
---|
212 | ENDIF |
---|
213 | |
---|
214 | ! |
---|
215 | !-- Determine output parameters |
---|
216 | SELECT CASE( mode ) |
---|
217 | |
---|
218 | CASE( 'min' ) |
---|
219 | |
---|
220 | value = fmin(1) |
---|
221 | value_ijk = fmin_ijk |
---|
222 | |
---|
223 | CASE( 'max' ) |
---|
224 | |
---|
225 | value = fmax(1) |
---|
226 | value_ijk = fmax_ijk |
---|
227 | |
---|
228 | CASE( 'minmax' ) |
---|
229 | |
---|
230 | value = fmin(1) |
---|
231 | value_ijk = fmin_ijk |
---|
232 | value1 = fmax(1) |
---|
233 | value1_ijk = fmax_ijk |
---|
234 | |
---|
235 | CASE( 'abs' ) |
---|
236 | |
---|
237 | value = fmax(1) |
---|
238 | value_ijk = fmax_ijk |
---|
239 | IF ( fmax_ijk(1) < 0 ) THEN |
---|
240 | value = -value |
---|
241 | value_ijk(1) = -value_ijk(1) - 10 |
---|
242 | ENDIF |
---|
243 | |
---|
244 | END SELECT |
---|
245 | |
---|
246 | ! |
---|
247 | !-- Limit index values to the range 0..nx, 0..ny |
---|
248 | IF ( value_ijk(3) == -1 ) value_ijk(3) = nx |
---|
249 | IF ( value_ijk(3) == nx+1 ) value_ijk(3) = 0 |
---|
250 | IF ( value_ijk(2) == -1 ) value_ijk(2) = ny |
---|
251 | IF ( value_ijk(2) == ny+1 ) value_ijk(2) = 0 |
---|
252 | |
---|
253 | |
---|
254 | END SUBROUTINE global_min_max |
---|