1 |
jmc |
1.3 |
C $Header: /u/gcmpack/MITgcm/eesupp/src/scatter_xz.F,v 1.2 2007/09/23 22:52:50 jmc Exp $ |
2 |
jmc |
1.2 |
C $Name: $ |
3 |
|
|
|
4 |
|
|
#include "CPP_EEOPTIONS.h" |
5 |
heimbach |
1.1 |
|
6 |
|
|
SUBROUTINE SCATTER_XZ( global, local, myThid ) |
7 |
|
|
C Scatter elements of a x-z 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 |
jmc |
1.2 |
INTEGER mythid |
15 |
heimbach |
1.1 |
Real*8 global(Nx) |
16 |
|
|
_RL local(1-OLx:sNx+OLx,nSx,nSy) |
17 |
|
|
|
18 |
|
|
INTEGER iG, i, bi, bj |
19 |
|
|
#ifdef ALLOW_USE_MPI |
20 |
|
|
|
21 |
|
|
_RL temp(1-OLx:sNx+OLx,nSx,nSy) |
22 |
|
|
|
23 |
|
|
INTEGER istatus(MPI_STATUS_SIZE), ierr |
24 |
|
|
INTEGER isource, itag, npe |
25 |
|
|
INTEGER lbuff |
26 |
|
|
#endif /* ALLOW_USE_MPI */ |
27 |
|
|
|
28 |
|
|
C-- Make everyone wait except for master thread. |
29 |
|
|
_BARRIER |
30 |
|
|
_BEGIN_MASTER( myThid ) |
31 |
|
|
|
32 |
|
|
#ifndef ALLOW_USE_MPI |
33 |
|
|
|
34 |
|
|
DO bj=1,nSy |
35 |
|
|
DO bi=1,nSx |
36 |
|
|
DO i=1,sNx |
37 |
|
|
iG = myXGlobalLo-1+(bi-1)*sNx+i |
38 |
|
|
local(i,bi,bj) = global(iG) |
39 |
|
|
ENDDO |
40 |
|
|
ENDDO |
41 |
|
|
ENDDO |
42 |
|
|
|
43 |
|
|
#else /* ALLOW_USE_MPI */ |
44 |
|
|
|
45 |
|
|
lbuff=(sNx+2*OLx)*nSx*nSy |
46 |
|
|
isource = 0 |
47 |
|
|
itag = 0 |
48 |
|
|
|
49 |
|
|
IF( mpiMyId .EQ. 0 ) THEN |
50 |
|
|
|
51 |
|
|
C-- Process 0 fills-in its local data |
52 |
|
|
npe = 0 |
53 |
|
|
DO bj=1,nSy |
54 |
|
|
DO bi=1,nSx |
55 |
|
|
DO i=1,sNx |
56 |
|
|
iG = mpi_myXGlobalLo(npe+1)-1+(bi-1)*sNx+i |
57 |
|
|
local(i,bi,bj) = global(iG) |
58 |
|
|
ENDDO |
59 |
|
|
ENDDO |
60 |
|
|
ENDDO |
61 |
|
|
|
62 |
|
|
C-- Process 0 sends local arrays to all other processes |
63 |
|
|
DO npe = 1, numberOfProcs-1 |
64 |
|
|
DO bj=1,nSy |
65 |
|
|
DO bi=1,nSx |
66 |
|
|
DO i=1,sNx |
67 |
|
|
iG = mpi_myXGlobalLo(npe+1)-1+(bi-1)*sNx+i |
68 |
|
|
temp(i,bi,bj) = global(iG) |
69 |
|
|
ENDDO |
70 |
|
|
ENDDO |
71 |
|
|
ENDDO |
72 |
|
|
CALL MPI_SEND (temp, lbuff, MPI_DOUBLE_PRECISION, |
73 |
|
|
& npe, itag, MPI_COMM_MODEL, ierr) |
74 |
|
|
ENDDO |
75 |
|
|
|
76 |
|
|
ELSE |
77 |
|
|
|
78 |
|
|
C-- All proceses except 0 receive local array from process 0 |
79 |
|
|
CALL MPI_RECV (local, lbuff, MPI_DOUBLE_PRECISION, |
80 |
|
|
& isource, itag, MPI_COMM_MODEL, istatus, ierr) |
81 |
|
|
|
82 |
|
|
ENDIF |
83 |
|
|
|
84 |
|
|
#endif /* ALLOW_USE_MPI */ |
85 |
|
|
|
86 |
|
|
_END_MASTER( myThid ) |
87 |
|
|
_BARRIER |
88 |
|
|
|
89 |
|
|
C-- Fill in edges. |
90 |
jmc |
1.3 |
CMM _EXCH_XY_RL( local, myThid ) |
91 |
heimbach |
1.1 |
|
92 |
|
|
RETURN |
93 |
|
|
END |