Changeset 683 for palm/trunk/SOURCE/exchange_horiz.f90
- Timestamp:
- Feb 9, 2011 2:25:15 PM (13 years ago)
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
palm/trunk/SOURCE/exchange_horiz.f90
r668 r683 4 4 ! Current revisions: 5 5 ! ----------------- 6 ! optional synchronous exchange (sendrecv) implemented, code partly reformatted 6 7 ! 7 8 ! Former revisions: … … 48 49 INTEGER, DIMENSION(MPI_STATUS_SIZE,4) :: wait_stat 49 50 #endif 50 INTEGER :: i,nbgp_local51 INTEGER :: i, nbgp_local 51 52 REAL, DIMENSION(nzb:nzt+1,nys-nbgp_local:nyn+nbgp_local, & 52 53 nxl-nbgp_local:nxr+nbgp_local) :: ar … … 54 55 CALL cpu_log( log_point_s(2), 'exchange_horiz', 'start' ) 55 56 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 58 63 ELSE 59 i = 064 i = 0 60 65 END IF 66 61 67 #if defined( __parallel ) 62 68 … … 74 80 ELSE 75 81 76 req = 082 IF ( synchronous_exchange ) THEN 77 83 ! 78 !-- Send left boundary, receive right one79 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 ) 83 89 ! 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 ) 85 95 96 ELSE 86 97 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 ) 91 111 92 CALL MPI_WAITALL( 4, req, wait_stat, ierr ) 112 CALL MPI_WAITALL( 4, req, wait_stat, ierr ) 113 114 ENDIF 93 115 94 116 ENDIF … … 106 128 ELSE 107 129 108 req = 0130 IF ( synchronous_exchange ) THEN 109 131 ! 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 ) 112 143 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 117 147 ! 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 124 163 125 164 ENDIF
Note: See TracChangeset
for help on using the changeset viewer.