source: CPL/oasis3-mct/branches/OASIS3-MCT_2.0_branch/lib/mct/mct/m_GeneralGridComms.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: 48.6 KB
Line 
1!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2!    Math and Computer Science Division, Argonne National Laboratory   !
3!-----------------------------------------------------------------------
4! CVS m_GeneralGridComms.F90,v 1.23 2004-04-21 22:16:33 jacob Exp
5! CVS MCT_2_8_0
6!BOP -------------------------------------------------------------------
7!
8! !MODULE: m_GeneralGridComms - Communications for the GeneralGrid type.
9!
10! !DESCRIPTION:
11!
12! In this module, we define communications methods specific to the
13! {\tt GeneralGrid} class (see the module {\tt m\_GeneralGrid} for more
14! information about this class and its methods).
15!
16! !INTERFACE:
17 module m_GeneralGridComms
18!
19! !USES:
20!
21      use m_GeneralGrid ! GeneralGrid class and its methods
22
23
24      implicit none
25
26      private   ! except
27
28      public :: gather          ! gather all local vectors to the root
29      public :: scatter         ! scatter from the root to all PEs
30      public :: bcast           ! bcast from root to all PEs
31      public :: send            ! Blocking SEND
32      public :: recv            ! Blocking RECEIVE
33
34    interface gather ; module procedure &
35              GM_gather_, &
36              GSM_gather_ 
37    end interface
38    interface scatter ; module procedure &
39              GM_scatter_, &
40              GSM_scatter_ 
41    end interface
42    interface bcast ; module procedure bcast_ ; end interface
43    interface send  ; module procedure send_  ; end interface
44    interface recv  ; module procedure recv_  ; end interface
45
46! !REVISION HISTORY:
47!       27Apr01 - J.W. Larson <larson@mcs.anl.gov> - Initial module/APIs
48!       07Jun01 - J.W. Larson <larson@mcs.anl.gov> - Added point-to-point
49!       27Mar02 - J.W. Larson <larson@mcs.anl.gov> - Overhaul of error
50!                 handling calls throughout this module.
51!       05Aug02 - E. Ong <eong@mcs.anl.gov> - Added buffer association
52!                 error checks to avoid making bad MPI calls
53!EOP ___________________________________________________________________
54
55  character(len=*),parameter :: myname='MCT::m_GeneralGridComms'
56
57 contains
58
59!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
60!    Math and Computer Science Division, Argonne National Laboratory   !
61!BOP -------------------------------------------------------------------
62!
63! !IROUTINE: send_ - Point-to-point blocking send for the GeneralGrid.
64!
65! !DESCRIPTION:  The point-to-point send routine {\tt send\_()} sends
66! the input {\tt GeneralGrid} argument {\tt iGGrid} to component
67! {\tt comp\_id}. 
68! The message is identified by the tag defined by the {\tt INTEGER}
69! argument {\tt TagBase}.  The value of {\tt TagBase} must match the
70! value used in the call to {\tt recv\_()} on process {\tt dest}.  The
71! success (failure) of this operation corresponds to a zero (nonzero)
72! value for the output {\tt INTEGER} flag {\tt status}.
73! The argument will be sent to the local root of the component.
74!
75! {\bf N.B.}:  One must avoid assigning elsewhere the MPI tag values
76! between {\tt TagBase} and {\tt TagBase+20}, inclusive.  This is
77! because {\tt send\_()} performs one send operation set up the header
78! transfer, up to five {\tt List\_send} operations (two {\tt MPI\_SEND}
79! calls in each), two send operations to transfer {\tt iGGrid\%descend(:)},
80! and finally the send of the {\tt AttrVect} component {\tt iGGrid\%data}
81! (which comprises eight {\tt MPI\_SEND} operations).
82!
83! !INTERFACE:
84
85 subroutine send_(iGGrid, comp_id, TagBase, status)
86
87!
88! !USES:
89!
90      use m_stdio
91      use m_die
92      use m_mpif90
93
94      use m_GeneralGrid, only : GeneralGrid
95      use m_GeneralGrid, only : GeneralGrid_init => init
96      use m_GeneralGrid, only : GeneralGrid_lsize => lsize
97
98      use m_MCTWorld, only : ComponentToWorldRank
99      use m_MCTWorld, only : ThisMCTWorld
100
101      use m_AttrVectComms,only : AttrVect_send => send
102
103      use m_List, only : List_send => send
104      use m_List, only : List_allocated => allocated
105
106      implicit none
107
108! !INPUT PARAMETERS:
109!
110      type(GeneralGrid), intent(in) :: iGGrid
111      integer,           intent(in) :: comp_id
112      integer,           intent(in) :: TagBase
113
114! !OUTPUT PARAMETERS:
115!
116      integer, optional, intent(out) :: status
117
118! !REVISION HISTORY:
119!       04Jun01 - J.W. Larson <larson@mcs.anl.gov> - API Specification.
120!       07Jun01 - J.W. Larson <larson@mcs.anl.gov> - Initial version.
121!       10Jun01 - J.W. Larson <larson@mcs.anl.gov> - Bug fixes--now works.
122!       11Jun01 - R. Jacob <jacob@mcs.anl.gov> use component id as input
123!                 argument.
124!       13Jun01 - J.W. Larson <larson@mcs.anl.gov> - Initialize status
125!                 (if present).
126!       15Feb02 - J.W. Larson <larson@mcs.anl.gov> - Made input argument
127!                 comm optional.
128!       13Jun02 - J.W. Larson <larson@mcs.anl.gov> - Removed the argument
129!                 comm.  This routine is now explicitly for intercomponent
130!                 communications only.
131!EOP ___________________________________________________________________
132
133  character(len=*),parameter :: myname_=myname//'::send_'
134
135  integer :: ierr
136  integer :: dest
137  logical :: HeaderAssoc(6)
138
139      ! Initialize status (if present)
140
141  if(present(status)) status = 0
142
143  dest = ComponentToWorldRank(0, comp_id, ThisMCTWorld)
144
145      ! Step 1. Check elements of the GeneralGrid header to see
146      ! which components of it are allocated.  Load the results
147      ! into HeaderAssoc(:), and send it to process dest.
148
149  HeaderAssoc(1) = List_allocated(iGGrid%coordinate_list)
150  HeaderAssoc(2) = List_allocated(iGGrid%coordinate_sort_order)
151  HeaderAssoc(3) = associated(iGGrid%descend)
152  HeaderAssoc(4) = List_allocated(iGGrid%weight_list)
153  HeaderAssoc(5) = List_allocated(iGGrid%other_list)
154  HeaderAssoc(6) = List_allocated(iGGrid%index_list)
155
156  call MPI_SEND(HeaderAssoc, 6, MP_LOGICAL, dest, TagBase, ThisMCTWorld%MCT_comm, ierr)
157  if(ierr /= 0) then
158     call MP_perr_die(myname_,':: MPI_SEND(HeaderAssoc...',ierr)
159  endif
160
161       ! Step 2.  If iGGrid%coordinate_list is defined, send it.
162
163  if(HeaderAssoc(1)) then
164    call List_send(iGGrid%coordinate_list, dest, TagBase+1, ThisMCTWorld%MCT_comm, ierr)
165    if(ierr /= 0) then
166       write(stderr,*) myname_,':: call List_send(iGGrid%coordinate_list...', &
167            'Error flag ierr = ',ierr
168       if(present(status)) then
169          status = ierr
170          return
171       else
172          call die(myname_,':: call List_send(iGGrid%coordinate_list...',ierr)
173       endif
174    endif
175  else  ! This constitutes an error, as a GeneralGrid must have coordinates
176
177     if(present(status)) then
178        write(stderr,*) myname_,':: Error.  GeneralGrid%coordinate_list undefined.'
179        status = -1
180        return
181     else
182        call die(myname_,'::  Error.  GeneralGrid%coordinate_list undefined.',-1)
183     endif
184
185  endif ! if(HeaderAssoc(1))...
186
187       ! Step 3.  If iGGrid%coordinate_sort_order is defined, send it.
188
189  if(HeaderAssoc(2)) then
190    call List_send(iGGrid%coordinate_sort_order, dest, TagBase+3, ThisMCTWorld%MCT_comm, ierr)
191    if(ierr /= 0) then
192       if(present(status)) then
193          write(stderr,*) myname_,':: call List_send(iGGrid%coordinate_sort_order...'
194          status = ierr
195          return
196       else
197          call die(myname_,':: call List_send(iGGrid%coordinate_sort_order...',ierr)
198       endif
199    endif
200
201  endif ! if(HeaderAssoc(2))...
202
203       ! Step 4.  If iGGrid%descend is allocated, determine its size,
204       ! send this size, and then send the elements of iGGrid%descend.
205     
206  if(HeaderAssoc(3)) then
207
208     if(size(iGGrid%descend)<=0) call die(myname_,'size(iGGrid%descend)<=0')
209
210     call MPI_SEND(size(iGGrid%descend), 1, MP_type(size(iGGrid%descend)), &
211                   dest, TagBase+5, ThisMCTWorld%MCT_comm, ierr)
212     if(ierr /= 0) then
213        call MP_perr_die(myname_,':: call MPI_SEND(size(iGGrid%descend)...',ierr)
214     endif
215
216     call MPI_SEND(iGGrid%descend, size(iGGrid%descend), MP_type(iGGrid%descend(1)), &
217                   dest, TagBase+6, ThisMCTWorld%MCT_comm, ierr)
218     if(ierr /= 0) then
219        call MP_perr_die(myname_,':: call MPI_SEND(iGGrid%descend...',ierr)
220     endif
221
222  endif ! if(HeaderAssoc(3))...
223
224       ! Step 5.  If iGGrid%weight_list is defined, send it.
225
226  if(HeaderAssoc(4)) then
227
228    call List_send(iGGrid%weight_list, dest, TagBase+7, ThisMCTWorld%MCT_comm, ierr)
229    if(ierr /= 0) then
230       if(present(status)) then
231          write(stderr,*) myname_,':: call List_send(iGGrid%weight_list...'
232          status = ierr
233          return
234       else
235          call die(myname_,':: call List_send(iGGrid%weight_list...',ierr)
236       endif
237    endif
238
239  endif ! if(HeaderAssoc(4))...
240
241       ! Step 6.  If iGGrid%other_list is defined, send it.
242
243  if(HeaderAssoc(5)) then
244
245    call List_send(iGGrid%other_list, dest, TagBase+9, ThisMCTWorld%MCT_comm, ierr)
246    if(ierr /= 0) then
247       if(present(status)) then
248          write(stderr,*) myname_,':: call List_send(iGGrid%other_list...'
249          status = ierr
250          return
251       else
252          call die(myname_,':: call List_send(iGGrid%other_list...',ierr)
253       endif
254    endif
255
256  endif ! if(HeaderAssoc(5))...
257
258       ! Step 7.  If iGGrid%index_list is defined, send it.
259
260  if(HeaderAssoc(6)) then
261
262    call List_send(iGGrid%index_list, dest, TagBase+11, ThisMCTWorld%MCT_comm, ierr)
263    if(ierr /= 0) then
264       if(present(status)) then
265          write(stderr,*) myname_,':: call List_send(iGGrid%index_list...'
266          status = ierr
267          return
268       else
269          call die(myname_,':: call List_send(iGGrid%index_list...',ierr)
270       endif
271    endif
272
273  else  ! This constitutes an error, as a GeneralGrid must at a minimum
274        ! contain the index GlobGridNum
275
276     if(present(status)) then
277        write(stderr,*) myname_,':: Error.  GeneralGrid%index_list undefined.'
278        status = -2
279        return
280     else
281        call die(myname_,'::  Error.  GeneralGrid%index_list undefined.',-2)
282     endif
283
284  endif ! if(HeaderAssoc(6))...
285
286       ! Step 8.  Finally, send the AttrVect iGGrid%data.
287
288  call AttrVect_send(iGGrid%data, dest, TagBase+13, ThisMCTWorld%MCT_comm, ierr)
289  if(ierr /= 0) then
290     if(present(status)) then
291        write(stderr,*) myname_,':: call AttrVect_send(iGGrid%data...'
292        status = ierr
293        return
294     else
295        call die(myname_,':: call AttrVect_send(iGGrid%data...',ierr)
296     endif
297  endif
298
299       ! The GeneralGrid send is now complete.
300
301 end subroutine send_
302
303!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
304!    Math and Computer Science Division, Argonne National Laboratory   !
305!BOP -------------------------------------------------------------------
306!
307! !IROUTINE: recv_ - Point-to-point blocking recv for the GeneralGrid.
308!
309! !DESCRIPTION:  The point-to-point receive routine {\tt recv\_()}
310! receives the output {\tt GeneralGrid} argument {\tt oGGrid} from component
311! {\tt comp\_id}.  The message is identified by the tag defined by the
312! {\tt INTEGER} argument {\tt TagBase}.  The value of {\tt TagBase} must
313! match the value used in the call to {\tt send\_()} on the other component.
314! The success (failure) of this operation corresponds to a zero (nonzero)
315! value for the output {\tt INTEGER} flag {\tt status}.
316!
317! {\bf N.B.}:  This routine assumes that the {\tt GeneralGrid} argument
318! {\tt oGGrid} is uninitialized on input; that is, all the {\tt List}
319! components are blank, the {\tt LOGICAL} array {\tt oGGrid\%descend} is
320! unallocated, and the {\tt AttrVect} component {\tt oGGrid\%data} is
321! uninitialized.  The {\tt GeneralGrid} {\tt oGGrid} represents allocated
322! memory.  When the user no longer needs {\tt oGGrid}, it should be
323! deallocated by invoking {\tt GeneralGrid\_clean()} (see
324! {\tt m\_GeneralGrid} for further details).
325!
326! {\bf N.B.}:  One must avoid assigning elsewhere the MPI tag values
327! between {\tt TagBase} and {\tt TagBase+20}, inclusive.  This is
328! because {\tt recv\_()} performs one receive operation set up the header
329! transfer, up to five {\tt List\_recv} operations (two {\tt MPI\_RECV}
330! calls in each), two receive operations to transfer {\tt iGGrid\%descend(:)},
331! and finally the receive of the {\tt AttrVect} component {\tt iGGrid\%data}
332! (which comprises eight {\tt MPI\_RECV} operations).
333!
334! !INTERFACE:
335
336 subroutine recv_(oGGrid, comp_id, TagBase, status)
337
338!
339! !USES:
340!
341      use m_stdio
342      use m_die
343      use m_mpif90
344
345      use m_GeneralGrid, only : GeneralGrid
346      use m_GeneralGrid, only : GeneralGrid_init => init
347      use m_GeneralGrid, only : GeneralGrid_lsize => lsize
348
349      use m_MCTWorld, only : ComponentToWorldRank
350      use m_MCTWorld, only : ThisMCTWorld
351
352      use m_AttrVectComms,only : AttrVect_recv => recv
353
354      use m_List,only : List_recv => recv
355      use m_List,only : List_nullify => nullify
356
357      implicit none
358
359! !INPUT PARAMETERS:
360!
361      integer,           intent(in) :: comp_id
362      integer,           intent(in) :: TagBase
363
364! !OUTPUT PARAMETERS:
365!
366      type(GeneralGrid), intent(out) :: oGGrid
367      integer, optional, intent(out) :: status
368
369! !REVISION HISTORY:
370!       04Jun01 - J.W. Larson <larson@mcs.anl.gov> - API Specification.
371!       07Jun01 - J.W. Larson <larson@mcs.anl.gov> - Initial version.
372!       10Jun01 - J.W. Larson <larson@mcs.anl.gov> - Bug fixes--now works.
373!       11Jun01 - R. Jacob <jacob@mcs.anl.gov> use component id as input
374!                 argument.
375!       13Jun01 - J.W. Larson <larson@mcs.anl.gov> - Initialize status
376!                 (if present).
377!       13Jun02 - J.W. Larson <larson@mcs.anl.gov> - Removed the argument
378!                 comm.  This routine is now explicitly for intercomponent
379!                 communications only.
380!EOP ___________________________________________________________________
381
382  character(len=*),parameter :: myname_=myname//'::recv_'
383
384  integer :: ierr
385  integer :: source
386  integer :: MPstatus(MP_STATUS_SIZE), DescendSize
387  logical :: HeaderAssoc(6)
388
389! for now, assume the components root is the source.
390  source = ComponentToWorldRank(0, comp_id, ThisMCTWorld)
391
392      ! Step 1. Receive the elements of the LOGICAL flag array
393      ! HeaderAssoc.  TRUE entries in this array correspond to
394      ! Check elements of the GeneralGrid header that are not
395      ! blank, and are being sent by process source.
396      !     
397      ! The significance of the entries of HeaderAssoc has been
398      ! defined in send_().  Here are the definitions of these
399      ! values:
400      !
401      !  HeaderAssoc(1) = List_allocated(oGGrid%coordinate_list)
402      !  HeaderAssoc(2) = List_allocated(oGGrid%coordinate_sort_order)
403      !  HeaderAssoc(3) = associated(oGGrid%descend)
404      !  HeaderAssoc(4) = List_allocated(oGGrid%weight_list)
405      !  HeaderAssoc(5) = List_allocated(oGGrid%other_list)
406      !  HeaderAssoc(6) = List_allocated(oGGrid%index_list)
407
408      ! Initialize status (if present)
409
410  if(present(status)) status = 0
411
412      ! Step 1.  Nullify oGGrid components, set HeaderAssoc(:) to .FALSE.,
413      ! then receive incoming HeaderAssoc(:) data
414
415  call List_nullify(oGGrid%coordinate_list)
416  call List_nullify(oGGrid%coordinate_sort_order)
417  call List_nullify(oGGrid%weight_list)
418  call List_nullify(oGGrid%other_list)
419  call List_nullify(oGGrid%index_list)
420  nullify(oGGrid%descend)
421
422  HeaderAssoc = .FALSE.
423
424  call MPI_RECV(HeaderAssoc, 6, MP_LOGICAL, source, TagBase, ThisMCTWorld%MCT_comm, MPstatus, ierr)
425  if(ierr /= 0) then
426     call MP_perr_die(myname_,':: MPI_RECV(HeaderAssoc...',ierr)
427  endif
428
429       ! Step 2.  If oGGrid%coordinate_list is defined, receive it.
430
431  if(HeaderAssoc(1)) then
432    call List_recv(oGGrid%coordinate_list, source, TagBase+1, ThisMCTWorld%MCT_comm, ierr)
433    if(ierr /= 0) then
434       if(present(status)) then
435          write(stderr,*) myname_,':: call List_recv(oGGrid%coordinate_list...'
436          status = ierr
437          return
438       else
439          call die(myname_,':: call List_recv(oGGrid%coordinate_list...',ierr)
440       endif
441    endif
442  else  ! This constitutes an error, as a GeneralGrid must have coordinates
443
444     if(present(status)) then
445        write(stderr,*) myname_,':: Error.  GeneralGrid%coordinate_list undefined.'
446        status = -1
447        return
448     else
449        call die(myname_,'::  Error.  GeneralGrid%coordinate_list undefined.',-1)
450     endif
451
452  endif ! if(HeaderAssoc(1))...
453
454       ! Step 3.  If oGGrid%coordinate_sort_order is defined, receive it.
455
456  if(HeaderAssoc(2)) then
457     call List_recv(oGGrid%coordinate_sort_order, source, TagBase+3, ThisMCTWorld%MCT_comm, ierr)
458     if(ierr /= 0) then
459        if(present(status)) then
460           write(stderr,*) myname_,':: Error calling ',&
461                'List_recv(oGGrid%coordinate_sort_order...'
462           status = ierr
463           return
464        else
465           call die(myname_,':: call List_recv(oGGrid%coordinate_sort_order...', ierr)
466        endif
467     endif
468  endif ! if(HeaderAssoc(2))...
469
470       ! Step 4.  If oGGrid%descend is allocated, determine its size,
471       ! receive this size, allocate oGGrid%descend, and then receive
472       ! the elements of oGGrid%descend.
473     
474  if(HeaderAssoc(3)) then
475
476     call MPI_RECV(DescendSize, 1, MP_type(DescendSize), &
477                   source, TagBase+5, ThisMCTWorld%MCT_comm, MPstatus, ierr)
478     if(ierr /= 0) then
479           call MP_perr_die(myname_,':: call MPI_RECV(size(oGGrid%descend)...',ierr)
480     endif
481
482     allocate(oGGrid%descend(DescendSize), stat=ierr)
483     if(ierr /= 0) then
484        if(present(status)) then
485           write(stderr,*) myname_,':: allocate(oGGrid%descend...'
486           status = ierr
487           return
488        else
489           call die(myname_,':: allocate(oGGrid%descend... failed.',ierr)
490        endif
491     endif
492
493     call MPI_RECV(oGGrid%descend, DescendSize, MP_type(oGGrid%descend(1)), &
494                   source, TagBase+6, ThisMCTWorld%MCT_comm, MPstatus, ierr)
495     if(ierr /= 0) then
496           call MP_perr_die(myname_,':: call MPI_RECV(oGGrid%descend...',ierr)
497     endif
498
499  endif ! if(HeaderAssoc(3))...
500
501       ! Step 5.  If oGGrid%weight_list is defined, receive it.
502
503  if(HeaderAssoc(4)) then
504
505    call List_recv(oGGrid%weight_list, source, TagBase+7, ThisMCTWorld%MCT_comm, ierr)
506    if(ierr /= 0) then
507       if(present(status)) then
508          write(stderr,*) myname_,':: call List_recv(oGGrid%weight_list...'
509          status = ierr
510          return
511       else
512          call die(myname_,':: call List_recv(oGGrid%weight_list...',ierr)
513       endif
514    endif
515
516  endif ! if(HeaderAssoc(4))...
517
518       ! Step 6.  If oGGrid%other_list is defined, receive it.
519
520  if(HeaderAssoc(5)) then
521
522    call List_recv(oGGrid%other_list, source, TagBase+9, ThisMCTWorld%MCT_comm, ierr)
523    if(ierr /= 0) then
524       if(present(status)) then
525          write(stderr,*) myname_,':: call List_recv(oGGrid%other_list...'
526          status = ierr
527          return
528       else
529          call die(myname_,':: call List_recv(oGGrid%other_list...',ierr)
530       endif
531    endif
532
533  endif ! if(HeaderAssoc(5))...
534
535       ! Step 7.  If oGGrid%index_list is defined, receive it.
536
537  if(HeaderAssoc(6)) then
538
539    call List_recv(oGGrid%index_list, source, TagBase+11, ThisMCTWorld%MCT_comm, ierr)
540    if(ierr /= 0) then
541       if(present(status)) then
542          write(stderr,*) myname_,':: call List_recv(oGGrid%index_list...'
543          status = ierr
544          return
545       else
546          call die(myname_,':: call List_recv(oGGrid%index_list...',ierr)
547       endif
548    endif
549
550  else  ! This constitutes an error, as a GeneralGrid must at a minimum
551        ! contain the index GlobGridNum
552
553     if(present(status)) then
554        write(stderr,*) myname_,':: Error.  GeneralGrid%index_list undefined.'
555        status = -2
556        return
557     else
558        call die(myname_,'::  Error.  GeneralGrid%index_list undefined.',-2)
559     endif
560
561  endif ! if(HeaderAssoc(6))...
562
563       ! Step 8.  Finally, receive the AttrVect oGGrid%data.
564
565  call AttrVect_recv(oGGrid%data, source, TagBase+13, ThisMCTWorld%MCT_comm, ierr)
566  if(ierr /= 0) then
567     if(present(status)) then
568        write(stderr,*) myname_,':: call AttrVect_recv(oGGrid%data...'
569        status = ierr
570        return
571     else
572        call die(myname_,':: call AttrVect_recv(oGGrid%data...',ierr)
573     endif
574  endif
575
576       ! The GeneralGrid receive is now complete.
577
578 end subroutine recv_
579
580!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
581!    Math and Computer Science Division, Argonne National Laboratory   !
582!BOP -------------------------------------------------------------------
583!
584! !IROUTINE: GM_gather_ - gather a GeneralGrid using input GlobalMap.
585!
586! !DESCRIPTION:  {\tt GM\_gather\_()} takes an input {\tt GeneralGrid}
587! argument {\tt iG} whose decomposition on the communicator associated
588! with the F90 handle {\tt comm} is described by the {\tt GlobalMap}
589! argument {\tt GMap}, and gathers it to the {\tt GeneralGrid} output
590! argument {\tt oG} on the {\tt root}.  The success (failure) of this
591! operation is reported as a zero (nonzero) value in the optional
592! {\tt INTEGER} output argument {\tt stat}.
593
594! {\bf N.B.}:  An important assumption made here is that the distributed
595! {\tt GeneralGrid} {\tt iG} has been initialized with the same
596! coordinate system, sort order, other real attributes, and the same
597! indexing attributes for all processes on {\tt comm}.
598!
599! {\bf N.B.}:  Once the gridpoint data of the {\tt GeneralGrid} are assembled
600! on the {\tt root}, they are stored in the order determined by the input
601! {\tt GlobalMap} {\tt GMap}.  The user may need to sorted these gathered
602! data to order them in  accordance with the {\tt coordinate\_sort\_order}
603! attribute of {\tt iG}.
604!
605! {\bf N.B.}:  The output {\tt GeneralGrid} {\tt oG} represents allocated
606! memory on the {\tt root}.  When the user no longer needs {\tt oG} it
607! should be deallocated using {\tt GeneralGrid\_clean()} to avoid a memory
608! leak
609!
610! !INTERFACE:
611!
612 subroutine GM_gather_(iG, oG, GMap, root, comm, stat)
613!
614! !USES:
615!
616      use m_stdio
617      use m_die
618      use m_mpif90
619
620      use m_GlobalMap, only : GlobalMap
621      use m_GlobalMap, only : GlobalMap_gsize => gsize
622
623      use m_GeneralGrid, only : GeneralGrid
624      use m_GeneralGrid, only : GeneralGrid_init => init
625
626      use m_AttrVectComms,only : AttrVect_Gather => gather
627
628      implicit none
629
630! !INPUT PARAMETERS:
631!
632      type(GeneralGrid), intent(in)  :: iG
633      type(GlobalMap),   intent(in)  :: GMap
634      integer,           intent(in)  :: root
635      integer,           intent(in)  :: comm
636
637! !OUTPUT PARAMETERS:
638!
639      type(GeneralGrid), intent(out) :: oG
640      integer, optional, intent(out) :: stat
641
642! !REVISION HISTORY:
643!       27Apr01 - J.W. Larson <larson@mcs.anl.gov> - API Specification.
644!       02May01 - J.W. Larson <larson@mcs.anl.gov> - Initial code.
645!       13Jun01 - J.W. Larson <larson@mcs.anl.gov> - Initialize stat
646!                 (if present).
647!EOP ___________________________________________________________________
648
649 character(len=*),parameter :: myname_=myname//'::GM_gather_'
650!Process ID
651 integer :: myID
652!Error flag
653 integer :: ierr
654!Number of points on the _Gathered_ grid:
655 integer :: length
656
657      ! Initialize stat (if present)
658
659  if(present(stat)) stat = 0
660
661       ! Which process am I?
662
663  call MPI_COMM_RANK(comm, myID, ierr)
664  if(ierr /= 0) then
665    call MP_perr_die(myname_,'call MPI_COMM_RANK()',ierr)
666  endif
667
668  if(myID == root) then ! prepare oG:
669
670       ! The length of the _gathered_ GeneralGrid oG is determined by
671       ! the GlobalMap function GlobalMap_gsize()
672
673     length = GlobalMap_gsize(GMap)
674
675       ! Initialize attributes of oG from iG
676     call copyGeneralGridHeader_(iG,oG) 
677
678  endif
679
680       ! Gather gridpoint data in iG%data to oG%data
681
682  call AttrVect_Gather(iG%data, oG%data, GMap, root, comm, ierr)
683
684  if(ierr /= 0) then
685     write(stderr,*) myname_,':: Error--call AttrVect_Gather() failed.', &
686          ' ierr = ',ierr
687     if(present(stat)) then
688        stat=ierr
689        return
690     else
691        call die(myname_,'call AttrVect_Gather(ig%data...',ierr)
692     endif
693  endif
694
695 end subroutine GM_gather_
696
697!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
698!    Math and Computer Science Division, Argonne National Laboratory   !
699!BOP -------------------------------------------------------------------
700!
701! !IROUTINE: GSM_gather_ - gather a GeneralGrid using input GlobalSegMap.
702!
703! !DESCRIPTION:  {\tt GMS\_gather\_()} takes an input {\tt GeneralGrid}
704! argument {\tt iG} whose decomposition on the communicator associated
705! with the F90 handle {\tt comm} is described by the {\tt GlobalSegMap}
706! argument {\tt GSMap}, and gathers it to the {\tt GeneralGrid} output
707! argument {\tt oG} on the {\tt root}.  The success (failure) of this
708! operation is reported as a zero (nonzero) value in the optional
709! {\tt INTEGER} output argument {\tt stat}.
710!
711! {\bf N.B.}:  An important assumption made here is that the distributed
712! {\tt GeneralGrid} {\tt iG} has been initialized with the same
713! coordinate system, sort order, other real attributes, and the same
714! indexing attributes for all processes on {\tt comm}.
715!
716! {\bf N.B.}:  Once the gridpoint data of the {\tt GeneralGrid} are assembled
717! on the {\tt root}, they are stored in the order determined by the input
718! {\tt GlobalSegMap} {\tt GSMap}.  The user may need to sorted these gathered
719! data to order them in  accordance with the {\tt coordinate\_sort\_order}
720! attribute of {\tt iG}.
721!
722! {\bf N.B.}:  The output {\tt GeneralGrid} {\tt oG} represents allocated
723! memory on the {\tt root}.  When the user no longer needs {\tt oG} it
724! should be deallocated using {\tt GeneralGrid\_clean()} to avoid a memory
725! leak
726!
727! !INTERFACE:
728
729 subroutine GSM_gather_(iG, oG, GSMap, root, comm, stat)
730!
731! !USES:
732!
733      use m_stdio
734      use m_die
735      use m_mpif90
736
737      use m_GlobalSegMap, only : GlobalSegMap
738      use m_GlobalSegMap, only : GlobalSegMap_lsize => lsize
739      use m_GlobalSegMap, only : GlobalSegMap_gsize => gsize
740
741      use m_GeneralGrid, only : GeneralGrid
742      use m_GeneralGrid, only : GeneralGrid_init => init
743      use m_GeneralGrid, only : GeneralGrid_lsize => lsize
744
745      use m_AttrVectComms,only : AttrVect_Gather => gather
746
747      implicit none
748
749! !INPUT PARAMETERS:
750!
751      type(GeneralGrid),  intent(in)  :: iG
752      type(GlobalSegMap), intent(in)  :: GSMap
753      integer,            intent(in)  :: root
754      integer,            intent(in)  :: comm
755
756! !OUTPUT PARAMETERS:
757!
758      type(GeneralGrid),  intent(out) :: oG
759      integer, optional,  intent(out) :: stat
760
761! !REVISION HISTORY:
762!       27Apr01 - J.W. Larson <larson@mcs.anl.gov> - API Specification.
763!       01May01 - J.W. Larson <larson@mcs.anl.gov> - Working Version.
764!       13Jun01 - J.W. Larson <larson@mcs.anl.gov> - Initialize stat
765!                 (if present).
766!EOP ___________________________________________________________________
767
768  character(len=*),parameter :: myname_=myname//'::GSM_gather_'
769
770!Process ID
771 integer :: myID
772!Error flag
773 integer :: ierr
774!Number of points on the _Gathered_ grid:
775 integer :: length
776
777      ! Initialize stat (if present)
778
779  if(present(stat)) stat = 0
780
781       ! Which process am I?
782
783  call MPI_COMM_RANK(comm, myID, ierr)
784  if(ierr /= 0) then
785    call MP_perr_die(myname_,'MPI_COMM_RANK()',ierr)
786  endif
787
788  if(myID == root) then ! prepare oG:
789
790       ! The length of the _gathered_ GeneralGrid oG is determined by
791       ! the GlobalMap function GlobalSegMap_gsize()
792
793     length = GlobalSegMap_gsize(GSMap)
794
795       ! Initialize attributes of oG from iG
796     call copyGeneralGridHeader_(iG,oG) 
797
798  endif
799
800       ! Gather gridpoint data in iG%data to oG%data
801
802  call AttrVect_Gather(iG%data, oG%data, GSMap, root, comm, ierr)
803  if(ierr /= 0) then
804     write(stderr,*) myname_,':: ERROR--call AttrVect_Gather() failed.', &
805          ' ierr = ',ierr
806    if(present(stat)) then
807       stat=ierr
808       return
809    else
810       call die(myname_)
811    endif
812  endif
813
814 end subroutine GSM_gather_
815
816!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
817!    Math and Computer Science Division, Argonne National Laboratory   !
818!BOP -------------------------------------------------------------------
819!
820! !IROUTINE: GM_scatter_ - scatter a GeneralGrid using input GlobalMap.
821!
822! !DESCRIPTION:  {\tt GM\_scatter\_()} takes an input {\tt GeneralGrid}
823! argument {\tt iG} (valid only on the {\tt root} process), and scatters
824! it to the distributed {\tt GeneralGrid} variable {\tt oG}.  The
825! {\tt GeneralGrid} {\tt oG} is distributed on the communicator
826! associated with the F90  handle {\tt comm} using the domain
827! decomposition described by the {\tt GlobalMap} argument {\tt GMap}.
828! The success (failure) of this operation is reported as a zero (nonzero)
829! value in the optional {\tt INTEGER} output argument {\tt stat}.
830!
831! {\bf N.B.}:  The output {\tt GeneralGrid} {\tt oG} represents allocated
832! memory on the {\tt root}.  When the user no longer needs {\tt oG} it
833! should be deallocated using {\tt GeneralGrid\_clean()} to avoid a memory
834! leak.
835!
836! !INTERFACE:
837
838 subroutine GM_scatter_(iG, oG, GMap, root, comm, stat)
839!
840! !USES:
841!
842      use m_stdio
843      use m_die
844      use m_mpif90
845
846      use m_GlobalMap, only : GlobalMap
847      use m_GlobalMap, only : GlobalMap_lsize => lsize
848      use m_GlobalMap, only : GlobalMap_gsize => gsize
849
850      use m_AttrVectComms, only : AttrVect_scatter => scatter
851
852      use m_GeneralGrid, only : GeneralGrid
853      use m_GeneralGrid, only : GeneralGrid_init => init
854      use m_GeneralGrid, only : GeneralGrid_lsize => lsize
855
856      implicit none
857
858! !INPUT PARAMETERS:
859!
860      type(GeneralGrid), intent(in)  :: iG
861      type(GlobalMap),   intent(in)  :: GMap
862      integer,           intent(in)  :: root
863      integer,           intent(in)  :: comm
864
865! !OUTPUT PARAMETERS:
866!
867      type(GeneralGrid), intent(out) :: oG
868      integer, optional, intent(out) :: stat
869
870! !REVISION HISTORY:
871!       27Apr01 - J.W. Larson <larson@mcs.anl.gov> - API Specification.
872!       04Jun01 - J.W. Larson <larson@mcs.anl.gov> - Changed comms model
873!                 to MPI-style (i.e. iG valid on root only).
874!       13Jun01 - J.W. Larson <larson@mcs.anl.gov> - Initialize stat
875!                 (if present).
876!EOP ___________________________________________________________________
877
878  character(len=*),parameter :: myname_=myname//'::GM_scatter_'
879
880  logical :: DescendAssoc
881  integer :: DescendSize
882  integer :: ierr, myID
883
884      ! Initialize status (if present)
885
886  if(present(stat)) stat = 0
887
888       ! Step 1.  Determine process ID number myID
889
890  call MPI_COMM_RANK(comm, myID, ierr)
891  if(ierr /= 0) then
892     call MP_perr_die(myname_,'MPI_COMM_RANK(comm...',ierr)
893  endif
894
895       ! Step 2.  On the root, initialize the List and LOGICAL
896       ! attributes of the GeneralGrid variable iG to oG.
897
898  if(myID == root) then
899     call copyGeneralGridHeader_(iG, oG)
900  endif
901
902       ! Step 3.  Broadcast from the root the List and LOGICAL
903       ! attributes of the GeneralGrid variable oG.
904
905  call bcastGeneralGridHeader_(oG, root, comm, ierr)
906  if(ierr /= 0) then
907     write(stderr,*) myname_,':: Error calling bcastGeneralGridHeader_().',&
908          ' ierr = ',ierr
909     if(present(stat)) then
910        stat = ierr
911        return
912     else
913        call die(myname_,'call bcastGeneralGridHeader_(oG...',ierr)
914     endif
915  endif
916
917
918       ! Step 4.  Using the GeneralMap GMap, scatter the AttrVect
919       ! portion of the input GeneralGrid iG to the GeneralGrid oG.
920
921  call AttrVect_scatter(iG%data, oG%data, GMap, root, comm, ierr)
922  if(ierr /= 0) then
923     write(stderr,*) myname_,':: Error calling AttrVect_scatter(iG%data...',&
924          ' ierr = ',ierr
925     if(present(stat)) then
926        stat = ierr
927        return
928     else
929        call die(myname_,'call AttrVect_scatter(iG%data...',ierr)
930     endif
931  endif
932
933       ! The GeneralGrid scatter is now complete.
934
935 end subroutine GM_scatter_
936
937!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
938!    Math and Computer Science Division, Argonne National Laboratory   !
939!BOP -------------------------------------------------------------------
940!
941! !IROUTINE: GSM_scatter_ - scatter a GeneralGrid using input GlobalSegMap.
942!
943! !DESCRIPTION:  {\tt GM\_scatter\_()} takes an input {\tt GeneralGrid}
944! argument {\tt iG} (valid only on the {\tt root} process), and scatters
945! it to the distributed {\tt GeneralGrid} variable {\tt oG}.  The
946! {\tt GeneralGrid} {\tt oG} is distributed on the communicator
947! associated with the F90  handle {\tt comm} using the domain
948! decomposition described by the {\tt GlobalSegMap} argument {\tt GSMap}.
949! The success (failure) of this operation is reported as a zero (nonzero)
950! value in the optional {\tt INTEGER} output argument {\tt stat}.
951!
952! {\bf N.B.}:  The output {\tt GeneralGrid} {\tt oG} represents allocated
953! memory on the {\tt root}.  When the user no longer needs {\tt oG} it
954! should be deallocated using {\tt GeneralGrid\_clean()} to avoid a memory
955! leak.
956!
957! !INTERFACE:
958
959 subroutine GSM_scatter_(iG, oG, GSMap, root, comm, stat)
960!
961! !USES:
962!
963      use m_stdio
964      use m_die
965      use m_mpif90
966
967      use m_GlobalSegMap, only : GlobalSegMap
968      use m_GlobalSegMap, only : GlobalSegMap_lsize => lsize
969      use m_GlobalSegMap, only : GlobalSegMap_gsize => gsize
970
971      use m_AttrVectComms, only : AttrVect_scatter => scatter
972
973      use m_GeneralGrid, only : GeneralGrid
974      use m_GeneralGrid, only : GeneralGrid_init => init
975      use m_GeneralGrid, only : GeneralGrid_lsize => lsize
976
977      implicit none
978
979! !INPUT PARAMETERS:
980!
981      type(GeneralGrid),  intent(in)  :: iG
982      type(GlobalSegMap), intent(in)  :: GSMap
983      integer,            intent(in)  :: root
984      integer,            intent(in)  :: comm
985
986! !OUTPUT PARAMETERS:
987!
988      type(GeneralGrid),  intent(out) :: oG
989      integer, optional,  intent(out) :: stat
990
991! !REVISION HISTORY:
992!       27Apr01 - J.W. Larson <larson@mcs.anl.gov> - API Specification.
993!       04Jun01 - J.W. Larson <larson@mcs.anl.gov> - Initial code.
994!       13Jun01 - J.W. Larson <larson@mcs.anl.gov> - Initialize stat
995!                 (if present).
996!EOP ___________________________________________________________________
997
998  character(len=*),parameter :: myname_=myname//'::GSM_scatter_'
999
1000  integer :: ierr, myID
1001
1002      ! Initialize stat (if present)
1003
1004  if(present(stat)) stat = 0
1005
1006       ! Step 1.  Determine process ID number myID
1007
1008  call MPI_COMM_RANK(comm, myID, ierr)
1009  if(ierr /= 0) then
1010     call MP_perr_die(myname_,'MPI_COMM_RANK(comm...',ierr)
1011  endif
1012
1013       ! Step 2.  On the root, initialize the List and LOGICAL
1014       ! attributes of the GeneralGrid variable iG to oG.
1015
1016  if(myID == root) then
1017     call copyGeneralGridHeader_(iG, oG)
1018  endif
1019
1020       ! Step 3.  Broadcast from the root the List and LOGICAL
1021       ! attributes of the GeneralGrid variable oG.
1022
1023  call bcastGeneralGridHeader_(oG, root, comm, ierr)
1024  if(ierr /= 0) then
1025     write(stderr,*) myname_,':: Error calling bcastGeneralGridHeader_(...',&
1026          ' ierr = ',ierr
1027     if(present(stat)) then
1028        stat = ierr
1029        return
1030     else
1031        call die(myname_,'bcastGeneralGridHeader_(oG...',ierr)
1032     endif
1033  endif
1034
1035       ! Step 4.  Using the GeneralSegMap GSMap, scatter the AttrVect
1036       ! portion of the input GeneralGrid iG to the GeneralGrid oG.
1037
1038  call AttrVect_scatter(iG%data, oG%data, GSMap, root, comm, ierr)
1039  if(ierr /= 0) then
1040     write(stderr,*) myname_,':: Error calling AttrVect_scatter(iG%data...',&
1041          ' ierr = ',ierr
1042     if(present(stat)) then
1043        stat = ierr
1044        return
1045     else
1046        call die(myname_,'call AttrVect_scatter(iG%data...',ierr)
1047     endif
1048  endif
1049
1050       ! The GeneralGrid scatter is now complete.
1051
1052 end subroutine GSM_scatter_
1053
1054!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1055!    Math and Computer Science Division, Argonne National Laboratory   !
1056!BOP -------------------------------------------------------------------
1057!
1058! !IROUTINE: bcast_ - Broadcast a GeneralGrid.
1059!
1060! !DESCRIPTION:  {\tt bcast\_()} takes an input {\tt GeneralGrid}
1061! argument {\tt ioG} (valid only on the {\tt root} process), and
1062! broadcasts it to all processes on the communicator associated with the
1063! F90  handle {\tt comm}.  The success (failure) of this operation is
1064! reported as a zero (nonzero) value in the optional {\tt INTEGER}
1065! output argument {\tt stat}.
1066!
1067! {\bf N.B.}:  On the non-root processes, the output {\tt GeneralGrid}
1068! {\tt ioG} represents allocated memory.  When the user no longer needs
1069! {\tt ioG} it should be deallocated by invoking {\tt GeneralGrid\_clean()}.
1070! Failure to do so risks a memory leak.
1071!
1072! !INTERFACE:
1073
1074 subroutine bcast_(ioG, root, comm, stat)
1075!
1076! !USES:
1077!
1078      use m_stdio
1079      use m_die
1080      use m_mpif90
1081
1082      use m_GlobalSegMap, only : GlobalSegMap
1083      use m_GlobalSegMap, only : GlobalSegMap_lsize => lsize
1084      use m_GlobalSegMap, only : GlobalSegMap_gsize => gsize
1085
1086      use m_GeneralGrid, only : GeneralGrid
1087      use m_GeneralGrid, only : GeneralGrid_init => init
1088      use m_GeneralGrid, only : GeneralGrid_lsize => lsize
1089
1090      use m_AttrVectComms,only : AttrVect_bcast => bcast
1091
1092      implicit none
1093
1094! !INPUT PARAMETERS:
1095!
1096      integer,           intent(in)    :: root
1097      integer,           intent(in)    :: comm
1098
1099! !INPUT/OUTPUT PARAMETERS:
1100!
1101      type(GeneralGrid), intent(inout) :: ioG
1102
1103! !OUTPUT PARAMETERS:
1104!
1105      integer, optional, intent(out)   :: stat
1106
1107! !REVISION HISTORY:
1108!       27Apr01 - J.W. Larson <larson@mcs.anl.gov> - API Specification.
1109!       02May01 - J.W. Larson <larson@mcs.anl.gov> - Initial version.
1110!       13Jun01 - J.W. Larson <larson@mcs.anl.gov> - Initialize stat
1111!                 (if present).
1112!EOP ___________________________________________________________________
1113
1114  character(len=*),parameter :: myname_=myname//'::bcast_'
1115
1116  integer :: ierr, myID
1117
1118      ! Initialize status (if present)
1119
1120  if(present(stat)) stat = 0
1121
1122       ! Step 1.  Determine process ID number myID
1123
1124  call MPI_COMM_RANK(comm, myID, ierr)
1125  if(ierr /= 0) then
1126     call MP_perr_die(myname_,'MPI_COMM_RANK(comm...',ierr)
1127  endif
1128
1129       ! Step 2.  Broadcast from the root the List and LOGICAL
1130       ! attributes of the GeneralGrid variable ioG.
1131
1132  call bcastGeneralGridHeader_(ioG, root, comm, ierr)
1133  if(ierr /= 0) then
1134     write(stderr,*) myname_,':: Error calling bcastGeneralGridHeader_(...',&
1135          ' ierr = ',ierr
1136     if(present(stat)) then
1137        stat = ierr
1138        return
1139     else
1140        call die(myname_)
1141     endif
1142  endif
1143
1144       ! Step 3.  Broadcast ioG%data from the root.
1145
1146  call AttrVect_bcast(ioG%data, root, comm, ierr)
1147  if(ierr /= 0) then
1148     write(stderr,*) myname_,':: Error calling AttrVect_scatter(iG%data...',&
1149          ' ierr = ',ierr
1150     if(present(stat)) then
1151        stat = ierr
1152        return
1153     else
1154        call die(myname_)
1155     endif
1156  endif
1157
1158       ! The GeneralGrid broadcast is now complete.
1159
1160 end subroutine bcast_
1161
1162!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1163!    Math and Computer Science Division, Argonne National Laboratory   !
1164!BOP -------------------------------------------------------------------
1165!
1166! !IROUTINE: bcastGeneralGridHeader_ - Broadcast the GeneralGrid Header.
1167!
1168! !DESCRIPTION:  This routine broadcasts the header information from
1169! the input {\tt GeneralGrid} argument {\tt ioGGrid} (on input valid
1170! on the {\tt root} only).  This broadcast is from the {\tt root} to
1171! all processes on the communicator associated with the fortran 90
1172! {\tt INTEGER} handle {\tt comm}.  The success (failure) of this operation
1173! corresponds to a zero (nonzero) value for the output {\tt INTEGER} flag
1174! {\tt stat}.
1175!
1176! The {\em header information} in a {\tt GeneralGrid} variable comprises
1177! all the non-{\tt AttrVect} components of the {\tt GeneralGrid}; that
1178! is, everything except the gridpoint coordinate, geometry, and index
1179! data stored in {\tt iGGrid\%data}.  This information includes:
1180! \begin{enumerate}
1181! \item The coordinates in {\tt iGGrid\%coordinate\_list}
1182! \item The coordinate sort order in {\tt iGGrid\%coordinate\_sort\_order}
1183! \item The area/volume weights in {\tt iGGrid\%weight\_list}
1184! \item Other {\tt REAL} geometric information in {\tt iGGrid\%other\_list}
1185! \item Indexing information in {\tt iGGrid\%index\_list}
1186! \item The {\tt LOGICAL} descending/ascending order sort flags in
1187! {\tt iGGrid\%descend(:)}.
1188! \end{enumerate}
1189!
1190! !INTERFACE:
1191
1192 subroutine bcastGeneralGridHeader_(ioGGrid, root, comm, stat)
1193!
1194! !USES:
1195!
1196      use m_stdio
1197      use m_die
1198      use m_mpif90
1199
1200      use m_GlobalSegMap, only : GlobalSegMap
1201      use m_GlobalSegMap, only : GlobalSegMap_lsize => lsize
1202      use m_GlobalSegMap, only : GlobalSegMap_gsize => gsize
1203
1204      use m_GeneralGrid, only : GeneralGrid
1205      use m_GeneralGrid, only : GeneralGrid_init => init
1206      use m_GeneralGrid, only : GeneralGrid_lsize => lsize
1207
1208      use m_List, only : List
1209      use m_List, only : List_allocated => allocated
1210      use m_List, only : List_nullify => nullify
1211      use m_List, only : List_bcast => bcast
1212
1213      implicit none
1214
1215! !INPUT PARAMETERS:
1216!
1217      integer,           intent(in)    :: root
1218      integer,           intent(in)    :: comm
1219
1220! !INPUT/OUTPUT PARAMETERS:
1221!
1222      type(GeneralGrid), intent(inout) :: ioGGrid
1223
1224! !OUTPUT PARAMETERS:
1225!
1226      integer, optional, intent(out)   :: stat
1227
1228! !REVISION HISTORY:
1229!       05Jun01 - J.W. Larson <larson@mcs.anl.gov> - Initial code.
1230!       13Jun01 - J.W. Larson <larson@mcs.anl.gov> - Initialize stat
1231!                 (if present).
1232!       05Aug02 - E. Ong <eong@mcs.anl.gov> - added association checking
1233!EOP ___________________________________________________________________
1234
1235  character(len=*),parameter :: myname_=myname//'::bcastGeneralGridHeader_'
1236
1237! Process ID
1238  integer :: myID
1239! Error flag
1240  integer :: ierr
1241! Size of array ioGGrid%descend(:)
1242  integer :: DescendSize
1243! Header-Assocation array
1244  logical :: HeaderAssoc(6)
1245
1246      ! Initialize stat (if present)
1247
1248  if(present(stat)) stat = 0
1249
1250       ! Determine process ID number myID
1251
1252  call MPI_COMM_RANK(comm, myID, ierr)
1253  if(ierr /= 0) then
1254     call MP_perr_die(myname_,'MPI_COMM_RANK(comm...',ierr)
1255  endif
1256
1257       ! Step 0.5. Check elements of the GeneralGrid header to see
1258       ! which components of it are allocated.  Load the results
1259       ! into HeaderAssoc(:), and broadcast it.
1260
1261  if(myID == root) then
1262
1263     HeaderAssoc(1) = List_allocated(ioGGrid%coordinate_list)
1264     HeaderAssoc(2) = List_allocated(ioGGrid%coordinate_sort_order)
1265     HeaderAssoc(3) = List_allocated(ioGGrid%weight_list)
1266     HeaderAssoc(4) = List_allocated(ioGGrid%other_list)
1267     HeaderAssoc(5) = List_allocated(ioGGrid%index_list)
1268     HeaderAssoc(6) = associated(ioGGrid%descend)
1269
1270  else
1271
1272     call List_nullify(ioGGrid%coordinate_list)
1273     call List_nullify(ioGGrid%coordinate_sort_order)
1274     call List_nullify(ioGGrid%weight_list)
1275     call List_nullify(ioGGrid%other_list)
1276     call List_nullify(ioGGrid%index_list)
1277     nullify(ioGGrid%descend) 
1278
1279  endif
1280
1281  call MPI_BCAST(HeaderAssoc,6,MP_LOGICAL,root,comm,ierr)
1282
1283       ! Step 1. Broadcast List attributes of the GeneralGrid.
1284
1285  if(HeaderAssoc(1)) then
1286     call List_bcast(ioGGrid%coordinate_list, root, comm, ierr)
1287     if(ierr /= 0) then
1288        write(stderr,*) myname_,'List_bcast(ioGGrid%coordinate_list... failed.',&
1289             ' ierr = ',ierr
1290        if(present(stat)) then
1291           stat = ierr
1292           return
1293        else
1294           call die(myname_)
1295        endif
1296     endif
1297  endif
1298
1299  if(HeaderAssoc(2)) then
1300     call List_bcast(ioGGrid%coordinate_sort_order, root, comm, ierr)
1301     if(ierr /= 0) then
1302        write(stderr,*) myname_,'List_bcast(ioGGrid%coordinate_sort_order... failed', &
1303             ' ierr = ',ierr
1304        if(present(stat)) then
1305           stat = ierr
1306           return
1307        else
1308           call die(myname_)
1309        endif
1310     endif
1311  endif
1312
1313  if(HeaderAssoc(3)) then
1314     call List_bcast(ioGGrid%weight_list, root, comm, ierr)
1315     if(ierr /= 0) then
1316        write(stderr,*) myname_,'List_bcast(ioGGrid%weight_list... failed',&
1317             ' ierr = ',ierr
1318        if(present(stat)) then
1319           stat = ierr
1320           return
1321        else
1322           call die(myname_)
1323        endif
1324     endif
1325  endif
1326
1327  if(HeaderAssoc(4)) then
1328     call List_bcast(ioGGrid%other_list, root, comm, ierr)
1329     if(ierr /= 0) then
1330        write(stderr,*) myname_,'List_bcast(ioGGrid%other_list... failed',&
1331             ' ierr = ',ierr
1332        if(present(stat)) then
1333           stat = ierr
1334           return
1335        else
1336           call die(myname_)
1337        endif
1338     endif
1339  endif
1340
1341  if(HeaderAssoc(5)) then
1342     call List_bcast(ioGGrid%index_list, root, comm, ierr)
1343     if(ierr /= 0) then
1344        write(stderr,*) myname_,'List_bcast(ioGGrid%index_list... failed',&
1345             ' ierr = ',ierr
1346        if(present(stat)) then
1347           stat = ierr
1348           return
1349        else
1350           call die(myname_)
1351        endif
1352     endif
1353  endif
1354
1355       ! If ioGGrid%descend is associated on the root, prepare and
1356       ! execute its broadcast
1357
1358  if(HeaderAssoc(6)) then
1359
1360       ! On the root, get the size of ioGGrid%descend(:)
1361
1362     if(myID == root) then
1363        DescendSize = size(ioGGrid%descend)
1364        if(DescendSize<=0) call die(myname_,'size(ioGGrid%descend)<=0')
1365     endif
1366
1367       ! Broadcast the size of ioGGrid%descend(:) from the root.
1368
1369     call MPI_BCAST(DescendSize, 1, MP_INTEGER, root, comm, ierr)
1370     if(ierr /= 0) then
1371        call MP_perr_die(myname_,'MPI_BCAST(DescendSize...',ierr)
1372     endif
1373
1374       ! Off the root, allocate ioGGrid%descend(:)
1375
1376     if(myID /= root) then
1377        allocate(ioGGrid%descend(DescendSize), stat=ierr)
1378        if(ierr /= 0) then
1379           write(stderr,*) myname_,':: ERROR in allocate(ioGGrid%descend...',&
1380                ' ierr = ',ierr
1381           call die(myname_)
1382        endif
1383     endif
1384 
1385      ! Finally, broadcast ioGGrid%descend(:) from the root
1386
1387     call MPI_BCAST(ioGGrid%descend, DescendSize, MP_LOGICAL, root, &
1388                    comm, ierr)
1389     if(ierr /= 0) then
1390        call MP_perr_die(myname_,'MPI_BCAST(ioGGrid%descend...',ierr)
1391     endif
1392
1393  endif
1394
1395       ! The broadcast of the GeneralGrid Header from the &
1396       ! root is complete.
1397
1398
1399 end subroutine bcastGeneralGridHeader_
1400
1401!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1402!    Math and Computer Science Division, Argonne National Laboratory   !
1403!BOP -------------------------------------------------------------------
1404!
1405! !IROUTINE: copyGeneralGridHeader_ - Copy the GeneralGrid Header.
1406!
1407! !DESCRIPTION:  This routine copies the header information from the
1408! input {\tt GeneralGrid} argument {\tt iGGrid} to the output
1409! {\tt GeneralGrid} argument {\tt oGGrid}.  The {\em header information}
1410! in a {\tt GeneralGrid} variable comprises all the non-{\tt AttrVect}
1411! components of the {\tt GeneralGrid}; that is, everything except the
1412! gridpoint coordinate, geometry, and index data stored in
1413! {\tt iGGrid\%data}.  This information includes:
1414! \begin{enumerate}
1415! \item The coordinates in {\tt iGGrid\%coordinate\_list}
1416! \item The coordinate sort order in {\tt iGGrid\%coordinate\_sort\_order}
1417! \item The area/volume weights in {\tt iGGrid\%weight\_list}
1418! \item Other {\tt REAL} geometric information in {\tt iGGrid\%other\_list}
1419! \item Indexing information in {\tt iGGrid\%index\_list}
1420! \item The {\tt LOGICAL} descending/ascending order sort flags in
1421! {\tt iGGrid\%descend(:)}.
1422! \end{enumerate}
1423!
1424! !INTERFACE:
1425
1426 subroutine copyGeneralGridHeader_(iGGrid, oGGrid)
1427!
1428! !USES:
1429!
1430      use m_stdio
1431      use m_die
1432
1433      use m_List, only : List
1434      use m_List, only : List_copy => copy
1435      use m_List, only : List_allocated => allocated
1436      use m_List, only : List_nullify => nullify
1437
1438      use m_GeneralGrid, only : GeneralGrid
1439
1440      implicit none
1441
1442! !INPUT PARAMETERS:
1443!
1444      type(GeneralGrid), intent(in)  :: iGGrid
1445
1446! !OUTPUT PARAMETERS:
1447!
1448      type(GeneralGrid), intent(out) :: oGGrid
1449
1450! !REVISION HISTORY:
1451!       05Jun01 - J.W. Larson <larson@mcs.anl.gov> - Initial code.
1452!       08Aug01 - E.T. Ong <eong@mcs.anl.gov> - changed list assignments(=)
1453!                 to list copy.
1454!       05Aug02 - E. Ong <eong@mcs.anl.gov> - added association checking
1455!EOP ___________________________________________________________________
1456
1457  character(len=*),parameter :: myname_=myname//'::copyGeneralGridHeader_'
1458
1459  logical :: DescendAssoc
1460  integer :: DescendSize, i, ierr
1461
1462       ! Step 1. Copy GeneralGrid List attributes from iGGrid
1463       ! to oGGrid.
1464
1465  call List_nullify(oGGrid%coordinate_list)
1466  call List_nullify(oGGrid%coordinate_sort_order)
1467  call List_nullify(oGGrid%weight_list)
1468  call List_nullify(oGGrid%other_list)
1469  call List_nullify(oGGrid%index_list)
1470  nullify(oGGrid%descend)
1471
1472  if(List_allocated(iGGrid%coordinate_list)) then
1473     call List_copy(oGGrid%coordinate_list,iGGrid%coordinate_list)
1474  endif
1475
1476  if(List_allocated(iGGrid%coordinate_sort_order)) then
1477     call List_copy(oGGrid%coordinate_sort_order,iGGrid%coordinate_sort_order)
1478  endif
1479
1480  if(List_allocated(iGGrid%weight_list)) then
1481     call List_copy(oGGrid%weight_list,iGGrid%weight_list)
1482  endif
1483
1484  if(List_allocated(iGGrid%other_list)) then
1485     call List_copy(oGGrid%other_list,iGGrid%other_list)
1486  endif
1487
1488  if(List_allocated(iGGrid%index_list)) then
1489     call List_copy(oGGrid%index_list,iGGrid%index_list)
1490  endif
1491
1492  DescendAssoc = associated(iGGrid%descend) 
1493  if(DescendAssoc) then
1494
1495     DescendSize = size(iGGrid%descend)
1496     allocate(oGGrid%descend(DescendSize), stat=ierr)
1497     if(ierr /= 0) then
1498        write(stderr,*) myname_,':: ERROR--allocate(iGGrid%descend(... failed.',&
1499             ' ierr = ', ierr, 'DescendSize = ', DescendSize
1500        call die(myname_)
1501     endif
1502     do i=1,DescendSize
1503        oGGrid%descend(i) = iGGrid%descend(i)
1504     end do
1505
1506  endif
1507
1508       ! The GeneralGrid header copy is now complete.
1509
1510 end subroutine copyGeneralGridHeader_
1511
1512 end module m_GeneralGridComms
1513
1514
1515
1516
1517
1518
1519
1520
1521
1522
1523
1524
1525
1526
1527
1528
1529
1530
1531
1532
1533
1534
1535
1536
Note: See TracBrowser for help on using the repository browser.