SUBROUTINE init_pegrid !------------------------------------------------------------------------------! ! Current revisions: ! ----------------- ! Bugfix: bc_lr/ns_cyc/dirrad/raddir replaced by bc_lr/ns, because variables ! are not yet set here; grid_level set to 0 ! ! ATTENTION: nnz_x undefined problem still has to be solved!!!!!!!! ! TEST OUTPUT (TO BE REMOVED) logging mpi2 ierr values ! ! Former revisions: ! ----------------- ! $Id: init_pegrid.f90 722 2011-04-11 06:21:09Z raasch $ ! ! 709 2011-03-30 09:31:40Z raasch ! formatting adjustments ! ! 707 2011-03-29 11:39:40Z raasch ! bc_lr/ns replaced by bc_lr/ns_cyc/dirrad/raddir ! ! 667 2010-12-23 12:06:00Z suehring/gryschka ! Moved determination of target_id's from init_coupling ! Determination of parameters needed for coupling (coupling_topology, ngp_a, ! ngp_o) with different grid/processor-topology in ocean and atmosphere ! Adaption of ngp_xy, ngp_y to a dynamic number of ghost points. ! The maximum_grid_level changed from 1 to 0. 0 is the normal grid, 1 to ! maximum_grid_level the grids for multigrid, in which 0 and 1 are normal grids. ! This distinction is due to reasons of data exchange and performance for the ! normal grid and grids in poismg. ! The definition of MPI-Vectors adapted to a dynamic numer of ghost points. ! New MPI-Vectors for data exchange between left and right boundaries added. ! This is due to reasons of performance (10% faster). ! ! 646 2010-12-15 13:03:52Z raasch ! lctit is now using a 2d decomposition by default ! ! 622 2010-12-10 08:08:13Z raasch ! optional barriers included in order to speed up collective operations ! ! 438 2010-02-01 04:32:43Z raasch ! 2d-decomposition is default for Cray-XT machines ! ! 274 2009-03-26 15:11:21Z heinze ! Output of messages replaced by message handling routine. ! ! 206 2008-10-13 14:59:11Z raasch ! Implementation of a MPI-1 coupling: added __parallel within the __mpi2 part ! 2d-decomposition is default on SGI-ICE systems ! ! 197 2008-09-16 15:29:03Z raasch ! multigrid levels are limited by subdomains if mg_switch_to_pe0_level = -1, ! nz is used instead nnz for calculating mg-levels ! Collect on PE0 horizontal index bounds from all other PEs, ! broadcast the id of the inflow PE (using the respective communicator) ! ! 114 2007-10-10 00:03:15Z raasch ! Allocation of wall flag arrays for multigrid solver ! ! 108 2007-08-24 15:10:38Z letzel ! Intercommunicator (comm_inter) and derived data type (type_xy) for ! coupled model runs created, assign coupling_mode_remote, ! indices nxlu and nysv are calculated (needed for non-cyclic boundary ! conditions) ! ! 82 2007-04-16 15:40:52Z raasch ! Cpp-directive lcmuk changed to intel_openmp_bug, setting of host on lcmuk by ! cpp-directive removed ! ! 75 2007-03-22 09:54:05Z raasch ! uxrp, vynp eliminated, ! dirichlet/neumann changed to dirichlet/radiation, etc., ! poisfft_init is only called if fft-solver is switched on ! ! RCS Log replace by Id keyword, revision history cleaned up ! ! Revision 1.28 2006/04/26 13:23:32 raasch ! lcmuk does not understand the !$ comment so a cpp-directive is required ! ! Revision 1.1 1997/07/24 11:15:09 raasch ! Initial revision ! ! ! Description: ! ------------ ! Determination of the virtual processor topology (if not prescribed by the ! user)and computation of the grid point number and array bounds of the local ! domains. !------------------------------------------------------------------------------! USE control_parameters USE fft_xy USE grid_variables USE indices USE pegrid USE poisfft_mod USE poisfft_hybrid_mod USE statistics USE transpose_indices IMPLICIT NONE INTEGER :: gathered_size, i, id_inflow_l, id_recycling_l, ind(5), j, k, & maximum_grid_level_l, mg_switch_to_pe0_level_l, mg_levels_x, & mg_levels_y, mg_levels_z, nnx_y, nnx_z, nny_x, nny_z, nnz_x, & nnz_y, numproc_sqr, nx_total, nxl_l, nxr_l, nyn_l, nys_l, & nzb_l, nzt_l, omp_get_num_threads, subdomain_size INTEGER, DIMENSION(:), ALLOCATABLE :: ind_all, nxlf, nxrf, nynf, nysf INTEGER, DIMENSION(2) :: pdims_remote LOGICAL :: found ! !-- Get the number of OpenMP threads !$OMP PARALLEL #if defined( __intel_openmp_bug ) threads_per_task = omp_get_num_threads() #else !$ threads_per_task = omp_get_num_threads() #endif !$OMP END PARALLEL #if defined( __parallel ) ! !-- Determine the processor topology or check it, if prescribed by the user IF ( npex == -1 .AND. npey == -1 ) THEN ! !-- Automatic determination of the topology !-- The default on SMP- and cluster-hosts is a 1d-decomposition along x IF ( host(1:3) == 'ibm' .OR. host(1:3) == 'nec' .OR. & ( host(1:2) == 'lc' .AND. host(3:5) /= 'sgi' .AND. & host(3:4) /= 'xt' .AND. host(3:5) /= 'tit' ) .OR. & host(1:3) == 'dec' ) THEN pdims(1) = numprocs pdims(2) = 1 ELSE numproc_sqr = SQRT( REAL( numprocs ) ) pdims(1) = MAX( numproc_sqr , 1 ) DO WHILE ( MOD( numprocs , pdims(1) ) /= 0 ) pdims(1) = pdims(1) - 1 ENDDO pdims(2) = numprocs / pdims(1) ENDIF ELSEIF ( npex /= -1 .AND. npey /= -1 ) THEN ! !-- Prescribed by user. Number of processors on the prescribed topology !-- must be equal to the number of PEs available to the job IF ( ( npex * npey ) /= numprocs ) THEN WRITE( message_string, * ) 'number of PEs of the prescribed ', & 'topology (', npex*npey,') does not match & the number of ', & 'PEs available to the job (', numprocs, ')' CALL message( 'init_pegrid', 'PA0221', 1, 2, 0, 6, 0 ) ENDIF pdims(1) = npex pdims(2) = npey ELSE ! !-- If the processor topology is prescribed by the user, the number of !-- PEs must be given in both directions message_string = 'if the processor topology is prescribed by the, ' // & ' user& both values of "npex" and "npey" must be given ' // & 'in the &NAMELIST-parameter file' CALL message( 'init_pegrid', 'PA0222', 1, 2, 0, 6, 0 ) ENDIF ! !-- The hybrid solver can only be used in case of a 1d-decomposition along x IF ( pdims(2) /= 1 .AND. psolver == 'poisfft_hybrid' ) THEN message_string = 'psolver = "poisfft_hybrid" can only be' // & '& used in case of a 1d-decomposition along x' CALL message( 'init_pegrid', 'PA0223', 1, 2, 0, 6, 0 ) ENDIF ! !-- For communication speedup, set barriers in front of collective !-- communications by default on SGI-type systems IF ( host(3:5) == 'sgi' ) collective_wait = .TRUE. ! !-- If necessary, set horizontal boundary conditions to non-cyclic IF ( bc_lr /= 'cyclic' ) cyclic(1) = .FALSE. IF ( bc_ns /= 'cyclic' ) cyclic(2) = .FALSE. ! !-- Create the virtual processor grid CALL MPI_CART_CREATE( comm_palm, ndim, pdims, cyclic, reorder, & comm2d, ierr ) CALL MPI_COMM_RANK( comm2d, myid, ierr ) WRITE (myid_char,'(''_'',I4.4)') myid CALL MPI_CART_COORDS( comm2d, myid, ndim, pcoord, ierr ) CALL MPI_CART_SHIFT( comm2d, 0, 1, pleft, pright, ierr ) CALL MPI_CART_SHIFT( comm2d, 1, 1, psouth, pnorth, ierr ) ! !-- Determine sub-topologies for transpositions !-- Transposition from z to x: remain_dims(1) = .TRUE. remain_dims(2) = .FALSE. CALL MPI_CART_SUB( comm2d, remain_dims, comm1dx, ierr ) CALL MPI_COMM_RANK( comm1dx, myidx, ierr ) ! !-- Transposition from x to y remain_dims(1) = .FALSE. remain_dims(2) = .TRUE. CALL MPI_CART_SUB( comm2d, remain_dims, comm1dy, ierr ) CALL MPI_COMM_RANK( comm1dy, myidy, ierr ) ! !-- Find a grid (used for array d) which will match the transposition demands IF ( grid_matching == 'strict' ) THEN nxa = nx; nya = ny; nza = nz ELSE found = .FALSE. xn: DO nxa = nx, 2*nx ! !-- Meet conditions for nx IF ( MOD( nxa+1, pdims(1) ) /= 0 .OR. & MOD( nxa+1, pdims(2) ) /= 0 ) CYCLE xn yn: DO nya = ny, 2*ny ! !-- Meet conditions for ny IF ( MOD( nya+1, pdims(2) ) /= 0 .OR. & MOD( nya+1, pdims(1) ) /= 0 ) CYCLE yn zn: DO nza = nz, 2*nz ! !-- Meet conditions for nz IF ( ( MOD( nza, pdims(1) ) /= 0 .AND. pdims(1) /= 1 .AND. & pdims(2) /= 1 ) .OR. & ( MOD( nza, pdims(2) ) /= 0 .AND. dt_dosp /= 9999999.9 & ) ) THEN CYCLE zn ELSE found = .TRUE. EXIT xn ENDIF ENDDO zn ENDDO yn ENDDO xn IF ( .NOT. found ) THEN message_string = 'no matching grid for transpositions found' CALL message( 'init_pegrid', 'PA0224', 1, 2, 0, 6, 0 ) ENDIF ENDIF ! !-- Calculate array bounds in x-direction for every PE. !-- The last PE along x may get less grid points than the others ALLOCATE( nxlf(0:pdims(1)-1), nxrf(0:pdims(1)-1), nynf(0:pdims(2)-1), & nysf(0:pdims(2)-1), nnx_pe(0:pdims(1)-1), nny_pe(0:pdims(2)-1) ) IF ( MOD( nxa+1 , pdims(1) ) /= 0 ) THEN WRITE( message_string, * ) 'x-direction: gridpoint number (',nx+1,') ',& 'is not an& integral divisor of the number ', & 'processors (', pdims(1),')' CALL message( 'init_pegrid', 'PA0225', 1, 2, 0, 6, 0 ) ELSE nnx = ( nxa + 1 ) / pdims(1) IF ( nnx*pdims(1) - ( nx + 1) > nnx ) THEN WRITE( message_string, * ) 'x-direction: nx does not match the', & 'requirements given by the number of PEs &used', & '& please use nx = ', nx - ( pdims(1) - ( nnx*pdims(1) & - ( nx + 1 ) ) ), ' instead of nx =', nx CALL message( 'init_pegrid', 'PA0226', 1, 2, 0, 6, 0 ) ENDIF ENDIF ! !-- Left and right array bounds, number of gridpoints DO i = 0, pdims(1)-1 nxlf(i) = i * nnx nxrf(i) = ( i + 1 ) * nnx - 1 nnx_pe(i) = MIN( nx, nxrf(i) ) - nxlf(i) + 1 ENDDO ! !-- Calculate array bounds in y-direction for every PE. IF ( MOD( nya+1 , pdims(2) ) /= 0 ) THEN WRITE( message_string, * ) 'y-direction: gridpoint number (',ny+1,') ', & 'is not an& integral divisor of the number of', & 'processors (', pdims(2),')' CALL message( 'init_pegrid', 'PA0227', 1, 2, 0, 6, 0 ) ELSE nny = ( nya + 1 ) / pdims(2) IF ( nny*pdims(2) - ( ny + 1) > nny ) THEN WRITE( message_string, * ) 'y-direction: ny does not match the', & 'requirements given by the number of PEs &used ', & '& please use ny = ', ny - ( pdims(2) - ( nnx*pdims(2) & - ( ny + 1 ) ) ), ' instead of ny =', ny CALL message( 'init_pegrid', 'PA0228', 1, 2, 0, 6, 0 ) ENDIF ENDIF ! !-- South and north array bounds DO j = 0, pdims(2)-1 nysf(j) = j * nny nynf(j) = ( j + 1 ) * nny - 1 nny_pe(j) = MIN( ny, nynf(j) ) - nysf(j) + 1 ENDDO ! !-- Local array bounds of the respective PEs nxl = nxlf(pcoord(1)) nxra = nxrf(pcoord(1)) nxr = MIN( nx, nxra ) nys = nysf(pcoord(2)) nyna = nynf(pcoord(2)) nyn = MIN( ny, nyna ) nzb = 0 nzta = nza nzt = MIN( nz, nzta ) nnz = nza ! !-- Set switches to define if the PE is situated at the border of the virtual !-- processor grid IF ( nxl == 0 ) left_border_pe = .TRUE. IF ( nxr == nx ) right_border_pe = .TRUE. IF ( nys == 0 ) south_border_pe = .TRUE. IF ( nyn == ny ) north_border_pe = .TRUE. ! !-- Calculate array bounds and gridpoint numbers for the transposed arrays !-- (needed in the pressure solver) !-- For the transposed arrays, cyclic boundaries as well as top and bottom !-- boundaries are omitted, because they are obstructive to the transposition ! !-- 1. transposition z --> x !-- This transposition is not neccessary in case of a 1d-decomposition along x, !-- except that the uptream-spline method is switched on IF ( pdims(2) /= 1 .OR. momentum_advec == 'ups-scheme' .OR. & scalar_advec == 'ups-scheme' ) THEN IF ( pdims(2) == 1 .AND. ( momentum_advec == 'ups-scheme' .OR. & scalar_advec == 'ups-scheme' ) ) THEN message_string = '1d-decomposition along x ' // & 'chosen but nz restrictions may occur' // & '& since ups-scheme is activated' CALL message( 'init_pegrid', 'PA0229', 0, 1, 0, 6, 0 ) ENDIF nys_x = nys nyn_xa = nyna nyn_x = nyn nny_x = nny IF ( MOD( nza , pdims(1) ) /= 0 ) THEN WRITE( message_string, * ) 'transposition z --> x:', & '&nz=',nz,' is not an integral divisior of pdims(1)=', & pdims(1) CALL message( 'init_pegrid', 'PA0230', 1, 2, 0, 6, 0 ) ENDIF nnz_x = nza / pdims(1) nzb_x = 1 + myidx * nnz_x nzt_xa = ( myidx + 1 ) * nnz_x nzt_x = MIN( nzt, nzt_xa ) sendrecvcount_zx = nnx * nny * nnz_x ELSE ! !--- Setting of dummy values because otherwise variables are undefined in !--- the next step x --> y !--- WARNING: This case has still to be clarified!!!!!!!!!!!! nnz_x = 1 nzb_x = 1 nzt_xa = 1 nzt_x = 1 nny_x = nny ENDIF ! !-- 2. transposition x --> y nnz_y = nnz_x nzb_y = nzb_x nzt_ya = nzt_xa nzt_y = nzt_x IF ( MOD( nxa+1 , pdims(2) ) /= 0 ) THEN WRITE( message_string, * ) 'transposition x --> y:', & '&nx+1=',nx+1,' is not an integral divisor of ',& 'pdims(2)=',pdims(2) CALL message( 'init_pegrid', 'PA0231', 1, 2, 0, 6, 0 ) ENDIF nnx_y = (nxa+1) / pdims(2) nxl_y = myidy * nnx_y nxr_ya = ( myidy + 1 ) * nnx_y - 1 nxr_y = MIN( nx, nxr_ya ) sendrecvcount_xy = nnx_y * nny_x * nnz_y ! !-- 3. transposition y --> z (ELSE: x --> y in case of 1D-decomposition !-- along x) IF ( pdims(2) /= 1 .OR. momentum_advec == 'ups-scheme' .OR. & scalar_advec == 'ups-scheme' ) THEN ! !-- y --> z !-- This transposition is not neccessary in case of a 1d-decomposition !-- along x, except that the uptream-spline method is switched on nnx_z = nnx_y nxl_z = nxl_y nxr_za = nxr_ya nxr_z = nxr_y IF ( MOD( nya+1 , pdims(1) ) /= 0 ) THEN WRITE( message_string, * ) 'transposition y --> z:', & '& ny+1=',ny+1,' is not an integral divisor of ',& 'pdims(1)=',pdims(1) CALL message( 'init_pegrid', 'PA0232', 1, 2, 0, 6, 0 ) ENDIF nny_z = (nya+1) / pdims(1) nys_z = myidx * nny_z nyn_za = ( myidx + 1 ) * nny_z - 1 nyn_z = MIN( ny, nyn_za ) sendrecvcount_yz = nnx_y * nny_z * nnz_y ELSE ! !-- x --> y. This condition must be fulfilled for a 1D-decomposition along x IF ( MOD( nya+1 , pdims(1) ) /= 0 ) THEN WRITE( message_string, * ) 'transposition x --> y:', & '& ny+1=',ny+1,' is not an integral divisor of ',& 'pdims(1)=',pdims(1) CALL message( 'init_pegrid', 'PA0233', 1, 2, 0, 6, 0 ) ENDIF ENDIF ! !-- Indices for direct transpositions z --> y (used for calculating spectra) IF ( dt_dosp /= 9999999.9 ) THEN IF ( MOD( nza, pdims(2) ) /= 0 ) THEN WRITE( message_string, * ) 'direct transposition z --> y (needed ', & 'for spectra):& nz=',nz,' is not an integral divisor of ',& 'pdims(2)=',pdims(2) CALL message( 'init_pegrid', 'PA0234', 1, 2, 0, 6, 0 ) ELSE nxl_yd = nxl nxr_yda = nxra nxr_yd = nxr nzb_yd = 1 + myidy * ( nza / pdims(2) ) nzt_yda = ( myidy + 1 ) * ( nza / pdims(2) ) nzt_yd = MIN( nzt, nzt_yda ) sendrecvcount_zyd = nnx * nny * ( nza / pdims(2) ) ENDIF ENDIF ! !-- Indices for direct transpositions y --> x (they are only possible in case !-- of a 1d-decomposition along x) IF ( pdims(2) == 1 ) THEN nny_x = nny / pdims(1) nys_x = myid * nny_x nyn_xa = ( myid + 1 ) * nny_x - 1 nyn_x = MIN( ny, nyn_xa ) nzb_x = 1 nzt_xa = nza nzt_x = nz sendrecvcount_xy = nnx * nny_x * nza ENDIF ! !-- Indices for direct transpositions x --> y (they are only possible in case !-- of a 1d-decomposition along y) IF ( pdims(1) == 1 ) THEN nnx_y = nnx / pdims(2) nxl_y = myid * nnx_y nxr_ya = ( myid + 1 ) * nnx_y - 1 nxr_y = MIN( nx, nxr_ya ) nzb_y = 1 nzt_ya = nza nzt_y = nz sendrecvcount_xy = nnx_y * nny * nza ENDIF ! !-- Arrays for storing the array bounds are needed any more DEALLOCATE( nxlf , nxrf , nynf , nysf ) ! !-- Collect index bounds from other PEs (to be written to restart file later) ALLOCATE( hor_index_bounds(4,0:numprocs-1) ) IF ( myid == 0 ) THEN hor_index_bounds(1,0) = nxl hor_index_bounds(2,0) = nxr hor_index_bounds(3,0) = nys hor_index_bounds(4,0) = nyn ! !-- Receive data from all other PEs DO i = 1, numprocs-1 CALL MPI_RECV( ibuf, 4, MPI_INTEGER, i, MPI_ANY_TAG, comm2d, status, & ierr ) hor_index_bounds(:,i) = ibuf(1:4) ENDDO ELSE ! !-- Send index bounds to PE0 ibuf(1) = nxl ibuf(2) = nxr ibuf(3) = nys ibuf(4) = nyn CALL MPI_SEND( ibuf, 4, MPI_INTEGER, 0, myid, comm2d, ierr ) ENDIF #if defined( __print ) ! !-- Control output IF ( myid == 0 ) THEN PRINT*, '*** processor topology ***' PRINT*, ' ' PRINT*, 'myid pcoord left right south north idx idy nxl: nxr',& &' nys: nyn' PRINT*, '------------------------------------------------------------',& &'-----------' WRITE (*,1000) 0, pcoord(1), pcoord(2), pleft, pright, psouth, pnorth, & myidx, myidy, nxl, nxr, nys, nyn 1000 FORMAT (I4,2X,'(',I3,',',I3,')',3X,I4,2X,I4,3X,I4,2X,I4,2X,I3,1X,I3, & 2(2X,I4,':',I4)) ! !-- Receive data from the other PEs DO i = 1,numprocs-1 CALL MPI_RECV( ibuf, 12, MPI_INTEGER, i, MPI_ANY_TAG, comm2d, status, & ierr ) WRITE (*,1000) i, ( ibuf(j) , j = 1,12 ) ENDDO ELSE ! !-- Send data to PE0 ibuf(1) = pcoord(1); ibuf(2) = pcoord(2); ibuf(3) = pleft ibuf(4) = pright; ibuf(5) = psouth; ibuf(6) = pnorth; ibuf(7) = myidx ibuf(8) = myidy; ibuf(9) = nxl; ibuf(10) = nxr; ibuf(11) = nys ibuf(12) = nyn CALL MPI_SEND( ibuf, 12, MPI_INTEGER, 0, myid, comm2d, ierr ) ENDIF #endif #if defined( __parallel ) #if defined( __mpi2 ) ! !-- In case of coupled runs, get the port name on PE0 of the atmosphere model !-- and pass it to PE0 of the ocean model IF ( myid == 0 ) THEN IF ( coupling_mode == 'atmosphere_to_ocean' ) THEN CALL MPI_OPEN_PORT( MPI_INFO_NULL, port_name, ierr ) CALL MPI_PUBLISH_NAME( 'palm_coupler', MPI_INFO_NULL, port_name, & ierr ) ! !-- Write a flag file for the ocean model and the other atmosphere !-- processes. !-- There seems to be a bug in MPICH2 which causes hanging processes !-- in case that execution of LOOKUP_NAME is continued too early !-- (i.e. before the port has been created) OPEN( 90, FILE='COUPLING_PORT_OPENED', FORM='FORMATTED' ) WRITE ( 90, '(''TRUE'')' ) CLOSE ( 90 ) ELSEIF ( coupling_mode == 'ocean_to_atmosphere' ) THEN ! !-- Continue only if the atmosphere model has created the port. !-- There seems to be a bug in MPICH2 which causes hanging processes !-- in case that execution of LOOKUP_NAME is continued too early !-- (i.e. before the port has been created) INQUIRE( FILE='COUPLING_PORT_OPENED', EXIST=found ) DO WHILE ( .NOT. found ) INQUIRE( FILE='COUPLING_PORT_OPENED', EXIST=found ) ENDDO CALL MPI_LOOKUP_NAME( 'palm_coupler', MPI_INFO_NULL, port_name, ierr ) ENDIF ENDIF ! !-- In case of coupled runs, establish the connection between the atmosphere !-- and the ocean model and define the intercommunicator (comm_inter) CALL MPI_BARRIER( comm2d, ierr ) IF ( coupling_mode == 'atmosphere_to_ocean' ) THEN CALL MPI_COMM_ACCEPT( port_name, MPI_INFO_NULL, 0, MPI_COMM_WORLD, & comm_inter, ierr ) coupling_mode_remote = 'ocean_to_atmosphere' ELSEIF ( coupling_mode == 'ocean_to_atmosphere' ) THEN CALL MPI_COMM_CONNECT( port_name, MPI_INFO_NULL, 0, MPI_COMM_WORLD, & comm_inter, ierr ) coupling_mode_remote = 'atmosphere_to_ocean' ENDIF #endif ! !-- Determine the number of ghost point layers IF ( scalar_advec == 'ws-scheme' .OR. momentum_advec == 'ws-scheme' ) THEN nbgp = 3 ELSE nbgp = 1 ENDIF ! !-- Create a new MPI derived datatype for the exchange of surface (xy) data, !-- which is needed for coupled atmosphere-ocean runs. !-- First, calculate number of grid points of an xy-plane. ngp_xy = ( nxr - nxl + 1 + 2 * nbgp ) * ( nyn - nys + 1 + 2 * nbgp ) CALL MPI_TYPE_VECTOR( ngp_xy, 1, nzt-nzb+2, MPI_REAL, type_xy, ierr ) CALL MPI_TYPE_COMMIT( type_xy, ierr ) IF ( TRIM( coupling_mode ) /= 'uncoupled' ) THEN ! !-- Pass the number of grid points of the atmosphere model to !-- the ocean model and vice versa IF ( coupling_mode == 'atmosphere_to_ocean' ) THEN nx_a = nx ny_a = ny IF ( myid == 0 ) THEN CALL MPI_SEND( nx_a, 1, MPI_INTEGER, numprocs, 1, comm_inter, & ierr ) CALL MPI_SEND( ny_a, 1, MPI_INTEGER, numprocs, 2, comm_inter, & ierr ) CALL MPI_SEND( pdims, 2, MPI_INTEGER, numprocs, 3, comm_inter, & ierr ) CALL MPI_RECV( nx_o, 1, MPI_INTEGER, numprocs, 4, comm_inter, & status, ierr ) CALL MPI_RECV( ny_o, 1, MPI_INTEGER, numprocs, 5, comm_inter, & status, ierr ) CALL MPI_RECV( pdims_remote, 2, MPI_INTEGER, numprocs, 6, & comm_inter, status, ierr ) ENDIF CALL MPI_BCAST( nx_o, 1, MPI_INTEGER, 0, comm2d, ierr ) CALL MPI_BCAST( ny_o, 1, MPI_INTEGER, 0, comm2d, ierr ) CALL MPI_BCAST( pdims_remote, 2, MPI_INTEGER, 0, comm2d, ierr ) ELSEIF ( coupling_mode == 'ocean_to_atmosphere' ) THEN nx_o = nx ny_o = ny IF ( myid == 0 ) THEN CALL MPI_RECV( nx_a, 1, MPI_INTEGER, 0, 1, comm_inter, status, & ierr ) CALL MPI_RECV( ny_a, 1, MPI_INTEGER, 0, 2, comm_inter, status, & ierr ) CALL MPI_RECV( pdims_remote, 2, MPI_INTEGER, 0, 3, comm_inter, & status, ierr ) CALL MPI_SEND( nx_o, 1, MPI_INTEGER, 0, 4, comm_inter, ierr ) CALL MPI_SEND( ny_o, 1, MPI_INTEGER, 0, 5, comm_inter, ierr ) CALL MPI_SEND( pdims, 2, MPI_INTEGER, 0, 6, comm_inter, ierr ) ENDIF CALL MPI_BCAST( nx_a, 1, MPI_INTEGER, 0, comm2d, ierr) CALL MPI_BCAST( ny_a, 1, MPI_INTEGER, 0, comm2d, ierr) CALL MPI_BCAST( pdims_remote, 2, MPI_INTEGER, 0, comm2d, ierr) ENDIF ngp_a = ( nx_a+1 + 2 * nbgp ) * ( ny_a+1 + 2 * nbgp ) ngp_o = ( nx_o+1 + 2 * nbgp ) * ( ny_o+1 + 2 * nbgp ) ! !-- Determine if the horizontal grid and the number of PEs in ocean and !-- atmosphere is same or not IF ( nx_o == nx_a .AND. ny_o == ny_a .AND. & pdims(1) == pdims_remote(1) .AND. pdims(2) == pdims_remote(2) ) & THEN coupling_topology = 0 ELSE coupling_topology = 1 ENDIF ! !-- Determine the target PEs for the exchange between ocean and !-- atmosphere (comm2d) IF ( coupling_topology == 0 ) THEN ! !-- In case of identical topologies, every atmosphere PE has exactly one !-- ocean PE counterpart and vice versa IF ( TRIM( coupling_mode ) == 'atmosphere_to_ocean' ) THEN target_id = myid + numprocs ELSE target_id = myid ENDIF ELSE ! !-- In case of nonequivalent topology in ocean and atmosphere only for !-- PE0 in ocean and PE0 in atmosphere a target_id is needed, since !-- data echxchange between ocean and atmosphere will be done only !-- between these PEs. IF ( myid == 0 ) THEN IF ( TRIM( coupling_mode ) == 'atmosphere_to_ocean' ) THEN target_id = numprocs ELSE target_id = 0 ENDIF ENDIF ENDIF ENDIF #endif #else ! !-- Array bounds when running on a single PE (respectively a non-parallel !-- machine) nxl = 0 nxr = nx nxra = nx nnx = nxr - nxl + 1 nys = 0 nyn = ny nyna = ny nny = nyn - nys + 1 nzb = 0 nzt = nz nzta = nz nnz = nz ALLOCATE( hor_index_bounds(4,0:0) ) hor_index_bounds(1,0) = nxl hor_index_bounds(2,0) = nxr hor_index_bounds(3,0) = nys hor_index_bounds(4,0) = nyn ! !-- Array bounds for the pressure solver (in the parallel code, these bounds !-- are the ones for the transposed arrays) nys_x = nys nyn_x = nyn nyn_xa = nyn nzb_x = nzb + 1 nzt_x = nzt nzt_xa = nzt nxl_y = nxl nxr_y = nxr nxr_ya = nxr nzb_y = nzb + 1 nzt_y = nzt nzt_ya = nzt nxl_z = nxl nxr_z = nxr nxr_za = nxr nys_z = nys nyn_z = nyn nyn_za = nyn #endif ! !-- Calculate number of grid levels necessary for the multigrid poisson solver !-- as well as the gridpoint indices on each level IF ( psolver == 'multigrid' ) THEN ! !-- First calculate number of possible grid levels for the subdomains mg_levels_x = 1 mg_levels_y = 1 mg_levels_z = 1 i = nnx DO WHILE ( MOD( i, 2 ) == 0 .AND. i /= 2 ) i = i / 2 mg_levels_x = mg_levels_x + 1 ENDDO j = nny DO WHILE ( MOD( j, 2 ) == 0 .AND. j /= 2 ) j = j / 2 mg_levels_y = mg_levels_y + 1 ENDDO k = nz ! do not use nnz because it might be > nz due to transposition ! requirements DO WHILE ( MOD( k, 2 ) == 0 .AND. k /= 2 ) k = k / 2 mg_levels_z = mg_levels_z + 1 ENDDO maximum_grid_level = MIN( mg_levels_x, mg_levels_y, mg_levels_z ) ! !-- Find out, if the total domain allows more levels. These additional !-- levels are identically processed on all PEs. IF ( numprocs > 1 .AND. mg_switch_to_pe0_level /= -1 ) THEN IF ( mg_levels_z > MIN( mg_levels_x, mg_levels_y ) ) THEN mg_switch_to_pe0_level_l = maximum_grid_level mg_levels_x = 1 mg_levels_y = 1 i = nx+1 DO WHILE ( MOD( i, 2 ) == 0 .AND. i /= 2 ) i = i / 2 mg_levels_x = mg_levels_x + 1 ENDDO j = ny+1 DO WHILE ( MOD( j, 2 ) == 0 .AND. j /= 2 ) j = j / 2 mg_levels_y = mg_levels_y + 1 ENDDO maximum_grid_level_l = MIN( mg_levels_x, mg_levels_y, mg_levels_z ) IF ( maximum_grid_level_l > mg_switch_to_pe0_level_l ) THEN mg_switch_to_pe0_level_l = maximum_grid_level_l - & mg_switch_to_pe0_level_l + 1 ELSE mg_switch_to_pe0_level_l = 0 ENDIF ELSE mg_switch_to_pe0_level_l = 0 maximum_grid_level_l = maximum_grid_level ENDIF ! !-- Use switch level calculated above only if it is not pre-defined !-- by user IF ( mg_switch_to_pe0_level == 0 ) THEN IF ( mg_switch_to_pe0_level_l /= 0 ) THEN mg_switch_to_pe0_level = mg_switch_to_pe0_level_l maximum_grid_level = maximum_grid_level_l ENDIF ELSE ! !-- Check pre-defined value and reset to default, if neccessary IF ( mg_switch_to_pe0_level < mg_switch_to_pe0_level_l .OR. & mg_switch_to_pe0_level >= maximum_grid_level_l ) THEN message_string = 'mg_switch_to_pe0_level ' // & 'out of range and reset to default (=0)' CALL message( 'init_pegrid', 'PA0235', 0, 1, 0, 6, 0 ) mg_switch_to_pe0_level = 0 ELSE ! !-- Use the largest number of possible levels anyway and recalculate !-- the switch level to this largest number of possible values maximum_grid_level = maximum_grid_level_l ENDIF ENDIF ENDIF ALLOCATE( grid_level_count(maximum_grid_level), & nxl_mg(maximum_grid_level), nxr_mg(maximum_grid_level), & nyn_mg(maximum_grid_level), nys_mg(maximum_grid_level), & nzt_mg(maximum_grid_level) ) grid_level_count = 0 nxl_l = nxl; nxr_l = nxr; nys_l = nys; nyn_l = nyn; nzt_l = nzt DO i = maximum_grid_level, 1 , -1 IF ( i == mg_switch_to_pe0_level ) THEN #if defined( __parallel ) ! !-- Save the grid size of the subdomain at the switch level, because !-- it is needed in poismg. ind(1) = nxl_l; ind(2) = nxr_l ind(3) = nys_l; ind(4) = nyn_l ind(5) = nzt_l ALLOCATE( ind_all(5*numprocs), mg_loc_ind(5,0:numprocs-1) ) CALL MPI_ALLGATHER( ind, 5, MPI_INTEGER, ind_all, 5, & MPI_INTEGER, comm2d, ierr ) DO j = 0, numprocs-1 DO k = 1, 5 mg_loc_ind(k,j) = ind_all(k+j*5) ENDDO ENDDO DEALLOCATE( ind_all ) ! !-- Calculate the grid size of the total domain nxr_l = ( nxr_l-nxl_l+1 ) * pdims(1) - 1 nxl_l = 0 nyn_l = ( nyn_l-nys_l+1 ) * pdims(2) - 1 nys_l = 0 ! !-- The size of this gathered array must not be larger than the !-- array tend, which is used in the multigrid scheme as a temporary !-- array subdomain_size = ( nxr - nxl + 3 ) * ( nyn - nys + 3 ) * & ( nzt - nzb + 2 ) gathered_size = ( nxr_l - nxl_l + 3 ) * ( nyn_l - nys_l + 3 ) * & ( nzt_l - nzb + 2 ) IF ( gathered_size > subdomain_size ) THEN message_string = 'not enough memory for storing ' // & 'gathered multigrid data on PE0' CALL message( 'init_pegrid', 'PA0236', 1, 2, 0, 6, 0 ) ENDIF #else message_string = 'multigrid gather/scatter impossible ' // & 'in non parallel mode' CALL message( 'init_pegrid', 'PA0237', 1, 2, 0, 6, 0 ) #endif ENDIF nxl_mg(i) = nxl_l nxr_mg(i) = nxr_l nys_mg(i) = nys_l nyn_mg(i) = nyn_l nzt_mg(i) = nzt_l nxl_l = nxl_l / 2 nxr_l = nxr_l / 2 nys_l = nys_l / 2 nyn_l = nyn_l / 2 nzt_l = nzt_l / 2 ENDDO ELSE maximum_grid_level = 0 ENDIF ! !-- Default level 0 tells exchange_horiz that all ghost planes have to be !-- exchanged. grid_level is adjusted in poismg, where only one ghost plane !-- is required. grid_level = 0 #if defined( __parallel ) ! !-- Gridpoint number for the exchange of ghost points (y-line for 2D-arrays) ngp_y = nyn - nys + 1 + 2 * nbgp ! !-- Define new MPI derived datatypes for the exchange of ghost points in !-- x- and y-direction for 2D-arrays (line) CALL MPI_TYPE_VECTOR( nxr-nxl+1+2*nbgp, nbgp, ngp_y, MPI_REAL, type_x, & ierr ) CALL MPI_TYPE_COMMIT( type_x, ierr ) CALL MPI_TYPE_VECTOR( nxr-nxl+1+2*nbgp, nbgp, ngp_y, MPI_INTEGER, & type_x_int, ierr ) CALL MPI_TYPE_COMMIT( type_x_int, ierr ) CALL MPI_TYPE_VECTOR( nbgp, ngp_y, ngp_y, MPI_REAL, type_y, ierr ) CALL MPI_TYPE_COMMIT( type_y, ierr ) CALL MPI_TYPE_VECTOR( nbgp, ngp_y, ngp_y, MPI_INTEGER, type_y_int, ierr ) CALL MPI_TYPE_COMMIT( type_y_int, ierr ) ! !-- Calculate gridpoint numbers for the exchange of ghost points along x !-- (yz-plane for 3D-arrays) and define MPI derived data type(s) for the !-- exchange of ghost points in y-direction (xz-plane). !-- Do these calculations for the model grid and (if necessary) also !-- for the coarser grid levels used in the multigrid method ALLOCATE ( ngp_yz(0:maximum_grid_level), type_xz(0:maximum_grid_level),& type_yz(0:maximum_grid_level) ) nxl_l = nxl; nxr_l = nxr; nys_l = nys; nyn_l = nyn; nzb_l = nzb; nzt_l = nzt ! !-- Discern between the model grid, which needs nbgp ghost points and !-- grid levels for the multigrid scheme. In the latter case only one !-- ghost point is necessary. !-- First definition of MPI-datatypes for exchange of ghost layers on normal !-- grid. The following loop is needed for data exchange in poismg.f90. ! !-- Determine number of grid points of yz-layer for exchange ngp_yz(0) = (nzt - nzb + 2) * (nyn - nys + 1 + 2 * nbgp) ! !-- Define an MPI-datatype for the exchange of left/right boundaries. !-- Although data are contiguous in physical memory (which does not !-- necessarily require an MPI-derived datatype), the data exchange between !-- left and right PE's using the MPI-derived type is 10% faster than without. CALL MPI_TYPE_VECTOR( nxr-nxl+1+2*nbgp, nbgp*(nzt-nzb+2), ngp_yz(0), & MPI_REAL, type_xz(0), ierr ) CALL MPI_TYPE_COMMIT( type_xz(0), ierr ) CALL MPI_TYPE_VECTOR( nbgp, ngp_yz(0), ngp_yz(0), MPI_REAL, type_yz(0), & ierr ) CALL MPI_TYPE_COMMIT( type_yz(0), ierr ) ! !-- Definition of MPI-datatypes for multigrid method (coarser level grids) IF ( psolver == 'multigrid' ) THEN ! !-- Definition of MPI-datatyoe as above, but only 1 ghost level is used DO i = maximum_grid_level, 1 , -1 ngp_yz(i) = (nzt_l - nzb_l + 2) * (nyn_l - nys_l + 3) CALL MPI_TYPE_VECTOR( nxr_l-nxl_l+3, nzt_l-nzb_l+2, ngp_yz(i), & MPI_REAL, type_xz(i), ierr ) CALL MPI_TYPE_COMMIT( type_xz(i), ierr ) CALL MPI_TYPE_VECTOR( 1, ngp_yz(i), ngp_yz(i), MPI_REAL, type_yz(i), & ierr ) CALL MPI_TYPE_COMMIT( type_yz(i), ierr ) nxl_l = nxl_l / 2 nxr_l = nxr_l / 2 nys_l = nys_l / 2 nyn_l = nyn_l / 2 nzt_l = nzt_l / 2 ENDDO ENDIF #endif #if defined( __parallel ) ! !-- Setting of flags for inflow/outflow conditions in case of non-cyclic !-- horizontal boundary conditions. IF ( pleft == MPI_PROC_NULL ) THEN IF ( bc_lr == 'dirichlet/radiation' ) THEN inflow_l = .TRUE. ELSEIF ( bc_lr == 'radiation/dirichlet' ) THEN outflow_l = .TRUE. ENDIF ENDIF IF ( pright == MPI_PROC_NULL ) THEN IF ( bc_lr == 'dirichlet/radiation' ) THEN outflow_r = .TRUE. ELSEIF ( bc_lr == 'radiation/dirichlet' ) THEN inflow_r = .TRUE. ENDIF ENDIF IF ( psouth == MPI_PROC_NULL ) THEN IF ( bc_ns == 'dirichlet/radiation' ) THEN outflow_s = .TRUE. ELSEIF ( bc_ns == 'radiation/dirichlet' ) THEN inflow_s = .TRUE. ENDIF ENDIF IF ( pnorth == MPI_PROC_NULL ) THEN IF ( bc_ns == 'dirichlet/radiation' ) THEN inflow_n = .TRUE. ELSEIF ( bc_ns == 'radiation/dirichlet' ) THEN outflow_n = .TRUE. ENDIF ENDIF ! !-- Broadcast the id of the inflow PE IF ( inflow_l ) THEN id_inflow_l = myidx ELSE id_inflow_l = 0 ENDIF IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) CALL MPI_ALLREDUCE( id_inflow_l, id_inflow, 1, MPI_INTEGER, MPI_SUM, & comm1dx, ierr ) ! !-- Broadcast the id of the recycling plane !-- WARNING: needs to be adjusted in case of inflows other than from left side! IF ( ( recycling_width / dx ) >= nxl .AND. & ( recycling_width / dx ) <= nxr ) THEN id_recycling_l = myidx ELSE id_recycling_l = 0 ENDIF IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) CALL MPI_ALLREDUCE( id_recycling_l, id_recycling, 1, MPI_INTEGER, MPI_SUM, & comm1dx, ierr ) #else IF ( bc_lr == 'dirichlet/radiation' ) THEN inflow_l = .TRUE. outflow_r = .TRUE. ELSEIF ( bc_lr == 'radiation/dirichlet' ) THEN outflow_l = .TRUE. inflow_r = .TRUE. ENDIF IF ( bc_ns == 'dirichlet/radiation' ) THEN inflow_n = .TRUE. outflow_s = .TRUE. ELSEIF ( bc_ns == 'radiation/dirichlet' ) THEN outflow_n = .TRUE. inflow_s = .TRUE. ENDIF #endif ! !-- At the outflow, u or v, respectively, have to be calculated for one more !-- grid point. IF ( outflow_l ) THEN nxlu = nxl + 1 ELSE nxlu = nxl ENDIF IF ( outflow_s ) THEN nysv = nys + 1 ELSE nysv = nys ENDIF IF ( psolver == 'poisfft_hybrid' ) THEN CALL poisfft_hybrid_ini ELSEIF ( psolver == 'poisfft' ) THEN CALL poisfft_init ENDIF ! !-- Allocate wall flag arrays used in the multigrid solver IF ( psolver == 'multigrid' ) THEN DO i = maximum_grid_level, 1, -1 SELECT CASE ( i ) CASE ( 1 ) ALLOCATE( wall_flags_1(nzb:nzt_mg(i)+1, & nys_mg(i)-1:nyn_mg(i)+1, & nxl_mg(i)-1:nxr_mg(i)+1) ) CASE ( 2 ) ALLOCATE( wall_flags_2(nzb:nzt_mg(i)+1, & nys_mg(i)-1:nyn_mg(i)+1, & nxl_mg(i)-1:nxr_mg(i)+1) ) CASE ( 3 ) ALLOCATE( wall_flags_3(nzb:nzt_mg(i)+1, & nys_mg(i)-1:nyn_mg(i)+1, & nxl_mg(i)-1:nxr_mg(i)+1) ) CASE ( 4 ) ALLOCATE( wall_flags_4(nzb:nzt_mg(i)+1, & nys_mg(i)-1:nyn_mg(i)+1, & nxl_mg(i)-1:nxr_mg(i)+1) ) CASE ( 5 ) ALLOCATE( wall_flags_5(nzb:nzt_mg(i)+1, & nys_mg(i)-1:nyn_mg(i)+1, & nxl_mg(i)-1:nxr_mg(i)+1) ) CASE ( 6 ) ALLOCATE( wall_flags_6(nzb:nzt_mg(i)+1, & nys_mg(i)-1:nyn_mg(i)+1, & nxl_mg(i)-1:nxr_mg(i)+1) ) CASE ( 7 ) ALLOCATE( wall_flags_7(nzb:nzt_mg(i)+1, & nys_mg(i)-1:nyn_mg(i)+1, & nxl_mg(i)-1:nxr_mg(i)+1) ) CASE ( 8 ) ALLOCATE( wall_flags_8(nzb:nzt_mg(i)+1, & nys_mg(i)-1:nyn_mg(i)+1, & nxl_mg(i)-1:nxr_mg(i)+1) ) CASE ( 9 ) ALLOCATE( wall_flags_9(nzb:nzt_mg(i)+1, & nys_mg(i)-1:nyn_mg(i)+1, & nxl_mg(i)-1:nxr_mg(i)+1) ) CASE ( 10 ) ALLOCATE( wall_flags_10(nzb:nzt_mg(i)+1, & nys_mg(i)-1:nyn_mg(i)+1, & nxl_mg(i)-1:nxr_mg(i)+1) ) CASE DEFAULT message_string = 'more than 10 multigrid levels' CALL message( 'init_pegrid', 'PA0238', 1, 2, 0, 6, 0 ) END SELECT ENDDO ENDIF END SUBROUTINE init_pegrid