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

Last change on this file since 590 was 482, checked in by raasch, 15 years ago

New:
---
compare_palm_logs is additionally compiled with mbuild -u (Makefile in trunk/UTIL)

make options (mopts) to be set by configuration file implemented (mrun, mbuild)

humidity=.T. is now usable for runs with topography. wall_humidityflux and
wall_scalarflux are the corresponding new parin arrays.
(check_parameters, init_3d_model, parin)

Large scale vertical motion (subsidence/ascent) can be applied to the
prognostic equation for the potential temperature. (check_parameters, header,
Makefile, modules, parin, prognostic_equations, read_var_list, subsidence,
write_var_list)

A simple method for installing and running palm (with limited features)
has been added. (Makefile, palm_simple_install, palm_simple_run)

Changed:


2d-decomposition is default for Cray-XT machines (init_pegrid)

var_ts is replaced by dots_max (modules,init_3d_model)

Every cloud droplet has now an own weighting factor and can be deleted due to
collisions. Condensation and collision of cloud droplets are adjusted
accordingly. (advec_particles)

Collision efficiency for large cloud droplets has changed according to table of
Rogers and Yau. (collision_efficiency)

Errors:


Bugfix for generating serial jobs (subjob)

Bugfix: index problem concerning gradient_level indices removed (header)

Dimension of array stat in cascade change to prevent type problems with
mpi2 libraries (poisfft_hybrid)

Loop was split to make runs reproducible when using ifort compiler.
(disturb_field)

Bugfix: exchange of ghost points for prho included (time_integration)

Bugfix: calculation of time-averaged surface heatfluxes (sum_up_3d_data)

Bugfix: calculation of precipitation_rate (calc_precipitation)

Bugfix: initial data assignments to some dvrp arrays changed due to error
messages from gfortran compiler (modules)

Bugfix: calculation of cloud droplet velocity (advec_particles)

Bugfix: transfer of particles at south/left edge (advec_particles)

Bugfix: calculation of collision_efficiency (collision_efficiency)

Bugfix: initialisation of var_mod (subsidence)

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