File indexing completed on 2018-03-02 18:36:09 UTC
view on githubraw file Latest commit add29e06 on 2018-01-31 20:35:05 UTC
463a0fddee Patr*0001 #include "CPP_OPTIONS.h"
0002
0003 SUBROUTINE GATHER_YZ( global, local, myThid )
0004
0005 IMPLICIT NONE
0006 #include "SIZE.h"
0007 #include "EEPARAMS.h"
0008 #include "EESUPPORT.h"
0009
0010
0011 INTEGER mythid
0012 Real*8 global(Ny)
0013 _RL local(1-OLy:sNy+OLy,nSx,nSy)
0014
0015 INTEGER jG, j, bi, bj
0016 #ifdef ALLOW_USE_MPI
0017
0018 _RL temp(1-OLy:sNy+OLy,nSx,nSy)
0019
0020 INTEGER istatus(MPI_STATUS_SIZE), ierr
0021 INTEGER lbuff, idest, itag, npe, ready_to_receive
0022 #endif /* ALLOW_USE_MPI */
0023
0024
0025 _BARRIER
0026 _BEGIN_MASTER( myThid )
0027
0028 #ifndef ALLOW_USE_MPI
0029
0030 DO bj=1,nSy
0031 DO bi=1,nSx
0032 DO j=1,sNy
0033 jG = myYGlobalLo-1+(bi-1)*sNy+j
0034 global(jG) = local(j,bi,bj)
0035 ENDDO
0036 ENDDO
0037 ENDDO
0038
0039 #else /* ALLOW_USE_MPI */
0040
0041 lbuff = (sNy+2*OLy)*nSx*nSy
0042 idest = 0
0043 itag = 0
0044 ready_to_receive = 0
0045
0046 IF( mpiMyId .EQ. 0 ) THEN
0047
0048
0049 npe = 0
0050 DO bj=1,nSy
0051 DO bi=1,nSx
0052 DO j=1,sNy
0053 jG = mpi_myYGlobalLo(npe+1)-1+(bi-1)*sNy+j
0054 global(jG) = local(j,bi,bj)
0055 ENDDO
0056 ENDDO
0057 ENDDO
0058
0059
0060 DO npe = 1, numberOfProcs-1
f5cbf7b96d Dimi*0061 #ifndef DISABLE_MPI_READY_TO_RECEIVE
463a0fddee Patr*0062 CALL MPI_SEND (ready_to_receive, 1, MPI_INTEGER,
0063 & npe, itag, MPI_COMM_MODEL, ierr)
f5cbf7b96d Dimi*0064 #endif
463a0fddee Patr*0065 CALL MPI_RECV (temp, lbuff, MPI_DOUBLE_PRECISION,
0066 & npe, itag, MPI_COMM_MODEL, istatus, ierr)
0067
0068
0069 DO bj=1,nSy
0070 DO bi=1,nSx
0071 DO j=1,sNy
0072 jG = mpi_myYGlobalLo(npe+1)-1+(bi-1)*sNy+j
0073 global(jG) = temp(j,bi,bj)
0074 ENDDO
0075 ENDDO
0076 ENDDO
0077 ENDDO
0078
0079 ELSE
0080
0081
f5cbf7b96d Dimi*0082 #ifndef DISABLE_MPI_READY_TO_RECEIVE
463a0fddee Patr*0083 CALL MPI_RECV (ready_to_receive, 1, MPI_INTEGER,
0084 & idest, itag, MPI_COMM_MODEL, istatus, ierr)
f5cbf7b96d Dimi*0085 #endif
463a0fddee Patr*0086 CALL MPI_SEND (local, lbuff, MPI_DOUBLE_PRECISION,
0087 & idest, itag, MPI_COMM_MODEL, ierr)
0088
0089 ENDIF
0090
0091 #endif /* ALLOW_USE_MPI */
0092
0093 _END_MASTER( myThid )
0094 _BARRIER
0095
0096 RETURN
0097 END