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

Last change on this file since 3780 was 3761, checked in by raasch, 6 years ago

unused variables removed, OpenACC directives re-formatted, statements added to avoid compiler warnings

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