1 |
jmc |
1.1 |
!======================================================================= |
2 |
|
|
subroutine coupsend_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,dest,tag,comm,ierr |
14 |
|
|
! Functions |
15 |
|
|
integer mitcplr_match_comp |
16 |
|
|
integer generate_tag |
17 |
|
|
! Local |
18 |
|
|
integer compind,numprocs |
19 |
|
|
integer i,j,ij,n,bibj |
20 |
|
|
integer Ni,Io,Nj,Jo |
21 |
|
|
! ------------------------------------------------------------------ |
22 |
|
|
|
23 |
|
|
! Establish who I am communicating with |
24 |
|
|
compind=mitcplr_match_comp( component ) |
25 |
|
|
if (compind.le.0) stop 'coupsend_r8tiles: Bad component id' |
26 |
|
|
comm=MPI_COMM_compcplr( compind ) |
27 |
|
|
numprocs=num_component_procs(compind) |
28 |
|
|
if (numprocs.lt.1) then |
29 |
|
|
write(LogUnit,*) 'coupsend_r8tiles: compind = ',compind |
30 |
|
|
stop 'coupsend_r8tiles: numprocs < 1' |
31 |
|
|
endif |
32 |
|
|
if (VERB) |
33 |
|
|
& write(LogUnit,*) 'coupsend_r8tiles: ',component_Name(compind) |
34 |
|
|
if (VERB) |
35 |
|
|
& write(LogUnit,*) 'coupsend_r8tiles: dataname=',dataname |
36 |
|
|
|
37 |
|
|
! Foreach component process |
38 |
|
|
do n=1,numprocs |
39 |
|
|
|
40 |
|
|
! Foreach tile on that process |
41 |
|
|
do bibj=1,component_num_tiles(n,compind) |
42 |
|
|
|
43 |
|
|
! Create header |
44 |
|
|
Io=component_tile_i0(bibj,n,compind) |
45 |
|
|
Jo=component_tile_j0(bibj,n,compind) |
46 |
|
|
Ni=component_tile_nx(bibj,n,compind) |
47 |
|
|
Nj=component_tile_ny(bibj,n,compind) |
48 |
|
|
r8buf(1)=float( Io ) |
49 |
|
|
r8buf(2)=float( Jo ) |
50 |
|
|
r8buf(3)=float( Ni ) |
51 |
|
|
r8buf(4)=float( Nj ) |
52 |
|
|
call mitcplr_char2real( dataname, r8buf(9) ) |
53 |
|
|
|
54 |
|
|
! Pack data |
55 |
|
|
do j=1,Nj |
56 |
|
|
do i=1,Ni |
57 |
|
|
ij=HEADER_SIZE+i+Ni*(j-1) |
58 |
|
|
r8buf(ij)=arr(Io+i-1,Jo+j-1) |
59 |
|
|
enddo |
60 |
|
|
enddo |
61 |
|
|
|
62 |
|
|
! Send message |
63 |
|
|
count=HEADER_SIZE+Ni*Nj |
64 |
|
|
dtype=MPI_DOUBLE_PRECISION |
65 |
|
|
tag=generate_tag(123,bibj,dataname) |
66 |
|
|
dest=rank_component_procs(n,compind) |
67 |
|
|
|
68 |
|
|
if (VERB) then |
69 |
|
|
write(LogUnit,*) 'coupsend_r8tiles: calling MPI_Send dest=', |
70 |
|
|
& dest,' proc=',n,'/',numprocs,' tile=',bibj |
71 |
|
|
call flush(LogUnit) |
72 |
|
|
endif |
73 |
|
|
call MPI_Send( r8buf, count, dtype, dest, tag, comm, ierr ) |
74 |
|
|
if (VERB) then |
75 |
|
|
write(LogUnit,*) 'coupsend_r8tiles: returned ierr=',ierr |
76 |
|
|
call flush(LogUnit) |
77 |
|
|
endif |
78 |
|
|
|
79 |
|
|
if (ierr.ne.0) then |
80 |
|
|
write(LogUnit,*) 'coupsend_r8tiles: rank(W,G)=', |
81 |
|
|
& my_rank_in_world,my_rank_in_global, |
82 |
|
|
& ' ierr=',ierr |
83 |
|
|
stop 'coupsend_r8tiles: MPI_Recv failed' |
84 |
|
|
endif |
85 |
|
|
|
86 |
|
|
enddo ! bibj |
87 |
|
|
|
88 |
|
|
enddo ! n |
89 |
|
|
|
90 |
|
|
! ------------------------------------------------------------------ |
91 |
|
|
return |
92 |
|
|
end |
93 |
|
|
!======================================================================= |