1 |
!======================================================================= |
2 |
subroutine comprecv_r8tiles( dataname, Ni,Oi,Nj,Oj,Nk,Tx,Ty, arr ) |
3 |
implicit none |
4 |
! Arguments |
5 |
character*(*) dataname |
6 |
integer Ni,Oi,Nj,Oj,Io,Jo,Nk,Tx,Ty |
7 |
real*8 arr(1-Oi:Ni+Oi,1-Oj:Nj+Oj,Nk,Tx,Ty) |
8 |
! Predefined constants/arrays |
9 |
#include "CPLR_SIG.h" |
10 |
! MPI variables |
11 |
#include "mpif.h" |
12 |
integer count,dtype,rank,tag,comm,ierr |
13 |
integer stat(MPI_STATUS_SIZE) |
14 |
! Functions |
15 |
integer generate_tag |
16 |
! Local |
17 |
integer i,j,ij,nx,ny,k,bibj,bi,bj |
18 |
character*(MAXLEN_COMP_NAME) recvdname |
19 |
! ------------------------------------------------------------------ |
20 |
|
21 |
if (HEADER_SIZE+Ni*Nj.gt.MAX_R8_BUFLEN) |
22 |
& stop 'comprecv_r8tiles: Nx*Ny too big' |
23 |
|
24 |
! Foreach tile which is non-blank |
25 |
do bibj=1,my_num_tiles |
26 |
|
27 |
bi=my_tile_bi(bibj) |
28 |
bj=my_tile_bj(bibj) |
29 |
|
30 |
! Receive message |
31 |
count=MAX_R8_BUFLEN |
32 |
dtype=MPI_DOUBLE_PRECISION |
33 |
tag=generate_tag(123,bibj,dataname) |
34 |
rank=my_coupler_rank |
35 |
comm=MPI_COMM_myglobal |
36 |
|
37 |
if (VERB) then |
38 |
write(LogUnit,*) 'comprecv_r8tiles: calling MPI_Recv rank=',rank |
39 |
write(LogUnit,*) 'comprecv_r8tiles: dataname=',dataname |
40 |
call flush(LogUnit) |
41 |
endif |
42 |
call MPI_Recv(r8buf, count, dtype, rank, tag, comm, stat, ierr) |
43 |
if (VERB) then |
44 |
write(LogUnit,*) 'comprecv_r8tiles: returned ierr=',ierr |
45 |
call flush(LogUnit) |
46 |
endif |
47 |
|
48 |
if (ierr.ne.0) then |
49 |
write(LogUnit,*) 'comprecv_r8tiles: rank(W,G)=', |
50 |
& my_rank_in_world,my_rank_in_global, |
51 |
& ' ierr=',ierr |
52 |
stop 'comprecv_r8tiles: MPI_Recv failed' |
53 |
endif |
54 |
|
55 |
! Extract buffer |
56 |
Io=int(0.5+r8buf(1)) |
57 |
Jo=int(0.5+r8buf(2)) |
58 |
nx=int(0.5+r8buf(3)) |
59 |
ny=int(0.5+r8buf(4)) |
60 |
call mitcplr_real2char( r8buf(9), recvdname ) |
61 |
|
62 |
if (Io.ne.my_tile_i0(bibj)) stop 'comprecv_r8tiles: bad Io' |
63 |
if (Jo.ne.my_tile_j0(bibj)) stop 'comprecv_r8tiles: bad Jo' |
64 |
if (nx.ne.my_tile_nx(bibj)) stop 'comprecv_r8tiles: bad nx' |
65 |
if (ny.ne.my_tile_ny(bibj)) stop 'comprecv_r8tiles: bad ny' |
66 |
if (recvdname .ne. dataname) then |
67 |
write(LogUnit,*) 'comprecv_r8tiles: recvdname = ',recvdname |
68 |
write(LogUnit,*) 'comprecv_r8tiles: dataname = ',dataname |
69 |
stop 'comprecv_r8tiles: recvdname != dataname' |
70 |
endif |
71 |
|
72 |
! Copy buffer to interior of tile |
73 |
k=1 |
74 |
do j=1,Nj |
75 |
do i=1,Ni |
76 |
ij=HEADER_SIZE+i+Ni*(j-1) |
77 |
arr(i,j,k,bi,bj)=r8buf(ij) |
78 |
enddo |
79 |
enddo |
80 |
|
81 |
enddo ! bibj |
82 |
|
83 |
! ------------------------------------------------------------------ |
84 |
return |
85 |
end |
86 |
!======================================================================= |