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

Last change on this file since 3884 was 3884, checked in by Giersch, 5 years ago

id_recycling is only calculated in case of tubulent inflow

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