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

Last change on this file since 3430 was 3347, checked in by suehring, 6 years ago

Offline nesting revised and separated from large_scale_forcing_mod; Time-dependent synthetic turbulence generator; bugfixes in USM/LSM radiation model and init_pegrid

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