File indexing completed on 2018-03-02 18:38:29 UTC
view on githubraw file Latest commit add29e06 on 2018-01-31 20:35:05 UTC
bf4be02920 Jean*0001
0002 subroutine coupsend_r4tiles( component, dataname, Nx, Ny, arr )
0003 implicit none
fbcf25275b Jean*0004
0005 #include "CPLR_SIG.h"
0006
0007 #include "mpif.h"
bf4be02920 Jean*0008
0009 character*(*) component
0010 character*(*) dataname
0011 integer Nx,Ny
0012 real*4 arr(Nx,Ny)
0013
0014 integer mitcplr_match_comp
0015 integer generate_tag
fbcf25275b Jean*0016 external mitcplr_match_comp
0017 external generate_tag
bf4be02920 Jean*0018
fbcf25275b Jean*0019 integer count,dtype,dest,tag,comm,ierr
bf4be02920 Jean*0020 integer compind,numprocs
0021 integer i,j,ij,n,bibj
0022 integer Ni,Io,Nj,Jo
0023
0024
0025
0026 compind=mitcplr_match_comp( component )
0027 if (compind.le.0) stop 'coupsend_r4tiles: Bad component id'
0028 comm=MPI_COMM_compcplr( compind )
0029 numprocs=num_component_procs(compind)
0030 if (numprocs.lt.1) then
0031 write(LogUnit,*) 'coupsend_r4tiles: compind = ',compind
0032 stop 'coupsend_r4tiles: numprocs < 1'
0033 endif
0034 if (VERB)
0035 & write(LogUnit,*) 'coupsend_r4tiles: ',component_Name(compind)
0036 if (VERB)
0037 & write(LogUnit,*) 'coupsend_r4tiles: dataname=',dataname
0038
0039
0040 do n=1,numprocs
0041
0042
0043 do bibj=1,component_num_tiles(n,compind)
0044
0045
0046 Io=component_tile_i0(bibj,n,compind)
0047 Jo=component_tile_j0(bibj,n,compind)
0048 Ni=component_tile_nx(bibj,n,compind)
0049 Nj=component_tile_ny(bibj,n,compind)
0050 r4buf(1)=float( Io )
0051 r4buf(2)=float( Jo )
0052 r4buf(3)=float( Ni )
0053 r4buf(4)=float( Nj )
0054 call mitcplr_char2real( dataname, r4buf(9) )
0055
0056
0057 do j=1,Nj
0058 do i=1,Ni
0059 ij=HEADER_SIZE+i+Ni*(j-1)
0060 r4buf(ij)=arr(Io+i-1,Jo+j-1)
0061 enddo
0062 enddo
0063
0064
0065 count=HEADER_SIZE+Ni*Nj
0066 dtype=MPI_REAL
0067 tag=generate_tag(123,bibj,dataname)
0068 dest=rank_component_procs(n,compind)
0069
0070 if (VERB) then
0071 write(LogUnit,*) 'coupsend_r4tiles: calling MPI_Send dest=',
0072 & dest,' proc=',n,'/',numprocs,' tile=',bibj
0073 call flush(LogUnit)
0074 endif
0075 call MPI_Send( r4buf, count, dtype, dest, tag, comm, ierr )
0076 if (VERB) then
0077 write(LogUnit,*) 'coupsend_r4tiles: returned ierr=',ierr
0078 call flush(LogUnit)
0079 endif
0080
0081 if (ierr.ne.0) then
0082 write(LogUnit,*) 'coupsend_r4tiles: rank(W,G)=',
0083 & my_rank_in_world,my_rank_in_global,
0084 & ' ierr=',ierr
fbcf25275b Jean*0085 stop 'coupsend_r4tiles: MPI_Send failed'
bf4be02920 Jean*0086 endif
0087
0088 enddo
0089
0090 enddo
0091
0092
0093 return
0094 end
0095