source: CONFIG/UNIFORM/v7/IPSLCM7/SOURCES/OASIS3-MCT/m_AttrVectComms.F90 @ 6554

Last change on this file since 6554 was 6346, checked in by cetlod, 16 months ago

IPSLCM7_WORK : 1st step of coupling DYNAMICO-LMDZ-ORCHIDEE-NEMO4.2 + Oasis3-mct.5.0

File size: 51.8 KB
Line 
1!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2!    Math and Computer Science Division, Argonne National Laboratory   !
3!-----------------------------------------------------------------------
4! CVS $Id$
5! CVS $Name$
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!CE
694! It has been suggested to comment out the test that triggers the error in
695! lib/mct/mct/m_AttrVectComms.F90, and that works in my case
696! Eric Maisonnave had exactly the same problem at the MetOffice and he did
697! exactly
698!  what you did. There is a ticket opened on the problem #2472 and we work on it
699!
700!  if(GlobalSegMap_haloed(GSMap)) then
701!     ierr = 1
702!     call die(myname_,"Input GlobalSegMap haloed--not allowed",ierr)
703!  endif
704!CE
705       ! Which process am I?
706
707  call MPI_COMM_RANK(comm, myID, ierr)
708
709  if(ierr /= 0) then
710        call MP_perr_die(myname_,':: call MPI_COMM_RANK()',ierr)
711  endif
712       ! How many processes are there on this communicator?
713
714  call MPI_COMM_SIZE(comm, NumProcs, ierr)
715
716  if(ierr /= 0) then
717     call MP_perr_die(myname_,':: call MPI_COMM_SIZE()',ierr)
718  endif
719
720       ! Processor Check: Do the processors on GSMap match those in comm?
721
722  if(MAXVAL(GSMap%pe_loc) > (NumProcs-1)) then
723     stat=2
724     write(stderr,*) myname_, &
725       ":: Procs in GSMap%pe_loc do not match procs in communicator ", &
726       NumProcs-1, MAXVAL(GSMap%pe_loc)
727     call die(myname_, &
728          "Procs in GSMap%pe_loc do not match procs in communicator",stat)
729  endif
730
731  if(myID == root) then
732
733       ! Allocate a precursor to a GlobalMap accordingly...
734
735     allocate(lns(0:NumProcs-1), stat=ierr)
736
737       ! And Load it...
738
739     lns(:)=0
740     do iseg=1,GSMap%ngseg
741        n = GSMap%pe_loc(iseg)
742        lns(n) = lns(n) + GSMap%length(iseg)
743     end do
744
745  else
746
747     allocate(lns(0)) ! This conforms to F90 standard for shaped arguments.
748
749  endif ! if(myID == root)
750
751       ! Determine the component id of GSMap:
752
753  comp_id = GlobalSegMap_comp_id(GSMap)
754
755       ! Create working GlobalMap workGMap (used for the gather):
756
757  call GlobalMap_init(workGMap, comp_id, lns, root, comm)
758
759       ! Gather the Data process-by-process to workV...
760       ! do not include stat argument; bypass an argument check in gm_gather.
761
762  call GM_gather_(iV, workV, workGMap, root, comm, stat)
763
764       ! On the root, initialize oV, and load the contents of
765       !workV into it...
766
767  if(myID == root) then
768
769! bug fix:  gstorage will be bigger than gssize if GSmap is
770! haloed.  But gstorage may be smaller than gsize if GSmap
771! is masked.  So take the maximum.  RLJ
772     gstorage = GlobalSegMap_GlobalStorage(GSMap)
773     gssize = GlobalSegMap_gsize(GSMap)
774     global_storage = MAX(gstorage,gssize)
775
776     call AttrVect_init(oV,iV,global_storage)
777     call AttrVect_zero(oV)
778
779     if (present(rdefault)) then
780       if (AttrVect_nRAttr(oV) > 0) oV%rAttr=rdefault
781     endif
782     if (present(idefault)) then
783       if (AttrVect_nIAttr(oV) > 0) oV%iAttr=idefault
784     endif
785
786       ! On the root, allocate current position index for
787       ! each process chunk:
788
789     allocate(current_pos(0:NumProcs-1), stat=ierr)
790
791     if(ierr /= 0) then
792        write(stderr,*) myname_,':: allocate(current_pos(..) failed,', &
793             'stat = ',ierr
794        if(present(stat)) then
795           stat=ierr
796        else
797           call die(myname_,'allocate(current_pos(..) failed.' )
798        endif
799     endif
800
801       ! Initialize current_pos(:) using GMap%displs(:)
802
803     do n=0,NumProcs-1
804        current_pos(n) = workGMap%displs(n) + 1
805     end do
806
807       ! Load each segment of iV into its appropriate segment
808       ! of workV:
809
810     ngseg = GlobalSegMap_ngseg(GSMap)
811
812     do n=1,ngseg
813
814       ! Determine which process owns segment n:
815
816        pe = GSMap%pe_loc(n)
817
818       ! Input map (lower/upper indicess) of segment of iV:
819
820        ilb = current_pos(pe)
821        iub = current_pos(pe) + GSMap%length(n) - 1
822
823       ! Output map of (lower/upper indicess) segment of workV:
824
825        olb = GSMap%start(n)
826        oub = GSMap%start(n) + GSMap%length(n) - 1
827
828       ! Increment current_pos(n) for next time:
829
830        current_pos(pe) = current_pos(pe) + GSMap%length(n)
831
832       ! Now we are equipped to do the copy:
833
834        do m=1,AttrVect_nIAttr(iV)
835           oV%iAttr(m,olb:oub) = workV%iAttr(m,ilb:iub)
836        end do
837
838        do m=1,AttrVect_nRAttr(iV)
839           oV%rAttr(m,olb:oub) = workV%rAttr(m,ilb:iub)
840        end do
841
842     end do ! do n=1,ngseg
843
844       ! Clean up current_pos, which was only allocated on the root
845
846     deallocate(current_pos, stat=ierr)
847     if(ierr /= 0) then
848        write(stderr,*) myname_,'error in deallocate(current_pos), stat=',ierr
849        if(present(stat)) then
850           stat=ierr
851        else
852           call die(myname_)
853        endif
854     endif
855  endif ! if(myID == root)
856
857       ! At this point, we are finished.  The data have been gathered
858       ! to oV
859
860       ! Finally, clean up allocated structures:
861
862  if(myID == root) call AttrVect_clean(workV)
863  call GlobalMap_clean(workGMap)
864
865  deallocate(lns, stat=ierr)
866
867  if(ierr /= 0) then
868    write(stderr,*) myname_,'error in deallocate(lns), stat=',ierr
869    if(present(stat)) then
870       stat=ierr
871    else
872       call die(myname_)
873    endif
874  endif
875
876 end subroutine GSM_gather_
877
878!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
879!    Math and Computer Science Division, Argonne National Laboratory   !
880!BOP -------------------------------------------------------------------
881!
882! !IROUTINE: GM_scatter_ - Scatter an AttrVect Using a GlobalMap
883!
884! !DESCRIPTION:
885! The routine {\tt GM\_scatter\_} takes an input {\tt AttrVect} type
886! {\tt iV} (valid only on the root), and scatters it to a distributed
887! {\tt AttrVect} {\tt oV}.  The input {\tt GlobalMap} argument
888! {\tt GMap} dictates how {\tt iV} is scattered to {\tt oV}.  The
889! success (failure) of this routine is reported in the zero (non-zero)
890! value of the output argument {\tt stat}.
891!
892! {\bf N.B.}:  The output {\tt AttrVect} argument {\tt oV} represents
893! dynamically allocated memory.  When it is no longer needed, it should
894! be deallocated by invoking {\tt AttrVect\_clean()} (see the module
895! {\tt m\_AttrVect} for more details).
896!
897! !INTERFACE:
898
899 subroutine GM_scatter_(iV, oV, GMap, root, comm, stat)
900!
901! !USES:
902!
903      use m_stdio
904      use m_die
905      use m_mpif90
906      use m_realkinds, only : FP
907
908      use m_List, only : List
909      use m_List, only : List_copy => copy
910      use m_List, only : List_bcast => bcast
911      use m_List, only : List_clean => clean
912      use m_List, only : List_nullify => nullify
913      use m_List, only : List_nitem => nitem
914
915      use m_GlobalMap, only : GlobalMap
916      use m_GlobalMap, only : GlobalMap_lsize => lsize
917      use m_GlobalMap, only : GlobalMap_gsize => gsize
918
919      use m_AttrVect, only : AttrVect
920      use m_AttrVect, only : AttrVect_init => init
921      use m_AttrVect, only : AttrVect_zero => zero
922      use m_AttrVect, only : AttrVect_lsize => lsize
923      use m_AttrVect, only : AttrVect_nIAttr => nIAttr
924      use m_AttrVect, only : AttrVect_nRAttr => nRAttr
925      use m_AttrVect, only : AttrVect_clean => clean
926
927      implicit none
928
929! !INPUT PARAMETERS:
930!
931      type(AttrVect),           intent(in)  :: iV
932      type(GlobalMap),          intent(in)  :: GMap
933      integer,                  intent(in)  :: root
934      integer,                  intent(in)  :: comm
935
936! !OUTPUT PARAMETERS:
937!
938      type(AttrVect),           intent(out) :: oV
939      integer,        optional, intent(out) :: stat
940
941! !REVISION HISTORY:
942! 21Apr98 - Jing Guo <guo@thunder> - initial prototype/prolog/code
943! 27Oct00 - J.W. Larson <larson@mcs.anl.gov> - relocated from
944!           m_AttrVect
945! 15Jan01 - J.W. Larson <larson@mcs.anl.gov> - renamed  GM_scatter_
946!  8Feb01 - J.W. Larson <larson@mcs.anl.gov> - add logic to prevent
947!           empty calls (i.e. no data in buffer) to MPI_SCATTERV()
948! 27Apr01 - R.L. Jacob <jacob@mcs.anl.gov> - small bug fix to
949!           integer attribute scatter
950!  9May01 - J.W. Larson <larson@mcs.anl.gov> - Re-vamped comms model
951!           to reflect MPI comms model for the scatter.  Tidied up
952!           the prologue, too.
953! 18May01 - R.L. Jacob <jacob@mcs.anl.gov> - use MP_Type function
954!           to determine type for mpi_scatterv
955!  8Aug01 - E.T. Ong <eong@mcs.anl.gov> - replace list assignment(=)
956!           with list copy to avoid compiler errors in pgf90.
957! 13Dec01 - E.T. Ong <eong@mcs.anl.gov> - allow scatter with an
958!           AttrVect containing only an iList or rList.
959!EOP ___________________________________________________________________
960
961  character(len=*),parameter :: myname_=myname//'::GM_scatter_'
962  integer :: nIA,nRA,niV,noV,ier
963  integer :: myID
964  integer :: mp_type_Av
965  type(List) :: iList, rList
966  type(AttrVect) :: nonRootAV
967
968  if(present(stat)) stat=0
969
970  call MP_comm_rank(comm,myID,ier)
971  if(ier /= 0) then
972    call MP_perr_die(myname_,'MP_comm_rank()',ier)
973  endif
974
975        ! Verify the input: a _gathered_ vector
976
977  if(myID == root) then
978
979     niV = GlobalMap_gsize(GMap)  ! the _gathered_ local size
980     noV = AttrVect_lsize(iV)     ! the length of the input AttrVect iV
981
982     if(niV /= noV) then
983        write(stderr,'(2a,i5,a,i8,a,i8)') myname_,      &
984             ': myID = ',myID,'.  Invalid input on root, gsize(GMap) =',&
985             niV,', lsize(iV) =',noV
986        if(present(stat)) then
987           stat=-1
988        else
989           call die(myname_)
990        endif
991     endif
992
993  endif
994
995        ! On the root, read the integer and real attribute
996        ! lists off of iV.
997
998  call List_nullify(iList)
999  call List_nullify(rList)
1000
1001  if(myID == root) then
1002
1003     ! Count the number of real and integer attributes
1004
1005     nIA = AttrVect_nIAttr(iV)  ! number of INTEGER attributes
1006     nRA = AttrVect_nRAttr(iV)  ! number of REAL attributes
1007
1008     if(nIA > 0) then
1009        call List_copy(iList,iV%iList)
1010     endif
1011
1012     if(nRA > 0) then
1013        call List_copy(rList,iV%rList)
1014     endif
1015
1016  endif
1017
1018        ! From the root, broadcast iList and rList
1019
1020  call MPI_BCAST(nIA,1,MP_INTEGER,root,comm,ier)
1021  if(ier /= 0) call MP_perr(myname_,'MPI_BCAST(nIA)',ier)
1022
1023  call MPI_BCAST(nRA,1,MP_INTEGER,root,comm,ier)
1024  if(ier /= 0) call MP_perr(myname_,'MPI_BCAST(nRA)',ier)
1025
1026  if(nIA>0) call List_bcast(iList, root, comm)
1027  if(nRA>0) call List_bcast(rList, root, comm)
1028
1029  noV = GlobalMap_lsize(GMap) ! the _scatterd_ local size
1030
1031        ! On all processes, use List data and noV to initialize oV
1032
1033  call AttrVect_init(oV, iList, rList, noV)
1034  call AttrVect_zero(oV)
1035
1036        ! Initialize a dummy AttrVect for non-root MPI calls
1037
1038  if(myID/=root) then
1039    call AttrVect_init(nonRootAV,oV,1)
1040    call AttrVect_zero(nonRootAV)
1041  endif
1042
1043
1044  if(nIA > 0) then
1045
1046     if(myID == root) then
1047
1048        call MPI_scatterv(iV%iAttr,GMap%counts*nIA, &
1049             GMap%displs*nIA,MP_INTEGER,oV%iAttr,   &
1050             noV*nIA,MP_INTEGER,root,comm,ier )
1051        if(ier /= 0) then
1052           call MP_perr_die(myname_,'MPI_scatterv(iAttr) on root',ier)
1053        endif
1054
1055     else
1056
1057        call MPI_scatterv(nonRootAV%iAttr,GMap%counts*nIA, &
1058             GMap%displs*nIA,MP_INTEGER,oV%iAttr,          &
1059             noV*nIA,MP_INTEGER,root,comm,ier )
1060        if(ier /= 0) then
1061           call MP_perr_die(myname_,'MPI_scatterv(iAttr) off root',ier)
1062        endif
1063
1064     endif   ! if(myID == root)
1065
1066     call List_clean(iList)
1067
1068  endif   ! if(nIA > 0)
1069
1070  mp_type_Av = MP_Type(1._FP)   ! set mpi type to same as AV%rAttr
1071
1072  if(nRA > 0) then
1073
1074     if(myID == root) then
1075
1076
1077        call MPI_scatterv(iV%rAttr,GMap%counts*nRA, &
1078             GMap%displs*nRA,mp_type_Av,oV%rAttr,   &
1079             noV*nRA,mp_type_Av,root,comm,ier )
1080        if(ier /= 0) then
1081           call MP_perr_die(myname_,'MPI_scatterv(rAttr) on root',ier)
1082        endif
1083
1084     else
1085
1086
1087        call MPI_scatterv(nonRootAV%rAttr,GMap%counts*nRA, &
1088             GMap%displs*nRA,mp_type_Av,oV%rAttr,          &
1089             noV*nRA,mp_type_Av,root,comm,ier )
1090        if(ier /= 0) then
1091           call MP_perr_die(myname_,'MPI_scatterv(rAttr) off root',ier)
1092        endif
1093
1094     endif
1095
1096     call List_clean(rList)
1097
1098  endif
1099
1100  if(myID /= root) then
1101     call AttrVect_clean(nonRootAV,ier)
1102     if(ier /= 0) then
1103        write(stderr,'(2a,i4)') myname_,        &
1104             ':: AttrVect_clean(nonRootAV) failed for non-root &
1105             &process: myID = ', myID
1106        call die(myname_,':: AttrVect_clean failed &
1107             &for nonRootAV off of root',ier)
1108     endif
1109  endif
1110
1111 end subroutine GM_scatter_
1112
1113!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1114!    Math and Computer Science Division, Argonne National Laboratory   !
1115!BOP -------------------------------------------------------------------
1116!
1117! !IROUTINE: GSM_scatter_ - Scatter an AttrVect using a GlobalSegMap
1118!
1119! !DESCRIPTION:
1120! The routine {\tt GSM\_scatter\_} takes an input {\tt AttrVect} type
1121! {\tt iV} (valid only on the root), and scatters it to a distributed
1122! {\tt AttrVect} {\tt oV}.  The input {\tt GlobalSegMap} argument
1123! {\tt GSMap} dictates how {\tt iV} is scattered to {\tt oV}.  The
1124! success (failure) of this routine is reported in the zero (non-zero)
1125! value of the output argument {\tt stat}.
1126!
1127! {\tt GSM\_scatter\_()} converts the problem of scattering data
1128! according to a {\tt GlobalSegMap} into the simpler problem of
1129! scattering data as specified by a {\tt GlobalMap}.  The {\tt GlobalMap}
1130! variable {\tt GMap} is created based on the local storage requirements
1131! for each distributed piece of {\tt iV}.  On the root, a complete
1132! (including halo points) copy of {\tt iV} is stored in
1133! the temporary {\tt AttrVect} variable {\tt workV} (the length of
1134! {\tt workV} is {\tt GlobalSegMap\_GlobalStorage(GSMap)}).  The
1135! variable {\tt workV} is segmented by process, and segments are
1136! copied into it by process, but ordered in the same order the segments
1137! appear in {\tt GSMap}.  Once {\tt workV} is loaded, the data are
1138! scattered to the output {\tt AttrVect} {\tt oV} by a call to the
1139! routine {\tt GM\_scatter\_()} defined in this module, with {\tt workV}
1140! and {\tt GMap} as the input arguments.
1141!
1142! {\bf N.B.:}  This algorithm  assumes that memory access times are much
1143! shorter than message-passing transmission times.
1144!
1145! {\bf N.B.}:  The output {\tt AttrVect} argument {\tt oV} represents
1146! dynamically allocated memory.  When it is no longer needed, it should
1147! be deallocated by invoking {\tt AttrVect\_clean()} (see the module
1148! {\tt m\_AttrVect} for more details).
1149!
1150! !INTERFACE:
1151
1152 subroutine GSM_scatter_(iV, oV, GSMap, root, comm, stat)
1153!
1154! !USES:
1155!
1156! Environment utilities from mpeu:
1157
1158      use m_stdio
1159      use m_die
1160      use m_mpif90
1161
1162      use m_List, only : List_nullify => nullify
1163
1164! GlobalSegMap and associated services:
1165      use m_GlobalSegMap, only : GlobalSegMap
1166      use m_GlobalSegMap, only : GlobalSegMap_comp_id => comp_id
1167      use m_GlobalSegMap, only : GlobalSegMap_ngseg => ngseg
1168      use m_GlobalSegMap, only : GlobalSegMap_lsize => lsize
1169      use m_GlobalSegMap, only : GlobalSegMap_gsize => gsize
1170      use m_GlobalSegMap, only : GlobalSegMap_GlobalStorage => GlobalStorage
1171! AttrVect and associated services:
1172      use m_AttrVect, only : AttrVect
1173      use m_AttrVect, only : AttrVect_init => init
1174      use m_AttrVect, only : AttrVect_zero => zero
1175      use m_AttrVect,  only : AttrVect_lsize => lsize
1176      use m_AttrVect, only : AttrVect_nIAttr => nIAttr
1177      use m_AttrVect, only : AttrVect_nRAttr => nRAttr
1178      use m_AttrVect, only : AttrVect_clean => clean
1179! GlobalMap and associated services:
1180      use m_GlobalMap, only : GlobalMap
1181      use m_GlobalMap, only : GlobalMap_init => init
1182      use m_GlobalMap, only : GlobalMap_clean => clean
1183
1184      implicit none
1185
1186! !INPUT PARAMETERS:
1187!
1188      type(AttrVect),            intent(in)  :: iV
1189      type(GlobalSegMap),        intent(in)  :: GSMap
1190      integer,                   intent(in)  :: root
1191      integer,                   intent(in)  :: comm
1192
1193! !OUTPUT PARAMETERS:
1194!
1195      type(AttrVect),            intent(out) :: oV
1196      integer,        optional,  intent(out) :: stat
1197
1198! !REVISION HISTORY:
1199! 15Jan01 - J.W. Larson <larson@mcs.anl.gov> - API specification.
1200!  8Feb01 - J.W. Larson <larson@mcs.anl.gov> - Initial code.
1201! 25Feb01 - J.W. Larson <larson@mcs.anl.gov> - Bug fix--replaced
1202!           call to GlobalSegMap_lsize with call to the new fcn.
1203!           GlobalSegMap_ProcessStorage().
1204! 26Apr01 - R.L. Jacob <jacob@mcs.anl.gov> - add use statement for
1205!           AttVect_clean
1206! 26Apr01 - J.W. Larson <larson@mcs.anl.gov> - bug fixes--data
1207!           misalignment in use of the GlobalMap to compute the
1208!           memory map into workV, and initialization of workV
1209!           on all processes.
1210!  9May01 - J.W. Larson <larson@mcs.anl.gov> - tidied up prologue
1211! 15May01 - Larson / Jacob <larson@mcs.anl.gov> - stopped initializing
1212!           workV on off-root processes (no longer necessary).
1213! 13Jun01 - J.W. Larson <larson@mcs.anl.gov> - Initialize stat
1214!           (if present).
1215! 20Jun01 - J.W. Larson <larson@mcs.anl.gov> - Fixed a subtle bug
1216!           appearing on AIX regarding the fact workV is uninitial-
1217!           ized on non-root processes.  This is fixed by nullifying
1218!           all the pointers in workV for non-root processes.
1219! 20Aug01 - E.T. Ong <eong@mcs.anl.gov> - Added argument check
1220!           for matching processors in gsmap and comm.
1221! 13Dec01 - E.T. Ong <eong@mcs.anl.gov> - got rid of restriction
1222!           GlobalStorage(GSMap)==AttrVect_lsize(AV) to allow for
1223!           GSMap to be haloed.
1224! 11Aug08 - R. Jacob <jacob@mcs.anl.gov> - remove call to ProcessStorage
1225!           and replace with faster algorithm provided by Pat Worley
1226!EOP ___________________________________________________________________
1227
1228  character(len=*),parameter :: myname_=myname//'::GSM_scatter_'
1229
1230! Temporary workspace AttrVect:
1231  type(AttrVect) :: workV
1232! Component ID and number of segments for GSMap:
1233  integer :: comp_id, ngseg, iseg
1234! Total length of GSMap segments laid end-to-end:
1235  integer :: global_storage
1236! Error Flag
1237  integer :: ierr
1238! Number of processes on communicator, and local rank:
1239  integer :: NumProcs, myID
1240! Total local storage on each pe according to GSMap:
1241  integer, dimension(:), allocatable :: lns
1242! Temporary GlobalMap used to scatter the segmented (by pe) data
1243  type(GlobalMap) :: GMap
1244! Loop counters and temporary indices:
1245  integer :: m, n, ilb, iub, olb, oub, pe
1246! workV segment tracking index array:
1247  integer, dimension(:), allocatable :: current_pos
1248
1249      ! Initialize stat (if present)
1250
1251  if(present(stat)) stat = 0
1252
1253       ! Which process am I?
1254
1255  call MPI_COMM_RANK(comm, myID, ierr)
1256
1257  if(ierr /= 0) then
1258     call MP_perr_die(myname_,'MPI_COMM_RANK',ierr)
1259  endif
1260
1261  if(myID == root) then
1262
1263     if(GSMap%gsize > AttrVect_lsize(iV)) then
1264        write(stderr,'(2a,i5,a,i8,a,i8)') myname_,      &
1265             ': myID = ',myID,'.  Invalid input, GSMap%gsize =',&
1266             GSMap%gsize, ', lsize(iV) =',AttrVect_lsize(iV)
1267        if(present(stat)) then
1268           stat=-1
1269        else
1270           call die(myname_)
1271        endif
1272     endif
1273
1274  endif
1275
1276       ! On the root, initialize a work AttrVect type of the
1277       ! above length, and with the same attribute lists as iV.
1278       ! on other processes, initialize workV only with the
1279       ! attribute information, but no storage.
1280
1281  if(myID == root) then
1282
1283     global_storage = GlobalSegMap_GlobalStorage(GSMap)
1284     call AttrVect_init(workV, iV, global_storage)
1285     call AttrVect_zero(workV)
1286
1287  else
1288       ! nullify workV just to be safe
1289
1290     call List_nullify(workV%iList)
1291     call List_nullify(workV%rList)
1292     nullify(workV%iAttr)
1293     nullify(workV%rAttr)
1294
1295  endif
1296
1297       ! Return to processing on the root to load workV:
1298
1299       ! How many processes are there on this communicator?
1300
1301  call MPI_COMM_SIZE(comm, NumProcs, ierr)
1302
1303  if(ierr /= 0) then
1304     call MP_perr_die(myname_,'MPI_COMM_SIZE',ierr)
1305  endif
1306
1307       ! Processor Check: Do the processors on GSMap match those in comm?
1308
1309  if(MAXVAL(GSMap%pe_loc) > (NumProcs-1)) then
1310     write(stderr,*) myname_, &
1311          ":: Procs in GSMap%pe_loc do not match procs in communicator ", &
1312          NumProcs-1, MAXVAL(GSMap%pe_loc)
1313     if(present(stat)) then
1314        stat=1
1315        return
1316     else
1317        call die(myname_)
1318     endif
1319  endif
1320
1321  if(myID == root) then
1322
1323       ! Allocate a precursor to a GlobalMap accordingly...
1324
1325     allocate(lns(0:NumProcs-1), stat=ierr)
1326     if(ierr /= 0) then
1327        write(stderr,*) myname_,':: allocate(lns...) failed, stat=',ierr
1328        if(present(stat)) then
1329           stat=ierr
1330        else
1331           call die(myname_,'allocate(lns)',ierr)
1332        endif
1333     endif
1334
1335       ! And Load it...
1336
1337     lns(:)=0
1338     do iseg=1,GSMap%ngseg
1339        n = GSMap%pe_loc(iseg)
1340        lns(n) = lns(n) + GSMap%length(iseg)
1341     end do
1342
1343  endif ! if(myID == root)
1344
1345        ! Non-root processes call GlobalMap_init with lns,
1346        ! although this argument is not used in the
1347        ! subroutine. Since it correspond to a dummy shaped array arguments
1348        ! in GlobslMap_init, the Fortran 90 standard dictates that the actual
1349        ! argument must contain complete shape information. Therefore,
1350        ! the array argument must be allocated on all processes.
1351
1352  if(myID /= root) then
1353
1354     allocate(lns(1),stat=ierr)
1355     if(ierr /= 0) then
1356        write(stderr,*) myname_,':: allocate(lns...) failed, stat=',ierr
1357        if(present(stat)) then
1358           stat=ierr
1359           return
1360        else
1361           call die(myname_,'allocate(lns(1))',ierr)
1362        endif
1363     endif
1364
1365  endif ! if(myID /= root)...
1366
1367       ! Create a GlobalMap describing the 1-D decomposition
1368       ! of workV:
1369
1370  comp_id = GlobalSegMap_comp_id(GSMap)
1371
1372  call GlobalMap_init(GMap, comp_id, lns, root, comm)
1373
1374       ! On the root, load workV:
1375
1376  if(myID == root) then
1377
1378       ! On the root, allocate current position index for
1379       ! each process chunk:
1380
1381     allocate(current_pos(0:NumProcs-1), stat=ierr)
1382     if(ierr /= 0) then
1383        write(stderr,*) myname_,':: allocate(current_pos..) failed, stat=', &
1384             ierr
1385        if(present(stat)) then
1386           stat=ierr
1387           return
1388        else
1389           call die(myname_,'allocate(current_pos)',ierr)
1390        endif
1391     endif
1392
1393       ! Initialize current_pos(:) using GMap%displs(:)
1394
1395     do n=0,NumProcs-1
1396        current_pos(n) = GMap%displs(n) + 1
1397     end do
1398
1399       ! Load each segment of iV into its appropriate segment
1400       ! of workV:
1401
1402     ngseg = GlobalSegMap_ngseg(GSMap)
1403
1404     do n=1,ngseg
1405
1406       ! Determine which process owns segment n:
1407
1408        pe = GSMap%pe_loc(n)
1409
1410       ! Input map (lower/upper indicess) of segment of iV:
1411
1412        ilb = GSMap%start(n)
1413        iub = GSMap%start(n) + GSMap%length(n) - 1
1414
1415       ! Output map of (lower/upper indicess) segment of workV:
1416
1417        olb = current_pos(pe)
1418        oub = current_pos(pe) + GSMap%length(n) - 1
1419
1420       ! Increment current_pos(n) for next time:
1421
1422        current_pos(pe) = current_pos(pe) + GSMap%length(n)
1423
1424       ! Now we are equipped to do the copy:
1425
1426        do m=1,AttrVect_nIAttr(iV)
1427           workV%iAttr(m,olb:oub) = iV%iAttr(m,ilb:iub)
1428        end do
1429
1430        do m=1,AttrVect_nRAttr(iV)
1431           workV%rAttr(m,olb:oub) = iV%rAttr(m,ilb:iub)
1432        end do
1433
1434     end do ! do n=1,ngseg
1435
1436       ! Clean up current_pos, which was only allocated on the root
1437
1438     deallocate(current_pos, stat=ierr)
1439     if(ierr /= 0) then
1440        write(stderr,*) myname_,':: deallocate(current_pos) failed.  ', &
1441             'stat = ',ierr
1442        if(present(stat)) then
1443           stat=ierr
1444           return
1445        else
1446           call die(myname_,'deallocate(current_pos)',ierr)
1447        endif
1448     endif
1449
1450  endif ! if(myID == root)
1451
1452       ! Now we are in business...we have:  1) an AttrVect laid out
1453       ! in contiguous segments, each segment corresponding to a
1454       ! process, and in the same order dictated by GSMap;
1455       ! 2) a GlobalMap telling us which segment of workV goes to
1456       ! which process.  Thus, we can us GM_scatter_() to achieve
1457       ! our goal.
1458
1459  call GM_scatter_(workV, oV, GMap, root, comm, ierr)
1460  if(ierr /= 0) then
1461     write(stderr,*) myname,':: ERROR in return from GM_scatter_(), ierr=',&
1462          ierr
1463     if(present(stat)) then
1464        stat = ierr
1465        return
1466     else
1467        call die(myname_,'ERROR returning from GM_scatter_()',ierr)
1468     endif
1469  endif
1470
1471       ! Finally, clean up allocated structures:
1472
1473  if(myID == root) then
1474     call AttrVect_clean(workV)
1475  endif
1476
1477  call GlobalMap_clean(GMap)
1478
1479  deallocate(lns, stat=ierr)
1480  if(ierr /= 0) then
1481     write(stderr,*) myname_,':: ERROR in deallocate(lns), ierr=',ierr
1482     if(present(stat)) then
1483        stat=ierr
1484        return
1485     else
1486        call die(myname_,'deallocate(lns)',ierr)
1487     endif
1488  endif
1489
1490 end subroutine GSM_scatter_
1491
1492!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1493!    Math and Computer Science Division, Argonne National Laboratory   !
1494!BOP -------------------------------------------------------------------
1495!
1496! !IROUTINE: bcast_ - Broadcast an AttrVect
1497!
1498! !DESCRIPTION:  This routine takes an {\tt AttrVect} argument {\tt aV}
1499! (at input, valid on the root only), and broadcasts it to all the
1500! processes associated with the communicator handle {\tt comm}.  The
1501! success (failure) of this routine is reported in the zero (non-zero)
1502! value of the output argument {\tt stat}.
1503!
1504! {\bf N.B.}:  The output (on non-root processes) {\tt AttrVect} argument
1505! {\tt aV} represents dynamically allocated memory.  When it is no longer
1506! needed, it should be deallocated by invoking {\tt AttrVect\_clean()}
1507! (see the module {\tt m\_AttrVect} for details).
1508!
1509! !INTERFACE:
1510
1511 subroutine bcast_(aV, root, comm, stat)
1512!
1513! !USES:
1514!
1515      use m_stdio
1516      use m_die
1517      use m_mpif90
1518      use m_String, only : String,bcast,char,String_clean
1519      use m_String, only : String_bcast => bcast
1520      use m_List, only : List_get => get
1521      use m_AttrVect, only : AttrVect
1522      use m_AttrVect, only : AttrVect_init => init
1523      use m_AttrVect, only : AttrVect_zero => zero
1524      use m_AttrVect, only : AttrVect_lsize => lsize
1525      use m_AttrVect, only : AttrVect_nIAttr => nIAttr
1526      use m_AttrVect, only : AttrVect_nRAttr => nRAttr
1527
1528      implicit none
1529
1530! !INPUT PARAMETERS:
1531!
1532      integer,                  intent(in)    :: root
1533      integer,                  intent(in)    :: comm
1534
1535! !INPUT/OUTPUT PARAMETERS:
1536!
1537      type(AttrVect),           intent(inout) :: aV ! (IN) on the root,
1538                                                    ! (OUT) elsewhere
1539
1540! !OUTPUT PARAMETERS:
1541!
1542      integer,        optional, intent(out)   :: stat
1543
1544! !REVISION HISTORY:
1545! 27Apr98 - Jing Guo <guo@thunder> - initial prototype/prologue/code
1546! 27Oct00 - J.W. Larson <larson@mcs.anl.gov> - relocated from
1547!           m_AttrVect
1548!  9May01 - J.W. Larson <larson@mcs.anl.gov> - tidied up prologue
1549! 18May01 - R.L. Jacob <jacob@mcs.anl.gov> - use MP_Type function
1550!           to determine type for bcast
1551! 19Dec01 - E.T. Ong <eong@mcs.anl.gov> - adjusted for case of AV with
1552!           only integer or real attribute
1553!EOP ___________________________________________________________________
1554
1555  character(len=*),parameter :: myname_=myname//'::bcast_'
1556  type(String) :: iLStr,rLStr
1557  integer :: nIA, nRA, lsize
1558  integer :: myID
1559  integer :: ier
1560  integer :: mp_Type_aV
1561
1562  if(present(stat)) stat=0
1563
1564  call MP_comm_rank(comm,myID,ier)
1565  if(ier /= 0) then
1566    call MP_perr_die(myname_,'MP_comm_rank()',ier)
1567  endif
1568
1569       ! Broadcaast to all PEs
1570
1571  if(myID == root) then
1572     nIA = AttrVect_nIAttr(aV)
1573     nRA = AttrVect_nRAttr(aV)
1574     lsize = AttrVect_lsize(aV)
1575  endif
1576
1577  call MPI_bcast(nIA,1,MP_INTEGER,root,comm,ier)
1578  if(ier /= 0) then
1579    call MP_perr_die(myname_,'MPI_bcast(nIA)',ier)
1580  endif
1581
1582  call MPI_bcast(nRA,1,MP_INTEGER,root,comm,ier)
1583  if(ier /= 0) then
1584    call MP_perr_die(myname_,'MPI_bcast(nRA)',ier)
1585  endif
1586
1587  call MPI_bcast(lsize,1,MP_INTEGER,root,comm,ier)
1588  if(ier /= 0) then
1589    call MP_perr_die(myname_,'MPI_bcast(lsize)',ier)
1590  endif
1591
1592        ! Convert the two Lists to two Strings
1593
1594  if(nIA>0) then
1595
1596     if(myID == root) call List_get(iLStr,aV%iList)
1597
1598     call String_bcast(iLStr,root,comm,stat=ier)        ! bcast.String()
1599
1600     if(ier /= 0) then
1601        write(stderr,*) myname_,'bcast.String(iLstr), ier=',ier
1602        if(present(stat)) then
1603           stat=ier
1604           return
1605        else
1606           call die(myname_,'String_bcast(iLStr) failed',ier)
1607        endif
1608     endif ! if(ier /= 0)...
1609
1610  endif ! if(nIA > 0)...
1611
1612
1613  if(nRA>0) then
1614
1615     if(myID == root) call List_get(rLStr,aV%rList)
1616
1617     call String_bcast(rLStr,root,comm,stat=ier)        ! bcast.String()
1618     if(ier /= 0) then
1619        write(stderr,*) myname_,'bcast.String(iLstr), ier=',ier
1620        if(present(stat)) then
1621           stat=ier
1622           return
1623        else
1624           call die(myname_,'String_bcast(iLStr) failed',ier)
1625        endif
1626     endif ! if(ier /= 0)...
1627
1628  endif ! if(nRA > 0)...
1629
1630  if(myID /= root) then
1631
1632     if( (nIA>0) .and. (nRA>0) ) then
1633        call AttrVect_init(aV,iList=char(iLStr),rList=char(rLStr), &
1634                           lsize=lsize)
1635     endif
1636
1637     if( (nIA>0) .and. (nRA<=0) ) then
1638        call AttrVect_init(aV,iList=char(iLStr),lsize=lsize)
1639     endif
1640
1641     if( (nIA<=0) .and. (nRA>0) ) then
1642        call AttrVect_init(aV,rList=char(rLStr),lsize=lsize)
1643     endif
1644
1645     if( (nIA<=0) .and. (nRA<=0) ) then
1646        write(stderr,*) myname_,':: Nonpositive numbers of both ',&
1647             'real AND integer attributes.  nIA =',nIA,' nRA=',nRA
1648        if(present(stat)) then
1649           stat = -1
1650           return
1651        else
1652           call die(myname_,'AV has not been initialized',-1)
1653        endif
1654     endif ! if((nIA<= 0) .and. (nRA<=0))...
1655
1656     call AttrVect_zero(aV)
1657
1658
1659  endif ! if(myID /= root)...
1660
1661  if(nIA > 0) then
1662
1663     mp_Type_aV=MP_Type(av%iAttr)
1664     call MPI_bcast(aV%iAttr,nIA*lsize,mp_Type_aV,root,comm,ier)
1665     if(ier /= 0) then
1666        call MP_perr_die(myname_,'MPI_bcast(iAttr) failed.',ier)
1667     endif
1668
1669     call String_clean(iLStr)
1670
1671  endif
1672
1673  if(nRA > 0) then
1674
1675     mp_Type_aV=MP_Type(av%rAttr)
1676     call MPI_bcast(aV%rAttr,nRA*lsize,mp_Type_aV,root,comm,ier)
1677     if(ier /= 0) then
1678        call MP_perr_die(myname_,'MPI_bcast(rAttr) failed.',ier)
1679     endif
1680
1681     call String_clean(rLStr)
1682
1683  endif
1684
1685 end subroutine bcast_
1686
1687 end module m_AttrVectComms
1688
1689
1690
Note: See TracBrowser for help on using the repository browser.