1 |
jmc |
1.1 |
!======================================================================= |
2 |
|
|
subroutine MITCOMPONENT_tile_register( ntx, nty, ireg, rreg ) |
3 |
|
|
implicit none |
4 |
|
|
|
5 |
|
|
! Arguments |
6 |
|
|
integer ntx,nty |
7 |
|
|
integer ireg(4,ntx,nty) |
8 |
|
|
real*4 rreg(4,ntx,nty) |
9 |
|
|
|
10 |
|
|
! MPI variables |
11 |
|
|
#include "mpif.h" |
12 |
|
|
integer myid, numprocs, ierr, rc |
13 |
|
|
|
14 |
|
|
! Predefined constants/arrays |
15 |
|
|
#include "CPLR_SIG.h" |
16 |
|
|
|
17 |
|
|
! Functions |
18 |
|
|
integer mitcplr_match_comp |
19 |
|
|
integer generate_tag |
20 |
|
|
|
21 |
|
|
! Local |
22 |
|
|
integer bi,bj,num_tiles,bibj |
23 |
|
|
integer count,datatype,dest,tag,comm |
24 |
|
|
|
25 |
|
|
! ------------------------------------------------------------------ |
26 |
|
|
|
27 |
|
|
num_tiles=0 |
28 |
|
|
do bj=1,nty |
29 |
|
|
do bi=1,ntx |
30 |
|
|
if ( ireg(1,bi,bj)*ireg(2,bi,bj).gt.0 ) then |
31 |
|
|
num_tiles=num_tiles+1 |
32 |
|
|
my_tile_nx(num_tiles)=ireg(1,bi,bj) |
33 |
|
|
my_tile_ny(num_tiles)=ireg(2,bi,bj) |
34 |
|
|
my_tile_i0(num_tiles)=ireg(3,bi,bj) |
35 |
|
|
my_tile_j0(num_tiles)=ireg(4,bi,bj) |
36 |
|
|
my_tile_bi(num_tiles)=bi |
37 |
|
|
my_tile_bj(num_tiles)=bj |
38 |
|
|
if (VERB) |
39 |
|
|
& write(LogUnit,*) 'MITCOMPONENT_tile_register: bi,bj=',bi,bj, |
40 |
|
|
& ' nx,ny=',ireg(1,bi,bj),ireg(2,bi,bj) |
41 |
|
|
else |
42 |
|
|
if (VERB) write(LogUnit,*) |
43 |
|
|
& 'MITCOMPONENT_tile_register: blank bi,bj=',bi,bj |
44 |
|
|
endif |
45 |
|
|
enddo |
46 |
|
|
enddo |
47 |
|
|
my_num_tiles=num_tiles |
48 |
|
|
if (VERB) write(LogUnit,*) |
49 |
|
|
& 'MITCOMPONENT_tile_register: num_tiles =',num_tiles |
50 |
|
|
|
51 |
|
|
if (num_tiles.lt.1) |
52 |
|
|
& stop 'MITCOMPONENT_tile_register: num_tiles < 1' |
53 |
|
|
if (num_tiles.gt.MAX_TILES) |
54 |
|
|
& stop 'MITCOMPONENT_tile_register: num_tiles > MAX_TILES' |
55 |
|
|
|
56 |
|
|
! Set up buffer |
57 |
|
|
ibuf(1)=num_tiles |
58 |
|
|
|
59 |
|
|
! Send message |
60 |
|
|
count=1 |
61 |
|
|
datatype=MPI_INTEGER |
62 |
|
|
dest=my_coupler_rank |
63 |
|
|
tag=generate_tag(112,my_rank_in_global,'Register Tiles') |
64 |
|
|
comm=MPI_COMM_myglobal |
65 |
|
|
|
66 |
|
|
call MPI_Send( ibuf, count, datatype, dest, tag, comm, ierr ) |
67 |
|
|
|
68 |
|
|
if (ierr.ne.0) then |
69 |
|
|
write(LogUnit,*) 'MITCOMPONENT_tile_register: rank(W,G,L)=', |
70 |
|
|
& my_rank_in_world,my_rank_in_global,my_rank_in_local, |
71 |
|
|
& ' ierr=',ierr |
72 |
|
|
stop 'MITCOMPONENT_tile_register: MPI_Send failed' |
73 |
|
|
endif |
74 |
|
|
|
75 |
|
|
do bibj=1,my_num_tiles |
76 |
|
|
|
77 |
|
|
! Set up buffer |
78 |
|
|
bi=my_tile_bi(bibj) |
79 |
|
|
bj=my_tile_bj(bibj) |
80 |
|
|
ibuf(1)=my_tile_nx(bibj) |
81 |
|
|
ibuf(2)=my_tile_ny(bibj) |
82 |
|
|
ibuf(3)=my_tile_i0(bibj) |
83 |
|
|
ibuf(4)=my_tile_j0(bibj) |
84 |
|
|
|
85 |
|
|
! Send message |
86 |
|
|
count=4 |
87 |
|
|
datatype=MPI_INTEGER |
88 |
|
|
dest=my_coupler_rank |
89 |
|
|
tag=generate_tag(113,bibj,'Register each tile') |
90 |
|
|
comm=MPI_COMM_myglobal |
91 |
|
|
|
92 |
|
|
call MPI_Send( ibuf, count, datatype, dest, tag, comm, ierr ) |
93 |
|
|
|
94 |
|
|
if (ierr.ne.0) then |
95 |
|
|
write(LogUnit,*) 'MITCOMPONENT_tile_register: rank(W,G,L)=', |
96 |
|
|
& my_rank_in_world,my_rank_in_global,my_rank_in_local, |
97 |
|
|
& ' ierr=',ierr |
98 |
|
|
stop 'MITCOMPONENT_tile_register: MPI_Send failed' |
99 |
|
|
endif |
100 |
|
|
|
101 |
|
|
enddo |
102 |
|
|
|
103 |
|
|
! ------------------------------------------------------------------ |
104 |
|
|
call flush(LogUnit) |
105 |
|
|
return |
106 |
|
|
end |
107 |
|
|
!======================================================================= |