source: CPL/oasis3-mct/branches/OASIS3-MCT_2.0_branch/lib/mct/mct/m_Transfer.F90 @ 4775

Last change on this file since 4775 was 4775, checked in by aclsce, 5 years ago
  • Imported oasis3-mct from Cerfacs svn server (not suppotred anymore).

The version has been extracted from https://oasis3mct.cerfacs.fr/svn/branches/OASIS3-MCT_2.0_branch/oasis3-mct@1818

File size: 24.0 KB
Line 
1!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2!    Math and Computer Science Division, Argonne National Laboratory   !
3!-----------------------------------------------------------------------
4! CVS m_Transfer.F90,v 1.13 2008-01-28 06:48:03 jacob Exp
5! CVS MCT_2_8_0
6!BOP -------------------------------------------------------------------
7!
8! !MODULE: m_Transfer - Routines for the MxN transfer of Attribute Vectors
9!
10! !DESCRIPTION:
11! This module provides routines for doing MxN transfer of data in an
12! Attribute Vector between two components on separate sets of MPI processes.
13! Uses the Router datatype.
14!
15! !SEE ALSO:
16! m_Rearranger
17
18! !INTERFACE:
19
20 module m_Transfer
21
22! !USES:
23  use m_MCTWorld, only : MCTWorld
24  use m_MCTWorld, only : ThisMCTWorld
25  use m_AttrVect, only : AttrVect
26  use m_AttrVect, only : nIAttr,nRAttr
27  use m_AttrVect, only : Permute, Unpermute
28  use m_AttrVect, only : AttrVect_init => init
29  use m_AttrVect, only : AttrVect_copy => copy 
30  use m_AttrVect, only : AttrVect_clean => clean
31  use m_AttrVect, only : lsize
32  use m_Router,   only : Router
33
34  use m_mpif90
35  use m_die
36  use m_stdio
37
38  implicit none
39
40  private ! except
41
42! !PUBLIC MEMBER FUNCTIONS:
43
44  public  :: isend
45  public  :: send
46  public  :: waitsend
47  public  :: irecv
48  public  :: recv
49  public  :: waitrecv
50
51
52  interface isend    ; module procedure isend_    ; end interface
53  interface send     ; module procedure send_     ; end interface
54  interface waitsend ; module procedure waitsend_ ; end interface
55  interface irecv    ; module procedure irecv_    ; end interface
56  interface recv     ; module procedure recv_     ; end interface
57  interface waitrecv ; module procedure waitrecv_ ; end interface
58
59! !DEFINED PARAMETERS:
60
61  integer,parameter                    :: DefaultTag = 600
62
63! !REVISION HISTORY:
64! 08Nov02 - R. Jacob <jacob@mcs.anl.gov> - make new module by combining
65!           MCT_Send, MCT_Recv and MCT_Recvsum
66! 11Nov02 - R. Jacob <jacob@mcs.anl.gov> - Remove MCT_Recvsum and use
67!           optional argument in recv_ to do the same thing.
68! 23Jul03 - R. Jacob <jacob@mcs.anl.gov> - Move buffers for data and
69!           MPI_Reqest and MPI_Status arrays to Router.  Use them.
70! 24Jul03 - R. Jacob <jacob@mcs.anl.gov> - Split send_ into isend_ and
71!           waitsend_.  Redefine send_.
72! 22Jan08 - R. Jacob <jacob@mcs.anl.gov> - Handle unordered GSMaps
73!EOP ___________________________________________________________________
74
75  character(len=*),parameter :: myname='MCT::m_Transfer'
76
77  contains
78
79!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
80!    Math and Computer Science Division, Argonne National Laboratory   !
81!BOP -------------------------------------------------------------------
82!
83! !IROUTINE: isend_ - Distributed non-blocking send of an Attribute Vector
84!
85! !DESCRIPTION:
86! Send the the data in the {\tt AttrVect} {\tt aV} to the
87! component specified in the {\tt Router} {\tt Rout}.  An error will
88! result if the size of the attribute vector does not match the size
89! parameter stored in the {\tt Router}.
90!
91! Requires a corresponding {\tt recv\_} or {\tt irecv\_} to be called on the other component.
92!
93! The optional argument {\tt Tag} can be used to set the tag value used in
94! the data transfer.  DefaultTag will be used otherwise. {\tt Tag} must be
95! the same in the matching {\tt recv\_} or {\tt irecv\_}.
96!
97! {\bf N.B.:} The {\tt AttrVect} argument in the corresponding
98! {\tt recv\_} call is assumed to have exactly the same attributes
99! in exactly the same order as {\tt aV}.
100!
101! !INTERFACE:
102
103 subroutine isend_(aVin, Rout, Tag)
104
105!
106! !USES:
107!
108      implicit none
109
110! !INPUT PARAMETERS:
111!
112
113      Type(AttrVect),target,intent(in)    :: aVin
114      Type(Router),         intent(inout) :: Rout
115      integer,optional,     intent(in) :: Tag
116
117! !REVISION HISTORY:
118! 07Feb01 - R. Jacob <jacob@mcs.anl.gov> - initial prototype
119! 08Feb01 - R. Jacob <jacob@mcs.anl.gov> - First working code
120! 18May01 - R. Jacob <jacob@mcs.anl.gov> - use MP_Type to determine type in mpi_send
121! 07Jun01 - R. Jacob <jacob@mcs.anl.gov> - remove logic to check "direction" of Router. 
122!           remove references to ThisMCTWorld%mylrank
123! 03Aug01 - E. Ong <eong@mcs.anl.gov> - Explicitly specify the starting address in mpi_send. 
124! 15Feb02 - R. Jacob <jacob@mcs.anl.gov> - Use MCT_comm
125! 26Mar02 - E. Ong <eong@mcs.anl.gov> - Apply faster copy order
126! 26Sep02 - R. Jacob <jacob@mcs.anl.gov> - Check Av against Router lAvsize
127! 05Nov02 - R. Jacob <jacob@mcs.anl.gov> - Remove iList, rList arguments.
128! 08Nov02 - R. Jacob <jacob@mcs.anl.gov> - MCT_Send is now send_ in m_Transfer
129! 11Nov02 - R. Jacob <jacob@mcs.anl.gov> - Use DefaultTag and add optional Tag argument
130! 25Jul03 - R. Jacob <jacob@mcs.anl.gov> - Split into isend_ and waitsend_
131! 22Jan08 - R. Jacob <jacob@mcs.anl.gov> - Handle unordered GSMaps by permuting before send.
132!           remove special case for sending one segment directly from Av which probably
133!           wasn't safe.
134!EOP ___________________________________________________________________
135
136  character(len=*),parameter :: myname_=myname//'::isend_'
137  integer ::    numi,numr,i,j,k,ier
138  integer ::    mycomp,othercomp
139  integer ::    AttrIndex,VectIndex,seg_start,seg_end
140  integer ::    proc,nseg,mytag
141  integer ::    mp_Type_rp1
142  logical ::    unordered
143  type(AttrVect),pointer :: Av
144  type(AttrVect),target  :: Avtmp
145
146!--------------------------------------------------------
147
148! Return if no one to send to
149  if(Rout%nprocs .eq. 0 ) RETURN
150
151! set up Av to send from
152  unordered = associated(Rout%permarr)
153  if (unordered) then
154     call AttrVect_init(Avtmp,Avin,lsize(Avin))
155     call AttrVect_copy(Avin,aVtmp)
156     call Permute(aVtmp,Rout%permarr)
157     Av => Avtmp
158  else
159     Av => Avin
160  endif
161
162!check Av size against Router
163!
164  if(lsize(aV) /= Rout%lAvsize) then
165    write(stderr,'(2a)') myname_, &
166    ' MCTERROR:  AV size not appropriate for this Router...exiting'
167    call die(myname_)
168  endif
169
170! get ids of components involved in this communication
171  mycomp=Rout%comp1id
172  othercomp=Rout%comp2id
173
174
175! find total number of real and integer vectors
176! for now, assume we are sending all of them
177  Rout%numiatt = nIAttr(aV)
178  Rout%numratt = nRAttr(aV)
179  numi = Rout%numiatt
180  numr = Rout%numratt
181
182!!!!!!!!!!!!!! IF SENDING INTEGER DATA
183  if(numi .ge. 1) then
184
185! allocate buffers to hold all outgoing data
186  do proc=1,Rout%nprocs
187     allocate(Rout%ip1(proc)%pi(Rout%locsize(proc)*numi),stat=ier)
188     if(ier/=0) call die(myname_,'allocate(Rout%ip1%pi)',ier)
189  enddo
190
191  endif
192
193!!!!!!!!!!!!!! IF SENDING REAL DATA
194  if(numr .ge. 1) then
195
196! allocate buffers to hold all outgoing data
197  do proc=1,Rout%nprocs
198     allocate(Rout%rp1(proc)%pr(Rout%locsize(proc)*numr),stat=ier)
199     if(ier/=0) call die(myname_,'allocate(Rout%rp1%pr)',ier)
200  enddo
201
202  mp_Type_rp1=MP_Type(Rout%rp1(1)%pr(1))
203
204  endif
205
206
207  ! Load data going to each processor
208  do proc = 1,Rout%nprocs
209   
210     j=1
211     k=1
212
213     ! load the correct pieces of the integer and real vectors
214     ! if Rout%num_segs(proc)=1, then this will do one loop
215     do nseg = 1,Rout%num_segs(proc)
216         seg_start = Rout%seg_starts(proc,nseg)
217         seg_end = seg_start + Rout%seg_lengths(proc,nseg)-1
218         do VectIndex = seg_start,seg_end
219            do AttrIndex = 1,numi
220                Rout%ip1(proc)%pi(j) = aV%iAttr(AttrIndex,VectIndex)
221                j=j+1
222            enddo
223            do AttrIndex = 1,numr
224                Rout%rp1(proc)%pr(k) = aV%rAttr(AttrIndex,VectIndex)
225                k=k+1
226            enddo
227         enddo
228     enddo
229
230
231
232     ! Send the integer data
233     if(numi .ge. 1) then
234
235        ! set tag
236        mytag = DefaultTag
237        if(present(Tag)) mytag=Tag
238
239
240        call MPI_ISEND(Rout%ip1(proc)%pi(1), &
241             Rout%locsize(proc)*numi,MP_INTEGER,Rout%pe_list(proc), &
242             mytag,ThisMCTWorld%MCT_comm,Rout%ireqs(proc),ier)
243           
244        if(ier /= 0) call MP_perr_die(myname_,'MPI_ISEND(ints)',ier)
245
246     endif
247
248     ! Send the real data
249     if(numr .ge. 1) then
250
251       ! set tag
252       mytag = DefaultTag + 1 
253       if(present(Tag)) mytag=Tag +1
254
255
256       call MPI_ISEND(Rout%rp1(proc)%pr(1), &
257            Rout%locsize(proc)*numr,mp_Type_rp1,Rout%pe_list(proc), &
258            mytag,ThisMCTWorld%MCT_comm,Rout%rreqs(proc),ier)
259
260
261       if(ier /= 0) call MP_perr_die(myname_,'MPI_ISEND(reals)',ier)
262
263    endif
264
265  enddo
266
267  if (unordered) then
268     call AttrVect_clean(aVtmp)
269     nullify(aV)
270  else
271     nullify(aV)
272  endif
273
274end subroutine isend_
275
276!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
277!    Math and Computer Science Division, Argonne National Laboratory   !
278!BOP -------------------------------------------------------------------
279!
280! !IROUTINE: waitsend_ - Wait for a distributed non-blocking send to complete
281!
282! !DESCRIPTION:
283! Wait for the data being sent with the {\tt Router} {\tt Rout} to complete.
284!
285! !INTERFACE:
286
287 subroutine waitsend_(Rout)
288
289!
290! !USES:
291!
292      implicit none
293
294! !INPUT PARAMETERS:
295!
296      Type(Router),         intent(inout) :: Rout
297
298! !REVISION HISTORY:
299! 24Jul03 - R. Jacob <jacob@mcs.anl.gov> - First working version is
300!           the wait part of original send_
301!EOP ___________________________________________________________________
302  character(len=*),parameter :: myname_=myname//'::waitsend_'
303  integer ::   proc,ier
304
305! Return if nothing to wait for
306  if(Rout%nprocs .eq. 0 ) RETURN
307
308  ! wait for all sends to complete
309  if(Rout%numiatt .ge. 1) then
310
311     call MPI_WAITALL(Rout%nprocs,Rout%ireqs,Rout%istatus,ier)
312     if(ier /= 0) call MP_perr_die(myname_,'MPI_WAITALL(ints)',ier)
313
314     do proc=1,Rout%nprocs
315        deallocate(Rout%ip1(proc)%pi,stat=ier)
316        if(ier/=0) call die(myname_,'deallocate(ip1%pi)',ier)
317     enddo
318
319  endif
320
321  if(Rout%numratt .ge. 1) then
322
323     call MPI_WAITALL(Rout%nprocs,Rout%rreqs,Rout%rstatus,ier)
324     if(ier /= 0) call MP_perr_die(myname_,'MPI_WAITALL(reals)',ier)
325
326     do proc=1,Rout%nprocs
327        deallocate(Rout%rp1(proc)%pr,stat=ier)
328        if(ier/=0) call die(myname_,'deallocate(rp1%pi)',ier)
329     enddo
330
331  endif
332
333
334end subroutine waitsend_
335
336!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
337!    Math and Computer Science Division, Argonne National Laboratory   !
338!BOP -------------------------------------------------------------------
339!
340! !IROUTINE: send_ - Distributed blocking send of an Attribute Vector
341!
342! !DESCRIPTION:
343! Send the the data in the {\tt AttrVect} {\tt aV} to the
344! component specified in the {\tt Router} {\tt Rout}.  An error will
345! result if the size of the attribute vector does not match the size
346! parameter stored in the {\tt Router}.
347!
348! Requires a corresponding {\tt recv\_} or {\tt irecv\_} to be called on the other
349! component.
350!
351! The optional argument {\tt Tag} can be used to set the tag value used in
352! the data transfer.  DefaultTag will be used otherwise. {\tt Tag} must be
353! the same in the matching {\tt recv\_} or {\tt irecv\_}.
354!
355! {\bf N.B.:} The {\tt AttrVect} argument in the corresponding
356! {\tt recv} call is assumed to have exactly the same attributes
357! in exactly the same order as {\tt aV}.
358!
359! !INTERFACE:
360
361 subroutine send_(aV, Rout, Tag)
362
363!
364! !USES:
365!
366      implicit none
367
368! !INPUT PARAMETERS:
369!
370
371      Type(AttrVect),       intent(in)    :: aV 
372      Type(Router),         intent(inout) :: Rout
373      integer,optional,     intent(in)    :: Tag
374
375! !REVISION HISTORY:
376! 24Jul03 - R. Jacob <jacob@mcs.anl.gov> - New version uses isend and waitsend
377!EOP ___________________________________________________________________
378  character(len=*),parameter :: myname_=myname//'::send_'
379
380  call isend_(aV,Rout,Tag)
381
382  call waitsend_(Rout)
383
384end subroutine send_
385
386!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
387!    Math and Computer Science Division, Argonne National Laboratory   !
388!BOP -------------------------------------------------------------------
389!
390! !IROUTINE: irecv_ - Distributed receive of an Attribute Vector
391!
392! !DESCRIPTION:
393! Recieve into the {\tt AttrVect} {\tt aV} the data coming from the
394! component specified in the {\tt Router} {\tt Rout}.  An error will
395! result if the size of the attribute vector does not match the size
396! parameter stored in the {\tt Router}.
397!
398! Requires a corresponding {\tt send\_} or {\tt isend\_} to be called
399! on the other component.
400!
401! The optional argument {\tt Tag} can be used to set the tag value used in
402! the data transfer.  DefaultTag will be used otherwise. {\tt Tag} must be
403! the same in the matching {\tt send\_} or {\tt isend\_}.
404!
405! If data for a grid point is coming from more than one process, {\tt recv\_}
406! will overwrite the duplicate values leaving the last received value
407! in the output aV.  If the optional argument {\tt Sum} is invoked, the output
408! will contain the sum of any duplicate values received for the same grid point.
409!
410! Will return as soon as MPI\_IRECV's are posted.  Call {\tt waitrecv\_} to
411! complete the receive operation.
412!
413! {\bf N.B.:} The {\tt AttrVect} argument in the corresponding
414! {\tt send\_} call is assumed to have exactly the same attributes
415! in exactly the same order as {\tt aV}.
416!
417! !INTERFACE:
418
419 subroutine irecv_(aV, Rout, Tag, Sum)
420!
421! !USES:
422!
423      implicit none
424
425! !INPUT/OUTPUT PARAMETERS:
426!
427      Type(AttrVect),       intent(inout) :: aV
428
429! !INPUT PARAMETERS:
430!
431      Type(Router),         intent(inout) :: Rout
432      integer,optional,     intent(in)    :: Tag
433      logical,optional,     intent(in)    :: Sum
434
435! !REVISION HISTORY:
436! 07Feb01 - R. Jacob <jacob@mcs.anl.gov> - initial prototype
437! 07Jun01 - R. Jacob <jacob@mcs.anl.gov> - remove logic to
438!           check "direction" of Router.  remove references
439!           to ThisMCTWorld%mylrank
440! 03Aug01 - E.T. Ong <eong@mcs.anl.gov> - explicity specify starting
441!           address in MPI_RECV
442! 27Nov01 - E.T. Ong <eong@mcs.anl.gov> - deallocated to prevent
443!           memory leaks
444! 15Feb02 - R. Jacob <jacob@mcs.anl.gov> - Use MCT_comm
445! 26Mar02 - E. Ong <eong@mcs.anl.gov> - Apply faster copy order.
446! 26Sep02 - R. Jacob <jacob@mcs.anl.gov> - Check Av against Router lAvsize
447! 08Nov02 - R. Jacob <jacob@mcs.anl.gov> - MCT_Recv is now recv_ in m_Transfer
448! 11Nov02 - R. Jacob <jacob@mcs.anl.gov> - Add optional Sum argument to
449!           tell recv_ to sum data for the same point received from multiple
450!           processors.  Replaces recvsum_ which had replaced MCT_Recvsum.
451!           Use DefaultTag and add optional Tag argument
452! 25Jul03 - R. Jacob <jacob@mcs.anl.gov> - break into irecv_ and waitrecv_
453!EOP ___________________________________________________________________
454
455  character(len=*),parameter :: myname_=myname//'::irecv_'
456  integer ::    numi,numr,i,j,k,ier
457  integer ::    mycomp,othercomp
458  integer ::    seg_start,seg_end
459  integer ::    proc,numprocs,nseg,mytag
460  integer ::    mp_Type_rp1
461  logical ::    DoSum
462
463!--------------------------------------------------------
464
465! Return if no one to receive from
466  if(Rout%nprocs .eq. 0 ) RETURN
467
468!check Av size against Router
469!
470  if(lsize(aV) /= Rout%lAvsize) then
471    write(stderr,'(2a)') myname_, &
472    ' MCTERROR:  AV size not appropriate for this Router...exiting'
473    call die(myname_)
474  endif
475
476  DoSum = .false.
477  if(present(Sum)) DoSum=Sum
478
479
480  mycomp=Rout%comp1id
481  othercomp=Rout%comp2id
482
483!  find total number of real and integer vectors
484! for now, assume we are receiving all of them
485  Rout%numiatt = nIAttr(aV)
486  Rout%numratt = nRAttr(aV)
487  numi = Rout%numiatt
488  numr = Rout%numratt
489
490!!!!!!!!!!!!!! IF RECEVING INTEGER DATA
491  if(numi .ge. 1) then
492
493! allocate buffers to hold all incoming data
494  do proc=1,Rout%nprocs
495     allocate(Rout%ip1(proc)%pi(Rout%locsize(proc)*numi),stat=ier)
496     if(ier/=0) call die(myname_,'allocate(Rout%ip1%pi)',ier)
497  enddo
498
499  endif
500
501!!!!!!!!!!!!!! IF RECEIVING REAL DATA
502  if(numr .ge. 1) then
503
504! allocate buffers to hold all incoming data
505  do proc=1,Rout%nprocs
506     allocate(Rout%rp1(proc)%pr(Rout%locsize(proc)*numr),stat=ier)
507     if(ier/=0) call die(myname_,'allocate(Rout%rp1%pr)',ier)
508  enddo
509
510  mp_Type_rp1=MP_Type(Rout%rp1(1)%pr(1))
511
512  endif
513
514  ! Post all MPI_IRECV
515  do proc=1,Rout%nprocs
516
517     ! receive the integer data
518     if(numi .ge. 1) then
519
520        ! set tag
521        mytag = DefaultTag
522        if(present(Tag)) mytag=Tag
523
524        if( Rout%num_segs(proc) > 1 .or. DoSum ) then
525
526           call MPI_IRECV(Rout%ip1(proc)%pi(1), &
527                Rout%locsize(proc)*numi,MP_INTEGER,Rout%pe_list(proc), &
528                mytag,ThisMCTWorld%MCT_comm,Rout%ireqs(proc),ier)
529
530        else
531
532           call MPI_IRECV(aV%iAttr(1,Rout%seg_starts(proc,1)), &
533                Rout%locsize(proc)*numi,MP_INTEGER,Rout%pe_list(proc), &
534                mytag,ThisMCTWorld%MCT_comm,Rout%ireqs(proc),ier)
535
536        endif
537
538        if(ier /= 0) call MP_perr_die(myname_,'MPI_IRECV(ints)',ier)
539
540     endif
541
542     ! receive the real data
543     if(numr .ge. 1) then
544
545        ! corresponding tag logic must be in send_
546        mytag = DefaultTag + 1
547        if(present(Tag)) mytag=Tag +1
548
549        if( Rout%num_segs(proc) > 1 .or. DoSum ) then
550
551           call MPI_IRECV(Rout%rp1(proc)%pr(1), &
552                Rout%locsize(proc)*numr,mp_Type_rp1,Rout%pe_list(proc), &
553                mytag,ThisMCTWorld%MCT_comm,Rout%rreqs(proc),ier)
554
555        else
556
557           call MPI_IRECV(aV%rAttr(1,Rout%seg_starts(proc,1)), &
558                Rout%locsize(proc)*numr,mp_Type_rp1,Rout%pe_list(proc), &
559                mytag,ThisMCTWorld%MCT_comm,Rout%rreqs(proc),ier)
560
561        endif
562
563        if(ier /= 0) call MP_perr_die(myname_,'MPI_IRECV(reals)',ier)
564
565     endif
566
567  enddo
568
569end subroutine irecv_
570
571!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
572!    Math and Computer Science Division, Argonne National Laboratory   !
573!BOP -------------------------------------------------------------------
574!
575! !IROUTINE: waitrecv_ - Wait for a distributed non-blocking recv to complete
576!
577! !DESCRIPTION:
578! Wait for the data being received with the {\tt Router} {\tt Rout} to complete.
579! When done, copy the data into the {\tt AttrVect} {\tt aV}.
580!
581! !INTERFACE:
582
583 subroutine waitrecv_(aV, Rout, Sum)
584
585!
586! !USES:
587!
588      implicit none
589
590! !INPUT/OUTPUT PARAMETERS:
591!
592      Type(AttrVect),       intent(inout) :: aV
593      Type(Router),         intent(inout) :: Rout
594
595! !INPUT PARAMETERS:
596!
597      logical,optional,     intent(in)    :: Sum
598
599
600! !REVISION HISTORY:
601! 25Jul03 - R. Jacob <jacob@mcs.anl.gov> - First working version is the wait
602!           and copy parts from old recv_.
603! 25Jan08 - R. Jacob <jacob@mcs.anl.gov> - Handle unordered GSMaps by
604!           applying permutation to received array.
605!EOP ___________________________________________________________________
606  character(len=*),parameter :: myname_=myname//'::waitrecv_'
607  integer ::   proc,ier,j,k,nseg
608  integer ::   AttrIndex,VectIndex,seg_start,seg_end
609  logical ::   DoSum
610  logical ::   unordered
611
612! Return if nothing to wait for
613  if(Rout%nprocs .eq. 0 ) RETURN
614
615!check Av size against Router
616!
617  if(lsize(aV) /= Rout%lAvsize) then
618    write(stderr,'(2a)') myname_, &
619    ' MCTERROR:  AV size not appropriate for this Router...exiting'
620    call die(myname_)
621  endif
622
623  unordered = associated(Rout%permarr)
624
625  DoSum = .false.
626  if(present(Sum)) DoSum=Sum
627
628  ! wait for all recieves to complete
629  if(Rout%numiatt .ge. 1)  then
630
631    call MPI_WAITALL(Rout%nprocs,Rout%ireqs,Rout%istatus,ier)
632    if(ier /= 0) call MP_perr_die(myname_,'MPI_WAITALL(ints)',ier)
633
634  endif
635
636  if(Rout%numratt .ge. 1) then
637
638    call MPI_WAITALL(Rout%nprocs,Rout%rreqs,Rout%rstatus,ier)
639    if(ier /= 0) call MP_perr_die(myname_,'MPI_WAITALL(reals)',ier)
640
641  endif
642
643  ! Load data which came from each processor
644  do proc=1,Rout%nprocs
645
646     if( (Rout%num_segs(proc) > 1) .or. DoSum ) then
647
648        j=1
649        k=1
650
651        if(DoSum) then
652        ! sum the correct pieces of the integer and real vectors
653          do nseg = 1,Rout%num_segs(proc)
654           seg_start = Rout%seg_starts(proc,nseg)
655           seg_end = seg_start + Rout%seg_lengths(proc,nseg)-1
656           do VectIndex = seg_start,seg_end
657              do AttrIndex = 1,Rout%numiatt
658                 aV%iAttr(AttrIndex,VectIndex)= &
659                 aV%iAttr(AttrIndex,VectIndex)+Rout%ip1(proc)%pi(j)
660                 j=j+1
661              enddo
662              do AttrIndex = 1,Rout%numratt
663                 aV%rAttr(AttrIndex,VectIndex)= &
664                 aV%rAttr(AttrIndex,VectIndex)+Rout%rp1(proc)%pr(k)
665                 k=k+1
666              enddo
667           enddo
668          enddo
669        else
670        ! load the correct pieces of the integer and real vectors
671          do nseg = 1,Rout%num_segs(proc)
672           seg_start = Rout%seg_starts(proc,nseg)
673           seg_end = seg_start + Rout%seg_lengths(proc,nseg)-1
674           do VectIndex = seg_start,seg_end
675              do AttrIndex = 1,Rout%numiatt
676                 aV%iAttr(AttrIndex,VectIndex)=Rout%ip1(proc)%pi(j)
677                 j=j+1
678              enddo
679              do AttrIndex = 1,Rout%numratt
680                 aV%rAttr(AttrIndex,VectIndex)=Rout%rp1(proc)%pr(k)
681                 k=k+1
682              enddo
683           enddo
684          enddo
685        endif
686       
687     endif
688
689  enddo
690
691!........................WAITANY METHOD................................
692!
693!....NOTE: Make status argument a 1-dimensional array
694!  ! Load data which came from each processor
695!  do numprocs = 1,Rout%nprocs
696!     ! Load the integer data
697!     if(Rout%numiatt .ge. 1) then
698!       call MPI_WAITANY(Rout%nprocs,Rout%ireqs,proc,Rout%istatus,ier)
699!       if(ier /= 0) call MP_perr_die(myname_,'MPI_WAITANY(ints)',ier)
700!       j=1
701!       ! load the correct pieces of the integer vectors
702!       do nseg = 1,Rout%num_segs(proc)
703!          seg_start = Rout%seg_starts(proc,nseg)
704!          seg_end = seg_start + Rout%seg_lengths(proc,nseg)-1
705!          do VectIndex = seg_start,seg_end
706!             do AttrIndex = 1,Rout%numiatt
707!                aV%iAttr(AttrIndex,VectIndex)=Rout%ip1(proc)%pi(j)
708!                j=j+1
709!             enddo
710!          enddo
711!       enddo
712!     endif
713!     ! Load the real data
714!     if(numr .ge. 1) then
715!       call MPI_WAITANY(Rout%nprocs,Rout%rreqs,proc,Rout%rstatus,ier)
716!       if(ier /= 0) call MP_perr_die(myname_,'MPI_WAITANY(reals)',ier)
717!       k=1
718!       ! load the correct pieces of the real vectors
719!       do nseg = 1,Rout%num_segs(proc)
720!          seg_start = Rout%seg_starts(proc,nseg)
721!          seg_end = seg_start + Rout%seg_lengths(proc,nseg)-1
722!          do VectIndex = seg_start,seg_end
723!             do AttrIndex = 1,numr
724!                aV%rAttr(AttrIndex,VectIndex)=Rout%rp1(proc)%pr(k)
725!                k=k+1
726!             enddo
727!          enddo
728!       enddo
729!     endif
730!  enddo
731!........................................................................
732
733  ! Deallocate all structures
734  if(Rout%numiatt .ge. 1) then
735
736     ! Deallocate the receive buffers
737     do proc=1,Rout%nprocs
738        deallocate(Rout%ip1(proc)%pi,stat=ier)
739        if(ier/=0) call die(myname_,'deallocate(Rout%ip1%pi)',ier)
740     enddo
741
742  endif
743
744  if(Rout%numratt .ge. 1) then
745
746     ! Deallocate the receive buffers
747     do proc=1,Rout%nprocs
748        deallocate(Rout%rp1(proc)%pr,stat=ier)
749        if(ier/=0) call die(myname_,'deallocate(Rout%rp1%pr)',ier)
750     enddo
751
752  endif
753
754  if (unordered) call Unpermute(aV,Rout%permarr)
755
756end subroutine waitrecv_
757
758!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
759!    Math and Computer Science Division, Argonne National Laboratory   !
760!BOP -------------------------------------------------------------------
761!
762! !IROUTINE: recv_ - Distributed receive of an Attribute Vector
763!
764! !DESCRIPTION:
765! Recieve into the {\tt AttrVect} {\tt aV} the data coming from the
766! component specified in the {\tt Router} {\tt Rout}.  An error will
767! result if the size of the attribute vector does not match the size
768! parameter stored in the {\tt Router}.
769!
770! Requires a corresponding {\tt send\_} or {\tt isend\_}to be called
771! on the other component.
772!
773! The optional argument {\tt Tag} can be used to set the tag value used in
774! the data transfer.  DefaultTag will be used otherwise. {\tt Tag} must be
775! the same in the matching {\tt send\_}
776!
777! If data for a grid point is coming from more than one process, {\tt recv\_}
778! will overwrite the duplicate values leaving the last received value
779! in the output aV.  If the optional argument {\tt Sum} is invoked, the output
780! will contain the sum of any duplicate values received for the same grid point.
781!
782! Will not return until all data has been received.
783!
784! {\bf N.B.:} The {\tt AttrVect} argument in the corresponding
785! {\tt send\_} call is assumed to have exactly the same attributes
786! in exactly the same order as {\tt aV}.
787!
788! !INTERFACE:
789
790 subroutine recv_(aV, Rout, Tag, Sum)
791!
792! !USES:
793!
794      implicit none
795
796! !INPUT/OUTPUT PARAMETERS:
797!
798      Type(AttrVect),       intent(inout) :: aV
799
800! !INPUT PARAMETERS:
801!
802      Type(Router),         intent(inout)    :: Rout
803      integer,optional,     intent(in)    :: Tag
804      logical,optional,     intent(in)    :: Sum
805
806! !REVISION HISTORY:
807! 25Jul03 - R. Jacob <jacob@mcs.anl.gov> - Rewrite using irecv and waitrecv
808!EOP ___________________________________________________________________
809  character(len=*),parameter :: myname_=myname//'::recv_'
810
811  call irecv_(aV,Rout,Tag,Sum)
812
813  call waitrecv_(aV,Rout,Sum)
814
815end subroutine recv_
816
817
818end module m_Transfer
Note: See TracBrowser for help on using the repository browser.