1 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
2 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
3 | !----------------------------------------------------------------------- |
---|
4 | ! CVS m_SparseMatrixDecomp.F90,v 1.15 2008-05-12 01:59:32 jacob Exp |
---|
5 | ! CVS MCT_2_8_0 |
---|
6 | !BOP ------------------------------------------------------------------- |
---|
7 | ! |
---|
8 | ! !MODULE: m_SparseMatrixDecomp -- Parallel sparse matrix decomposition. |
---|
9 | ! |
---|
10 | ! !DESCRIPTION: |
---|
11 | ! The {\tt SparseMatrix} datatype provides sparse matrix storage for |
---|
12 | ! the parallel matrix-vector multiplication ${\bf y} = {\bf M} {\bf x}$. |
---|
13 | ! This module provides services to create decompositions for the |
---|
14 | ! {\tt SparseMatrix}. The matrix decompositions available are row |
---|
15 | ! and column decompositions. They are generated by invoking the |
---|
16 | ! appropriate routine in this module, and passing the corresponding |
---|
17 | ! {\em vector} decomposition. For a row (column) decomposition, one |
---|
18 | ! invokes the routine {\tt ByRow()} ({\tt ByColumn()}), passing the |
---|
19 | ! domain decomposition for the vector {\bf y} ({\bf x}). |
---|
20 | ! |
---|
21 | ! !INTERFACE: |
---|
22 | |
---|
23 | module m_SparseMatrixDecomp |
---|
24 | |
---|
25 | private ! except |
---|
26 | |
---|
27 | ! !PUBLIC MEMBER FUNCTIONS: |
---|
28 | ! |
---|
29 | public :: ByColumn |
---|
30 | public :: ByRow |
---|
31 | |
---|
32 | |
---|
33 | interface ByColumn ; module procedure & |
---|
34 | ByColumnGSMap_ |
---|
35 | end interface |
---|
36 | |
---|
37 | interface ByRow ; module procedure & |
---|
38 | ByRowGSMap_ |
---|
39 | end interface |
---|
40 | |
---|
41 | ! !REVISION HISTORY: |
---|
42 | ! 13Apr01 - J.W. Larson <larson@mcs.anl.gov> - initial prototype |
---|
43 | ! and API specifications. |
---|
44 | ! 03Aug01 - E. Ong <eong@mcs.anl.gov> - in ByRowGSMap and ByColumnGSMap, |
---|
45 | ! call GlobalSegMap_init on non-root processes with actual |
---|
46 | ! shaped arguments to satisfy Fortran 90 standard. See |
---|
47 | ! comments in ByRowGSMap/ByColumnGSMap. |
---|
48 | !EOP ___________________________________________________________________ |
---|
49 | |
---|
50 | character(len=*),parameter :: myname='MCT::m_SparseMatrixDecomp' |
---|
51 | |
---|
52 | contains |
---|
53 | |
---|
54 | !------------------------------------------------------------------------- |
---|
55 | ! Math + Computer Science Division / Argonne National Laboratory ! |
---|
56 | !------------------------------------------------------------------------- |
---|
57 | !BOP |
---|
58 | ! |
---|
59 | ! !IROUTINE: ByColumnGSMap_ - Generate Row-based GlobalSegMap for SparseMatrix |
---|
60 | ! |
---|
61 | ! !INTERFACE: |
---|
62 | |
---|
63 | subroutine ByColumnGSMap_(xGSMap, sMat, sMGSMap, root, comm) |
---|
64 | ! |
---|
65 | ! !USES: |
---|
66 | ! |
---|
67 | use m_die, only: MP_perr_die,die |
---|
68 | |
---|
69 | use m_List, only: List |
---|
70 | use m_List, only: List_init => init |
---|
71 | use m_List, only: List_clean => clean |
---|
72 | |
---|
73 | use m_AttrVect, only: AttrVect |
---|
74 | use m_AttrVect, only: AttrVect_init => init |
---|
75 | use m_AttrVect, only: AttrVect_zero => zero |
---|
76 | use m_AttrVect, only: AttrVect_lsize => lsize |
---|
77 | use m_AttrVect, only: AttrVect_indexIA => indexIA |
---|
78 | use m_AttrVect, only: AttrVect_copy => copy |
---|
79 | use m_AttrVect, only: AttrVect_clean => clean |
---|
80 | |
---|
81 | use m_AttrVectComms, only: AttrVect_scatter => scatter |
---|
82 | use m_AttrVectComms, only: AttrVect_gather => gather |
---|
83 | |
---|
84 | use m_GlobalMap, only : GlobalMap |
---|
85 | use m_GlobalMap, only : GlobalMap_init => init |
---|
86 | use m_GlobalMap, only : GlobalMap_clean => clean |
---|
87 | |
---|
88 | use m_GlobalSegMap, only: GlobalSegMap |
---|
89 | use m_GlobalSegMap, only: GlobalSegMap_init => init |
---|
90 | use m_GlobalSegMap, only: GlobalSegMap_peLocs => peLocs |
---|
91 | use m_GlobalSegMap, only: GlobalSegMap_comp_id => comp_id |
---|
92 | |
---|
93 | use m_SparseMatrix, only: SparseMatrix |
---|
94 | use m_SparseMatrix, only: SparseMatrix_lsize => lsize |
---|
95 | use m_SparseMatrix, only: SparseMatrix_SortPermute => SortPermute |
---|
96 | |
---|
97 | implicit none |
---|
98 | |
---|
99 | ! !INPUT PARAMETERS: |
---|
100 | ! |
---|
101 | type(GlobalSegMap), intent(in) :: xGSMap |
---|
102 | integer, intent(in) :: root |
---|
103 | integer, intent(in) :: comm |
---|
104 | |
---|
105 | ! !INPUT/OUTPUT PARAMETERS: |
---|
106 | ! |
---|
107 | type(SparseMatrix), intent(inout) :: sMat |
---|
108 | |
---|
109 | ! !OUTPUT PARAMETERS: |
---|
110 | ! |
---|
111 | type(GlobalSegMap), intent(out) :: sMGSMap |
---|
112 | |
---|
113 | ! !DESCRIPTION: This routine is invoked from all processes on the |
---|
114 | ! communicator {\tt comm} to create from an input {\tt SparseMatrix} |
---|
115 | ! {\tt sMat} (valid only on the {\tt root} process) and an input |
---|
116 | ! {\bf x}-vector decomposition described by the {\tt GlobalSegMap} |
---|
117 | ! argument {\tt xGSMap} (valid at least on the {\tt root}) to create |
---|
118 | ! an output {\tt GlobalSegMap} decomposition of the matrix elements |
---|
119 | ! {\tt sMGSMap}, which is valid on all processes on the communicator. |
---|
120 | ! This matrix {\tt GlobalSegMap} describes the corresponding column |
---|
121 | ! decomposition of {\tt sMat}. |
---|
122 | ! |
---|
123 | ! {\bf N.B.}: The argument {\tt sMat} is returned sorted in lexicographic |
---|
124 | ! order by column and row. |
---|
125 | ! |
---|
126 | ! !REVISION HISTORY: |
---|
127 | ! |
---|
128 | ! 13Apr01 - J.W. Larson <larson@mcs.anl.gov> - initial API spec. |
---|
129 | ! 26Apr01 - R.L. Jacob <jacob@mcs.anl.gov> - add use statements for |
---|
130 | ! GlobalSegMap_init and GSMap_peLocs. |
---|
131 | ! Add gsize argument required to GSMap_peLocs. |
---|
132 | ! Add underscore to ComputeSegments call so it matches |
---|
133 | ! the subroutine decleration. |
---|
134 | ! change attribute on starts,lengths, and pe_locs to |
---|
135 | ! pointer to match GSMap_init. |
---|
136 | ! add use m_die statement |
---|
137 | ! 26Apr01 - J.W. Larson <larson@mcs.anl.gov> - fixed major logic bug |
---|
138 | ! that had all processes executing some operations that |
---|
139 | ! should only occur on the root. |
---|
140 | ! 09Jul03 - E.T. Ong <eong@mcs.anl.gov> - call pe_locs in parallel. |
---|
141 | ! reduce the serial sort from gcol:grow to just gcol. |
---|
142 | !EOP |
---|
143 | !------------------------------------------------------------------------- |
---|
144 | |
---|
145 | character(len=*),parameter :: myname_=myname//'ByColumnGSMap_' |
---|
146 | ! Process ID number |
---|
147 | integer :: myID, mySIZE |
---|
148 | ! Attributes for the output GlobalSegMap |
---|
149 | integer :: gsize, comp_id, ngseg |
---|
150 | ! Temporary array for identifying each matrix element column and |
---|
151 | ! process ID destination |
---|
152 | type(AttrVect) :: gcol |
---|
153 | type(AttrVect) :: dist_gcol |
---|
154 | type(AttrVect) :: element_pe_locs |
---|
155 | type(AttrVect) :: dist_element_pe_locs |
---|
156 | ! Index variables for the AttrVects |
---|
157 | integer :: dist_gsize |
---|
158 | integer :: gcol_index |
---|
159 | integer :: element_pe_locs_index |
---|
160 | ! Temporary array for initializing GlobalMap Decomposition |
---|
161 | integer,dimension(:), allocatable :: counts |
---|
162 | ! GlobalMap for setting up decomposition to call pe_locs |
---|
163 | type(GlobalMap) :: dist_GMap |
---|
164 | ! Temporary arrays for matrix GlobalSegMap attributes |
---|
165 | integer, dimension(:), pointer :: starts, lengths, pe_locs |
---|
166 | ! List storage for sorting keys |
---|
167 | type(List) :: sort_keys |
---|
168 | ! Error flag |
---|
169 | integer :: ierr |
---|
170 | ! Loop index |
---|
171 | integer :: i |
---|
172 | |
---|
173 | ! Determine process id number myID |
---|
174 | |
---|
175 | call MPI_COMM_RANK(comm, myID, ierr) |
---|
176 | if(ierr /= 0) then |
---|
177 | call MP_perr_die(myname_,'call MPI_COMM_RANK(...',ierr) |
---|
178 | endif |
---|
179 | |
---|
180 | ! Determine the number of processors in communicator |
---|
181 | |
---|
182 | call MPI_COMM_SIZE(comm, mySIZE, ierr) |
---|
183 | if(ierr /= 0) then |
---|
184 | call MP_perr_die(myname_,'call MPI_COMM_SIZE(...',ierr) |
---|
185 | endif |
---|
186 | |
---|
187 | ! Allocate space for GlobalMap length information |
---|
188 | |
---|
189 | allocate(counts(0:mySIZE-1),stat=ierr) |
---|
190 | if(ierr/=0) call die(myname_,"allocate(counts)",ierr) |
---|
191 | |
---|
192 | ! First step: a lot of prep work on the root only: |
---|
193 | |
---|
194 | if(myID == root) then |
---|
195 | |
---|
196 | ! Sort the matrix entries in sMat by column. |
---|
197 | ! First, create the key list... |
---|
198 | |
---|
199 | call List_init(sort_keys,'gcol') |
---|
200 | |
---|
201 | ! Now perform the sort/permute... |
---|
202 | |
---|
203 | call SparseMatrix_SortPermute(sMat, sort_keys) |
---|
204 | |
---|
205 | call List_clean(sort_keys) |
---|
206 | |
---|
207 | ! The global size of matrix GlobalSegMap is the number nonzero |
---|
208 | ! elements in sMat. |
---|
209 | |
---|
210 | gsize = SparseMatrix_lsize(sMat) |
---|
211 | |
---|
212 | ! Allocate storage space for matrix element column indices and |
---|
213 | ! process ID destinations |
---|
214 | |
---|
215 | call AttrVect_init(aV=gcol, iList="gcol", lsize=gsize) |
---|
216 | |
---|
217 | ! Extract global column information and place in array gCol |
---|
218 | |
---|
219 | call AttrVect_copy(aVin=sMat%data, aVout=gcol, iList="gcol") |
---|
220 | |
---|
221 | ! Setup GlobalMap decomposition lengths: |
---|
222 | |
---|
223 | do i=0,mySIZE-1 |
---|
224 | counts(i) = gsize/mySIZE |
---|
225 | enddo |
---|
226 | counts(mySIZE-1) = counts(mySIZE-1) + mod(gsize,mySIZE) |
---|
227 | |
---|
228 | endif |
---|
229 | |
---|
230 | ! Initialize GlobalMap so that we can scatter the global row |
---|
231 | ! information. The GlobalMap will inherit the component ID |
---|
232 | ! from xGSMap |
---|
233 | |
---|
234 | comp_id = GlobalSegMap_comp_id(xGSMap) |
---|
235 | |
---|
236 | call GlobalMap_init(GMap=dist_GMap, comp_id=comp_id, lns=counts, & |
---|
237 | root=root, comm=comm) |
---|
238 | |
---|
239 | call AttrVect_scatter(iV=gcol, oV=dist_gcol, GMap=dist_GMap, & |
---|
240 | root=root, comm=comm) |
---|
241 | |
---|
242 | ! Similarly, we want to scatter the element_pe_locs using the |
---|
243 | ! same decomposition |
---|
244 | |
---|
245 | dist_gsize = AttrVect_lsize(dist_gcol) |
---|
246 | |
---|
247 | call AttrVect_init(aV=dist_element_pe_locs, iList="element_pe_locs", & |
---|
248 | lsize=dist_gsize) |
---|
249 | call AttrVect_zero(dist_element_pe_locs) |
---|
250 | |
---|
251 | ! Compute process ID destination for each matrix element, |
---|
252 | ! and store in the AttrVect element_pe_locs |
---|
253 | |
---|
254 | gcol_index = AttrVect_indexIA(dist_gcol,"gcol", dieWith=myname_) |
---|
255 | element_pe_locs_index = AttrVect_indexIA(dist_element_pe_locs, & |
---|
256 | "element_pe_locs", dieWith=myname_) |
---|
257 | |
---|
258 | call GlobalSegMap_peLocs(xGSMap, dist_gsize, & |
---|
259 | dist_gcol%iAttr(gcol_index,1:dist_gsize), & |
---|
260 | dist_element_pe_locs%iAttr(element_pe_locs_index,1:dist_gsize)) |
---|
261 | |
---|
262 | call AttrVect_gather(iV=dist_element_pe_locs, oV=element_pe_locs, & |
---|
263 | GMap=dist_GMap, root=root, comm=comm) |
---|
264 | |
---|
265 | ! Back to the root operations |
---|
266 | |
---|
267 | if(myID == root) then |
---|
268 | |
---|
269 | ! Sanity check: Is the globalsize of sMat the same as the |
---|
270 | ! gathered size of element_pe_locs? |
---|
271 | |
---|
272 | if(gsize /= AttrVect_lsize(element_pe_locs)) then |
---|
273 | call die(myname_,"gsize /= AttrVect_lsize(element_pe_locs) & |
---|
274 | & on root process") |
---|
275 | endif |
---|
276 | |
---|
277 | ! Using the entries of gCol and element_pe_locs, build the |
---|
278 | ! output GlobalSegMap attribute arrays starts(:), lengths(:), |
---|
279 | ! and pe_locs(:) |
---|
280 | |
---|
281 | gcol_index = AttrVect_indexIA(gcol,"gcol", dieWith=myname_) |
---|
282 | element_pe_locs_index = AttrVect_indexIA(element_pe_locs, & |
---|
283 | "element_pe_locs", dieWith=myname_) |
---|
284 | |
---|
285 | call ComputeSegments_(element_pe_locs%iAttr(element_pe_locs_index, & |
---|
286 | 1:gsize), & |
---|
287 | gcol%iAttr(gcol_index,1:gsize), & |
---|
288 | gsize, ngseg, starts, lengths, pe_locs) |
---|
289 | ! Clean up on the root |
---|
290 | |
---|
291 | call AttrVect_clean(gcol) |
---|
292 | call AttrVect_clean(element_pe_locs) |
---|
293 | |
---|
294 | endif ! if(myID == root) |
---|
295 | |
---|
296 | ! Non-root processes call GlobalSegMap_init with root_start, |
---|
297 | ! root_length, and root_pe_loc, although these arguments are |
---|
298 | ! not used in the subroutine. Since these correspond to dummy |
---|
299 | ! shaped array arguments in initr_, the Fortran 90 standard |
---|
300 | ! dictates that the actual arguments must contain complete shape |
---|
301 | ! information. Therefore, these array arguments must be |
---|
302 | ! allocated on all processes. |
---|
303 | |
---|
304 | if(myID /= root) then |
---|
305 | allocate(starts(0),lengths(0),pe_locs(0),stat=ierr) |
---|
306 | if(ierr /= 0) then |
---|
307 | call die(myname_,'non-root allocate(starts...',ierr) |
---|
308 | endif |
---|
309 | endif |
---|
310 | |
---|
311 | ! Using this local data on the root, create the SparseMatrix |
---|
312 | ! GlobalSegMap sMGSMap (which will be valid on all processes |
---|
313 | ! on the communicator: |
---|
314 | |
---|
315 | call GlobalSegMap_init(sMGSMap, ngseg, starts, lengths, pe_locs, & |
---|
316 | root, comm, comp_id, gsize) |
---|
317 | |
---|
318 | ! Clean up |
---|
319 | |
---|
320 | call GlobalMap_clean(dist_GMap) |
---|
321 | call AttrVect_clean(dist_gcol) |
---|
322 | call AttrVect_clean(dist_element_pe_locs) |
---|
323 | |
---|
324 | deallocate(starts, lengths, pe_locs, counts, stat=ierr) |
---|
325 | if(ierr /= 0) then |
---|
326 | call die(myname_,'deallocate(starts...',ierr) |
---|
327 | endif |
---|
328 | |
---|
329 | |
---|
330 | end subroutine ByColumnGSMap_ |
---|
331 | |
---|
332 | !------------------------------------------------------------------------- |
---|
333 | ! Math + Computer Science Division / Argonne National Laboratory ! |
---|
334 | !------------------------------------------------------------------------- |
---|
335 | !BOP |
---|
336 | ! |
---|
337 | ! !IROUTINE: ByRowGSMap_ - Generate Row-based GlobalSegMap for SparseMatrix |
---|
338 | ! |
---|
339 | ! !INTERFACE: |
---|
340 | |
---|
341 | subroutine ByRowGSMap_(yGSMap, sMat, sMGSMap, root, comm) |
---|
342 | ! |
---|
343 | ! !USES: |
---|
344 | ! |
---|
345 | |
---|
346 | use m_die, only: MP_perr_die,die |
---|
347 | |
---|
348 | use m_List, only: List |
---|
349 | use m_List, only: List_init => init |
---|
350 | use m_List, only: List_clean => clean |
---|
351 | |
---|
352 | use m_AttrVect, only: AttrVect |
---|
353 | use m_AttrVect, only: AttrVect_init => init |
---|
354 | use m_AttrVect, only: AttrVect_lsize => lsize |
---|
355 | use m_AttrVect, only: AttrVect_indexIA => indexIA |
---|
356 | use m_AttrVect, only: AttrVect_copy => copy |
---|
357 | use m_AttrVect, only: AttrVect_clean => clean |
---|
358 | use m_AttrVect, only: AttrVect_zero => zero |
---|
359 | |
---|
360 | use m_AttrVectComms, only: AttrVect_scatter => scatter |
---|
361 | use m_AttrVectComms, only: AttrVect_gather => gather |
---|
362 | |
---|
363 | use m_GlobalMap, only : GlobalMap |
---|
364 | use m_GlobalMap, only : GlobalMap_init => init |
---|
365 | use m_GlobalMap, only : GlobalMap_clean => clean |
---|
366 | |
---|
367 | use m_GlobalSegMap, only: GlobalSegMap |
---|
368 | use m_GlobalSegMap, only: GlobalSegMap_init => init |
---|
369 | use m_GlobalSegMap, only: GlobalSegMap_peLocs => peLocs |
---|
370 | use m_GlobalSegMap, only: GlobalSegMap_comp_id => comp_id |
---|
371 | |
---|
372 | use m_SparseMatrix, only: SparseMatrix |
---|
373 | use m_SparseMatrix, only: SparseMatrix_lsize => lsize |
---|
374 | use m_SparseMatrix, only: SparseMatrix_SortPermute => SortPermute |
---|
375 | |
---|
376 | implicit none |
---|
377 | |
---|
378 | ! !INPUT PARAMETERS: |
---|
379 | ! |
---|
380 | type(GlobalSegMap), intent(in) :: yGSMap |
---|
381 | integer, intent(in) :: root |
---|
382 | integer, intent(in) :: comm |
---|
383 | |
---|
384 | ! !INPUT/OUTPUT PARAMETERS: |
---|
385 | ! |
---|
386 | type(SparseMatrix), intent(inout) :: sMat |
---|
387 | |
---|
388 | ! !OUTPUT PARAMETERS: |
---|
389 | ! |
---|
390 | type(GlobalSegMap), intent(out) :: sMGSMap |
---|
391 | |
---|
392 | ! !DESCRIPTION: This routine is invoked from all processes on the |
---|
393 | ! communicator {\tt comm} to create from an input {\tt SparseMatrix} |
---|
394 | ! {\tt sMat} (valid only on the {\tt root} process) and an input |
---|
395 | ! {\bf y}-vector decomposition described by the {\tt GlobalSegMap} |
---|
396 | ! argument {\tt yGSMap} (valid at least on the {\tt root}) to create |
---|
397 | ! an output {\tt GlobalSegMap} decomposition of the matrix elements |
---|
398 | ! {\tt sMGSMap}, which is valid on all processes on the communicator. |
---|
399 | ! This matrix {\tt GlobalSegMap} describes the corresponding row |
---|
400 | ! decomposition of {\tt sMat}. |
---|
401 | ! |
---|
402 | ! {\bf N.B.}: The argument {\tt sMat} is returned sorted in lexicographic |
---|
403 | ! order by row and column. |
---|
404 | ! |
---|
405 | ! !REVISION HISTORY: |
---|
406 | ! |
---|
407 | ! 13Apr01 - J.W. Larson <larson@mcs.anl.gov> - initial API spec. |
---|
408 | ! 26Apr01 - R.L. Jacob <jacob@mcs.anl.gov> - add use statements for |
---|
409 | ! GlobalSegMap_init and GSMap_peLocs. |
---|
410 | ! Add gsize argument required to GSMap_peLocs. |
---|
411 | ! Add underscore to ComputeSegments call so it matches |
---|
412 | ! the subroutine decleration. |
---|
413 | ! change attribute on starts,lengths, and pe_locs to |
---|
414 | ! pointer to match GSMap_init. |
---|
415 | ! 26Apr01 - J.W. Larson <larson@mcs.anl.gov> - fixed major logic bug |
---|
416 | ! that had all processes executing some operations that |
---|
417 | ! should only occur on the root. |
---|
418 | ! 09Jun03 - E.T. Ong <eong@mcs.anl.gov> - call peLocs in parallel. |
---|
419 | ! reduce the serial sort from grow:gcol to just grow. |
---|
420 | !EOP |
---|
421 | !------------------------------------------------------------------------- |
---|
422 | |
---|
423 | character(len=*),parameter :: myname_=myname//'ByRowGSMap_' |
---|
424 | ! Process ID number and communicator size |
---|
425 | integer :: myID, mySIZE |
---|
426 | ! Attributes for the output GlobalSegMap |
---|
427 | integer :: gsize, comp_id, ngseg |
---|
428 | ! Temporary array for identifying each matrix element row and |
---|
429 | ! process ID destination |
---|
430 | type(AttrVect) :: grow |
---|
431 | type(AttrVect) :: dist_grow |
---|
432 | type(AttrVect) :: element_pe_locs |
---|
433 | type(AttrVect) :: dist_element_pe_locs |
---|
434 | ! Index variables for AttrVects |
---|
435 | integer :: dist_gsize |
---|
436 | integer :: grow_index |
---|
437 | integer :: element_pe_locs_index |
---|
438 | ! Temporary array for initializing GlobalMap Decomposition |
---|
439 | integer,dimension(:), allocatable :: counts |
---|
440 | ! GlobalMap for setting up decomposition to call pe_locs |
---|
441 | type(GlobalMap) :: dist_GMap |
---|
442 | ! Temporary arrays for matrix GlobalSegMap attributes |
---|
443 | integer, dimension(:), pointer :: starts, lengths, pe_locs |
---|
444 | ! List storage for sorting keys |
---|
445 | type(List) :: sort_keys |
---|
446 | ! Error flag |
---|
447 | integer :: ierr |
---|
448 | ! Loop index |
---|
449 | integer :: i |
---|
450 | |
---|
451 | ! Determine process id number myID |
---|
452 | |
---|
453 | call MPI_COMM_RANK(comm, myID, ierr) |
---|
454 | if(ierr /= 0) then |
---|
455 | call MP_perr_die(myname_,'call MPI_COMM_RANK(...',ierr) |
---|
456 | endif |
---|
457 | |
---|
458 | ! Determine the number of processors in communicator |
---|
459 | |
---|
460 | call MPI_COMM_SIZE(comm, mySIZE, ierr) |
---|
461 | if(ierr /= 0) then |
---|
462 | call MP_perr_die(myname_,'call MPI_COMM_SIZE(...',ierr) |
---|
463 | endif |
---|
464 | |
---|
465 | ! Allocate space for GlobalMap length information |
---|
466 | |
---|
467 | allocate(counts(0:mySIZE-1),stat=ierr) |
---|
468 | if(ierr/=0) call die(myname_,"allocate(counts)",ierr) |
---|
469 | |
---|
470 | ! First step: a lot of prep work on the root only: |
---|
471 | |
---|
472 | if(myID == root) then |
---|
473 | |
---|
474 | ! Sort the matrix entries in sMat by row. |
---|
475 | ! First, create the key list... |
---|
476 | |
---|
477 | call List_init(sort_keys,'grow') |
---|
478 | |
---|
479 | ! Now perform the sort/permute... |
---|
480 | |
---|
481 | call SparseMatrix_SortPermute(sMat, sort_keys) |
---|
482 | |
---|
483 | call List_clean(sort_keys) |
---|
484 | |
---|
485 | ! The global size of matrix GlobalSegMap is the number of rows. |
---|
486 | |
---|
487 | gsize = SparseMatrix_lsize(sMat) |
---|
488 | |
---|
489 | ! Allocate storage space for matrix element row indices and |
---|
490 | ! process ID destinations |
---|
491 | |
---|
492 | call AttrVect_init(aV=grow, iList="grow", lsize=gsize) |
---|
493 | |
---|
494 | ! Extract global row information and place in AttrVect grow |
---|
495 | |
---|
496 | call AttrVect_copy(aVin=sMat%data, aVout=grow, iList="grow") |
---|
497 | |
---|
498 | ! Setup GlobalMap decomposition lengths: |
---|
499 | ! Give any extra points to the last process |
---|
500 | |
---|
501 | do i=0,mySIZE-1 |
---|
502 | counts(i) = gsize/mySIZE |
---|
503 | enddo |
---|
504 | counts(mySIZE-1) = counts(mySIZE-1) + mod(gsize,mySIZE) |
---|
505 | |
---|
506 | endif |
---|
507 | |
---|
508 | ! Initialize GlobalMap and scatter the global row information. |
---|
509 | ! The GlobalMap will inherit the component ID from yGSMap |
---|
510 | |
---|
511 | comp_id = GlobalSegMap_comp_id(yGSMap) |
---|
512 | |
---|
513 | call GlobalMap_init(GMap=dist_GMap, comp_id=comp_id, lns=counts, & |
---|
514 | root=root, comm=comm) |
---|
515 | |
---|
516 | call AttrVect_scatter(iV=grow, oV=dist_grow, GMap=dist_GMap, & |
---|
517 | root=root, comm=comm) |
---|
518 | |
---|
519 | ! Similarly, we want to scatter the element_pe_locs using the |
---|
520 | ! same decomposition |
---|
521 | |
---|
522 | dist_gsize = AttrVect_lsize(dist_grow) |
---|
523 | |
---|
524 | call AttrVect_init(aV=dist_element_pe_locs, iList="element_pe_locs", & |
---|
525 | lsize=dist_gsize) |
---|
526 | call AttrVect_zero(dist_element_pe_locs) |
---|
527 | |
---|
528 | ! Compute process ID destination for each matrix element, |
---|
529 | ! and store in the AttrVect element_pe_locs |
---|
530 | |
---|
531 | grow_index = AttrVect_indexIA(dist_grow,"grow", dieWith=myname_) |
---|
532 | element_pe_locs_index = AttrVect_indexIA(dist_element_pe_locs, & |
---|
533 | "element_pe_locs", dieWith=myname_) |
---|
534 | |
---|
535 | call GlobalSegMap_peLocs(yGSMap, dist_gsize, & |
---|
536 | dist_grow%iAttr(grow_index,1:dist_gsize), & |
---|
537 | dist_element_pe_locs%iAttr(element_pe_locs_index,1:dist_gsize)) |
---|
538 | |
---|
539 | ! Gather element_pe_locs on root so that we can call compute_segments |
---|
540 | |
---|
541 | call AttrVect_gather(iV=dist_element_pe_locs, oV=element_pe_locs, & |
---|
542 | GMap=dist_GMap, root=root, comm=comm) |
---|
543 | |
---|
544 | ! Back to the root operations |
---|
545 | |
---|
546 | if(myID == root) then |
---|
547 | |
---|
548 | ! Sanity check: Is the globalsize of sMat the same as the |
---|
549 | ! gathered size of element_pe_locs? |
---|
550 | |
---|
551 | if(gsize /= AttrVect_lsize(element_pe_locs)) then |
---|
552 | call die(myname_,"gsize /= AttrVect_lsize(element_pe_locs) & |
---|
553 | & on root process") |
---|
554 | endif |
---|
555 | |
---|
556 | ! Using the entries of grow and element_pe_locs, build the |
---|
557 | ! output GlobalSegMap attribute arrays starts(:), lengths(:), |
---|
558 | ! and pe_locs(:) |
---|
559 | |
---|
560 | grow_index = AttrVect_indexIA(grow,"grow", dieWith=myname_) |
---|
561 | element_pe_locs_index = AttrVect_indexIA(element_pe_locs, & |
---|
562 | "element_pe_locs", dieWith=myname_) |
---|
563 | |
---|
564 | call ComputeSegments_(element_pe_locs%iAttr(element_pe_locs_index, & |
---|
565 | 1:gsize), & |
---|
566 | grow%iAttr(grow_index,1:gsize), & |
---|
567 | gsize, ngseg, starts, lengths, pe_locs) |
---|
568 | |
---|
569 | ! Clean up on the root |
---|
570 | |
---|
571 | call AttrVect_clean(grow) |
---|
572 | call AttrVect_clean(element_pe_locs) |
---|
573 | |
---|
574 | endif ! if(myID == root) |
---|
575 | |
---|
576 | ! Non-root processes call GlobalSegMap_init with root_start, |
---|
577 | ! root_length, and root_pe_loc, although these arguments are |
---|
578 | ! not used in the subroutine. Since these correspond to dummy |
---|
579 | ! shaped array arguments in initr_, the Fortran 90 standard |
---|
580 | ! dictates that the actual arguments must contain complete shape |
---|
581 | ! information. Therefore, these array arguments must be |
---|
582 | ! allocated on all processes. |
---|
583 | |
---|
584 | if(myID /= root) then |
---|
585 | allocate(starts(0),lengths(0),pe_locs(0),stat=ierr) |
---|
586 | if(ierr /= 0) then |
---|
587 | call die(myname_,'non-root allocate(starts...',ierr) |
---|
588 | endif |
---|
589 | endif |
---|
590 | |
---|
591 | ! Using this local data on the root, create the SparseMatrix |
---|
592 | ! GlobalSegMap sMGSMap (which will be valid on all processes |
---|
593 | ! on the communicator. The GlobalSegMap will inherit the |
---|
594 | ! component ID from yGSMap |
---|
595 | |
---|
596 | call GlobalSegMap_init(sMGSMap, ngseg, starts, lengths, pe_locs, & |
---|
597 | root, comm, comp_id, gsize) |
---|
598 | |
---|
599 | ! Clean up: |
---|
600 | |
---|
601 | call GlobalMap_clean(dist_GMap) |
---|
602 | call AttrVect_clean(dist_grow) |
---|
603 | call AttrVect_clean(dist_element_pe_locs) |
---|
604 | |
---|
605 | deallocate(starts, lengths, pe_locs, counts, stat=ierr) |
---|
606 | if(ierr /= 0) then |
---|
607 | call die(myname_,'deallocate(starts...',ierr) |
---|
608 | endif |
---|
609 | |
---|
610 | |
---|
611 | end subroutine ByRowGSMap_ |
---|
612 | |
---|
613 | !------------------------------------------------------------------------- |
---|
614 | ! Math + Computer Science Division / Argonne National Laboratory ! |
---|
615 | !------------------------------------------------------------------------- |
---|
616 | !BOP |
---|
617 | ! |
---|
618 | ! !IROUTINE: ComputeSegments_ - Create segments from list data. |
---|
619 | ! |
---|
620 | ! !INTERFACE: |
---|
621 | |
---|
622 | subroutine ComputeSegments_(element_pe_locs, elements, num_elements, & |
---|
623 | nsegs, seg_starts, seg_lengths, seg_pe_locs) |
---|
624 | ! |
---|
625 | ! !USES: |
---|
626 | ! |
---|
627 | |
---|
628 | use m_die, only: die |
---|
629 | |
---|
630 | implicit none |
---|
631 | |
---|
632 | ! !INPUT PARAMETERS: |
---|
633 | ! |
---|
634 | integer, dimension(:), intent(in) :: element_pe_locs |
---|
635 | integer, dimension(:), intent(in) :: elements |
---|
636 | integer, intent(in) :: num_elements |
---|
637 | |
---|
638 | ! !OUTPUT PARAMETERS: |
---|
639 | ! |
---|
640 | integer, intent(out) :: nsegs |
---|
641 | integer, dimension(:), pointer :: seg_starts |
---|
642 | integer, dimension(:), pointer :: seg_lengths |
---|
643 | integer, dimension(:), pointer :: seg_pe_locs |
---|
644 | |
---|
645 | ! !DESCRIPTION: This routine examins an input list of {\tt num\_elements} |
---|
646 | ! process ID locations stored in the array {\tt element\_pe\_locs}, counts |
---|
647 | ! the number of contiguous segments {\tt nsegs}, and returns the segment |
---|
648 | ! start index, length, and process ID location in the arrays {\tt seg\_starts(:)}, |
---|
649 | ! {\tt seg\_lengths(:)}, and {\tt seg\_pe\_locs(:)}, respectively. |
---|
650 | ! |
---|
651 | ! {\bf N.B.}: The argument {\tt sMat} is returned sorted in lexicographic |
---|
652 | ! order by row and column. |
---|
653 | ! |
---|
654 | ! !REVISION HISTORY: |
---|
655 | ! |
---|
656 | ! 18Apr01 - J.W. Larson <larson@mcs.anl.gov> - initial version. |
---|
657 | ! 28Aug01 - M.J. Zavislak <zavislak@mcs.anl.gov> |
---|
658 | ! Changed first sanity check to get size(element_pe_locs) |
---|
659 | ! instead of size(elements) |
---|
660 | !EOP |
---|
661 | !------------------------------------------------------------------------- |
---|
662 | character(len=*),parameter :: myname_=myname//'ComputeSegments_' |
---|
663 | |
---|
664 | integer :: i, ierr, iseg |
---|
665 | |
---|
666 | ! Input argument sanity checks: |
---|
667 | |
---|
668 | if(size(element_pe_locs) < num_elements) then |
---|
669 | call die(myname_,'input argument array element_pe_locs too small', & |
---|
670 | num_elements-size(element_pe_locs)) |
---|
671 | endif |
---|
672 | |
---|
673 | if(size(elements) < num_elements) then |
---|
674 | call die(myname_,'input argument array elements too small', & |
---|
675 | num_elements-size(elements)) |
---|
676 | endif |
---|
677 | |
---|
678 | ! First pass: how many segments? |
---|
679 | |
---|
680 | do i=1,num_elements |
---|
681 | |
---|
682 | if(i == 1) then ! bootstrap segment count |
---|
683 | |
---|
684 | nsegs = 1 |
---|
685 | |
---|
686 | else ! usual point/segment processing |
---|
687 | |
---|
688 | ! New segment? If so, increment nsegs. |
---|
689 | |
---|
690 | if((elements(i) > elements(i-1) + 1) .or. & |
---|
691 | (element_pe_locs(i) /= element_pe_locs(i-1))) then ! new segment |
---|
692 | nsegs = nsegs + 1 |
---|
693 | endif |
---|
694 | |
---|
695 | endif ! if(i == 1) block |
---|
696 | |
---|
697 | end do ! do i=1,num_elements |
---|
698 | |
---|
699 | allocate(seg_starts(nsegs), seg_lengths(nsegs), seg_pe_locs(nsegs), & |
---|
700 | stat=ierr) |
---|
701 | |
---|
702 | if(ierr /= 0) then |
---|
703 | call die(myname_,'allocate(seg_starts...',ierr) |
---|
704 | endif |
---|
705 | |
---|
706 | ! Second pass: fill in segment data. |
---|
707 | |
---|
708 | ! NOTE: Structure of this loop was changed from a for loop |
---|
709 | ! to avoid a faulty vectorization on the SUPER-UX compiler |
---|
710 | |
---|
711 | i=1 |
---|
712 | ASSIGN_LOOP: do |
---|
713 | |
---|
714 | if(i == 1) then ! bootstrap first segment info. |
---|
715 | |
---|
716 | iseg = 1 |
---|
717 | seg_starts(iseg) = 1 |
---|
718 | seg_lengths(iseg) = 1 |
---|
719 | seg_pe_locs(iseg) = element_pe_locs(iseg) |
---|
720 | |
---|
721 | else ! do usual point/segment processing |
---|
722 | |
---|
723 | ! New segment? This happens if 1) elements(i) > elements(i-1) + 1, or |
---|
724 | ! 2) element_pe_locs(i) /= element_pe_locs(i-1). |
---|
725 | |
---|
726 | if((elements(i) > elements(i-1) + 1) .or. & |
---|
727 | (element_pe_locs(i) /= element_pe_locs(i-1))) then ! new segment |
---|
728 | |
---|
729 | ! Initialize new segment |
---|
730 | iseg = iseg + 1 |
---|
731 | seg_starts(iseg) = i |
---|
732 | seg_lengths(iseg) = 1 |
---|
733 | seg_pe_locs(iseg) = element_pe_locs(i) |
---|
734 | |
---|
735 | else |
---|
736 | |
---|
737 | ! Increment current segment length |
---|
738 | seg_lengths(iseg) = seg_lengths(iseg) + 1 |
---|
739 | |
---|
740 | endif ! If new segment block |
---|
741 | |
---|
742 | endif ! if(i == 1) block |
---|
743 | |
---|
744 | ! Prepare index i for the next loop around; |
---|
745 | if(i>=num_elements) EXIT |
---|
746 | i = i + 1 |
---|
747 | |
---|
748 | end do ASSIGN_LOOP |
---|
749 | |
---|
750 | if(iseg /= nsegs) then |
---|
751 | call die(myname_,'segment number difference',iseg-nsegs) |
---|
752 | endif |
---|
753 | |
---|
754 | end subroutine ComputeSegments_ |
---|
755 | |
---|
756 | end module m_SparseMatrixDecomp |
---|