1 |
C $Header: /u/u0/gcmpack/MITgcm/eesupp/src/exch_rx_recv_get_y.template,v 1.6 2005/11/09 17:22:08 cnh Exp $ |
2 |
C $Name: $ |
3 |
#include "CPP_EEOPTIONS.h" |
4 |
|
5 |
CBOP |
6 |
C !ROUTINE: EXCH_RX_RECV_GET_Y |
7 |
|
8 |
C !INTERFACE: |
9 |
SUBROUTINE EXCH_RX_RECV_GET_Y( array, |
10 |
I myOLw, myOLe, myOLs, myOLn, myNz, |
11 |
I exchWidthX, exchWidthY, |
12 |
I theSimulationMode, theCornerMode, myThid ) |
13 |
IMPLICIT NONE |
14 |
|
15 |
C !DESCRIPTION: |
16 |
C *==========================================================* |
17 |
C | SUBROUTINE RECV_GET_Y |
18 |
C | o "Send" or "put" Y edges for RX array. |
19 |
C *==========================================================* |
20 |
C | Routine that invokes actual message passing send or |
21 |
C | direct "put" of data to update X faces of an XY[R] array. |
22 |
C *==========================================================* |
23 |
|
24 |
C !USES: |
25 |
C == Global variables == |
26 |
#include "SIZE.h" |
27 |
#include "EEPARAMS.h" |
28 |
#include "EESUPPORT.h" |
29 |
#include "EXCH.h" |
30 |
|
31 |
C !INPUT/OUTPUT PARAMETERS: |
32 |
C == Routine arguments == |
33 |
C array :: Array with edges to exchange. |
34 |
C myOLw :: West, East, North and South overlap region sizes. |
35 |
C myOLe |
36 |
C myOLn |
37 |
C myOLs |
38 |
C exchWidthX :: Width of data region exchanged. |
39 |
C exchWidthY |
40 |
C theSimulationMode :: Forward or reverse mode exchange ( provides |
41 |
C support for adjoint integration of code. ) |
42 |
C theCornerMode :: Flag indicating whether corner updates are |
43 |
C needed. |
44 |
C myThid :: Thread number of this instance of S/R EXCH... |
45 |
C eBl :: Edge buffer level |
46 |
INTEGER myOLw |
47 |
INTEGER myOLe |
48 |
INTEGER myOLs |
49 |
INTEGER myOLn |
50 |
INTEGER myNz |
51 |
_RX array(1-myOLw:sNx+myOLe, |
52 |
& 1-myOLs:sNy+myOLn, |
53 |
& myNZ, nSx, nSy) |
54 |
INTEGER exchWidthX |
55 |
INTEGER exchWidthY |
56 |
INTEGER theSimulationMode |
57 |
INTEGER theCornerMode |
58 |
INTEGER myThid |
59 |
|
60 |
C !LOCAL VARIABLES: |
61 |
C == Local variables == |
62 |
C I, J, K, iMin, iMax, iB :: Loop counters and extents |
63 |
C bi, bj |
64 |
C biS, bjS :: South tile indices |
65 |
C biN, bjN :: North tile indices |
66 |
C eBl :: Current exchange buffer level |
67 |
C theProc, theTag, theType, :: Variables used in message building |
68 |
C theSize |
69 |
C southCommMode :: Working variables holding type |
70 |
C northCommMode of communication a particular |
71 |
C tile face uses. |
72 |
C spinCount :: Exchange statistics counter |
73 |
C mpiStatus :: MPI error code |
74 |
INTEGER I, J, K, iMin, iMax, jMin, jMax, iB, iB0 |
75 |
INTEGER bi, bj, biS, bjS, biN, bjN |
76 |
INTEGER eBl |
77 |
INTEGER southCommMode |
78 |
INTEGER northCommMode |
79 |
INTEGER spinCount |
80 |
#ifdef ALLOW_USE_MPI |
81 |
INTEGER theProc, theTag, theType, theSize |
82 |
INTEGER mpiStatus(MPI_STATUS_SIZE,4), mpiRc |
83 |
#endif |
84 |
CEOP |
85 |
|
86 |
INTEGER myBxLoSave(MAX_NO_THREADS) |
87 |
INTEGER myBxHiSave(MAX_NO_THREADS) |
88 |
INTEGER myByLoSave(MAX_NO_THREADS) |
89 |
INTEGER myByHiSave(MAX_NO_THREADS) |
90 |
LOGICAL doingSingleThreadedComms |
91 |
|
92 |
doingSingleThreadedComms = .FALSE. |
93 |
#ifdef ALLOW_USE_MPI |
94 |
#ifndef ALWAYS_USE_MPI |
95 |
IF ( usingMPI ) THEN |
96 |
#endif |
97 |
C Set default behavior to have MPI comms done by a single thread. |
98 |
C Most MPI implementations don't support concurrent comms from |
99 |
C several threads. |
100 |
IF ( nThreads .GT. 1 ) THEN |
101 |
_BARRIER |
102 |
_BEGIN_MASTER( myThid ) |
103 |
DO I=1,nThreads |
104 |
myBxLoSave(I) = myBxLo(I) |
105 |
myBxHiSave(I) = myBxHi(I) |
106 |
myByLoSave(I) = myByLo(I) |
107 |
myByHiSave(I) = myByHi(I) |
108 |
ENDDO |
109 |
C Comment out loop below and myB[xy][Lo|Hi](1) settings below |
110 |
C if you want to get multi-threaded MPI comms. |
111 |
DO I=1,nThreads |
112 |
myBxLo(I) = 0 |
113 |
myBxHi(I) = -1 |
114 |
myByLo(I) = 0 |
115 |
myByHi(I) = -1 |
116 |
ENDDO |
117 |
myBxLo(1) = 1 |
118 |
myBxHi(1) = nSx |
119 |
myByLo(1) = 1 |
120 |
myByHi(1) = nSy |
121 |
doingSingleThreadedComms = .TRUE. |
122 |
_END_MASTER( myThid ) |
123 |
_BARRIER |
124 |
ENDIF |
125 |
#ifndef ALWAYS_USE_MPI |
126 |
ENDIF |
127 |
#endif |
128 |
#endif |
129 |
C-- Under a "put" scenario we |
130 |
C-- i. set completetion signal for buffer we put into. |
131 |
C-- ii. wait for completetion signal indicating data has been put in |
132 |
C-- our buffer. |
133 |
C-- Under a messaging mode we "receive" the message. |
134 |
C-- Under a "get" scenario we |
135 |
C-- i. Check that the data is ready. |
136 |
C-- ii. Read the data. |
137 |
C-- iii. Set data read flag + memory sync. |
138 |
|
139 |
|
140 |
DO bj=myByLo(myThid),myByHi(myThid) |
141 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
142 |
ebL = exchangeBufLevel(1,bi,bj) |
143 |
southCommMode = _tileCommModeS(bi,bj) |
144 |
northCommMode = _tileCommModeN(bi,bj) |
145 |
biN = _tileBiN(bi,bj) |
146 |
bjN = _tileBjN(bi,bj) |
147 |
biS = _tileBiS(bi,bj) |
148 |
bjS = _tileBjS(bi,bj) |
149 |
IF ( southCommMode .EQ. COMM_MSG ) THEN |
150 |
#ifdef ALLOW_USE_MPI |
151 |
#ifndef ALWAYS_USE_MPI |
152 |
IF ( usingMPI ) THEN |
153 |
#endif |
154 |
theProc = tilePidS(bi,bj) |
155 |
theTag = _tileTagRecvS(bi,bj) |
156 |
theType = _MPI_TYPE_RX |
157 |
theSize = sNx*exchWidthY*myNz |
158 |
IF ( theCornerMode .EQ. EXCH_UPDATE_CORNERS ) |
159 |
& theSize = (sNx+2*exchWidthX)*exchWidthY*myNz |
160 |
CALL MPI_Recv( southRecvBuf_RX(1,eBl,bi,bj), theSize, theType, |
161 |
& theProc, theTag, MPI_COMM_MODEL, |
162 |
& mpiStatus, mpiRc ) |
163 |
#ifndef ALWAYS_USE_MPI |
164 |
ENDIF |
165 |
#endif |
166 |
#endif /* ALLOW_USE_MPI */ |
167 |
ENDIF |
168 |
IF ( northCommMode .EQ. COMM_MSG ) THEN |
169 |
#ifdef ALLOW_USE_MPI |
170 |
#ifndef ALWAYS_USE_MPI |
171 |
IF ( usingMPI ) THEN |
172 |
#endif |
173 |
theProc = tilePidN(bi,bj) |
174 |
theTag = _tileTagRecvN(bi,bj) |
175 |
theType = _MPI_TYPE_RX |
176 |
theSize = sNx*exchWidthY*myNz |
177 |
IF ( theCornerMode .EQ. EXCH_UPDATE_CORNERS ) |
178 |
& theSize = (sNx+2*exchWidthX)*exchWidthY*myNz |
179 |
CALL MPI_Recv( northRecvBuf_RX(1,eBl,bi,bj), theSize, theType, |
180 |
& theProc, theTag, MPI_COMM_MODEL, |
181 |
& mpiStatus, mpiRc ) |
182 |
#ifndef ALWAYS_USE_MPI |
183 |
ENDIF |
184 |
#endif |
185 |
#endif /* ALLOW_USE_MPI */ |
186 |
ENDIF |
187 |
ENDDO |
188 |
ENDDO |
189 |
|
190 |
C-- Wait for buffers I am going read to be ready. |
191 |
IF ( exchUsesBarrier ) THEN |
192 |
C o On some machines ( T90 ) use system barrier rather than spinning. |
193 |
CALL BARRIER( myThid ) |
194 |
ELSE |
195 |
C o Spin waiting for completetion flag. This avoids a global-lock |
196 |
C i.e. we only lock waiting for data that we need. |
197 |
DO bj=myByLo(myThid),myByHi(myThid) |
198 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
199 |
ebL = exchangeBufLevel(1,bi,bj) |
200 |
southCommMode = _tileCommModeS(bi,bj) |
201 |
northCommMode = _tileCommModeN(bi,bj) |
202 |
spinCount = 0 |
203 |
10 CONTINUE |
204 |
CALL FOOL_THE_COMPILER( spinCount ) |
205 |
spinCount = spinCount+1 |
206 |
C IF ( myThid .EQ. 1 .AND. spinCount .GT. _EXCH_SPIN_LIMIT ) THEN |
207 |
C STOP ' S/R EXCH_RECV_GET_Y: spinCount .GT. _EXCH_SPIN_LIMIT' |
208 |
C ENDIF |
209 |
IF ( southRecvAck(eBl,bi,bj) .EQ. 0. ) GOTO 10 |
210 |
IF ( northRecvAck(eBl,bi,bj) .EQ. 0. ) GOTO 10 |
211 |
C Clear requests |
212 |
southRecvAck(eBl,bi,bj) = 0. |
213 |
northRecvAck(eBl,bi,bj) = 0. |
214 |
C Update statistics |
215 |
IF ( exchCollectStatistics ) THEN |
216 |
exchRecvYExchCount(1,bi,bj) = exchRecvYExchCount(1,bi,bj)+1 |
217 |
exchRecvYSpinCount(1,bi,bj) = |
218 |
& exchRecvYSpinCount(1,bi,bj)+spinCount |
219 |
exchRecvYSpinMax(1,bi,bj) = |
220 |
& MAX(exchRecvYSpinMax(1,bi,bj),spinCount) |
221 |
exchRecvYSpinMin(1,bi,bj) = |
222 |
& MIN(exchRecvYSpinMin(1,bi,bj),spinCount) |
223 |
ENDIF |
224 |
|
225 |
|
226 |
IF ( exchNReqsY(1,bi,bj) .GT. 0 ) THEN |
227 |
#ifdef ALLOW_USE_MPI |
228 |
#ifndef ALWAYS_USE_MPI |
229 |
IF ( usingMPI ) THEN |
230 |
#endif |
231 |
CALL MPI_Waitall( exchNReqsY(1,bi,bj), exchReqIdY(1,1,bi,bj), |
232 |
& mpiStatus, mpiRC ) |
233 |
#ifndef ALWAYS_USE_MPI |
234 |
ENDIF |
235 |
#endif |
236 |
#endif /* ALLOW_USE_MPI */ |
237 |
ENDIF |
238 |
C Clear outstanding requests counter |
239 |
exchNReqsY(1,bi,bj) = 0 |
240 |
ENDDO |
241 |
ENDDO |
242 |
ENDIF |
243 |
|
244 |
C-- Read from the buffers |
245 |
DO bj=myByLo(myThid),myByHi(myThid) |
246 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
247 |
|
248 |
ebL = exchangeBufLevel(1,bi,bj) |
249 |
biN = _tileBiN(bi,bj) |
250 |
bjN = _tileBjN(bi,bj) |
251 |
biS = _tileBiS(bi,bj) |
252 |
bjS = _tileBjS(bi,bj) |
253 |
southCommMode = _tileCommModeS(bi,bj) |
254 |
northCommMode = _tileCommModeN(bi,bj) |
255 |
IF ( theCornerMode .EQ. EXCH_UPDATE_CORNERS ) THEN |
256 |
iMin = 1-exchWidthX |
257 |
iMax = sNx+exchWidthX |
258 |
ELSE |
259 |
iMin = 1 |
260 |
iMax = sNx |
261 |
ENDIF |
262 |
IF ( _theSimulationMode .EQ. FORWARD_SIMULATION ) THEN |
263 |
jMin = sNy+1 |
264 |
jMax = sNy+exchWidthY |
265 |
iB0 = 0 |
266 |
IF ( northCommMode .EQ. COMM_PUT |
267 |
& .OR. northCommMode .EQ. COMM_MSG ) THEN |
268 |
iB = 0 |
269 |
DO K=1,myNz |
270 |
DO J=jMin,jMax |
271 |
DO I=iMin,iMax |
272 |
iB = iB + 1 |
273 |
array(I,J,K,bi,bj) = northRecvBuf_RX(iB,eBl,bi,bj) |
274 |
ENDDO |
275 |
ENDDO |
276 |
ENDDO |
277 |
ELSEIF ( northCommMode .EQ. COMM_GET ) THEN |
278 |
DO K=1,myNz |
279 |
iB = iB0 |
280 |
DO J=jMin,jMax |
281 |
iB = iB+1 |
282 |
DO I=iMin,iMax |
283 |
array(I,J,K,bi,bj) = array(I,iB,K,biN,bjN) |
284 |
ENDDO |
285 |
ENDDO |
286 |
ENDDO |
287 |
ENDIF |
288 |
ELSEIF ( _theSimulationMode .EQ. REVERSE_SIMULATION ) THEN |
289 |
jMin = sNy-exchWidthY+1 |
290 |
jMax = sNy |
291 |
iB0 = 1-exchWidthY-1 |
292 |
IF ( northCommMode .EQ. COMM_PUT |
293 |
& .OR. northCommMode .EQ. COMM_MSG ) THEN |
294 |
iB = 0 |
295 |
DO K=1,myNz |
296 |
DO J=jMin,jMax |
297 |
DO I=iMin,iMax |
298 |
iB = iB + 1 |
299 |
array(I,J,K,bi,bj) = |
300 |
& array(I,J,K,bi,bj)+northRecvBuf_RX(iB,eBl,bi,bj) |
301 |
ENDDO |
302 |
ENDDO |
303 |
ENDDO |
304 |
ELSEIF ( northCommMode .EQ. COMM_GET ) THEN |
305 |
DO K=1,myNz |
306 |
iB = iB0 |
307 |
DO J=jMin,jMax |
308 |
iB = iB+1 |
309 |
DO I=iMin,iMax |
310 |
array(I,J,K,bi,bj) = |
311 |
& array(I,J,K,bi,bj)+array(I,iB,K,biN,bjN) |
312 |
ENDDO |
313 |
ENDDO |
314 |
ENDDO |
315 |
ENDIF |
316 |
ENDIF |
317 |
|
318 |
IF ( _theSimulationMode .EQ. FORWARD_SIMULATION ) THEN |
319 |
jMin = 1-exchWidthY |
320 |
jMax = 0 |
321 |
iB0 = sNy-exchWidthY |
322 |
IF ( southCommMode .EQ. COMM_PUT |
323 |
& .OR. southCommMode .EQ. COMM_MSG ) THEN |
324 |
iB = 0 |
325 |
DO K=1,myNz |
326 |
DO J=jMin,jMax |
327 |
DO I=iMin,iMax |
328 |
iB = iB + 1 |
329 |
array(I,J,K,bi,bj) = southRecvBuf_RX(iB,eBl,bi,bj) |
330 |
ENDDO |
331 |
ENDDO |
332 |
ENDDO |
333 |
ELSEIF ( southCommMode .EQ. COMM_GET ) THEN |
334 |
DO K=1,myNz |
335 |
iB = iB0 |
336 |
DO J=jMin,jMax |
337 |
iB = iB+1 |
338 |
DO I=iMin,iMax |
339 |
array(I,J,K,bi,bj) = array(I,iB,K,biS,bjS) |
340 |
ENDDO |
341 |
ENDDO |
342 |
ENDDO |
343 |
ENDIF |
344 |
ELSEIF ( _theSimulationMode .EQ. REVERSE_SIMULATION ) THEN |
345 |
jMin = 1 |
346 |
jMax = 1+exchWidthY-1 |
347 |
iB0 = sNy |
348 |
IF ( southCommMode .EQ. COMM_PUT |
349 |
& .OR. southCommMode .EQ. COMM_MSG ) THEN |
350 |
iB = 0 |
351 |
DO K=1,myNz |
352 |
DO J=jMin,jMax |
353 |
DO I=iMin,iMax |
354 |
iB = iB + 1 |
355 |
array(I,J,K,bi,bj) = |
356 |
& array(I,J,K,bi,bj)+southRecvBuf_RX(iB,eBl,bi,bj) |
357 |
ENDDO |
358 |
ENDDO |
359 |
ENDDO |
360 |
ELSEIF ( southCommMode .EQ. COMM_GET ) THEN |
361 |
DO K=1,myNz |
362 |
iB = iB0 |
363 |
DO J=jMin,jMax |
364 |
iB = iB+1 |
365 |
DO I=iMin,iMax |
366 |
array(I,J,K,bi,bj) = |
367 |
& array(I,J,K,bi,bj)+array(I,iB,K,biS,bjS) |
368 |
ENDDO |
369 |
ENDDO |
370 |
ENDDO |
371 |
ENDIF |
372 |
ENDIF |
373 |
ENDDO |
374 |
ENDDO |
375 |
|
376 |
_BARRIER |
377 |
IF ( doingSingleThreadedComms ) THEN |
378 |
C Restore saved settings that were stored to allow |
379 |
C single thred comms. |
380 |
_BEGIN_MASTER(myThid) |
381 |
DO I=1,nThreads |
382 |
myBxLo(I) = myBxLoSave(I) |
383 |
myBxHi(I) = myBxHiSave(I) |
384 |
myByLo(I) = myByLoSave(I) |
385 |
myByHi(I) = myByHiSave(I) |
386 |
ENDDO |
387 |
_END_MASTER(myThid) |
388 |
ENDIF |
389 |
_BARRIER |
390 |
|
391 |
RETURN |
392 |
END |