source: CPL/oasis3-mct/branches/OASIS3-MCT_2.0_branch/lib/mct/mct/m_AttrVectComms.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: 51.8 KB
Line 
1!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2!    Math and Computer Science Division, Argonne National Laboratory   !
3!-----------------------------------------------------------------------
4! CVS m_AttrVectComms.F90,v 1.39 2009-03-05 16:41:47 jacob Exp
5! CVS MCT_2_8_0
6!BOP -------------------------------------------------------------------
7!
8! !MODULE: m_AttrVectComms - MPI Communications Methods for the AttrVect
9!
10! !DESCRIPTION:
11!
12! This module defines the communications methods for the {\tt AttrVect}
13! datatype (see the module {\tt m\_AttrVect} for more information about
14! this class and its methods).  MCT's communications are implemented
15! in terms of the Message Passing Interface (MPI) standard, and we have
16! as best as possible, made the interfaces to these routines appear as
17! similar as possible to the corresponding MPI routines.  For the
18! { \tt AttrVect}, we supply {\em blocking} point-to-point send and
19! receive operations.  We also supply the following collective
20! operations: broadcast, gather, and scatter.  The gather and scatter
21! operations rely on domain decomposition descriptors that are defined
22! elsewhere in MCT:  the {\tt GlobalMap}, which is a one-dimensional
23! decomposition (see the MCT module {\tt m\_GlobalMap} for more details);
24! and the {\tt GlobalSegMap}, which is a segmented decomposition capable
25! of supporting multidimensional domain decompositions (see the MCT module
26! {\tt m\_GlobalSegMap} for more details).
27!
28! !INTERFACE:
29 module m_AttrVectComms
30!
31! !USES:
32!
33      use m_AttrVect ! AttrVect class and its methods
34
35      implicit none
36
37      private   ! except
38
39      public :: gather          ! gather all local vectors to the root
40      public :: scatter         ! scatter from the root to all PEs
41      public :: bcast           ! bcast from root to all PEs
42      public :: send            ! send an AttrVect
43      public :: recv            ! receive an AttrVect
44
45    interface gather ; module procedure &
46              GM_gather_, &
47              GSM_gather_ 
48    end interface
49    interface scatter ; module procedure &
50              GM_scatter_, &
51              GSM_scatter_ 
52    end interface
53    interface bcast  ; module procedure bcast_  ; end interface
54    interface send  ; module procedure send_  ; end interface
55    interface recv  ; module procedure recv_  ; end interface
56
57! !REVISION HISTORY:
58! 27Oct00 - J.W. Larson <larson@mcs.anl.gov> - relocated routines
59!           from m_AttrVect to create this module.
60! 15Jan01 - J.W. Larson <larson@mcs.anl.gov> - Added APIs for
61!           GSM_gather_() and GSM_scatter_().
62!  9May01 - J.W. Larson <larson@mcs.anl.gov> - Modified GM_scatter_
63!           so its communication model agrees with MPI_scatter().
64!           Also tidied up prologues in all module routines.
65!  7Jun01 - J.W. Larson <larson@mcs.anl.gov> - Added send()
66!           and recv().
67!  3Aug01 - E.T. Ong <eong@mcs.anl.gov> - in GSM_scatter, call
68!           GlobalMap_init with actual shaped array to satisfy
69!           Fortran 90 standard. See comment in subroutine.
70! 23Aug01 - E.T. Ong <eong@mcs.anl.gov> - replaced assignment(=)
71!           with copy for list type to avoid compiler bugs in pgf90.
72!           Added more error checking in gsm scatter. Fixed minor bugs
73!          in gsm and gm gather.
74! 13Dec01 - E.T. Ong <eong@mcs.anl.gov> - GSM_scatter, allow users
75!           to scatter with a haloed GSMap. Fixed some bugs in
76!           GM_scatter.
77! 19Dec01 - E.T. Ong <eong@mcs.anl.gov> - allow bcast of an AttrVect
78!           with only an integer or real attribute.
79! 27Mar02 - J.W. Larson <larson@mcs.anl.gov> - Corrected usage of
80!           m_die routines throughout this module.
81!EOP ___________________________________________________________________
82
83  character(len=*),parameter :: myname='MCT::m_AttrVectComms'
84
85 contains
86
87!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
88!    Math and Computer Science Division, Argonne National Laboratory   !
89!BOP -------------------------------------------------------------------
90!
91! !IROUTINE: send_ - Point-to-point Send of an AttrVect
92!
93! !DESCRIPTION:  This routine takes an input {\tt AttrVect} argument
94! {\tt inAV} and sends it to processor {\tt dest} on the communicator
95! associated with the Fortran {\tt INTEGER} MPI communicator handle
96! {\tt comm}.  The overalll message is tagged by the input {\tt INTEGER}
97! argument {\tt TagBase}.  The success (failure) of this operation is
98! reported in the zero (nonzero) optional output argument {\tt status}.
99!
100! {\bf N.B.}:  One must avoid assigning elsewhere the MPI tag values
101! between {\tt TagBase} and {\tt TagBase+7}, inclusive.  This is
102! because {\tt send\_()} performs the send of the {\tt AttrVect} as
103! a series of eight send operations.
104!
105! !INTERFACE:
106
107 subroutine send_(inAV, dest, TagBase, comm, status)
108!
109! !USES:
110!
111      use m_stdio
112      use m_mpif90
113      use m_die
114
115      use m_List, only : List
116      use m_List, only : List_allocated => allocated
117      use m_List, only : List_nitem => nitem
118      use m_List, only : List_send => send
119
120      use m_AttrVect, only : AttrVect
121      use m_AttrVect, only : AttrVect_lsize => lsize
122
123      implicit none
124
125! !INPUT PARAMETERS:
126!
127      type(AttrVect),     intent(in)  :: inAV
128      integer,            intent(in)  :: dest
129      integer,            intent(in)  :: TagBase
130      integer,            intent(in)  :: comm
131
132! !OUTPUT PARAMETERS:
133!
134      integer, optional,  intent(out) :: status
135
136! !REVISION HISTORY:
137!  7Jun01 - J.W. Larson - initial version.
138! 13Jun01 - J.W. Larson <larson@mcs.anl.gov> - Initialize status
139!           (if present).
140!EOP ___________________________________________________________________
141
142 character(len=*),parameter :: myname_=myname//'::send_'
143
144 logical :: ListAssoc(2)
145 integer :: ierr
146 integer :: AVlength
147
148      ! Initialize status (if present)
149
150  if(present(status)) status = 0
151
152
153       ! Step 1. Are inAV%iList and inAV%rList filled?  Store
154       ! the answers in the LOGICAL array ListAssoc and send.
155
156  ListAssoc(1) = List_allocated(inAV%iList) 
157  ListAssoc(2) = List_allocated(inAV%rList) 
158
159  if(.NOT. (ListAssoc(1).or.ListAssoc(2)) ) then
160     call die(myname_,"inAV has not been initialized")
161  endif
162     
163  call MPI_SEND(ListAssoc, 2, MP_LOGICAL, dest, TagBase, comm, ierr)
164  if(ierr /= 0) then
165     call MP_perr_die(myname_,':: MPI_SEND(ListAssoc...',ierr)
166  endif
167
168
169       ! Step 2. Send non-blank inAV%iList and inAV%rList.
170
171  if(ListAssoc(1)) then
172    call List_send(inAV%iList, dest, TagBase+1, comm, ierr)
173    if(ierr /= 0) then
174       if(present(status)) then
175          write(stderr,*) myname_,':: call List_send(inAV%iList...'
176          status = ierr
177          return
178       else
179          call die(myname_,':: call List_send(inAV%iList...',ierr)
180       endif
181    endif
182  endif
183
184  if(ListAssoc(2)) then
185    call List_send(inAV%rList, dest, TagBase+3, comm, ierr)
186    if(ierr /= 0) then
187       if(present(status)) then
188          write(stderr,*) myname_,':: call List_send(inAV%rList...'
189          status = ierr
190          return
191       else
192          call die(myname_,':: call List_send(inAV%rList...',ierr)
193       endif
194    endif
195  endif
196
197       ! Step 3. Determine and send the lengths of inAV%iAttr(:,:)
198       ! and inAV%rAttr(:,:).
199
200  AVlength = AttrVect_lsize(inAV)
201
202  if(AVlength<=0) then
203     call die(myname_,"Size of inAV <= 0",AVLength)
204  endif
205
206  call MPI_SEND(AVlength, 1, MP_type(AVlength), dest, TagBase+5, &
207                comm, ierr)
208  if(ierr /= 0) then
209     call MP_perr_die(myname_,':: call MPI_SEND(AVlength...',ierr)
210  endif
211
212       ! Step 4. If AVlength > 0, we may have INTEGER and REAL
213       ! data to send.  Send as needed.
214
215  if(AVlength > 0) then
216
217     if(ListAssoc(1)) then
218
219       ! Send the INTEGER data stored in inAV%iAttr(:,:)
220
221        call MPI_SEND(inAV%iAttr(1,1), AVlength*List_nitem(inAV%iList), &
222                      MP_type(inAV%iAttr(1,1)), dest, TagBase+6,   &
223                      comm, ierr)
224        if(ierr /= 0) then
225           call MP_perr_die(myname_,':: call MPI_SEND(inAV%iAttr...',ierr)
226        endif
227
228     endif ! if(associated(inAV%rList))
229
230     if(ListAssoc(2)) then
231
232       ! Send the REAL data stored in inAV%rAttr(:,:)
233
234        call MPI_SEND(inAV%rAttr(1,1), AVlength*List_nitem(inAV%rList), &
235                      MP_type(inAV%rAttr(1,1)), dest, TagBase+7,   &
236                      comm, ierr)
237        if(ierr /= 0) then
238           call MP_perr_die(myname_,':: call MPI_SEND(inAV%rAttr...',ierr)
239        endif
240
241     endif ! if(associated(inAV%rList))
242
243  endif ! if (AVlength > 0)
244
245 end subroutine send_
246
247!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
248!    Math and Computer Science Division, Argonne National Laboratory   !
249!BOP -------------------------------------------------------------------
250!
251! !IROUTINE: recv_ - Point-to-point Receive of an AttrVect
252!
253! !DESCRIPTION:  This routine receives the output {\tt AttrVect} argument
254! {\tt outAV} from processor {\tt source} on the communicator associated
255! with the Fortran {\tt INTEGER} MPI communicator handle {\tt comm}.  The
256! overall message is tagged by the input {\tt INTEGER} argument
257! {\tt TagBase}.  The success (failure) of this operation is reported in
258! the zero (nonzero) optional output argument {\tt status}.
259!
260! {\bf N.B.}:  One must avoid assigning elsewhere the MPI tag values
261! between {\tt TagBase} and {\tt TagBase+7}, inclusive.  This is
262! because {\tt recv\_()} performs the receive of the {\tt AttrVect} as
263! a series of eight receive operations.
264!
265! !INTERFACE:
266
267 subroutine recv_(outAV, dest, TagBase, comm, status)
268!
269! !USES:
270!
271      use m_stdio
272      use m_mpif90
273      use m_die
274
275      use m_List, only : List
276      use m_List, only : List_nitem => nitem
277      use m_List, only : List_recv => recv
278
279      use m_AttrVect, only : AttrVect
280
281      implicit none
282
283! !INPUT PARAMETERS:
284!
285      integer,            intent(in)  :: dest
286      integer,            intent(in)  :: TagBase
287      integer,            intent(in)  :: comm
288
289! !OUTPUT PARAMETERS:
290!
291      type(AttrVect),     intent(out) :: outAV
292      integer, optional,  intent(out) :: status
293
294! !REVISION HISTORY:
295!  7Jun01 - J.W. Larson - initial working version.
296! 13Jun01 - J.W. Larson <larson@mcs.anl.gov> - Initialize status
297!           (if present).
298!EOP ___________________________________________________________________
299
300 character(len=*),parameter :: myname_=myname//'::recv_'
301
302 logical :: ListAssoc(2)
303 integer :: ierr
304 integer :: AVlength
305 integer :: MPstatus(MP_STATUS_SIZE)
306
307      ! Initialize status (if present)
308
309  if(present(status)) status = 0
310
311
312       ! Step 1. Are outAV%iList and outAV%rList filled?  TRUE
313       ! entries in the LOGICAL array ListAssoc(:) correspond
314       ! to Non-blank Lists...that is:
315       !
316       ! ListAssoc(1) = .TRUE. <==> associated(outAV%iList%bf)
317       ! ListAssoc(2) = .TRUE. <==> associated(outAV%rList%bf)
318
319  call MPI_RECV(ListAssoc, 2, MP_LOGICAL, dest, TagBase, comm, &
320                MPstatus, ierr)
321  if(ierr /= 0) then
322     call MP_perr_die(myname_,':: MPI_RECV(ListAssoc...',ierr)
323  endif
324
325
326       ! Step 2. Receive non-blank outAV%iList and outAV%rList.
327
328  if(ListAssoc(1)) then
329    call List_recv(outAV%iList, dest, TagBase+1, comm, ierr)
330    if(ierr /= 0) then
331       if(present(status)) then
332          write(stderr,*) myname_,':: call List_recv(outAV%iList...'
333          status = ierr
334          return
335       else
336          call die(myname_,':: call List_recv(outAV%iList...',ierr)
337       endif
338    endif
339  endif
340
341  if(ListAssoc(2)) then
342    call List_recv(outAV%rList, dest, TagBase+3, comm, ierr)
343    if(ierr /= 0) then
344       if(present(status)) then
345          write(stderr,*) myname_,':: call List_recv(outAV%rList...'
346          status = ierr
347          return
348       else
349          call die(myname_,':: call List_recv(outAV%rList...',ierr)
350       endif
351    endif
352  endif
353
354       ! Step 3. Receive the lengths of outAV%iAttr(:,:) and outAV%rAttr(:,:).
355
356  call MPI_RECV(AVlength, 1, MP_type(AVlength), dest, TagBase+5, &
357                comm, MPstatus, ierr)
358  if(ierr /= 0) then
359     call MP_perr_die(myname_,':: call MPI_RECV(AVlength...',ierr)
360  endif
361
362       ! Step 4. If AVlength > 0, we may have to receive INTEGER
363       ! and/or REAL data.  Receive as needed.
364
365  if(AVlength > 0) then
366
367     if(ListAssoc(1)) then
368
369       ! Allocate outAV%iAttr(:,:)
370
371        allocate(outAV%iAttr(List_nitem(outAV%iList),AVlength), stat=ierr)
372        if(ierr/=0) call die(myname_,"allocate(outAV%iAttr)",ierr)
373
374       ! Receive the INTEGER data to outAV%iAttr(:,:)
375
376        call MPI_RECV(outAV%iAttr(1,1), AVlength*List_nitem(outAV%iList), &
377                      MP_type(outAV%iAttr(1,1)), dest, TagBase+6,   &
378                      comm, MPstatus, ierr)
379        if(ierr /= 0) then
380           call MP_perr_die(myname_,':: call MPI_RECV(outAV%iAttr...',ierr)
381        endif
382
383     endif ! if(associated(outAV%rList))
384
385     if(ListAssoc(2)) then
386
387       ! Allocate outAV%rAttr(:,:)
388
389        allocate(outAV%rAttr(List_nitem(outAV%rList),AVlength), stat=ierr)
390        if(ierr/=0) call die(myname_,"allocate(outAV%rAttr)",ierr)
391
392       ! Receive the REAL data to outAV%rAttr(:,:)
393
394        call MPI_RECV(outAV%rAttr(1,1), AVlength*List_nitem(outAV%rList), &
395                      MP_type(outAV%rAttr(1,1)), dest, TagBase+7,   &
396                      comm, MPstatus, ierr)
397        if(ierr /= 0) then
398           call MP_perr_die(myname_,':: call MPI_RECV(outAV%rAttr...',ierr)
399        endif
400
401     endif ! if(associated(outAV%rList))
402
403  endif ! if (AVlength > 0)
404
405 end subroutine recv_
406
407!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
408!    Math and Computer Science Division, Argonne National Laboratory   !
409!BOP -------------------------------------------------------------------
410!
411! !IROUTINE: GM_gather_ - Gather an AttrVect Distributed by a GlobalMap
412!
413! !DESCRIPTION:
414! This routine gathers a {\em distributed} {\tt AttrVect} {\tt iV} to
415! the {\tt root} process, and returns it in the output {\tt AttrVect}
416! argument {\tt oV}.  The decomposition of {\tt iV} is described by
417! the input {\tt GlobalMap} argument {\tt GMap}.  The input {\tt INTEGER}
418! argument {\tt comm} is the Fortran integer MPI communicator handle.
419! The success (failure) of this operation corresponds to a zero (nonzero)
420! value of the optional output {\tt INTEGER} argument {\tt stat}.
421!
422! !INTERFACE:
423
424 subroutine GM_gather_(iV, oV, GMap, root, comm, stat)
425!
426! !USES:
427!
428      use m_stdio
429      use m_die
430      use m_mpif90
431      use m_realkinds, only : FP
432      use m_GlobalMap, only : GlobalMap
433      use m_GlobalMap, only : GlobalMap_lsize => lsize
434      use m_GlobalMap, only : GlobalMap_gsize => gsize
435      use m_AttrVect, only : AttrVect
436      use m_AttrVect, only : AttrVect_init => init
437      use m_AttrVect, only : AttrVect_zero => zero 
438      use m_AttrVect, only : AttrVect_lsize => lsize
439      use m_AttrVect, only : AttrVect_nIAttr => nIAttr
440      use m_AttrVect, only : AttrVect_nRAttr => nRAttr
441      use m_AttrVect, only : AttrVect_clean => clean
442      use m_FcComms,  only : fc_gatherv_int, fc_gatherv_fp
443
444      implicit none
445
446! !INPUT PARAMETERS:
447!
448      type(AttrVect),           intent(in)  :: iV
449      type(GlobalMap),          intent(in)  :: GMap
450      integer,                  intent(in)  :: root
451      integer,                  intent(in)  :: comm
452
453! !OUTPUT PARAMETERS:
454!
455      type(AttrVect),           intent(out) :: oV
456      integer,        optional, intent(out) :: stat
457
458! !REVISION HISTORY:
459! 15Apr98 - Jing Guo <guo@thunder> - initial prototype/prolog/code
460! 27Oct00 - J.W. Larson <larson@mcs.anl.gov> - relocated from
461!           m_AttrVect
462! 15Jan01 - J.W. Larson <larson@mcs.anl.gov> - renamed GM_gather_
463!  9May01 - J.W. Larson <larson@mcs.anl.gov> - tidied up prologue
464! 18May01 - R.L. Jacob <jacob@mcs.anl.gov> - use MP_Type function
465!           to determine type for mpi_gatherv
466! 31Jan09 - P.H. Worley <worleyph@ornl.gov> - replaced call to
467!           MPI_gatherv with call to flow controlled gather routines
468!EOP ___________________________________________________________________
469
470  character(len=*),parameter :: myname_=myname//'::GM_gather_'
471  integer :: nIA,nRA,niV,noV,ier
472  integer :: myID
473  integer :: mp_type_Av
474  type(AttrVect) :: nonRootAV
475
476  if(present(stat)) stat=0
477
478  call MP_comm_rank(comm, myID, ier)
479  if(ier /= 0) then
480     call MP_perr_die(myname_,':: call MP_COMM_RANK()',ier)
481  endif
482
483        ! Verify the input: a _scatterd_ vector
484
485  niV=GlobalMap_lsize(GMap)
486  noV=AttrVect_lsize(iV)
487
488  if(niV /= noV) then
489     write(stderr,'(2a,i4,a,i4,a,i4)') myname_, &
490          ': invalid input, lsize(GMap) =',niV, &
491          ', lsize(iV) =',noV, 'myID =', myID
492     if(.not.present(stat)) call die(myname_)
493     stat=-1
494     return
495  endif
496
497  noV=GlobalMap_gsize(GMap) ! the gathered local size, as for the output
498
499  if(myID == root) then
500     call AttrVect_init(oV,iV,noV)
501     call AttrVect_zero(oV)
502  else
503     call AttrVect_init(nonRootAV,iV,1)
504     call AttrVect_zero(nonRootAV)
505  endif
506
507  niV=GlobalMap_lsize(GMap) ! the scattered local size, as for the input
508
509  nIA=AttrVect_nIAttr(iV)       ! number of INTEGER attributes
510  nRA=AttrVect_nRAttr(iV)       ! number of REAL attributes
511
512  mp_type_Av = MP_Type(1._FP)   ! set mpi type to same as AV%rAttr
513
514  if(nIA > 0) then
515     
516     if(myID == root) then
517
518        call fc_gatherv_int(iV%iAttr,niV*nIA,MP_INTEGER,                &
519             oV%iAttr,GMap%counts*nIA,GMap%displs*nIA,             &
520             MP_INTEGER,root,comm)
521
522     else
523       
524        call fc_gatherv_int(iV%iAttr,niV*nIA,MP_INTEGER,                &
525             nonRootAV%iAttr,GMap%counts*nIA,GMap%displs*nIA,      &
526             MP_INTEGER,root,comm)
527
528     endif  ! if(myID == root)
529       
530  endif  ! if(nIA > 0)
531
532  if(nRA > 0) then
533
534     if(myID == root) then
535
536        call fc_gatherv_fp(iV%rAttr,niV*nRA,mp_type_Av,         &
537             oV%rAttr,GMap%counts*nRA,GMap%displs*nRA,             & 
538             mp_type_Av,root,comm)
539
540     else
541
542        call fc_gatherv_fp(iV%rAttr,niV*nRA,mp_type_Av,         &
543             nonRootAV%rAttr,GMap%counts*nRA,GMap%displs*nRA,      &
544             mp_type_Av,root,comm)
545
546     endif  ! if(myID == root)
547
548  endif  ! if(nRA > 0)
549
550
551
552  if(myID /= root) then
553     call AttrVect_clean(nonRootAV,ier)
554     if(ier /= 0) then
555        write(stderr,'(2a,i4)') myname_,        &
556             ':: AttrVect_clean(nonRootAV) failed for non-root &
557             &process: myID = ', myID
558        call die(myname_,':: AttrVect_clean failed &
559             &for nonRootAV off of root',ier)
560     endif
561  endif
562
563 end subroutine GM_gather_
564
565!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
566!    Math and Computer Science Division, Argonne National Laboratory   !
567!BOP -------------------------------------------------------------------
568!
569! !IROUTINE: GSM_gather_ - Gather an AttrVect Distributed by a GlobalSegMap
570!
571! !DESCRIPTION:
572! The routine {\tt GSM\_gather\_()} takes a distributed input
573! {\tt AttrVect} argument {\tt iV}, whose decomposition is described
574! by the input {\tt GlobalSegMap} argument {\tt GSMap}, and gathers
575! it to the output {\tt AttrVect} argument {\tt oV}.  The gathered
576! {\tt AttrVect} {\tt oV} is valid only on the root process specified
577! by the input argument {\tt root}.  The communicator used to gather
578! the data is specified by the argument {\tt comm}.  The success (failure)
579! is reported in the zero (non-zero) value of the output argument
580! {\tt stat}.
581!
582! {\tt GSM\_gather\_()} converts the problem of gathering data
583! according to a {\tt GlobalSegMap} into the simpler problem of
584! gathering data as specified by a {\tt GlobalMap}.  The {\tt GlobalMap}
585! variable {\tt GMap} is created based on the local storage requirements
586! for each distributed piece of {\tt iV}.  On the root, a complete
587! (including halo points) gathered copy of {\tt iV} is collected into
588! the temporary {\tt AttrVect} variable {\tt workV} (the length of
589! {\tt workV} is the larger of {\tt GlobalSegMap\_GlobalStorage(GSMap)} or
590! {\tt GlobalSegMap\_GlobalSize(GSMap)}).  The
591! variable {\tt workV} is segmented by process, and segments are
592! copied into it by process, but ordered in the same order the segments
593! appear in {\tt GSMap}.  Once {\tt workV} is loaded, the data are
594! copied segment-by-segment to their appropriate locations in the output
595! {\tt AttrVect} {\tt oV}.
596!
597! !INTERFACE:
598
599 subroutine GSM_gather_(iV, oV, GSMap, root, comm, stat, rdefault, idefault)
600!
601! !USES:
602!
603! Message-passing environment utilities (mpeu) modules:
604      use m_stdio
605      use m_die
606      use m_mpif90
607      use m_realkinds, only: FP
608! GlobalSegMap and associated services:
609      use m_GlobalSegMap, only : GlobalSegMap
610      use m_GlobalSegMap, only : GlobalSegMap_comp_id => comp_id
611      use m_GlobalSegMap, only : GlobalSegMap_ngseg => ngseg
612      use m_GlobalSegMap, only : GlobalSegMap_lsize => lsize
613      use m_GlobalSegMap, only : GlobalSegMap_gsize => gsize
614      use m_GlobalSegMap, only : GlobalSegMap_haloed => haloed
615      use m_GlobalSegMap, only : GlobalSegMap_GlobalStorage => GlobalStorage
616! AttrVect and associated services:
617      use m_AttrVect, only : AttrVect
618      use m_AttrVect, only : AttrVect_init => init
619      use m_AttrVect, only : AttrVect_zero => zero
620      use m_AttrVect, only : AttrVect_lsize => lsize
621      use m_AttrVect, only : AttrVect_nIAttr => nIAttr
622      use m_AttrVect, only : AttrVect_nRAttr => nRAttr
623      use m_AttrVect, only : AttrVect_clean => clean 
624! GlobalMap and associated services:
625      use m_GlobalMap, only : GlobalMap
626      use m_GlobalMap, only : GlobalMap_init => init
627      use m_GlobalMap, only : GlobalMap_clean => clean
628
629      implicit none
630
631! !INPUT PARAMETERS:
632!
633      type(AttrVect),            intent(in)  :: iV
634      type(GlobalSegMap),        intent(in)  :: GSMap
635      integer,                   intent(in)  :: root
636      integer,                   intent(in)  :: comm
637      real(FP), optional,        intent(in)  :: rdefault
638      integer, optional,         intent(in)  :: idefault
639
640! !OUTPUT PARAMETERS:
641!
642      type(AttrVect),            intent(out) :: oV
643      integer,        optional,  intent(out) :: stat
644
645! !REVISION HISTORY:
646! 15Jan01 - J.W. Larson <larson@mcs.anl.gov> - API specification.
647! 25Feb01 - J.W. Larson <larson@mcs.anl.gov> - Prototype code.
648! 26Apr01 - R.L. Jacob <jacob@mcs.anl.gov> - add use statement for
649!           AttVect_clean
650!  9May01 - J.W. Larson <larson@mcs.anl.gov> - tidied up prologue
651! 13Jun01 - J.W. Larson <larson@mcs.anl.gov> - Initialize stat
652!           (if present).
653! 20Aug01 - E.T. Ong <eong@mcs.anl.gov> - Added error checking for
654!           matching processors in gsmap and comm. Corrected
655!           current_pos assignment.
656! 23Nov01 - R. Jacob <jacob@mcs.anl.gov> - zero the oV before copying in
657!           gathered data.
658! 27Jul07 - R. Loy <rloy@mcs.anl.gov> - add Tony's suggested improvement
659!           for a default value in the output AV
660! 11Aug08 - R. Jacob <jacob@mcs.anl.gov> - add Pat Worley's faster way
661!           to initialize lns
662!EOP ___________________________________________________________________
663
664  character(len=*),parameter :: myname_=myname//'::GSM_gather_'
665
666! Temporary workspace AttrVect:
667  type(AttrVect) :: workV
668! Component ID and number of segments for GSMap:
669  integer :: comp_id, ngseg, iseg
670! Total length of GSMap segments laid end-to-end:
671  integer :: global_storage
672! Error Flag
673  integer :: ierr
674! Number of processes on communicator, and local rank:
675  integer :: NumProcs, myID
676! Total local storage on each pe according to GSMap:
677  integer, dimension(:), allocatable :: lns
678! Temporary GlobalMap used to scatter the segmented (by pe) data
679  type(GlobalMap) :: workGMap
680! Loop counters and temporary indices:
681  integer :: m, n, ilb, iub, olb, oub, pe
682! workV segment tracking index array:
683  integer, dimension(:), allocatable :: current_pos
684! workV sizes
685  integer :: gssize, gstorage
686
687      ! Initialize stat (if present)
688
689  if(present(stat)) stat = 0
690
691       ! Initial Check:  If GSMap contains halo points, die
692       
693  if(GlobalSegMap_haloed(GSMap)) then
694     ierr = 1
695     call die(myname_,"Input GlobalSegMap haloed--not allowed",ierr)
696  endif
697
698       ! Which process am I?
699
700  call MPI_COMM_RANK(comm, myID, ierr)
701
702  if(ierr /= 0) then
703        call MP_perr_die(myname_,':: call MPI_COMM_RANK()',ierr)
704  endif
705       ! How many processes are there on this communicator?
706
707  call MPI_COMM_SIZE(comm, NumProcs, ierr)
708
709  if(ierr /= 0) then
710     call MP_perr_die(myname_,':: call MPI_COMM_SIZE()',ierr)
711  endif
712
713       ! Processor Check: Do the processors on GSMap match those in comm?
714
715  if(MAXVAL(GSMap%pe_loc) > (NumProcs-1)) then
716     stat=2
717     write(stderr,*) myname_, &
718       ":: Procs in GSMap%pe_loc do not match procs in communicator ", &
719       NumProcs-1, MAXVAL(GSMap%pe_loc)
720     call die(myname_, &
721          "Procs in GSMap%pe_loc do not match procs in communicator",stat)
722  endif
723
724  if(myID == root) then
725
726       ! Allocate a precursor to a GlobalMap accordingly...
727
728     allocate(lns(0:NumProcs-1), stat=ierr)
729
730       ! And Load it...
731
732     lns(:)=0
733     do iseg=1,GSMap%ngseg
734        n = GSMap%pe_loc(iseg)
735        lns(n) = lns(n) + GSMap%length(iseg)
736     end do
737
738  else
739
740     allocate(lns(0)) ! This conforms to F90 standard for shaped arguments.
741
742  endif ! if(myID == root)
743
744       ! Determine the component id of GSMap:
745
746  comp_id = GlobalSegMap_comp_id(GSMap)
747
748       ! Create working GlobalMap workGMap (used for the gather):
749
750  call GlobalMap_init(workGMap, comp_id, lns, root, comm)
751
752       ! Gather the Data process-by-process to workV...
753       ! do not include stat argument; bypass an argument check in gm_gather.
754
755  call GM_gather_(iV, workV, workGMap, root, comm, stat)
756
757       ! On the root, initialize oV, and load the contents of
758       !workV into it...
759
760  if(myID == root) then
761
762! bug fix:  gstorage will be bigger than gssize if GSmap is
763! haloed.  But gstorage may be smaller than gsize if GSmap
764! is masked.  So take the maximum.  RLJ
765     gstorage = GlobalSegMap_GlobalStorage(GSMap)
766     gssize = GlobalSegMap_gsize(GSMap)
767     global_storage = MAX(gstorage,gssize)
768
769     call AttrVect_init(oV,iV,global_storage)
770     call AttrVect_zero(oV)
771
772     if (present(rdefault)) then
773       if (AttrVect_nRAttr(oV) > 0) oV%rAttr=rdefault
774     endif
775     if (present(idefault)) then
776       if (AttrVect_nIAttr(oV) > 0) oV%iAttr=idefault
777     endif
778
779       ! On the root, allocate current position index for
780       ! each process chunk:
781
782     allocate(current_pos(0:NumProcs-1), stat=ierr)
783
784     if(ierr /= 0) then
785        write(stderr,*) myname_,':: allocate(current_pos(..) failed,', &
786             'stat = ',ierr
787        if(present(stat)) then
788           stat=ierr
789        else
790           call die(myname_,'allocate(current_pos(..) failed.' )
791        endif
792     endif
793
794       ! Initialize current_pos(:) using GMap%displs(:)
795
796     do n=0,NumProcs-1
797        current_pos(n) = workGMap%displs(n) + 1
798     end do
799
800       ! Load each segment of iV into its appropriate segment
801       ! of workV:
802     
803     ngseg = GlobalSegMap_ngseg(GSMap)
804
805     do n=1,ngseg
806
807       ! Determine which process owns segment n:
808
809        pe = GSMap%pe_loc(n)
810
811       ! Input map (lower/upper indicess) of segment of iV:
812
813        ilb = current_pos(pe)
814        iub = current_pos(pe) + GSMap%length(n) - 1
815
816       ! Output map of (lower/upper indicess) segment of workV:
817
818        olb = GSMap%start(n)
819        oub = GSMap%start(n) + GSMap%length(n) - 1
820
821       ! Increment current_pos(n) for next time:
822
823        current_pos(pe) = current_pos(pe) + GSMap%length(n)
824
825       ! Now we are equipped to do the copy:
826
827        do m=1,AttrVect_nIAttr(iV)
828           oV%iAttr(m,olb:oub) = workV%iAttr(m,ilb:iub)
829        end do
830
831        do m=1,AttrVect_nRAttr(iV)
832           oV%rAttr(m,olb:oub) = workV%rAttr(m,ilb:iub)
833        end do
834
835     end do ! do n=1,ngseg
836
837       ! Clean up current_pos, which was only allocated on the root
838
839     deallocate(current_pos, stat=ierr)
840     if(ierr /= 0) then
841        write(stderr,*) myname_,'error in deallocate(current_pos), stat=',ierr
842        if(present(stat)) then
843           stat=ierr
844        else
845           call die(myname_)
846        endif
847     endif
848  endif ! if(myID == root)
849
850       ! At this point, we are finished.  The data have been gathered
851       ! to oV
852
853       ! Finally, clean up allocated structures:
854
855  if(myID == root) call AttrVect_clean(workV)
856  call GlobalMap_clean(workGMap)
857
858  deallocate(lns, stat=ierr)
859
860  if(ierr /= 0) then
861    write(stderr,*) myname_,'error in deallocate(lns), stat=',ierr
862    if(present(stat)) then
863       stat=ierr
864    else
865       call die(myname_)
866    endif
867  endif
868
869 end subroutine GSM_gather_
870
871!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
872!    Math and Computer Science Division, Argonne National Laboratory   !
873!BOP -------------------------------------------------------------------
874!
875! !IROUTINE: GM_scatter_ - Scatter an AttrVect Using a GlobalMap
876!
877! !DESCRIPTION:
878! The routine {\tt GM\_scatter\_} takes an input {\tt AttrVect} type
879! {\tt iV} (valid only on the root), and scatters it to a distributed
880! {\tt AttrVect} {\tt oV}.  The input {\tt GlobalMap} argument
881! {\tt GMap} dictates how {\tt iV} is scattered to {\tt oV}.  The
882! success (failure) of this routine is reported in the zero (non-zero)
883! value of the output argument {\tt stat}.
884!
885! {\bf N.B.}:  The output {\tt AttrVect} argument {\tt oV} represents
886! dynamically allocated memory.  When it is no longer needed, it should
887! be deallocated by invoking {\tt AttrVect\_clean()} (see the module
888! {\tt m\_AttrVect} for more details).
889!
890! !INTERFACE:
891
892 subroutine GM_scatter_(iV, oV, GMap, root, comm, stat)
893!
894! !USES:
895!
896      use m_stdio
897      use m_die
898      use m_mpif90
899      use m_realkinds, only : FP
900
901      use m_List, only : List
902      use m_List, only : List_copy => copy
903      use m_List, only : List_bcast => bcast
904      use m_List, only : List_clean => clean
905      use m_List, only : List_nullify => nullify
906      use m_List, only : List_nitem => nitem
907
908      use m_GlobalMap, only : GlobalMap
909      use m_GlobalMap, only : GlobalMap_lsize => lsize
910      use m_GlobalMap, only : GlobalMap_gsize => gsize
911
912      use m_AttrVect, only : AttrVect
913      use m_AttrVect, only : AttrVect_init => init
914      use m_AttrVect, only : AttrVect_zero => zero
915      use m_AttrVect, only : AttrVect_lsize => lsize
916      use m_AttrVect, only : AttrVect_nIAttr => nIAttr
917      use m_AttrVect, only : AttrVect_nRAttr => nRAttr
918      use m_AttrVect, only : AttrVect_clean => clean
919
920      implicit none
921
922! !INPUT PARAMETERS:
923!
924      type(AttrVect),           intent(in)  :: iV
925      type(GlobalMap),          intent(in)  :: GMap
926      integer,                  intent(in)  :: root
927      integer,                  intent(in)  :: comm
928
929! !OUTPUT PARAMETERS:
930!
931      type(AttrVect),           intent(out) :: oV
932      integer,        optional, intent(out) :: stat
933
934! !REVISION HISTORY:
935! 21Apr98 - Jing Guo <guo@thunder> - initial prototype/prolog/code
936! 27Oct00 - J.W. Larson <larson@mcs.anl.gov> - relocated from
937!           m_AttrVect
938! 15Jan01 - J.W. Larson <larson@mcs.anl.gov> - renamed  GM_scatter_
939!  8Feb01 - J.W. Larson <larson@mcs.anl.gov> - add logic to prevent
940!           empty calls (i.e. no data in buffer) to MPI_SCATTERV()
941! 27Apr01 - R.L. Jacob <jacob@mcs.anl.gov> - small bug fix to
942!           integer attribute scatter
943!  9May01 - J.W. Larson <larson@mcs.anl.gov> - Re-vamped comms model
944!           to reflect MPI comms model for the scatter.  Tidied up
945!           the prologue, too.
946! 18May01 - R.L. Jacob <jacob@mcs.anl.gov> - use MP_Type function
947!           to determine type for mpi_scatterv
948!  8Aug01 - E.T. Ong <eong@mcs.anl.gov> - replace list assignment(=)
949!           with list copy to avoid compiler errors in pgf90.
950! 13Dec01 - E.T. Ong <eong@mcs.anl.gov> - allow scatter with an
951!           AttrVect containing only an iList or rList.
952!EOP ___________________________________________________________________
953
954  character(len=*),parameter :: myname_=myname//'::GM_scatter_'
955  integer :: nIA,nRA,niV,noV,ier
956  integer :: myID
957  integer :: mp_type_Av
958  type(List) :: iList, rList
959  type(AttrVect) :: nonRootAV
960
961  if(present(stat)) stat=0
962
963  call MP_comm_rank(comm,myID,ier)
964  if(ier /= 0) then
965    call MP_perr_die(myname_,'MP_comm_rank()',ier)
966  endif
967
968        ! Verify the input: a _gathered_ vector
969
970  if(myID == root) then
971
972     niV = GlobalMap_gsize(GMap)  ! the _gathered_ local size
973     noV = AttrVect_lsize(iV)     ! the length of the input AttrVect iV
974
975     if(niV /= noV) then
976        write(stderr,'(2a,i5,a,i8,a,i8)') myname_,      &
977             ': myID = ',myID,'.  Invalid input on root, gsize(GMap) =',&
978             niV,', lsize(iV) =',noV
979        if(present(stat)) then
980           stat=-1
981        else
982           call die(myname_)
983        endif
984     endif
985
986  endif
987
988        ! On the root, read the integer and real attribute
989        ! lists off of iV.
990
991  call List_nullify(iList)
992  call List_nullify(rList)
993
994  if(myID == root) then
995
996     ! Count the number of real and integer attributes
997
998     nIA = AttrVect_nIAttr(iV)  ! number of INTEGER attributes
999     nRA = AttrVect_nRAttr(iV)  ! number of REAL attributes
1000
1001     if(nIA > 0) then
1002        call List_copy(iList,iV%iList)
1003     endif
1004
1005     if(nRA > 0) then
1006        call List_copy(rList,iV%rList)
1007     endif
1008
1009  endif
1010 
1011        ! From the root, broadcast iList and rList
1012
1013  call MPI_BCAST(nIA,1,MP_INTEGER,root,comm,ier)
1014  if(ier /= 0) call MP_perr(myname_,'MPI_BCAST(nIA)',ier)
1015
1016  call MPI_BCAST(nRA,1,MP_INTEGER,root,comm,ier)
1017  if(ier /= 0) call MP_perr(myname_,'MPI_BCAST(nRA)',ier)
1018
1019  if(nIA>0) call List_bcast(iList, root, comm)
1020  if(nRA>0) call List_bcast(rList, root, comm)
1021
1022  noV = GlobalMap_lsize(GMap) ! the _scatterd_ local size
1023
1024        ! On all processes, use List data and noV to initialize oV
1025
1026  call AttrVect_init(oV, iList, rList, noV)
1027  call AttrVect_zero(oV)
1028
1029        ! Initialize a dummy AttrVect for non-root MPI calls
1030
1031  if(myID/=root) then
1032    call AttrVect_init(nonRootAV,oV,1)
1033    call AttrVect_zero(nonRootAV)
1034  endif
1035
1036
1037  if(nIA > 0) then
1038
1039     if(myID == root) then
1040
1041        call MPI_scatterv(iV%iAttr(1,1),GMap%counts*nIA,        &
1042             GMap%displs*nIA,MP_INTEGER,oV%iAttr(1,1),          &
1043             noV*nIA,MP_INTEGER,root,comm,ier )
1044        if(ier /= 0) then
1045           call MP_perr_die(myname_,'MPI_scatterv(iAttr) on root',ier)
1046        endif
1047
1048     else
1049
1050        call MPI_scatterv(nonRootAV%iAttr(1,1),GMap%counts*nIA, &
1051             GMap%displs*nIA,MP_INTEGER,oV%iAttr(1,1),          &
1052             noV*nIA,MP_INTEGER,root,comm,ier )
1053        if(ier /= 0) then
1054           call MP_perr_die(myname_,'MPI_scatterv(iAttr) off root',ier)
1055        endif
1056
1057     endif   ! if(myID == root)
1058
1059     call List_clean(iList)
1060
1061  endif   ! if(nIA > 0)
1062
1063  mp_type_Av = MP_Type(1._FP)   ! set mpi type to same as AV%rAttr
1064
1065  if(nRA > 0) then
1066
1067     if(myID == root) then
1068
1069
1070        call MPI_scatterv(iV%rAttr(1,1),GMap%counts*nRA,        &
1071             GMap%displs*nRA,mp_type_Av,oV%rAttr(1,1),          &
1072             noV*nRA,mp_type_Av,root,comm,ier )
1073        if(ier /= 0) then
1074           call MP_perr_die(myname_,'MPI_scatterv(rAttr) on root',ier)
1075        endif
1076
1077     else
1078
1079
1080        call MPI_scatterv(nonRootAV%rAttr,GMap%counts*nRA,      &
1081             GMap%displs*nRA,mp_type_Av,oV%rAttr,          &
1082             noV*nRA,mp_type_Av,root,comm,ier )
1083        if(ier /= 0) then
1084           call MP_perr_die(myname_,'MPI_scatterv(rAttr) off root',ier)
1085        endif
1086
1087     endif
1088
1089     call List_clean(rList)
1090
1091  endif
1092
1093  if(myID /= root) then
1094     call AttrVect_clean(nonRootAV,ier)
1095     if(ier /= 0) then
1096        write(stderr,'(2a,i4)') myname_,        &
1097             ':: AttrVect_clean(nonRootAV) failed for non-root &
1098             &process: myID = ', myID
1099        call die(myname_,':: AttrVect_clean failed &
1100             &for nonRootAV off of root',ier)
1101     endif
1102  endif
1103
1104 end subroutine GM_scatter_
1105
1106!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1107!    Math and Computer Science Division, Argonne National Laboratory   !
1108!BOP -------------------------------------------------------------------
1109!
1110! !IROUTINE: GSM_scatter_ - Scatter an AttrVect using a GlobalSegMap
1111!
1112! !DESCRIPTION:
1113! The routine {\tt GSM\_scatter\_} takes an input {\tt AttrVect} type
1114! {\tt iV} (valid only on the root), and scatters it to a distributed
1115! {\tt AttrVect} {\tt oV}.  The input {\tt GlobalSegMap} argument
1116! {\tt GSMap} dictates how {\tt iV} is scattered to {\tt oV}.  The
1117! success (failure) of this routine is reported in the zero (non-zero)
1118! value of the output argument {\tt stat}.
1119!
1120! {\tt GSM\_scatter\_()} converts the problem of scattering data
1121! according to a {\tt GlobalSegMap} into the simpler problem of
1122! scattering data as specified by a {\tt GlobalMap}.  The {\tt GlobalMap}
1123! variable {\tt GMap} is created based on the local storage requirements
1124! for each distributed piece of {\tt iV}.  On the root, a complete
1125! (including halo points) copy of {\tt iV} is stored in
1126! the temporary {\tt AttrVect} variable {\tt workV} (the length of
1127! {\tt workV} is {\tt GlobalSegMap\_GlobalStorage(GSMap)}).  The
1128! variable {\tt workV} is segmented by process, and segments are
1129! copied into it by process, but ordered in the same order the segments
1130! appear in {\tt GSMap}.  Once {\tt workV} is loaded, the data are
1131! scattered to the output {\tt AttrVect} {\tt oV} by a call to the
1132! routine {\tt GM\_scatter\_()} defined in this module, with {\tt workV}
1133! and {\tt GMap} as the input arguments.
1134!
1135! {\bf N.B.:}  This algorithm  assumes that memory access times are much
1136! shorter than message-passing transmission times.
1137!
1138! {\bf N.B.}:  The output {\tt AttrVect} argument {\tt oV} represents
1139! dynamically allocated memory.  When it is no longer needed, it should
1140! be deallocated by invoking {\tt AttrVect\_clean()} (see the module
1141! {\tt m\_AttrVect} for more details).
1142!
1143! !INTERFACE:
1144
1145 subroutine GSM_scatter_(iV, oV, GSMap, root, comm, stat)
1146!
1147! !USES:
1148!
1149! Environment utilities from mpeu:
1150
1151      use m_stdio
1152      use m_die
1153      use m_mpif90
1154
1155      use m_List, only : List_nullify => nullify
1156
1157! GlobalSegMap and associated services:
1158      use m_GlobalSegMap, only : GlobalSegMap
1159      use m_GlobalSegMap, only : GlobalSegMap_comp_id => comp_id
1160      use m_GlobalSegMap, only : GlobalSegMap_ngseg => ngseg
1161      use m_GlobalSegMap, only : GlobalSegMap_lsize => lsize
1162      use m_GlobalSegMap, only : GlobalSegMap_gsize => gsize
1163      use m_GlobalSegMap, only : GlobalSegMap_GlobalStorage => GlobalStorage
1164! AttrVect and associated services:
1165      use m_AttrVect, only : AttrVect
1166      use m_AttrVect, only : AttrVect_init => init
1167      use m_AttrVect, only : AttrVect_zero => zero
1168      use m_AttrVect,  only : AttrVect_lsize => lsize
1169      use m_AttrVect, only : AttrVect_nIAttr => nIAttr
1170      use m_AttrVect, only : AttrVect_nRAttr => nRAttr
1171      use m_AttrVect, only : AttrVect_clean => clean 
1172! GlobalMap and associated services:
1173      use m_GlobalMap, only : GlobalMap
1174      use m_GlobalMap, only : GlobalMap_init => init
1175      use m_GlobalMap, only : GlobalMap_clean => clean
1176
1177      implicit none
1178
1179! !INPUT PARAMETERS:
1180!
1181      type(AttrVect),            intent(in)  :: iV
1182      type(GlobalSegMap),        intent(in)  :: GSMap
1183      integer,                   intent(in)  :: root
1184      integer,                   intent(in)  :: comm
1185
1186! !OUTPUT PARAMETERS:
1187!
1188      type(AttrVect),            intent(out) :: oV
1189      integer,        optional,  intent(out) :: stat
1190
1191! !REVISION HISTORY:
1192! 15Jan01 - J.W. Larson <larson@mcs.anl.gov> - API specification.
1193!  8Feb01 - J.W. Larson <larson@mcs.anl.gov> - Initial code.
1194! 25Feb01 - J.W. Larson <larson@mcs.anl.gov> - Bug fix--replaced
1195!           call to GlobalSegMap_lsize with call to the new fcn.
1196!           GlobalSegMap_ProcessStorage().
1197! 26Apr01 - R.L. Jacob <jacob@mcs.anl.gov> - add use statement for
1198!           AttVect_clean
1199! 26Apr01 - J.W. Larson <larson@mcs.anl.gov> - bug fixes--data
1200!           misalignment in use of the GlobalMap to compute the
1201!           memory map into workV, and initialization of workV
1202!           on all processes.
1203!  9May01 - J.W. Larson <larson@mcs.anl.gov> - tidied up prologue
1204! 15May01 - Larson / Jacob <larson@mcs.anl.gov> - stopped initializing
1205!           workV on off-root processes (no longer necessary).
1206! 13Jun01 - J.W. Larson <larson@mcs.anl.gov> - Initialize stat
1207!           (if present).
1208! 20Jun01 - J.W. Larson <larson@mcs.anl.gov> - Fixed a subtle bug
1209!           appearing on AIX regarding the fact workV is uninitial-
1210!           ized on non-root processes.  This is fixed by nullifying
1211!           all the pointers in workV for non-root processes.
1212! 20Aug01 - E.T. Ong <eong@mcs.anl.gov> - Added argument check
1213!           for matching processors in gsmap and comm.
1214! 13Dec01 - E.T. Ong <eong@mcs.anl.gov> - got rid of restriction
1215!           GlobalStorage(GSMap)==AttrVect_lsize(AV) to allow for
1216!           GSMap to be haloed.
1217! 11Aug08 - R. Jacob <jacob@mcs.anl.gov> - remove call to ProcessStorage
1218!           and replace with faster algorithm provided by Pat Worley
1219!EOP ___________________________________________________________________
1220
1221  character(len=*),parameter :: myname_=myname//'::GSM_scatter_'
1222
1223! Temporary workspace AttrVect:
1224  type(AttrVect) :: workV
1225! Component ID and number of segments for GSMap:
1226  integer :: comp_id, ngseg, iseg
1227! Total length of GSMap segments laid end-to-end:
1228  integer :: global_storage
1229! Error Flag
1230  integer :: ierr
1231! Number of processes on communicator, and local rank:
1232  integer :: NumProcs, myID
1233! Total local storage on each pe according to GSMap:
1234  integer, dimension(:), allocatable :: lns
1235! Temporary GlobalMap used to scatter the segmented (by pe) data
1236  type(GlobalMap) :: GMap
1237! Loop counters and temporary indices:
1238  integer :: m, n, ilb, iub, olb, oub, pe
1239! workV segment tracking index array:
1240  integer, dimension(:), allocatable :: current_pos
1241
1242      ! Initialize stat (if present)
1243
1244  if(present(stat)) stat = 0
1245
1246       ! Which process am I?
1247
1248  call MPI_COMM_RANK(comm, myID, ierr)
1249
1250  if(ierr /= 0) then
1251     call MP_perr_die(myname_,'MPI_COMM_RANK',ierr)
1252  endif
1253
1254  if(myID == root) then
1255     
1256     if(GSMap%gsize > AttrVect_lsize(iV)) then
1257        write(stderr,'(2a,i5,a,i8,a,i8)') myname_,      &
1258             ': myID = ',myID,'.  Invalid input, GSMap%gsize =',&
1259             GSMap%gsize, ', lsize(iV) =',AttrVect_lsize(iV)
1260        if(present(stat)) then
1261           stat=-1
1262        else
1263           call die(myname_)
1264        endif
1265     endif
1266
1267  endif     
1268
1269       ! On the root, initialize a work AttrVect type of the
1270       ! above length, and with the same attribute lists as iV.
1271       ! on other processes, initialize workV only with the
1272       ! attribute information, but no storage.
1273
1274  if(myID == root) then
1275
1276     global_storage = GlobalSegMap_GlobalStorage(GSMap)
1277     call AttrVect_init(workV, iV, global_storage)
1278     call AttrVect_zero(workV)
1279
1280  else
1281       ! nullify workV just to be safe
1282 
1283     call List_nullify(workV%iList)
1284     call List_nullify(workV%rList)
1285     nullify(workV%iAttr)
1286     nullify(workV%rAttr)
1287
1288  endif
1289
1290       ! Return to processing on the root to load workV:
1291
1292       ! How many processes are there on this communicator?
1293
1294  call MPI_COMM_SIZE(comm, NumProcs, ierr)
1295
1296  if(ierr /= 0) then
1297     call MP_perr_die(myname_,'MPI_COMM_SIZE',ierr)
1298  endif
1299
1300       ! Processor Check: Do the processors on GSMap match those in comm?
1301
1302  if(MAXVAL(GSMap%pe_loc) > (NumProcs-1)) then
1303     write(stderr,*) myname_, &
1304          ":: Procs in GSMap%pe_loc do not match procs in communicator ", &
1305          NumProcs-1, MAXVAL(GSMap%pe_loc)
1306     if(present(stat)) then
1307        stat=1
1308        return
1309     else
1310        call die(myname_)
1311     endif
1312  endif
1313
1314  if(myID == root) then
1315
1316       ! Allocate a precursor to a GlobalMap accordingly...
1317
1318     allocate(lns(0:NumProcs-1), stat=ierr)
1319     if(ierr /= 0) then
1320        write(stderr,*) myname_,':: allocate(lns...) failed, stat=',ierr
1321        if(present(stat)) then
1322           stat=ierr
1323        else
1324           call die(myname_,'allocate(lns)',ierr)
1325        endif
1326     endif
1327
1328       ! And Load it...
1329
1330     lns(:)=0
1331     do iseg=1,GSMap%ngseg
1332        n = GSMap%pe_loc(iseg)
1333        lns(n) = lns(n) + GSMap%length(iseg)
1334     end do
1335
1336  endif ! if(myID == root)
1337
1338        ! Non-root processes call GlobalMap_init with lns,
1339        ! although this argument is not used in the
1340        ! subroutine. Since it correspond to a dummy shaped array arguments
1341        ! in GlobslMap_init, the Fortran 90 standard dictates that the actual
1342        ! argument must contain complete shape information. Therefore,
1343        ! the array argument must be allocated on all processes.
1344
1345  if(myID /= root) then
1346
1347     allocate(lns(1),stat=ierr)
1348     if(ierr /= 0) then
1349        write(stderr,*) myname_,':: allocate(lns...) failed, stat=',ierr
1350        if(present(stat)) then
1351           stat=ierr
1352           return
1353        else
1354           call die(myname_,'allocate(lns(1))',ierr)
1355        endif
1356     endif
1357
1358  endif ! if(myID /= root)...
1359
1360       ! Create a GlobalMap describing the 1-D decomposition
1361       ! of workV:
1362
1363  comp_id = GlobalSegMap_comp_id(GSMap)
1364
1365  call GlobalMap_init(GMap, comp_id, lns, root, comm)
1366
1367       ! On the root, load workV:
1368
1369  if(myID == root) then
1370
1371       ! On the root, allocate current position index for
1372       ! each process chunk:
1373
1374     allocate(current_pos(0:NumProcs-1), stat=ierr)
1375     if(ierr /= 0) then
1376        write(stderr,*) myname_,':: allocate(current_pos..) failed, stat=', &
1377             ierr
1378        if(present(stat)) then
1379           stat=ierr
1380           return
1381        else
1382           call die(myname_,'allocate(current_pos)',ierr)
1383        endif
1384     endif
1385
1386       ! Initialize current_pos(:) using GMap%displs(:)
1387
1388     do n=0,NumProcs-1
1389        current_pos(n) = GMap%displs(n) + 1
1390     end do
1391
1392       ! Load each segment of iV into its appropriate segment
1393       ! of workV:
1394
1395     ngseg = GlobalSegMap_ngseg(GSMap)
1396
1397     do n=1,ngseg
1398
1399       ! Determine which process owns segment n:
1400
1401        pe = GSMap%pe_loc(n)
1402
1403       ! Input map (lower/upper indicess) of segment of iV:
1404
1405        ilb = GSMap%start(n)
1406        iub = GSMap%start(n) + GSMap%length(n) - 1
1407
1408       ! Output map of (lower/upper indicess) segment of workV:
1409
1410        olb = current_pos(pe)
1411        oub = current_pos(pe) + GSMap%length(n) - 1
1412
1413       ! Increment current_pos(n) for next time:
1414
1415        current_pos(pe) = current_pos(pe) + GSMap%length(n)
1416
1417       ! Now we are equipped to do the copy:
1418
1419        do m=1,AttrVect_nIAttr(iV)
1420           workV%iAttr(m,olb:oub) = iV%iAttr(m,ilb:iub)
1421        end do
1422
1423        do m=1,AttrVect_nRAttr(iV)
1424           workV%rAttr(m,olb:oub) = iV%rAttr(m,ilb:iub)
1425        end do
1426
1427     end do ! do n=1,ngseg
1428
1429       ! Clean up current_pos, which was only allocated on the root
1430
1431     deallocate(current_pos, stat=ierr)
1432     if(ierr /= 0) then
1433        write(stderr,*) myname_,':: deallocate(current_pos) failed.  ', &
1434             'stat = ',ierr
1435        if(present(stat)) then
1436           stat=ierr
1437           return
1438        else
1439           call die(myname_,'deallocate(current_pos)',ierr)
1440        endif
1441     endif
1442
1443  endif ! if(myID == root)
1444
1445       ! Now we are in business...we have:  1) an AttrVect laid out
1446       ! in contiguous segments, each segment corresponding to a
1447       ! process, and in the same order dictated by GSMap;
1448       ! 2) a GlobalMap telling us which segment of workV goes to
1449       ! which process.  Thus, we can us GM_scatter_() to achieve
1450       ! our goal.
1451
1452  call GM_scatter_(workV, oV, GMap, root, comm, ierr)
1453  if(ierr /= 0) then
1454     write(stderr,*) myname,':: ERROR in return from GM_scatter_(), ierr=',&
1455          ierr
1456     if(present(stat)) then
1457        stat = ierr
1458        return
1459     else
1460        call die(myname_,'ERROR returning from GM_scatter_()',ierr)
1461     endif
1462  endif
1463
1464       ! Finally, clean up allocated structures:
1465 
1466  if(myID == root) then
1467     call AttrVect_clean(workV)
1468  endif
1469
1470  call GlobalMap_clean(GMap)
1471
1472  deallocate(lns, stat=ierr)
1473  if(ierr /= 0) then
1474     write(stderr,*) myname_,':: ERROR in deallocate(lns), ierr=',ierr
1475     if(present(stat)) then
1476        stat=ierr
1477        return
1478     else
1479        call die(myname_,'deallocate(lns)',ierr)
1480     endif
1481  endif
1482
1483 end subroutine GSM_scatter_
1484
1485!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1486!    Math and Computer Science Division, Argonne National Laboratory   !
1487!BOP -------------------------------------------------------------------
1488!
1489! !IROUTINE: bcast_ - Broadcast an AttrVect
1490!
1491! !DESCRIPTION:  This routine takes an {\tt AttrVect} argument {\tt aV}
1492! (at input, valid on the root only), and broadcasts it to all the
1493! processes associated with the communicator handle {\tt comm}.  The
1494! success (failure) of this routine is reported in the zero (non-zero)
1495! value of the output argument {\tt stat}.
1496!
1497! {\bf N.B.}:  The output (on non-root processes) {\tt AttrVect} argument
1498! {\tt aV} represents dynamically allocated memory.  When it is no longer
1499! needed, it should be deallocated by invoking {\tt AttrVect\_clean()}
1500! (see the module {\tt m\_AttrVect} for details).
1501!
1502! !INTERFACE:
1503
1504 subroutine bcast_(aV, root, comm, stat)
1505!
1506! !USES:
1507!
1508      use m_stdio
1509      use m_die
1510      use m_mpif90
1511      use m_String, only : String,bcast,char,String_clean
1512      use m_String, only : String_bcast => bcast
1513      use m_List, only : List_get => get
1514      use m_AttrVect, only : AttrVect
1515      use m_AttrVect, only : AttrVect_init => init
1516      use m_AttrVect, only : AttrVect_zero => zero
1517      use m_AttrVect, only : AttrVect_lsize => lsize
1518      use m_AttrVect, only : AttrVect_nIAttr => nIAttr
1519      use m_AttrVect, only : AttrVect_nRAttr => nRAttr
1520
1521      implicit none
1522
1523! !INPUT PARAMETERS:
1524!
1525      integer,                  intent(in)    :: root
1526      integer,                  intent(in)    :: comm
1527
1528! !INPUT/OUTPUT PARAMETERS:
1529!
1530      type(AttrVect),           intent(inout) :: aV ! (IN) on the root,
1531                                                    ! (OUT) elsewhere
1532
1533! !OUTPUT PARAMETERS:
1534!
1535      integer,        optional, intent(out)   :: stat
1536
1537! !REVISION HISTORY:
1538! 27Apr98 - Jing Guo <guo@thunder> - initial prototype/prologue/code
1539! 27Oct00 - J.W. Larson <larson@mcs.anl.gov> - relocated from
1540!           m_AttrVect
1541!  9May01 - J.W. Larson <larson@mcs.anl.gov> - tidied up prologue
1542! 18May01 - R.L. Jacob <jacob@mcs.anl.gov> - use MP_Type function
1543!           to determine type for bcast
1544! 19Dec01 - E.T. Ong <eong@mcs.anl.gov> - adjusted for case of AV with
1545!           only integer or real attribute
1546!EOP ___________________________________________________________________
1547
1548  character(len=*),parameter :: myname_=myname//'::bcast_'
1549  type(String) :: iLStr,rLStr
1550  integer :: nIA, nRA, lsize
1551  integer :: myID
1552  integer :: ier
1553  integer :: mp_Type_aV
1554
1555  if(present(stat)) stat=0
1556
1557  call MP_comm_rank(comm,myID,ier)
1558  if(ier /= 0) then
1559    call MP_perr_die(myname_,'MP_comm_rank()',ier)
1560  endif
1561
1562       ! Broadcaast to all PEs
1563
1564  if(myID == root) then
1565     nIA = AttrVect_nIAttr(aV)
1566     nRA = AttrVect_nRAttr(aV)
1567     lsize = AttrVect_lsize(aV)
1568  endif
1569
1570  call MPI_bcast(nIA,1,MP_INTEGER,root,comm,ier)
1571  if(ier /= 0) then
1572    call MP_perr_die(myname_,'MPI_bcast(nIA)',ier)
1573  endif
1574
1575  call MPI_bcast(nRA,1,MP_INTEGER,root,comm,ier)
1576  if(ier /= 0) then
1577    call MP_perr_die(myname_,'MPI_bcast(nRA)',ier)
1578  endif
1579
1580  call MPI_bcast(lsize,1,MP_INTEGER,root,comm,ier)
1581  if(ier /= 0) then
1582    call MP_perr_die(myname_,'MPI_bcast(lsize)',ier)
1583  endif
1584
1585        ! Convert the two Lists to two Strings
1586
1587  if(nIA>0) then
1588
1589     if(myID == root) call List_get(iLStr,aV%iList)
1590
1591     call String_bcast(iLStr,root,comm,stat=ier)        ! bcast.String()
1592
1593     if(ier /= 0) then
1594        write(stderr,*) myname_,'bcast.String(iLstr), ier=',ier
1595        if(present(stat)) then
1596           stat=ier
1597           return
1598        else
1599           call die(myname_,'String_bcast(iLStr) failed',ier)
1600        endif
1601     endif ! if(ier /= 0)...
1602
1603  endif ! if(nIA > 0)...
1604
1605
1606  if(nRA>0) then
1607
1608     if(myID == root) call List_get(rLStr,aV%rList)
1609
1610     call String_bcast(rLStr,root,comm,stat=ier)        ! bcast.String()
1611     if(ier /= 0) then
1612        write(stderr,*) myname_,'bcast.String(iLstr), ier=',ier
1613        if(present(stat)) then
1614           stat=ier
1615           return
1616        else
1617           call die(myname_,'String_bcast(iLStr) failed',ier)
1618        endif
1619     endif ! if(ier /= 0)...
1620
1621  endif ! if(nRA > 0)...
1622
1623  if(myID /= root) then
1624     
1625     if( (nIA>0) .and. (nRA>0) ) then
1626        call AttrVect_init(aV,iList=char(iLStr),rList=char(rLStr), &
1627                           lsize=lsize)
1628     endif
1629
1630     if( (nIA>0) .and. (nRA<=0) ) then
1631        call AttrVect_init(aV,iList=char(iLStr),lsize=lsize)
1632     endif
1633
1634     if( (nIA<=0) .and. (nRA>0) ) then
1635        call AttrVect_init(aV,rList=char(rLStr),lsize=lsize)
1636     endif
1637
1638     if( (nIA<=0) .and. (nRA<=0) ) then
1639        write(stderr,*) myname_,':: Nonpositive numbers of both ',&
1640             'real AND integer attributes.  nIA =',nIA,' nRA=',nRA
1641        if(present(stat)) then
1642           stat = -1
1643           return
1644        else
1645           call die(myname_,'AV has not been initialized',-1)
1646        endif
1647     endif ! if((nIA<= 0) .and. (nRA<=0))...
1648
1649     call AttrVect_zero(aV)
1650
1651
1652  endif ! if(myID /= root)...
1653
1654  if(nIA > 0) then
1655
1656     mp_Type_aV=MP_Type(av%iAttr)
1657     call MPI_bcast(aV%iAttr,nIA*lsize,mp_Type_aV,root,comm,ier)
1658     if(ier /= 0) then
1659        call MP_perr_die(myname_,'MPI_bcast(iAttr) failed.',ier)
1660     endif
1661
1662     call String_clean(iLStr)
1663
1664  endif
1665
1666  if(nRA > 0) then
1667
1668     mp_Type_aV=MP_Type(av%rAttr)
1669     call MPI_bcast(aV%rAttr,nRA*lsize,mp_Type_aV,root,comm,ier)
1670     if(ier /= 0) then
1671        call MP_perr_die(myname_,'MPI_bcast(rAttr) failed.',ier)
1672     endif
1673
1674     call String_clean(rLStr)
1675
1676  endif
1677
1678 end subroutine bcast_
1679
1680 end module m_AttrVectComms
1681
1682
1683
Note: See TracBrowser for help on using the repository browser.