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

Last change on this file since 4560 was 4461, checked in by raasch, 5 years ago

extensions to allow usage of alternative communicators in exchange_horiz

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