source: CPL/oasis3-mct/branches/OASIS3-MCT_2.0_branch/lib/mct/mct/m_GlobalSegMapComms.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: 18.4 KB
Line 
1!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2!    Math and Computer Science Division, Argonne National Laboratory   !
3!-----------------------------------------------------------------------
4! CVS m_GlobalSegMapComms.F90,v 1.6 2004-04-23 23:58:47 jacob Exp
5! CVS MCT_2_8_0
6!BOP -------------------------------------------------------------------
7!
8! !MODULE: m_GlobalSegMapComms - GlobalSegMap Communications Support
9!
10! !DESCRIPTION:
11!
12! This module provides communications support for the {\tt GlobalSegMap}
13! datatype.  Both blocking and non-blocking point-to-point communications
14! are provided for send (analogues to {\tt MPI\_SEND()/MPI\_ISEND()})
15! A receive and broadcast method is also supplied.
16!
17! !INTERFACE:
18
19 module m_GlobalSegMapComms
20
21      implicit none
22
23      private   ! except
24
25! !PUBLIC MEMBER FUNCTIONS:
26
27      public :: send
28      public :: recv
29      public :: isend
30      public :: bcast
31
32      interface bcast ; module procedure bcast_ ; end interface
33      interface send  ; module procedure send_  ; end interface
34      interface recv  ; module procedure recv_  ; end interface
35      interface isend ; module procedure isend_ ; end interface
36
37! !REVISION HISTORY:
38! 11Aug03 - J.W. Larson <larson@mcs.anl.gov> - initial version
39!
40!EOP ___________________________________________________________________
41
42  character(len=*),parameter :: myname='MCT::m_GlobalSegMapComms'
43
44 contains
45
46!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
47!    Math and Computer Science Division, Argonne National Laboratory   !
48!BOP -------------------------------------------------------------------
49!
50! !IROUTINE: send_ - Point-to-point blocking Send of a GlobalSegMap
51!
52! !DESCRIPTION:
53! This routine performs a blocking send of a {\tt GlobalSegMap} (the
54! input argument {\tt outgoingGSMap}) to the root processor on component
55! {\tt comp\_id}. The input {\tt INTEGER} argument {\tt TagBase}
56! is used to generate tags for the messages associated with this operation;
57! there are six messages involved, so the user should avoid using tag
58! values {\tt TagBase} and {\tt TagBase + 5}.  All six messages are blocking.
59! The success (failure) of this operation is reported in the zero
60! (non-zero) value of the optional {\tt INTEGER} output variable {\tt status}.
61!
62! !INTERFACE:
63
64 subroutine send_(outgoingGSMap, comp_id, TagBase, status)
65
66!
67! !USES:
68!
69      use m_mpif90
70      use m_die, only : MP_perr_die,die
71      use m_stdio
72
73      use m_GlobalSegMap, only : GlobalSegMap
74      use m_GlobalSegMap, only : GlobalSegMap_ngseg => ngseg
75      use m_GlobalSegMap, only : GlobalSegMap_comp_id => comp_ID
76      use m_GlobalSegMap, only : GlobalSegMap_gsize => gsize
77
78      use m_MCTWorld, only : ComponentToWorldRank
79      use m_MCTWorld, only : ThisMCTWorld
80
81      implicit none
82
83! !INPUT PARAMETERS:
84
85  type(GlobalSegMap),    intent(IN)  :: outgoingGSMap
86  integer,               intent(IN)  :: comp_id
87  integer,               intent(IN)  :: TagBase
88
89! !OUTPUT PARAMETERS:
90
91  integer, optional,     intent(OUT) :: status
92
93! !REVISION HISTORY:
94! 13Aug03 - J.W. Larson <larson@mcs.anl.gov> - API and initial version.
95! 26Aug03 - R. Jacob <jacob@mcs.anl.gov> - use same method as isend_
96! 05Mar04 - R. Jacob <jacob@mcs.anl.gov> - match new isend_ method.
97!EOP ___________________________________________________________________
98  character(len=*),parameter :: myname_=myname//'::send_'
99
100  integer :: ierr
101  integer :: destID
102  integer :: nsegs
103
104  if(present(status)) status = 0 ! the success value
105
106  destID = ComponentToWorldRank(0, comp_id, ThisMCTWorld)
107
108       ! Next, send the buffer size to destID so it can prepare a
109       ! receive buffer of the correct size.
110  nsegs = GlobalSegMap_ngseg(outgoingGSMap)
111
112  call MPI_SEND(outgoingGSMap%comp_id, 1, MP_Type(outgoingGSMap%comp_id), destID, &
113                TagBase, ThisMCTWorld%MCT_comm, ierr)
114  if(ierr /= 0) then
115     call MP_perr_die(myname_, 'Send compid failed',ierr)
116  endif
117
118  call MPI_SEND(outgoingGSMap%ngseg, 1, MP_Type(outgoingGSMap%ngseg), destID, &
119                TagBase+1, ThisMCTWorld%MCT_comm, ierr)
120  if(ierr /= 0) then
121     call MP_perr_die(myname_, 'Send ngseg failed',ierr)
122  endif
123
124  call MPI_SEND(outgoingGSMap%gsize, 1, MP_Type(outgoingGSMap%gsize), destID, &
125                TagBase+2, ThisMCTWorld%MCT_comm, ierr)
126  if(ierr /= 0) then
127     call MP_perr_die(myname_, 'Send gsize failed',ierr)
128  endif
129
130
131       ! Send segment information data (3 messages)
132
133  call MPI_SEND(outgoingGSMap%start, nsegs, &
134                MP_Type(outgoingGSMap%start(1)), &
135                destID, TagBase+3, ThisMCTWorld%MCT_comm, ierr)
136  if(ierr /= 0) then
137     call MP_perr_die(myname_, 'Send outgoingGSMap%start failed',ierr)
138  endif
139
140  call MPI_SEND(outgoingGSMap%length, nsegs, &
141                MP_Type(outgoingGSMap%length(1)), &
142                destID, TagBase+4, ThisMCTWorld%MCT_comm, ierr)
143  if(ierr /= 0) then
144     call MP_perr_die(myname_, 'Send outgoingGSMap%length failed',ierr)
145  endif
146
147  call MPI_SEND(outgoingGSMap%pe_loc, nsegs, &
148                MP_Type(outgoingGSMap%pe_loc(1)), &
149                destID, TagBase+5, ThisMCTWorld%MCT_comm, ierr)
150  if(ierr /= 0) then
151     call MP_perr_die(myname_, 'Send outgoingGSMap%pe_loc failed',ierr)
152  endif
153
154 end subroutine send_
155
156!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
157!    Math and Computer Science Division, Argonne National Laboratory   !
158!BOP -------------------------------------------------------------------
159!
160! !IROUTINE: isend_ - Point-to-point Non-blocking Send of a GlobalSegMap
161!
162! !DESCRIPTION:
163! This routine performs a non-blocking send of a {\tt GlobalSegMap} (the
164! input argument {\tt outgoingGSMap}) to the root processor on component
165! {\tt comp\_id}  The input {\tt INTEGER} argument {\tt TagBase}
166! is used to generate tags for the messages associated with this operation;
167! there are six messages involved, so the user should avoid using tag
168! values {\tt TagBase} and {\tt TagBase + 5}.  All six messages are non-
169! blocking, and the request handles for them are returned in the output
170! {\tt INTEGER} array {\tt reqHandle}, which can be checked for completion
171! using any of MPI's wait functions.  The success (failure) of
172! this operation is reported in the zero (non-zero) value of the optional
173! {\tt INTEGER} output variable {\tt status}.
174!
175! {\bf N.B.}:  Data is sent directly out of {\tt outgoingGSMap} so it
176! must not be deleted until the send has completed.
177!
178! {\bf N.B.}:  The array {\tt reqHandle} represents allocated memory that
179! must be deallocated when it is no longer needed.  Failure to do so will
180! create a memory leak.
181!
182! !INTERFACE:
183
184 subroutine isend_(outgoingGSMap, comp_id, TagBase, reqHandle, status)
185
186!
187! !USES:
188!
189      use m_mpif90
190      use m_die, only : MP_perr_die,die
191      use m_stdio
192
193      use m_GlobalSegMap, only : GlobalSegMap
194      use m_GlobalSegMap, only : GlobalSegMap_ngseg => ngseg
195
196      use m_MCTWorld, only : ComponentToWorldRank
197      use m_MCTWorld, only : ThisMCTWorld
198
199      implicit none
200
201! !INPUT PARAMETERS:
202
203  type(GlobalSegMap),    intent(IN)  :: outgoingGSMap
204  integer,               intent(IN)  :: comp_id
205  integer,               intent(IN)  :: TagBase
206
207! !OUTPUT PARAMETERS:
208
209  integer, dimension(:), pointer     :: reqHandle
210  integer, optional,     intent(OUT) :: status
211
212! !REVISION HISTORY:
213! 13Aug03 - J.W. Larson <larson@mcs.anl.gov> - API and initial version.
214! 05Mar04 - R. Jacob <jacob@mcs.anl.gov> - Send everything directly out
215!           of input GSMap.  Don't use a SendBuffer.
216!
217!EOP ___________________________________________________________________
218
219  character(len=*),parameter :: myname_=myname//'::isend_'
220
221  integer :: ierr,destID,nsegs
222
223  if(present(status)) status = 0 ! the success value
224
225  destID = ComponentToWorldRank(0, comp_id, ThisMCTWorld)
226
227  allocate(reqHandle(6), stat=ierr)
228  if(ierr /= 0) then
229     write(stderr,'(2a,i8)') myname_, &
230          'FATAL--allocation of send buffer failed with ierr=',ierr
231     call die(myname_)
232  endif
233
234       ! Next, send the buffer size to destID so it can prepare a
235       ! receive buffer of the correct size (3 messages).
236  nsegs = GlobalSegMap_ngseg(outgoingGSMap)
237
238  call MPI_ISEND(outgoingGSMap%comp_id, 1, MP_Type(outgoingGSMap%comp_id), destID, &
239                TagBase, ThisMCTWorld%MCT_comm, reqHandle(1), ierr)
240  if(ierr /= 0) then
241     call MP_perr_die(myname_, 'Send compid failed',ierr)
242  endif
243
244  call MPI_ISEND(outgoingGSMap%ngseg, 1, MP_Type(outgoingGSMap%ngseg), destID, &
245                TagBase+1, ThisMCTWorld%MCT_comm, reqHandle(2), ierr)
246  if(ierr /= 0) then
247     call MP_perr_die(myname_, 'Send ngseg failed',ierr)
248  endif
249
250  call MPI_ISEND(outgoingGSMap%gsize, 1, MP_Type(outgoingGSMap%gsize), destID, &
251                TagBase+2, ThisMCTWorld%MCT_comm, reqHandle(3), ierr)
252  if(ierr /= 0) then
253     call MP_perr_die(myname_, 'Send gsize failed',ierr)
254  endif
255
256       ! Send segment information data (3 messages)
257
258  call MPI_ISEND(outgoingGSMap%start, nsegs, &
259                MP_Type(outgoingGSMap%start(1)), &
260                destID, TagBase+3, ThisMCTWorld%MCT_comm, reqHandle(4), ierr)
261  if(ierr /= 0) then
262     call MP_perr_die(myname_, 'Send outgoingGSMap%start failed',ierr)
263  endif
264
265  call MPI_ISEND(outgoingGSMap%length, nsegs, &
266                MP_Type(outgoingGSMap%length(1)), &
267                destID, TagBase+4, ThisMCTWorld%MCT_comm, reqHandle(5), ierr)
268  if(ierr /= 0) then
269     call MP_perr_die(myname_, 'Send outgoingGSMap%length failed',ierr)
270  endif
271
272  call MPI_ISEND(outgoingGSMap%pe_loc, nsegs, &
273                MP_Type(outgoingGSMap%pe_loc(1)), &
274                destID, TagBase+5, ThisMCTWorld%MCT_comm, reqHandle(6), ierr)
275  if(ierr /= 0) then
276     call MP_perr_die(myname_, 'Send outgoingGSMap%pe_loc failed',ierr)
277  endif
278
279 end subroutine isend_
280
281!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
282!    Math and Computer Science Division, Argonne National Laboratory   !
283!BOP -------------------------------------------------------------------
284!
285! !IROUTINE: recv_ - Point-to-point blocking Receive of a GlobalSegMap
286!
287! !DESCRIPTION:
288! This routine performs a blocking receive of a {\tt GlobalSegMap} (the
289! input argument {\tt outgoingGSMap}) from the root processor on component
290! {\tt comp\_id}. The input {\tt INTEGER} argument {\tt TagBase}
291! is used to generate tags for the messages associated with this operation;
292! there are six messages involved, so the user should avoid using tag
293! values {\tt TagBase} and {\tt TagBase + 5}. The success (failure) of this
294! operation is reported in the zero (non-zero) value of the optional {\tt INTEGER}
295! output variable {\tt status}.
296!
297! !INTERFACE:
298
299 subroutine recv_(incomingGSMap, comp_id, TagBase, status)
300
301!
302! !USES:
303!
304      use m_mpif90
305      use m_die, only : MP_perr_die, die
306      use m_stdio
307
308      use m_GlobalSegMap, only : GlobalSegMap
309      use m_GlobalSegMap, only : GlobalSegMap_init => init
310
311      use m_MCTWorld, only : ComponentToWorldRank
312      use m_MCTWorld, only : ThisMCTWorld
313
314      implicit none
315
316! !INPUT PARAMETERS:
317
318  integer,               intent(IN)  :: comp_id
319  integer,               intent(IN)  :: TagBase
320
321! !OUTPUT PARAMETERS:
322
323  type(GlobalSegMap),    intent(OUT) :: incomingGSMap
324  integer, optional,     intent(OUT) :: status
325
326! !REVISION HISTORY:
327! 13Aug03 - J.W. Larson <larson@mcs.anl.gov> - API and initial version.
328! 25Aug03 - R.Jacob <larson@mcs.anl.gov> - rename to recv_.
329!EOP ___________________________________________________________________
330
331  character(len=*),parameter :: myname_=myname//'::recv_'
332
333  integer :: ierr,sourceID
334  integer :: MPstatus(MP_STATUS_SIZE)
335  integer :: RecvBuffer(3)
336
337  if(present(status)) status = 0 ! the success value
338
339  sourceID = ComponentToWorldRank(0, comp_id, ThisMCTWorld)
340
341       ! Receive the GlobalSegMap's basic constants:  component id,
342       ! grid size, and number of segments.  The number of segments
343       ! is needed to construct the arrays into which segment
344       ! information will be received.  Thus, this receive blocks.
345
346  call MPI_RECV(RecvBuffer(1), 1, MP_Type(RecvBuffer(1)), sourceID, &
347                TagBase, ThisMCTWorld%MCT_comm, MPstatus, ierr)
348  if(ierr /= 0) then
349     call MP_perr_die(myname_, 'Receive of compid failed',ierr)
350  endif
351  call MPI_RECV(RecvBuffer(2), 1, MP_Type(RecvBuffer(2)), sourceID, &
352                TagBase+1, ThisMCTWorld%MCT_comm, MPstatus, ierr)
353  if(ierr /= 0) then
354     call MP_perr_die(myname_, 'Receive of ngseg failed',ierr)
355  endif
356  call MPI_RECV(RecvBuffer(3), 1, MP_Type(RecvBuffer(3)), sourceID, &
357                TagBase+2, ThisMCTWorld%MCT_comm, MPstatus, ierr)
358  if(ierr /= 0) then
359     call MP_perr_die(myname_, 'Receive of gsize failed',ierr)
360  endif
361
362       ! Create Empty GlobaSegMap into which segment information
363       ! will be received
364
365  call GlobalSegMap_init(incomingGSMap, RecvBuffer(1), RecvBuffer(2), &
366                         RecvBuffer(3))
367
368       ! Receive segment information data (3 messages)
369
370  call MPI_RECV(incomingGSMap%start, RecvBuffer(2), &
371                MP_Type(incomingGSMap%start(1)), &
372                sourceID, TagBase+3, ThisMCTWorld%MCT_comm, MPstatus, ierr)
373  if(ierr /= 0) then
374     call MP_perr_die(myname_, 'Recv incomingGSMap%start failed',ierr)
375  endif
376
377  call MPI_RECV(incomingGSMap%length, RecvBuffer(2), &
378                MP_Type(incomingGSMap%length(1)), &
379                sourceID, TagBase+4, ThisMCTWorld%MCT_comm, MPstatus, ierr)
380  if(ierr /= 0) then
381     call MP_perr_die(myname_, 'Recv incomingGSMap%length failed',ierr)
382  endif
383
384  call MPI_RECV(incomingGSMap%pe_loc, RecvBuffer(2), &
385                MP_Type(incomingGSMap%pe_loc(1)), &
386                sourceID, TagBase+5, ThisMCTWorld%MCT_comm, MPstatus, ierr)
387  if(ierr /= 0) then
388     call MP_perr_die(myname_, 'Recv incomingGSMap%pe_loc failed',ierr)
389  endif
390
391 end subroutine recv_
392
393!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
394!    Math and Computer Science Division, Argonne National Laboratory   !
395!BOP -------------------------------------------------------------------
396!
397! !IROUTINE: bcast_ - broadcast a GlobalSegMap object
398!
399! !DESCRIPTION:
400!
401! The routine {\tt bcast\_()} takes the input/output {\em GlobalSegMap}
402! argument {\tt GSMap} (on input valid only on the {\tt root} process,
403! on output valid on all processes) and broadcasts it to all processes
404! on the communicator associated with the F90 handle {\tt comm}.  The
405! success (failure) of this operation is returned as a zero (non-zero)
406! value of the optional output {\tt INTEGER} argument {\tt status}.
407!
408! !INTERFACE:
409
410 subroutine bcast_(GSMap, root, comm, status)
411
412!
413! !USES:
414!
415      use m_mpif90
416      use m_die, only : MP_perr_die,die
417      use m_stdio
418
419      use m_GlobalSegMap, only : GlobalSegMap
420
421      implicit none
422
423! !INPUT PARAMETERS:
424
425      integer,            intent(in)     :: root
426      integer,            intent(in)     :: comm
427
428! !INPUT/OUTPUT PARAMETERS:
429
430      type(GlobalSegMap), intent(inout)  :: GSMap  ! Output GlobalSegMap
431
432! !OUTPUT PARAMETERS:
433
434      integer, optional,  intent(out)    :: status ! global vector size
435
436! !REVISION HISTORY:
437! 17Oct01 - J.W. Larson <larson@mcs.anl.gov> - Initial version.
438! 11Aug03 - J.W. Larson <larson@mcs.anl.gov> - Relocated from original
439!           location in m_GlobalSegMap.
440!EOP ___________________________________________________________________
441
442  character(len=*),parameter :: myname_=myname//'::bcast_'
443
444  integer :: myID, ierr, n
445  integer, dimension(:), allocatable :: IntBuffer
446
447       ! Step One:  which process am I?
448
449  call MP_COMM_RANK(comm, myID, ierr)
450  if(ierr /= 0) call MP_perr_die(myname_,'MP_comm_rank()',ierr)
451
452       ! Step Two:  Broadcast the scalar bits of the GlobalSegMap from
453       ! the root.
454
455  allocate(IntBuffer(3), stat=ierr) ! allocate buffer space (all PEs)
456  if(ierr /= 0) then
457    if(.not. present(status)) then
458       call die(myname_,'allocate(IntBuffer)',ierr)
459    else
460       write(stderr,*) myname_,':: error during allocate(IntBuffer)'
461       status = 2
462       return
463    endif
464  endif
465
466  if(myID == root) then ! pack the buffer
467     IntBuffer(1) = GSMap%comp_id
468     IntBuffer(2) = GSMap%ngseg
469     IntBuffer(3) = GSMap%gsize
470  endif
471
472  call MPI_BCAST(IntBuffer, 3, MP_type(IntBuffer(1)), root, comm, ierr)
473  if(ierr /= 0) call MP_perr_die(myname_,'MPI_BCAST(IntBuffer)',ierr)
474
475  if(myID /= root) then ! unpack from buffer to GSMap
476     GSMap%comp_id = IntBuffer(1)
477     GSMap%ngseg = IntBuffer(2)
478     GSMap%gsize = IntBuffer(3)
479  endif
480
481  deallocate(IntBuffer, stat=ierr) ! deallocate buffer space
482  if(ierr /= 0) then
483    if(.not. present(status)) then
484       call die(myname_,'deallocate(IntBuffer)',ierr)
485    else
486       write(stderr,*) myname_,':: error during deallocate(IntBuffer)'
487       status = 4
488       return
489    endif
490  endif
491
492       ! Step Three:  Broadcast the vector bits of GSMap from the root.
493       ! Pack them into one big array to save latency costs associated
494       ! with multiple broadcasts.
495
496  allocate(IntBuffer(3*GSMap%ngseg), stat=ierr) ! allocate buffer space (all PEs)
497  if(ierr /= 0) then
498    if(.not. present(status)) then
499       call die(myname_,'second allocate(IntBuffer)',ierr)
500    else
501       write(stderr,*) myname_,':: error during second allocate(IntBuffer)'
502       status = 5
503       return
504    endif
505  endif
506
507  if(myID == root) then ! pack outgoing broadcast buffer
508     do n=1,GSMap%ngseg
509        IntBuffer(n) = GSMap%start(n)
510        IntBuffer(GSMap%ngseg+n) = GSMap%length(n)
511        IntBuffer(2*GSMap%ngseg+n) = GSMap%pe_loc(n)
512     end do
513  endif
514
515  call MPI_BCAST(IntBuffer, 3*GSMap%ngseg, MP_Type(IntBuffer(1)), root, comm, ierr)
516  if(ierr /= 0) call MP_perr_die(myname_,'Error in second MPI_BCAST(IntBuffer)',ierr)
517
518  if(myID /= root) then ! Allocate GSMap%start, GSMap%length,...and fill them
519
520     allocate(GSMap%start(GSMap%ngseg), GSMap%length(GSMap%ngseg), &
521              GSMap%pe_loc(GSMap%ngseg), stat=ierr)
522     if(ierr /= 0) then
523        if(.not. present(status)) then
524           call die(myname_,'off-root allocate(GSMap%start...)',ierr)
525        else
526           write(stderr,*) myname_,':: error during off-root allocate(GSMap%start...)'
527           status = 7
528           return
529        endif
530     endif
531
532     do n=1,GSMap%ngseg ! unpack the buffer into the GlobalSegMap
533        GSMap%start(n) = IntBuffer(n)
534        GSMap%length(n) = IntBuffer(GSMap%ngseg+n)
535        GSMap%pe_loc(n) = IntBuffer(2*GSMap%ngseg+n)
536     end do
537
538  endif
539
540       ! Clean up buffer space:
541
542  deallocate(IntBuffer, stat=ierr) 
543  if(ierr /= 0) then
544    if(.not. present(status)) then
545       call die(myname_,'second deallocate(IntBuffer)',ierr)
546    else
547       write(stderr,*) myname_,':: error during second deallocate(IntBuffer)'
548       status = 8
549       return
550    endif
551  endif
552
553 end subroutine bcast_
554
555 end module m_GlobalSegMapComms
Note: See TracBrowser for help on using the repository browser.