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

Last change on this file since 4218 was 4182, checked in by scharf, 5 years ago
  • corrected "Former revisions" section
  • minor formatting in "Former revisions" section
  • added "Author" section
  • Property svn:keywords set to Id
File size: 10.2 KB
Line 
1!> @file global_min_max.f90
2!------------------------------------------------------------------------------!
3! This file is part of the PALM model system.
4!
5! PALM is free software: you can redistribute it and/or modify it under the
6! terms of the GNU General Public License as published by the Free Software
7! Foundation, either version 3 of the License, or (at your option) any later
8! version.
9!
10! PALM is distributed in the hope that it will be useful, but WITHOUT ANY
11! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR
12! A PARTICULAR PURPOSE.  See the GNU General Public License for more details.
13!
14! You should have received a copy of the GNU General Public License along with
15! PALM. If not, see <http://www.gnu.org/licenses/>.
16!
17! Copyright 1997-2019 Leibniz Universitaet Hannover
18!------------------------------------------------------------------------------!
19!
20! Current revisions:
21! ------------------
22!
23!
24! Former revisions:
25! -----------------
26! $Id: global_min_max.f90 4182 2019-08-22 15:20:23Z knoop $
27! Corrected "Former revisions" section
28!
29! 3655 2019-01-07 16:51:22Z knoop
30! Corrected "Former revisions" section
31!
32! Revision 1.1  1997/07/24 11:14:03  raasch
33! Initial revision
34!
35!
36! Description:
37! ------------
38!> Determine the array minimum/maximum and the corresponding indices.
39!------------------------------------------------------------------------------!
40 SUBROUTINE global_min_max( i1, i2, j1, j2, k1, k2, ar, mode, offset, value, &
41                            value_ijk, value1, value1_ijk )
42 
43
44    USE indices,                                                               &
45        ONLY:  nbgp, ny, nx
46       
47    USE kinds
48   
49    USE pegrid
50
51    IMPLICIT NONE
52
53    CHARACTER (LEN=*) ::  mode  !<
54
55    INTEGER(iwp) ::  i              !<
56    INTEGER(iwp) ::  i1             !<
57    INTEGER(iwp) ::  i2             !<
58    INTEGER(iwp) ::  id_fmax        !<
59    INTEGER(iwp) ::  id_fmin        !<
60    INTEGER(iwp) ::  j              !<
61    INTEGER(iwp) ::  j1             !<
62    INTEGER(iwp) ::  j2             !<
63    INTEGER(iwp) ::  k              !<
64    INTEGER(iwp) ::  k1             !<
65    INTEGER(iwp) ::  k2             !<
66    INTEGER(iwp) ::  fmax_ijk(3)    !<
67    INTEGER(iwp) ::  fmax_ijk_l(3)  !<
68    INTEGER(iwp) ::  fmin_ijk(3)    !<
69    INTEGER(iwp) ::  fmin_ijk_l(3)  !<
70    INTEGER(iwp) ::  value_ijk(3)   !<
71   
72    INTEGER(iwp), OPTIONAL ::  value1_ijk(3)  !<
73   
74    REAL(wp) ::  offset                 !<
75    REAL(wp) ::  value                  !<
76    REAL(wp) ::  ar(i1:i2,j1:j2,k1:k2)  !<
77   
78#if defined( __ibm )
79    REAL(sp) ::  fmax(2)    !<
80    REAL(sp) ::  fmax_l(2)  !<
81    REAL(sp) ::  fmin(2)    !<
82    REAL(sp) ::  fmin_l(2)  !<
83             ! on 32bit-machines MPI_2REAL must not be replaced
84             ! by MPI_2DOUBLE_PRECISION
85#else
86    REAL(wp) ::  fmax(2)    !<
87    REAL(wp) ::  fmax_l(2)  !<
88    REAL(wp) ::  fmin(2)    !<
89    REAL(wp) ::  fmin_l(2)  !<
90#endif
91    REAL(wp), OPTIONAL ::  value1  !<
92
93
94!
95!-- Determine array minimum
96    IF ( mode == 'min'  .OR.  mode == 'minmax' )  THEN
97
98!
99!--    Determine the local minimum
100       fmin_ijk_l = MINLOC( ar )
101       fmin_ijk_l(1) = i1 + fmin_ijk_l(1) - 1 ! MINLOC assumes lowerbound = 1
102       fmin_ijk_l(2) = j1 + fmin_ijk_l(2) - nbgp
103       fmin_ijk_l(3) = k1 + fmin_ijk_l(3) - nbgp
104       fmin_l(1)  = ar(fmin_ijk_l(1),fmin_ijk_l(2),fmin_ijk_l(3))
105
106#if defined( __parallel )
107       fmin_l(2)  = myid
108       IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
109       CALL MPI_ALLREDUCE( fmin_l, fmin, 1, MPI_2REAL, MPI_MINLOC, comm2d, &
110                           ierr )
111
112!
113!--    Determine the global minimum. Result stored on PE0.
114       id_fmin = fmin(2)
115       IF ( id_fmin /= 0 )  THEN
116          IF ( myid == 0 )  THEN
117             CALL MPI_RECV( fmin_ijk, 3, MPI_INTEGER, id_fmin, 0, comm2d, &
118                            status, ierr )
119          ELSEIF ( myid == id_fmin )  THEN
120             CALL MPI_SEND( fmin_ijk_l, 3, MPI_INTEGER, 0, 0, comm2d, ierr )
121          ENDIF
122       ELSE
123          fmin_ijk = fmin_ijk_l
124       ENDIF
125!
126!--    Send the indices of the just determined array minimum to other PEs
127       CALL MPI_BCAST( fmin_ijk, 3, MPI_INTEGER, 0, comm2d, ierr )
128#else
129       fmin(1)  = fmin_l(1)
130       fmin_ijk = fmin_ijk_l
131#endif
132
133    ENDIF
134
135!
136!-- Determine array maximum
137    IF ( mode == 'max'  .OR.  mode == 'minmax' )  THEN
138
139!
140!--    Determine the local maximum
141       fmax_ijk_l = MAXLOC( ar )
142       fmax_ijk_l(1) = i1 + fmax_ijk_l(1) - 1 ! MAXLOC assumes lowerbound = 1
143       fmax_ijk_l(2) = j1 + fmax_ijk_l(2) - nbgp
144       fmax_ijk_l(3) = k1 + fmax_ijk_l(3) - nbgp
145       fmax_l(1) = ar(fmax_ijk_l(1),fmax_ijk_l(2),fmax_ijk_l(3))
146
147#if defined( __parallel )
148       fmax_l(2)  = myid
149       IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
150       CALL MPI_ALLREDUCE( fmax_l, fmax, 1, MPI_2REAL, MPI_MAXLOC, comm2d, &
151                           ierr )
152
153!
154!--    Determine the global maximum. Result stored on PE0.
155       id_fmax = fmax(2)
156       IF ( id_fmax /= 0 )  THEN
157          IF ( myid == 0 )  THEN
158             CALL MPI_RECV( fmax_ijk, 3, MPI_INTEGER, id_fmax, 0, comm2d, &
159                            status, ierr )
160          ELSEIF ( myid == id_fmax )  THEN
161             CALL MPI_SEND( fmax_ijk_l, 3, MPI_INTEGER, 0, 0, comm2d, ierr )
162          ENDIF
163       ELSE
164          fmax_ijk = fmax_ijk_l
165       ENDIF
166!
167!--    send the indices of the just determined array maximum to other PEs
168       CALL MPI_BCAST( fmax_ijk, 3, MPI_INTEGER, 0, comm2d, ierr )
169#else
170       fmax(1)  = fmax_l(1)
171       fmax_ijk = fmax_ijk_l
172#endif
173
174    ENDIF
175
176!
177!-- Determine absolute array maximum
178    IF ( mode == 'abs' )  THEN
179
180!
181!--    Determine the local absolut maximum
182       fmax_l(1)     = 0.0_wp
183       fmax_ijk_l(1) =  i1
184       fmax_ijk_l(2) =  j1
185       fmax_ijk_l(3) =  k1
186       DO  k = k1, k2
187          DO  j = j1, j2
188             DO  i = i1, i2
189                IF ( ABS( ar(i,j,k) ) > fmax_l(1) )  THEN
190                   fmax_l(1) = ABS( ar(i,j,k) )
191                   fmax_ijk_l(1) = i
192                   fmax_ijk_l(2) = j
193                   fmax_ijk_l(3) = k
194                ENDIF
195             ENDDO
196          ENDDO
197       ENDDO
198
199!
200!--    Set a flag in case that the determined value is negative.
201!--    A constant offset has to be subtracted in order to handle the special
202!--    case i=0 correctly
203       IF ( ar(fmax_ijk_l(1),fmax_ijk_l(2),fmax_ijk_l(3)) < 0.0_wp )  THEN
204          fmax_ijk_l(1) = -fmax_ijk_l(1) - 10
205       ENDIF
206
207#if defined( __parallel )
208       fmax_l(2)  = myid
209       IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
210       CALL MPI_ALLREDUCE( fmax_l, fmax, 1, MPI_2REAL, MPI_MAXLOC, comm2d, &
211                           ierr )
212
213!
214!--    Determine the global absolut maximum. Result stored on PE0.
215       id_fmax = fmax(2)
216       IF ( id_fmax /= 0 )  THEN
217          IF ( myid == 0 )  THEN
218             CALL MPI_RECV( fmax_ijk, 3, MPI_INTEGER, id_fmax, 0, comm2d, &
219                            status, ierr )
220          ELSEIF ( myid == id_fmax )  THEN
221             CALL MPI_SEND( fmax_ijk_l, 3, MPI_INTEGER, 0, 0, comm2d, ierr )
222          ENDIF
223       ELSE
224          fmax_ijk = fmax_ijk_l
225       ENDIF
226!
227!--    Send the indices of the just determined absolut maximum to other PEs
228       CALL MPI_BCAST( fmax_ijk, 3, MPI_INTEGER, 0, comm2d, ierr )
229#else
230       fmax(1)  = fmax_l(1)
231       fmax_ijk = fmax_ijk_l
232#endif
233
234    ENDIF
235
236!
237!-- Determine absolute maximum of ( array - offset )
238    IF ( mode == 'absoff' )  THEN
239
240!
241!--    Determine the local absolut maximum
242       fmax_l(1)     = 0.0_wp
243       fmax_ijk_l(1) =  i1
244       fmax_ijk_l(2) =  j1
245       fmax_ijk_l(3) =  k1
246       DO  k = k1, k2
247          DO  j = j1, j2
248!
249!--          Attention: the lowest gridpoint is excluded here, because there
250!--          ---------  is no advection at nzb=0 and mode 'absoff' is only
251!--                     used for calculating u,v extrema for CFL-criteria
252             DO  i = i1+1, i2
253                IF ( ABS( ar(i,j,k) - offset ) > fmax_l(1) )  THEN
254                   fmax_l(1) = ABS( ar(i,j,k) - offset )
255                   fmax_ijk_l(1) = i
256                   fmax_ijk_l(2) = j
257                   fmax_ijk_l(3) = k
258                ENDIF
259             ENDDO
260          ENDDO
261       ENDDO
262
263!
264!--    Set a flag in case that the determined value is negative.
265!--    A constant offset has to be subtracted in order to handle the special
266!--    case i=0 correctly
267       IF ( ar(fmax_ijk_l(1),fmax_ijk_l(2),fmax_ijk_l(3)) < 0.0_wp )  THEN
268          fmax_ijk_l(1) = -fmax_ijk_l(1) - 10
269       ENDIF
270
271#if defined( __parallel )
272       fmax_l(2)  = myid
273       IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
274       CALL MPI_ALLREDUCE( fmax_l, fmax, 1, MPI_2REAL, MPI_MAXLOC, comm2d, &
275                           ierr )
276
277!
278!--    Determine the global absolut maximum. Result stored on PE0.
279       id_fmax = fmax(2)
280       IF ( id_fmax /= 0 )  THEN
281          IF ( myid == 0 )  THEN
282             CALL MPI_RECV( fmax_ijk, 3, MPI_INTEGER, id_fmax, 0, comm2d, &
283                            status, ierr )
284          ELSEIF ( myid == id_fmax )  THEN
285             CALL MPI_SEND( fmax_ijk_l, 3, MPI_INTEGER, 0, 0, comm2d, ierr )
286          ENDIF
287       ELSE
288          fmax_ijk = fmax_ijk_l
289       ENDIF
290!
291!--    Send the indices of the just determined absolut maximum to other PEs
292       CALL MPI_BCAST( fmax_ijk, 3, MPI_INTEGER, 0, comm2d, ierr )
293#else
294       fmax(1)  = fmax_l(1)
295       fmax_ijk = fmax_ijk_l
296#endif
297
298    ENDIF
299
300!
301!-- Determine output parameters
302    SELECT CASE( mode )
303
304       CASE( 'min' )
305
306          value     = fmin(1)
307          value_ijk = fmin_ijk
308
309       CASE( 'max' )
310
311          value     = fmax(1)
312          value_ijk = fmax_ijk
313
314       CASE( 'minmax' )
315
316          value      = fmin(1)
317          value_ijk  = fmin_ijk
318          value1     = fmax(1)
319          value1_ijk = fmax_ijk
320
321       CASE( 'abs', 'absoff' )
322
323          value     = fmax(1)
324          value_ijk = fmax_ijk
325          IF ( fmax_ijk(1) < 0 )  THEN
326             value        = -value
327             value_ijk(1) = -value_ijk(1) - 10         !???
328          ENDIF
329
330    END SELECT
331
332!
333!-- Limit index values to the range 0..nx, 0..ny
334    IF ( value_ijk(3) < 0  ) value_ijk(3) = nx +1 + value_ijk(3)
335    IF ( value_ijk(3) > nx ) value_ijk(3) = value_ijk(3) - (nx+1)
336    IF ( value_ijk(2) < 0  ) value_ijk(2) = ny +1 + value_ijk(2)
337    IF ( value_ijk(2) > ny ) value_ijk(2) = value_ijk(2) - (ny+1)
338
339
340 END SUBROUTINE global_min_max
Note: See TracBrowser for help on using the repository browser.