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

Last change on this file since 4402 was 4360, checked in by suehring, 5 years ago

Bugfix in output of time-averaged plant-canopy quanities; Output of plant-canopy data only where tall canopy is defined; land-surface model: fix wrong location strings; tests: update urban test case; all source code files: copyright update

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