ccm_scatter

CCM scatter

Routine:

ccm_scatter

Purpose:

Send data from one task to all other tasks in a parallel application. By default, the zeroth task in the parallel application sends the data. Xout is overwritten on all tasks.

Minimal calling sequence:

call ccm_scatter(xin,xout)

Required Arguments:

xin :: integer, real, double precision, complex, logical, character array,intent (inout)
The values in xin will be scattered from one task to all other tasks in a parallel application. The "root" or zeroth task in the parallel application sends the data.
xout :: integer, real, double precision, complex, logical, character scaler or array,intent (inout)
The values in xin will be scattered from one task into xout on all tasks in a parallel application. Xout is overwritten on all tasks.

Call with all Optional Arguments:

call ccm_scatter(xin,xout,root,the_err)
root :: integer,intent (in)
The number of the task that sends the data, defaults to the root task of the parallel application.
the_err :: integer, intent (out)
Error code 0 = success, != 0 failure.
See Specifying Optional Arguments for the syntax for using optional arguments.

Notes:

Xin and xout can be any intrinsic Fortran data type: integer, real, double precision, complex, logical or, character. Xout can be a scaler or an array. If Xout is an array then it must be the same size in all tasks. The entire xin array will be scattered with the first portion of the array going to the zeroth task and the last portion going to the last task. Each task will receive the same number of values.

First example:


program ccm_scatter_x1
    use ccm
    implicit none
    integer :: myid,numprocs,igot,i
    integer ,allocatable :: to_send(:)
    real local_time,global_time
    call ccm_init(myid,numprocs)
    allocate(to_send(0:numprocs-1))
    do i=0,numprocs-1
       to_send(i)=i*i
    enddo
    call ccm_scatter(to_send,igot)
    write(*,*)"for i= ",myid," igot= ",igot
    call ccm_close()
end program

Example output on 4 processors


[ccm_host:~/ccm/source]% ccm_scatter_x1
 for i=   0  igot=   0
 for i=   2  igot=   4
 for i=   1  igot=   1
 for i=   3  igot=   9
[ccm_host:~/ccm/source] % 

The call to ccm_init initializes the communication package. The zeroth task fills an array with i**2. These values are then scattered to the tasks, each getting a single value equal to i**2. The call to ccm_close closes the communication package.

Second Example:


program ccm_scatter_x2
	use ccm
	implicit none
	real :: i40(40),i3d(5,4,2),got(8)
	real,allocatable :: rs(:,:)
	integer:: order(3),new_shape(3)
	integer:: myid,numprocs,i,j,k
	do i=1,40
		i40(i)=i
	enddo
	i3d=reshape(i40,shape(i3d))
	call ccm_init(myid,numprocs)
	if(myid .eq. 0)then
		write(*,*)"3d array of bounds ",shape(i3d)
		do k=1,2
			write(*,*)"plane ",k
			do i=1,5
				write(*,"(4f5.0)")(i3d(i,j,k),j=1,4)
			enddo
		enddo
	endif
!
	call ccm_barrier(0.95,.true.)
	allocate(rs(2,5))
	order=(/1,2,3/)
	new_shape=(/5,4,2/)
	call ccm_scatter(reshape(i3d,shape=new_shape,order=order),rs)
	if(myid .eq. 0)write(*,"(/3i3,5x,3i3)")order,new_shape
	do i=1,2
		if(i.eq. 1)then
			write(*,"(i4,(5f5.0))")myid,(rs(i,j),j=1,5)
		else
			write(*,"(4x,(5f5.0))")(rs(i,j),j=1,5)
		endif
	enddo
	deallocate(rs)
	call ccm_barrier(0.95,.true.)
	allocate(rs(10,1))
	order=(/3,1,2/)
	new_shape=(/2,5,4/)
	call ccm_scatter(reshape(i3d,shape=new_shape,order=order),rs)
	if(myid .eq. 0)write(*,"(/3i3,5x,3i3)")order,new_shape
	write(*,"(i4,(10f5.0))")myid,(rs(i,1),i=1,10)
	call ccm_barrier(0.95,.true.)
	order=(/3,1,2/)
	new_shape=(/5,4,2/)
	call ccm_scatter(reshape(i3d,shape=new_shape,order=order),rs)
	if(myid .eq. 0)write(*,"(/3i3,5x,3i3)")order,new_shape
	write(*,"(i4,(10f5.0))")myid,(rs(i,1),i=1,10)
	deallocate(rs)
	call ccm_barrier(0.95,.true.)
	allocate(rs(5,2))
	order=(/2,3,1/)
	new_shape=(/5,4,2/)
	call ccm_scatter(reshape(i3d,shape=new_shape,order=order),rs)
	if(myid .eq. 0)write(*,"(/3i3,5x,3i3)")order,new_shape
	do j=1,2
		if(j.eq. 1)then
			write(*,"(i4,(5f5.0))")myid,(rs(i,j),i=1,5)
		else
			write(*,"(4x,(5f5.0))")(rs(i,j),i=1,5)
		endif
	enddo
	call ccm_barrier(0.95,.true.)
	deallocate(rs)
	allocate(rs(5,2))
	order=(/2,1,3/)
	new_shape=(/5,4,2/)
	call ccm_scatter(reshape(i3d,shape=new_shape,order=order),rs)
	if(myid .eq. 0)write(*,"(/3i3,5x,3i3)")order,new_shape
	do j=1,2
		if(j.eq. 1)then
			write(*,"(i4,(5f5.0))")myid,(rs(i,j),i=1,5)
		else
			write(*,"(4x,(5f5.0))")(rs(i,j),i=1,5)
		endif
	enddo
	deallocate(rs)
	call ccm_close()
