1 |
C $Header: /u/gcmpack/MITgcm/eesupp/src/scatter_2d.F,v 1.8 2008/11/26 08:41:57 mlosch Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "CPP_EEOPTIONS.h" |
5 |
|
6 |
SUBROUTINE SCATTER_2D( global, local, myThid ) |
7 |
C Scatter elements of a 2-D array from mpi process 0 to all processes. |
8 |
IMPLICIT NONE |
9 |
#include "SIZE.h" |
10 |
#include "EEPARAMS.h" |
11 |
#include "EESUPPORT.h" |
12 |
C mythid - thread number for this instance of the routine. |
13 |
C global,local - working arrays used to transfer 2-D fields |
14 |
INTEGER mythid |
15 |
Real*8 global(Nx,Ny) |
16 |
_RL local(1:sNx,1:sNy,nSx,nSy) |
17 |
|
18 |
INTEGER iG,jG, i, j, bi, bj |
19 |
#ifdef ALLOW_USE_MPI |
20 |
_RL temp(1:sNx,1:sNy,nSx,nSy) |
21 |
INTEGER istatus(MPI_STATUS_SIZE), ierr |
22 |
INTEGER isource, itag, npe |
23 |
INTEGER lbuff |
24 |
#endif /* ALLOW_USE_MPI */ |
25 |
|
26 |
C-- Make everyone wait except for master thread. |
27 |
_BARRIER |
28 |
_BEGIN_MASTER( myThid ) |
29 |
|
30 |
#ifndef ALLOW_USE_MPI |
31 |
|
32 |
DO bj=1,nSy |
33 |
DO bi=1,nSx |
34 |
DO j=1,sNy |
35 |
DO i=1,sNx |
36 |
iG = myXGlobalLo-1+(bi-1)*sNx+i |
37 |
jG = myYGlobalLo-1+(bj-1)*sNy+j |
38 |
local(i,j,bi,bj) = global(iG,jG) |
39 |
ENDDO |
40 |
ENDDO |
41 |
ENDDO |
42 |
ENDDO |
43 |
|
44 |
#else /* ALLOW_USE_MPI */ |
45 |
|
46 |
lbuff=sNx*nSx*sNy*nSy |
47 |
isource = 0 |
48 |
itag = 0 |
49 |
|
50 |
IF( mpiMyId .EQ. 0 ) THEN |
51 |
|
52 |
C-- Process 0 fills-in its local data |
53 |
npe = 0 |
54 |
DO bj=1,nSy |
55 |
DO bi=1,nSx |
56 |
DO j=1,sNy |
57 |
#ifdef TARGET_NEC_SX |
58 |
!cdir novector |
59 |
#endif |
60 |
DO i=1,sNx |
61 |
iG = mpi_myXGlobalLo(npe+1)-1+(bi-1)*sNx+i |
62 |
jG = mpi_myYGlobalLo(npe+1)-1+(bj-1)*sNy+j |
63 |
local(i,j,bi,bj) = global(iG,jG) |
64 |
ENDDO |
65 |
ENDDO |
66 |
ENDDO |
67 |
ENDDO |
68 |
|
69 |
C-- Process 0 sends local arrays to all other processes |
70 |
DO npe = 1, numberOfProcs-1 |
71 |
DO bj=1,nSy |
72 |
DO bi=1,nSx |
73 |
DO j=1,sNy |
74 |
#ifdef TARGET_NEC_SX |
75 |
!cdir novector |
76 |
#endif |
77 |
DO i=1,sNx |
78 |
iG = mpi_myXGlobalLo(npe+1)-1+(bi-1)*sNx+i |
79 |
jG = mpi_myYGlobalLo(npe+1)-1+(bj-1)*sNy+j |
80 |
temp(i,j,bi,bj) = global(iG,jG) |
81 |
ENDDO |
82 |
ENDDO |
83 |
ENDDO |
84 |
ENDDO |
85 |
CALL MPI_SEND (temp, lbuff, MPI_DOUBLE_PRECISION, |
86 |
& npe, itag, MPI_COMM_MODEL, ierr) |
87 |
ENDDO |
88 |
|
89 |
ELSE |
90 |
|
91 |
C-- All proceses except 0 receive local array from process 0 |
92 |
CALL MPI_RECV (local, lbuff, MPI_DOUBLE_PRECISION, |
93 |
& isource, itag, MPI_COMM_MODEL, istatus, ierr) |
94 |
|
95 |
ENDIF |
96 |
|
97 |
#endif /* ALLOW_USE_MPI */ |
98 |
|
99 |
_END_MASTER( myThid ) |
100 |
_BARRIER |
101 |
|
102 |
RETURN |
103 |
END |