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