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

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

in a nested run spend 3 ghost points also for 2nd-order pw-scheme

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