/[MITgcm]/MITgcm/eesupp/src/exch_rx_recv_get_x.template
ViewVC logotype

Annotation of /MITgcm/eesupp/src/exch_rx_recv_get_x.template

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.11 - (hide annotations) (download)
Fri Apr 4 20:18:34 2008 UTC (16 years, 2 months ago) by utke
Branch: MAIN
CVS Tags: checkpoint60, checkpoint61, checkpoint59q, checkpoint59p, checkpoint59r
Changes since 1.10: +8 -6 lines
fix requestid passing for ampi calls

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

  ViewVC Help
Powered by ViewVC 1.1.22