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

Last change on this file since 480 was 438, checked in by raasch, 14 years ago

2d-decomposition is default on Cray-XT machines

  • Property svn:keywords set to Id
File size: 34.0 KB
Line 
1 SUBROUTINE init_pegrid
2
3!------------------------------------------------------------------------------!
4! Current revisions:
5! -----------------
6!
7! 2d-decomposition is default for Cray-XT machines
8! ATTENTION: nnz_x undefined problem still has to be solved!!!!!!!!
9! TEST OUTPUT (TO BE REMOVED) logging mpi2 ierr values
10!
11! Former revisions:
12! -----------------
13! $Id: init_pegrid.f90 438 2010-02-01 04:32:43Z raasch $
14!
15! 274 2009-03-26 15:11:21Z heinze
16! Output of messages replaced by message handling routine.
17!
18! 206 2008-10-13 14:59:11Z raasch
19! Implementation of a MPI-1 coupling: added __parallel within the __mpi2 part
20! 2d-decomposition is default on SGI-ICE systems
21!
22! 197 2008-09-16 15:29:03Z raasch
23! multigrid levels are limited by subdomains if mg_switch_to_pe0_level = -1,
24! nz is used instead nnz for calculating mg-levels
25! Collect on PE0 horizontal index bounds from all other PEs,
26! broadcast the id of the inflow PE (using the respective communicator)
27!
28! 114 2007-10-10 00:03:15Z raasch
29! Allocation of wall flag arrays for multigrid solver
30!
31! 108 2007-08-24 15:10:38Z letzel
32! Intercommunicator (comm_inter) and derived data type (type_xy) for
33! coupled model runs created, assign coupling_mode_remote,
34! indices nxlu and nysv are calculated (needed for non-cyclic boundary
35! conditions)
36!
37! 82 2007-04-16 15:40:52Z raasch
38! Cpp-directive lcmuk changed to intel_openmp_bug, setting of host on lcmuk by
39! cpp-directive removed
40!
41! 75 2007-03-22 09:54:05Z raasch
42! uxrp, vynp eliminated,
43! dirichlet/neumann changed to dirichlet/radiation, etc.,
44! poisfft_init is only called if fft-solver is switched on
45!
46! RCS Log replace by Id keyword, revision history cleaned up
47!
48! Revision 1.28  2006/04/26 13:23:32  raasch
49! lcmuk does not understand the !$ comment so a cpp-directive is required
50!
51! Revision 1.1  1997/07/24 11:15:09  raasch
52! Initial revision
53!
54!
55! Description:
56! ------------
57! Determination of the virtual processor topology (if not prescribed by the
58! user)and computation of the grid point number and array bounds of the local
59! domains.
60!------------------------------------------------------------------------------!
61
62    USE control_parameters
63    USE fft_xy
64    USE grid_variables
65    USE indices
66    USE pegrid
67    USE poisfft_mod
68    USE poisfft_hybrid_mod
69    USE statistics
70    USE transpose_indices
71
72
73    IMPLICIT NONE
74
75    INTEGER ::  gathered_size, i, id_inflow_l, id_recycling_l, ind(5), j, k, &
76                maximum_grid_level_l, mg_switch_to_pe0_level_l, mg_levels_x, &
77                mg_levels_y, mg_levels_z, nnx_y, nnx_z, nny_x, nny_z, nnz_x, &
78                nnz_y, numproc_sqr, nx_total, nxl_l, nxr_l, nyn_l, nys_l,    &
79                nzb_l, nzt_l, omp_get_num_threads, subdomain_size
80
81    INTEGER, DIMENSION(:), ALLOCATABLE ::  ind_all, nxlf, nxrf, nynf, nysf
82
83    LOGICAL ::  found
84
85!
86!-- Get the number of OpenMP threads
87    !$OMP PARALLEL
88#if defined( __intel_openmp_bug )
89    threads_per_task = omp_get_num_threads()
90#else
91!$  threads_per_task = omp_get_num_threads()
92#endif
93    !$OMP END PARALLEL
94
95
96#if defined( __parallel )
97!
98!-- Determine the processor topology or check it, if prescribed by the user
99    IF ( npex == -1  .AND.  npey == -1 )  THEN
100
101!
102!--    Automatic determination of the topology
103!--    The default on SMP- and cluster-hosts is a 1d-decomposition along x
104       IF ( host(1:3) == 'ibm'  .OR.  host(1:3) == 'nec'      .OR. &
105            ( host(1:2) == 'lc'  .AND.  host(3:5) /= 'sgi'  .AND.  &
106              host(3:4) /= 'xt' )                             .OR. &
107             host(1:3) == 'dec' )  THEN
108
109          pdims(1) = numprocs
110          pdims(2) = 1
111
112       ELSE
113
114          numproc_sqr = SQRT( REAL( numprocs ) )
115          pdims(1)    = MAX( numproc_sqr , 1 )
116          DO  WHILE ( MOD( numprocs , pdims(1) ) /= 0 )
117             pdims(1) = pdims(1) - 1
118          ENDDO
119          pdims(2) = numprocs / pdims(1)
120
121       ENDIF
122
123    ELSEIF ( npex /= -1  .AND.  npey /= -1 )  THEN
124
125!
126!--    Prescribed by user. Number of processors on the prescribed topology
127!--    must be equal to the number of PEs available to the job
128       IF ( ( npex * npey ) /= numprocs )  THEN
129          WRITE( message_string, * ) 'number of PEs of the prescribed ',      & 
130                 'topology (', npex*npey,') does not match & the number of ', & 
131                 'PEs available to the job (', numprocs, ')'
132          CALL message( 'init_pegrid', 'PA0221', 1, 2, 0, 6, 0 )
133       ENDIF
134       pdims(1) = npex
135       pdims(2) = npey
136
137    ELSE
138!
139!--    If the processor topology is prescribed by the user, the number of
140!--    PEs must be given in both directions
141       message_string = 'if the processor topology is prescribed by the, ' //  &
142                   ' user& both values of "npex" and "npey" must be given ' // &
143                   'in the &NAMELIST-parameter file'
144       CALL message( 'init_pegrid', 'PA0222', 1, 2, 0, 6, 0 )
145
146    ENDIF
147
148!
149!-- The hybrid solver can only be used in case of a 1d-decomposition along x
150    IF ( pdims(2) /= 1  .AND.  psolver == 'poisfft_hybrid' )  THEN
151       message_string = 'psolver = "poisfft_hybrid" can only be' // &
152                        '& used in case of a 1d-decomposition along x'
153       CALL message( 'init_pegrid', 'PA0223', 1, 2, 0, 6, 0 )
154    ENDIF
155
156!
157!-- If necessary, set horizontal boundary conditions to non-cyclic
158    IF ( bc_lr /= 'cyclic' )  cyclic(1) = .FALSE.
159    IF ( bc_ns /= 'cyclic' )  cyclic(2) = .FALSE.
160
161!
162!-- Create the virtual processor grid
163    CALL MPI_CART_CREATE( comm_palm, ndim, pdims, cyclic, reorder, &
164                          comm2d, ierr )
165    CALL MPI_COMM_RANK( comm2d, myid, ierr )
166    WRITE (myid_char,'(''_'',I4.4)')  myid
167
168    CALL MPI_CART_COORDS( comm2d, myid, ndim, pcoord, ierr )
169    CALL MPI_CART_SHIFT( comm2d, 0, 1, pleft, pright, ierr )
170    CALL MPI_CART_SHIFT( comm2d, 1, 1, psouth, pnorth, ierr )
171
172!
173!-- Determine sub-topologies for transpositions
174!-- Transposition from z to x:
175    remain_dims(1) = .TRUE.
176    remain_dims(2) = .FALSE.
177    CALL MPI_CART_SUB( comm2d, remain_dims, comm1dx, ierr )
178    CALL MPI_COMM_RANK( comm1dx, myidx, ierr )
179!
180!-- Transposition from x to y
181    remain_dims(1) = .FALSE.
182    remain_dims(2) = .TRUE.
183    CALL MPI_CART_SUB( comm2d, remain_dims, comm1dy, ierr )
184    CALL MPI_COMM_RANK( comm1dy, myidy, ierr )
185
186
187!
188!-- Find a grid (used for array d) which will match the transposition demands
189    IF ( grid_matching == 'strict' )  THEN
190
191       nxa = nx;  nya = ny;  nza = nz
192
193    ELSE
194
195       found = .FALSE.
196   xn: DO  nxa = nx, 2*nx
197!
198!--       Meet conditions for nx
199          IF ( MOD( nxa+1, pdims(1) ) /= 0 .OR. &
200               MOD( nxa+1, pdims(2) ) /= 0 )  CYCLE xn
201
202      yn: DO  nya = ny, 2*ny
203!
204!--          Meet conditions for ny
205             IF ( MOD( nya+1, pdims(2) ) /= 0 .OR. &
206                  MOD( nya+1, pdims(1) ) /= 0 )  CYCLE yn
207
208
209         zn: DO  nza = nz, 2*nz
210!
211!--             Meet conditions for nz
212                IF ( ( MOD( nza, pdims(1) ) /= 0  .AND.  pdims(1) /= 1  .AND. &
213                       pdims(2) /= 1 )  .OR.                                  &
214                     ( MOD( nza, pdims(2) ) /= 0  .AND.  dt_dosp /= 9999999.9 &
215                     ) )  THEN
216                   CYCLE zn
217                ELSE
218                   found = .TRUE.
219                   EXIT xn
220                ENDIF
221
222             ENDDO zn
223
224          ENDDO yn
225
226       ENDDO xn
227
228       IF ( .NOT. found )  THEN
229          message_string = 'no matching grid for transpositions found'
230          CALL message( 'init_pegrid', 'PA0224', 1, 2, 0, 6, 0 )
231       ENDIF
232
233    ENDIF
234
235!
236!-- Calculate array bounds in x-direction for every PE.
237!-- The last PE along x may get less grid points than the others
238    ALLOCATE( nxlf(0:pdims(1)-1), nxrf(0:pdims(1)-1), nynf(0:pdims(2)-1), &
239              nysf(0:pdims(2)-1), nnx_pe(0:pdims(1)-1), nny_pe(0:pdims(2)-1) )
240
241    IF ( MOD( nxa+1 , pdims(1) ) /= 0 )  THEN
242       WRITE( message_string, * ) 'x-direction: gridpoint number (',nx+1,') ',&
243                               'is not an& integral divisor of the number ',  &
244                               'processors (', pdims(1),')'
245       CALL message( 'init_pegrid', 'PA0225', 1, 2, 0, 6, 0 )
246    ELSE
247       nnx  = ( nxa + 1 ) / pdims(1)
248       IF ( nnx*pdims(1) - ( nx + 1) > nnx )  THEN
249          WRITE( message_string, * ) 'x-direction: nx does not match the',    & 
250                       'requirements given by the number of PEs &used',       &
251                       '& please use nx = ', nx - ( pdims(1) - ( nnx*pdims(1) &
252                                   - ( nx + 1 ) ) ), ' instead of nx =', nx
253          CALL message( 'init_pegrid', 'PA0226', 1, 2, 0, 6, 0 )
254       ENDIF
255    ENDIF   
256
257!
258!-- Left and right array bounds, number of gridpoints
259    DO  i = 0, pdims(1)-1
260       nxlf(i)   = i * nnx
261       nxrf(i)   = ( i + 1 ) * nnx - 1
262       nnx_pe(i) = MIN( nx, nxrf(i) ) - nxlf(i) + 1
263    ENDDO
264
265!
266!-- Calculate array bounds in y-direction for every PE.
267    IF ( MOD( nya+1 , pdims(2) ) /= 0 )  THEN
268       WRITE( message_string, * ) 'y-direction: gridpoint number (',ny+1,') ', &
269                           'is not an& integral divisor of the number of',     &
270                           'processors (', pdims(2),')'
271       CALL message( 'init_pegrid', 'PA0227', 1, 2, 0, 6, 0 )
272    ELSE
273       nny  = ( nya + 1 ) / pdims(2)
274       IF ( nny*pdims(2) - ( ny + 1) > nny )  THEN
275          WRITE( message_string, * ) 'y-direction: ny does not match the',    &
276                       'requirements given by the number of PEs &used ',      &
277                       '& please use ny = ', ny - ( pdims(2) - ( nnx*pdims(2) &
278                                     - ( ny + 1 ) ) ), ' instead of ny =', ny
279          CALL message( 'init_pegrid', 'PA0228', 1, 2, 0, 6, 0 )
280       ENDIF
281    ENDIF   
282
283!
284!-- South and north array bounds
285    DO  j = 0, pdims(2)-1
286       nysf(j)   = j * nny
287       nynf(j)   = ( j + 1 ) * nny - 1
288       nny_pe(j) = MIN( ny, nynf(j) ) - nysf(j) + 1
289    ENDDO
290
291!
292!-- Local array bounds of the respective PEs
293    nxl  = nxlf(pcoord(1))
294    nxra = nxrf(pcoord(1))
295    nxr  = MIN( nx, nxra )
296    nys  = nysf(pcoord(2))
297    nyna = nynf(pcoord(2))
298    nyn  = MIN( ny, nyna )
299    nzb  = 0
300    nzta = nza
301    nzt  = MIN( nz, nzta )
302    nnz  = nza
303
304!
305!-- Calculate array bounds and gridpoint numbers for the transposed arrays
306!-- (needed in the pressure solver)
307!-- For the transposed arrays, cyclic boundaries as well as top and bottom
308!-- boundaries are omitted, because they are obstructive to the transposition
309
310!
311!-- 1. transposition  z --> x
312!-- This transposition is not neccessary in case of a 1d-decomposition along x,
313!-- except that the uptream-spline method is switched on
314    IF ( pdims(2) /= 1  .OR.  momentum_advec == 'ups-scheme'  .OR. &
315         scalar_advec == 'ups-scheme' )  THEN
316
317       IF ( pdims(2) == 1  .AND. ( momentum_advec == 'ups-scheme'  .OR. &
318            scalar_advec == 'ups-scheme' ) )  THEN
319          message_string = '1d-decomposition along x ' // &
320                           'chosen but nz restrictions may occur' // &
321                           '& since ups-scheme is activated'
322          CALL message( 'init_pegrid', 'PA0229', 0, 1, 0, 6, 0 )
323       ENDIF
324       nys_x  = nys
325       nyn_xa = nyna
326       nyn_x  = nyn
327       nny_x  = nny
328       IF ( MOD( nza , pdims(1) ) /= 0 )  THEN
329          WRITE( message_string, * ) 'transposition z --> x:',                &
330                       '&nz=',nz,' is not an integral divisior of pdims(1)=', &
331                                                                   pdims(1)
332          CALL message( 'init_pegrid', 'PA0230', 1, 2, 0, 6, 0 )
333       ENDIF
334       nnz_x  = nza / pdims(1)
335       nzb_x  = 1 + myidx * nnz_x
336       nzt_xa = ( myidx + 1 ) * nnz_x
337       nzt_x  = MIN( nzt, nzt_xa )
338
339       sendrecvcount_zx = nnx * nny * nnz_x
340
341    ELSE
342!
343!---   Setting of dummy values because otherwise variables are undefined in
344!---   the next step  x --> y
345!---   WARNING: This case has still to be clarified!!!!!!!!!!!!
346       nnz_x  = 1
347       nzb_x  = 1
348       nzt_xa = 1
349       nzt_x  = 1
350       nny_x  = nny
351
352    ENDIF
353
354!
355!-- 2. transposition  x --> y
356    nnz_y  = nnz_x
357    nzb_y  = nzb_x
358    nzt_ya = nzt_xa
359    nzt_y  = nzt_x
360    IF ( MOD( nxa+1 , pdims(2) ) /= 0 )  THEN
361       WRITE( message_string, * ) 'transposition x --> y:',                &
362                         '&nx+1=',nx+1,' is not an integral divisor of ',&
363                         'pdims(2)=',pdims(2)
364       CALL message( 'init_pegrid', 'PA0231', 1, 2, 0, 6, 0 )
365    ENDIF
366    nnx_y = (nxa+1) / pdims(2)
367    nxl_y = myidy * nnx_y
368    nxr_ya = ( myidy + 1 ) * nnx_y - 1
369    nxr_y  = MIN( nx, nxr_ya )
370
371    sendrecvcount_xy = nnx_y * nny_x * nnz_y
372
373!
374!-- 3. transposition  y --> z  (ELSE:  x --> y  in case of 1D-decomposition
375!-- along x)
376    IF ( pdims(2) /= 1  .OR.  momentum_advec == 'ups-scheme'  .OR. &
377         scalar_advec == 'ups-scheme' )  THEN
378!
379!--    y --> z
380!--    This transposition is not neccessary in case of a 1d-decomposition
381!--    along x, except that the uptream-spline method is switched on
382       nnx_z  = nnx_y
383       nxl_z  = nxl_y
384       nxr_za = nxr_ya
385       nxr_z  = nxr_y
386       IF ( MOD( nya+1 , pdims(1) ) /= 0 )  THEN
387          WRITE( message_string, * ) 'transposition y --> z:',            &
388                            '& ny+1=',ny+1,' is not an integral divisor of ',&
389                            'pdims(1)=',pdims(1)
390          CALL message( 'init_pegrid', 'PA0232', 1, 2, 0, 6, 0 )
391       ENDIF
392       nny_z  = (nya+1) / pdims(1)
393       nys_z  = myidx * nny_z
394       nyn_za = ( myidx + 1 ) * nny_z - 1
395       nyn_z  = MIN( ny, nyn_za )
396
397       sendrecvcount_yz = nnx_y * nny_z * nnz_y
398
399    ELSE
400!
401!--    x --> y. This condition must be fulfilled for a 1D-decomposition along x
402       IF ( MOD( nya+1 , pdims(1) ) /= 0 )  THEN
403          WRITE( message_string, * ) 'transposition x --> y:',               &
404                            '& ny+1=',ny+1,' is not an integral divisor of ',&
405                            'pdims(1)=',pdims(1)
406          CALL message( 'init_pegrid', 'PA0233', 1, 2, 0, 6, 0 )
407       ENDIF
408
409    ENDIF
410
411!
412!-- Indices for direct transpositions z --> y (used for calculating spectra)
413    IF ( dt_dosp /= 9999999.9 )  THEN
414       IF ( MOD( nza, pdims(2) ) /= 0 )  THEN
415          WRITE( message_string, * ) 'direct transposition z --> y (needed ', &
416                    'for spectra):& nz=',nz,' is not an integral divisor of ',&
417                    'pdims(2)=',pdims(2)
418          CALL message( 'init_pegrid', 'PA0234', 1, 2, 0, 6, 0 )
419       ELSE
420          nxl_yd  = nxl
421          nxr_yda = nxra
422          nxr_yd  = nxr
423          nzb_yd  = 1 + myidy * ( nza / pdims(2) )
424          nzt_yda = ( myidy + 1 ) * ( nza / pdims(2) )
425          nzt_yd  = MIN( nzt, nzt_yda )
426
427          sendrecvcount_zyd = nnx * nny * ( nza / pdims(2) )
428       ENDIF
429    ENDIF
430
431!
432!-- Indices for direct transpositions y --> x (they are only possible in case
433!-- of a 1d-decomposition along x)
434    IF ( pdims(2) == 1 )  THEN
435       nny_x  = nny / pdims(1)
436       nys_x  = myid * nny_x
437       nyn_xa = ( myid + 1 ) * nny_x - 1
438       nyn_x  = MIN( ny, nyn_xa )
439       nzb_x  = 1
440       nzt_xa = nza
441       nzt_x  = nz
442       sendrecvcount_xy = nnx * nny_x * nza
443    ENDIF
444
445!
446!-- Indices for direct transpositions x --> y (they are only possible in case
447!-- of a 1d-decomposition along y)
448    IF ( pdims(1) == 1 )  THEN
449       nnx_y  = nnx / pdims(2)
450       nxl_y  = myid * nnx_y
451       nxr_ya = ( myid + 1 ) * nnx_y - 1
452       nxr_y  = MIN( nx, nxr_ya )
453       nzb_y  = 1
454       nzt_ya = nza
455       nzt_y  = nz
456       sendrecvcount_xy = nnx_y * nny * nza
457    ENDIF
458
459!
460!-- Arrays for storing the array bounds are needed any more
461    DEALLOCATE( nxlf , nxrf , nynf , nysf )
462
463!
464!-- Collect index bounds from other PEs (to be written to restart file later)
465    ALLOCATE( hor_index_bounds(4,0:numprocs-1) )
466
467    IF ( myid == 0 )  THEN
468
469       hor_index_bounds(1,0) = nxl
470       hor_index_bounds(2,0) = nxr
471       hor_index_bounds(3,0) = nys
472       hor_index_bounds(4,0) = nyn
473
474!
475!--    Receive data from all other PEs
476       DO  i = 1, numprocs-1
477          CALL MPI_RECV( ibuf, 4, MPI_INTEGER, i, MPI_ANY_TAG, comm2d, status, &
478                         ierr )
479          hor_index_bounds(:,i) = ibuf(1:4)
480       ENDDO
481
482    ELSE
483!
484!--    Send index bounds to PE0
485       ibuf(1) = nxl
486       ibuf(2) = nxr
487       ibuf(3) = nys
488       ibuf(4) = nyn
489       CALL MPI_SEND( ibuf, 4, MPI_INTEGER, 0, myid, comm2d, ierr )
490
491    ENDIF
492
493#if defined( __print )
494!
495!-- Control output
496    IF ( myid == 0 )  THEN
497       PRINT*, '*** processor topology ***'
498       PRINT*, ' '
499       PRINT*, 'myid   pcoord    left right  south north  idx idy   nxl: nxr',&
500               &'   nys: nyn'
501       PRINT*, '------------------------------------------------------------',&
502               &'-----------'
503       WRITE (*,1000)  0, pcoord(1), pcoord(2), pleft, pright, psouth, pnorth, &
504                       myidx, myidy, nxl, nxr, nys, nyn
5051000   FORMAT (I4,2X,'(',I3,',',I3,')',3X,I4,2X,I4,3X,I4,2X,I4,2X,I3,1X,I3, &
506               2(2X,I4,':',I4))
507
508!
509!--    Receive data from the other PEs
510       DO  i = 1,numprocs-1
511          CALL MPI_RECV( ibuf, 12, MPI_INTEGER, i, MPI_ANY_TAG, comm2d, status, &
512                         ierr )
513          WRITE (*,1000)  i, ( ibuf(j) , j = 1,12 )
514       ENDDO
515    ELSE
516
517!
518!--    Send data to PE0
519       ibuf(1) = pcoord(1); ibuf(2) = pcoord(2); ibuf(3) = pleft
520       ibuf(4) = pright; ibuf(5) = psouth; ibuf(6) = pnorth; ibuf(7) = myidx
521       ibuf(8) = myidy; ibuf(9) = nxl; ibuf(10) = nxr; ibuf(11) = nys
522       ibuf(12) = nyn
523       CALL MPI_SEND( ibuf, 12, MPI_INTEGER, 0, myid, comm2d, ierr )       
524    ENDIF
525#endif
526
527#if defined( __parallel )
528#if defined( __mpi2 )
529!
530!-- In case of coupled runs, get the port name on PE0 of the atmosphere model
531!-- and pass it to PE0 of the ocean model
532    IF ( myid == 0 )  THEN
533
534       IF ( coupling_mode == 'atmosphere_to_ocean' )  THEN
535
536          CALL MPI_OPEN_PORT( MPI_INFO_NULL, port_name, ierr )
537!
538!--       TEST OUTPUT (TO BE REMOVED)
539          WRITE(9,*)  TRIM( coupling_mode ),  &
540               ', ierr after MPI_OPEN_PORT: ', ierr
541          CALL LOCAL_FLUSH( 9 )
542
543          CALL MPI_PUBLISH_NAME( 'palm_coupler', MPI_INFO_NULL, port_name, &
544                                 ierr )
545!
546!--       TEST OUTPUT (TO BE REMOVED)
547          WRITE(9,*)  TRIM( coupling_mode ),  &
548               ', ierr after MPI_PUBLISH_NAME: ', ierr
549          CALL LOCAL_FLUSH( 9 )
550
551!
552!--       Write a flag file for the ocean model and the other atmosphere
553!--       processes.
554!--       There seems to be a bug in MPICH2 which causes hanging processes
555!--       in case that execution of LOOKUP_NAME is continued too early
556!--       (i.e. before the port has been created)
557          OPEN( 90, FILE='COUPLING_PORT_OPENED', FORM='FORMATTED' )
558          WRITE ( 90, '(''TRUE'')' )
559          CLOSE ( 90 )
560
561       ELSEIF ( coupling_mode == 'ocean_to_atmosphere' )  THEN
562
563!
564!--       Continue only if the atmosphere model has created the port.
565!--       There seems to be a bug in MPICH2 which causes hanging processes
566!--       in case that execution of LOOKUP_NAME is continued too early
567!--       (i.e. before the port has been created)
568          INQUIRE( FILE='COUPLING_PORT_OPENED', EXIST=found )
569          DO WHILE ( .NOT. found )
570             INQUIRE( FILE='COUPLING_PORT_OPENED', EXIST=found )
571          ENDDO
572
573          CALL MPI_LOOKUP_NAME( 'palm_coupler', MPI_INFO_NULL, port_name, ierr )
574!
575!--       TEST OUTPUT (TO BE REMOVED)
576          WRITE(9,*)  TRIM( coupling_mode ),  &
577               ', ierr after MPI_LOOKUP_NAME: ', ierr
578          CALL LOCAL_FLUSH( 9 )
579
580
581       ENDIF
582
583    ENDIF
584
585!
586!-- In case of coupled runs, establish the connection between the atmosphere
587!-- and the ocean model and define the intercommunicator (comm_inter)
588    CALL MPI_BARRIER( comm2d, ierr )
589    IF ( coupling_mode == 'atmosphere_to_ocean' )  THEN
590
591       PRINT*, '... before COMM_ACCEPT'
592       CALL MPI_COMM_ACCEPT( port_name, MPI_INFO_NULL, 0, MPI_COMM_WORLD, &
593                             comm_inter, ierr )
594       PRINT*, '--- ierr = ', ierr
595       PRINT*, '--- comm_inter atmosphere = ', comm_inter
596
597       coupling_mode_remote = 'ocean_to_atmosphere'
598
599    ELSEIF ( coupling_mode == 'ocean_to_atmosphere' )  THEN
600
601       IF ( myid == 0 )  PRINT*, '*** read: ', port_name, '  ierr = ', ierr
602       PRINT*, '... before COMM_CONNECT'
603       CALL MPI_COMM_CONNECT( port_name, MPI_INFO_NULL, 0, MPI_COMM_WORLD, &
604                              comm_inter, ierr )
605       PRINT*, '--- ierr = ', ierr
606       PRINT*, '--- comm_inter ocean      = ', comm_inter
607
608       coupling_mode_remote = 'atmosphere_to_ocean'
609
610    ENDIF
611#endif
612
613!
614!-- In case of coupled runs, create a new MPI derived datatype for the
615!-- exchange of surface (xy) data .
616!-- Gridpoint number for the exchange of ghost points (xy-plane)
617    ngp_xy  = ( nxr - nxl + 3 ) * ( nyn - nys + 3 )
618
619!
620!-- Define a new MPI derived datatype for the exchange of ghost points in
621!-- y-direction for 2D-arrays (line)
622    CALL MPI_TYPE_VECTOR( ngp_xy, 1, nzt-nzb+2, MPI_REAL, type_xy, ierr )
623    CALL MPI_TYPE_COMMIT( type_xy, ierr )
624#endif
625
626#else
627
628!
629!-- Array bounds when running on a single PE (respectively a non-parallel
630!-- machine)
631    nxl  = 0
632    nxr  = nx
633    nxra = nx
634    nnx  = nxr - nxl + 1
635    nys  = 0
636    nyn  = ny
637    nyna = ny
638    nny  = nyn - nys + 1
639    nzb  = 0
640    nzt  = nz
641    nzta = nz
642    nnz  = nz
643
644    ALLOCATE( hor_index_bounds(4,0:0) )
645    hor_index_bounds(1,0) = nxl
646    hor_index_bounds(2,0) = nxr
647    hor_index_bounds(3,0) = nys
648    hor_index_bounds(4,0) = nyn
649
650!
651!-- Array bounds for the pressure solver (in the parallel code, these bounds
652!-- are the ones for the transposed arrays)
653    nys_x  = nys
654    nyn_x  = nyn
655    nyn_xa = nyn
656    nzb_x  = nzb + 1
657    nzt_x  = nzt
658    nzt_xa = nzt
659
660    nxl_y  = nxl
661    nxr_y  = nxr
662    nxr_ya = nxr
663    nzb_y  = nzb + 1
664    nzt_y  = nzt
665    nzt_ya = nzt
666
667    nxl_z  = nxl
668    nxr_z  = nxr
669    nxr_za = nxr
670    nys_z  = nys
671    nyn_z  = nyn
672    nyn_za = nyn
673
674#endif
675
676!
677!-- Calculate number of grid levels necessary for the multigrid poisson solver
678!-- as well as the gridpoint indices on each level
679    IF ( psolver == 'multigrid' )  THEN
680
681!
682!--    First calculate number of possible grid levels for the subdomains
683       mg_levels_x = 1
684       mg_levels_y = 1
685       mg_levels_z = 1
686
687       i = nnx
688       DO WHILE ( MOD( i, 2 ) == 0  .AND.  i /= 2 )
689          i = i / 2
690          mg_levels_x = mg_levels_x + 1
691       ENDDO
692
693       j = nny
694       DO WHILE ( MOD( j, 2 ) == 0  .AND.  j /= 2 )
695          j = j / 2
696          mg_levels_y = mg_levels_y + 1
697       ENDDO
698
699       k = nz    ! do not use nnz because it might be > nz due to transposition
700                 ! requirements
701       DO WHILE ( MOD( k, 2 ) == 0  .AND.  k /= 2 )
702          k = k / 2
703          mg_levels_z = mg_levels_z + 1
704       ENDDO
705
706       maximum_grid_level = MIN( mg_levels_x, mg_levels_y, mg_levels_z )
707
708!
709!--    Find out, if the total domain allows more levels. These additional
710!--    levels are processed on PE0 only.
711       IF ( numprocs > 1  .AND.  mg_switch_to_pe0_level /= -1 )  THEN
712          IF ( mg_levels_z > MIN( mg_levels_x, mg_levels_y ) )  THEN
713             mg_switch_to_pe0_level_l = maximum_grid_level
714
715             mg_levels_x = 1
716             mg_levels_y = 1
717
718             i = nx+1
719             DO WHILE ( MOD( i, 2 ) == 0  .AND.  i /= 2 )
720                i = i / 2
721                mg_levels_x = mg_levels_x + 1
722             ENDDO
723
724             j = ny+1
725             DO WHILE ( MOD( j, 2 ) == 0  .AND.  j /= 2 )
726                j = j / 2
727                mg_levels_y = mg_levels_y + 1
728             ENDDO
729
730             maximum_grid_level_l = MIN( mg_levels_x, mg_levels_y, mg_levels_z )
731
732             IF ( maximum_grid_level_l > mg_switch_to_pe0_level_l )  THEN
733                mg_switch_to_pe0_level_l = maximum_grid_level_l - &
734                                           mg_switch_to_pe0_level_l + 1
735             ELSE
736                mg_switch_to_pe0_level_l = 0
737             ENDIF
738          ELSE
739             mg_switch_to_pe0_level_l = 0
740             maximum_grid_level_l = maximum_grid_level
741          ENDIF
742
743!
744!--       Use switch level calculated above only if it is not pre-defined
745!--       by user
746          IF ( mg_switch_to_pe0_level == 0 )  THEN
747
748             IF ( mg_switch_to_pe0_level_l /= 0 )  THEN
749                mg_switch_to_pe0_level = mg_switch_to_pe0_level_l
750                maximum_grid_level     = maximum_grid_level_l
751             ENDIF
752
753          ELSE
754!
755!--          Check pre-defined value and reset to default, if neccessary
756             IF ( mg_switch_to_pe0_level < mg_switch_to_pe0_level_l  .OR.  &
757                  mg_switch_to_pe0_level >= maximum_grid_level_l )  THEN
758                message_string = 'mg_switch_to_pe0_level ' // &
759                                 'out of range and reset to default (=0)'
760                CALL message( 'init_pegrid', 'PA0235', 0, 1, 0, 6, 0 )
761                mg_switch_to_pe0_level = 0
762             ELSE
763!
764!--             Use the largest number of possible levels anyway and recalculate
765!--             the switch level to this largest number of possible values
766                maximum_grid_level = maximum_grid_level_l
767
768             ENDIF
769          ENDIF
770
771       ENDIF
772
773       ALLOCATE( grid_level_count(maximum_grid_level),                   &
774                 nxl_mg(maximum_grid_level), nxr_mg(maximum_grid_level), &
775                 nyn_mg(maximum_grid_level), nys_mg(maximum_grid_level), &
776                 nzt_mg(maximum_grid_level) )
777
778       grid_level_count = 0
779       nxl_l = nxl; nxr_l = nxr; nys_l = nys; nyn_l = nyn; nzt_l = nzt
780
781       DO  i = maximum_grid_level, 1 , -1
782
783          IF ( i == mg_switch_to_pe0_level )  THEN
784#if defined( __parallel )
785!
786!--          Save the grid size of the subdomain at the switch level, because
787!--          it is needed in poismg.
788!--          Array bounds of the local subdomain grids are gathered on PE0
789             ind(1) = nxl_l; ind(2) = nxr_l
790             ind(3) = nys_l; ind(4) = nyn_l
791             ind(5) = nzt_l
792             ALLOCATE( ind_all(5*numprocs), mg_loc_ind(5,0:numprocs-1) )
793             CALL MPI_ALLGATHER( ind, 5, MPI_INTEGER, ind_all, 5, &
794                                 MPI_INTEGER, comm2d, ierr )
795             DO  j = 0, numprocs-1
796                DO  k = 1, 5
797                   mg_loc_ind(k,j) = ind_all(k+j*5)
798                ENDDO
799             ENDDO
800             DEALLOCATE( ind_all )
801!
802!--          Calculate the grid size of the total domain gathered on PE0
803             nxr_l = ( nxr_l-nxl_l+1 ) * pdims(1) - 1
804             nxl_l = 0
805             nyn_l = ( nyn_l-nys_l+1 ) * pdims(2) - 1
806             nys_l = 0
807!
808!--          The size of this gathered array must not be larger than the
809!--          array tend, which is used in the multigrid scheme as a temporary
810!--          array
811             subdomain_size = ( nxr - nxl + 3 )     * ( nyn - nys + 3 )     * &
812                              ( nzt - nzb + 2 )
813             gathered_size  = ( nxr_l - nxl_l + 3 ) * ( nyn_l - nys_l + 3 ) * &
814                              ( nzt_l - nzb + 2 )
815
816             IF ( gathered_size > subdomain_size )  THEN
817                message_string = 'not enough memory for storing ' // &
818                                 'gathered multigrid data on PE0'
819                CALL message( 'init_pegrid', 'PA0236', 1, 2, 0, 6, 0 )
820             ENDIF
821#else
822             message_string = 'multigrid gather/scatter impossible ' // &
823                          'in non parallel mode'
824             CALL message( 'init_pegrid', 'PA0237', 1, 2, 0, 6, 0 )
825#endif
826          ENDIF
827
828          nxl_mg(i) = nxl_l
829          nxr_mg(i) = nxr_l
830          nys_mg(i) = nys_l
831          nyn_mg(i) = nyn_l
832          nzt_mg(i) = nzt_l
833
834          nxl_l = nxl_l / 2 
835          nxr_l = nxr_l / 2
836          nys_l = nys_l / 2 
837          nyn_l = nyn_l / 2 
838          nzt_l = nzt_l / 2 
839       ENDDO
840
841    ELSE
842
843       maximum_grid_level = 1
844
845    ENDIF
846
847    grid_level = maximum_grid_level
848
849#if defined( __parallel )
850!
851!-- Gridpoint number for the exchange of ghost points (y-line for 2D-arrays)
852    ngp_y  = nyn - nys + 1
853
854!
855!-- Define a new MPI derived datatype for the exchange of ghost points in
856!-- y-direction for 2D-arrays (line)
857    CALL MPI_TYPE_VECTOR( nxr-nxl+3, 1, ngp_y+2, MPI_REAL, type_x, ierr )
858    CALL MPI_TYPE_COMMIT( type_x, ierr )
859    CALL MPI_TYPE_VECTOR( nxr-nxl+3, 1, ngp_y+2, MPI_INTEGER, type_x_int, ierr )
860    CALL MPI_TYPE_COMMIT( type_x_int, ierr )
861
862!
863!-- Calculate gridpoint numbers for the exchange of ghost points along x
864!-- (yz-plane for 3D-arrays) and define MPI derived data type(s) for the
865!-- exchange of ghost points in y-direction (xz-plane).
866!-- Do these calculations for the model grid and (if necessary) also
867!-- for the coarser grid levels used in the multigrid method
868    ALLOCATE ( ngp_yz(maximum_grid_level), type_xz(maximum_grid_level) )
869
870    nxl_l = nxl; nxr_l = nxr; nys_l = nys; nyn_l = nyn; nzb_l = nzb; nzt_l = nzt
871         
872    DO i = maximum_grid_level, 1 , -1
873       ngp_yz(i) = (nzt_l - nzb_l + 2) * (nyn_l - nys_l + 3)
874
875       CALL MPI_TYPE_VECTOR( nxr_l-nxl_l+3, nzt_l-nzb_l+2, ngp_yz(i), &
876                             MPI_REAL, type_xz(i), ierr )
877       CALL MPI_TYPE_COMMIT( type_xz(i), ierr )
878
879       nxl_l = nxl_l / 2 
880       nxr_l = nxr_l / 2
881       nys_l = nys_l / 2 
882       nyn_l = nyn_l / 2 
883       nzt_l = nzt_l / 2 
884    ENDDO
885#endif
886
887#if defined( __parallel )
888!
889!-- Setting of flags for inflow/outflow conditions in case of non-cyclic
890!-- horizontal boundary conditions.
891    IF ( pleft == MPI_PROC_NULL )  THEN
892       IF ( bc_lr == 'dirichlet/radiation' )  THEN
893          inflow_l  = .TRUE.
894       ELSEIF ( bc_lr == 'radiation/dirichlet' )  THEN
895          outflow_l = .TRUE.
896       ENDIF
897    ENDIF
898
899    IF ( pright == MPI_PROC_NULL )  THEN
900       IF ( bc_lr == 'dirichlet/radiation' )  THEN
901          outflow_r = .TRUE.
902       ELSEIF ( bc_lr == 'radiation/dirichlet' )  THEN
903          inflow_r  = .TRUE.
904       ENDIF
905    ENDIF
906
907    IF ( psouth == MPI_PROC_NULL )  THEN
908       IF ( bc_ns == 'dirichlet/radiation' )  THEN
909          outflow_s = .TRUE.
910       ELSEIF ( bc_ns == 'radiation/dirichlet' )  THEN
911          inflow_s  = .TRUE.
912       ENDIF
913    ENDIF
914
915    IF ( pnorth == MPI_PROC_NULL )  THEN
916       IF ( bc_ns == 'dirichlet/radiation' )  THEN
917          inflow_n  = .TRUE.
918       ELSEIF ( bc_ns == 'radiation/dirichlet' )  THEN
919          outflow_n = .TRUE.
920       ENDIF
921    ENDIF
922
923!
924!-- Broadcast the id of the inflow PE
925    IF ( inflow_l )  THEN
926       id_inflow_l = myidx
927    ELSE
928       id_inflow_l = 0
929    ENDIF
930    CALL MPI_ALLREDUCE( id_inflow_l, id_inflow, 1, MPI_INTEGER, MPI_SUM, &
931                        comm1dx, ierr )
932
933!
934!-- Broadcast the id of the recycling plane
935!-- WARNING: needs to be adjusted in case of inflows other than from left side!
936    IF ( ( recycling_width / dx ) >= nxl  .AND.  ( recycling_width / dx ) <= nxr ) &
937    THEN
938       id_recycling_l = myidx
939    ELSE
940       id_recycling_l = 0
941    ENDIF
942    CALL MPI_ALLREDUCE( id_recycling_l, id_recycling, 1, MPI_INTEGER, MPI_SUM, &
943                        comm1dx, ierr )
944
945#else
946    IF ( bc_lr == 'dirichlet/radiation' )  THEN
947       inflow_l  = .TRUE.
948       outflow_r = .TRUE.
949    ELSEIF ( bc_lr == 'radiation/dirichlet' )  THEN
950       outflow_l = .TRUE.
951       inflow_r  = .TRUE.
952    ENDIF
953
954    IF ( bc_ns == 'dirichlet/radiation' )  THEN
955       inflow_n  = .TRUE.
956       outflow_s = .TRUE.
957    ELSEIF ( bc_ns == 'radiation/dirichlet' )  THEN
958       outflow_n = .TRUE.
959       inflow_s  = .TRUE.
960    ENDIF
961#endif
962!
963!-- At the outflow, u or v, respectively, have to be calculated for one more
964!-- grid point.
965    IF ( outflow_l )  THEN
966       nxlu = nxl + 1
967    ELSE
968       nxlu = nxl
969    ENDIF
970    IF ( outflow_s )  THEN
971       nysv = nys + 1
972    ELSE
973       nysv = nys
974    ENDIF
975
976    IF ( psolver == 'poisfft_hybrid' )  THEN
977       CALL poisfft_hybrid_ini
978    ELSEIF ( psolver == 'poisfft' )  THEN
979       CALL poisfft_init
980    ENDIF
981
982!
983!-- Allocate wall flag arrays used in the multigrid solver
984    IF ( psolver == 'multigrid' )  THEN
985
986       DO  i = maximum_grid_level, 1, -1
987
988           SELECT CASE ( i )
989
990              CASE ( 1 )
991                 ALLOCATE( wall_flags_1(nzb:nzt_mg(i)+1,         &
992                                        nys_mg(i)-1:nyn_mg(i)+1, &
993                                        nxl_mg(i)-1:nxr_mg(i)+1) )
994
995              CASE ( 2 )
996                 ALLOCATE( wall_flags_2(nzb:nzt_mg(i)+1,         &
997                                        nys_mg(i)-1:nyn_mg(i)+1, &
998                                        nxl_mg(i)-1:nxr_mg(i)+1) )
999
1000              CASE ( 3 )
1001                 ALLOCATE( wall_flags_3(nzb:nzt_mg(i)+1,         &
1002                                        nys_mg(i)-1:nyn_mg(i)+1, &
1003                                        nxl_mg(i)-1:nxr_mg(i)+1) )
1004
1005              CASE ( 4 )
1006                 ALLOCATE( wall_flags_4(nzb:nzt_mg(i)+1,         &
1007                                        nys_mg(i)-1:nyn_mg(i)+1, &
1008                                        nxl_mg(i)-1:nxr_mg(i)+1) )
1009
1010              CASE ( 5 )
1011                 ALLOCATE( wall_flags_5(nzb:nzt_mg(i)+1,         &
1012                                        nys_mg(i)-1:nyn_mg(i)+1, &
1013                                        nxl_mg(i)-1:nxr_mg(i)+1) )
1014
1015              CASE ( 6 )
1016                 ALLOCATE( wall_flags_6(nzb:nzt_mg(i)+1,         &
1017                                        nys_mg(i)-1:nyn_mg(i)+1, &
1018                                        nxl_mg(i)-1:nxr_mg(i)+1) )
1019
1020              CASE ( 7 )
1021                 ALLOCATE( wall_flags_7(nzb:nzt_mg(i)+1,         &
1022                                        nys_mg(i)-1:nyn_mg(i)+1, &
1023                                        nxl_mg(i)-1:nxr_mg(i)+1) )
1024
1025              CASE ( 8 )
1026                 ALLOCATE( wall_flags_8(nzb:nzt_mg(i)+1,         &
1027                                        nys_mg(i)-1:nyn_mg(i)+1, &
1028                                        nxl_mg(i)-1:nxr_mg(i)+1) )
1029
1030              CASE ( 9 )
1031                 ALLOCATE( wall_flags_9(nzb:nzt_mg(i)+1,         &
1032                                        nys_mg(i)-1:nyn_mg(i)+1, &
1033                                        nxl_mg(i)-1:nxr_mg(i)+1) )
1034
1035              CASE ( 10 )
1036                 ALLOCATE( wall_flags_10(nzb:nzt_mg(i)+1,        &
1037                                        nys_mg(i)-1:nyn_mg(i)+1, &
1038                                        nxl_mg(i)-1:nxr_mg(i)+1) )
1039
1040              CASE DEFAULT
1041                 message_string = 'more than 10 multigrid levels'
1042                 CALL message( 'init_pegrid', 'PA0238', 1, 2, 0, 6, 0 )
1043
1044          END SELECT
1045
1046       ENDDO
1047
1048    ENDIF
1049
1050 END SUBROUTINE init_pegrid
Note: See TracBrowser for help on using the repository browser.