/[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.9 - (hide annotations) (download)
Tue Mar 18 21:34:01 2008 UTC (17 years, 4 months ago) by utke
Branch: MAIN
Changes since 1.8: +42 -1 lines
aMPI prototype

1 utke 1.9 C $Header: /u/gcmpack/MITgcm/eesupp/src/exch_rx_recv_get_x.template,v 1.8 2008/02/20 20:19:00 jmc 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     INTEGER theProc, theTag, theType, theSize
80     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     CALL ampi_recv_RX(
163     & westRecvBuf_RX(1,eBl,bi,bj) ,
164     & theSize ,
165     & theType ,
166     & theProc ,
167     & theTag ,
168     & MPI_COMM_MODEL ,
169     & mpiStatus ,
170     & mpiRc )
171     # endif /* ALLOW_AUTODIFF_OPENAD */
172 adcroft 1.1 #ifndef ALWAYS_USE_MPI
173     ENDIF
174     #endif
175     #endif /* ALLOW_USE_MPI */
176     ENDIF
177     IF ( eastCommMode .EQ. COMM_MSG ) THEN
178     #ifdef ALLOW_USE_MPI
179     #ifndef ALWAYS_USE_MPI
180     IF ( usingMPI ) THEN
181     #endif
182     theProc = tilePidE(bi,bj)
183     theTag = _tileTagRecvE(bi,bj)
184 dimitri 1.3 theType = _MPI_TYPE_RX
185 adcroft 1.1 theSize = sNy*exchWidthX*myNz
186 utke 1.9 # ifndef ALLOW_AUTODIFF_OPENAD
187 adcroft 1.1 CALL MPI_Recv( eastRecvBuf_RX(1,eBl,bi,bj), theSize, theType,
188     & theProc, theTag, MPI_COMM_MODEL,
189     & mpiStatus, mpiRc )
190 utke 1.9 # else
191     CALL ampi_recv_RX(
192     & eastRecvBuf_RX(1,eBl,bi,bj) ,
193     & theSize ,
194     & theType ,
195     & theProc ,
196     & theTag ,
197     & MPI_COMM_MODEL ,
198     & mpiStatus ,
199     & mpiRc )
200     # endif /* ALLOW_AUTODIFF_OPENAD */
201 adcroft 1.1 #ifndef ALWAYS_USE_MPI
202     ENDIF
203     #endif
204     #endif /* ALLOW_USE_MPI */
205     ENDIF
206     ENDDO
207     ENDDO
208    
209     C-- Wait for buffers I am going read to be ready.
210     IF ( exchUsesBarrier ) THEN
211     C o On some machines ( T90 ) use system barrier rather than spinning.
212     CALL BARRIER( myThid )
213     ELSE
214     C o Spin waiting for completetion flag. This avoids a global-lock
215     C i.e. we only lock waiting for data that we need.
216     DO bj=myByLo(myThid),myByHi(myThid)
217     DO bi=myBxLo(myThid),myBxHi(myThid)
218     spinCount = 0
219     ebL = exchangeBufLevel(1,bi,bj)
220     westCommMode = _tileCommModeW(bi,bj)
221     eastCommMode = _tileCommModeE(bi,bj)
222 utke 1.9 # ifndef ALLOW_AUTODIFF_OPENAD
223 adcroft 1.1 10 CONTINUE
224 jmc 1.4 CALL FOOL_THE_COMPILER( spinCount )
225 adcroft 1.1 spinCount = spinCount+1
226     C IF ( myThid .EQ. 1 .AND. spinCount .GT. _EXCH_SPIN_LIMIT ) THEN
227     C WRITE(*,*) ' eBl = ', ebl
228     C STOP ' S/R EXCH_RECV_GET_X: spinCount .GT. _EXCH_SPIN_LIMIT'
229     C ENDIF
230 jmc 1.8 IF ( westRecvAck(eBl,bi,bj) .EQ. 0 ) GOTO 10
231     IF ( eastRecvAck(eBl,bi,bj) .EQ. 0 ) GOTO 10
232 utke 1.9 # else
233     do while ((westRecvAck(eBl,bi,bj) .EQ. 0.
234     & .or.
235     & eastRecvAck(eBl,bi,bj) .EQ. 0. ))
236     CALL FOOL_THE_COMPILER( spinCount )
237     spinCount = spinCount+1
238     end do
239     # endif /* ALLOW_AUTODIFF_OPENAD */
240 adcroft 1.1 C Clear outstanding requests
241 jmc 1.8 westRecvAck(eBl,bi,bj) = 0
242     eastRecvAck(eBl,bi,bj) = 0
243 adcroft 1.1
244     IF ( exchNReqsX(1,bi,bj) .GT. 0 ) THEN
245     #ifdef ALLOW_USE_MPI
246     #ifndef ALWAYS_USE_MPI
247     IF ( usingMPI ) THEN
248     #endif
249 utke 1.9 # ifndef ALLOW_AUTODIFF_OPENAD
250 adcroft 1.1 CALL MPI_Waitall( exchNReqsX(1,bi,bj), exchReqIdX(1,1,bi,bj),
251     & mpiStatus, mpiRC )
252 utke 1.9 # else
253     CALL ampi_waitall(
254     & exchNReqsX(1,bi,bj),
255     & exchReqIdX(1,1,bi,bj),
256     & mpiStatus,
257     & mpiRC )
258     # endif /* ALLOW_AUTODIFF_OPENAD */
259 adcroft 1.1 #ifndef ALWAYS_USE_MPI
260     ENDIF
261     #endif
262     #endif /* ALLOW_USE_MPI */
263     ENDIF
264     C Clear outstanding requests counter
265     exchNReqsX(1,bi,bj) = 0
266     C Update statistics
267     IF ( exchCollectStatistics ) THEN
268     exchRecvXExchCount(1,bi,bj) = exchRecvXExchCount(1,bi,bj)+1
269     exchRecvXSpinCount(1,bi,bj) =
270     & exchRecvXSpinCount(1,bi,bj)+spinCount
271     exchRecvXSpinMax(1,bi,bj) =
272     & MAX(exchRecvXSpinMax(1,bi,bj),spinCount)
273     exchRecvXSpinMin(1,bi,bj) =
274     & MIN(exchRecvXSpinMin(1,bi,bj),spinCount)
275     ENDIF
276     ENDDO
277     ENDDO
278     ENDIF
279    
280     C-- Read from the buffers
281     DO bj=myByLo(myThid),myByHi(myThid)
282     DO bi=myBxLo(myThid),myBxHi(myThid)
283    
284     ebL = exchangeBufLevel(1,bi,bj)
285     biE = _tileBiE(bi,bj)
286     bjE = _tileBjE(bi,bj)
287     biW = _tileBiW(bi,bj)
288     bjW = _tileBjW(bi,bj)
289     westCommMode = _tileCommModeW(bi,bj)
290     eastCommMode = _tileCommModeE(bi,bj)
291     IF ( _theSimulationMode .EQ. FORWARD_SIMULATION ) THEN
292     iMin = sNx+1
293     iMax = sNx+exchWidthX
294     iB0 = 0
295     IF ( eastCommMode .EQ. COMM_PUT
296     & .OR. eastCommMode .EQ. COMM_MSG ) THEN
297     iB = 0
298     DO K=1,myNz
299     DO J=1,sNy
300     DO I=iMin,iMax
301     iB = iB + 1
302     array(I,J,K,bi,bj) = eastRecvBuf_RX(iB,eBl,bi,bj)
303     ENDDO
304     ENDDO
305     ENDDO
306     ELSEIF ( eastCommMode .EQ. COMM_GET ) THEN
307     DO K=1,myNz
308     DO J=1,sNy
309     iB = iB0
310     DO I=iMin,iMax
311     iB = iB+1
312     array(I,J,K,bi,bj) = array(iB,J,K,biE,bjE)
313     ENDDO
314     ENDDO
315     ENDDO
316     ENDIF
317     ELSEIF ( _theSimulationMode .EQ. REVERSE_SIMULATION ) THEN
318     iMin = sNx-exchWidthX+1
319     iMax = sNx
320     iB0 = 1-exchWidthX-1
321     IF ( eastCommMode .EQ. COMM_PUT
322     & .OR. eastCommMode .EQ. COMM_MSG ) THEN
323     iB = 0
324     DO K=1,myNz
325     DO J=1,sNy
326     DO I=iMin,iMax
327     iB = iB + 1
328     array(I,J,K,bi,bj) =
329     & array(I,J,K,bi,bj)+eastRecvBuf_RX(iB,eBl,bi,bj)
330     ENDDO
331     ENDDO
332     ENDDO
333     ELSEIF ( eastCommMode .EQ. COMM_GET ) THEN
334     DO K=1,myNz
335     DO J=1,sNy
336     iB = iB0
337     DO I=iMin,iMax
338     iB = iB+1
339     array(I,J,K,bi,bj) =
340     & array(I,J,K,bi,bj)+array(iB,J,K,biE,bjE)
341     ENDDO
342     ENDDO
343     ENDDO
344     ENDIF
345     ENDIF
346     IF ( _theSimulationMode .EQ. FORWARD_SIMULATION ) THEN
347     iMin = 1-exchWidthX
348     iMax = 0
349     iB0 = sNx-exchWidthX
350     IF ( westCommMode .EQ. COMM_PUT
351     & .OR. westCommMode .EQ. COMM_MSG ) THEN
352     iB = 0
353     DO K=1,myNz
354     DO J=1,sNy
355     DO I=iMin,iMax
356     iB = iB + 1
357     array(I,J,K,bi,bj) = westRecvBuf_RX(iB,eBl,bi,bj)
358     ENDDO
359     ENDDO
360     ENDDO
361     ELSEIF ( westCommMode .EQ. COMM_GET ) THEN
362     DO K=1,myNz
363     DO J=1,sNy
364     iB = iB0
365     DO I=iMin,iMax
366     iB = iB+1
367     array(I,J,K,bi,bj) = array(iB,J,K,biW,bjW)
368     ENDDO
369     ENDDO
370     ENDDO
371     ENDIF
372     ELSEIF ( _theSimulationMode .EQ. REVERSE_SIMULATION ) THEN
373     iMin = 1
374     iMax = 1+exchWidthX-1
375     iB0 = sNx
376     IF ( westCommMode .EQ. COMM_PUT
377     & .OR. westCommMode .EQ. COMM_MSG ) THEN
378     iB = 0
379     DO K=1,myNz
380     DO J=1,sNy
381     DO I=iMin,iMax
382     iB = iB + 1
383     array(I,J,K,bi,bj) =
384     & array(I,J,K,bi,bj)+westRecvBuf_RX(iB,eBl,bi,bj)
385     ENDDO
386     ENDDO
387     ENDDO
388     ELSEIF ( westCommMode .EQ. COMM_GET ) THEN
389     DO K=1,myNz
390     DO J=1,sNy
391     iB = iB0
392     DO I=iMin,iMax
393     iB = iB+1
394     array(I,J,K,bi,bj) =
395     & array(I,J,K,bi,bj)+array(iB,J,K,biW,bjW)
396     ENDDO
397     ENDDO
398     ENDDO
399     ENDIF
400     ENDIF
401    
402     ENDDO
403     ENDDO
404    
405 cnh 1.7 _BARRIER
406 cnh 1.6 IF ( doingSingleThreadedComms ) THEN
407     C Restore saved settings that were stored to allow
408     C single thred comms.
409     _BEGIN_MASTER(myThid)
410     DO I=1,nThreads
411     myBxLo(I) = myBxLoSave(I)
412     myBxHi(I) = myBxHiSave(I)
413     myByLo(I) = myByLoSave(I)
414     myByHi(I) = myByHiSave(I)
415     ENDDO
416     _END_MASTER(myThid)
417     ENDIF
418 cnh 1.7 _BARRIER
419 cnh 1.5
420 adcroft 1.1 RETURN
421     END

  ViewVC Help
Powered by ViewVC 1.1.22