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

Last change on this file since 4444 was 4444, checked in by raasch, 4 years ago

bugfix: cpp-directives for serial mode added

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