Parent Directory
|
Revision Log
|
Revision Graph
|
Patch
--- MITgcm_contrib/cg2d_bench/exch.F 2006/05/12 21:58:05 1.1
+++ MITgcm_contrib/cg2d_bench/exch.F 2006/05/12 22:22:11 1.2
@@ -1,3 +1,4 @@
+C $Id: exch.F,v 1.2 2006/05/12 22:22:11 ce107 Exp $
SUBROUTINE EXCH_XY_R8( arr )
C == Global variables ==
@@ -7,17 +8,20 @@
#ifdef ALLOW_MPI
#include "mpif.h"
+#include "MPI_INFO.h"
#endif
-#include "MPI_INFO.h"
#include "JAM_INFO.h"
C == Routine arguments ==
- Real*8 arr(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
+ Real arr(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
C == Local variables ==
INTEGER I, J
INTEGER northProc, southProc
+#ifdef DECOMP2D
+ INTEGER eastProc, westProc
+#endif
INTEGER farProc1, farProc2
INTEGER toPid, fromPid
INTEGER rc
@@ -25,9 +29,12 @@
INTEGER myFourWayRank, exchangePhase
#ifdef ALLOW_MPI
+ INTEGER mpiReq(8)
+ INTEGER mpiStat(MPI_STATUS_SIZE,8)
INTEGER mpiStatus(MPI_STATUS_SIZE)
#endif
+#ifndef DECOMP2D
C East-west halo update (without corners)
DO J=1,sNy
DO I=1,OLx
@@ -35,7 +42,8 @@
arr(sNx+I,J) = arr(1+I-1 ,J)
ENDDO
ENDDO
-
+#endif
+
C Phase 1 pairing
C | 0 | ---> | 1 |
C | 0 | <--- | 1 |
@@ -50,6 +58,58 @@
C
#ifdef USE_MPI_EXCH
+#ifdef DECOMP2D
+C East-West exchanges
+#ifdef USE_SNDRCV
+C Send West, receive from East
+ CALL MPI_Sendrecv(arr(1,1), 1, ewslice, mpi_westId, 100,
+ $ arr(sNx+1,1), 1, ewslice, mpi_eastId, 100, comm_use,
+ $ mpiStatus, rc)
+C Send East, receive from West
+ CALL MPI_Sendrecv(arr(sNx-OLx+1,1), 1, ewslice, mpi_eastId, 200,
+ $ arr(1-OLx,1), 1, ewslice, mpi_westId, 200, comm_use,
+ $ mpiStatus,rc)
+
+C North-South exchanges
+
+C Send South, receive from North
+ CALL MPI_Sendrecv(arr(1,1), 1, nsslice, mpi_southId, 300,
+ $ arr(1,sNy+1), 1, nsslice, mpi_northId, 300, comm_use,
+ $ mpiStatus, rc)
+C Send North, receive from South
+ CALL MPI_Sendrecv(arr(1,sNy-OLy+1), 1, nsslice, mpi_northId, 400,
+ $ arr(1,1-OLy), 1, nsslice, mpi_southId, 400, comm_use,
+ $ mpiStatus,rc)
+#else
+C Send West, receive from East
+ CALL MPI_Isend(arr(1,1), 1, ewslice, mpi_westId, 100,
+ $ comm_use, mpiReq(1), rc)
+ CALL MPI_Irecv(arr(sNx+1,1), 1, ewslice, mpi_eastId, 100,
+ $ comm_use, mpiReq(2), rc)
+C Send East, receive from West
+ CALL MPI_Isend(arr(sNx-OLx+1,1), 1, ewslice, mpi_eastId, 200,
+ $ comm_use, mpiReq(3), rc)
+ CALL MPI_Irecv(arr(1-OLx,1), 1, ewslice, mpi_westId, 200,
+ $ comm_use, mpiReq(4),rc)
+
+C North-South exchanges
+
+C Send South, receive from North
+ CALL MPI_Isend(arr(1,1), 1, nsslice, mpi_southId, 300,
+ $ comm_use, mpiReq(5), rc)
+ CALL MPI_Irecv(arr(1,sNy+1), 1, nsslice, mpi_northId, 300,
+ $ comm_use, mpiReq(6), rc)
+C Send North, receive from South
+ CALL MPI_Isend(arr(1,sNy-OLy+1), 1, nsslice, mpi_northId, 400,
+ $ comm_use, mpiReq(7), rc)
+ CALL MPI_Irecv(arr(1,1-OLy), 1, nsslice, mpi_southId, 400,
+ $ comm_use, mpiReq(8),rc)
+
+ CALL MPI_Waitall(8, mpiReq, mpiStat, rc)
+
+#endif
+
+#else
C North-south halo update (without corners)
C Put my edges into a buffers
IF ( MOD(myProcId,2) .EQ. 0 ) THEN
@@ -76,21 +136,22 @@
ENDIF
C Even-odd pairs
IF ( farProc1 .NE. myProcId ) THEN
- CALL MPI_Sendrecv_replace(exchBuf1,sNx,MPI_DOUBLE_PRECISION,
+ CALL MPI_Sendrecv_replace(exchBuf1,sNx,_MPI_TYPE_REAL,
& farProc1,0,
& farProc1,MPI_ANY_TAG,
- & MPI_COMM_WORLD,mpiStatus,
+ & comm_use,mpiStatus,
& rc)
ENDIF
C Odd-even pairs
IF ( farProc2 .NE. myProcId ) THEN
- CALL MPI_Sendrecv_replace(exchBuf2,sNx,MPI_DOUBLE_PRECISION,
+ CALL MPI_Sendrecv_replace(exchBuf2,sNx,_MPI_TYPE_REAL,
& farProc2,0,
& farProc2,MPI_ANY_TAG,
- & MPI_COMM_WORLD,mpiStatus,
+ & comm_use,mpiStatus,
& rc)
ENDIF
#endif
+#endif
#ifdef USE_JAM_EXCH
@@ -106,7 +167,8 @@
farProc1 = northProc
farProc2 = southProc
IF ( farProc1 .NE. myProcId ) THEN
- CALL JAM_EXCHANGE(farProc1,arr(I,sNy),arr(I,sNy+1),sNx*8,jam_exchKey)
+ CALL JAM_EXCHANGE(farProc1,arr(I,sNy),arr(I,sNy+1),sNx*8
+ $ ,jam_exchKey)
jam_exchKey = jam_exchKey+1
ENDIF
@@ -150,7 +212,8 @@
21 CONTINUE
IF ( farProc2 .NE. myProcId ) THEN
- CALL JAM_EXCHANGE(farProc2,arr(I,sNy),arr(I,sNy+1),sNx*8,jam_exchKey)
+ CALL JAM_EXCHANGE(farProc2,arr(I,sNy),arr(I,sNy+1),sNx*8
+ $ ,jam_exchKey)
jam_exchKey = jam_exchKey+1
ENDIF
@@ -159,7 +222,7 @@
ENDIF
#endif
-#ifdef USE_MPI_EXCH
+#if defined(USE_MPI_EXCH) && !defined(DECOMP2D)
C Fill overlap regions from the buffers
IF ( MOD(myProcId,2) .EQ. 0 ) THEN
DO I=1,sNx
| ViewVC Help | |
| Powered by ViewVC 1.1.22 |