1 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
2 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
3 | !----------------------------------------------------------------------- |
---|
4 | ! CVS twocmp.seqUnvn.F90,v 1.6 2007-12-19 17:13:17 rloy Exp |
---|
5 | ! CVS MCT_2_8_0 |
---|
6 | !BOP ------------------------------------------------------------------- |
---|
7 | ! |
---|
8 | ! !ROUTINE: twocomponentUneven.sequential |
---|
9 | ! |
---|
10 | ! !DESCRIPTION: Provide a simple example of using MCT to connect two components |
---|
11 | ! In this case the models are running sequentialy but the second model |
---|
12 | ! is only running on 1 processor. |
---|
13 | ! |
---|
14 | ! !INTERFACE: |
---|
15 | ! |
---|
16 | program twosequn |
---|
17 | ! |
---|
18 | ! !USES: |
---|
19 | ! |
---|
20 | !--- Get only the things needed from MCT |
---|
21 | use m_MCTWorld,only: MCTWorld_init => init |
---|
22 | |
---|
23 | use m_GlobalSegMap,only: GlobalSegMap |
---|
24 | use m_GlobalSegMap,only: MCT_GSMap_init => init |
---|
25 | use m_GlobalSegMap,only: MCT_GSMap_lsize => lsize |
---|
26 | |
---|
27 | use m_AttrVect,only : AttrVect |
---|
28 | use m_AttrVect,only : MCT_AtrVt_init => init |
---|
29 | use m_AttrVect,only : MCT_AtrVt_zero => zero |
---|
30 | use m_AttrVect,only : MCT_AtrVt_lsize => lsize |
---|
31 | use m_AttrVect,only : MCT_AtrVt_indexRA => indexRA |
---|
32 | use m_AttrVect,only : MCT_AtrVt_importRA => importRAttr |
---|
33 | |
---|
34 | use m_Rearranger,only: Rearranger |
---|
35 | use m_Rearranger,only: MCT_Rearranger_init => init |
---|
36 | use m_Rearranger,only: MCT_Rearrange => Rearrange |
---|
37 | |
---|
38 | implicit none |
---|
39 | |
---|
40 | include 'mpif.h' |
---|
41 | |
---|
42 | integer,parameter :: ngx = 6 ! points in x-direction |
---|
43 | integer,parameter :: ngy = 4 ! points in y-direction |
---|
44 | |
---|
45 | integer ier,world_group,model2_group,myrank2,myrank3 |
---|
46 | integer,dimension(:),pointer :: myids,mycomms,peloc2 |
---|
47 | integer,dimension(:,:),pointer :: GlobalId |
---|
48 | integer :: comm1,comm2,asize,mysize,i,myproc |
---|
49 | integer :: commsize |
---|
50 | integer,dimension(1) :: start1,length1,ranks |
---|
51 | integer,dimension(:),allocatable :: start2,length2 |
---|
52 | !----------------------------------------------------------------------- |
---|
53 | ! The Main program. |
---|
54 | ! We are implementing a single-executable, sequential-execution system. |
---|
55 | ! Because its sequential, communication occurs through the main using |
---|
56 | ! arguments. The second component is only running on 1 processor |
---|
57 | |
---|
58 | type(GlobalSegMap) :: GSmap1,GSmap2 |
---|
59 | type(AttrVect) :: av1,av2 |
---|
60 | type(Rearranger) :: Rearr |
---|
61 | |
---|
62 | call MPI_init(ier) |
---|
63 | |
---|
64 | call mpi_comm_size(MPI_COMM_WORLD, mysize,ier) |
---|
65 | if(mysize .gt. 12) then |
---|
66 | write(6,*)"Must run on less than 12 processors" |
---|
67 | stop |
---|
68 | endif |
---|
69 | call mpi_comm_rank(MPI_COMM_WORLD, myproc,ier) |
---|
70 | |
---|
71 | ! the first model is running on all the processors so give |
---|
72 | ! it a dubplicate of MPI_COMM_WORLD for its communicator |
---|
73 | call mpi_comm_dup(MPI_COMM_WORLD,comm1,ier) |
---|
74 | |
---|
75 | ! the second model is only running on one processor |
---|
76 | ! so use mpi_groups methods to define its communicator |
---|
77 | call mpi_comm_group(MPI_COMM_WORLD,world_group,ier) |
---|
78 | |
---|
79 | ! need a communicator that only has the first processor |
---|
80 | ranks(1)=0 |
---|
81 | ! define the group |
---|
82 | call mpi_group_incl(world_group,1,ranks,model2_group,ier) |
---|
83 | ! now define the communicator |
---|
84 | ! first initialize it |
---|
85 | comm2=MPI_COMM_NULL |
---|
86 | call mpi_comm_create(MPI_COMM_WORLD,model2_group,comm2,ier) |
---|
87 | |
---|
88 | ! don't need the groups anymore |
---|
89 | call mpi_group_free(world_group,ier) |
---|
90 | call mpi_group_free(model2_group,ier) |
---|
91 | |
---|
92 | ! allocate arrays for the ids and comms |
---|
93 | allocate(myids(2),mycomms(2)) |
---|
94 | |
---|
95 | ! Set the arrays to their values. |
---|
96 | myids(1)=1 |
---|
97 | myids(2)=2 |
---|
98 | mycomms(1)=comm1 |
---|
99 | mycomms(2)=comm2 |
---|
100 | |
---|
101 | ! now call the initm_ version of MCTWorld_init |
---|
102 | call MCTWorld_init(2,MPI_COMM_WORLD,mycomms,myids) |
---|
103 | |
---|
104 | |
---|
105 | ! first gsmap is the grid decomposed in one dimension |
---|
106 | ! there is 1 segment per processor |
---|
107 | length1(1)= (ngx * ngy)/mysize |
---|
108 | start1(1)= myproc * length1(1) + 1 |
---|
109 | |
---|
110 | write(6,*)'gsmap1', myproc,length1(1),start1(1) |
---|
111 | call MCT_GSMap_init(GSMap1,start1,length1,0,comm1,1) |
---|
112 | |
---|
113 | ! second gsmap is the grid on one processor |
---|
114 | |
---|
115 | ! for GSMap init to work, the size of the start and length arrays |
---|
116 | ! must equal the number of local segments. So I must allocate |
---|
117 | ! size zero arrays on the other processors. |
---|
118 | if(myproc .eq. 0) then |
---|
119 | allocate(start2(1),length2(1)) |
---|
120 | length2(1) = ngx*ngy |
---|
121 | start2(1) = 1 |
---|
122 | else |
---|
123 | allocate(start2(0),length2(0)) |
---|
124 | endif |
---|
125 | |
---|
126 | call MCT_GSMap_init(GSMap2,start2,length2,0,comm1,2) |
---|
127 | write(6,*)'gsmap2', myproc,GSMap2%ngseg,GSmap2%gsize,GSmap2%start(1), & |
---|
128 | GSmap2%pe_loc(1),GSmap2%length(1) |
---|
129 | |
---|
130 | |
---|
131 | ! initialize an Av on each GSMap |
---|
132 | call MCT_AtrVt_init(av1,rList="field1:field2",lsize=MCT_GSMap_lsize(GSMap1,comm1)) |
---|
133 | |
---|
134 | ! Use comm1 because lsize of GSMap2 on comm1 will return 0 on non-root processors. |
---|
135 | ! We need av2 to be full-sized on proc 0 and 0 size on other processors. |
---|
136 | call MCT_AtrVt_init(av2,rList="field1:field2",lsize=MCT_GSMap_lsize(GSMap2,comm1)) |
---|
137 | |
---|
138 | |
---|
139 | ! create a rearranger. Use the communicator which contains all processors |
---|
140 | ! involved in the rearrangement, comm1 |
---|
141 | call MCT_Rearranger_init(GSMap1,GSMap2,comm1,Rearr) |
---|
142 | |
---|
143 | !-------------end of initialization steps |
---|
144 | |
---|
145 | |
---|
146 | ! Start up model1 which fills av1 with data. |
---|
147 | call model1(comm1,av1) |
---|
148 | |
---|
149 | ! print out Av data |
---|
150 | do i=1,MCT_AtrVt_lsize(av1) |
---|
151 | write(6,*) "model 1 data", myproc,i,av1%rAttr(1,i),av1%rAttr(2,i) |
---|
152 | enddo |
---|
153 | |
---|
154 | ! rearrange data from model1 so that model2 can use it. |
---|
155 | call MCT_Rearrange(av1,av2,Rearr) |
---|
156 | |
---|
157 | ! pass data to model2 (which will print it out) |
---|
158 | ! model2 should only run on one processor. |
---|
159 | if(myproc .eq. 0) then |
---|
160 | call model2(comm2,av2) |
---|
161 | endif |
---|
162 | |
---|
163 | |
---|
164 | ! all done |
---|
165 | call MPI_Barrier(MPI_COMM_WORLD,ier) |
---|
166 | if (myproc==0) write(6,*) 'All Done' |
---|
167 | |
---|
168 | call mpi_finalize(ier) |
---|
169 | |
---|
170 | contains |
---|
171 | |
---|
172 | !----------------------------------------------------------------------- |
---|
173 | !----------------------------------------------------------------------- |
---|
174 | ! !ROUTINE: |
---|
175 | subroutine model1(comm1,mod1av) ! the first model |
---|
176 | |
---|
177 | implicit none |
---|
178 | |
---|
179 | integer :: comm1,mysize,ier,asize,myproc |
---|
180 | integer :: fieldindx,avsize,i |
---|
181 | integer,dimension(1) :: start,length |
---|
182 | real,pointer :: testarray(:) |
---|
183 | |
---|
184 | type(GlobalSegMap) :: GSmap |
---|
185 | type(AttrVect) :: mod1av |
---|
186 | !--------------------------- |
---|
187 | |
---|
188 | ! find local rank and size |
---|
189 | call mpi_comm_size(comm1,mysize,ier) |
---|
190 | call mpi_comm_rank(comm1,myproc,ier) |
---|
191 | write(6,*)"model1 myproc,mysize",myproc,mysize |
---|
192 | |
---|
193 | |
---|
194 | avsize = MCT_AtrVt_lsize(mod1av) |
---|
195 | write(6,*)"model 1 myproc, av size", myproc,avsize |
---|
196 | |
---|
197 | ! Fill Av with some data |
---|
198 | ! fill first attribute the direct way |
---|
199 | fieldindx = MCT_AtrVt_indexRA(mod1av,"field1") |
---|
200 | do i=1,avsize |
---|
201 | mod1av%rAttr(fieldindx,i) = float(i+ 20*myproc) |
---|
202 | enddo |
---|
203 | |
---|
204 | ! fill second attribute using Av import function |
---|
205 | allocate(testarray(avsize)) |
---|
206 | do i=1,avsize |
---|
207 | testarray(i)= cos((float(i+ 20*myproc)/24.) * 3.14) |
---|
208 | enddo |
---|
209 | call MCT_AtrVt_importRA(mod1av,"field2",testarray) |
---|
210 | |
---|
211 | |
---|
212 | end subroutine model1 |
---|
213 | |
---|
214 | !----------------------------------------------------------------------- |
---|
215 | !----------------------------------------------------------------------- |
---|
216 | ! !ROUTINE: |
---|
217 | subroutine model2(comm2,mod2av) |
---|
218 | |
---|
219 | implicit none |
---|
220 | |
---|
221 | integer :: comm2,mysize,ier,asize,myproc |
---|
222 | integer :: i |
---|
223 | type(AttrVect) :: mod2av |
---|
224 | !--------------------------- |
---|
225 | |
---|
226 | ! find local rank and size |
---|
227 | call mpi_comm_size(comm2,mysize,ier) |
---|
228 | call mpi_comm_rank(comm2,myproc,ier) |
---|
229 | write(6,*)"model2 myproc,mysize",myproc,mysize |
---|
230 | |
---|
231 | asize = MCT_AtrVt_lsize(mod2av) |
---|
232 | write(6,*)"model 2 myproc, av size", myproc,asize |
---|
233 | |
---|
234 | ! print out Av data |
---|
235 | do i=1,asize |
---|
236 | write(6,*) "model 2 data after", myproc,i,mod2av%rAttr(1,i),mod2av%rAttr(2,i) |
---|
237 | enddo |
---|
238 | |
---|
239 | |
---|
240 | end subroutine model2 |
---|
241 | |
---|
242 | end |
---|