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

Last change on this file since 622 was 622, checked in by raasch, 13 years ago

New:
---

Optional barriers included in order to speed up collective operations
MPI_ALLTOALL and MPI_ALLREDUCE. This feature is controlled with new initial
parameter collective_wait. Default is .FALSE, but .TRUE. on SGI-type
systems. (advec_particles, advec_s_bc, buoyancy, check_for_restart,
cpu_statistics, data_output_2d, data_output_ptseries, flow_statistics,
global_min_max, inflow_turbulence, init_3d_model, init_particles, init_pegrid,
init_slope, parin, pres, poismg, set_particle_attributes, timestep,
read_var_list, user_statistics, write_compressed, write_var_list)

Adjustments for Kyushu Univ. (lcrte, ibmku). Concerning hybrid
(MPI/openMP) runs, the number of openMP threads per MPI tasks can now
be given as an argument to mrun-option -O. (mbuild, mrun, subjob)

Changed:


Initialization of the module command changed for SGI-ICE/lcsgi (mbuild, subjob)

Errors:


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