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

Last change on this file since 273 was 254, checked in by heinze, 16 years ago

Output of messages replaced by message handling routine.

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