source: palm/tags/release-3.4a/SOURCE/init_pegrid.f90 @ 4003

Last change on this file since 4003 was 139, checked in by raasch, 16 years ago

New:
---

Plant canopy model of Watanabe (2004,BLM 112,307-341) added.
It can be switched on by the inipar parameter plant_canopy.
The inipar parameter canopy_mode can be used to prescribe a
plant canopy type. The default case is a homogeneous plant
canopy. Heterogeneous distributions of the leaf area
density and the canopy drag coefficient can be defined in the
new routine user_init_plant_canopy (user_interface).
The inipar parameters lad_surface, lad_vertical_gradient and
lad_vertical_gradient_level can be used in order to
prescribe the vertical profile of leaf area density. The
inipar parameter drag_coefficient determines the canopy
drag coefficient.
Finally, the inipar parameter pch_index determines the
index of the upper boundary of the plant canopy.

Allow new case bc_uv_t = 'dirichlet_0' for channel flow.

For unknown variables (CASE DEFAULT) call new subroutine user_data_output_dvrp

Pressure boundary conditions for vertical walls added to the multigrid solver.
They are applied using new wall flag arrays (wall_flags_..) which are defined
for each grid level. New argument gls added to routine user_init_grid
(user_interface).

Frequence of sorting particles can be controlled with new particles_par
parameter dt_sort_particles. Sorting is moved from the SGS timestep loop in
advec_particles after the end of this loop.

advec_particles, check_parameters, data_output_dvrp, header, init_3d_model, init_grid, init_particles, init_pegrid, modules, package_parin, parin, plant_canopy_model, read_var_list, read_3d_binary, user_interface, write_var_list, write_3d_binary

Changed:


Redefine initial nzb_local as the actual total size of topography (later the
extent of topography in nzb_local is reduced by 1dx at the E topography walls
and by 1dy at the N topography walls to form the basis for nzb_s_inner);
for consistency redefine 'single_building' case.

Vertical profiles now based on nzb_s_inner; they are divided by
ngp_2dh_s_inner (scalars, procucts of scalars) and ngp_2dh (staggered velocity
components and their products, procucts of scalars and velocity components),
respectively.

Allow two instead of one digit to specify isosurface and slicer variables.

Status of 3D-volume NetCDF data file only depends on switch netcdf_64bit_3d (check_open)

prognostic_equations include the respective wall_*flux in the parameter list of
calls of diffusion_s. Same as before, only the values of wall_heatflux(0:4)
can be assigned. At present, wall_humidityflux, wall_qflux, wall_salinityflux,
and wall_scalarflux are kept zero. diffusion_s uses the respective wall_*flux
instead of wall_heatflux. This update serves two purposes:

  • it avoids errors in calculations with humidity/scalar/salinity and prescribed

non-zero wall_heatflux,

  • it prepares PALM for a possible assignment of wall fluxes of

humidity/scalar/salinity in a future release.

buoyancy, check_open, data_output_dvrp, diffusion_s, diffusivities, flow_statistics, header, init_3d_model, init_dvrp, init_grid, modules, prognostic_equations

Errors:


Bugfix: summation of sums_l_l in diffusivities.

Several bugfixes in the ocean part: Initial density rho is calculated
(init_ocean). Error in initializing u_init and v_init removed
(check_parameters). Calculation of density flux now starts from
nzb+1 (production_e).

Bugfix: pleft/pright changed to pnorth/psouth in sendrecv of particle tail
numbers along y, small bugfixes in the SGS part (advec_particles)

Bugfix: model_string needed a default value (combine_plot_fields)

Bugfix: wavenumber calculation for even nx in routines maketri (poisfft)

Bugfix: assignment of fluxes at walls

Bugfix: absolute value of f must be used when calculating the Blackadar mixing length (init_1d_model)

advec_particles, check_parameters, combine_plot_fields, diffusion_s, diffusivities, init_ocean, init_1d_model, poisfft, production_e

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