1 |
!======================================================================= |
2 |
subroutine couprecv_r8tiles( component, dataname, Nx, Ny, arr ) |
3 |
implicit none |
4 |
! Arguments |
5 |
character*(*) component |
6 |
character*(*) dataname |
7 |
integer Nx,Ny |
8 |
real*8 arr(Nx,Ny) |
9 |
! Predefined constants/arrays |
10 |
#include "CPLR_SIG.h" |
11 |
! MPI variables |
12 |
#include "mpif.h" |
13 |
integer count,dtype,rank,tag,comm,ierr |
14 |
integer stat(MPI_STATUS_SIZE) |
15 |
! Functions |
16 |
integer mitcplr_match_comp |
17 |
integer generate_tag |
18 |
! Local |
19 |
integer compind,numprocs |
20 |
integer i,j,ij,n,bibj |
21 |
integer Ni,Io,Nj,Jo |
22 |
character*(MAXLEN_COMP_NAME) recvdname |
23 |
! ------------------------------------------------------------------ |
24 |
|
25 |
! Establish who I am communicating with |
26 |
compind=mitcplr_match_comp( component ) |
27 |
if (compind.le.0) stop 'couprecv_r8tiles: Bad component id' |
28 |
comm=MPI_COMM_compcplr( compind ) |
29 |
numprocs=num_component_procs(compind) |
30 |
if (numprocs.lt.1) then |
31 |
write(LogUnit,*) 'couprecv_r8tiles: compind = ',compind |
32 |
stop 'couprecv_r8tiles: numprocs < 1' |
33 |
endif |
34 |
if (VERB) |
35 |
& write(LogUnit,*) 'couprecv_r8tiles: ',component_Name(compind) |
36 |
if (VERB) |
37 |
& write(LogUnit,*) 'couprecv_r8tiles: dataname=',dataname |
38 |
|
39 |
! Foreach component process |
40 |
do n=1,numprocs |
41 |
|
42 |
! Foreach tile within process |
43 |
do bibj=1,component_num_tiles(n,compind) |
44 |
|
45 |
! Receive message |
46 |
count=MAX_R8_BUFLEN |
47 |
dtype=MPI_DOUBLE_PRECISION |
48 |
tag=generate_tag(103,bibj,dataname) |
49 |
rank=rank_component_procs(n,compind) |
50 |
|
51 |
if (VERB) then |
52 |
write(LogUnit,*) |
53 |
& 'couprecv_r8tiles: calling MPI_Recv rank=',rank, |
54 |
& ' proc=',n,'/',numprocs,' tile=',bibj |
55 |
call flush(LogUnit) |
56 |
endif |
57 |
call MPI_Recv(r8buf, count, dtype, rank, tag, comm, stat, ierr) |
58 |
if (VERB) then |
59 |
write(LogUnit,*) 'couprecv_r8tiles: returned ierr=',ierr |
60 |
call flush(LogUnit) |
61 |
endif |
62 |
|
63 |
if (ierr.ne.0) then |
64 |
write(LogUnit,*) 'couprecv_r8tiles: rank(W,G)=', |
65 |
& my_rank_in_world,my_rank_in_global, |
66 |
& ' ierr=',ierr |
67 |
stop 'couprecv_r8tiles: MPI_Recv failed' |
68 |
endif |
69 |
|
70 |
! Extract header |
71 |
Io=int(0.5+r8buf(1)) |
72 |
Ni=int(0.5+r8buf(2)) |
73 |
Jo=int(0.5+r8buf(3)) |
74 |
Nj=int(0.5+r8buf(4)) |
75 |
|
76 |
if (Io+Ni-1.gt.Nx .or. Io.lt.1) then |
77 |
write(LogUnit,*) 'couprecv_r8tiles: Io,Ni=',Io,Ni |
78 |
stop 'couprecv_r8tiles: Incompatible header/target array' |
79 |
endif |
80 |
if (Jo+Nj-1.gt.Ny .or. Jo.lt.1) then |
81 |
write(LogUnit,*) 'couprecv_r8tiles: Jo,Nj=',Jo,Nj |
82 |
stop 'couprecv_r8tiles: Incompatible header/target array' |
83 |
endif |
84 |
if (Io.ne.component_tile_i0(bibj,n,compind)) |
85 |
& stop 'couprecv_r8tiles: Io != component_tile_i0' |
86 |
if (Jo.ne.component_tile_j0(bibj,n,compind)) |
87 |
& stop 'couprecv_r8tiles: Jo != component_tile_j0' |
88 |
if (Ni.ne.component_tile_nx(bibj,n,compind)) |
89 |
& stop 'couprecv_r8tiles: Ni != component_tile_nx' |
90 |
if (Nj.ne.component_tile_ny(bibj,n,compind)) |
91 |
& stop 'couprecv_r8tiles: Nj != component_tile_ny' |
92 |
|
93 |
call mitcplr_real2char( r8buf(9), recvdname ) |
94 |
if (recvdname .ne. dataname) then |
95 |
write(LogUnit,*) 'couprecv_r8tiles: recvdname = ',recvdname |
96 |
write(LogUnit,*) 'couprecv_r8tiles: dataname = ',dataname |
97 |
stop 'couprecv_r8tiles: recvdname != dataname' |
98 |
endif |
99 |
|
100 |
! Extract data |
101 |
do j=1,Nj |
102 |
do i=1,Ni |
103 |
ij=HEADER_SIZE+i+Ni*(j-1) |
104 |
arr(Io+i-1,Jo+j-1)=r8buf(ij) |
105 |
enddo |
106 |
enddo |
107 |
|
108 |
enddo ! bibj |
109 |
|
110 |
enddo ! n |
111 |
|
112 |
! ------------------------------------------------------------------ |
113 |
return |
114 |
end |
115 |
!======================================================================= |