end program

Example output on 4 processors


[ccm_home:~/ccm/source] % ccm_scatter_x2
 3d array of bounds   5  4  2
 plane   1
   1.   6.  11.  16.
   2.   7.  12.  17.
   3.   8.  13.  18.
   4.   9.  14.  19.
   5.  10.  15.  20.
 plane   2
  21.  26.  31.  36.
  22.  27.  32.  37.
  23.  28.  33.  38.
  24.  29.  34.  39.
  25.  30.  35.  40.

  1  2  3       5  4  2
   0   1.   3.   5.   7.   9.
       2.   4.   6.   8.  10.
   1  11.  13.  15.  17.  19.
      12.  14.  16.  18.  20.
   2  21.  23.  25.  27.  29.
      22.  24.  26.  28.  30.
   3  31.  33.  35.  37.  39.
      32.  34.  36.  38.  40.

  3  1  2       2  5  4
   0   1.   5.   9.  13.  17.  21.  25.  29.  33.  37.
   1   2.   6.  10.  14.  18.  22.  26.  30.  34.  38.
   2   3.   7.  11.  15.  19.  23.  27.  31.  35.  39.
   3   4.   8.  12.  16.  20.  24.  28.  32.  36.  40.

  3  1  2       5  4  2
   0   1.   3.   5.   7.   9.  11.  13.  15.  17.  19.
   1  21.  23.  25.  27.  29.  31.  33.  35.  37.  39.
   2   2.   4.   6.   8.  10.  12.  14.  16.  18.  20.
   3  22.  24.  26.  28.  30.  32.  34.  36.  38.  40.

  2  3  1       5  4  2
   0   1.   9.  17.  25.  33.
       2.  10.  18.  26.  34.
   1   3.  11.  19.  27.  35.
       4.  12.  20.  28.  36.
   2   5.  13.  21.  29.  37.
       6.  14.  22.  30.  38.
   3   7.  15.  23.  31.  39.
       8.  16.  24.  32.  40.

  2  1  3       5  4  2
   0   1.   5.   9.  13.  17.
       2.   6.  10.  14.  18.
   1   3.   7.  11.  15.  19.
       4.   8.  12.  16.  20.
   2  21.  25.  29.  33.  37.
      22.  26.  30.  34.  38.
   3  23.  27.  31.  35.  39.
      24.  28.  32.  36.  40.
[ccm_home:~/ccm/source] % 

The call to ccm_init initializes the communication package. The program starts with the root task holding a 40 element 3d array that contains the values 1-40. Scatter is used along the Fortran 90 reshape routine to produce various spreads of the data to processors. Task zero prints the order and shape that are input to the reshape function. All tasks print the received data. This example shows that many different distributions are possible using only the scatter and reshape routines. The type of spread that is best is dependent on the numerical algorithm. Note: you will only see this exact output on a machine that preserves ordering of output values. This is not practical on a real machine. The values shown above were rearranged to the "proper" order.



Error conditions:


If the error checking level is set to ccm_checksize the following error condition(s) are checked: These conditions may lead to a Underling communications error if not detected.
Back to API and user's guide