1 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
2 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
3 | !----------------------------------------------------------------------- |
---|
4 | ! CVS m_SparseMatrixToMaps.F90,v 1.13 2007-05-11 02:09:48 rloy Exp |
---|
5 | ! CVS MCT_2_8_0 |
---|
6 | !BOP ------------------------------------------------------------------- |
---|
7 | ! |
---|
8 | ! !MODULE: m_SparseMatrixToMaps -- Maps from the Sparse Matrix |
---|
9 | ! |
---|
10 | ! !DESCRIPTION: |
---|
11 | ! The {\tt SparseMatrix} provides consolidated (on one process) or |
---|
12 | ! distributed sparse matrix storage for the operation |
---|
13 | ! ${\bf y} = {\bf M} {\bf x}$, where {\bf x} and {\bf y} are vectors, |
---|
14 | ! and {\bf M} is a matrix. In performing parallel matrix-vector |
---|
15 | ! multiplication, one has numerous options regarding the decomposition |
---|
16 | ! of the matrix {\bf M}, and the vectors {\bf y} and {\bf x}. |
---|
17 | ! This module provides services to generate mct mapping components---the |
---|
18 | ! {\tt GlobalMap} and {\tt GlobalSegMap} for the vectors {\bf y} and/or |
---|
19 | ! {\bf x} based on the decomposition of the sparse matrix {\bf M}. |
---|
20 | ! |
---|
21 | ! !INTERFACE: |
---|
22 | |
---|
23 | module m_SparseMatrixToMaps |
---|
24 | ! |
---|
25 | ! !USES: |
---|
26 | ! |
---|
27 | use m_SparseMatrix, only : SparseMatrix |
---|
28 | |
---|
29 | implicit none |
---|
30 | |
---|
31 | private ! except |
---|
32 | |
---|
33 | public :: SparseMatrixToXGlobalSegMap |
---|
34 | public :: SparseMatrixToYGlobalSegMap |
---|
35 | |
---|
36 | interface SparseMatrixToXGlobalSegMap ; module procedure & |
---|
37 | SparseMatrixToXGlobalSegMap_ |
---|
38 | end interface |
---|
39 | |
---|
40 | interface SparseMatrixToYGlobalSegMap ; module procedure & |
---|
41 | SparseMatrixToYGlobalSegMap_ |
---|
42 | end interface |
---|
43 | |
---|
44 | ! !REVISION HISTORY: |
---|
45 | ! 13Apr01 - J.W. Larson <larson@mcs.anl.gov> - initial prototype |
---|
46 | ! and API specifications. |
---|
47 | !EOP ___________________________________________________________________ |
---|
48 | |
---|
49 | character(len=*),parameter :: myname='MCT::m_SparseMatrixToMaps' |
---|
50 | |
---|
51 | contains |
---|
52 | |
---|
53 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
54 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
55 | !BOP ------------------------------------------------------------------- |
---|
56 | ! |
---|
57 | ! !IROUTINE: SparseMatrixToXGlobalSegMap_ - Generate X GlobalSegmap. |
---|
58 | ! |
---|
59 | ! !DESCRIPTION: Given an input {\tt SparseMatrix} argument {\tt sMat}, |
---|
60 | ! this routine generates an output {\tt GlobalSegMap} variable |
---|
61 | ! {\tt xGSMap}, which describes the domain decomposition of the vector |
---|
62 | ! {\bf x} in the distributed matrix-vector multiplication |
---|
63 | ! $${\bf y} = {\bf M} {\bf x}.$$ |
---|
64 | ! |
---|
65 | ! !INTERFACE: |
---|
66 | |
---|
67 | subroutine SparseMatrixToXGlobalSegMap_(sMat, xGSMap, root, comm, comp_id) |
---|
68 | ! |
---|
69 | ! !USES: |
---|
70 | ! |
---|
71 | use m_stdio, only : stderr |
---|
72 | use m_die, only : die |
---|
73 | use m_mpif90 |
---|
74 | |
---|
75 | use m_List, only : List |
---|
76 | use m_List, only : List_init => init |
---|
77 | use m_List, only : List_clean => clean |
---|
78 | |
---|
79 | use m_SparseMatrix, only : SparseMatrix |
---|
80 | use m_SparseMatrix, only : SparseMatrix_nCols => nCols |
---|
81 | use m_SparseMatrix, only : SparseMatrix_lsize => lsize |
---|
82 | use m_SparseMatrix, only : SparseMatrix_indexIA => indexIA |
---|
83 | use m_SparseMatrix, only : SparseMatrix_SortPermute => SortPermute |
---|
84 | |
---|
85 | use m_GlobalSegMap, only : GlobalSegMap |
---|
86 | use m_GlobalSegMap, only : GlobalSegMap_init => init |
---|
87 | |
---|
88 | implicit none |
---|
89 | |
---|
90 | ! !INPUT PARAMETERS: |
---|
91 | ! |
---|
92 | integer, intent(in) :: root ! communicator root |
---|
93 | integer, intent(in) :: comm ! communicator handle |
---|
94 | integer, intent(in) :: comp_id ! component id |
---|
95 | |
---|
96 | ! !INPUT/OUTPUT PARAMETERS: |
---|
97 | ! |
---|
98 | type(SparseMatrix), intent(inout) :: sMat ! input SparseMatrix |
---|
99 | |
---|
100 | ! !OUTPUT PARAMETERS: |
---|
101 | ! |
---|
102 | type(GlobalSegMap), intent(out) :: xGSMap ! segmented decomposition |
---|
103 | ! for x |
---|
104 | ! !REVISION HISTORY: |
---|
105 | ! 13Apr01 - J.W. Larson <larson@mcs.anl.gov> - API specification. |
---|
106 | ! 25Apr01 - J.W. Larson <larson@mcs.anl.gov> - First version. |
---|
107 | ! 27Apr01 - J.W. Larson <larson@mcs.anl.gov> - Bug fix--intent of |
---|
108 | ! argument sMat changed from (IN) to (INOUT) |
---|
109 | ! 27Apr01 - R.L. Jacob <jacob@mcs.anl.gov> - bug fix-- add use |
---|
110 | ! statement for SortPermute |
---|
111 | ! 01May01 - R.L. Jacob <jacob@mcs.anl.gov> - make comp_id an |
---|
112 | ! input argument |
---|
113 | !EOP ___________________________________________________________________ |
---|
114 | |
---|
115 | character(len=*),parameter :: myname_=myname//'::SparseMatrixToXGlobalSegMap_' |
---|
116 | |
---|
117 | ! SparseMatrix attributes: |
---|
118 | integer :: lsize |
---|
119 | ! GlobalSegMap input attributes: |
---|
120 | integer :: gsize, ngseg |
---|
121 | integer, dimension(:), pointer :: starts, lengths |
---|
122 | ! Temporary array for identifying each matrix element column and |
---|
123 | ! process ID destination |
---|
124 | integer, dimension(:), allocatable :: gCol, element_pe_locs |
---|
125 | ! Index to identify the gcol attribute in sMat: |
---|
126 | integer :: igCol |
---|
127 | ! Matrix element sorting keys list: |
---|
128 | type(List) :: sort_keys |
---|
129 | ! Loop index and error flag: |
---|
130 | integer :: i, ierr |
---|
131 | |
---|
132 | ! Determine he local number of matrix elements lsize |
---|
133 | |
---|
134 | lsize = SparseMatrix_lsize(sMat) |
---|
135 | |
---|
136 | ! The value of gsize is taken from the number of columns in sMat: |
---|
137 | |
---|
138 | gsize = SparseMatrix_nCols(sMat) |
---|
139 | |
---|
140 | ! Sort SparseMatrix entries by global column index gcol, then |
---|
141 | ! global row index. |
---|
142 | |
---|
143 | ! Create Sort keys list |
---|
144 | |
---|
145 | call List_init(sort_keys,'gcol:grow') |
---|
146 | |
---|
147 | ! Sort and permute the entries of sMat into lexicographic order |
---|
148 | ! by global column, then global row. |
---|
149 | |
---|
150 | call SparseMatrix_SortPermute(sMat, sort_keys) |
---|
151 | |
---|
152 | ! Clean up sort keys list |
---|
153 | |
---|
154 | call List_clean(sort_keys) |
---|
155 | |
---|
156 | ! Allocate storage space for matrix element column indices and |
---|
157 | ! process ID destinations |
---|
158 | |
---|
159 | allocate(gCol(lsize), stat=ierr) |
---|
160 | |
---|
161 | if(ierr /= 0) then |
---|
162 | call die(myname_,'allocate(gCol...',ierr) |
---|
163 | endif |
---|
164 | |
---|
165 | ! Extract global column information and place in array gCol |
---|
166 | |
---|
167 | igCol = SparseMatrix_indexIA(sMat, 'gcol', dieWith=myname_) |
---|
168 | |
---|
169 | do i=1, lsize |
---|
170 | gCol(i) = sMat%data%iAttr(igCol,i) |
---|
171 | end do |
---|
172 | |
---|
173 | ! Scan sorted entries of gCol to count segments (ngseg), and |
---|
174 | ! their starting indices and lengths (returned in the arrays |
---|
175 | ! starts(:) and lengths(:), respectively) |
---|
176 | |
---|
177 | call ComputeSegments_(gCol, lsize, ngseg, starts, lengths) |
---|
178 | |
---|
179 | ! Now we have sufficient data to call the GlobalSegMap |
---|
180 | ! initialization using distributed data: |
---|
181 | |
---|
182 | call GlobalSegMap_init(xGSMap, starts, lengths, root, comm, & |
---|
183 | comp_id, gsize=gsize) |
---|
184 | |
---|
185 | ! clean up temporary arrays gCol(:), starts(:) and lengths(:), |
---|
186 | ! (the latter two were allocated in the call to the routine |
---|
187 | ! ComputeSegments_()) |
---|
188 | |
---|
189 | deallocate(gCol, starts, lengths, stat=ierr) |
---|
190 | |
---|
191 | if(ierr /= 0) then |
---|
192 | call die(myname_,'deallocate(gCol...',ierr) |
---|
193 | endif |
---|
194 | |
---|
195 | end subroutine SparseMatrixToXGlobalSegMap_ |
---|
196 | |
---|
197 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
198 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
199 | !BOP ------------------------------------------------------------------- |
---|
200 | ! |
---|
201 | ! !IROUTINE: SparseMatrixToYGlobalSegMap_ - Generate Y GlobalSegmap. |
---|
202 | ! |
---|
203 | ! !DESCRIPTION: Given an input {\tt SparseMatrix} argument {\tt sMat}, |
---|
204 | ! this routine generates an output {\tt GlobalSegMap} variable |
---|
205 | ! {\tt yGSMap}, which describes the domain decomposition of the vector |
---|
206 | ! {\bf y} in the distributed matrix-vector multiplication |
---|
207 | ! ${\bf y} = {\bf M} {\bf x}$. |
---|
208 | ! |
---|
209 | ! !INTERFACE: |
---|
210 | |
---|
211 | subroutine SparseMatrixToYGlobalSegMap_(sMat, yGSMap, root, comm, comp_id) |
---|
212 | ! |
---|
213 | ! !USES: |
---|
214 | ! |
---|
215 | use m_stdio, only : stderr |
---|
216 | use m_die, only : die |
---|
217 | |
---|
218 | use m_List, only : List |
---|
219 | use m_List, only : List_init => init |
---|
220 | use m_List, only : List_clean => clean |
---|
221 | |
---|
222 | use m_SparseMatrix, only : SparseMatrix |
---|
223 | use m_SparseMatrix, only : SparseMatrix_nRows => nRows |
---|
224 | use m_SparseMatrix, only : SparseMatrix_lsize => lsize |
---|
225 | use m_SparseMatrix, only : SparseMatrix_indexIA => indexIA |
---|
226 | use m_SparseMatrix, only : SparseMatrix_SortPermute => SortPermute |
---|
227 | |
---|
228 | use m_GlobalSegMap, only : GlobalSegMap |
---|
229 | use m_GlobalSegMap, only : GlobalSegMap_init => init |
---|
230 | |
---|
231 | implicit none |
---|
232 | |
---|
233 | ! !INPUT PARAMETERS: |
---|
234 | ! |
---|
235 | integer, intent(in) :: root ! communicator root |
---|
236 | integer, intent(in) :: comm ! communicator handle |
---|
237 | integer, intent(in) :: comp_id ! component id |
---|
238 | |
---|
239 | ! !INPUT/OUTPUT PARAMETERS: |
---|
240 | ! |
---|
241 | type(SparseMatrix), intent(inout) :: sMat ! input SparseMatrix |
---|
242 | |
---|
243 | ! !OUTPUT PARAMETERS: |
---|
244 | ! |
---|
245 | type(GlobalSegMap), intent(out) :: yGSMap ! segmented decomposition |
---|
246 | ! for y |
---|
247 | ! !REVISION HISTORY: |
---|
248 | ! 13Apr01 - J.W. Larson <larson@mcs.anl.gov> - API specification. |
---|
249 | ! 25Apr01 - J.W. Larson <larson@mcs.anl.gov> - initial code. |
---|
250 | ! 27Apr01 - J.W. Larson <larson@mcs.anl.gov> - Bug fix--intent of |
---|
251 | ! argument sMat changed from (IN) to (INOUT) |
---|
252 | ! 27Apr01 - R.L. Jacob <jacob@mcs.anl.gov> - bug fix-- add use |
---|
253 | ! statement for SortPermute |
---|
254 | ! 01May01 - R.L. Jacob <jacob@mcs.anl.gov> - make comp_id an |
---|
255 | ! input argument |
---|
256 | ! 07May02 - J.W. Larson <larson@mcs.anl.gov> - Changed interface to |
---|
257 | ! make it consistent with SparseMatrixToXGlobalSegMap_(). |
---|
258 | !EOP ___________________________________________________________________ |
---|
259 | |
---|
260 | character(len=*),parameter :: myname_=myname//'::SparseMatrixToYGlobalSegMap_' |
---|
261 | |
---|
262 | ! SparseMatrix attributes: |
---|
263 | integer :: lsize |
---|
264 | ! GlobalSegMap input attributes: |
---|
265 | integer :: gsize, ngseg |
---|
266 | integer, dimension(:), pointer :: starts, lengths |
---|
267 | ! Temporary array for identifying each matrix element column and |
---|
268 | ! process ID destination |
---|
269 | integer, dimension(:), allocatable :: gRow, element_pe_locs |
---|
270 | ! Index to identify the gRow attribute in sMat: |
---|
271 | integer :: igRow |
---|
272 | ! Matrix element sorting keys list: |
---|
273 | type(List) :: sort_keys |
---|
274 | ! Loop index and error flag: |
---|
275 | integer :: i, ierr |
---|
276 | |
---|
277 | ! Determine he local number of matrix elements lsize |
---|
278 | |
---|
279 | lsize = SparseMatrix_lsize(sMat) |
---|
280 | |
---|
281 | ! The value of gsize is taken from the number of columns in sMat: |
---|
282 | |
---|
283 | gsize = SparseMatrix_nRows(sMat) |
---|
284 | |
---|
285 | ! Sort SparseMatrix entries by global column index grow, then |
---|
286 | ! global row index. |
---|
287 | |
---|
288 | ! Create Sort keys list |
---|
289 | |
---|
290 | call List_init(sort_keys,'grow:gcol') |
---|
291 | |
---|
292 | ! Sort and permute the entries of sMat into lexicographic order |
---|
293 | ! by global column, then global row. |
---|
294 | |
---|
295 | call SparseMatrix_SortPermute(sMat, sort_keys) |
---|
296 | |
---|
297 | ! Clean up sort keys list |
---|
298 | |
---|
299 | call List_clean(sort_keys) |
---|
300 | |
---|
301 | ! Allocate storage space for matrix element column indices and |
---|
302 | ! process ID destinations |
---|
303 | |
---|
304 | allocate(gRow(lsize), stat=ierr) |
---|
305 | |
---|
306 | if(ierr /= 0) then |
---|
307 | call die(myname_,'allocate(gRow...',ierr) |
---|
308 | endif |
---|
309 | |
---|
310 | ! Extract global column information and place in array gRow |
---|
311 | |
---|
312 | igRow = SparseMatrix_indexIA(sMat,'grow', dieWith=myname_) |
---|
313 | |
---|
314 | do i=1, lsize |
---|
315 | gRow(i) = sMat%data%iAttr(igRow,i) |
---|
316 | end do |
---|
317 | |
---|
318 | ! Scan sorted entries of gRow to count segments (ngseg), and |
---|
319 | ! their starting indices and lengths (returned in the arrays |
---|
320 | ! starts(:) and lengths(:), respectively) |
---|
321 | |
---|
322 | call ComputeSegments_(gRow, lsize, ngseg, starts, lengths) |
---|
323 | |
---|
324 | ! Now we have sufficient data to call the GlobalSegMap |
---|
325 | ! initialization using distributed data: |
---|
326 | |
---|
327 | call GlobalSegMap_init(yGSMap, starts, lengths, root, comm, & |
---|
328 | comp_id, gsize=gsize) |
---|
329 | |
---|
330 | ! clean up temporary arrays gRow(:), starts(:) and lengths(:), |
---|
331 | ! (the latter two were allocated in the call to the routine |
---|
332 | ! ComputeSegments_()) |
---|
333 | |
---|
334 | deallocate(gRow, starts, lengths, stat=ierr) |
---|
335 | |
---|
336 | if(ierr /= 0) then |
---|
337 | call die(myname_,'deallocate(gRow...',ierr) |
---|
338 | endif |
---|
339 | |
---|
340 | end subroutine SparseMatrixToYGlobalSegMap_ |
---|
341 | |
---|
342 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
343 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
344 | !BOP ------------------------------------------------------------------- |
---|
345 | ! |
---|
346 | ! !IROUTINE: CreateSegments_ - Generate segment information. |
---|
347 | ! |
---|
348 | ! !DESCRIPTION: This routine examines an input {\tt INTEGER} list of |
---|
349 | ! numbers {\tt indices} (of length {\tt num\_indices}), determines the |
---|
350 | ! number of segments of consecutive numbers (or runs) {\tt nsegs}. The |
---|
351 | ! starting indices for each run, and their lengths are returned in the |
---|
352 | ! {\tt INTEGER} arrays {\tt starts(:)} and {\tt lengths(:)}, respectively. |
---|
353 | ! |
---|
354 | ! !INTERFACE: |
---|
355 | |
---|
356 | subroutine ComputeSegments_(indices, num_indices, nsegs, starts, lengths) |
---|
357 | |
---|
358 | ! |
---|
359 | ! !USES: |
---|
360 | ! |
---|
361 | use m_stdio, only : stderr |
---|
362 | use m_die, only : die |
---|
363 | |
---|
364 | implicit none |
---|
365 | ! |
---|
366 | ! !INPUT PARAMETERS: |
---|
367 | ! |
---|
368 | |
---|
369 | integer, dimension(:), intent(in) :: indices |
---|
370 | integer, intent(in) :: num_indices |
---|
371 | ! |
---|
372 | ! !OUTPUT PARAMETERS: |
---|
373 | ! |
---|
374 | integer, intent(out) :: nsegs |
---|
375 | integer, dimension(:), pointer :: starts |
---|
376 | integer, dimension(:), pointer :: lengths |
---|
377 | |
---|
378 | |
---|
379 | ! !REVISION HISTORY: |
---|
380 | ! 19Apr01 - J.W. Larson <larson@mcs.anl.gov> - API specification. |
---|
381 | ! 25Apr01 - J.W. Larson <larson@mcs.anl.gov> - Initial code. |
---|
382 | ! 27Apr01 - J.W. Larson <larson@mcs.anl.gov> - Bug fix--error in |
---|
383 | ! computation of segment starts/lengths. |
---|
384 | ! 27Nov01 - E.T. Ong <eong@mcs.anl.gov> - Bug fix--initialize |
---|
385 | ! nsegs=0 in case num_indices=0. |
---|
386 | !EOP ___________________________________________________________________ |
---|
387 | |
---|
388 | character(len=*),parameter :: myname_=myname//'::ComputeSegments_' |
---|
389 | |
---|
390 | integer :: i, ierr |
---|
391 | |
---|
392 | ! First pass: count the segments |
---|
393 | |
---|
394 | nsegs = 0 |
---|
395 | |
---|
396 | do i=1,num_indices |
---|
397 | |
---|
398 | if(i == 1) then ! bootstrap segment counting process |
---|
399 | |
---|
400 | nsegs = 1 |
---|
401 | |
---|
402 | else |
---|
403 | |
---|
404 | if(indices(i) > indices(i-1) + 1) then ! new segment |
---|
405 | nsegs = nsegs + 1 |
---|
406 | endif |
---|
407 | |
---|
408 | endif ! if(i==1) |
---|
409 | |
---|
410 | end do ! do i=1, num_indices |
---|
411 | |
---|
412 | ! Allocate storage space for starts(:) and lengths(:) |
---|
413 | |
---|
414 | allocate(starts(nsegs), lengths(nsegs), stat=ierr) |
---|
415 | |
---|
416 | if(ierr /= 0) then |
---|
417 | call die(myname_,'allocate(starts...',ierr) |
---|
418 | endif |
---|
419 | |
---|
420 | ! Second pass: compute segment start/length info |
---|
421 | |
---|
422 | do i=1,num_indices |
---|
423 | |
---|
424 | select case(i) |
---|
425 | case(1) ! bootstrap segment counting process |
---|
426 | nsegs = 1 |
---|
427 | starts(nsegs) = indices(i) |
---|
428 | ! rml patch |
---|
429 | lengths(nsegs) = 1 |
---|
430 | case default |
---|
431 | |
---|
432 | if(i == num_indices) then ! last point |
---|
433 | if(indices(i) > indices(i-1) + 1) then ! new segment with 1 pt. |
---|
434 | ! first, close the books on the penultimate segment: |
---|
435 | lengths(nsegs) = indices(i-1) - starts(nsegs) + 1 |
---|
436 | nsegs = nsegs + 1 |
---|
437 | starts(nsegs) = indices(i) |
---|
438 | lengths(nsegs) = 1 ! (just one point) |
---|
439 | else |
---|
440 | lengths(nsegs) = indices(i) - starts(nsegs) + 1 |
---|
441 | endif |
---|
442 | else |
---|
443 | if(indices(i) > indices(i-1) + 1) then ! new segment |
---|
444 | lengths(nsegs) = indices(i-1) - starts(nsegs) + 1 |
---|
445 | nsegs = nsegs + 1 |
---|
446 | starts(nsegs) = indices(i) |
---|
447 | endif |
---|
448 | endif |
---|
449 | |
---|
450 | end select ! select case(i) |
---|
451 | |
---|
452 | end do ! do i=1, num_indices |
---|
453 | |
---|
454 | end subroutine ComputeSegments_ |
---|
455 | |
---|
456 | end module m_SparseMatrixToMaps |
---|