1 |
C $Header: $ |
2 |
C $Name: $ |
3 |
|
4 |
!======================================================================= |
5 |
subroutine MITCOUPLER_tile_register( compName, nnx, nny ) |
6 |
implicit none |
7 |
|
8 |
! Arguments |
9 |
character*(*) compName |
10 |
integer nnx, nny |
11 |
|
12 |
! MPI variables |
13 |
#include "mpif.h" |
14 |
|
15 |
! Predefined constants/arrays |
16 |
#include "CPLR_SIG.h" |
17 |
|
18 |
! Functions |
19 |
integer mitcplr_match_comp |
20 |
integer generate_tag |
21 |
|
22 |
! Local |
23 |
integer n,numprocs |
24 |
integer comm |
25 |
integer compind,count,dtype,tag,rank |
26 |
integer ierr, rc |
27 |
integer stat(MPI_STATUS_SIZE) |
28 |
integer j,numtiles |
29 |
|
30 |
! ------------------------------------------------------------------ |
31 |
|
32 |
! Establish who I am communicating with |
33 |
compind=mitcplr_match_comp( compName ) |
34 |
if (compind.le.0) stop 'MITCOUPLER_tile_register: Bad component' |
35 |
comm=MPI_COMM_compcplr( compind ) |
36 |
numprocs=num_component_procs(compind) |
37 |
if (numprocs.lt.1) then |
38 |
write(LogUnit,*) 'MITCOUPLER_tile_register: compind = ',compind |
39 |
stop 'MITCOUPLER_tile_register: numprocs < 1' |
40 |
endif |
41 |
|
42 |
! Foreach component process |
43 |
do n=1,numprocs |
44 |
|
45 |
! Receive message |
46 |
count=MAX_IBUF |
47 |
dtype=MPI_INTEGER |
48 |
tag=generate_tag(112,n,'Register Tiles') |
49 |
rank=rank_component_procs(n,compind) |
50 |
|
51 |
call MPI_Recv(ibuf, count, dtype, rank, tag, comm, stat, ierr) |
52 |
|
53 |
if (ierr.ne.0) then |
54 |
write(LogUnit,*) 'MITCOUPLER_tile_register: rank(W,G)=', |
55 |
& my_rank_in_world,my_rank_in_global, |
56 |
& ' ierr=',ierr |
57 |
stop 'MITCOUPLER_tile_register: MPI_Recv failed' |
58 |
endif |
59 |
|
60 |
numtiles=ibuf(1) |
61 |
if (numtiles.lt.1 .or. numtiles.gt.MAX_TILES) then |
62 |
write(LogUnit,*) 'MITCOUPLER_tile_register: #tiles = ',numtiles |
63 |
stop 'MITCOUPLER_tile_register: invalid value for numtiles' |
64 |
endif |
65 |
component_num_tiles(n,compind)=numtiles |
66 |
|
67 |
do j=1,numtiles |
68 |
|
69 |
! Receive message |
70 |
count=MAX_IBUF |
71 |
dtype=MPI_INTEGER |
72 |
tag=generate_tag(113,j,'Register each tile') |
73 |
rank=rank_component_procs(n,compind) |
74 |
|
75 |
call MPI_Recv(ibuf, count, dtype, rank, tag, comm, stat, ierr) |
76 |
|
77 |
if (ierr.ne.0) then |
78 |
write(LogUnit,*) 'MITCOUPLER_tile_register: rank(W,G)=', |
79 |
& my_rank_in_world,my_rank_in_global, |
80 |
& ' ierr=',ierr |
81 |
stop 'MITCOUPLER_tile_register: MPI_Recv failed' |
82 |
endif |
83 |
|
84 |
component_tile_nx(j,n,compind)=ibuf(1) |
85 |
component_tile_ny(j,n,compind)=ibuf(2) |
86 |
component_tile_i0(j,n,compind)=ibuf(3) |
87 |
component_tile_j0(j,n,compind)=ibuf(4) |
88 |
|
89 |
if (VERB) then |
90 |
write(LogUnit,*) 'MITCOUPLER_tile_register:', |
91 |
& my_rank_in_world,my_rank_in_global, |
92 |
& ' proc,tile,nx,ny,i0,j0 = ',n,j,ibuf(1),ibuf(2),ibuf(3),ibuf(4) |
93 |
endif |
94 |
|
95 |
enddo ! j |
96 |
write(LogUnit,*) 'MITCOUPLER_tile_register:', |
97 |
& my_rank_in_world,my_rank_in_global, |
98 |
& ' rank = ',rank, |
99 |
& ' num_tiles = ',numtiles |
100 |
|
101 |
enddo ! n |
102 |
|
103 |
! ------------------------------------------------------------------ |
104 |
call flush(LogUnit) |
105 |
return |
106 |
end |
107 |
!======================================================================= |