source: palm/trunk/SOURCE/init_pegrid.f90 @ 1491

Last change on this file since 1491 was 1469, checked in by maronga, 10 years ago

last commit documented

  • Property svn:keywords set to Id
File size: 41.2 KB
Line 
1 SUBROUTINE init_pegrid
2
3!--------------------------------------------------------------------------------!
4! This file is part of PALM.
5!
6! PALM is free software: you can redistribute it and/or modify it under the terms
7! of the GNU General Public License as published by the Free Software Foundation,
8! either version 3 of the License, or (at your option) any later 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-2014 Leibniz Universitaet Hannover
18!--------------------------------------------------------------------------------!
19!
20! Current revisions:
21! ------------------
22!
23!
24! Former revisions:
25! -----------------
26! $Id: init_pegrid.f90 1469 2014-09-24 14:09:56Z heinze $
27!
28! 1468 2014-09-24 14:06:57Z maronga
29! Adapted for use on up to 6-digit processor cores
30!
31! 1435 2014-07-21 10:37:02Z keck
32! bugfix: added missing parameter coupling_mode_remote to ONLY-attribute
33!
34! 1402 2014-05-09 14:25:13Z raasch
35! location messages modified
36!
37! 1384 2014-05-02 14:31:06Z raasch
38! location messages added
39!
40! 1353 2014-04-08 15:21:23Z heinze
41! REAL constants provided with KIND-attribute
42!
43! 1322 2014-03-20 16:38:49Z raasch
44! REAL functions provided with KIND-attribute
45!
46! 1320 2014-03-20 08:40:49Z raasch
47! ONLY-attribute added to USE-statements,
48! kind-parameters added to all INTEGER and REAL declaration statements,
49! kinds are defined in new module kinds,
50! revision history before 2012 removed,
51! comment fields (!:) to be used for variable explanations added to
52! all variable declaration statements
53!
54! 1304 2014-03-12 10:29:42Z raasch
55! bugfix: single core MPI runs missed some settings of transpose indices
56!
57! 1212 2013-08-15 08:46:27Z raasch
58! error message for poisfft_hybrid removed
59!
60! 1159 2013-05-21 11:58:22Z fricke
61! dirichlet/neumann and neumann/dirichlet removed
62!
63! 1139 2013-04-18 07:25:03Z raasch
64! bugfix for calculating the id of the PE carrying the recycling plane
65!
66! 1111 2013-03-08 23:54:10Z raasch
67! initialization of poisfft moved to module poisfft
68!
69! 1092 2013-02-02 11:24:22Z raasch
70! unused variables removed
71!
72! 1056 2012-11-16 15:28:04Z raasch
73! Indices for arrays n.._mg start from zero due to definition of arrays f2 and
74! p2 as automatic arrays in recursive subroutine next_mg_level
75!
76! 1041 2012-11-06 02:36:29Z raasch
77! a 2d virtual processor topology is used by default for all machines
78!
79! 1036 2012-10-22 13:43:42Z raasch
80! code put under GPL (PALM 3.9)
81!
82! 1003 2012-09-14 14:35:53Z raasch
83! subdomains must have identical size (grid matching = "match" removed)
84!
85! 1001 2012-09-13 14:08:46Z raasch
86! all actions concerning upstream-spline-method removed
87!
88! 978 2012-08-09 08:28:32Z fricke
89! dirichlet/neumann and neumann/dirichlet added
90! nxlu and nysv are also calculated for inflow boundary
91!
92! 809 2012-01-30 13:32:58Z maronga
93! Bugfix: replaced .AND. and .NOT. with && and ! in the preprocessor directives
94!
95! 807 2012-01-25 11:53:51Z maronga
96! New cpp directive "__check" implemented which is used by check_namelist_files
97!
98! Revision 1.1  1997/07/24 11:15:09  raasch
99! Initial revision
100!
101!
102! Description:
103! ------------
104! Determination of the virtual processor topology (if not prescribed by the
105! user)and computation of the grid point number and array bounds of the local
106! domains.
107!------------------------------------------------------------------------------!
108
109    USE control_parameters,                                                    &
110        ONLY:  bc_lr, bc_ns, coupling_mode, coupling_mode_remote,              &
111               coupling_topology, dt_dosp, gathered_size, grid_level,          &
112               grid_level_count, host, inflow_l, inflow_n, inflow_r, inflow_s, &
113               io_blocks, io_group, maximum_grid_level,                        &
114               maximum_parallel_io_streams, message_string,                    &
115               mg_switch_to_pe0_level, momentum_advec, psolver, outflow_l,     &
116               outflow_n, outflow_r, outflow_s, recycling_width, scalar_advec, &
117               subdomain_size 
118
119    USE grid_variables,                                                        &
120        ONLY:  dx
121       
122    USE indices,                                                               &
123        ONLY:  mg_loc_ind, nbgp, nnx, nny, nnz, nx, nx_a, nx_o, nxl, nxl_mg,   &
124               nxlu, nxr, nxr_mg, ny, ny_a, ny_o, nyn, nyn_mg, nys, nys_mg,    &
125               nysv, nz, nzb, nzt, nzt_mg, wall_flags_1, wall_flags_2,         &
126               wall_flags_3, wall_flags_4, wall_flags_5, wall_flags_6,         &
127               wall_flags_7, wall_flags_8, wall_flags_9, wall_flags_10
128
129    USE kinds
130     
131    USE pegrid
132 
133    USE transpose_indices,                                                     &
134        ONLY:  nxl_y, nxl_yd, nxl_z, nxr_y, nxr_yd, nxr_z, nyn_x, nyn_z, nys_x,&
135               nys_z, nzb_x, nzb_y, nzb_yd, nzt_x, nzt_yd, nzt_y
136
137    IMPLICIT NONE
138
139    INTEGER(iwp) ::  i                        !:
140    INTEGER(iwp) ::  id_inflow_l              !:
141    INTEGER(iwp) ::  id_recycling_l           !:
142    INTEGER(iwp) ::  ind(5)                   !:
143    INTEGER(iwp) ::  j                        !:
144    INTEGER(iwp) ::  k                        !:
145    INTEGER(iwp) ::  maximum_grid_level_l     !:
146    INTEGER(iwp) ::  mg_levels_x              !:
147    INTEGER(iwp) ::  mg_levels_y              !:
148    INTEGER(iwp) ::  mg_levels_z              !:
149    INTEGER(iwp) ::  mg_switch_to_pe0_level_l !:
150    INTEGER(iwp) ::  nnx_y                    !:
151    INTEGER(iwp) ::  nnx_z                    !:
152    INTEGER(iwp) ::  nny_x                    !:
153    INTEGER(iwp) ::  nny_z                    !:
154    INTEGER(iwp) ::  nnz_x                    !:
155    INTEGER(iwp) ::  nnz_y                    !:
156    INTEGER(iwp) ::  numproc_sqr              !:
157    INTEGER(iwp) ::  nxl_l                    !:
158    INTEGER(iwp) ::  nxr_l                    !:
159    INTEGER(iwp) ::  nyn_l                    !:
160    INTEGER(iwp) ::  nys_l                    !:
161    INTEGER(iwp) ::  nzb_l                    !:
162    INTEGER(iwp) ::  nzt_l                    !:
163    INTEGER(iwp) ::  omp_get_num_threads      !:
164
165    INTEGER(iwp), DIMENSION(:), ALLOCATABLE ::  ind_all !:
166    INTEGER(iwp), DIMENSION(:), ALLOCATABLE ::  nxlf    !:
167    INTEGER(iwp), DIMENSION(:), ALLOCATABLE ::  nxrf    !:
168    INTEGER(iwp), DIMENSION(:), ALLOCATABLE ::  nynf    !:
169    INTEGER(iwp), DIMENSION(:), ALLOCATABLE ::  nysf    !:
170
171    INTEGER(iwp), DIMENSION(2) :: pdims_remote          !:
172
173#if defined( __mpi2 )
174    LOGICAL ::  found                                   !:
175#endif
176
177!
178!-- Get the number of OpenMP threads
179    !$OMP PARALLEL
180#if defined( __intel_openmp_bug )
181    threads_per_task = omp_get_num_threads()
182#else
183!$  threads_per_task = omp_get_num_threads()
184#endif
185    !$OMP END PARALLEL
186
187
188#if defined( __parallel )
189
190    CALL location_message( 'creating virtual PE grids + MPI derived data types', &
191                           .FALSE. )
192!
193!-- Determine the processor topology or check it, if prescribed by the user
194    IF ( npex == -1  .AND.  npey == -1 )  THEN
195
196!
197!--    Automatic determination of the topology
198       numproc_sqr = SQRT( REAL( numprocs, KIND=wp ) )
199       pdims(1)    = MAX( numproc_sqr , 1 )
200       DO  WHILE ( MOD( numprocs , pdims(1) ) /= 0 )
201          pdims(1) = pdims(1) - 1
202       ENDDO
203       pdims(2) = numprocs / pdims(1)
204
205    ELSEIF ( npex /= -1  .AND.  npey /= -1 )  THEN
206
207!
208!--    Prescribed by user. Number of processors on the prescribed topology
209!--    must be equal to the number of PEs available to the job
210       IF ( ( npex * npey ) /= numprocs )  THEN
211          WRITE( message_string, * ) 'number of PEs of the prescribed ',      & 
212                 'topology (', npex*npey,') does not match & the number of ', & 
213                 'PEs available to the job (', numprocs, ')'
214          CALL message( 'init_pegrid', 'PA0221', 1, 2, 0, 6, 0 )
215       ENDIF
216       pdims(1) = npex
217       pdims(2) = npey
218
219    ELSE
220!
221!--    If the processor topology is prescribed by the user, the number of
222!--    PEs must be given in both directions
223       message_string = 'if the processor topology is prescribed by the, ' //  &
224                   ' user& both values of "npex" and "npey" must be given ' // &
225                   'in the &NAMELIST-parameter file'
226       CALL message( 'init_pegrid', 'PA0222', 1, 2, 0, 6, 0 )
227
228    ENDIF
229
230!
231!-- For communication speedup, set barriers in front of collective
232!-- communications by default on SGI-type systems
233    IF ( host(3:5) == 'sgi' )  collective_wait = .TRUE.
234
235!
236!-- If necessary, set horizontal boundary conditions to non-cyclic
237    IF ( bc_lr /= 'cyclic' )  cyclic(1) = .FALSE.
238    IF ( bc_ns /= 'cyclic' )  cyclic(2) = .FALSE.
239
240
241#if ! defined( __check)
242!
243!-- Create the virtual processor grid
244    CALL MPI_CART_CREATE( comm_palm, ndim, pdims, cyclic, reorder, &
245                          comm2d, ierr )
246    CALL MPI_COMM_RANK( comm2d, myid, ierr )
247    WRITE (myid_char,'(''_'',I6.6)')  myid
248
249    CALL MPI_CART_COORDS( comm2d, myid, ndim, pcoord, ierr )
250    CALL MPI_CART_SHIFT( comm2d, 0, 1, pleft, pright, ierr )
251    CALL MPI_CART_SHIFT( comm2d, 1, 1, psouth, pnorth, ierr )
252
253!
254!-- Determine sub-topologies for transpositions
255!-- Transposition from z to x:
256    remain_dims(1) = .TRUE.
257    remain_dims(2) = .FALSE.
258    CALL MPI_CART_SUB( comm2d, remain_dims, comm1dx, ierr )
259    CALL MPI_COMM_RANK( comm1dx, myidx, ierr )
260!
261!-- Transposition from x to y
262    remain_dims(1) = .FALSE.
263    remain_dims(2) = .TRUE.
264    CALL MPI_CART_SUB( comm2d, remain_dims, comm1dy, ierr )
265    CALL MPI_COMM_RANK( comm1dy, myidy, ierr )
266
267#endif
268
269!
270!-- Calculate array bounds along x-direction for every PE.
271    ALLOCATE( nxlf(0:pdims(1)-1), nxrf(0:pdims(1)-1), nynf(0:pdims(2)-1), &
272              nysf(0:pdims(2)-1) )
273
274    IF ( MOD( nx+1 , pdims(1) ) /= 0 )  THEN
275       WRITE( message_string, * ) 'x-direction: gridpoint number (',nx+1,') ',&
276                               'is not an& integral divisor of the number ',  &
277                               'processors (', pdims(1),')'
278       CALL message( 'init_pegrid', 'PA0225', 1, 2, 0, 6, 0 )
279    ELSE
280       nnx  = ( nx + 1 ) / pdims(1)
281       IF ( nnx*pdims(1) - ( nx + 1) > nnx )  THEN
282          WRITE( message_string, * ) 'x-direction: nx does not match the',    & 
283                       'requirements given by the number of PEs &used',       &
284                       '& please use nx = ', nx - ( pdims(1) - ( nnx*pdims(1) &
285                                   - ( nx + 1 ) ) ), ' instead of nx =', nx
286          CALL message( 'init_pegrid', 'PA0226', 1, 2, 0, 6, 0 )
287       ENDIF
288    ENDIF   
289
290!
291!-- Left and right array bounds, number of gridpoints
292    DO  i = 0, pdims(1)-1
293       nxlf(i)   = i * nnx
294       nxrf(i)   = ( i + 1 ) * nnx - 1
295    ENDDO
296
297!
298!-- Calculate array bounds in y-direction for every PE.
299    IF ( MOD( ny+1 , pdims(2) ) /= 0 )  THEN
300       WRITE( message_string, * ) 'y-direction: gridpoint number (',ny+1,') ', &
301                           'is not an& integral divisor of the number of',     &
302                           'processors (', pdims(2),')'
303       CALL message( 'init_pegrid', 'PA0227', 1, 2, 0, 6, 0 )
304    ELSE
305       nny  = ( ny + 1 ) / pdims(2)
306       IF ( nny*pdims(2) - ( ny + 1) > nny )  THEN
307          WRITE( message_string, * ) 'y-direction: ny does not match the',    &
308                       'requirements given by the number of PEs &used ',      &
309                       '& please use ny = ', ny - ( pdims(2) - ( nnx*pdims(2) &
310                                     - ( ny + 1 ) ) ), ' instead of ny =', ny
311          CALL message( 'init_pegrid', 'PA0228', 1, 2, 0, 6, 0 )
312       ENDIF
313    ENDIF   
314
315!
316!-- South and north array bounds
317    DO  j = 0, pdims(2)-1
318       nysf(j)   = j * nny
319       nynf(j)   = ( j + 1 ) * nny - 1
320    ENDDO
321
322!
323!-- Local array bounds of the respective PEs
324    nxl = nxlf(pcoord(1))
325    nxr = nxrf(pcoord(1))
326    nys = nysf(pcoord(2))
327    nyn = nynf(pcoord(2))
328    nzb = 0
329    nzt = nz
330    nnz = nz
331
332!
333!-- Set switches to define if the PE is situated at the border of the virtual
334!-- processor grid
335    IF ( nxl == 0 )   left_border_pe  = .TRUE.
336    IF ( nxr == nx )  right_border_pe = .TRUE.
337    IF ( nys == 0 )   south_border_pe = .TRUE.
338    IF ( nyn == ny )  north_border_pe = .TRUE.
339
340!
341!-- Calculate array bounds and gridpoint numbers for the transposed arrays
342!-- (needed in the pressure solver)
343!-- For the transposed arrays, cyclic boundaries as well as top and bottom
344!-- boundaries are omitted, because they are obstructive to the transposition
345
346!
347!-- 1. transposition  z --> x
348!-- This transposition is not neccessary in case of a 1d-decomposition along x
349    nys_x = nys
350    nyn_x = nyn
351    nny_x = nny
352    nnz_x = nz / pdims(1)
353    nzb_x = 1 + myidx * nnz_x
354    nzt_x = ( myidx + 1 ) * nnz_x
355    sendrecvcount_zx = nnx * nny * nnz_x
356
357    IF ( pdims(2) /= 1 )  THEN
358       IF ( MOD( nz , pdims(1) ) /= 0 )  THEN
359          WRITE( message_string, * ) 'transposition z --> x:',                &
360                       '&nz=',nz,' is not an integral divisior of pdims(1)=', &
361                                                                   pdims(1)
362          CALL message( 'init_pegrid', 'PA0230', 1, 2, 0, 6, 0 )
363       ENDIF
364    ENDIF
365
366!
367!-- 2. transposition  x --> y
368    nnz_y = nnz_x
369    nzb_y = nzb_x
370    nzt_y = nzt_x
371    IF ( MOD( nx+1 , pdims(2) ) /= 0 )  THEN
372       WRITE( message_string, * ) 'transposition x --> y:',                &
373                         '&nx+1=',nx+1,' is not an integral divisor of ',&
374                         'pdims(2)=',pdims(2)
375       CALL message( 'init_pegrid', 'PA0231', 1, 2, 0, 6, 0 )
376    ENDIF
377    nnx_y = (nx+1) / pdims(2)
378    nxl_y = myidy * nnx_y
379    nxr_y = ( myidy + 1 ) * nnx_y - 1
380    sendrecvcount_xy = nnx_y * nny_x * nnz_y
381
382!
383!-- 3. transposition  y --> z  (ELSE:  x --> y  in case of 1D-decomposition
384!-- along x)
385    nnx_z = nnx_y
386    nxl_z = nxl_y
387    nxr_z = nxr_y
388    nny_z = (ny+1) / pdims(1)
389    nys_z = myidx * nny_z
390    nyn_z = ( myidx + 1 ) * nny_z - 1
391    sendrecvcount_yz = nnx_y * nny_z * nnz_y
392
393    IF ( pdims(2) /= 1 )  THEN
394!
395!--    y --> z
396!--    This transposition is not neccessary in case of a 1d-decomposition
397!--    along x, except that the uptream-spline method is switched on
398       IF ( MOD( ny+1 , pdims(1) ) /= 0 )  THEN
399          WRITE( message_string, * ) 'transposition y --> z:',            &
400                            '& ny+1=',ny+1,' is not an integral divisor of ',&
401                            'pdims(1)=',pdims(1)
402          CALL message( 'init_pegrid', 'PA0232', 1, 2, 0, 6, 0 )
403       ENDIF
404
405    ELSE
406!
407!--    x --> y. This condition must be fulfilled for a 1D-decomposition along x
408       IF ( MOD( ny+1 , pdims(1) ) /= 0 )  THEN
409          WRITE( message_string, * ) 'transposition x --> y:',               &
410                            '& ny+1=',ny+1,' is not an integral divisor of ',&
411                            'pdims(1)=',pdims(1)
412          CALL message( 'init_pegrid', 'PA0233', 1, 2, 0, 6, 0 )
413       ENDIF
414
415    ENDIF
416
417!
418!-- Indices for direct transpositions z --> y (used for calculating spectra)
419    IF ( dt_dosp /= 9999999.9_wp )  THEN
420       IF ( MOD( nz, pdims(2) ) /= 0 )  THEN
421          WRITE( message_string, * ) 'direct transposition z --> y (needed ', &
422                    'for spectra):& nz=',nz,' is not an integral divisor of ',&
423                    'pdims(2)=',pdims(2)
424          CALL message( 'init_pegrid', 'PA0234', 1, 2, 0, 6, 0 )
425       ELSE
426          nxl_yd = nxl
427          nxr_yd = nxr
428          nzb_yd = 1 + myidy * ( nz / pdims(2) )
429          nzt_yd = ( myidy + 1 ) * ( nz / pdims(2) )
430          sendrecvcount_zyd = nnx * nny * ( nz / pdims(2) )
431       ENDIF
432    ENDIF
433
434!
435!-- Indices for direct transpositions y --> x (they are only possible in case
436!-- of a 1d-decomposition along x)
437    IF ( pdims(2) == 1 )  THEN
438       nny_x = nny / pdims(1)
439       nys_x = myid * nny_x
440       nyn_x = ( myid + 1 ) * nny_x - 1
441       nzb_x = 1
442       nzt_x = nz
443       sendrecvcount_xy = nnx * nny_x * nz
444    ENDIF
445
446!
447!-- Indices for direct transpositions x --> y (they are only possible in case
448!-- of a 1d-decomposition along y)
449    IF ( pdims(1) == 1 )  THEN
450       nnx_y = nnx / pdims(2)
451       nxl_y = myid * nnx_y
452       nxr_y = ( myid + 1 ) * nnx_y - 1
453       nzb_y = 1
454       nzt_y = nz
455       sendrecvcount_xy = nnx_y * nny * nz
456    ENDIF
457
458!
459!-- Arrays for storing the array bounds are needed any more
460    DEALLOCATE( nxlf , nxrf , nynf , nysf )
461
462
463#if ! defined( __check)
464!
465!-- Collect index bounds from other PEs (to be written to restart file later)
466    ALLOCATE( hor_index_bounds(4,0:numprocs-1) )
467
468    IF ( myid == 0 )  THEN
469
470       hor_index_bounds(1,0) = nxl
471       hor_index_bounds(2,0) = nxr
472       hor_index_bounds(3,0) = nys
473       hor_index_bounds(4,0) = nyn
474
475!
476!--    Receive data from all other PEs
477       DO  i = 1, numprocs-1
478          CALL MPI_RECV( ibuf, 4, MPI_INTEGER, i, MPI_ANY_TAG, comm2d, status, &
479                         ierr )
480          hor_index_bounds(:,i) = ibuf(1:4)
481       ENDDO
482
483    ELSE
484!
485!--    Send index bounds to PE0
486       ibuf(1) = nxl
487       ibuf(2) = nxr
488       ibuf(3) = nys
489       ibuf(4) = nyn
490       CALL MPI_SEND( ibuf, 4, MPI_INTEGER, 0, myid, comm2d, ierr )
491
492    ENDIF
493
494#endif
495
496#if defined( __print )
497!
498!-- Control output
499    IF ( myid == 0 )  THEN
500       PRINT*, '*** processor topology ***'
501       PRINT*, ' '
502       PRINT*, 'myid   pcoord    left right  south north  idx idy   nxl: nxr',&
503               &'   nys: nyn'
504       PRINT*, '------------------------------------------------------------',&
505               &'-----------'
506       WRITE (*,1000)  0, pcoord(1), pcoord(2), pleft, pright, psouth, pnorth, &
507                       myidx, myidy, nxl, nxr, nys, nyn
5081000   FORMAT (I4,2X,'(',I3,',',I3,')',3X,I4,2X,I4,3X,I4,2X,I4,2X,I3,1X,I3, &
509               2(2X,I4,':',I4))
510
511!
512!--    Receive data from the other PEs
513       DO  i = 1,numprocs-1
514          CALL MPI_RECV( ibuf, 12, MPI_INTEGER, i, MPI_ANY_TAG, comm2d, status, &
515                         ierr )
516          WRITE (*,1000)  i, ( ibuf(j) , j = 1,12 )
517       ENDDO
518    ELSE
519
520!
521!--    Send data to PE0
522       ibuf(1) = pcoord(1); ibuf(2) = pcoord(2); ibuf(3) = pleft
523       ibuf(4) = pright; ibuf(5) = psouth; ibuf(6) = pnorth; ibuf(7) = myidx
524       ibuf(8) = myidy; ibuf(9) = nxl; ibuf(10) = nxr; ibuf(11) = nys
525       ibuf(12) = nyn
526       CALL MPI_SEND( ibuf, 12, MPI_INTEGER, 0, myid, comm2d, ierr )       
527    ENDIF
528#endif
529
530#if defined( __parallel ) && ! defined( __check)
531#if defined( __mpi2 )
532!
533!-- In case of coupled runs, get the port name on PE0 of the atmosphere model
534!-- and pass it to PE0 of the ocean model
535    IF ( myid == 0 )  THEN
536
537       IF ( coupling_mode == 'atmosphere_to_ocean' )  THEN
538
539          CALL MPI_OPEN_PORT( MPI_INFO_NULL, port_name, ierr )
540
541          CALL MPI_PUBLISH_NAME( 'palm_coupler', MPI_INFO_NULL, port_name, &
542                                 ierr )
543
544!
545!--       Write a flag file for the ocean model and the other atmosphere
546!--       processes.
547!--       There seems to be a bug in MPICH2 which causes hanging processes
548!--       in case that execution of LOOKUP_NAME is continued too early
549!--       (i.e. before the port has been created)
550          OPEN( 90, FILE='COUPLING_PORT_OPENED', FORM='FORMATTED' )
551          WRITE ( 90, '(''TRUE'')' )
552          CLOSE ( 90 )
553
554       ELSEIF ( coupling_mode == 'ocean_to_atmosphere' )  THEN
555
556!
557!--       Continue only if the atmosphere model has created the port.
558!--       There seems to be a bug in MPICH2 which causes hanging processes
559!--       in case that execution of LOOKUP_NAME is continued too early
560!--       (i.e. before the port has been created)
561          INQUIRE( FILE='COUPLING_PORT_OPENED', EXIST=found )
562          DO WHILE ( .NOT. found )
563             INQUIRE( FILE='COUPLING_PORT_OPENED', EXIST=found )
564          ENDDO
565
566          CALL MPI_LOOKUP_NAME( 'palm_coupler', MPI_INFO_NULL, port_name, ierr )
567
568       ENDIF
569
570    ENDIF
571
572!
573!-- In case of coupled runs, establish the connection between the atmosphere
574!-- and the ocean model and define the intercommunicator (comm_inter)
575    CALL MPI_BARRIER( comm2d, ierr )
576    IF ( coupling_mode == 'atmosphere_to_ocean' )  THEN
577
578       CALL MPI_COMM_ACCEPT( port_name, MPI_INFO_NULL, 0, MPI_COMM_WORLD, &
579                             comm_inter, ierr )
580       coupling_mode_remote = 'ocean_to_atmosphere'
581
582    ELSEIF ( coupling_mode == 'ocean_to_atmosphere' )  THEN
583
584       CALL MPI_COMM_CONNECT( port_name, MPI_INFO_NULL, 0, MPI_COMM_WORLD, &
585                              comm_inter, ierr )
586       coupling_mode_remote = 'atmosphere_to_ocean'
587
588    ENDIF
589#endif
590
591!
592!-- Determine the number of ghost point layers
593    IF ( scalar_advec == 'ws-scheme' .OR. momentum_advec == 'ws-scheme' )  THEN
594       nbgp = 3
595    ELSE
596       nbgp = 1
597    ENDIF 
598
599!
600!-- Create a new MPI derived datatype for the exchange of surface (xy) data,
601!-- which is needed for coupled atmosphere-ocean runs.
602!-- First, calculate number of grid points of an xy-plane.
603    ngp_xy  = ( nxr - nxl + 1 + 2 * nbgp ) * ( nyn - nys + 1 + 2 * nbgp )
604    CALL MPI_TYPE_VECTOR( ngp_xy, 1, nzt-nzb+2, MPI_REAL, type_xy, ierr )
605    CALL MPI_TYPE_COMMIT( type_xy, ierr )
606
607    IF ( TRIM( coupling_mode ) /= 'uncoupled' )  THEN
608   
609!
610!--    Pass the number of grid points of the atmosphere model to
611!--    the ocean model and vice versa
612       IF ( coupling_mode == 'atmosphere_to_ocean' )  THEN
613
614          nx_a = nx
615          ny_a = ny
616
617          IF ( myid == 0 )  THEN
618
619             CALL MPI_SEND( nx_a, 1, MPI_INTEGER, numprocs, 1, comm_inter,  &
620                            ierr )
621             CALL MPI_SEND( ny_a, 1, MPI_INTEGER, numprocs, 2, comm_inter,  &
622                            ierr )
623             CALL MPI_SEND( pdims, 2, MPI_INTEGER, numprocs, 3, comm_inter, &
624                            ierr )
625             CALL MPI_RECV( nx_o, 1, MPI_INTEGER, numprocs, 4, comm_inter,  &
626                            status, ierr )
627             CALL MPI_RECV( ny_o, 1, MPI_INTEGER, numprocs, 5, comm_inter,  &
628                            status, ierr )
629             CALL MPI_RECV( pdims_remote, 2, MPI_INTEGER, numprocs, 6,      &
630                            comm_inter, status, ierr )
631          ENDIF
632
633          CALL MPI_BCAST( nx_o, 1, MPI_INTEGER, 0, comm2d, ierr )
634          CALL MPI_BCAST( ny_o, 1, MPI_INTEGER, 0, comm2d, ierr ) 
635          CALL MPI_BCAST( pdims_remote, 2, MPI_INTEGER, 0, comm2d, ierr )
636       
637       ELSEIF ( coupling_mode == 'ocean_to_atmosphere' )  THEN
638
639          nx_o = nx
640          ny_o = ny 
641
642          IF ( myid == 0 ) THEN
643
644             CALL MPI_RECV( nx_a, 1, MPI_INTEGER, 0, 1, comm_inter, status, &
645                            ierr )
646             CALL MPI_RECV( ny_a, 1, MPI_INTEGER, 0, 2, comm_inter, status, &
647                            ierr )
648             CALL MPI_RECV( pdims_remote, 2, MPI_INTEGER, 0, 3, comm_inter, &
649                            status, ierr )
650             CALL MPI_SEND( nx_o, 1, MPI_INTEGER, 0, 4, comm_inter, ierr )
651             CALL MPI_SEND( ny_o, 1, MPI_INTEGER, 0, 5, comm_inter, ierr )
652             CALL MPI_SEND( pdims, 2, MPI_INTEGER, 0, 6, comm_inter, ierr )
653          ENDIF
654
655          CALL MPI_BCAST( nx_a, 1, MPI_INTEGER, 0, comm2d, ierr)
656          CALL MPI_BCAST( ny_a, 1, MPI_INTEGER, 0, comm2d, ierr) 
657          CALL MPI_BCAST( pdims_remote, 2, MPI_INTEGER, 0, comm2d, ierr) 
658
659       ENDIF
660 
661       ngp_a = ( nx_a+1 + 2 * nbgp ) * ( ny_a+1 + 2 * nbgp )
662       ngp_o = ( nx_o+1 + 2 * nbgp ) * ( ny_o+1 + 2 * nbgp )
663
664!
665!--    Determine if the horizontal grid and the number of PEs in ocean and
666!--    atmosphere is same or not
667       IF ( nx_o == nx_a  .AND.  ny_o == ny_a  .AND.  &
668            pdims(1) == pdims_remote(1) .AND. pdims(2) == pdims_remote(2) ) &
669       THEN
670          coupling_topology = 0
671       ELSE
672          coupling_topology = 1
673       ENDIF 
674
675!
676!--    Determine the target PEs for the exchange between ocean and
677!--    atmosphere (comm2d)
678       IF ( coupling_topology == 0 )  THEN
679!
680!--       In case of identical topologies, every atmosphere PE has exactly one
681!--       ocean PE counterpart and vice versa
682          IF ( TRIM( coupling_mode ) == 'atmosphere_to_ocean' ) THEN
683             target_id = myid + numprocs
684          ELSE
685             target_id = myid 
686          ENDIF
687
688       ELSE
689!
690!--       In case of nonequivalent topology in ocean and atmosphere only for
691!--       PE0 in ocean and PE0 in atmosphere a target_id is needed, since
692!--       data echxchange between ocean and atmosphere will be done only
693!--       between these PEs.   
694          IF ( myid == 0 )  THEN
695
696             IF ( TRIM( coupling_mode ) == 'atmosphere_to_ocean' )  THEN
697                target_id = numprocs 
698             ELSE
699                target_id = 0
700             ENDIF
701
702          ENDIF
703
704       ENDIF
705
706    ENDIF
707
708
709#endif
710
711#else
712
713!
714!-- Array bounds when running on a single PE (respectively a non-parallel
715!-- machine)
716    nxl = 0
717    nxr = nx
718    nnx = nxr - nxl + 1
719    nys = 0
720    nyn = ny
721    nny = nyn - nys + 1
722    nzb = 0
723    nzt = nz
724    nnz = nz
725
726    ALLOCATE( hor_index_bounds(4,0:0) )
727    hor_index_bounds(1,0) = nxl
728    hor_index_bounds(2,0) = nxr
729    hor_index_bounds(3,0) = nys
730    hor_index_bounds(4,0) = nyn
731
732!
733!-- Array bounds for the pressure solver (in the parallel code, these bounds
734!-- are the ones for the transposed arrays)
735    nys_x = nys
736    nyn_x = nyn
737    nzb_x = nzb + 1
738    nzt_x = nzt
739
740    nxl_y = nxl
741    nxr_y = nxr
742    nzb_y = nzb + 1
743    nzt_y = nzt
744
745    nxl_z = nxl
746    nxr_z = nxr
747    nys_z = nys
748    nyn_z = nyn
749
750#endif
751
752!
753!-- Calculate number of grid levels necessary for the multigrid poisson solver
754!-- as well as the gridpoint indices on each level
755    IF ( psolver == 'multigrid' )  THEN
756
757!
758!--    First calculate number of possible grid levels for the subdomains
759       mg_levels_x = 1
760       mg_levels_y = 1
761       mg_levels_z = 1
762
763       i = nnx
764       DO WHILE ( MOD( i, 2 ) == 0  .AND.  i /= 2 )
765          i = i / 2
766          mg_levels_x = mg_levels_x + 1
767       ENDDO
768
769       j = nny
770       DO WHILE ( MOD( j, 2 ) == 0  .AND.  j /= 2 )
771          j = j / 2
772          mg_levels_y = mg_levels_y + 1
773       ENDDO
774
775       k = nz    ! do not use nnz because it might be > nz due to transposition
776                 ! requirements
777       DO WHILE ( MOD( k, 2 ) == 0  .AND.  k /= 2 )
778          k = k / 2
779          mg_levels_z = mg_levels_z + 1
780       ENDDO
781
782       maximum_grid_level = MIN( mg_levels_x, mg_levels_y, mg_levels_z )
783
784!
785!--    Find out, if the total domain allows more levels. These additional
786!--    levels are identically processed on all PEs.
787       IF ( numprocs > 1  .AND.  mg_switch_to_pe0_level /= -1 )  THEN
788
789          IF ( mg_levels_z > MIN( mg_levels_x, mg_levels_y ) )  THEN
790
791             mg_switch_to_pe0_level_l = maximum_grid_level
792
793             mg_levels_x = 1
794             mg_levels_y = 1
795
796             i = nx+1
797             DO WHILE ( MOD( i, 2 ) == 0  .AND.  i /= 2 )
798                i = i / 2
799                mg_levels_x = mg_levels_x + 1
800             ENDDO
801
802             j = ny+1
803             DO WHILE ( MOD( j, 2 ) == 0  .AND.  j /= 2 )
804                j = j / 2
805                mg_levels_y = mg_levels_y + 1
806             ENDDO
807
808             maximum_grid_level_l = MIN( mg_levels_x, mg_levels_y, mg_levels_z )
809
810             IF ( maximum_grid_level_l > mg_switch_to_pe0_level_l )  THEN
811                mg_switch_to_pe0_level_l = maximum_grid_level_l - &
812                                           mg_switch_to_pe0_level_l + 1
813             ELSE
814                mg_switch_to_pe0_level_l = 0
815             ENDIF
816
817          ELSE
818             mg_switch_to_pe0_level_l = 0
819             maximum_grid_level_l = maximum_grid_level
820
821          ENDIF
822
823!
824!--       Use switch level calculated above only if it is not pre-defined
825!--       by user
826          IF ( mg_switch_to_pe0_level == 0 )  THEN
827             IF ( mg_switch_to_pe0_level_l /= 0 )  THEN
828                mg_switch_to_pe0_level = mg_switch_to_pe0_level_l
829                maximum_grid_level     = maximum_grid_level_l
830             ENDIF
831
832          ELSE
833!
834!--          Check pre-defined value and reset to default, if neccessary
835             IF ( mg_switch_to_pe0_level < mg_switch_to_pe0_level_l  .OR.  &
836                  mg_switch_to_pe0_level >= maximum_grid_level_l )  THEN
837                message_string = 'mg_switch_to_pe0_level ' // &
838                                 'out of range and reset to default (=0)'
839                CALL message( 'init_pegrid', 'PA0235', 0, 1, 0, 6, 0 )
840                mg_switch_to_pe0_level = 0
841             ELSE
842!
843!--             Use the largest number of possible levels anyway and recalculate
844!--             the switch level to this largest number of possible values
845                maximum_grid_level = maximum_grid_level_l
846
847             ENDIF
848
849          ENDIF
850
851       ENDIF
852
853       ALLOCATE( grid_level_count(maximum_grid_level),                       &
854                 nxl_mg(0:maximum_grid_level), nxr_mg(0:maximum_grid_level), &
855                 nyn_mg(0:maximum_grid_level), nys_mg(0:maximum_grid_level), &
856                 nzt_mg(0:maximum_grid_level) )
857
858       grid_level_count = 0
859!
860!--    Index zero required as dummy due to definition of arrays f2 and p2 in
861!--    recursive subroutine next_mg_level
862       nxl_mg(0) = 0; nxr_mg(0) = 0; nyn_mg(0) = 0; nys_mg(0) = 0; nzt_mg(0) = 0
863
864       nxl_l = nxl; nxr_l = nxr; nys_l = nys; nyn_l = nyn; nzt_l = nzt
865
866       DO  i = maximum_grid_level, 1 , -1
867
868          IF ( i == mg_switch_to_pe0_level )  THEN
869#if defined( __parallel ) && ! defined( __check )
870!
871!--          Save the grid size of the subdomain at the switch level, because
872!--          it is needed in poismg.
873             ind(1) = nxl_l; ind(2) = nxr_l
874             ind(3) = nys_l; ind(4) = nyn_l
875             ind(5) = nzt_l
876             ALLOCATE( ind_all(5*numprocs), mg_loc_ind(5,0:numprocs-1) )
877             CALL MPI_ALLGATHER( ind, 5, MPI_INTEGER, ind_all, 5, &
878                                 MPI_INTEGER, comm2d, ierr )
879             DO  j = 0, numprocs-1
880                DO  k = 1, 5
881                   mg_loc_ind(k,j) = ind_all(k+j*5)
882                ENDDO
883             ENDDO
884             DEALLOCATE( ind_all )
885!
886!--          Calculate the grid size of the total domain
887             nxr_l = ( nxr_l-nxl_l+1 ) * pdims(1) - 1
888             nxl_l = 0
889             nyn_l = ( nyn_l-nys_l+1 ) * pdims(2) - 1
890             nys_l = 0
891!
892!--          The size of this gathered array must not be larger than the
893!--          array tend, which is used in the multigrid scheme as a temporary
894!--          array. Therefore the subdomain size of an PE is calculated and
895!--          the size of the gathered grid. These values are used in 
896!--          routines pres and poismg
897             subdomain_size = ( nxr - nxl + 2 * nbgp + 1 ) * &
898                              ( nyn - nys + 2 * nbgp + 1 ) * ( nzt - nzb + 2 )
899             gathered_size  = ( nxr_l - nxl_l + 3 ) * ( nyn_l - nys_l + 3 ) * &
900                              ( nzt_l - nzb + 2 )
901
902#elif ! defined ( __parallel )
903             message_string = 'multigrid gather/scatter impossible ' // &
904                          'in non parallel mode'
905             CALL message( 'init_pegrid', 'PA0237', 1, 2, 0, 6, 0 )
906#endif
907          ENDIF
908
909          nxl_mg(i) = nxl_l
910          nxr_mg(i) = nxr_l
911          nys_mg(i) = nys_l
912          nyn_mg(i) = nyn_l
913          nzt_mg(i) = nzt_l
914
915          nxl_l = nxl_l / 2 
916          nxr_l = nxr_l / 2
917          nys_l = nys_l / 2 
918          nyn_l = nyn_l / 2 
919          nzt_l = nzt_l / 2 
920
921       ENDDO
922
923!
924!--    Temporary problem: Currently calculation of maxerror iin routine poismg crashes
925!--    if grid data are collected on PE0 already on the finest grid level.
926!--    To be solved later.
927       IF ( maximum_grid_level == mg_switch_to_pe0_level )  THEN
928          message_string = 'grid coarsening on subdomain level cannot be performed'
929          CALL message( 'poismg', 'PA0236', 1, 2, 0, 6, 0 )
930       ENDIF
931
932    ELSE
933
934       maximum_grid_level = 0
935
936    ENDIF
937
938!
939!-- Default level 0 tells exchange_horiz that all ghost planes have to be
940!-- exchanged. grid_level is adjusted in poismg, where only one ghost plane
941!-- is required.
942    grid_level = 0
943
944#if defined( __parallel ) && ! defined ( __check )
945!
946!-- Gridpoint number for the exchange of ghost points (y-line for 2D-arrays)
947    ngp_y  = nyn - nys + 1 + 2 * nbgp
948
949!
950!-- Define new MPI derived datatypes for the exchange of ghost points in
951!-- x- and y-direction for 2D-arrays (line)
952    CALL MPI_TYPE_VECTOR( nxr-nxl+1+2*nbgp, nbgp, ngp_y, MPI_REAL, type_x, &
953                          ierr )
954    CALL MPI_TYPE_COMMIT( type_x, ierr )
955    CALL MPI_TYPE_VECTOR( nxr-nxl+1+2*nbgp, nbgp, ngp_y, MPI_INTEGER, &
956                          type_x_int, ierr )
957    CALL MPI_TYPE_COMMIT( type_x_int, ierr )
958
959    CALL MPI_TYPE_VECTOR( nbgp, ngp_y, ngp_y, MPI_REAL, type_y, ierr )
960    CALL MPI_TYPE_COMMIT( type_y, ierr )
961    CALL MPI_TYPE_VECTOR( nbgp, ngp_y, ngp_y, MPI_INTEGER, type_y_int, ierr )
962    CALL MPI_TYPE_COMMIT( type_y_int, ierr )
963
964
965!
966!-- Calculate gridpoint numbers for the exchange of ghost points along x
967!-- (yz-plane for 3D-arrays) and define MPI derived data type(s) for the
968!-- exchange of ghost points in y-direction (xz-plane).
969!-- Do these calculations for the model grid and (if necessary) also
970!-- for the coarser grid levels used in the multigrid method
971    ALLOCATE ( ngp_yz(0:maximum_grid_level), type_xz(0:maximum_grid_level),&
972               type_yz(0:maximum_grid_level) )
973
974    nxl_l = nxl; nxr_l = nxr; nys_l = nys; nyn_l = nyn; nzb_l = nzb; nzt_l = nzt
975
976!
977!-- Discern between the model grid, which needs nbgp ghost points and
978!-- grid levels for the multigrid scheme. In the latter case only one
979!-- ghost point is necessary.
980!-- First definition of MPI-datatypes for exchange of ghost layers on normal
981!-- grid. The following loop is needed for data exchange in poismg.f90.
982!
983!-- Determine number of grid points of yz-layer for exchange
984    ngp_yz(0) = (nzt - nzb + 2) * (nyn - nys + 1 + 2 * nbgp)
985
986!
987!-- Define an MPI-datatype for the exchange of left/right boundaries.
988!-- Although data are contiguous in physical memory (which does not
989!-- necessarily require an MPI-derived datatype), the data exchange between
990!-- left and right PE's using the MPI-derived type is 10% faster than without.
991    CALL MPI_TYPE_VECTOR( nxr-nxl+1+2*nbgp, nbgp*(nzt-nzb+2), ngp_yz(0), &
992                          MPI_REAL, type_xz(0), ierr )
993    CALL MPI_TYPE_COMMIT( type_xz(0), ierr )
994
995    CALL MPI_TYPE_VECTOR( nbgp, ngp_yz(0), ngp_yz(0), MPI_REAL, type_yz(0), &
996                          ierr ) 
997    CALL MPI_TYPE_COMMIT( type_yz(0), ierr )
998
999!
1000!-- Definition of MPI-datatypes for multigrid method (coarser level grids)
1001    IF ( psolver == 'multigrid' )  THEN
1002!   
1003!--    Definition of MPI-datatyoe as above, but only 1 ghost level is used
1004       DO  i = maximum_grid_level, 1 , -1
1005
1006          ngp_yz(i) = (nzt_l - nzb_l + 2) * (nyn_l - nys_l + 3)
1007
1008          CALL MPI_TYPE_VECTOR( nxr_l-nxl_l+3, nzt_l-nzb_l+2, ngp_yz(i), &
1009                                MPI_REAL, type_xz(i), ierr )
1010          CALL MPI_TYPE_COMMIT( type_xz(i), ierr )
1011
1012          CALL MPI_TYPE_VECTOR( 1, ngp_yz(i), ngp_yz(i), MPI_REAL, type_yz(i), &
1013                                ierr )
1014          CALL MPI_TYPE_COMMIT( type_yz(i), ierr )
1015
1016          nxl_l = nxl_l / 2
1017          nxr_l = nxr_l / 2
1018          nys_l = nys_l / 2
1019          nyn_l = nyn_l / 2
1020          nzt_l = nzt_l / 2
1021
1022       ENDDO
1023
1024    ENDIF
1025#endif
1026
1027#if defined( __parallel ) && ! defined ( __check )
1028!
1029!-- Setting of flags for inflow/outflow conditions in case of non-cyclic
1030!-- horizontal boundary conditions.
1031    IF ( pleft == MPI_PROC_NULL )  THEN
1032       IF ( bc_lr == 'dirichlet/radiation' )  THEN
1033          inflow_l  = .TRUE.
1034       ELSEIF ( bc_lr == 'radiation/dirichlet' )  THEN
1035          outflow_l = .TRUE.
1036       ENDIF
1037    ENDIF
1038
1039    IF ( pright == MPI_PROC_NULL )  THEN
1040       IF ( bc_lr == 'dirichlet/radiation' )  THEN
1041          outflow_r = .TRUE.
1042       ELSEIF ( bc_lr == 'radiation/dirichlet' )  THEN
1043          inflow_r  = .TRUE.
1044       ENDIF
1045    ENDIF
1046
1047    IF ( psouth == MPI_PROC_NULL )  THEN
1048       IF ( bc_ns == 'dirichlet/radiation' )  THEN
1049          outflow_s = .TRUE.
1050       ELSEIF ( bc_ns == 'radiation/dirichlet' )  THEN
1051          inflow_s  = .TRUE.
1052       ENDIF
1053    ENDIF
1054
1055    IF ( pnorth == MPI_PROC_NULL )  THEN
1056       IF ( bc_ns == 'dirichlet/radiation' )  THEN
1057          inflow_n  = .TRUE.
1058       ELSEIF ( bc_ns == 'radiation/dirichlet' )  THEN
1059          outflow_n = .TRUE.
1060       ENDIF
1061    ENDIF
1062
1063!
1064!-- Broadcast the id of the inflow PE
1065    IF ( inflow_l )  THEN
1066       id_inflow_l = myidx
1067    ELSE
1068       id_inflow_l = 0
1069    ENDIF
1070    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
1071    CALL MPI_ALLREDUCE( id_inflow_l, id_inflow, 1, MPI_INTEGER, MPI_SUM, &
1072                        comm1dx, ierr )
1073
1074!
1075!-- Broadcast the id of the recycling plane
1076!-- WARNING: needs to be adjusted in case of inflows other than from left side!
1077    IF ( NINT( recycling_width / dx ) >= nxl  .AND. &
1078         NINT( recycling_width / dx ) <= nxr )  THEN
1079       id_recycling_l = myidx
1080    ELSE
1081       id_recycling_l = 0
1082    ENDIF
1083    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
1084    CALL MPI_ALLREDUCE( id_recycling_l, id_recycling, 1, MPI_INTEGER, MPI_SUM, &
1085                        comm1dx, ierr )
1086
1087    CALL location_message( 'finished', .TRUE. )
1088
1089#elif ! defined ( __parallel )
1090    IF ( bc_lr == 'dirichlet/radiation' )  THEN
1091       inflow_l  = .TRUE.
1092       outflow_r = .TRUE.
1093    ELSEIF ( bc_lr == 'radiation/dirichlet' )  THEN
1094       outflow_l = .TRUE.
1095       inflow_r  = .TRUE.
1096    ENDIF
1097
1098    IF ( bc_ns == 'dirichlet/radiation' )  THEN
1099       inflow_n  = .TRUE.
1100       outflow_s = .TRUE.
1101    ELSEIF ( bc_ns == 'radiation/dirichlet' )  THEN
1102       outflow_n = .TRUE.
1103       inflow_s  = .TRUE.
1104    ENDIF
1105#endif
1106
1107!
1108!-- At the inflow or outflow, u or v, respectively, have to be calculated for
1109!-- one more grid point.
1110    IF ( inflow_l .OR. outflow_l )  THEN
1111       nxlu = nxl + 1
1112    ELSE
1113       nxlu = nxl
1114    ENDIF
1115    IF ( inflow_s .OR. outflow_s )  THEN
1116       nysv = nys + 1
1117    ELSE
1118       nysv = nys
1119    ENDIF
1120
1121!
1122!-- Allocate wall flag arrays used in the multigrid solver
1123    IF ( psolver == 'multigrid' )  THEN
1124
1125       DO  i = maximum_grid_level, 1, -1
1126
1127           SELECT CASE ( i )
1128
1129              CASE ( 1 )
1130                 ALLOCATE( wall_flags_1(nzb:nzt_mg(i)+1,         &
1131                                        nys_mg(i)-1:nyn_mg(i)+1, &
1132                                        nxl_mg(i)-1:nxr_mg(i)+1) )
1133
1134              CASE ( 2 )
1135                 ALLOCATE( wall_flags_2(nzb:nzt_mg(i)+1,         &
1136                                        nys_mg(i)-1:nyn_mg(i)+1, &
1137                                        nxl_mg(i)-1:nxr_mg(i)+1) )
1138
1139              CASE ( 3 )
1140                 ALLOCATE( wall_flags_3(nzb:nzt_mg(i)+1,         &
1141                                        nys_mg(i)-1:nyn_mg(i)+1, &
1142                                        nxl_mg(i)-1:nxr_mg(i)+1) )
1143
1144              CASE ( 4 )
1145                 ALLOCATE( wall_flags_4(nzb:nzt_mg(i)+1,         &
1146                                        nys_mg(i)-1:nyn_mg(i)+1, &
1147                                        nxl_mg(i)-1:nxr_mg(i)+1) )
1148
1149              CASE ( 5 )
1150                 ALLOCATE( wall_flags_5(nzb:nzt_mg(i)+1,         &
1151                                        nys_mg(i)-1:nyn_mg(i)+1, &
1152                                        nxl_mg(i)-1:nxr_mg(i)+1) )
1153
1154              CASE ( 6 )
1155                 ALLOCATE( wall_flags_6(nzb:nzt_mg(i)+1,         &
1156                                        nys_mg(i)-1:nyn_mg(i)+1, &
1157                                        nxl_mg(i)-1:nxr_mg(i)+1) )
1158
1159              CASE ( 7 )
1160                 ALLOCATE( wall_flags_7(nzb:nzt_mg(i)+1,         &
1161                                        nys_mg(i)-1:nyn_mg(i)+1, &
1162                                        nxl_mg(i)-1:nxr_mg(i)+1) )
1163
1164              CASE ( 8 )
1165                 ALLOCATE( wall_flags_8(nzb:nzt_mg(i)+1,         &
1166                                        nys_mg(i)-1:nyn_mg(i)+1, &
1167                                        nxl_mg(i)-1:nxr_mg(i)+1) )
1168
1169              CASE ( 9 )
1170                 ALLOCATE( wall_flags_9(nzb:nzt_mg(i)+1,         &
1171                                        nys_mg(i)-1:nyn_mg(i)+1, &
1172                                        nxl_mg(i)-1:nxr_mg(i)+1) )
1173
1174              CASE ( 10 )
1175                 ALLOCATE( wall_flags_10(nzb:nzt_mg(i)+1,        &
1176                                        nys_mg(i)-1:nyn_mg(i)+1, &
1177                                        nxl_mg(i)-1:nxr_mg(i)+1) )
1178
1179              CASE DEFAULT
1180                 message_string = 'more than 10 multigrid levels'
1181                 CALL message( 'init_pegrid', 'PA0238', 1, 2, 0, 6, 0 )
1182
1183          END SELECT
1184
1185       ENDDO
1186
1187    ENDIF
1188
1189!
1190!-- Calculate the number of groups into which parallel I/O is split.
1191!-- The default for files which are opened by all PEs (or where each
1192!-- PE opens his own independent file) is, that all PEs are doing input/output
1193!-- in parallel at the same time. This might cause performance or even more
1194!-- severe problems depending on the configuration of the underlying file
1195!-- system.
1196!-- First, set the default:
1197    IF ( maximum_parallel_io_streams == -1  .OR. &
1198         maximum_parallel_io_streams > numprocs )  THEN
1199       maximum_parallel_io_streams = numprocs
1200    ENDIF
1201
1202!
1203!-- Now calculate the number of io_blocks and the io_group to which the
1204!-- respective PE belongs. I/O of the groups is done in serial, but in parallel
1205!-- for all PEs belonging to the same group. A preliminary setting with myid
1206!-- based on MPI_COMM_WORLD has been done in parin.
1207    io_blocks = numprocs / maximum_parallel_io_streams
1208    io_group  = MOD( myid+1, io_blocks )
1209   
1210
1211 END SUBROUTINE init_pegrid
Note: See TracBrowser for help on using the repository browser.