source: palm/trunk/SOURCE/global_min_max.f90 @ 671

Last change on this file since 671 was 668, checked in by suehring, 14 years ago

last commit documented

  • Property svn:keywords set to Id
File size: 7.0 KB
RevLine 
[1]1 SUBROUTINE global_min_max( i1, i2, j1, j2, k1, k2, ar, mode, value, &
2                            value_ijk, value1, value1_ijk )
3
[623]4!------------------------------------------------------------------------------!
[484]5! Current revisions:
[1]6! -----------------
7!
8! Former revisions:
9! -----------------
[3]10! $Id: global_min_max.f90 668 2010-12-23 13:22:58Z heinze $
[623]11!
[668]12! 667 2010-12-23 12:06:00Z suehring/gryschka
13! Adapting of the index arrays, because MINLOC assumes lowerbound at 1 and not
14! at nbgp.
15!
[623]16! 622 2010-12-10 08:08:13Z raasch
17! optional barriers included in order to speed up collective operations
18!
19! Feb. 2007
[3]20! RCS Log replace by Id keyword, revision history cleaned up
21!
[1]22! Revision 1.11  2003/04/16 12:56:58  raasch
23! Index values of the extrema are limited to the range 0..nx, 0..ny
24!
25! Revision 1.1  1997/07/24 11:14:03  raasch
26! Initial revision
27!
28!
29! Description:
30! ------------
31! Determine the array minimum/maximum and the corresponding indices.
[623]32!------------------------------------------------------------------------------!
[1]33
34    USE indices
35    USE pegrid
36
37    IMPLICIT NONE
38
39    CHARACTER (LEN=*) ::  mode
40
41    INTEGER           ::  i, i1, i2, id_fmax, id_fmin, j, j1, j2, k, k1, k2, &
42                          fmax_ijk(3), fmax_ijk_l(3), fmin_ijk(3), &
43                          fmin_ijk_l(3), value_ijk(3)
44    INTEGER, OPTIONAL ::  value1_ijk(3)
45    REAL              ::  value, &
46                          ar(i1:i2,j1:j2,k1:k2)
47#if defined( __ibm )
48    REAL (KIND=4)     ::  fmax(2), fmax_l(2), fmin(2), fmin_l(2)  ! on 32bit-
49                          ! machines MPI_2REAL must not be replaced by
50                          ! MPI_2DOUBLE_PRECISION
51#else
52    REAL              ::  fmax(2), fmax_l(2), fmin(2), fmin_l(2)
53#endif
54    REAL, OPTIONAL    ::  value1
55
56
57!
58!-- Determine array minimum
59    IF ( mode == 'min'  .OR.  mode == 'minmax' )  THEN
60
61!
62!--    Determine the local minimum
63       fmin_ijk_l = MINLOC( ar )
[667]64       fmin_ijk_l(1) = i1 + fmin_ijk_l(1) - nbgp    ! MINLOC assumes lowerbound = 1
65       fmin_ijk_l(2) = j1 + fmin_ijk_l(2) - nbgp
[1]66       fmin_ijk_l(3) = k1 + fmin_ijk_l(3) - 1
67       fmin_l(1)  = ar(fmin_ijk_l(1),fmin_ijk_l(2),fmin_ijk_l(3))
68
69#if defined( __parallel )
70       fmin_l(2)  = myid
[622]71       IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
[623]72       CALL MPI_ALLREDUCE( fmin_l, fmin, 1, MPI_2REAL, MPI_MINLOC, comm2d, &
73                           ierr )
[1]74
75!
76!--    Determine the global minimum. Result stored on PE0.
77       id_fmin = fmin(2)
78       IF ( id_fmin /= 0 )  THEN
79          IF ( myid == 0 )  THEN
80             CALL MPI_RECV( fmin_ijk, 3, MPI_INTEGER, id_fmin, 0, comm2d, &
81                            status, ierr )
82          ELSEIF ( myid == id_fmin )  THEN
83             CALL MPI_SEND( fmin_ijk_l, 3, MPI_INTEGER, 0, 0, comm2d, ierr )
84          ENDIF
85       ELSE
86          fmin_ijk = fmin_ijk_l
87       ENDIF
88!
89!--    Send the indices of the just determined array minimum to other PEs
90       CALL MPI_BCAST( fmin_ijk, 3, MPI_INTEGER, 0, comm2d, ierr )
91#else
92       fmin(1)  = fmin_l(1)
93       fmin_ijk = fmin_ijk_l
94#endif
95
96    ENDIF
97
98!
99!-- Determine array maximum
100    IF ( mode == 'max'  .OR.  mode == 'minmax' )  THEN
101
102!
103!--    Determine the local maximum
104       fmax_ijk_l = MAXLOC( ar )
[667]105       fmax_ijk_l(1) = i1 + fmax_ijk_l(1) - nbgp    ! MAXLOC assumes lowerbound = 1
106       fmax_ijk_l(2) = j1 + fmax_ijk_l(2) - nbgp
[1]107       fmax_ijk_l(3) = k1 + fmax_ijk_l(3) - 1
108       fmax_l(1) = ar(fmax_ijk_l(1),fmax_ijk_l(2),fmax_ijk_l(3))
109
110#if defined( __parallel )
111       fmax_l(2)  = myid
[622]112       IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
[623]113       CALL MPI_ALLREDUCE( fmax_l, fmax, 1, MPI_2REAL, MPI_MAXLOC, comm2d, &
114                           ierr )
[1]115
116!
117!--    Determine the global maximum. Result stored on PE0.
118       id_fmax = fmax(2)
119       IF ( id_fmax /= 0 )  THEN
120          IF ( myid == 0 )  THEN
121             CALL MPI_RECV( fmax_ijk, 3, MPI_INTEGER, id_fmax, 0, comm2d, &
122                            status, ierr )
123          ELSEIF ( myid == id_fmax )  THEN
124             CALL MPI_SEND( fmax_ijk_l, 3, MPI_INTEGER, 0, 0, comm2d, ierr )
125          ENDIF
126       ELSE
127          fmax_ijk = fmax_ijk_l
128       ENDIF
129!
130!--    send the indices of the just determined array maximum to other PEs
131       CALL MPI_BCAST( fmax_ijk, 3, MPI_INTEGER, 0, comm2d, ierr )
132#else
133       fmax(1)  = fmax_l(1)
134       fmax_ijk = fmax_ijk_l
135#endif
136
137    ENDIF
138
139!
140!-- Determine absolute array maximum
141    IF ( mode == 'abs' )  THEN
142
143!
144!--    Determine the local absolut maximum
145       fmax_l(1)     = 0.0
146       fmax_ijk_l(1) =  i1
147       fmax_ijk_l(2) =  j1
148       fmax_ijk_l(3) =  k1
149       DO  k = k1, k2
150          DO  j = j1, j2
151             DO  i = i1, i2
152                IF ( ABS( ar(i,j,k) ) > fmax_l(1) )  THEN
153                   fmax_l(1) = ABS( ar(i,j,k) )
154                   fmax_ijk_l(1) = i
155                   fmax_ijk_l(2) = j
156                   fmax_ijk_l(3) = k
157                ENDIF
158             ENDDO
159          ENDDO
160       ENDDO
161
162!
163!--    Set a flag in case that the determined value is negative.
164!--    A constant offset has to be subtracted in order to handle the special
165!--    case i=0 correctly
166       IF ( ar(fmax_ijk_l(1),fmax_ijk_l(2),fmax_ijk_l(3)) < 0.0 )  THEN
167          fmax_ijk_l(1) = -fmax_ijk_l(1) - 10
168       ENDIF
169
170#if defined( __parallel )
171       fmax_l(2)  = myid
[622]172       IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
[1]173       CALL MPI_ALLREDUCE( fmax_l, fmax, 1, MPI_2REAL, MPI_MAXLOC, comm2d, &
174                           ierr )
175
176!
177!--    Determine the global absolut maximum. Result stored on PE0.
178       id_fmax = fmax(2)
179       IF ( id_fmax /= 0 )  THEN
180          IF ( myid == 0 )  THEN
181             CALL MPI_RECV( fmax_ijk, 3, MPI_INTEGER, id_fmax, 0, comm2d, &
182                            status, ierr )
183          ELSEIF ( myid == id_fmax )  THEN
184             CALL MPI_SEND( fmax_ijk_l, 3, MPI_INTEGER, 0, 0, comm2d, ierr )
185          ENDIF
186       ELSE
187          fmax_ijk = fmax_ijk_l
188       ENDIF
189!
190!--    Send the indices of the just determined absolut maximum to other PEs
191       CALL MPI_BCAST( fmax_ijk, 3, MPI_INTEGER, 0, comm2d, ierr )
192#else
193       fmax(1)  = fmax_l(1)
194       fmax_ijk = fmax_ijk_l
195#endif
196
197    ENDIF
198
199!
200!-- Determine output parameters
201    SELECT CASE( mode )
202
203       CASE( 'min' )
204
205          value     = fmin(1)
206          value_ijk = fmin_ijk
207
208       CASE( 'max' )
209
210          value     = fmax(1)
211          value_ijk = fmax_ijk
212
213       CASE( 'minmax' )
214
215          value      = fmin(1)
216          value_ijk  = fmin_ijk
217          value1     = fmax(1)
218          value1_ijk = fmax_ijk
219
220       CASE( 'abs' )
221
222          value     = fmax(1)
223          value_ijk = fmax_ijk
224          IF ( fmax_ijk(1) < 0 )  THEN
225             value        = -value
[667]226             value_ijk(1) = -value_ijk(1) - 10         !???
[1]227          ENDIF
228
229    END SELECT
230
231!
232!-- Limit index values to the range 0..nx, 0..ny
[667]233    IF ( value_ijk(3) < 0  ) value_ijk(3) = nx +1 + value_ijk(3)
234    IF ( value_ijk(3) > nx ) value_ijk(3) = value_ijk(3) - (nx+1)
235    IF ( value_ijk(2) < 0  ) value_ijk(2) = ny +1 + value_ijk(2)
236    IF ( value_ijk(2) > ny ) value_ijk(2) = value_ijk(2) - (ny+1)
[1]237
238
239 END SUBROUTINE global_min_max
Note: See TracBrowser for help on using the repository browser.