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

Last change on this file since 4180 was 4180, checked in by scharf, 5 years ago

removed comments in 'Former revisions' section that are older than 01.01.2019

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