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

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

Check if grid coarsening is possible on subdomain, in order to avoid that multigrid approach effectively reduces to a Gauss-Seidel scheme.

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