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