Changeset 14367
- Timestamp:
- 2021-02-02T08:51:42+01:00 (4 years ago)
- Location:
- NEMO/branches/2021/dev_r14312_MPI_Interface/src/OCE/LBC
- Files:
-
- 3 edited
Legend:
- Unmodified
- Added
- Removed
-
NEMO/branches/2021/dev_r14312_MPI_Interface/src/OCE/LBC/lbc_lnk_neicoll_generic.h90
r14366 r14367 17 17 INTEGER :: ip0j, ip1j, im0j, im1j 18 18 INTEGER :: ishti, ishtj, ishti2, ishtj2 19 INTEGER :: isz s, iszr19 INTEGER :: iszS, iszR 20 20 INTEGER :: ierr 21 21 INTEGER :: ihls, idx … … 26 26 INTEGER, DIMENSION(8) :: isizej, ishtSj, ishtRj, ishtPj 27 27 INTEGER, DIMENSION(8) :: ifill, iszall 28 INTEGER, DIMENSION(:), ALLOCATABLE :: i counts, icountr! number of elements to be sent/received29 INTEGER, DIMENSION(:), ALLOCATABLE :: i displs, idisplr! displacement in halos arrays28 INTEGER, DIMENSION(:), ALLOCATABLE :: iScnt, iRcnt ! number of elements to be sent/received 29 INTEGER, DIMENSION(:), ALLOCATABLE :: iSdpl, iRdpl ! displacement in halos arrays 30 30 LOGICAL, DIMENSION(8) :: llsend, llrecv 31 31 REAL(PRECISION) :: zland 32 REAL(PRECISION), DIMENSION(:), ALLOCATABLE :: zsnd, zrcv ! halos arrays 33 LOGICAL :: ll4only ! default: 8 neighbourgs 32 LOGICAL :: ll4only ! default: 8 neighbourgs 34 33 LOGICAL :: ll_IdoNFold 35 34 !!---------------------------------------------------------------------- … … 140 139 ! 141 140 ! Allocate local temporary arrays to be sent/received. 142 isz s= COUNT( llsend )143 isz r= COUNT( llrecv )144 ALLOCATE( i counts(iszs), icountr(iszr), idispls(iszs), idisplr(iszr) ) ! ok if iszs = 0 or iszr= 0141 iszS = COUNT( llsend ) 142 iszR = COUNT( llrecv ) 143 ALLOCATE( iScnt(iszS), iRcnt(iszR), iSdpl(iszS), iRdpl(iszR) ) ! ok if iszS = 0 or iszR = 0 145 144 iszall(:) = isizei(:) * isizej(:) * ipk * ipl * ipf 146 i counts(:) = PACK( iszall, mask = llsend ) ! ok if mask = .false.147 i countr(:) = PACK( iszall, mask = llrecv )148 i displs(1) = 0149 DO jn = 2,isz s150 i displs(jn) = idispls(jn-1) + icounts(jn-1) ! with _alltoallv: in units of sendtype151 END DO 152 i displr(1) = 0153 DO jn = 2,isz r154 i displr(jn) = idisplr(jn-1) + icountr(jn-1) ! with _alltoallv: in units of sendtype145 iScnt(:) = PACK( iszall, mask = llsend ) ! ok if mask = .false. 146 iRcnt(:) = PACK( iszall, mask = llrecv ) 147 iSdpl(1) = 0 148 DO jn = 2,iszS 149 iSdpl(jn) = iSdpl(jn-1) + iScnt(jn-1) ! with _alltoallv: in units of sendtype 150 END DO 151 iRdpl(1) = 0 152 DO jn = 2,iszR 153 iRdpl(jn) = iRdpl(jn-1) + iRcnt(jn-1) ! with _alltoallv: in units of sendtype 155 154 END DO 156 155 157 ALLOCATE( zsnd(SUM(icounts)), zrcv(SUM(icountr)) ) 156 ! Allocate buffer arrays to be sent/received if needed 157 iszS = SUM(iszall, mask = llsend) ! send buffer size 158 IF( ALLOCATED(BUFFSND) ) THEN 159 IF( SIZE(BUFFSND) < iszS ) DEALLOCATE(BUFFSND) ! send buffer is too small 160 ENDIF 161 IF( .NOT. ALLOCATED(BUFFSND) ) ALLOCATE( BUFFSND(iszS) ) 162 iszR = SUM(iszall, mask = llrecv) ! recv buffer size 163 IF( ALLOCATED(BUFFRCV) ) THEN 164 IF( SIZE(BUFFRCV) < iszR ) DEALLOCATE(BUFFRCV) ! recv buffer is too small 165 ENDIF 166 IF( .NOT. ALLOCATED(BUFFRCV) ) ALLOCATE( BUFFRCV(iszR) ) 158 167 159 168 ! fill sending buffer with ptab(jf)%pt4d … … 164 173 ishtj = ishtSj(jn) 165 174 DO jf = 1, ipf ; DO jl = 1, ipl ; DO jk = 1, ipk ; DO jj = 1,isizej(jn) ; DO ji = 1,isizei(jn) 166 zsnd(idx) = ptab(jf)%pt4d(ishti+ji,ishtj+jj,jk,jl)175 BUFFSND(idx) = ptab(jf)%pt4d(ishti+ji,ishtj+jj,jk,jl) 167 176 idx = idx + 1 168 177 END DO ; END DO ; END DO ; END DO ; END DO … … 175 184 ! 176 185 IF( ln_timing ) CALL tic_tac(.TRUE.) 177 CALL mpi_neighbor_alltoallv ( zsnd, icounts, idispls, MPI_TYPE, zrcv, icountr, idisplr, MPI_TYPE, impi_nc, ierr)186 CALL mpi_neighbor_alltoallv (BUFFSND, iScnt, iSdpl, MPI_TYPE, BUFFRCV, iRcnt, iRdpl, MPI_TYPE, impi_nc, ierr) 178 187 IF( ln_timing ) CALL tic_tac(.FALSE.) 179 188 ! … … 190 199 CASE ( jpfillmpi ) ! fill with data received by MPI 191 200 DO jf = 1, ipf ; DO jl = 1, ipl ; DO jk = 1, ipk ; DO jj = 1,isizej(jn) ; DO ji = 1,isizei(jn) 192 ptab(jf)%pt4d(ishti+ji,ishtj+jj,jk,jl) = zrcv(idx)201 ptab(jf)%pt4d(ishti+ji,ishtj+jj,jk,jl) = BUFFRCV(idx) 193 202 idx = idx + 1 194 203 END DO ; END DO ; END DO ; END DO ; END DO … … 212 221 END DO 213 222 214 DEALLOCATE( icounts, icountr, idispls, idisplr, zsnd, zrcv ) 223 DEALLOCATE( iScnt, iRcnt, iSdpl, iRdpl ) 224 IF( iszS > jpi*jpj ) DEALLOCATE(BUFFSND) ! blocking Send -> can directly deallocate 225 IF( iszR > jpi*jpj ) DEALLOCATE(BUFFRCV) ! blocking Recv -> can directly deallocate 215 226 216 227 ! potential "indirect self-periodicity" for the corners -
NEMO/branches/2021/dev_r14312_MPI_Interface/src/OCE/LBC/lbc_lnk_pt2pt_generic.h90
r14366 r14367 1 1 2 2 SUBROUTINE lbc_lnk_pt2pt_/**/PRECISION( cdname, ptab, cd_nat, psgn, kfld, kfillmode, pfillval, khls, lsend, lrecv, ld4only ) 3 3 CHARACTER(len=*) , INTENT(in ) :: cdname ! name of the calling subroutine … … 18 18 INTEGER :: ishti, ishtj, ishti2, ishtj2 19 19 INTEGER :: ifill_nfd, icomm, ierr 20 INTEGER :: ihls, idxs, idxr 20 INTEGER :: ihls, idxs, idxr, iszS, iszR 21 21 INTEGER, DIMENSION(4) :: iwewe, issnn 22 22 INTEGER, DIMENSION(8) :: isizei, ishtSi, ishtRi, ishtPi … … 26 26 INTEGER, DIMENSION(8) :: iStag, iRtag ! Send and Recv mpi_tag id 27 27 REAL(PRECISION) :: zland 28 REAL(PRECISION), DIMENSION(:), ALLOCATABLE :: zsnd, zrcv ! halos arrays29 28 LOGICAL, DIMENSION(8) :: llsend, llrecv 30 29 LOGICAL :: ll4only ! default: 8 neighbourgs … … 137 136 ! -------------------------------- ! 138 137 ! 139 ireq(:) = MPI_REQUEST_NULL ! default definition when no communication is done. understood by mpi_waitall140 138 iStag = (/ 1, 2, 3, 4, 5, 6, 7, 8 /) ! any value but each one must be different 141 139 ! define iRtag with the corresponding iStag, e.g. data received at west where sent at east. … … 153 151 END DO 154 152 155 ! Allocate local temporary arrays to be sent/received. 156 ALLOCATE( zsnd( SUM(iszall, mask = llsend) ), zrcv( SUM(iszall, mask = llrecv) ) ) 153 ! Allocate buffer arrays to be sent/received if needed 154 iszS = SUM(iszall, mask = llsend) ! send buffer size 155 IF( ALLOCATED(BUFFSND) ) THEN 156 CALL mpi_waitall(8, nreq_p2p, MPI_STATUSES_IGNORE, ierr) ! wait for Isend from the PREVIOUS call 157 IF( SIZE(BUFFSND) < iszS ) DEALLOCATE(BUFFSND) ! send buffer is too small 158 ENDIF 159 IF( .NOT. ALLOCATED(BUFFSND) ) ALLOCATE( BUFFSND(iszS) ) 160 iszR = SUM(iszall, mask = llrecv) ! recv buffer size 161 IF( ALLOCATED(BUFFRCV) ) THEN 162 IF( SIZE(BUFFRCV) < iszR ) DEALLOCATE(BUFFRCV) ! recv buffer is too small 163 ENDIF 164 IF( .NOT. ALLOCATED(BUFFRCV) ) ALLOCATE( BUFFRCV(iszR) ) 165 ! 166 ! default definition when no communication is done. understood by mpi_waitall 167 nreq_p2p(:) = MPI_REQUEST_NULL ! WARNING: Must be done after the call to mpi_waitall just above 157 168 ! 158 169 ! ----------------------------------------------- ! … … 165 176 ishtj = ishtSj(jn) 166 177 DO jf = 1, ipf ; DO jl = 1, ipl ; DO jk = 1, ipk ; DO jj = 1,isizej(jn) ; DO ji = 1,isizei(jn) 167 zsnd(idxs) = ptab(jf)%pt4d(ishti+ji,ishtj+jj,jk,jl)178 BUFFSND(idxs) = ptab(jf)%pt4d(ishti+ji,ishtj+jj,jk,jl) 168 179 idxs = idxs + 1 169 180 END DO ; END DO ; END DO ; END DO ; END DO … … 171 182 IF( ln_timing ) CALL tic_tac(.TRUE.) 172 183 ! non-blocking send of the west/east side using local buffer 173 CALL MPI_ISEND( zsnd(ishtS(jn)+1), iszall(jn), MPI_TYPE, mpiSnei(ihls,jn), iStag(jn), icomm, ireq(jn), ierr )184 CALL MPI_ISEND( BUFFSND(ishtS(jn)+1), iszall(jn), MPI_TYPE, mpiSnei(ihls,jn), iStag(jn), icomm, nreq_p2p(jn), ierr ) 174 185 IF( ln_timing ) CALL tic_tac(.FALSE.) 175 186 #endif … … 190 201 IF( ln_timing ) CALL tic_tac(.TRUE.) 191 202 ! ! blocking receive of the west/east halo in local temporary arrays 192 CALL MPI_RECV( zrcv(ishtR(jn)+1), iszall(jn), MPI_TYPE, mpiRnei(ihls,jn), iRtag(jn), icomm, MPI_STATUS_IGNORE, ierr )193 IF( ln_timing ) CALL tic_tac(.FALSE.) 194 #endif 195 DO jf = 1, ipf ; DO jl = 1, ipl ; DO jk = 1, ipk ; DO jj = 1,isizej(jn) ; DO ji = 1,isizei(jn) 196 ptab(jf)%pt4d(ishti+ji,ishtj+jj,jk,jl) = zrcv(idxr)203 CALL MPI_RECV( BUFFRCV(ishtR(jn)+1), iszall(jn), MPI_TYPE, mpiRnei(ihls,jn), iRtag(jn), icomm, MPI_STATUS_IGNORE, ierr ) 204 IF( ln_timing ) CALL tic_tac(.FALSE.) 205 #endif 206 DO jf = 1, ipf ; DO jl = 1, ipl ; DO jk = 1, ipk ; DO jj = 1,isizej(jn) ; DO ji = 1,isizei(jn) 207 ptab(jf)%pt4d(ishti+ji,ishtj+jj,jk,jl) = BUFFRCV(idxr) 197 208 idxr = idxr + 1 198 209 END DO ; END DO ; END DO ; END DO ; END DO … … 237 248 ishtj = ishtSj(jn) 238 249 DO jf = 1, ipf ; DO jl = 1, ipl ; DO jk = 1, ipk ; DO jj = 1,isizej(jn) ; DO ji = 1,isizei(jn) 239 zsnd(idxs) = ptab(jf)%pt4d(ishti+ji,ishtj+jj,jk,jl)250 BUFFSND(idxs) = ptab(jf)%pt4d(ishti+ji,ishtj+jj,jk,jl) 240 251 idxs = idxs + 1 241 252 END DO ; END DO ; END DO ; END DO ; END DO … … 243 254 IF( ln_timing ) CALL tic_tac(.TRUE.) 244 255 ! non-blocking send of the south/north side using local buffer 245 CALL MPI_ISEND( zsnd(ishtS(jn)+1), iszall(jn), MPI_TYPE, mpiSnei(ihls,jn), iStag(jn), icomm, ireq(jn), ierr )256 CALL MPI_ISEND( BUFFSND(ishtS(jn)+1), iszall(jn), MPI_TYPE, mpiSnei(ihls,jn), iStag(jn), icomm, nreq_p2p(jn), ierr ) 246 257 IF( ln_timing ) CALL tic_tac(.FALSE.) 247 258 #endif … … 262 273 IF( ln_timing ) CALL tic_tac(.TRUE.) 263 274 ! ! blocking receive of the south/north halo in local temporary arrays 264 CALL MPI_RECV( zrcv(ishtR(jn)+1), iszall(jn), MPI_TYPE, mpiRnei(ihls,jn), iRtag(jn), icomm, MPI_STATUS_IGNORE, ierr )265 IF( ln_timing ) CALL tic_tac(.FALSE.) 266 #endif 267 DO jf = 1, ipf ; DO jl = 1, ipl ; DO jk = 1, ipk ; DO jj = 1,isizej(jn) ; DO ji = 1,isizei(jn) 268 ptab(jf)%pt4d(ishti+ji,ishtj+jj,jk,jl) = zrcv(idxr)275 CALL MPI_RECV( BUFFRCV(ishtR(jn)+1), iszall(jn), MPI_TYPE, mpiRnei(ihls,jn), iRtag(jn), icomm, MPI_STATUS_IGNORE, ierr ) 276 IF( ln_timing ) CALL tic_tac(.FALSE.) 277 #endif 278 DO jf = 1, ipf ; DO jl = 1, ipl ; DO jk = 1, ipk ; DO jj = 1,isizej(jn) ; DO ji = 1,isizei(jn) 279 ptab(jf)%pt4d(ishti+ji,ishtj+jj,jk,jl) = BUFFRCV(idxr) 269 280 idxr = idxr + 1 270 281 END DO ; END DO ; END DO ; END DO ; END DO … … 298 309 ishtj = ishtSj(jn) 299 310 DO jf = 1, ipf ; DO jl = 1, ipl ; DO jk = 1, ipk ; DO jj = 1,isizej(jn) ; DO ji = 1,isizei(jn) 300 zsnd(idxs) = ptab(jf)%pt4d(ishti+ji,ishtj+jj,jk,jl)311 BUFFSND(idxs) = ptab(jf)%pt4d(ishti+ji,ishtj+jj,jk,jl) 301 312 idxs = idxs + 1 302 313 END DO ; END DO ; END DO ; END DO ; END DO … … 304 315 IF( ln_timing ) CALL tic_tac(.TRUE.) 305 316 ! non-blocking send of the corners using local buffer 306 CALL MPI_ISEND( zsnd(ishtS(jn)+1), iszall(jn), MPI_TYPE, mpiSnei(ihls,jn), iStag(jn), icomm, ireq(jn), ierr )317 CALL MPI_ISEND( BUFFSND(ishtS(jn)+1), iszall(jn), MPI_TYPE, mpiSnei(ihls,jn), iStag(jn), icomm, nreq_p2p(jn), ierr ) 307 318 IF( ln_timing ) CALL tic_tac(.FALSE.) 308 319 #endif … … 314 325 IF( ln_timing ) CALL tic_tac(.TRUE.) 315 326 ! blocking receive of the corner halo in local temporary arrays 316 CALL MPI_RECV( zrcv(ishtR(jn)+1), iszall(jn), MPI_TYPE, mpiRnei(ihls,jn), iRtag(jn), icomm, MPI_STATUS_IGNORE, ierr )327 CALL MPI_RECV( BUFFRCV(ishtR(jn)+1), iszall(jn), MPI_TYPE, mpiRnei(ihls,jn), iRtag(jn), icomm, MPI_STATUS_IGNORE, ierr ) 317 328 IF( ln_timing ) CALL tic_tac(.FALSE.) 318 329 #endif … … 320 331 ishtj = ishtRj(jn) 321 332 DO jf = 1, ipf ; DO jl = 1, ipl ; DO jk = 1, ipk ; DO jj = 1,isizej(jn) ; DO ji = 1,isizei(jn) 322 ptab(jf)%pt4d(ishti+ji,ishtj+jj,jk,jl) = zrcv(idxr)333 ptab(jf)%pt4d(ishti+ji,ishtj+jj,jk,jl) = BUFFRCV(idxr) 323 334 idxr = idxr + 1 324 335 END DO ; END DO ; END DO ; END DO ; END DO … … 328 339 ! -------------------------------------------- ! 329 340 ! 9. deallocate local temporary arrays ! 341 ! if they areg larger than jpi*jpj ! <- arbitrary max size... 330 342 ! -------------------------------------------- ! 331 343 ! 332 CALL mpi_waitall(8, ireq, MPI_STATUSES_IGNORE, ierr) 333 DEALLOCATE( zsnd, zrcv ) 344 IF( iszR > jpi*jpj ) DEALLOCATE(BUFFRCV) ! blocking receive -> can directly deallocate 345 IF( iszS > jpi*jpj ) THEN 346 CALL mpi_waitall(8, nreq_p2p, MPI_STATUSES_IGNORE, ierr) ! must wait before deallocate send buffer 347 DEALLOCATE(BUFFSND) 348 ENDIF 334 349 ! 335 350 END SUBROUTINE lbc_lnk_pt2pt_/**/PRECISION -
NEMO/branches/2021/dev_r14312_MPI_Interface/src/OCE/LBC/lbclnk.F90
r14366 r14367 50 50 PUBLIC lbc_lnk_icb ! iceberg lateral boundary conditions 51 51 52 REAL(dp), DIMENSION(:), ALLOCATABLE :: buffsnd_dp, buffrcv_dp ! MPI send/recv buffers 53 REAL(sp), DIMENSION(:), ALLOCATABLE :: buffsnd_sp, buffrcv_sp ! 54 INTEGER, DIMENSION(8) :: nreq_p2p ! request id for MPI_Isend in point-2-point communication 55 52 56 !! * Substitutions 53 57 !!# include "do_loop_substitute.h90" … … 126 130 #define PRECISION sp 127 131 # define MPI_TYPE MPI_REAL 132 # define BUFFSND buffsnd_sp 133 # define BUFFRCV buffrcv_sp 128 134 # include "lbc_lnk_pt2pt_generic.h90" 129 135 # include "lbc_lnk_neicoll_generic.h90" 130 136 # undef MPI_TYPE 137 # undef BUFFSND 138 # undef BUFFRCV 131 139 #undef PRECISION 132 140 !! … … 135 143 #define PRECISION dp 136 144 # define MPI_TYPE MPI_DOUBLE_PRECISION 145 # define BUFFSND buffsnd_dp 146 # define BUFFRCV buffrcv_dp 137 147 # include "lbc_lnk_pt2pt_generic.h90" 138 148 # include "lbc_lnk_neicoll_generic.h90" 139 149 # undef MPI_TYPE 150 # undef BUFFSND 151 # undef BUFFRCV 140 152 #undef PRECISION 141 153
Note: See TracChangeset
for help on using the changeset viewer.