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_r8tiles( dataname,Ni,Oi,Nj,Oj,Nk,Ti,Tj, arr )
0003 implicit none
0004
0005 character*(*) dataname
0006 integer Ni,Oi,Nj,Oj,Nk,Ti,Tj
0007 real*8 arr(1-Oi:Ni+Oi,1-Oj:Nj+Oj,Nk,Ti,Tj)
0008
0009 #include "CPLR_SIG.h"
0010
0011 #include "mpif.h"
0012 integer count,datatype,dest,tag,comm,ierr
0013
0014 integer generate_tag
0015
0016 integer i,j,ij,bi,bj,k,bibj
0017
0018
0019 if (HEADER_SIZE+Ni*Nj.gt.MAX_R8_BUFLEN)
0020 & stop 'compsend_r8tiles: Nx*Ny too big'
0021
0022
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_r8tiles: ni != my_tile_nx'
0030 if (Nj.ne.my_tile_ny(bibj))
0031 & stop 'compsend_r8tiles: nj != my_tile_ny'
0032
0033
0034 r8buf(1)=float( my_tile_i0(bibj) )
0035 r8buf(2)=float( my_tile_nx(bibj) )
0036 r8buf(3)=float( my_tile_j0(bibj) )
0037 r8buf(4)=float( my_tile_ny(bibj) )
e5266ce10c Jean*0038 call mitcplr_char2dbl( dataname, r8buf(9) )
bf4be02920 Jean*0039
0040
0041 k=1
0042 do j=1,Nj
0043 do i=1,Ni
0044 ij=HEADER_SIZE+i+Ni*(j-1)
0045 r8buf(ij)=arr(i,j,k,bi,bj)
0046 enddo
0047 enddo
0048
0049
0050 count=HEADER_SIZE+Ni*Nj
0051 datatype=MPI_DOUBLE_PRECISION
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_r8tiles: calling MPI_Send dest=',dest,' tile=',bibj
0059 write(LogUnit,*) 'compsend_r8tiles: dataname=',dataname
0060 call flush(LogUnit)
0061 endif
0062 call MPI_Send( r8buf, count, datatype, dest, tag, comm, ierr )
0063 if (VERB) then
0064 write(LogUnit,*) 'compsend_r8tiles: returned ierr=',ierr
0065 call flush(LogUnit)
0066 endif
0067
0068 if (ierr.ne.0) then
0069 write(LogUnit,*) 'compsend_r8tiles: rank(W,G,L)=',
0070 & my_rank_in_world,my_rank_in_global,my_rank_in_local,
0071 & ' ierr=',ierr
0072 stop 'compsend_r8tiles: MPI_Send failed'
0073 endif
0074
0075 enddo
0076
0077
0078 return
0079 end
0080