C $Header: /u/gcmpack/MITgcm/pkg/compon_communic/coupsend_r4tiles.F,v 1.3 2015/11/04 17:04:21 jmc Exp $
C $Name:  $

!=======================================================================
      subroutine COUPSEND_R4TILES( component, dataname, Nx, Ny, arr )
      implicit none
! Predefined constants/arrays
#include "CPLR_SIG.h"
! MPI variables
#include "mpif.h"
! Arguments
      character*(*) component
      character*(*) dataname
      integer Nx,Ny
      real*4 arr(Nx,Ny)
! Functions
      integer mitcplr_match_comp
      integer generate_tag
      external 
      external 
! Local
      integer count,dtype,dest,tag,comm,ierr
      integer compind,numprocs
      integer i,j,ij,n,bibj
      integer Ni,Io,Nj,Jo
!     ------------------------------------------------------------------

! Establish who I am communicating with
      compind=mitcplr_match_comp( component )
      if (compind.le.0) stop 'coupsend_r4tiles: Bad component id'
      comm=MPI_COMM_compcplr( compind )
      numprocs=num_component_procs(compind)
      if (numprocs.lt.1) then
       write(LogUnit,*) 'coupsend_r4tiles: compind = ',compind
       stop 'coupsend_r4tiles: numprocs < 1'
      endif
      if (VERB)
     &  write(LogUnit,*) 'coupsend_r4tiles: ',component_Name(compind)
      if (VERB)
     &  write(LogUnit,*) 'coupsend_r4tiles: dataname=',dataname

! Foreach component process
      do n=1,numprocs

! Foreach tile on that process
       do bibj=1,component_num_tiles(n,compind)

! Create header
        Io=component_tile_i0(bibj,n,compind)
        Jo=component_tile_j0(bibj,n,compind)
        Ni=component_tile_nx(bibj,n,compind)
        Nj=component_tile_ny(bibj,n,compind)
        r4buf(1)=float( Io )
        r4buf(2)=float( Jo )
        r4buf(3)=float( Ni )
        r4buf(4)=float( Nj )
        call MITCPLR_CHAR2REAL( dataname, r4buf(9) )

! Pack data
        do j=1,Nj
         do i=1,Ni
          ij=HEADER_SIZE+i+Ni*(j-1)
          r4buf(ij)=arr(Io+i-1,Jo+j-1)
         enddo
        enddo

! Send message
        count=HEADER_SIZE+Ni*Nj
        dtype=MPI_REAL
        tag=generate_tag(123,bibj,dataname)
        dest=rank_component_procs(n,compind)

        if (VERB) then
         write(LogUnit,*) 'coupsend_r4tiles: calling MPI_Send dest=',
     &     dest,' proc=',n,'/',numprocs,' tile=',bibj
         call FLUSH(LogUnit)
        endif
        call MPI_SEND( r4buf, count, dtype, dest, tag, comm, ierr )
        if (VERB) then
         write(LogUnit,*) 'coupsend_r4tiles: returned ierr=',ierr
         call FLUSH(LogUnit)
        endif

        if (ierr.ne.0) then
         write(LogUnit,*) 'coupsend_r4tiles: rank(W,G)=',
     &            my_rank_in_world,my_rank_in_global,
     &            ' ierr=',ierr
         stop 'coupsend_r4tiles: MPI_Send failed'
        endif

       enddo ! bibj

      enddo ! n

!     ------------------------------------------------------------------
      return
      end


!=======================================================================