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

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