Ignore:
Timestamp:
Feb 9, 2011 2:25:15 PM (13 years ago)
Author:
raasch
Message:

New:
---

optional exchange of ghost points in synchronous mode via MPI_SENDRCV,
steered by d3par parameter synchronous_exchange
(cpu_statistics, exchange_horiz, modules, parin)

openMP-parallelization of pressure solver (fft-method) for 2d-domain-decomposition
(poisfft, transpose)

Changed:


Errors:


mpt bugfix for netCDF4 usage (mrun)

File:
1 edited

Legend:

Unmodified
Added
Removed
  • palm/trunk/SOURCE/exchange_horiz.f90

    r668 r683  
    44! Current revisions:
    55! -----------------
     6! optional synchronous exchange (sendrecv) implemented, code partly reformatted
    67!
    78! Former revisions:
     
    4849    INTEGER, DIMENSION(MPI_STATUS_SIZE,4) ::  wait_stat
    4950#endif
    50     INTEGER :: i,nbgp_local
     51    INTEGER ::  i, nbgp_local
    5152    REAL, DIMENSION(nzb:nzt+1,nys-nbgp_local:nyn+nbgp_local, &
    5253                    nxl-nbgp_local:nxr+nbgp_local) ::  ar
     
    5455    CALL cpu_log( log_point_s(2), 'exchange_horiz', 'start' )
    5556
    56     IF ( exchange_mg == .TRUE. ) THEN
    57       i = grid_level
     57!
     58!-- In the Poisson multigrid solver arrays with coarser grids are used.
     59!-- Set i appropriately, because the coarser grids have different
     60!-- MPI datatypes type_xz, type_yz.
     61    IF ( exchange_mg == .TRUE. )  THEN
     62       i = grid_level
    5863    ELSE
    59       i = 0
     64       i = 0
    6065    END IF
     66
    6167#if defined( __parallel )
    6268
     
    7480    ELSE
    7581
    76        req = 0
     82       IF ( synchronous_exchange )  THEN
    7783!
    78 !--    Send left boundary, receive right one
    79        CALL MPI_ISEND(ar(nzb,nys-nbgp_local,nxl),1,type_yz(i),pleft,0,comm2d,&
    80                       req(1),ierr)
    81        CALL MPI_IRECV(ar(nzb,nys-nbgp_local,nxr+1),1,type_yz(i),pright,0,&
    82                      comm2d,req(2),ierr)
     84!--       Send left boundary, receive right one (synchronous)
     85          CALL MPI_SENDRECV(                                                   &
     86                       ar(nzb,nys-nbgp_local,nxl),   1, type_yz(i), pleft,  0, &
     87                       ar(nzb,nys-nbgp_local,nxr+1), 1, type_yz(i), pright, 0, &
     88                       comm2d, status, ierr )
    8389!
    84 !--    Send right boundary, receive left one
     90!--       Send right boundary, receive left one (synchronous)
     91          CALL MPI_SENDRECV(                                                   &
     92            ar(nzb,nys-nbgp_local,nxr+1-nbgp_local), 1, type_yz(i), pright, 1, &
     93            ar(nzb,nys-nbgp_local,nxl-nbgp_local),   1, type_yz(i), pleft,  1, &
     94                       comm2d, status, ierr )
    8595
     96       ELSE
    8697
    87        CALL MPI_ISEND(ar(nzb,nys-nbgp_local,nxr+1-nbgp_local),1,type_yz(i),pright, 1,  &
    88                       comm2d, req(3), ierr )
    89        CALL MPI_IRECV(ar(nzb,nys-nbgp_local,nxl-nbgp_local),1,type_yz(i),pleft,1,&
    90                       comm2d,req(4), ierr)
     98          req = 0
     99!
     100!--       Send left boundary, receive right one (asynchronous)
     101          CALL MPI_ISEND( ar(nzb,nys-nbgp_local,nxl),   1, type_yz(i), pleft,  &
     102                          0, comm2d, req(1), ierr )
     103          CALL MPI_IRECV( ar(nzb,nys-nbgp_local,nxr+1), 1, type_yz(i), pright, &
     104                          0, comm2d, req(2), ierr )
     105!
     106!--       Send right boundary, receive left one (asynchronous)
     107          CALL MPI_ISEND( ar(nzb,nys-nbgp_local,nxr+1-nbgp_local), 1,          &
     108                          type_yz(i), pright, 1, comm2d, req(3), ierr )
     109          CALL MPI_IRECV( ar(nzb,nys-nbgp_local,nxl-nbgp_local),   1,          &
     110                          type_yz(i), pleft,  1, comm2d, req(4), ierr )
    91111
    92        CALL MPI_WAITALL( 4, req, wait_stat, ierr )
     112          CALL MPI_WAITALL( 4, req, wait_stat, ierr )
     113
     114       ENDIF
    93115
    94116    ENDIF
     
    106128    ELSE
    107129
    108        req = 0
     130       IF ( synchronous_exchange )  THEN
    109131!
    110 !--    Send front boundary, receive rear one
    111 !--    MPI_ISEND initial send adress changed, type_xz() is sendet nbgp times
     132!--       Send front boundary, receive rear one (synchronous)
     133          CALL MPI_SENDRECV(                                                   &
     134                       ar(nzb,nys,nxl-nbgp_local),   1, type_xz(i), psouth, 0, &
     135                       ar(nzb,nyn+1,nxl-nbgp_local), 1, type_xz(i), pnorth, 0, &
     136                       comm2d, status, ierr )
     137!
     138!--       Send rear boundary, receive front one (synchronous)
     139          CALL MPI_SENDRECV(                                                   &
     140            ar(nzb,nyn-nbgp_local+1,nxl-nbgp_local), 1, type_xz(i), pnorth, 1, &
     141            ar(nzb,nys-nbgp_local,nxl-nbgp_local),   1, type_xz(i), psouth, 1, &
     142            comm2d, status, ierr )
    112143
    113        CALL MPI_ISEND( ar(nzb,nys,nxl-nbgp_local),1, type_xz(i), psouth, 0, &
    114                        comm2d, req(1), ierr )
    115        CALL MPI_IRECV( ar(nzb,nyn+1,nxl-nbgp_local),1, type_xz(i), pnorth, 0, &
    116                        comm2d, req(2), ierr )
     144       ELSE
     145
     146          req = 0
    117147!
    118 !--    Send rear boundary, receive front one
    119        CALL MPI_ISEND( ar(nzb,nyn-nbgp_local+1,nxl-nbgp_local),1, type_xz(i), pnorth, 1, &
    120                        comm2d, req(3), ierr )
    121        CALL MPI_IRECV( ar(nzb,nys-nbgp_local,nxl-nbgp_local),1, type_xz(i), psouth, 1, &
    122                        comm2d, req(4), ierr )
    123        call MPI_WAITALL( 4, req, wait_stat, ierr )
     148!--       Send front boundary, receive rear one (asynchronous)
     149          CALL MPI_ISEND( ar(nzb,nys,nxl-nbgp_local),   1, type_xz(i), psouth, &
     150                          0, comm2d, req(1), ierr )
     151          CALL MPI_IRECV( ar(nzb,nyn+1,nxl-nbgp_local), 1, type_xz(i), pnorth, &
     152                          0, comm2d, req(2), ierr )
     153!
     154!--       Send rear boundary, receive front one (asynchronous)
     155          CALL MPI_ISEND( ar(nzb,nyn-nbgp_local+1,nxl-nbgp_local), 1,          &
     156                          type_xz(i), pnorth, 1, comm2d, req(3), ierr )
     157          CALL MPI_IRECV( ar(nzb,nys-nbgp_local,nxl-nbgp_local),   1,          &
     158                          type_xz(i), psouth, 1, comm2d, req(4), ierr )
     159
     160          CALL MPI_WAITALL( 4, req, wait_stat, ierr )
     161
     162       ENDIF
    124163
    125164    ENDIF
Note: See TracChangeset for help on using the changeset viewer.