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