File indexing completed on 2018-03-02 18:38:30 UTC
view on githubraw file Latest commit add29e06 on 2018-01-31 20:35:05 UTC
bf4be02920 Jean*0001
44feb30e11 Jean*0002 subroutine MITCOMPONENT_tile_register( num_tiles, iReg )
bf4be02920 Jean*0003 implicit none
0004
0005
0006 #include "mpif.h"
44feb30e11 Jean*0007
bf4be02920 Jean*0008
0009
0010 #include "CPLR_SIG.h"
0011
44feb30e11 Jean*0012
0013 integer num_tiles
0014 integer iReg(6,num_tiles)
0015
bf4be02920 Jean*0016
0017 integer mitcplr_match_comp
0018 integer generate_tag
44feb30e11 Jean*0019 external mitcplr_match_comp
0020 external generate_tag
bf4be02920 Jean*0021
0022
44feb30e11 Jean*0023
0024 integer n
0025 integer count, datatype, dest, tag, comm, ierr
0026 integer ibuf(MAX_IBUF)
bf4be02920 Jean*0027
0028
0029
44feb30e11 Jean*0030 write(LogUnit,'(A,I6,A,2I4)')
0031 & 'MITCOMPONENT_tile_register (pId=', my_rank_in_local,
0032 & '): Starts ; num_tiles=', num_tiles
bf4be02920 Jean*0033
0034 if (num_tiles.lt.1)
44feb30e11 Jean*0035 & STOP 'MITCOMPONENT_tile_register: num_tiles < 1'
bf4be02920 Jean*0036 if (num_tiles.gt.MAX_TILES)
44feb30e11 Jean*0037 & STOP 'MITCOMPONENT_tile_register: num_tiles > MAX_TILES'
0038
0039 my_num_tiles = num_tiles
0040 do n=1,num_tiles
0041 my_tile_bi(n) = iReg(1,n)
0042 my_tile_bj(n) = iReg(2,n)
0043 my_tile_nx(n) = iReg(3,n)
0044 my_tile_ny(n) = iReg(4,n)
0045 my_tile_i0(n) = iReg(5,n)
0046 my_tile_j0(n) = iReg(6,n)
0047 write(LogUnit,'(A,I5,A,2I4,A,2I5,A,2I8)')
0048 & ' tile #', n,
0049 & ' ; bi,bj=', iReg(1,n), iReg(2,n),
0050 & ' ; Ni,Nj=', iReg(3,n), iReg(4,n),
0051 & ' ; Io,Jo=', iReg(5,n), iReg(6,n)
0052 enddo
bf4be02920 Jean*0053
0054
44feb30e11 Jean*0055 ibuf(1) = num_tiles
bf4be02920 Jean*0056
0057
0058 count=1
0059 datatype=MPI_INTEGER
0060 dest=my_coupler_rank
0061 tag=generate_tag(112,my_rank_in_global,'Register Tiles')
0062 comm=MPI_COMM_myglobal
0063
0064 call MPI_Send( ibuf, count, datatype, dest, tag, comm, ierr )
951926fb9b Jean*0065
bf4be02920 Jean*0066 if (ierr.ne.0) then
0067 write(LogUnit,*) 'MITCOMPONENT_tile_register: rank(W,G,L)=',
0068 & my_rank_in_world,my_rank_in_global,my_rank_in_local,
0069 & ' ierr=',ierr
44feb30e11 Jean*0070 STOP 'MITCOMPONENT_tile_register: MPI_Send failed'
bf4be02920 Jean*0071 endif
0072
44feb30e11 Jean*0073 do n=1,my_num_tiles
bf4be02920 Jean*0074
0075
44feb30e11 Jean*0076
0077
0078 ibuf(1) = my_tile_nx(n)
0079 ibuf(2) = my_tile_ny(n)
0080 ibuf(3) = my_tile_i0(n)
0081 ibuf(4) = my_tile_j0(n)
bf4be02920 Jean*0082
0083
0084 count=4
0085 datatype=MPI_INTEGER
0086 dest=my_coupler_rank
44feb30e11 Jean*0087 tag=generate_tag(113,n,'Register each tile')
bf4be02920 Jean*0088 comm=MPI_COMM_myglobal
0089
0090 call MPI_Send( ibuf, count, datatype, dest, tag, comm, ierr )
951926fb9b Jean*0091
bf4be02920 Jean*0092 if (ierr.ne.0) then
0093 write(LogUnit,*) 'MITCOMPONENT_tile_register: rank(W,G,L)=',
0094 & my_rank_in_world,my_rank_in_global,my_rank_in_local,
0095 & ' ierr=',ierr
44feb30e11 Jean*0096 STOP 'MITCOMPONENT_tile_register: MPI_Send failed'
bf4be02920 Jean*0097 endif
0098
0099 enddo
0100
44feb30e11 Jean*0101 write(LogUnit,'(A,I6,A,2I4)')
0102 & 'MITCOMPONENT_tile_register (pId=', my_rank_in_local, '): done'
0103
bf4be02920 Jean*0104
0105 call flush(LogUnit)
0106 return
0107 end
0108