Back to home page

MITgcm

 
 

    


File indexing completed on 2018-03-02 18:38:28 UTC

view on githubraw file Latest commit add29e06 on 2018-01-31 20:35:05 UTC
bf4be02920 Jean*0001 !=======================================================================
                0002       subroutine compsend_r4tiles( dataname,Ni,Oi,Nj,Oj,Nk,Ti,Tj, arr )
                0003       implicit none
                0004 ! Arguments
                0005       character*(*) dataname
                0006       integer Ni,Oi,Nj,Oj,Nk,Ti,Tj
                0007       real*4 arr(1-Oi:Ni+Oi,1-Oj:Nj+Oj,Nk,Ti,Tj)
                0008 ! Predefined constants/arrays
                0009 #include "CPLR_SIG.h"
                0010 ! MPI variables
                0011 #include "mpif.h"
                0012       integer count,datatype,dest,tag,comm,ierr
                0013 ! Functions
                0014       integer generate_tag
                0015 ! Local
                0016       integer i,j,ij,bi,bj,k,bibj
                0017 !     ------------------------------------------------------------------
                0018 
                0019       if (HEADER_SIZE+Ni*Nj.gt.MAX_R4_BUFLEN)
                0020      &    stop 'compsend_r4tiles: Nx*Ny too big'
                0021 
                0022 ! Foreach tile which is non-blank
                0023       do bibj=1,my_num_tiles
                0024 
                0025        bi=my_tile_bi(bibj)
                0026        bj=my_tile_bj(bibj)
                0027 
                0028        if (Ni.ne.my_tile_nx(bibj))
                0029      &     stop 'compsend_r4tiles: ni != my_tile_nx'
                0030        if (Nj.ne.my_tile_ny(bibj))
                0031      &     stop 'compsend_r4tiles: nj != my_tile_ny'
                0032 
                0033 ! Set up buffer
                0034        r4buf(1)=float( my_tile_i0(bibj) )
                0035        r4buf(2)=float( my_tile_nx(bibj) )
                0036        r4buf(3)=float( my_tile_j0(bibj) )
                0037        r4buf(4)=float( my_tile_ny(bibj) )
                0038        call mitcplr_char2real( dataname, r4buf(9) )
                0039 
                0040 ! Copy interior of tile to buffer
                0041        k=1
                0042        do j=1,Nj
                0043         do i=1,Ni
                0044          ij=HEADER_SIZE+i+Ni*(j-1)
                0045          r4buf(ij)=arr(i,j,k,bi,bj)
                0046         enddo
                0047        enddo
                0048 
                0049 ! Send message
                0050        count=HEADER_SIZE+Ni*Nj
                0051        datatype=MPI_REAL
                0052        dest=my_coupler_rank
                0053        tag=generate_tag(103,bibj,dataname)
                0054        comm=MPI_COMM_myglobal
                0055 
                0056        if (VERB) then
                0057         write(LogUnit,*)
                0058      &   'compsend_r4tiles: calling MPI_Send dest=',dest,' tile=',bibj
                0059         write(LogUnit,*) 'compsend_r4tiles: dataname=',dataname
                0060         call flush(LogUnit)
                0061        endif
                0062        call MPI_Send( r4buf, count, datatype, dest, tag, comm, ierr )
                0063        if (VERB) then
                0064         write(LogUnit,*) 'compsend_r4tiles: returned ierr=',ierr
                0065         call flush(LogUnit)
                0066        endif
                0067 
                0068        if (ierr.ne.0) then
                0069         write(LogUnit,*) 'compsend_r4tiles: rank(W,G,L)=',
                0070      &            my_rank_in_world,my_rank_in_global,my_rank_in_local,
                0071      &            ' ierr=',ierr
                0072         stop 'compsend_r4tiles: MPI_Send failed'
                0073        endif
                0074 
                0075       enddo ! bibj
                0076 
                0077 !     ------------------------------------------------------------------
                0078       return
                0079       end
                0080 !=======================================================================