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

Last change on this file since 1833 was 1833, checked in by raasch, 8 years ago

spectrum renamed spactra_par and further modularized, POINTER-attributes added in coupler-routines to avoid gfortran error messages

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