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

Last change on this file since 1929 was 1923, checked in by boeske, 9 years ago

last commit documented

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