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 couprecv_r8( component, dataname, Nx, Ny, arr )
0003 implicit none
0004
0005 character*(*) component
0006 character*(*) dataname
0007 integer Nx,Ny
0008 real*8 arr(Nx,Ny)
0009
0010 #include "CPLR_SIG.h"
0011
0012 #include "mpif.h"
0013 integer count,dtype,rank,tag,comm,ierr
0014 integer stat(MPI_STATUS_SIZE)
0015
0016 integer mitcplr_match_comp
0017 integer generate_tag
0018
0019 integer compind,numprocs
0020 integer i,j,ij,n
0021 integer Ni,Io,Nj,Jo
0022 character*(MAXLEN_COMP_NAME) recvdname
0023
0024
0025
0026 compind=mitcplr_match_comp( component )
0027 if (compind.le.0) stop 'couprecv_r8: Bad component id'
0028 comm=MPI_COMM_compcplr( compind )
0029 numprocs=num_component_procs(compind)
0030 if (numprocs.lt.1) then
0031 write(LogUnit,*) 'couprecv_r8: compind = ',compind
0032 stop 'couprecv_r8: numprocs < 1'
0033 endif
0034 if (VERB)
0035 & write(LogUnit,*) 'couprecv_r8: ',component_Name(compind)
0036 if (VERB)
0037 & write(LogUnit,*) 'couprecv_r8: dataname=',dataname
0038
0039
0040 do n=1,numprocs
0041
0042
0043 count=HEADER_SIZE+MAX_R8_BUFLEN
0044 dtype=MPI_DOUBLE_PRECISION
0045 tag=generate_tag(102,n,dataname)
0046 rank=rank_component_procs(n,compind)
0047 if (VERB) then
0048 write(LogUnit,*) 'couprecv_r8: calling MPI_Recv rank=',rank,
0049 & ' proc=',n,'/',numprocs
0050 call flush(LogUnit)
0051 endif
0052 call MPI_Recv(r8buf, count, dtype, rank, tag, comm, stat, ierr)
0053 if (VERB) then
0054 write(LogUnit,*) 'couprecv_r8: returned ierr=',ierr
0055 call flush(LogUnit)
0056 endif
0057
0058 if (ierr.ne.0) then
0059 write(LogUnit,*) 'couprecv_r8tiles: rank(W,G)=',
0060 & my_rank_in_world,my_rank_in_global,
0061 & ' ierr=',ierr
0062 stop 'couprecv_r8: MPI_Recv failed'
0063 endif
0064
0065
0066 Io=int(0.5+r8buf(1))
0067 Ni=int(0.5+r8buf(2))
0068 Jo=int(0.5+r8buf(3))
0069 Nj=int(0.5+r8buf(4))
0070
0071 if (Io+Ni-1.gt.Nx .or. Io.lt.1) then
0072 write(LogUnit,*) 'couprecv_r8tiles: Io,Ni=',Io,Ni
0073 stop 'couprecv_r8: Incompatible header/target array'
0074 endif
0075 if (Jo+Nj-1.gt.Ny .or. Jo.lt.1) then
0076 write(LogUnit,*) 'couprecv_r8tiles: Jo,Nj=',Jo,Nj
0077 stop 'couprecv_r8: Incompatible header/target array'
0078 endif
0079
e5266ce10c Jean*0080 call mitcplr_dbl2char( r8buf(9), recvdname )
bf4be02920 Jean*0081 if (recvdname .ne. dataname) then
0082 write(LogUnit,*) 'couprecv_r8: recvdname = ',recvdname
0083 write(LogUnit,*) 'couprecv_r8: dataname = ',dataname
0084 stop 'couprecv_r8: recvdname != dataname'
0085 endif
0086
0087
0088 do j=1,Nj
0089 do i=1,Ni
0090 ij=HEADER_SIZE+i+Ni*(j-1)
0091 arr(Io+i-1,Jo+j-1)=r8buf(ij)
0092 enddo
0093 enddo
0094
0095 enddo
0096
0097
0098 return
0099 end
0100