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

Last change on this file since 1804 was 1804, checked in by maronga, 8 years ago

removed parameter file check. update of mrungui for compilation with qt5

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