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

Last change on this file since 2178 was 2178, checked in by hellstea, 7 years ago

Nesting bugfixes

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