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

Last change on this file since 364 was 274, checked in by heinze, 16 years ago

Indentation of the message calls corrected

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