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

Last change on this file since 1965 was 1965, checked in by hellstea, 8 years ago

last commit documented

  • Property svn:keywords set to Id
File size: 44.0 KB
RevLine 
[1682]1!> @file init_pegrid.f90
[1036]2!--------------------------------------------------------------------------------!
3! This file is part of PALM.
4!
5! PALM is free software: you can redistribute it and/or modify it under the terms
6! of the GNU General Public License as published by the Free Software Foundation,
7! either version 3 of the License, or (at your option) any later version.
8!
9! PALM is distributed in the hope that it will be useful, but WITHOUT ANY
10! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR
11! A PARTICULAR PURPOSE.  See the GNU General Public License for more details.
12!
13! You should have received a copy of the GNU General Public License along with
14! PALM. If not, see <http://www.gnu.org/licenses/>.
15!
[1818]16! Copyright 1997-2016 Leibniz Universitaet Hannover
[1036]17!--------------------------------------------------------------------------------!
18!
[254]19! Current revisions:
[1322]20! ------------------
[1834]21!
[1965]22!
[1321]23! Former revisions:
24! -----------------
25! $Id: init_pegrid.f90 1965 2016-07-14 15:36:04Z hellstea $
26!
[1965]27! 1964 2016-07-14 15:35:18Z hellstea
28! Bugfix: erroneous setting of nest_bound_l/r/s/n = .TRUE. for vertical nesting mode removed.
29!
[1933]30! 1923 2016-05-31 16:37:07Z boeske
31! Initial version of purely vertical nesting introduced.
32!
[1923]33! 1922 2016-05-31 16:36:08Z boeske
34! Bugfix: array transposition checks restricted to cases if a fourier
35! transform is used , removed unused variable nnx_z
36!
[1834]37! 1833 2016-04-07 14:23:03Z raasch
38! spectra related variables moved to spectra_mod
39!
[1816]40! 1815 2016-04-06 13:49:59Z raasch
41! cpp-directives for intel openmp bug removed
42!
[1805]43! 1804 2016-04-05 16:30:18Z maronga
44! Removed code for parameter file check (__check)
45!
[1780]46! 1779 2016-03-03 08:01:28Z raasch
47! changes regarding nested domain removed: virtual PE grid will be automatically
48! calculated for nested runs too
49!
[1765]50! 1764 2016-02-28 12:45:19Z raasch
51! cpp-statements for nesting removed
52!
[1763]53! 1762 2016-02-25 12:31:13Z hellstea
54! Introduction of nested domain feature
55!
[1683]56! 1682 2015-10-07 23:56:08Z knoop
57! Code annotations made doxygen readable
58!
[1678]59! 1677 2015-10-02 13:25:23Z boeske
60! New MPI-data types for exchange of 3D integer arrays.
61!
[1576]62! 1575 2015-03-27 09:56:27Z raasch
63! adjustments for psolver-queries, calculation of ngp_xz added
64!
[1566]65! 1565 2015-03-09 20:59:31Z suehring
66! Refine if-clause for setting nbgp.
67!
[1558]68! 1557 2015-03-05 16:43:04Z suehring
69! Adjustment for monotonic limiter
70!
[1469]71! 1468 2014-09-24 14:06:57Z maronga
72! Adapted for use on up to 6-digit processor cores
73!
[1436]74! 1435 2014-07-21 10:37:02Z keck
75! bugfix: added missing parameter coupling_mode_remote to ONLY-attribute
76!
[1403]77! 1402 2014-05-09 14:25:13Z raasch
78! location messages modified
79!
[1385]80! 1384 2014-05-02 14:31:06Z raasch
81! location messages added
82!
[1354]83! 1353 2014-04-08 15:21:23Z heinze
84! REAL constants provided with KIND-attribute
85!
[1323]86! 1322 2014-03-20 16:38:49Z raasch
87! REAL functions provided with KIND-attribute
88!
[1321]89! 1320 2014-03-20 08:40:49Z raasch
[1320]90! ONLY-attribute added to USE-statements,
91! kind-parameters added to all INTEGER and REAL declaration statements,
92! kinds are defined in new module kinds,
93! revision history before 2012 removed,
94! comment fields (!:) to be used for variable explanations added to
95! all variable declaration statements
[760]96!
[1305]97! 1304 2014-03-12 10:29:42Z raasch
98! bugfix: single core MPI runs missed some settings of transpose indices
99!
[1213]100! 1212 2013-08-15 08:46:27Z raasch
101! error message for poisfft_hybrid removed
102!
[1160]103! 1159 2013-05-21 11:58:22Z fricke
104! dirichlet/neumann and neumann/dirichlet removed
105!
[1140]106! 1139 2013-04-18 07:25:03Z raasch
107! bugfix for calculating the id of the PE carrying the recycling plane
108!
[1112]109! 1111 2013-03-08 23:54:10Z raasch
110! initialization of poisfft moved to module poisfft
111!
[1093]112! 1092 2013-02-02 11:24:22Z raasch
113! unused variables removed
114!
[1057]115! 1056 2012-11-16 15:28:04Z raasch
116! Indices for arrays n.._mg start from zero due to definition of arrays f2 and
117! p2 as automatic arrays in recursive subroutine next_mg_level
118!
[1042]119! 1041 2012-11-06 02:36:29Z raasch
120! a 2d virtual processor topology is used by default for all machines
121!
[1037]122! 1036 2012-10-22 13:43:42Z raasch
123! code put under GPL (PALM 3.9)
124!
[1004]125! 1003 2012-09-14 14:35:53Z raasch
126! subdomains must have identical size (grid matching = "match" removed)
127!
[1002]128! 1001 2012-09-13 14:08:46Z raasch
129! all actions concerning upstream-spline-method removed
130!
[979]131! 978 2012-08-09 08:28:32Z fricke
132! dirichlet/neumann and neumann/dirichlet added
133! nxlu and nysv are also calculated for inflow boundary
134!
[810]135! 809 2012-01-30 13:32:58Z maronga
136! Bugfix: replaced .AND. and .NOT. with && and ! in the preprocessor directives
137!
[808]138! 807 2012-01-25 11:53:51Z maronga
139! New cpp directive "__check" implemented which is used by check_namelist_files
140!
[1]141! Revision 1.1  1997/07/24 11:15:09  raasch
142! Initial revision
143!
144!
145! Description:
146! ------------
[1682]147!> Determination of the virtual processor topology (if not prescribed by the
148!> user)and computation of the grid point number and array bounds of the local
149!> domains.
[1]150!------------------------------------------------------------------------------!
[1682]151 SUBROUTINE init_pegrid
152 
[1]153
[1320]154    USE control_parameters,                                                    &
[1435]155        ONLY:  bc_lr, bc_ns, coupling_mode, coupling_mode_remote,              &
[1833]156               coupling_topology, gathered_size, grid_level,                   &
[1435]157               grid_level_count, host, inflow_l, inflow_n, inflow_r, inflow_s, &
158               io_blocks, io_group, maximum_grid_level,                        &
159               maximum_parallel_io_streams, message_string,                    &
[1762]160               mg_switch_to_pe0_level, momentum_advec, nest_bound_l,           &
[1933]161               nest_bound_n, nest_bound_r, nest_bound_s, nest_domain, neutral, &
162               psolver, outflow_l, outflow_n, outflow_r, outflow_s,            &
163               recycling_width, scalar_advec, subdomain_size 
[1]164
[1320]165    USE grid_variables,                                                        &
166        ONLY:  dx
167       
168    USE indices,                                                               &
169        ONLY:  mg_loc_ind, nbgp, nnx, nny, nnz, nx, nx_a, nx_o, nxl, nxl_mg,   &
170               nxlu, nxr, nxr_mg, ny, ny_a, ny_o, nyn, nyn_mg, nys, nys_mg,    &
171               nysv, nz, nzb, nzt, nzt_mg, wall_flags_1, wall_flags_2,         &
172               wall_flags_3, wall_flags_4, wall_flags_5, wall_flags_6,         &
173               wall_flags_7, wall_flags_8, wall_flags_9, wall_flags_10
[1]174
[1320]175    USE kinds
176     
177    USE pegrid
[1933]178   
179    USE pmc_interface,                                                         &   
180        ONLY:  nesting_mode
181   
[1833]182    USE spectra_mod,                                                           &
[1922]183        ONLY:  calculate_spectra, dt_dosp
[1833]184
[1320]185    USE transpose_indices,                                                     &
186        ONLY:  nxl_y, nxl_yd, nxl_z, nxr_y, nxr_yd, nxr_z, nyn_x, nyn_z, nys_x,&
187               nys_z, nzb_x, nzb_y, nzb_yd, nzt_x, nzt_yd, nzt_y
[667]188
[1]189    IMPLICIT NONE
190
[1682]191    INTEGER(iwp) ::  i                        !<
192    INTEGER(iwp) ::  id_inflow_l              !<
193    INTEGER(iwp) ::  id_recycling_l           !<
194    INTEGER(iwp) ::  ind(5)                   !<
195    INTEGER(iwp) ::  j                        !<
196    INTEGER(iwp) ::  k                        !<
197    INTEGER(iwp) ::  maximum_grid_level_l     !<
198    INTEGER(iwp) ::  mg_levels_x              !<
199    INTEGER(iwp) ::  mg_levels_y              !<
200    INTEGER(iwp) ::  mg_levels_z              !<
201    INTEGER(iwp) ::  mg_switch_to_pe0_level_l !<
202    INTEGER(iwp) ::  nnx_y                    !<
203    INTEGER(iwp) ::  nnx_z                    !<
204    INTEGER(iwp) ::  nny_x                    !<
205    INTEGER(iwp) ::  nny_z                    !<
206    INTEGER(iwp) ::  nnz_x                    !<
207    INTEGER(iwp) ::  nnz_y                    !<
208    INTEGER(iwp) ::  numproc_sqr              !<
209    INTEGER(iwp) ::  nxl_l                    !<
210    INTEGER(iwp) ::  nxr_l                    !<
211    INTEGER(iwp) ::  nyn_l                    !<
212    INTEGER(iwp) ::  nys_l                    !<
213    INTEGER(iwp) ::  nzb_l                    !<
214    INTEGER(iwp) ::  nzt_l                    !<
215    INTEGER(iwp) ::  omp_get_num_threads      !<
[1]216
[1682]217    INTEGER(iwp), DIMENSION(:), ALLOCATABLE ::  ind_all !<
218    INTEGER(iwp), DIMENSION(:), ALLOCATABLE ::  nxlf    !<
219    INTEGER(iwp), DIMENSION(:), ALLOCATABLE ::  nxrf    !<
220    INTEGER(iwp), DIMENSION(:), ALLOCATABLE ::  nynf    !<
221    INTEGER(iwp), DIMENSION(:), ALLOCATABLE ::  nysf    !<
[1]222
[1682]223    INTEGER(iwp), DIMENSION(2) :: pdims_remote          !<
[667]224
[1092]225#if defined( __mpi2 )
[1682]226    LOGICAL ::  found                                   !<
[1092]227#endif
[1]228
229!
230!-- Get the number of OpenMP threads
231    !$OMP PARALLEL
232!$  threads_per_task = omp_get_num_threads()
233    !$OMP END PARALLEL
234
235
236#if defined( __parallel )
[667]237
[1402]238    CALL location_message( 'creating virtual PE grids + MPI derived data types', &
239                           .FALSE. )
[1764]240
[1]241!
[1764]242!--    Determine the processor topology or check it, if prescribed by the user
[1779]243    IF ( npex == -1  .AND.  npey == -1 )  THEN
[1]244
245!
[1764]246!--       Automatic determination of the topology
[1779]247       numproc_sqr = SQRT( REAL( numprocs, KIND=wp ) )
248       pdims(1)    = MAX( numproc_sqr , 1 )
249       DO  WHILE ( MOD( numprocs , pdims(1) ) /= 0 )
250          pdims(1) = pdims(1) - 1
251       ENDDO
252       pdims(2) = numprocs / pdims(1)
[1]253
[1779]254    ELSEIF ( npex /= -1  .AND.  npey /= -1 )  THEN
[1]255
256!
[1779]257!--    Prescribed by user. Number of processors on the prescribed topology
258!--    must be equal to the number of PEs available to the job
259       IF ( ( npex * npey ) /= numprocs )  THEN
260          WRITE( message_string, * ) 'number of PEs of the prescribed ',   &
261              'topology (', npex*npey,') does not match & the number of ', &
262              'PEs available to the job (', numprocs, ')'
263          CALL message( 'init_pegrid', 'PA0221', 1, 2, 0, 6, 0 )
264       ENDIF
265       pdims(1) = npex
266       pdims(2) = npey
[1]267
[1779]268    ELSE
[1]269!
[1779]270!--    If the processor topology is prescribed by the user, the number of
271!--    PEs must be given in both directions
272       message_string = 'if the processor topology is prescribed by th' //  &
273                'e user& both values of "npex" and "npey" must be given' // &
274                ' in the &NAMELIST-parameter file'
275       CALL message( 'init_pegrid', 'PA0222', 1, 2, 0, 6, 0 )
[1]276
277    ENDIF
278
279!
[622]280!-- For communication speedup, set barriers in front of collective
281!-- communications by default on SGI-type systems
282    IF ( host(3:5) == 'sgi' )  collective_wait = .TRUE.
283
284!
[1]285!-- If necessary, set horizontal boundary conditions to non-cyclic
[722]286    IF ( bc_lr /= 'cyclic' )  cyclic(1) = .FALSE.
287    IF ( bc_ns /= 'cyclic' )  cyclic(2) = .FALSE.
[1]288
[807]289
[1]290!
291!-- Create the virtual processor grid
292    CALL MPI_CART_CREATE( comm_palm, ndim, pdims, cyclic, reorder, &
293                          comm2d, ierr )
294    CALL MPI_COMM_RANK( comm2d, myid, ierr )
[1468]295    WRITE (myid_char,'(''_'',I6.6)')  myid
[1]296
297    CALL MPI_CART_COORDS( comm2d, myid, ndim, pcoord, ierr )
298    CALL MPI_CART_SHIFT( comm2d, 0, 1, pleft, pright, ierr )
299    CALL MPI_CART_SHIFT( comm2d, 1, 1, psouth, pnorth, ierr )
300
301!
302!-- Determine sub-topologies for transpositions
303!-- Transposition from z to x:
304    remain_dims(1) = .TRUE.
305    remain_dims(2) = .FALSE.
306    CALL MPI_CART_SUB( comm2d, remain_dims, comm1dx, ierr )
307    CALL MPI_COMM_RANK( comm1dx, myidx, ierr )
308!
309!-- Transposition from x to y
310    remain_dims(1) = .FALSE.
311    remain_dims(2) = .TRUE.
312    CALL MPI_CART_SUB( comm2d, remain_dims, comm1dy, ierr )
313    CALL MPI_COMM_RANK( comm1dy, myidy, ierr )
314
315
316!
[1003]317!-- Calculate array bounds along x-direction for every PE.
[1]318    ALLOCATE( nxlf(0:pdims(1)-1), nxrf(0:pdims(1)-1), nynf(0:pdims(2)-1), &
[1003]319              nysf(0:pdims(2)-1) )
[1]320
[1003]321    IF ( MOD( nx+1 , pdims(1) ) /= 0 )  THEN
[274]322       WRITE( message_string, * ) 'x-direction: gridpoint number (',nx+1,') ',&
323                               'is not an& integral divisor of the number ',  &
324                               'processors (', pdims(1),')'
[254]325       CALL message( 'init_pegrid', 'PA0225', 1, 2, 0, 6, 0 )
[1]326    ELSE
[1003]327       nnx  = ( nx + 1 ) / pdims(1)
[1]328       IF ( nnx*pdims(1) - ( nx + 1) > nnx )  THEN
[274]329          WRITE( message_string, * ) 'x-direction: nx does not match the',    & 
330                       'requirements given by the number of PEs &used',       &
331                       '& please use nx = ', nx - ( pdims(1) - ( nnx*pdims(1) &
332                                   - ( nx + 1 ) ) ), ' instead of nx =', nx
[254]333          CALL message( 'init_pegrid', 'PA0226', 1, 2, 0, 6, 0 )
[1]334       ENDIF
335    ENDIF   
336
337!
338!-- Left and right array bounds, number of gridpoints
339    DO  i = 0, pdims(1)-1
340       nxlf(i)   = i * nnx
341       nxrf(i)   = ( i + 1 ) * nnx - 1
342    ENDDO
343
344!
345!-- Calculate array bounds in y-direction for every PE.
[1003]346    IF ( MOD( ny+1 , pdims(2) ) /= 0 )  THEN
[274]347       WRITE( message_string, * ) 'y-direction: gridpoint number (',ny+1,') ', &
348                           'is not an& integral divisor of the number of',     &
349                           'processors (', pdims(2),')'
[254]350       CALL message( 'init_pegrid', 'PA0227', 1, 2, 0, 6, 0 )
[1]351    ELSE
[1003]352       nny  = ( ny + 1 ) / pdims(2)
[1]353       IF ( nny*pdims(2) - ( ny + 1) > nny )  THEN
[274]354          WRITE( message_string, * ) 'y-direction: ny does not match the',    &
355                       'requirements given by the number of PEs &used ',      &
356                       '& please use ny = ', ny - ( pdims(2) - ( nnx*pdims(2) &
[254]357                                     - ( ny + 1 ) ) ), ' instead of ny =', ny
358          CALL message( 'init_pegrid', 'PA0228', 1, 2, 0, 6, 0 )
[1]359       ENDIF
360    ENDIF   
361
362!
363!-- South and north array bounds
364    DO  j = 0, pdims(2)-1
365       nysf(j)   = j * nny
366       nynf(j)   = ( j + 1 ) * nny - 1
367    ENDDO
368
369!
370!-- Local array bounds of the respective PEs
[1003]371    nxl = nxlf(pcoord(1))
372    nxr = nxrf(pcoord(1))
373    nys = nysf(pcoord(2))
374    nyn = nynf(pcoord(2))
375    nzb = 0
376    nzt = nz
377    nnz = nz
[1]378
379!
[707]380!-- Set switches to define if the PE is situated at the border of the virtual
381!-- processor grid
382    IF ( nxl == 0 )   left_border_pe  = .TRUE.
383    IF ( nxr == nx )  right_border_pe = .TRUE.
384    IF ( nys == 0 )   south_border_pe = .TRUE.
385    IF ( nyn == ny )  north_border_pe = .TRUE.
386
387!
[1]388!-- Calculate array bounds and gridpoint numbers for the transposed arrays
389!-- (needed in the pressure solver)
390!-- For the transposed arrays, cyclic boundaries as well as top and bottom
391!-- boundaries are omitted, because they are obstructive to the transposition
392
393!
394!-- 1. transposition  z --> x
[1001]395!-- This transposition is not neccessary in case of a 1d-decomposition along x
[1922]396    IF ( psolver == 'poisfft'  .OR.  calculate_spectra )  THEN
[1304]397
[1922]398       IF ( pdims(2) /= 1 )  THEN
399          IF ( MOD( nz , pdims(1) ) /= 0 )  THEN
400             WRITE( message_string, * ) 'transposition z --> x:',              &
401                       '&nz=',nz,' is not an integral divisior of pdims(1)=',  &
[274]402                                                                   pdims(1)
[1922]403             CALL message( 'init_pegrid', 'PA0230', 1, 2, 0, 6, 0 )
404          ENDIF
[1]405       ENDIF
[1922]406
407       nys_x = nys
408       nyn_x = nyn
409       nny_x = nny
410       nnz_x = nz / pdims(1)
411       nzb_x = 1 + myidx * nnz_x
412       nzt_x = ( myidx + 1 ) * nnz_x
413       sendrecvcount_zx = nnx * nny * nnz_x
414
[1]415    ENDIF
416
[1922]417
418    IF ( psolver == 'poisfft' )  THEN 
[1]419!
[1922]420!--    2. transposition  x --> y
421       IF ( MOD( nx+1 , pdims(2) ) /= 0 )  THEN
422          WRITE( message_string, * ) 'transposition x --> y:',                 &
423                            '&nx+1=',nx+1,' is not an integral divisor of ',   &
424                            'pdims(2)=',pdims(2)
425          CALL message( 'init_pegrid', 'PA0231', 1, 2, 0, 6, 0 )
426       ENDIF
[1]427
[1922]428       nnz_y = nnz_x
429       nzb_y = nzb_x
430       nzt_y = nzt_x
431       nnx_y = (nx+1) / pdims(2)
432       nxl_y = myidy * nnx_y
433       nxr_y = ( myidy + 1 ) * nnx_y - 1
434       sendrecvcount_xy = nnx_y * nny_x * nnz_y
[1]435!
[1922]436!--    3. transposition  y --> z 
437!--    (ELSE:  x --> y  in case of 1D-decomposition along x)
438       nxl_z = nxl_y
439       nxr_z = nxr_y
440       nny_z = (ny+1) / pdims(1)
441       nys_z = myidx * nny_z
442       nyn_z = ( myidx + 1 ) * nny_z - 1
443       sendrecvcount_yz = nnx_y * nny_z * nnz_y
[1304]444
[1922]445       IF ( pdims(2) /= 1 )  THEN
[1]446!
[1922]447!--       y --> z
448!--       This transposition is not neccessary in case of a 1d-decomposition
449!--       along x, except that the uptream-spline method is switched on
450          IF ( MOD( ny+1 , pdims(1) ) /= 0 )  THEN
451             WRITE( message_string, * ) 'transposition y --> z:',              &
452                               '& ny+1=',ny+1,' is not an integral divisor of',&
453                               ' pdims(1)=',pdims(1)
454             CALL message( 'init_pegrid', 'PA0232', 1, 2, 0, 6, 0 )
455          ENDIF
[1]456
[1922]457       ELSE
[1]458!
[1922]459!--       x --> y
460!--       This condition must be fulfilled for a 1D-decomposition along x
461          IF ( MOD( ny+1 , pdims(1) ) /= 0 )  THEN
462             WRITE( message_string, * ) 'transposition x --> y:',              &
463                               '& ny+1=',ny+1,' is not an integral divisor of',&
464                               ' pdims(1)=',pdims(1)
465             CALL message( 'init_pegrid', 'PA0233', 1, 2, 0, 6, 0 )
466          ENDIF
467
[1]468       ENDIF
469
470    ENDIF
471
472!
473!-- Indices for direct transpositions z --> y (used for calculating spectra)
[1922]474    IF ( calculate_spectra )  THEN
[1003]475       IF ( MOD( nz, pdims(2) ) /= 0 )  THEN
[1922]476          WRITE( message_string, * ) 'direct transposition z --> y (needed ',  &
477                    'for spectra):& nz=',nz,' is not an integral divisor of ', &
[274]478                    'pdims(2)=',pdims(2)
[254]479          CALL message( 'init_pegrid', 'PA0234', 1, 2, 0, 6, 0 )
[1]480       ELSE
[1003]481          nxl_yd = nxl
482          nxr_yd = nxr
483          nzb_yd = 1 + myidy * ( nz / pdims(2) )
484          nzt_yd = ( myidy + 1 ) * ( nz / pdims(2) )
485          sendrecvcount_zyd = nnx * nny * ( nz / pdims(2) )
[1]486       ENDIF
487    ENDIF
488
[1922]489    IF ( psolver == 'poisfft'  .OR.  calculate_spectra )  THEN
[1]490!
[1922]491!--    Indices for direct transpositions y --> x
492!--    (they are only possible in case of a 1d-decomposition along x)
493       IF ( pdims(2) == 1 )  THEN
494          nny_x = nny / pdims(1)
495          nys_x = myid * nny_x
496          nyn_x = ( myid + 1 ) * nny_x - 1
497          nzb_x = 1
498          nzt_x = nz
499          sendrecvcount_xy = nnx * nny_x * nz
500       ENDIF
501
[1]502    ENDIF
503
[1922]504    IF ( psolver == 'poisfft' )  THEN
[1]505!
[1922]506!--    Indices for direct transpositions x --> y
507!--    (they are only possible in case of a 1d-decomposition along y)
508       IF ( pdims(1) == 1 )  THEN
509          nnx_y = nnx / pdims(2)
510          nxl_y = myid * nnx_y
511          nxr_y = ( myid + 1 ) * nnx_y - 1
512          nzb_y = 1
513          nzt_y = nz
514          sendrecvcount_xy = nnx_y * nny * nz
515       ENDIF
516
[1]517    ENDIF
518
519!
520!-- Arrays for storing the array bounds are needed any more
521    DEALLOCATE( nxlf , nxrf , nynf , nysf )
522
[807]523
[145]524!
525!-- Collect index bounds from other PEs (to be written to restart file later)
526    ALLOCATE( hor_index_bounds(4,0:numprocs-1) )
527
528    IF ( myid == 0 )  THEN
529
530       hor_index_bounds(1,0) = nxl
531       hor_index_bounds(2,0) = nxr
532       hor_index_bounds(3,0) = nys
533       hor_index_bounds(4,0) = nyn
534
535!
536!--    Receive data from all other PEs
537       DO  i = 1, numprocs-1
538          CALL MPI_RECV( ibuf, 4, MPI_INTEGER, i, MPI_ANY_TAG, comm2d, status, &
539                         ierr )
540          hor_index_bounds(:,i) = ibuf(1:4)
541       ENDDO
542
543    ELSE
544!
545!--    Send index bounds to PE0
546       ibuf(1) = nxl
547       ibuf(2) = nxr
548       ibuf(3) = nys
549       ibuf(4) = nyn
550       CALL MPI_SEND( ibuf, 4, MPI_INTEGER, 0, myid, comm2d, ierr )
551
552    ENDIF
553
[807]554
[1]555#if defined( __print )
556!
557!-- Control output
558    IF ( myid == 0 )  THEN
559       PRINT*, '*** processor topology ***'
560       PRINT*, ' '
561       PRINT*, 'myid   pcoord    left right  south north  idx idy   nxl: nxr',&
562               &'   nys: nyn'
563       PRINT*, '------------------------------------------------------------',&
564               &'-----------'
565       WRITE (*,1000)  0, pcoord(1), pcoord(2), pleft, pright, psouth, pnorth, &
566                       myidx, myidy, nxl, nxr, nys, nyn
5671000   FORMAT (I4,2X,'(',I3,',',I3,')',3X,I4,2X,I4,3X,I4,2X,I4,2X,I3,1X,I3, &
568               2(2X,I4,':',I4))
569
570!
[108]571!--    Receive data from the other PEs
[1]572       DO  i = 1,numprocs-1
573          CALL MPI_RECV( ibuf, 12, MPI_INTEGER, i, MPI_ANY_TAG, comm2d, status, &
574                         ierr )
575          WRITE (*,1000)  i, ( ibuf(j) , j = 1,12 )
576       ENDDO
577    ELSE
578
579!
580!--    Send data to PE0
581       ibuf(1) = pcoord(1); ibuf(2) = pcoord(2); ibuf(3) = pleft
582       ibuf(4) = pright; ibuf(5) = psouth; ibuf(6) = pnorth; ibuf(7) = myidx
583       ibuf(8) = myidy; ibuf(9) = nxl; ibuf(10) = nxr; ibuf(11) = nys
584       ibuf(12) = nyn
585       CALL MPI_SEND( ibuf, 12, MPI_INTEGER, 0, myid, comm2d, ierr )       
586    ENDIF
587#endif
588
[1804]589#if defined( __parallel )
[102]590#if defined( __mpi2 )
591!
592!-- In case of coupled runs, get the port name on PE0 of the atmosphere model
593!-- and pass it to PE0 of the ocean model
594    IF ( myid == 0 )  THEN
595
596       IF ( coupling_mode == 'atmosphere_to_ocean' )  THEN
597
598          CALL MPI_OPEN_PORT( MPI_INFO_NULL, port_name, ierr )
[108]599
[102]600          CALL MPI_PUBLISH_NAME( 'palm_coupler', MPI_INFO_NULL, port_name, &
601                                 ierr )
[108]602
603!
[104]604!--       Write a flag file for the ocean model and the other atmosphere
605!--       processes.
606!--       There seems to be a bug in MPICH2 which causes hanging processes
607!--       in case that execution of LOOKUP_NAME is continued too early
608!--       (i.e. before the port has been created)
609          OPEN( 90, FILE='COUPLING_PORT_OPENED', FORM='FORMATTED' )
610          WRITE ( 90, '(''TRUE'')' )
611          CLOSE ( 90 )
[102]612
613       ELSEIF ( coupling_mode == 'ocean_to_atmosphere' )  THEN
614
[104]615!
616!--       Continue only if the atmosphere model has created the port.
617!--       There seems to be a bug in MPICH2 which causes hanging processes
618!--       in case that execution of LOOKUP_NAME is continued too early
619!--       (i.e. before the port has been created)
620          INQUIRE( FILE='COUPLING_PORT_OPENED', EXIST=found )
621          DO WHILE ( .NOT. found )
622             INQUIRE( FILE='COUPLING_PORT_OPENED', EXIST=found )
623          ENDDO
624
[102]625          CALL MPI_LOOKUP_NAME( 'palm_coupler', MPI_INFO_NULL, port_name, ierr )
626
627       ENDIF
628
629    ENDIF
630
631!
632!-- In case of coupled runs, establish the connection between the atmosphere
633!-- and the ocean model and define the intercommunicator (comm_inter)
634    CALL MPI_BARRIER( comm2d, ierr )
635    IF ( coupling_mode == 'atmosphere_to_ocean' )  THEN
636
637       CALL MPI_COMM_ACCEPT( port_name, MPI_INFO_NULL, 0, MPI_COMM_WORLD, &
638                             comm_inter, ierr )
[108]639       coupling_mode_remote = 'ocean_to_atmosphere'
640
[102]641    ELSEIF ( coupling_mode == 'ocean_to_atmosphere' )  THEN
642
643       CALL MPI_COMM_CONNECT( port_name, MPI_INFO_NULL, 0, MPI_COMM_WORLD, &
644                              comm_inter, ierr )
[108]645       coupling_mode_remote = 'atmosphere_to_ocean'
646
[102]647    ENDIF
[206]648#endif
[102]649
[667]650!
[709]651!-- Determine the number of ghost point layers
[1565]652    IF ( ( scalar_advec == 'ws-scheme' .AND. .NOT. neutral ) .OR.             &
[1557]653         scalar_advec == 'ws-scheme-mono' .OR.                                &
654         momentum_advec == 'ws-scheme' )  THEN
[667]655       nbgp = 3
656    ELSE
657       nbgp = 1
[709]658    ENDIF 
[667]659
[102]660!
[709]661!-- Create a new MPI derived datatype for the exchange of surface (xy) data,
662!-- which is needed for coupled atmosphere-ocean runs.
663!-- First, calculate number of grid points of an xy-plane.
[667]664    ngp_xy  = ( nxr - nxl + 1 + 2 * nbgp ) * ( nyn - nys + 1 + 2 * nbgp )
[102]665    CALL MPI_TYPE_VECTOR( ngp_xy, 1, nzt-nzb+2, MPI_REAL, type_xy, ierr )
666    CALL MPI_TYPE_COMMIT( type_xy, ierr )
[667]667
[709]668    IF ( TRIM( coupling_mode ) /= 'uncoupled' )  THEN
[667]669   
670!
671!--    Pass the number of grid points of the atmosphere model to
672!--    the ocean model and vice versa
673       IF ( coupling_mode == 'atmosphere_to_ocean' )  THEN
674
675          nx_a = nx
676          ny_a = ny
677
[709]678          IF ( myid == 0 )  THEN
679
680             CALL MPI_SEND( nx_a, 1, MPI_INTEGER, numprocs, 1, comm_inter,  &
681                            ierr )
682             CALL MPI_SEND( ny_a, 1, MPI_INTEGER, numprocs, 2, comm_inter,  &
683                            ierr )
684             CALL MPI_SEND( pdims, 2, MPI_INTEGER, numprocs, 3, comm_inter, &
685                            ierr )
686             CALL MPI_RECV( nx_o, 1, MPI_INTEGER, numprocs, 4, comm_inter,  &
687                            status, ierr )
688             CALL MPI_RECV( ny_o, 1, MPI_INTEGER, numprocs, 5, comm_inter,  &
689                            status, ierr )
690             CALL MPI_RECV( pdims_remote, 2, MPI_INTEGER, numprocs, 6,      &
[667]691                            comm_inter, status, ierr )
692          ENDIF
693
[709]694          CALL MPI_BCAST( nx_o, 1, MPI_INTEGER, 0, comm2d, ierr )
695          CALL MPI_BCAST( ny_o, 1, MPI_INTEGER, 0, comm2d, ierr ) 
696          CALL MPI_BCAST( pdims_remote, 2, MPI_INTEGER, 0, comm2d, ierr )
[667]697       
698       ELSEIF ( coupling_mode == 'ocean_to_atmosphere' )  THEN
699
700          nx_o = nx
701          ny_o = ny 
702
703          IF ( myid == 0 ) THEN
[709]704
705             CALL MPI_RECV( nx_a, 1, MPI_INTEGER, 0, 1, comm_inter, status, &
706                            ierr )
707             CALL MPI_RECV( ny_a, 1, MPI_INTEGER, 0, 2, comm_inter, status, &
708                            ierr )
709             CALL MPI_RECV( pdims_remote, 2, MPI_INTEGER, 0, 3, comm_inter, &
710                            status, ierr )
711             CALL MPI_SEND( nx_o, 1, MPI_INTEGER, 0, 4, comm_inter, ierr )
712             CALL MPI_SEND( ny_o, 1, MPI_INTEGER, 0, 5, comm_inter, ierr )
713             CALL MPI_SEND( pdims, 2, MPI_INTEGER, 0, 6, comm_inter, ierr )
[667]714          ENDIF
715
716          CALL MPI_BCAST( nx_a, 1, MPI_INTEGER, 0, comm2d, ierr)
717          CALL MPI_BCAST( ny_a, 1, MPI_INTEGER, 0, comm2d, ierr) 
718          CALL MPI_BCAST( pdims_remote, 2, MPI_INTEGER, 0, comm2d, ierr) 
719
720       ENDIF
721 
[709]722       ngp_a = ( nx_a+1 + 2 * nbgp ) * ( ny_a+1 + 2 * nbgp )
723       ngp_o = ( nx_o+1 + 2 * nbgp ) * ( ny_o+1 + 2 * nbgp )
[667]724
725!
[709]726!--    Determine if the horizontal grid and the number of PEs in ocean and
727!--    atmosphere is same or not
728       IF ( nx_o == nx_a  .AND.  ny_o == ny_a  .AND.  &
[667]729            pdims(1) == pdims_remote(1) .AND. pdims(2) == pdims_remote(2) ) &
730       THEN
731          coupling_topology = 0
732       ELSE
733          coupling_topology = 1
734       ENDIF 
735
736!
737!--    Determine the target PEs for the exchange between ocean and
738!--    atmosphere (comm2d)
[709]739       IF ( coupling_topology == 0 )  THEN
740!
741!--       In case of identical topologies, every atmosphere PE has exactly one
742!--       ocean PE counterpart and vice versa
743          IF ( TRIM( coupling_mode ) == 'atmosphere_to_ocean' ) THEN
[667]744             target_id = myid + numprocs
745          ELSE
746             target_id = myid 
747          ENDIF
748
749       ELSE
750!
751!--       In case of nonequivalent topology in ocean and atmosphere only for
752!--       PE0 in ocean and PE0 in atmosphere a target_id is needed, since
[709]753!--       data echxchange between ocean and atmosphere will be done only
754!--       between these PEs.   
755          IF ( myid == 0 )  THEN
756
757             IF ( TRIM( coupling_mode ) == 'atmosphere_to_ocean' )  THEN
[667]758                target_id = numprocs 
759             ELSE
760                target_id = 0
761             ENDIF
[709]762
[667]763          ENDIF
[709]764
[667]765       ENDIF
766
767    ENDIF
768
769
[102]770#endif
771
[1]772#else
773
774!
775!-- Array bounds when running on a single PE (respectively a non-parallel
776!-- machine)
[1003]777    nxl = 0
778    nxr = nx
779    nnx = nxr - nxl + 1
780    nys = 0
781    nyn = ny
782    nny = nyn - nys + 1
783    nzb = 0
784    nzt = nz
785    nnz = nz
[1]786
[145]787    ALLOCATE( hor_index_bounds(4,0:0) )
788    hor_index_bounds(1,0) = nxl
789    hor_index_bounds(2,0) = nxr
790    hor_index_bounds(3,0) = nys
791    hor_index_bounds(4,0) = nyn
792
[1]793!
794!-- Array bounds for the pressure solver (in the parallel code, these bounds
795!-- are the ones for the transposed arrays)
[1003]796    nys_x = nys
797    nyn_x = nyn
798    nzb_x = nzb + 1
799    nzt_x = nzt
[1]800
[1003]801    nxl_y = nxl
802    nxr_y = nxr
803    nzb_y = nzb + 1
804    nzt_y = nzt
[1]805
[1003]806    nxl_z = nxl
807    nxr_z = nxr
808    nys_z = nys
809    nyn_z = nyn
[1]810
811#endif
812
813!
814!-- Calculate number of grid levels necessary for the multigrid poisson solver
815!-- as well as the gridpoint indices on each level
[1575]816    IF ( psolver(1:9) == 'multigrid' )  THEN
[1]817
818!
819!--    First calculate number of possible grid levels for the subdomains
820       mg_levels_x = 1
821       mg_levels_y = 1
822       mg_levels_z = 1
823
824       i = nnx
825       DO WHILE ( MOD( i, 2 ) == 0  .AND.  i /= 2 )
826          i = i / 2
827          mg_levels_x = mg_levels_x + 1
828       ENDDO
829
830       j = nny
831       DO WHILE ( MOD( j, 2 ) == 0  .AND.  j /= 2 )
832          j = j / 2
833          mg_levels_y = mg_levels_y + 1
834       ENDDO
835
[181]836       k = nz    ! do not use nnz because it might be > nz due to transposition
837                 ! requirements
[1]838       DO WHILE ( MOD( k, 2 ) == 0  .AND.  k /= 2 )
839          k = k / 2
840          mg_levels_z = mg_levels_z + 1
841       ENDDO
842
843       maximum_grid_level = MIN( mg_levels_x, mg_levels_y, mg_levels_z )
844
845!
846!--    Find out, if the total domain allows more levels. These additional
[709]847!--    levels are identically processed on all PEs.
[197]848       IF ( numprocs > 1  .AND.  mg_switch_to_pe0_level /= -1 )  THEN
[709]849
[1]850          IF ( mg_levels_z > MIN( mg_levels_x, mg_levels_y ) )  THEN
[709]851
[1]852             mg_switch_to_pe0_level_l = maximum_grid_level
853
854             mg_levels_x = 1
855             mg_levels_y = 1
856
857             i = nx+1
858             DO WHILE ( MOD( i, 2 ) == 0  .AND.  i /= 2 )
859                i = i / 2
860                mg_levels_x = mg_levels_x + 1
861             ENDDO
862
863             j = ny+1
864             DO WHILE ( MOD( j, 2 ) == 0  .AND.  j /= 2 )
865                j = j / 2
866                mg_levels_y = mg_levels_y + 1
867             ENDDO
868
869             maximum_grid_level_l = MIN( mg_levels_x, mg_levels_y, mg_levels_z )
870
871             IF ( maximum_grid_level_l > mg_switch_to_pe0_level_l )  THEN
872                mg_switch_to_pe0_level_l = maximum_grid_level_l - &
873                                           mg_switch_to_pe0_level_l + 1
874             ELSE
875                mg_switch_to_pe0_level_l = 0
876             ENDIF
[709]877
[1]878          ELSE
879             mg_switch_to_pe0_level_l = 0
880             maximum_grid_level_l = maximum_grid_level
[709]881
[1]882          ENDIF
883
884!
885!--       Use switch level calculated above only if it is not pre-defined
886!--       by user
887          IF ( mg_switch_to_pe0_level == 0 )  THEN
888             IF ( mg_switch_to_pe0_level_l /= 0 )  THEN
889                mg_switch_to_pe0_level = mg_switch_to_pe0_level_l
890                maximum_grid_level     = maximum_grid_level_l
891             ENDIF
892
893          ELSE
894!
895!--          Check pre-defined value and reset to default, if neccessary
896             IF ( mg_switch_to_pe0_level < mg_switch_to_pe0_level_l  .OR.  &
897                  mg_switch_to_pe0_level >= maximum_grid_level_l )  THEN
[254]898                message_string = 'mg_switch_to_pe0_level ' // &
899                                 'out of range and reset to default (=0)'
900                CALL message( 'init_pegrid', 'PA0235', 0, 1, 0, 6, 0 )
[1]901                mg_switch_to_pe0_level = 0
902             ELSE
903!
904!--             Use the largest number of possible levels anyway and recalculate
905!--             the switch level to this largest number of possible values
906                maximum_grid_level = maximum_grid_level_l
907
908             ENDIF
[709]909
[1]910          ENDIF
911
912       ENDIF
913
[1056]914       ALLOCATE( grid_level_count(maximum_grid_level),                       &
915                 nxl_mg(0:maximum_grid_level), nxr_mg(0:maximum_grid_level), &
916                 nyn_mg(0:maximum_grid_level), nys_mg(0:maximum_grid_level), &
917                 nzt_mg(0:maximum_grid_level) )
[1]918
919       grid_level_count = 0
[1056]920!
921!--    Index zero required as dummy due to definition of arrays f2 and p2 in
922!--    recursive subroutine next_mg_level
923       nxl_mg(0) = 0; nxr_mg(0) = 0; nyn_mg(0) = 0; nys_mg(0) = 0; nzt_mg(0) = 0
[778]924
[1]925       nxl_l = nxl; nxr_l = nxr; nys_l = nys; nyn_l = nyn; nzt_l = nzt
926
927       DO  i = maximum_grid_level, 1 , -1
928
929          IF ( i == mg_switch_to_pe0_level )  THEN
[1804]930#if defined( __parallel )
[1]931!
932!--          Save the grid size of the subdomain at the switch level, because
933!--          it is needed in poismg.
934             ind(1) = nxl_l; ind(2) = nxr_l
935             ind(3) = nys_l; ind(4) = nyn_l
936             ind(5) = nzt_l
937             ALLOCATE( ind_all(5*numprocs), mg_loc_ind(5,0:numprocs-1) )
938             CALL MPI_ALLGATHER( ind, 5, MPI_INTEGER, ind_all, 5, &
939                                 MPI_INTEGER, comm2d, ierr )
940             DO  j = 0, numprocs-1
941                DO  k = 1, 5
942                   mg_loc_ind(k,j) = ind_all(k+j*5)
943                ENDDO
944             ENDDO
945             DEALLOCATE( ind_all )
946!
[709]947!--          Calculate the grid size of the total domain
[1]948             nxr_l = ( nxr_l-nxl_l+1 ) * pdims(1) - 1
949             nxl_l = 0
950             nyn_l = ( nyn_l-nys_l+1 ) * pdims(2) - 1
951             nys_l = 0
952!
953!--          The size of this gathered array must not be larger than the
954!--          array tend, which is used in the multigrid scheme as a temporary
[778]955!--          array. Therefore the subdomain size of an PE is calculated and
956!--          the size of the gathered grid. These values are used in 
957!--          routines pres and poismg
958             subdomain_size = ( nxr - nxl + 2 * nbgp + 1 ) * &
959                              ( nyn - nys + 2 * nbgp + 1 ) * ( nzt - nzb + 2 )
[1]960             gathered_size  = ( nxr_l - nxl_l + 3 ) * ( nyn_l - nys_l + 3 ) * &
961                              ( nzt_l - nzb + 2 )
962
[1804]963#else
[254]964             message_string = 'multigrid gather/scatter impossible ' // &
[1]965                          'in non parallel mode'
[254]966             CALL message( 'init_pegrid', 'PA0237', 1, 2, 0, 6, 0 )
[1]967#endif
968          ENDIF
969
970          nxl_mg(i) = nxl_l
971          nxr_mg(i) = nxr_l
972          nys_mg(i) = nys_l
973          nyn_mg(i) = nyn_l
974          nzt_mg(i) = nzt_l
975
976          nxl_l = nxl_l / 2 
977          nxr_l = nxr_l / 2
978          nys_l = nys_l / 2 
979          nyn_l = nyn_l / 2 
980          nzt_l = nzt_l / 2 
[778]981
[1]982       ENDDO
983
[780]984!
985!--    Temporary problem: Currently calculation of maxerror iin routine poismg crashes
986!--    if grid data are collected on PE0 already on the finest grid level.
987!--    To be solved later.
988       IF ( maximum_grid_level == mg_switch_to_pe0_level )  THEN
989          message_string = 'grid coarsening on subdomain level cannot be performed'
990          CALL message( 'poismg', 'PA0236', 1, 2, 0, 6, 0 )
991       ENDIF
992
[1]993    ELSE
994
[667]995       maximum_grid_level = 0
[1]996
997    ENDIF
998
[722]999!
1000!-- Default level 0 tells exchange_horiz that all ghost planes have to be
1001!-- exchanged. grid_level is adjusted in poismg, where only one ghost plane
1002!-- is required.
1003    grid_level = 0
[1]1004
[1804]1005#if defined( __parallel )
[1]1006!
1007!-- Gridpoint number for the exchange of ghost points (y-line for 2D-arrays)
[667]1008    ngp_y  = nyn - nys + 1 + 2 * nbgp
[1]1009
1010!
[709]1011!-- Define new MPI derived datatypes for the exchange of ghost points in
1012!-- x- and y-direction for 2D-arrays (line)
1013    CALL MPI_TYPE_VECTOR( nxr-nxl+1+2*nbgp, nbgp, ngp_y, MPI_REAL, type_x, &
1014                          ierr )
[1]1015    CALL MPI_TYPE_COMMIT( type_x, ierr )
[709]1016    CALL MPI_TYPE_VECTOR( nxr-nxl+1+2*nbgp, nbgp, ngp_y, MPI_INTEGER, &
1017                          type_x_int, ierr )
[1]1018    CALL MPI_TYPE_COMMIT( type_x_int, ierr )
1019
[667]1020    CALL MPI_TYPE_VECTOR( nbgp, ngp_y, ngp_y, MPI_REAL, type_y, ierr )
1021    CALL MPI_TYPE_COMMIT( type_y, ierr )
1022    CALL MPI_TYPE_VECTOR( nbgp, ngp_y, ngp_y, MPI_INTEGER, type_y_int, ierr )
1023    CALL MPI_TYPE_COMMIT( type_y_int, ierr )
1024
1025
[1]1026!
1027!-- Calculate gridpoint numbers for the exchange of ghost points along x
1028!-- (yz-plane for 3D-arrays) and define MPI derived data type(s) for the
1029!-- exchange of ghost points in y-direction (xz-plane).
1030!-- Do these calculations for the model grid and (if necessary) also
1031!-- for the coarser grid levels used in the multigrid method
[1575]1032    ALLOCATE ( ngp_xz(0:maximum_grid_level), ngp_yz(0:maximum_grid_level),     &
1033               type_xz(0:maximum_grid_level), type_yz(0:maximum_grid_level) )
[1]1034
1035    nxl_l = nxl; nxr_l = nxr; nys_l = nys; nyn_l = nyn; nzb_l = nzb; nzt_l = nzt
[709]1036
[667]1037!
1038!-- Discern between the model grid, which needs nbgp ghost points and
1039!-- grid levels for the multigrid scheme. In the latter case only one
1040!-- ghost point is necessary.
[709]1041!-- First definition of MPI-datatypes for exchange of ghost layers on normal
[667]1042!-- grid. The following loop is needed for data exchange in poismg.f90.
1043!
1044!-- Determine number of grid points of yz-layer for exchange
1045    ngp_yz(0) = (nzt - nzb + 2) * (nyn - nys + 1 + 2 * nbgp)
[709]1046
[667]1047!
[709]1048!-- Define an MPI-datatype for the exchange of left/right boundaries.
1049!-- Although data are contiguous in physical memory (which does not
1050!-- necessarily require an MPI-derived datatype), the data exchange between
1051!-- left and right PE's using the MPI-derived type is 10% faster than without.
[667]1052    CALL MPI_TYPE_VECTOR( nxr-nxl+1+2*nbgp, nbgp*(nzt-nzb+2), ngp_yz(0), &
[709]1053                          MPI_REAL, type_xz(0), ierr )
[667]1054    CALL MPI_TYPE_COMMIT( type_xz(0), ierr )
[1]1055
[709]1056    CALL MPI_TYPE_VECTOR( nbgp, ngp_yz(0), ngp_yz(0), MPI_REAL, type_yz(0), &
1057                          ierr ) 
[667]1058    CALL MPI_TYPE_COMMIT( type_yz(0), ierr )
[709]1059
[667]1060!
[709]1061!-- Definition of MPI-datatypes for multigrid method (coarser level grids)
[1575]1062    IF ( psolver(1:9) == 'multigrid' )  THEN
[667]1063!   
[709]1064!--    Definition of MPI-datatyoe as above, but only 1 ghost level is used
1065       DO  i = maximum_grid_level, 1 , -1
1066
[1575]1067          ngp_xz(i) = (nzt_l - nzb_l + 2) * (nxr_l - nxl_l + 3)
[667]1068          ngp_yz(i) = (nzt_l - nzb_l + 2) * (nyn_l - nys_l + 3)
1069
1070          CALL MPI_TYPE_VECTOR( nxr_l-nxl_l+3, nzt_l-nzb_l+2, ngp_yz(i), &
[709]1071                                MPI_REAL, type_xz(i), ierr )
[667]1072          CALL MPI_TYPE_COMMIT( type_xz(i), ierr )
[1]1073
[709]1074          CALL MPI_TYPE_VECTOR( 1, ngp_yz(i), ngp_yz(i), MPI_REAL, type_yz(i), &
1075                                ierr )
[667]1076          CALL MPI_TYPE_COMMIT( type_yz(i), ierr )
1077
1078          nxl_l = nxl_l / 2
1079          nxr_l = nxr_l / 2
1080          nys_l = nys_l / 2
1081          nyn_l = nyn_l / 2
1082          nzt_l = nzt_l / 2
[709]1083
[667]1084       ENDDO
[709]1085
1086    ENDIF
[1677]1087!
1088!-- Define data types for exchange of 3D Integer arrays.
1089    ngp_yz_int = (nzt - nzb + 2) * (nyn - nys + 1 + 2 * nbgp)
1090
1091    CALL MPI_TYPE_VECTOR( nxr-nxl+1+2*nbgp, nbgp*(nzt-nzb+2), ngp_yz_int, &
1092                          MPI_INTEGER, type_xz_int, ierr )
1093    CALL MPI_TYPE_COMMIT( type_xz_int, ierr )
1094
1095    CALL MPI_TYPE_VECTOR( nbgp, ngp_yz_int, ngp_yz_int, MPI_INTEGER, type_yz_int, &
1096                          ierr )
1097    CALL MPI_TYPE_COMMIT( type_yz_int, ierr )
1098
[1]1099#endif
1100
[1804]1101#if defined( __parallel )
[1]1102!
[1933]1103!-- Setting of flags for inflow/outflow/nesting conditions.
[1964]1104    IF ( pleft == MPI_PROC_NULL )  THEN
1105       IF ( bc_lr == 'dirichlet/radiation' )  THEN
1106          inflow_l  = .TRUE.
1107       ELSEIF ( bc_lr == 'radiation/dirichlet' )  THEN
1108          outflow_l = .TRUE.
1109       ELSEIF ( bc_lr == 'nested' )  THEN
[1762]1110          nest_bound_l = .TRUE.
[1]1111       ENDIF
1112    ENDIF
[1933]1113 
[1964]1114    IF ( pright == MPI_PROC_NULL )  THEN
1115       IF ( bc_lr == 'dirichlet/radiation' )  THEN
1116          outflow_r = .TRUE.
1117       ELSEIF ( bc_lr == 'radiation/dirichlet' )  THEN
1118          inflow_r  = .TRUE.
1119       ELSEIF ( bc_lr == 'nested' )  THEN
[1762]1120          nest_bound_r = .TRUE.
[1]1121       ENDIF
1122    ENDIF
1123
[1964]1124    IF ( psouth == MPI_PROC_NULL )  THEN
1125       IF ( bc_ns == 'dirichlet/radiation' )  THEN
1126          outflow_s = .TRUE.
1127       ELSEIF ( bc_ns == 'radiation/dirichlet' )  THEN
1128          inflow_s  = .TRUE.
1129       ELSEIF ( bc_ns == 'nested' )  THEN
[1762]1130          nest_bound_s = .TRUE.
[1]1131       ENDIF
1132    ENDIF
1133
[1964]1134    IF ( pnorth == MPI_PROC_NULL )  THEN
1135       IF ( bc_ns == 'dirichlet/radiation' )  THEN
1136          inflow_n  = .TRUE.
1137       ELSEIF ( bc_ns == 'radiation/dirichlet' )  THEN
1138          outflow_n = .TRUE.
1139       ELSEIF ( bc_ns == 'nested' )  THEN
[1762]1140          nest_bound_n = .TRUE.
[1]1141       ENDIF
1142    ENDIF
[1933]1143       
[151]1144!
1145!-- Broadcast the id of the inflow PE
1146    IF ( inflow_l )  THEN
[163]1147       id_inflow_l = myidx
[151]1148    ELSE
1149       id_inflow_l = 0
1150    ENDIF
[622]1151    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
[151]1152    CALL MPI_ALLREDUCE( id_inflow_l, id_inflow, 1, MPI_INTEGER, MPI_SUM, &
1153                        comm1dx, ierr )
1154
[163]1155!
1156!-- Broadcast the id of the recycling plane
1157!-- WARNING: needs to be adjusted in case of inflows other than from left side!
[1139]1158    IF ( NINT( recycling_width / dx ) >= nxl  .AND. &
1159         NINT( recycling_width / dx ) <= nxr )  THEN
[163]1160       id_recycling_l = myidx
1161    ELSE
1162       id_recycling_l = 0
1163    ENDIF
[622]1164    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
[163]1165    CALL MPI_ALLREDUCE( id_recycling_l, id_recycling, 1, MPI_INTEGER, MPI_SUM, &
1166                        comm1dx, ierr )
1167
[1402]1168    CALL location_message( 'finished', .TRUE. )
[1384]1169
[1804]1170#else
[1159]1171    IF ( bc_lr == 'dirichlet/radiation' )  THEN
[1]1172       inflow_l  = .TRUE.
1173       outflow_r = .TRUE.
[1159]1174    ELSEIF ( bc_lr == 'radiation/dirichlet' )  THEN
[1]1175       outflow_l = .TRUE.
1176       inflow_r  = .TRUE.
1177    ENDIF
1178
[1159]1179    IF ( bc_ns == 'dirichlet/radiation' )  THEN
[1]1180       inflow_n  = .TRUE.
1181       outflow_s = .TRUE.
[1159]1182    ELSEIF ( bc_ns == 'radiation/dirichlet' )  THEN
[1]1183       outflow_n = .TRUE.
1184       inflow_s  = .TRUE.
1185    ENDIF
1186#endif
[807]1187
[106]1188!
[978]1189!-- At the inflow or outflow, u or v, respectively, have to be calculated for
1190!-- one more grid point.
[1762]1191    IF ( inflow_l .OR. outflow_l .OR. nest_bound_l )  THEN
[106]1192       nxlu = nxl + 1
1193    ELSE
1194       nxlu = nxl
1195    ENDIF
[1762]1196    IF ( inflow_s .OR. outflow_s .OR. nest_bound_s )  THEN
[106]1197       nysv = nys + 1
1198    ELSE
1199       nysv = nys
1200    ENDIF
[1]1201
[114]1202!
1203!-- Allocate wall flag arrays used in the multigrid solver
[1575]1204    IF ( psolver(1:9) == 'multigrid' )  THEN
[114]1205
1206       DO  i = maximum_grid_level, 1, -1
1207
1208           SELECT CASE ( i )
1209
1210              CASE ( 1 )
1211                 ALLOCATE( wall_flags_1(nzb:nzt_mg(i)+1,         &
1212                                        nys_mg(i)-1:nyn_mg(i)+1, &
1213                                        nxl_mg(i)-1:nxr_mg(i)+1) )
1214
1215              CASE ( 2 )
1216                 ALLOCATE( wall_flags_2(nzb:nzt_mg(i)+1,         &
1217                                        nys_mg(i)-1:nyn_mg(i)+1, &
1218                                        nxl_mg(i)-1:nxr_mg(i)+1) )
1219
1220              CASE ( 3 )
1221                 ALLOCATE( wall_flags_3(nzb:nzt_mg(i)+1,         &
1222                                        nys_mg(i)-1:nyn_mg(i)+1, &
1223                                        nxl_mg(i)-1:nxr_mg(i)+1) )
1224
1225              CASE ( 4 )
1226                 ALLOCATE( wall_flags_4(nzb:nzt_mg(i)+1,         &
1227                                        nys_mg(i)-1:nyn_mg(i)+1, &
1228                                        nxl_mg(i)-1:nxr_mg(i)+1) )
1229
1230              CASE ( 5 )
1231                 ALLOCATE( wall_flags_5(nzb:nzt_mg(i)+1,         &
1232                                        nys_mg(i)-1:nyn_mg(i)+1, &
1233                                        nxl_mg(i)-1:nxr_mg(i)+1) )
1234
1235              CASE ( 6 )
1236                 ALLOCATE( wall_flags_6(nzb:nzt_mg(i)+1,         &
1237                                        nys_mg(i)-1:nyn_mg(i)+1, &
1238                                        nxl_mg(i)-1:nxr_mg(i)+1) )
1239
1240              CASE ( 7 )
1241                 ALLOCATE( wall_flags_7(nzb:nzt_mg(i)+1,         &
1242                                        nys_mg(i)-1:nyn_mg(i)+1, &
1243                                        nxl_mg(i)-1:nxr_mg(i)+1) )
1244
1245              CASE ( 8 )
1246                 ALLOCATE( wall_flags_8(nzb:nzt_mg(i)+1,         &
1247                                        nys_mg(i)-1:nyn_mg(i)+1, &
1248                                        nxl_mg(i)-1:nxr_mg(i)+1) )
1249
1250              CASE ( 9 )
1251                 ALLOCATE( wall_flags_9(nzb:nzt_mg(i)+1,         &
1252                                        nys_mg(i)-1:nyn_mg(i)+1, &
1253                                        nxl_mg(i)-1:nxr_mg(i)+1) )
1254
1255              CASE ( 10 )
1256                 ALLOCATE( wall_flags_10(nzb:nzt_mg(i)+1,        &
1257                                        nys_mg(i)-1:nyn_mg(i)+1, &
1258                                        nxl_mg(i)-1:nxr_mg(i)+1) )
1259
1260              CASE DEFAULT
[254]1261                 message_string = 'more than 10 multigrid levels'
1262                 CALL message( 'init_pegrid', 'PA0238', 1, 2, 0, 6, 0 )
[114]1263
1264          END SELECT
1265
1266       ENDDO
1267
1268    ENDIF
1269
[759]1270!
1271!-- Calculate the number of groups into which parallel I/O is split.
1272!-- The default for files which are opened by all PEs (or where each
1273!-- PE opens his own independent file) is, that all PEs are doing input/output
1274!-- in parallel at the same time. This might cause performance or even more
1275!-- severe problems depending on the configuration of the underlying file
1276!-- system.
1277!-- First, set the default:
1278    IF ( maximum_parallel_io_streams == -1  .OR. &
1279         maximum_parallel_io_streams > numprocs )  THEN
1280       maximum_parallel_io_streams = numprocs
1281    ENDIF
1282
1283!
1284!-- Now calculate the number of io_blocks and the io_group to which the
1285!-- respective PE belongs. I/O of the groups is done in serial, but in parallel
1286!-- for all PEs belonging to the same group. A preliminary setting with myid
1287!-- based on MPI_COMM_WORLD has been done in parin.
1288    io_blocks = numprocs / maximum_parallel_io_streams
1289    io_group  = MOD( myid+1, io_blocks )
1290   
1291
[1]1292 END SUBROUTINE init_pegrid
Note: See TracBrowser for help on using the repository browser.