[1] | 1 | SUBROUTINE init_pegrid |
---|
| 2 | |
---|
| 3 | !------------------------------------------------------------------------------! |
---|
| 4 | ! Actual revisions: |
---|
| 5 | ! ----------------- |
---|
[145] | 6 | ! Collect on PE0 horizontal index bounds from all other PEs |
---|
[108] | 7 | ! TEST OUTPUT (TO BE REMOVED) logging mpi2 ierr values |
---|
[77] | 8 | ! |
---|
| 9 | ! Former revisions: |
---|
| 10 | ! ----------------- |
---|
| 11 | ! $Id: init_pegrid.f90 145 2008-01-09 08:17:38Z letzel $ |
---|
| 12 | ! |
---|
[139] | 13 | ! 114 2007-10-10 00:03:15Z raasch |
---|
| 14 | ! Allocation of wall flag arrays for multigrid solver |
---|
| 15 | ! |
---|
[110] | 16 | ! 108 2007-08-24 15:10:38Z letzel |
---|
| 17 | ! Intercommunicator (comm_inter) and derived data type (type_xy) for |
---|
| 18 | ! coupled model runs created, assign coupling_mode_remote, |
---|
| 19 | ! indices nxlu and nysv are calculated (needed for non-cyclic boundary |
---|
| 20 | ! conditions) |
---|
| 21 | ! |
---|
[83] | 22 | ! 82 2007-04-16 15:40:52Z raasch |
---|
| 23 | ! Cpp-directive lcmuk changed to intel_openmp_bug, setting of host on lcmuk by |
---|
| 24 | ! cpp-directive removed |
---|
| 25 | ! |
---|
[77] | 26 | ! 75 2007-03-22 09:54:05Z raasch |
---|
[73] | 27 | ! uxrp, vynp eliminated, |
---|
[75] | 28 | ! dirichlet/neumann changed to dirichlet/radiation, etc., |
---|
| 29 | ! poisfft_init is only called if fft-solver is switched on |
---|
[1] | 30 | ! |
---|
[3] | 31 | ! RCS Log replace by Id keyword, revision history cleaned up |
---|
| 32 | ! |
---|
[1] | 33 | ! Revision 1.28 2006/04/26 13:23:32 raasch |
---|
| 34 | ! lcmuk does not understand the !$ comment so a cpp-directive is required |
---|
| 35 | ! |
---|
| 36 | ! Revision 1.1 1997/07/24 11:15:09 raasch |
---|
| 37 | ! Initial revision |
---|
| 38 | ! |
---|
| 39 | ! |
---|
| 40 | ! Description: |
---|
| 41 | ! ------------ |
---|
| 42 | ! Determination of the virtual processor topology (if not prescribed by the |
---|
| 43 | ! user)and computation of the grid point number and array bounds of the local |
---|
| 44 | ! domains. |
---|
| 45 | !------------------------------------------------------------------------------! |
---|
| 46 | |
---|
| 47 | USE control_parameters |
---|
| 48 | USE fft_xy |
---|
| 49 | USE indices |
---|
| 50 | USE pegrid |
---|
| 51 | USE poisfft_mod |
---|
| 52 | USE poisfft_hybrid_mod |
---|
| 53 | USE statistics |
---|
| 54 | USE transpose_indices |
---|
| 55 | |
---|
| 56 | |
---|
| 57 | IMPLICIT NONE |
---|
| 58 | |
---|
| 59 | INTEGER :: gathered_size, i, ind(5), j, k, maximum_grid_level_l, & |
---|
| 60 | mg_switch_to_pe0_level_l, mg_levels_x, mg_levels_y, & |
---|
| 61 | mg_levels_z, nnx_y, nnx_z, nny_x, nny_z, nnz_x, nnz_y, & |
---|
| 62 | numproc_sqr, nx_total, nxl_l, nxr_l, nyn_l, nys_l, nzb_l, & |
---|
| 63 | nzt_l, omp_get_num_threads, subdomain_size |
---|
| 64 | |
---|
| 65 | INTEGER, DIMENSION(:), ALLOCATABLE :: ind_all, nxlf, nxrf, nynf, nysf |
---|
| 66 | |
---|
| 67 | LOGICAL :: found |
---|
| 68 | |
---|
| 69 | ! |
---|
| 70 | !-- Get the number of OpenMP threads |
---|
| 71 | !$OMP PARALLEL |
---|
[82] | 72 | #if defined( __intel_openmp_bug ) |
---|
[1] | 73 | threads_per_task = omp_get_num_threads() |
---|
| 74 | #else |
---|
| 75 | !$ threads_per_task = omp_get_num_threads() |
---|
| 76 | #endif |
---|
| 77 | !$OMP END PARALLEL |
---|
| 78 | |
---|
| 79 | |
---|
| 80 | #if defined( __parallel ) |
---|
| 81 | ! |
---|
| 82 | !-- Determine the processor topology or check it, if prescribed by the user |
---|
| 83 | IF ( npex == -1 .AND. npey == -1 ) THEN |
---|
| 84 | |
---|
| 85 | ! |
---|
| 86 | !-- Automatic determination of the topology |
---|
| 87 | !-- The default on SMP- and cluster-hosts is a 1d-decomposition along x |
---|
| 88 | IF ( host(1:3) == 'ibm' .OR. host(1:3) == 'nec' .OR. & |
---|
| 89 | host(1:2) == 'lc' .OR. host(1:3) == 'dec' ) THEN |
---|
| 90 | |
---|
| 91 | pdims(1) = numprocs |
---|
| 92 | pdims(2) = 1 |
---|
| 93 | |
---|
| 94 | ELSE |
---|
| 95 | |
---|
| 96 | numproc_sqr = SQRT( REAL( numprocs ) ) |
---|
| 97 | pdims(1) = MAX( numproc_sqr , 1 ) |
---|
| 98 | DO WHILE ( MOD( numprocs , pdims(1) ) /= 0 ) |
---|
| 99 | pdims(1) = pdims(1) - 1 |
---|
| 100 | ENDDO |
---|
| 101 | pdims(2) = numprocs / pdims(1) |
---|
| 102 | |
---|
| 103 | ENDIF |
---|
| 104 | |
---|
| 105 | ELSEIF ( npex /= -1 .AND. npey /= -1 ) THEN |
---|
| 106 | |
---|
| 107 | ! |
---|
| 108 | !-- Prescribed by user. Number of processors on the prescribed topology |
---|
| 109 | !-- must be equal to the number of PEs available to the job |
---|
| 110 | IF ( ( npex * npey ) /= numprocs ) THEN |
---|
| 111 | PRINT*, '+++ init_pegrid:' |
---|
| 112 | PRINT*, ' number of PEs of the prescribed topology (', npex*npey, & |
---|
| 113 | ') does not match the number of PEs available to the ', & |
---|
| 114 | 'job (', numprocs, ')' |
---|
| 115 | CALL local_stop |
---|
| 116 | ENDIF |
---|
| 117 | pdims(1) = npex |
---|
| 118 | pdims(2) = npey |
---|
| 119 | |
---|
| 120 | ELSE |
---|
| 121 | ! |
---|
| 122 | !-- If the processor topology is prescribed by the user, the number of |
---|
| 123 | !-- PEs must be given in both directions |
---|
| 124 | PRINT*, '+++ init_pegrid:' |
---|
| 125 | PRINT*, ' if the processor topology is prescribed by the user, ', & |
---|
| 126 | 'both values of "npex" and "npey" must be given in the ', & |
---|
| 127 | 'NAMELIST-parameter file' |
---|
| 128 | CALL local_stop |
---|
| 129 | |
---|
| 130 | ENDIF |
---|
| 131 | |
---|
| 132 | ! |
---|
| 133 | !-- The hybrid solver can only be used in case of a 1d-decomposition along x |
---|
| 134 | IF ( pdims(2) /= 1 .AND. psolver == 'poisfft_hybrid' ) THEN |
---|
| 135 | IF ( myid == 0 ) THEN |
---|
| 136 | PRINT*, '*** init_pegrid: psolver = "poisfft_hybrid" can only be' |
---|
| 137 | PRINT*, ' used in case of a 1d-decomposition along x' |
---|
| 138 | ENDIF |
---|
| 139 | ENDIF |
---|
| 140 | |
---|
| 141 | ! |
---|
| 142 | !-- If necessary, set horizontal boundary conditions to non-cyclic |
---|
| 143 | IF ( bc_lr /= 'cyclic' ) cyclic(1) = .FALSE. |
---|
| 144 | IF ( bc_ns /= 'cyclic' ) cyclic(2) = .FALSE. |
---|
| 145 | |
---|
| 146 | ! |
---|
| 147 | !-- Create the virtual processor grid |
---|
| 148 | CALL MPI_CART_CREATE( comm_palm, ndim, pdims, cyclic, reorder, & |
---|
| 149 | comm2d, ierr ) |
---|
| 150 | CALL MPI_COMM_RANK( comm2d, myid, ierr ) |
---|
| 151 | WRITE (myid_char,'(''_'',I4.4)') myid |
---|
| 152 | |
---|
| 153 | CALL MPI_CART_COORDS( comm2d, myid, ndim, pcoord, ierr ) |
---|
| 154 | CALL MPI_CART_SHIFT( comm2d, 0, 1, pleft, pright, ierr ) |
---|
| 155 | CALL MPI_CART_SHIFT( comm2d, 1, 1, psouth, pnorth, ierr ) |
---|
| 156 | |
---|
| 157 | ! |
---|
| 158 | !-- Determine sub-topologies for transpositions |
---|
| 159 | !-- Transposition from z to x: |
---|
| 160 | remain_dims(1) = .TRUE. |
---|
| 161 | remain_dims(2) = .FALSE. |
---|
| 162 | CALL MPI_CART_SUB( comm2d, remain_dims, comm1dx, ierr ) |
---|
| 163 | CALL MPI_COMM_RANK( comm1dx, myidx, ierr ) |
---|
| 164 | ! |
---|
| 165 | !-- Transposition from x to y |
---|
| 166 | remain_dims(1) = .FALSE. |
---|
| 167 | remain_dims(2) = .TRUE. |
---|
| 168 | CALL MPI_CART_SUB( comm2d, remain_dims, comm1dy, ierr ) |
---|
| 169 | CALL MPI_COMM_RANK( comm1dy, myidy, ierr ) |
---|
| 170 | |
---|
| 171 | |
---|
| 172 | ! |
---|
| 173 | !-- Find a grid (used for array d) which will match the transposition demands |
---|
| 174 | IF ( grid_matching == 'strict' ) THEN |
---|
| 175 | |
---|
| 176 | nxa = nx; nya = ny; nza = nz |
---|
| 177 | |
---|
| 178 | ELSE |
---|
| 179 | |
---|
| 180 | found = .FALSE. |
---|
| 181 | xn: DO nxa = nx, 2*nx |
---|
| 182 | ! |
---|
| 183 | !-- Meet conditions for nx |
---|
| 184 | IF ( MOD( nxa+1, pdims(1) ) /= 0 .OR. & |
---|
| 185 | MOD( nxa+1, pdims(2) ) /= 0 ) CYCLE xn |
---|
| 186 | |
---|
| 187 | yn: DO nya = ny, 2*ny |
---|
| 188 | ! |
---|
| 189 | !-- Meet conditions for ny |
---|
| 190 | IF ( MOD( nya+1, pdims(2) ) /= 0 .OR. & |
---|
| 191 | MOD( nya+1, pdims(1) ) /= 0 ) CYCLE yn |
---|
| 192 | |
---|
| 193 | |
---|
| 194 | zn: DO nza = nz, 2*nz |
---|
| 195 | ! |
---|
| 196 | !-- Meet conditions for nz |
---|
| 197 | IF ( ( MOD( nza, pdims(1) ) /= 0 .AND. pdims(1) /= 1 .AND. & |
---|
| 198 | pdims(2) /= 1 ) .OR. & |
---|
| 199 | ( MOD( nza, pdims(2) ) /= 0 .AND. dt_dosp /= 9999999.9 & |
---|
| 200 | ) ) THEN |
---|
| 201 | CYCLE zn |
---|
| 202 | ELSE |
---|
| 203 | found = .TRUE. |
---|
| 204 | EXIT xn |
---|
| 205 | ENDIF |
---|
| 206 | |
---|
| 207 | ENDDO zn |
---|
| 208 | |
---|
| 209 | ENDDO yn |
---|
| 210 | |
---|
| 211 | ENDDO xn |
---|
| 212 | |
---|
| 213 | IF ( .NOT. found ) THEN |
---|
| 214 | IF ( myid == 0 ) THEN |
---|
| 215 | PRINT*,'+++ init_pegrid: no matching grid for transpositions found' |
---|
| 216 | ENDIF |
---|
| 217 | CALL local_stop |
---|
| 218 | ENDIF |
---|
| 219 | |
---|
| 220 | ENDIF |
---|
| 221 | |
---|
| 222 | ! |
---|
| 223 | !-- Calculate array bounds in x-direction for every PE. |
---|
| 224 | !-- The last PE along x may get less grid points than the others |
---|
| 225 | ALLOCATE( nxlf(0:pdims(1)-1), nxrf(0:pdims(1)-1), nynf(0:pdims(2)-1), & |
---|
| 226 | nysf(0:pdims(2)-1), nnx_pe(0:pdims(1)-1), nny_pe(0:pdims(2)-1) ) |
---|
| 227 | |
---|
| 228 | IF ( MOD( nxa+1 , pdims(1) ) /= 0 ) THEN |
---|
| 229 | IF ( myid == 0 ) THEN |
---|
| 230 | PRINT*,'+++ x-direction: gridpoint number (',nx+1,') is not an' |
---|
| 231 | PRINT*,' integral divisor of the number of proces', & |
---|
| 232 | &'sors (', pdims(1),')' |
---|
| 233 | ENDIF |
---|
| 234 | CALL local_stop |
---|
| 235 | ELSE |
---|
| 236 | nnx = ( nxa + 1 ) / pdims(1) |
---|
| 237 | IF ( nnx*pdims(1) - ( nx + 1) > nnx ) THEN |
---|
| 238 | IF ( myid == 0 ) THEN |
---|
| 239 | PRINT*,'+++ x-direction: nx does not match the requirements ', & |
---|
| 240 | 'given by the number of PEs' |
---|
| 241 | PRINT*,' used' |
---|
| 242 | PRINT*,' please use nx = ', nx - ( pdims(1) - ( nnx*pdims(1) & |
---|
| 243 | - ( nx + 1 ) ) ), ' instead of nx =', nx |
---|
| 244 | ENDIF |
---|
| 245 | CALL local_stop |
---|
| 246 | ENDIF |
---|
| 247 | ENDIF |
---|
| 248 | |
---|
| 249 | ! |
---|
| 250 | !-- Left and right array bounds, number of gridpoints |
---|
| 251 | DO i = 0, pdims(1)-1 |
---|
| 252 | nxlf(i) = i * nnx |
---|
| 253 | nxrf(i) = ( i + 1 ) * nnx - 1 |
---|
| 254 | nnx_pe(i) = MIN( nx, nxrf(i) ) - nxlf(i) + 1 |
---|
| 255 | ENDDO |
---|
| 256 | |
---|
| 257 | ! |
---|
| 258 | !-- Calculate array bounds in y-direction for every PE. |
---|
| 259 | IF ( MOD( nya+1 , pdims(2) ) /= 0 ) THEN |
---|
| 260 | IF ( myid == 0 ) THEN |
---|
| 261 | PRINT*,'+++ y-direction: gridpoint number (',ny+1,') is not an' |
---|
| 262 | PRINT*,' integral divisor of the number of proces', & |
---|
| 263 | &'sors (', pdims(2),')' |
---|
| 264 | ENDIF |
---|
| 265 | CALL local_stop |
---|
| 266 | ELSE |
---|
| 267 | nny = ( nya + 1 ) / pdims(2) |
---|
| 268 | IF ( nny*pdims(2) - ( ny + 1) > nny ) THEN |
---|
| 269 | IF ( myid == 0 ) THEN |
---|
| 270 | PRINT*,'+++ x-direction: nx does not match the requirements ', & |
---|
| 271 | 'given by the number of PEs' |
---|
| 272 | PRINT*,' used' |
---|
| 273 | PRINT*,' please use nx = ', nx - ( pdims(1) - ( nnx*pdims(1) & |
---|
| 274 | - ( nx + 1 ) ) ), ' instead of nx =', nx |
---|
| 275 | ENDIF |
---|
| 276 | CALL local_stop |
---|
| 277 | ENDIF |
---|
| 278 | ENDIF |
---|
| 279 | |
---|
| 280 | ! |
---|
| 281 | !-- South and north array bounds |
---|
| 282 | DO j = 0, pdims(2)-1 |
---|
| 283 | nysf(j) = j * nny |
---|
| 284 | nynf(j) = ( j + 1 ) * nny - 1 |
---|
| 285 | nny_pe(j) = MIN( ny, nynf(j) ) - nysf(j) + 1 |
---|
| 286 | ENDDO |
---|
| 287 | |
---|
| 288 | ! |
---|
| 289 | !-- Local array bounds of the respective PEs |
---|
| 290 | nxl = nxlf(pcoord(1)) |
---|
| 291 | nxra = nxrf(pcoord(1)) |
---|
| 292 | nxr = MIN( nx, nxra ) |
---|
| 293 | nys = nysf(pcoord(2)) |
---|
| 294 | nyna = nynf(pcoord(2)) |
---|
| 295 | nyn = MIN( ny, nyna ) |
---|
| 296 | nzb = 0 |
---|
| 297 | nzta = nza |
---|
| 298 | nzt = MIN( nz, nzta ) |
---|
| 299 | nnz = nza |
---|
| 300 | |
---|
| 301 | ! |
---|
| 302 | !-- Calculate array bounds and gridpoint numbers for the transposed arrays |
---|
| 303 | !-- (needed in the pressure solver) |
---|
| 304 | !-- For the transposed arrays, cyclic boundaries as well as top and bottom |
---|
| 305 | !-- boundaries are omitted, because they are obstructive to the transposition |
---|
| 306 | |
---|
| 307 | ! |
---|
| 308 | !-- 1. transposition z --> x |
---|
| 309 | !-- This transposition is not neccessary in case of a 1d-decomposition along x, |
---|
| 310 | !-- except that the uptream-spline method is switched on |
---|
| 311 | IF ( pdims(2) /= 1 .OR. momentum_advec == 'ups-scheme' .OR. & |
---|
| 312 | scalar_advec == 'ups-scheme' ) THEN |
---|
| 313 | |
---|
| 314 | IF ( pdims(2) == 1 .AND. ( momentum_advec == 'ups-scheme' .OR. & |
---|
| 315 | scalar_advec == 'ups-scheme' ) ) THEN |
---|
| 316 | IF ( myid == 0 ) THEN |
---|
| 317 | PRINT*,'+++ WARNING: init_pegrid: 1d-decomposition along x ', & |
---|
| 318 | &'chosen but nz restrictions may occur' |
---|
| 319 | PRINT*,' since ups-scheme is activated' |
---|
| 320 | ENDIF |
---|
| 321 | ENDIF |
---|
| 322 | nys_x = nys |
---|
| 323 | nyn_xa = nyna |
---|
| 324 | nyn_x = nyn |
---|
| 325 | nny_x = nny |
---|
| 326 | IF ( MOD( nza , pdims(1) ) /= 0 ) THEN |
---|
| 327 | IF ( myid == 0 ) THEN |
---|
| 328 | PRINT*,'+++ transposition z --> x:' |
---|
| 329 | PRINT*,' nz=',nz,' is not an integral divisior of pdims(1)=', & |
---|
| 330 | &pdims(1) |
---|
| 331 | ENDIF |
---|
| 332 | CALL local_stop |
---|
| 333 | ENDIF |
---|
| 334 | nnz_x = nza / pdims(1) |
---|
| 335 | nzb_x = 1 + myidx * nnz_x |
---|
| 336 | nzt_xa = ( myidx + 1 ) * nnz_x |
---|
| 337 | nzt_x = MIN( nzt, nzt_xa ) |
---|
| 338 | |
---|
| 339 | sendrecvcount_zx = nnx * nny * nnz_x |
---|
| 340 | |
---|
| 341 | ENDIF |
---|
| 342 | |
---|
| 343 | ! |
---|
| 344 | !-- 2. transposition x --> y |
---|
| 345 | nnz_y = nnz_x |
---|
| 346 | nzb_y = nzb_x |
---|
| 347 | nzt_ya = nzt_xa |
---|
| 348 | nzt_y = nzt_x |
---|
| 349 | IF ( MOD( nxa+1 , pdims(2) ) /= 0 ) THEN |
---|
| 350 | IF ( myid == 0 ) THEN |
---|
| 351 | PRINT*,'+++ transposition x --> y:' |
---|
| 352 | PRINT*,' nx+1=',nx+1,' is not an integral divisor of ',& |
---|
| 353 | &'pdims(2)=',pdims(2) |
---|
| 354 | ENDIF |
---|
| 355 | CALL local_stop |
---|
| 356 | ENDIF |
---|
| 357 | nnx_y = (nxa+1) / pdims(2) |
---|
| 358 | nxl_y = myidy * nnx_y |
---|
| 359 | nxr_ya = ( myidy + 1 ) * nnx_y - 1 |
---|
| 360 | nxr_y = MIN( nx, nxr_ya ) |
---|
| 361 | |
---|
| 362 | sendrecvcount_xy = nnx_y * nny_x * nnz_y |
---|
| 363 | |
---|
| 364 | ! |
---|
| 365 | !-- 3. transposition y --> z (ELSE: x --> y in case of 1D-decomposition |
---|
| 366 | !-- along x) |
---|
| 367 | IF ( pdims(2) /= 1 .OR. momentum_advec == 'ups-scheme' .OR. & |
---|
| 368 | scalar_advec == 'ups-scheme' ) THEN |
---|
| 369 | ! |
---|
| 370 | !-- y --> z |
---|
| 371 | !-- This transposition is not neccessary in case of a 1d-decomposition |
---|
| 372 | !-- along x, except that the uptream-spline method is switched on |
---|
| 373 | nnx_z = nnx_y |
---|
| 374 | nxl_z = nxl_y |
---|
| 375 | nxr_za = nxr_ya |
---|
| 376 | nxr_z = nxr_y |
---|
| 377 | IF ( MOD( nya+1 , pdims(1) ) /= 0 ) THEN |
---|
| 378 | IF ( myid == 0 ) THEN |
---|
| 379 | PRINT*,'+++ Transposition y --> z:' |
---|
| 380 | PRINT*,' ny+1=',ny+1,' is not an integral divisor of ',& |
---|
| 381 | &'pdims(1)=',pdims(1) |
---|
| 382 | ENDIF |
---|
| 383 | CALL local_stop |
---|
| 384 | ENDIF |
---|
| 385 | nny_z = (nya+1) / pdims(1) |
---|
| 386 | nys_z = myidx * nny_z |
---|
| 387 | nyn_za = ( myidx + 1 ) * nny_z - 1 |
---|
| 388 | nyn_z = MIN( ny, nyn_za ) |
---|
| 389 | |
---|
| 390 | sendrecvcount_yz = nnx_y * nny_z * nnz_y |
---|
| 391 | |
---|
| 392 | ELSE |
---|
| 393 | ! |
---|
| 394 | !-- x --> y. This condition must be fulfilled for a 1D-decomposition along x |
---|
| 395 | IF ( MOD( nya+1 , pdims(1) ) /= 0 ) THEN |
---|
| 396 | IF ( myid == 0 ) THEN |
---|
| 397 | PRINT*,'+++ Transposition x --> y:' |
---|
| 398 | PRINT*,' ny+1=',ny+1,' is not an integral divisor of ',& |
---|
| 399 | &'pdims(1)=',pdims(1) |
---|
| 400 | ENDIF |
---|
| 401 | CALL local_stop |
---|
| 402 | ENDIF |
---|
| 403 | |
---|
| 404 | ENDIF |
---|
| 405 | |
---|
| 406 | ! |
---|
| 407 | !-- Indices for direct transpositions z --> y (used for calculating spectra) |
---|
| 408 | IF ( dt_dosp /= 9999999.9 ) THEN |
---|
| 409 | IF ( MOD( nza, pdims(2) ) /= 0 ) THEN |
---|
| 410 | IF ( myid == 0 ) THEN |
---|
| 411 | PRINT*,'+++ Direct transposition z --> y (needed for spectra):' |
---|
| 412 | PRINT*,' nz=',nz,' is not an integral divisor of ',& |
---|
| 413 | &'pdims(2)=',pdims(2) |
---|
| 414 | ENDIF |
---|
| 415 | CALL local_stop |
---|
| 416 | ELSE |
---|
| 417 | nxl_yd = nxl |
---|
| 418 | nxr_yda = nxra |
---|
| 419 | nxr_yd = nxr |
---|
| 420 | nzb_yd = 1 + myidy * ( nza / pdims(2) ) |
---|
| 421 | nzt_yda = ( myidy + 1 ) * ( nza / pdims(2) ) |
---|
| 422 | nzt_yd = MIN( nzt, nzt_yda ) |
---|
| 423 | |
---|
| 424 | sendrecvcount_zyd = nnx * nny * ( nza / pdims(2) ) |
---|
| 425 | ENDIF |
---|
| 426 | ENDIF |
---|
| 427 | |
---|
| 428 | ! |
---|
| 429 | !-- Indices for direct transpositions y --> x (they are only possible in case |
---|
| 430 | !-- of a 1d-decomposition along x) |
---|
| 431 | IF ( pdims(2) == 1 ) THEN |
---|
| 432 | nny_x = nny / pdims(1) |
---|
| 433 | nys_x = myid * nny_x |
---|
| 434 | nyn_xa = ( myid + 1 ) * nny_x - 1 |
---|
| 435 | nyn_x = MIN( ny, nyn_xa ) |
---|
| 436 | nzb_x = 1 |
---|
| 437 | nzt_xa = nza |
---|
| 438 | nzt_x = nz |
---|
| 439 | sendrecvcount_xy = nnx * nny_x * nza |
---|
| 440 | ENDIF |
---|
| 441 | |
---|
| 442 | ! |
---|
| 443 | !-- Indices for direct transpositions x --> y (they are only possible in case |
---|
| 444 | !-- of a 1d-decomposition along y) |
---|
| 445 | IF ( pdims(1) == 1 ) THEN |
---|
| 446 | nnx_y = nnx / pdims(2) |
---|
| 447 | nxl_y = myid * nnx_y |
---|
| 448 | nxr_ya = ( myid + 1 ) * nnx_y - 1 |
---|
| 449 | nxr_y = MIN( nx, nxr_ya ) |
---|
| 450 | nzb_y = 1 |
---|
| 451 | nzt_ya = nza |
---|
| 452 | nzt_y = nz |
---|
| 453 | sendrecvcount_xy = nnx_y * nny * nza |
---|
| 454 | ENDIF |
---|
| 455 | |
---|
| 456 | ! |
---|
| 457 | !-- Arrays for storing the array bounds are needed any more |
---|
| 458 | DEALLOCATE( nxlf , nxrf , nynf , nysf ) |
---|
| 459 | |
---|
[145] | 460 | ! |
---|
| 461 | !-- Collect index bounds from other PEs (to be written to restart file later) |
---|
| 462 | ALLOCATE( hor_index_bounds(4,0:numprocs-1) ) |
---|
| 463 | |
---|
| 464 | IF ( myid == 0 ) THEN |
---|
| 465 | |
---|
| 466 | hor_index_bounds(1,0) = nxl |
---|
| 467 | hor_index_bounds(2,0) = nxr |
---|
| 468 | hor_index_bounds(3,0) = nys |
---|
| 469 | hor_index_bounds(4,0) = nyn |
---|
| 470 | |
---|
| 471 | ! |
---|
| 472 | !-- Receive data from all other PEs |
---|
| 473 | DO i = 1, numprocs-1 |
---|
| 474 | CALL MPI_RECV( ibuf, 4, MPI_INTEGER, i, MPI_ANY_TAG, comm2d, status, & |
---|
| 475 | ierr ) |
---|
| 476 | hor_index_bounds(:,i) = ibuf(1:4) |
---|
| 477 | ENDDO |
---|
| 478 | |
---|
| 479 | ELSE |
---|
| 480 | ! |
---|
| 481 | !-- Send index bounds to PE0 |
---|
| 482 | ibuf(1) = nxl |
---|
| 483 | ibuf(2) = nxr |
---|
| 484 | ibuf(3) = nys |
---|
| 485 | ibuf(4) = nyn |
---|
| 486 | CALL MPI_SEND( ibuf, 4, MPI_INTEGER, 0, myid, comm2d, ierr ) |
---|
| 487 | |
---|
| 488 | ENDIF |
---|
| 489 | |
---|
[1] | 490 | #if defined( __print ) |
---|
| 491 | ! |
---|
| 492 | !-- Control output |
---|
| 493 | IF ( myid == 0 ) THEN |
---|
| 494 | PRINT*, '*** processor topology ***' |
---|
| 495 | PRINT*, ' ' |
---|
| 496 | PRINT*, 'myid pcoord left right south north idx idy nxl: nxr',& |
---|
| 497 | &' nys: nyn' |
---|
| 498 | PRINT*, '------------------------------------------------------------',& |
---|
| 499 | &'-----------' |
---|
| 500 | WRITE (*,1000) 0, pcoord(1), pcoord(2), pleft, pright, psouth, pnorth, & |
---|
| 501 | myidx, myidy, nxl, nxr, nys, nyn |
---|
| 502 | 1000 FORMAT (I4,2X,'(',I3,',',I3,')',3X,I4,2X,I4,3X,I4,2X,I4,2X,I3,1X,I3, & |
---|
| 503 | 2(2X,I4,':',I4)) |
---|
| 504 | |
---|
| 505 | ! |
---|
[108] | 506 | !-- Receive data from the other PEs |
---|
[1] | 507 | DO i = 1,numprocs-1 |
---|
| 508 | CALL MPI_RECV( ibuf, 12, MPI_INTEGER, i, MPI_ANY_TAG, comm2d, status, & |
---|
| 509 | ierr ) |
---|
| 510 | WRITE (*,1000) i, ( ibuf(j) , j = 1,12 ) |
---|
| 511 | ENDDO |
---|
| 512 | ELSE |
---|
| 513 | |
---|
| 514 | ! |
---|
| 515 | !-- Send data to PE0 |
---|
| 516 | ibuf(1) = pcoord(1); ibuf(2) = pcoord(2); ibuf(3) = pleft |
---|
| 517 | ibuf(4) = pright; ibuf(5) = psouth; ibuf(6) = pnorth; ibuf(7) = myidx |
---|
| 518 | ibuf(8) = myidy; ibuf(9) = nxl; ibuf(10) = nxr; ibuf(11) = nys |
---|
| 519 | ibuf(12) = nyn |
---|
| 520 | CALL MPI_SEND( ibuf, 12, MPI_INTEGER, 0, myid, comm2d, ierr ) |
---|
| 521 | ENDIF |
---|
| 522 | #endif |
---|
| 523 | |
---|
[102] | 524 | #if defined( __mpi2 ) |
---|
| 525 | ! |
---|
| 526 | !-- In case of coupled runs, get the port name on PE0 of the atmosphere model |
---|
| 527 | !-- and pass it to PE0 of the ocean model |
---|
| 528 | IF ( myid == 0 ) THEN |
---|
| 529 | |
---|
| 530 | IF ( coupling_mode == 'atmosphere_to_ocean' ) THEN |
---|
| 531 | |
---|
| 532 | CALL MPI_OPEN_PORT( MPI_INFO_NULL, port_name, ierr ) |
---|
[108] | 533 | ! |
---|
| 534 | !-- TEST OUTPUT (TO BE REMOVED) |
---|
| 535 | WRITE(9,*) TRIM( coupling_mode ), & |
---|
| 536 | ', ierr after MPI_OPEN_PORT: ', ierr |
---|
| 537 | CALL LOCAL_FLUSH( 9 ) |
---|
| 538 | |
---|
[102] | 539 | CALL MPI_PUBLISH_NAME( 'palm_coupler', MPI_INFO_NULL, port_name, & |
---|
| 540 | ierr ) |
---|
[104] | 541 | ! |
---|
[108] | 542 | !-- TEST OUTPUT (TO BE REMOVED) |
---|
| 543 | WRITE(9,*) TRIM( coupling_mode ), & |
---|
| 544 | ', ierr after MPI_PUBLISH_NAME: ', ierr |
---|
| 545 | CALL LOCAL_FLUSH( 9 ) |
---|
| 546 | |
---|
| 547 | ! |
---|
[104] | 548 | !-- Write a flag file for the ocean model and the other atmosphere |
---|
| 549 | !-- processes. |
---|
| 550 | !-- There seems to be a bug in MPICH2 which causes hanging processes |
---|
| 551 | !-- in case that execution of LOOKUP_NAME is continued too early |
---|
| 552 | !-- (i.e. before the port has been created) |
---|
| 553 | OPEN( 90, FILE='COUPLING_PORT_OPENED', FORM='FORMATTED' ) |
---|
| 554 | WRITE ( 90, '(''TRUE'')' ) |
---|
| 555 | CLOSE ( 90 ) |
---|
[102] | 556 | |
---|
| 557 | ELSEIF ( coupling_mode == 'ocean_to_atmosphere' ) THEN |
---|
| 558 | |
---|
[104] | 559 | ! |
---|
| 560 | !-- Continue only if the atmosphere model has created the port. |
---|
| 561 | !-- There seems to be a bug in MPICH2 which causes hanging processes |
---|
| 562 | !-- in case that execution of LOOKUP_NAME is continued too early |
---|
| 563 | !-- (i.e. before the port has been created) |
---|
| 564 | INQUIRE( FILE='COUPLING_PORT_OPENED', EXIST=found ) |
---|
| 565 | DO WHILE ( .NOT. found ) |
---|
| 566 | INQUIRE( FILE='COUPLING_PORT_OPENED', EXIST=found ) |
---|
| 567 | ENDDO |
---|
| 568 | |
---|
[102] | 569 | CALL MPI_LOOKUP_NAME( 'palm_coupler', MPI_INFO_NULL, port_name, ierr ) |
---|
[108] | 570 | ! |
---|
| 571 | !-- TEST OUTPUT (TO BE REMOVED) |
---|
| 572 | WRITE(9,*) TRIM( coupling_mode ), & |
---|
| 573 | ', ierr after MPI_LOOKUP_NAME: ', ierr |
---|
| 574 | CALL LOCAL_FLUSH( 9 ) |
---|
[102] | 575 | |
---|
[108] | 576 | |
---|
[102] | 577 | ENDIF |
---|
| 578 | |
---|
| 579 | ENDIF |
---|
| 580 | |
---|
| 581 | ! |
---|
| 582 | !-- In case of coupled runs, establish the connection between the atmosphere |
---|
| 583 | !-- and the ocean model and define the intercommunicator (comm_inter) |
---|
| 584 | CALL MPI_BARRIER( comm2d, ierr ) |
---|
| 585 | IF ( coupling_mode == 'atmosphere_to_ocean' ) THEN |
---|
| 586 | |
---|
| 587 | print*, '... before COMM_ACCEPT' |
---|
| 588 | CALL MPI_COMM_ACCEPT( port_name, MPI_INFO_NULL, 0, MPI_COMM_WORLD, & |
---|
| 589 | comm_inter, ierr ) |
---|
| 590 | print*, '--- ierr = ', ierr |
---|
| 591 | print*, '--- comm_inter atmosphere = ', comm_inter |
---|
| 592 | |
---|
[108] | 593 | coupling_mode_remote = 'ocean_to_atmosphere' |
---|
| 594 | |
---|
[102] | 595 | ELSEIF ( coupling_mode == 'ocean_to_atmosphere' ) THEN |
---|
| 596 | |
---|
| 597 | IF ( myid == 0 ) PRINT*, '*** read: ', port_name, ' ierr = ', ierr |
---|
| 598 | print*, '... before COMM_CONNECT' |
---|
| 599 | CALL MPI_COMM_CONNECT( port_name, MPI_INFO_NULL, 0, MPI_COMM_WORLD, & |
---|
| 600 | comm_inter, ierr ) |
---|
| 601 | print*, '--- ierr = ', ierr |
---|
| 602 | print*, '--- comm_inter ocean = ', comm_inter |
---|
| 603 | |
---|
[108] | 604 | coupling_mode_remote = 'atmosphere_to_ocean' |
---|
| 605 | |
---|
[102] | 606 | ENDIF |
---|
| 607 | |
---|
| 608 | ! |
---|
| 609 | !-- In case of coupled runs, create a new MPI derived datatype for the |
---|
| 610 | !-- exchange of surface (xy) data . |
---|
| 611 | !-- Gridpoint number for the exchange of ghost points (xy-plane) |
---|
| 612 | ngp_xy = ( nxr - nxl + 3 ) * ( nyn - nys + 3 ) |
---|
| 613 | |
---|
| 614 | ! |
---|
| 615 | !-- Define a new MPI derived datatype for the exchange of ghost points in |
---|
| 616 | !-- y-direction for 2D-arrays (line) |
---|
| 617 | CALL MPI_TYPE_VECTOR( ngp_xy, 1, nzt-nzb+2, MPI_REAL, type_xy, ierr ) |
---|
| 618 | CALL MPI_TYPE_COMMIT( type_xy, ierr ) |
---|
| 619 | #endif |
---|
| 620 | |
---|
[1] | 621 | #else |
---|
| 622 | |
---|
| 623 | ! |
---|
| 624 | !-- Array bounds when running on a single PE (respectively a non-parallel |
---|
| 625 | !-- machine) |
---|
| 626 | nxl = 0 |
---|
| 627 | nxr = nx |
---|
| 628 | nxra = nx |
---|
| 629 | nnx = nxr - nxl + 1 |
---|
| 630 | nys = 0 |
---|
| 631 | nyn = ny |
---|
| 632 | nyna = ny |
---|
| 633 | nny = nyn - nys + 1 |
---|
| 634 | nzb = 0 |
---|
| 635 | nzt = nz |
---|
| 636 | nzta = nz |
---|
| 637 | nnz = nz |
---|
| 638 | |
---|
[145] | 639 | ALLOCATE( hor_index_bounds(4,0:0) ) |
---|
| 640 | hor_index_bounds(1,0) = nxl |
---|
| 641 | hor_index_bounds(2,0) = nxr |
---|
| 642 | hor_index_bounds(3,0) = nys |
---|
| 643 | hor_index_bounds(4,0) = nyn |
---|
| 644 | |
---|
[1] | 645 | ! |
---|
| 646 | !-- Array bounds for the pressure solver (in the parallel code, these bounds |
---|
| 647 | !-- are the ones for the transposed arrays) |
---|
| 648 | nys_x = nys |
---|
| 649 | nyn_x = nyn |
---|
| 650 | nyn_xa = nyn |
---|
| 651 | nzb_x = nzb + 1 |
---|
| 652 | nzt_x = nzt |
---|
| 653 | nzt_xa = nzt |
---|
| 654 | |
---|
| 655 | nxl_y = nxl |
---|
| 656 | nxr_y = nxr |
---|
| 657 | nxr_ya = nxr |
---|
| 658 | nzb_y = nzb + 1 |
---|
| 659 | nzt_y = nzt |
---|
| 660 | nzt_ya = nzt |
---|
| 661 | |
---|
| 662 | nxl_z = nxl |
---|
| 663 | nxr_z = nxr |
---|
| 664 | nxr_za = nxr |
---|
| 665 | nys_z = nys |
---|
| 666 | nyn_z = nyn |
---|
| 667 | nyn_za = nyn |
---|
| 668 | |
---|
| 669 | #endif |
---|
| 670 | |
---|
| 671 | ! |
---|
| 672 | !-- Calculate number of grid levels necessary for the multigrid poisson solver |
---|
| 673 | !-- as well as the gridpoint indices on each level |
---|
| 674 | IF ( psolver == 'multigrid' ) THEN |
---|
| 675 | |
---|
| 676 | ! |
---|
| 677 | !-- First calculate number of possible grid levels for the subdomains |
---|
| 678 | mg_levels_x = 1 |
---|
| 679 | mg_levels_y = 1 |
---|
| 680 | mg_levels_z = 1 |
---|
| 681 | |
---|
| 682 | i = nnx |
---|
| 683 | DO WHILE ( MOD( i, 2 ) == 0 .AND. i /= 2 ) |
---|
| 684 | i = i / 2 |
---|
| 685 | mg_levels_x = mg_levels_x + 1 |
---|
| 686 | ENDDO |
---|
| 687 | |
---|
| 688 | j = nny |
---|
| 689 | DO WHILE ( MOD( j, 2 ) == 0 .AND. j /= 2 ) |
---|
| 690 | j = j / 2 |
---|
| 691 | mg_levels_y = mg_levels_y + 1 |
---|
| 692 | ENDDO |
---|
| 693 | |
---|
| 694 | k = nnz |
---|
| 695 | DO WHILE ( MOD( k, 2 ) == 0 .AND. k /= 2 ) |
---|
| 696 | k = k / 2 |
---|
| 697 | mg_levels_z = mg_levels_z + 1 |
---|
| 698 | ENDDO |
---|
| 699 | |
---|
| 700 | maximum_grid_level = MIN( mg_levels_x, mg_levels_y, mg_levels_z ) |
---|
| 701 | |
---|
| 702 | ! |
---|
| 703 | !-- Find out, if the total domain allows more levels. These additional |
---|
| 704 | !-- levels are processed on PE0 only. |
---|
| 705 | IF ( numprocs > 1 ) THEN |
---|
| 706 | IF ( mg_levels_z > MIN( mg_levels_x, mg_levels_y ) ) THEN |
---|
| 707 | mg_switch_to_pe0_level_l = maximum_grid_level |
---|
| 708 | |
---|
| 709 | mg_levels_x = 1 |
---|
| 710 | mg_levels_y = 1 |
---|
| 711 | |
---|
| 712 | i = nx+1 |
---|
| 713 | DO WHILE ( MOD( i, 2 ) == 0 .AND. i /= 2 ) |
---|
| 714 | i = i / 2 |
---|
| 715 | mg_levels_x = mg_levels_x + 1 |
---|
| 716 | ENDDO |
---|
| 717 | |
---|
| 718 | j = ny+1 |
---|
| 719 | DO WHILE ( MOD( j, 2 ) == 0 .AND. j /= 2 ) |
---|
| 720 | j = j / 2 |
---|
| 721 | mg_levels_y = mg_levels_y + 1 |
---|
| 722 | ENDDO |
---|
| 723 | |
---|
| 724 | maximum_grid_level_l = MIN( mg_levels_x, mg_levels_y, mg_levels_z ) |
---|
| 725 | |
---|
| 726 | IF ( maximum_grid_level_l > mg_switch_to_pe0_level_l ) THEN |
---|
| 727 | mg_switch_to_pe0_level_l = maximum_grid_level_l - & |
---|
| 728 | mg_switch_to_pe0_level_l + 1 |
---|
| 729 | ELSE |
---|
| 730 | mg_switch_to_pe0_level_l = 0 |
---|
| 731 | ENDIF |
---|
| 732 | ELSE |
---|
| 733 | mg_switch_to_pe0_level_l = 0 |
---|
| 734 | maximum_grid_level_l = maximum_grid_level |
---|
| 735 | ENDIF |
---|
| 736 | |
---|
| 737 | ! |
---|
| 738 | !-- Use switch level calculated above only if it is not pre-defined |
---|
| 739 | !-- by user |
---|
| 740 | IF ( mg_switch_to_pe0_level == 0 ) THEN |
---|
| 741 | |
---|
| 742 | IF ( mg_switch_to_pe0_level_l /= 0 ) THEN |
---|
| 743 | mg_switch_to_pe0_level = mg_switch_to_pe0_level_l |
---|
| 744 | maximum_grid_level = maximum_grid_level_l |
---|
| 745 | ENDIF |
---|
| 746 | |
---|
| 747 | ELSE |
---|
| 748 | ! |
---|
| 749 | !-- Check pre-defined value and reset to default, if neccessary |
---|
| 750 | IF ( mg_switch_to_pe0_level < mg_switch_to_pe0_level_l .OR. & |
---|
| 751 | mg_switch_to_pe0_level >= maximum_grid_level_l ) THEN |
---|
| 752 | IF ( myid == 0 ) THEN |
---|
| 753 | PRINT*, '+++ WARNING init_pegrid: mg_switch_to_pe0_level ', & |
---|
| 754 | 'out of range and reset to default (=0)' |
---|
| 755 | ENDIF |
---|
| 756 | mg_switch_to_pe0_level = 0 |
---|
| 757 | ELSE |
---|
| 758 | ! |
---|
| 759 | !-- Use the largest number of possible levels anyway and recalculate |
---|
| 760 | !-- the switch level to this largest number of possible values |
---|
| 761 | maximum_grid_level = maximum_grid_level_l |
---|
| 762 | |
---|
| 763 | ENDIF |
---|
| 764 | ENDIF |
---|
| 765 | |
---|
| 766 | ENDIF |
---|
| 767 | |
---|
| 768 | ALLOCATE( grid_level_count(maximum_grid_level), & |
---|
| 769 | nxl_mg(maximum_grid_level), nxr_mg(maximum_grid_level), & |
---|
| 770 | nyn_mg(maximum_grid_level), nys_mg(maximum_grid_level), & |
---|
| 771 | nzt_mg(maximum_grid_level) ) |
---|
| 772 | |
---|
| 773 | grid_level_count = 0 |
---|
| 774 | nxl_l = nxl; nxr_l = nxr; nys_l = nys; nyn_l = nyn; nzt_l = nzt |
---|
| 775 | |
---|
| 776 | DO i = maximum_grid_level, 1 , -1 |
---|
| 777 | |
---|
| 778 | IF ( i == mg_switch_to_pe0_level ) THEN |
---|
| 779 | #if defined( __parallel ) |
---|
| 780 | ! |
---|
| 781 | !-- Save the grid size of the subdomain at the switch level, because |
---|
| 782 | !-- it is needed in poismg. |
---|
| 783 | !-- Array bounds of the local subdomain grids are gathered on PE0 |
---|
| 784 | ind(1) = nxl_l; ind(2) = nxr_l |
---|
| 785 | ind(3) = nys_l; ind(4) = nyn_l |
---|
| 786 | ind(5) = nzt_l |
---|
| 787 | ALLOCATE( ind_all(5*numprocs), mg_loc_ind(5,0:numprocs-1) ) |
---|
| 788 | CALL MPI_ALLGATHER( ind, 5, MPI_INTEGER, ind_all, 5, & |
---|
| 789 | MPI_INTEGER, comm2d, ierr ) |
---|
| 790 | DO j = 0, numprocs-1 |
---|
| 791 | DO k = 1, 5 |
---|
| 792 | mg_loc_ind(k,j) = ind_all(k+j*5) |
---|
| 793 | ENDDO |
---|
| 794 | ENDDO |
---|
| 795 | DEALLOCATE( ind_all ) |
---|
| 796 | ! |
---|
| 797 | !-- Calculate the grid size of the total domain gathered on PE0 |
---|
| 798 | nxr_l = ( nxr_l-nxl_l+1 ) * pdims(1) - 1 |
---|
| 799 | nxl_l = 0 |
---|
| 800 | nyn_l = ( nyn_l-nys_l+1 ) * pdims(2) - 1 |
---|
| 801 | nys_l = 0 |
---|
| 802 | ! |
---|
| 803 | !-- The size of this gathered array must not be larger than the |
---|
| 804 | !-- array tend, which is used in the multigrid scheme as a temporary |
---|
| 805 | !-- array |
---|
| 806 | subdomain_size = ( nxr - nxl + 3 ) * ( nyn - nys + 3 ) * & |
---|
| 807 | ( nzt - nzb + 2 ) |
---|
| 808 | gathered_size = ( nxr_l - nxl_l + 3 ) * ( nyn_l - nys_l + 3 ) * & |
---|
| 809 | ( nzt_l - nzb + 2 ) |
---|
| 810 | |
---|
| 811 | IF ( gathered_size > subdomain_size ) THEN |
---|
| 812 | IF ( myid == 0 ) THEN |
---|
| 813 | PRINT*, '+++ init_pegrid: not enough memory for storing ', & |
---|
| 814 | 'gathered multigrid data on PE0' |
---|
| 815 | ENDIF |
---|
| 816 | CALL local_stop |
---|
| 817 | ENDIF |
---|
| 818 | #else |
---|
| 819 | PRINT*, '+++ init_pegrid: multigrid gather/scatter impossible ', & |
---|
| 820 | 'in non parallel mode' |
---|
| 821 | CALL local_stop |
---|
| 822 | #endif |
---|
| 823 | ENDIF |
---|
| 824 | |
---|
| 825 | nxl_mg(i) = nxl_l |
---|
| 826 | nxr_mg(i) = nxr_l |
---|
| 827 | nys_mg(i) = nys_l |
---|
| 828 | nyn_mg(i) = nyn_l |
---|
| 829 | nzt_mg(i) = nzt_l |
---|
| 830 | |
---|
| 831 | nxl_l = nxl_l / 2 |
---|
| 832 | nxr_l = nxr_l / 2 |
---|
| 833 | nys_l = nys_l / 2 |
---|
| 834 | nyn_l = nyn_l / 2 |
---|
| 835 | nzt_l = nzt_l / 2 |
---|
| 836 | ENDDO |
---|
| 837 | |
---|
| 838 | ELSE |
---|
| 839 | |
---|
| 840 | maximum_grid_level = 1 |
---|
| 841 | |
---|
| 842 | ENDIF |
---|
| 843 | |
---|
| 844 | grid_level = maximum_grid_level |
---|
| 845 | |
---|
| 846 | #if defined( __parallel ) |
---|
| 847 | ! |
---|
| 848 | !-- Gridpoint number for the exchange of ghost points (y-line for 2D-arrays) |
---|
| 849 | ngp_y = nyn - nys + 1 |
---|
| 850 | |
---|
| 851 | ! |
---|
| 852 | !-- Define a new MPI derived datatype for the exchange of ghost points in |
---|
| 853 | !-- y-direction for 2D-arrays (line) |
---|
| 854 | CALL MPI_TYPE_VECTOR( nxr-nxl+3, 1, ngp_y+2, MPI_REAL, type_x, ierr ) |
---|
| 855 | CALL MPI_TYPE_COMMIT( type_x, ierr ) |
---|
| 856 | CALL MPI_TYPE_VECTOR( nxr-nxl+3, 1, ngp_y+2, MPI_INTEGER, type_x_int, ierr ) |
---|
| 857 | CALL MPI_TYPE_COMMIT( type_x_int, ierr ) |
---|
| 858 | |
---|
| 859 | ! |
---|
| 860 | !-- Calculate gridpoint numbers for the exchange of ghost points along x |
---|
| 861 | !-- (yz-plane for 3D-arrays) and define MPI derived data type(s) for the |
---|
| 862 | !-- exchange of ghost points in y-direction (xz-plane). |
---|
| 863 | !-- Do these calculations for the model grid and (if necessary) also |
---|
| 864 | !-- for the coarser grid levels used in the multigrid method |
---|
| 865 | ALLOCATE ( ngp_yz(maximum_grid_level), type_xz(maximum_grid_level) ) |
---|
| 866 | |
---|
| 867 | nxl_l = nxl; nxr_l = nxr; nys_l = nys; nyn_l = nyn; nzb_l = nzb; nzt_l = nzt |
---|
| 868 | |
---|
| 869 | DO i = maximum_grid_level, 1 , -1 |
---|
| 870 | ngp_yz(i) = (nzt_l - nzb_l + 2) * (nyn_l - nys_l + 3) |
---|
| 871 | |
---|
| 872 | CALL MPI_TYPE_VECTOR( nxr_l-nxl_l+3, nzt_l-nzb_l+2, ngp_yz(i), & |
---|
| 873 | MPI_REAL, type_xz(i), ierr ) |
---|
| 874 | CALL MPI_TYPE_COMMIT( type_xz(i), ierr ) |
---|
| 875 | |
---|
| 876 | nxl_l = nxl_l / 2 |
---|
| 877 | nxr_l = nxr_l / 2 |
---|
| 878 | nys_l = nys_l / 2 |
---|
| 879 | nyn_l = nyn_l / 2 |
---|
| 880 | nzt_l = nzt_l / 2 |
---|
| 881 | ENDDO |
---|
| 882 | #endif |
---|
| 883 | |
---|
| 884 | #if defined( __parallel ) |
---|
| 885 | ! |
---|
| 886 | !-- Setting of flags for inflow/outflow conditions in case of non-cyclic |
---|
[106] | 887 | !-- horizontal boundary conditions. |
---|
[1] | 888 | IF ( pleft == MPI_PROC_NULL ) THEN |
---|
[73] | 889 | IF ( bc_lr == 'dirichlet/radiation' ) THEN |
---|
[1] | 890 | inflow_l = .TRUE. |
---|
[73] | 891 | ELSEIF ( bc_lr == 'radiation/dirichlet' ) THEN |
---|
[1] | 892 | outflow_l = .TRUE. |
---|
| 893 | ENDIF |
---|
| 894 | ENDIF |
---|
| 895 | |
---|
| 896 | IF ( pright == MPI_PROC_NULL ) THEN |
---|
[73] | 897 | IF ( bc_lr == 'dirichlet/radiation' ) THEN |
---|
[1] | 898 | outflow_r = .TRUE. |
---|
[73] | 899 | ELSEIF ( bc_lr == 'radiation/dirichlet' ) THEN |
---|
[1] | 900 | inflow_r = .TRUE. |
---|
| 901 | ENDIF |
---|
| 902 | ENDIF |
---|
| 903 | |
---|
| 904 | IF ( psouth == MPI_PROC_NULL ) THEN |
---|
[73] | 905 | IF ( bc_ns == 'dirichlet/radiation' ) THEN |
---|
[1] | 906 | outflow_s = .TRUE. |
---|
[73] | 907 | ELSEIF ( bc_ns == 'radiation/dirichlet' ) THEN |
---|
[1] | 908 | inflow_s = .TRUE. |
---|
| 909 | ENDIF |
---|
| 910 | ENDIF |
---|
| 911 | |
---|
| 912 | IF ( pnorth == MPI_PROC_NULL ) THEN |
---|
[73] | 913 | IF ( bc_ns == 'dirichlet/radiation' ) THEN |
---|
[1] | 914 | inflow_n = .TRUE. |
---|
[73] | 915 | ELSEIF ( bc_ns == 'radiation/dirichlet' ) THEN |
---|
[1] | 916 | outflow_n = .TRUE. |
---|
| 917 | ENDIF |
---|
| 918 | ENDIF |
---|
| 919 | |
---|
| 920 | #else |
---|
[73] | 921 | IF ( bc_lr == 'dirichlet/radiation' ) THEN |
---|
[1] | 922 | inflow_l = .TRUE. |
---|
| 923 | outflow_r = .TRUE. |
---|
[73] | 924 | ELSEIF ( bc_lr == 'radiation/dirichlet' ) THEN |
---|
[1] | 925 | outflow_l = .TRUE. |
---|
| 926 | inflow_r = .TRUE. |
---|
| 927 | ENDIF |
---|
| 928 | |
---|
[73] | 929 | IF ( bc_ns == 'dirichlet/radiation' ) THEN |
---|
[1] | 930 | inflow_n = .TRUE. |
---|
| 931 | outflow_s = .TRUE. |
---|
[73] | 932 | ELSEIF ( bc_ns == 'radiation/dirichlet' ) THEN |
---|
[1] | 933 | outflow_n = .TRUE. |
---|
| 934 | inflow_s = .TRUE. |
---|
| 935 | ENDIF |
---|
| 936 | #endif |
---|
[106] | 937 | ! |
---|
[110] | 938 | !-- At the outflow, u or v, respectively, have to be calculated for one more |
---|
| 939 | !-- grid point. |
---|
[106] | 940 | IF ( outflow_l ) THEN |
---|
| 941 | nxlu = nxl + 1 |
---|
| 942 | ELSE |
---|
| 943 | nxlu = nxl |
---|
| 944 | ENDIF |
---|
| 945 | IF ( outflow_s ) THEN |
---|
| 946 | nysv = nys + 1 |
---|
| 947 | ELSE |
---|
| 948 | nysv = nys |
---|
| 949 | ENDIF |
---|
[1] | 950 | |
---|
| 951 | IF ( psolver == 'poisfft_hybrid' ) THEN |
---|
| 952 | CALL poisfft_hybrid_ini |
---|
[75] | 953 | ELSEIF ( psolver == 'poisfft' ) THEN |
---|
[1] | 954 | CALL poisfft_init |
---|
| 955 | ENDIF |
---|
| 956 | |
---|
[114] | 957 | ! |
---|
| 958 | !-- Allocate wall flag arrays used in the multigrid solver |
---|
| 959 | IF ( psolver == 'multigrid' ) THEN |
---|
| 960 | |
---|
| 961 | DO i = maximum_grid_level, 1, -1 |
---|
| 962 | |
---|
| 963 | SELECT CASE ( i ) |
---|
| 964 | |
---|
| 965 | CASE ( 1 ) |
---|
| 966 | ALLOCATE( wall_flags_1(nzb:nzt_mg(i)+1, & |
---|
| 967 | nys_mg(i)-1:nyn_mg(i)+1, & |
---|
| 968 | nxl_mg(i)-1:nxr_mg(i)+1) ) |
---|
| 969 | |
---|
| 970 | CASE ( 2 ) |
---|
| 971 | ALLOCATE( wall_flags_2(nzb:nzt_mg(i)+1, & |
---|
| 972 | nys_mg(i)-1:nyn_mg(i)+1, & |
---|
| 973 | nxl_mg(i)-1:nxr_mg(i)+1) ) |
---|
| 974 | |
---|
| 975 | CASE ( 3 ) |
---|
| 976 | ALLOCATE( wall_flags_3(nzb:nzt_mg(i)+1, & |
---|
| 977 | nys_mg(i)-1:nyn_mg(i)+1, & |
---|
| 978 | nxl_mg(i)-1:nxr_mg(i)+1) ) |
---|
| 979 | |
---|
| 980 | CASE ( 4 ) |
---|
| 981 | ALLOCATE( wall_flags_4(nzb:nzt_mg(i)+1, & |
---|
| 982 | nys_mg(i)-1:nyn_mg(i)+1, & |
---|
| 983 | nxl_mg(i)-1:nxr_mg(i)+1) ) |
---|
| 984 | |
---|
| 985 | CASE ( 5 ) |
---|
| 986 | ALLOCATE( wall_flags_5(nzb:nzt_mg(i)+1, & |
---|
| 987 | nys_mg(i)-1:nyn_mg(i)+1, & |
---|
| 988 | nxl_mg(i)-1:nxr_mg(i)+1) ) |
---|
| 989 | |
---|
| 990 | CASE ( 6 ) |
---|
| 991 | ALLOCATE( wall_flags_6(nzb:nzt_mg(i)+1, & |
---|
| 992 | nys_mg(i)-1:nyn_mg(i)+1, & |
---|
| 993 | nxl_mg(i)-1:nxr_mg(i)+1) ) |
---|
| 994 | |
---|
| 995 | CASE ( 7 ) |
---|
| 996 | ALLOCATE( wall_flags_7(nzb:nzt_mg(i)+1, & |
---|
| 997 | nys_mg(i)-1:nyn_mg(i)+1, & |
---|
| 998 | nxl_mg(i)-1:nxr_mg(i)+1) ) |
---|
| 999 | |
---|
| 1000 | CASE ( 8 ) |
---|
| 1001 | ALLOCATE( wall_flags_8(nzb:nzt_mg(i)+1, & |
---|
| 1002 | nys_mg(i)-1:nyn_mg(i)+1, & |
---|
| 1003 | nxl_mg(i)-1:nxr_mg(i)+1) ) |
---|
| 1004 | |
---|
| 1005 | CASE ( 9 ) |
---|
| 1006 | ALLOCATE( wall_flags_9(nzb:nzt_mg(i)+1, & |
---|
| 1007 | nys_mg(i)-1:nyn_mg(i)+1, & |
---|
| 1008 | nxl_mg(i)-1:nxr_mg(i)+1) ) |
---|
| 1009 | |
---|
| 1010 | CASE ( 10 ) |
---|
| 1011 | ALLOCATE( wall_flags_10(nzb:nzt_mg(i)+1, & |
---|
| 1012 | nys_mg(i)-1:nyn_mg(i)+1, & |
---|
| 1013 | nxl_mg(i)-1:nxr_mg(i)+1) ) |
---|
| 1014 | |
---|
| 1015 | CASE DEFAULT |
---|
| 1016 | IF ( myid == 0 ) PRINT*, '+++ init_pegrid: more than 10 ', & |
---|
| 1017 | ' multigrid levels' |
---|
| 1018 | CALL local_stop |
---|
| 1019 | |
---|
| 1020 | END SELECT |
---|
| 1021 | |
---|
| 1022 | ENDDO |
---|
| 1023 | |
---|
| 1024 | ENDIF |
---|
| 1025 | |
---|
[1] | 1026 | END SUBROUTINE init_pegrid |
---|