Actual source code: iscoloring.c
1: /*$Id: iscoloring.c,v 1.69 2001/04/04 14:43:59 bsmith Exp $*/
3: #include "petscsys.h" /*I "petscsys.h" I*/
4: #include "petscis.h" /*I "petscis.h" I*/
6: /*@C
7: ISColoringDestroy - Destroys a coloring context.
9: Collective on ISColoring
11: Input Parameter:
12: . iscoloring - the coloring context
14: Level: advanced
16: .seealso: ISColoringView(), MatGetColoring()
17: @*/
18: int ISColoringDestroy(ISColoring iscoloring)
19: {
20: int i,ierr;
25: if (iscoloring->is) {
26: for (i=0; i<iscoloring->n; i++) {
27: ISDestroy(iscoloring->is[i]);
28: }
29: PetscFree(iscoloring->is);
30: }
31: if (iscoloring->colors) {
32: PetscFree(iscoloring->colors);
33: }
34: PetscCommDestroy_Private(&iscoloring->comm);
35: PetscFree(iscoloring);
36: return(0);
37: }
39: /*@C
40: ISColoringView - Views a coloring context.
42: Collective on ISColoring
44: Input Parameters:
45: + iscoloring - the coloring context
46: - viewer - the viewer
48: Level: advanced
50: .seealso: ISColoringDestroy(), ISColoringGetIS(), MatGetColoring()
51: @*/
52: int ISColoringView(ISColoring iscoloring,PetscViewer viewer)
53: {
54: int i,ierr;
55: PetscTruth isascii;
56: IS *is;
60: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(iscoloring->comm);
63: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
64: if (isascii) {
65: MPI_Comm comm;
66: int rank;
67: PetscObjectGetComm((PetscObject)viewer,&comm);
68: MPI_Comm_rank(comm,&rank);
69: PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Number of colors %dn",rank,iscoloring->n);
70: PetscViewerFlush(viewer);
71: } else {
72: SETERRQ1(1,"Viewer type %s not supported for ISColoring",((PetscObject)viewer)->type_name);
73: }
75: ISColoringGetIS(iscoloring,PETSC_IGNORE,&is);
76: for (i=0; i<iscoloring->n; i++) {
77: ISView(iscoloring->is[i],viewer);
78: }
79: ISColoringRestoreIS(iscoloring,&is);
80: return(0);
81: }
83: /*@C
84: ISColoringGetIS - Extracts index sets from the coloring context
86: Collective on ISColoring
88: Input Parameter:
89: . iscoloring - the coloring context
91: Output Parameters:
92: + nn - number of index sets in the coloring context
93: - is - array of index sets
95: Level: advanced
97: .seealso: ISColoringRestoreIS(), ISColoringView()
98: @*/
99: int ISColoringGetIS(ISColoring iscoloring,int *nn,IS *isis[])
100: {
106: if (nn) *nn = iscoloring->n;
107: if (isis) {
108: if (!iscoloring->is) {
109: int *mcolors,**ii,nc = iscoloring->n,i,base, n = iscoloring->N;
110: int *colors = iscoloring->colors;
111: IS *is;
113:
114: /* generate the lists of nodes for each color */
115: PetscMalloc((nc+1)*sizeof(int),&mcolors);
116: PetscMemzero(mcolors,nc*sizeof(int));
117: for (i=0; i<n; i++) {
118: mcolors[colors[i]]++;
119: }
121: PetscMalloc((nc+1)*sizeof(int*),&ii);
122: PetscMalloc((n+1)*sizeof(int),&ii[0]);
123: for (i=1; i<nc; i++) {
124: ii[i] = ii[i-1] + mcolors[i-1];
125: }
126:
127: MPI_Scan(&iscoloring->n,&base,1,MPI_INT,MPI_SUM,iscoloring->comm);
128: base -= iscoloring->n;
129: PetscMemzero(mcolors,nc*sizeof(int));
130: for (i=0; i<n; i++) {
131: ii[colors[i]][mcolors[colors[i]]++] = i + base;
132: }
133: PetscMalloc((nc+1)*sizeof(IS),&is);
134: for (i=0; i<nc; i++) {
135: ISCreateGeneral(iscoloring->comm,mcolors[i],ii[i],is+i);
136: }
138: iscoloring->is = is;
139: PetscFree(ii[0]);
140: PetscFree(ii);
141: PetscFree(mcolors);
142: }
143: *isis = iscoloring->is;
144: }
146: return(0);
147: }
149: /*@C
150: ISColoringGetIS - Restores the index sets extracted from the coloring context
152: Collective on ISColoring
154: Input Parameter:
155: + iscoloring - the coloring context
156: - is - array of index sets
158: Level: advanced
160: .seealso: ISColoringGetIS(), ISColoringView()
161: @*/
162: int ISColoringRestoreIS(ISColoring iscoloring,IS *is[])
163: {
166:
167: /* currently nothing is done here */
169: return(0);
170: }
173: /*@C
174: ISColoringCreate - Generates an ISColoring context from lists (provided
175: by each processor) of colors for each node.
177: Collective on MPI_Comm
179: Input Parameters:
180: + comm - communicator for the processors creating the coloring
181: . n - number of nodes on this processor
182: - colors - array containing the colors for this processor, color
183: numbers begin at 0. In C/C++ this array must have been obtained with PetscMalloc()
184: and should NOT be freed (The ISColoringDestroy() will free it).
186: Output Parameter:
187: . iscoloring - the resulting coloring data structure
189: Options Database Key:
190: . -is_coloring_view - Activates ISColoringView()
192: Level: advanced
194: .seealso: MatColoringCreate(), ISColoringView(), ISColoringDestroy()
195: @*/
196: int ISColoringCreate(MPI_Comm comm,int n,const int colors[],ISColoring *iscoloring)
197: {
198: int ierr,size,rank,base,top,tag,nc,ncwork,i;
199: PetscTruth flg;
200: MPI_Status status;
203: PetscNew(struct _p_ISColoring,iscoloring);
204: PetscCommDuplicate_Private(comm,&(*iscoloring)->comm,&tag);
205: comm = (*iscoloring)->comm;
207: /* compute the number of the first node on my processor */
208: MPI_Comm_size(comm,&size);
210: /* should use MPI_Scan() */
211: MPI_Comm_rank(comm,&rank);
212: if (!rank) {
213: base = 0;
214: top = n;
215: } else {
216: MPI_Recv(&base,1,MPI_INT,rank-1,tag,comm,&status);
217: top = base+n;
218: }
219: if (rank < size-1) {
220: MPI_Send(&top,1,MPI_INT,rank+1,tag,comm);
221: }
223: /* compute the total number of colors */
224: ncwork = 0;
225: for (i=0; i<n; i++) {
226: if (ncwork < colors[i]) ncwork = colors[i];
227: }
228: ncwork++;
229: MPI_Allreduce(&ncwork,&nc,1,MPI_INT,MPI_MAX,comm);
230: (*iscoloring)->n = nc;
231: (*iscoloring)->is = 0;
232: (*iscoloring)->colors = (int *)colors;
233: (*iscoloring)->N = n;
235: PetscOptionsHasName(PETSC_NULL,"-is_coloring_view",&flg);
236: if (flg) {
237: ISColoringView(*iscoloring,PETSC_VIEWER_STDOUT_((*iscoloring)->comm));
238: }
240: return(0);
241: }
243: /*@C
244: ISPartitioningToNumbering - Takes an ISPartitioning and on each processor
245: generates an IS that contains a new global node number for each index based
246: on the partitioing.
248: Collective on IS
250: Input Parameters
251: . partitioning - a partitioning as generated by MatPartitioningApply()
253: Output Parameter:
254: . is - on each processor the index set that defines the global numbers
255: (in the new numbering) for all the nodes currently (before the partitioning)
256: on that processor
258: Level: advanced
260: .seealso: MatPartitioningCreate(), AOCreateBasic(), ISPartioningCount()
262: @*/
263: int ISPartitioningToNumbering(IS part,IS *is)
264: {
265: MPI_Comm comm;
266: int i,ierr,size,*indices,np,n,*starts,*sums,*lsizes,*newi;
269: PetscObjectGetComm((PetscObject)part,&comm);
270: MPI_Comm_size(comm,&size);
272: /* count the number of partitions, make sure <= size */
273: ISGetLocalSize(part,&n);
274: ISGetIndices(part,&indices);
275: np = 0;
276: for (i=0; i<n; i++) {
277: np = PetscMax(np,indices[i]);
278: }
279: if (np >= size) {
280: SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Number of partitions %d larger than number of processors %d",np,size);
281: }
283: /*
284: lsizes - number of elements of each partition on this particular processor
285: sums - total number of "previous" nodes for any particular partition
286: starts - global number of first element in each partition on this processor
287: */
288: ierr = PetscMalloc(3*size*sizeof(int),&lsizes);
289: starts = lsizes + size;
290: sums = starts + size;
291: ierr = PetscMemzero(lsizes,size*sizeof(int));
292: for (i=0; i<n; i++) {
293: lsizes[indices[i]]++;
294: }
295: MPI_Allreduce(lsizes,sums,size,MPI_INT,MPI_SUM,comm);
296: MPI_Scan(lsizes,starts,size,MPI_INT,MPI_SUM,comm);
297: for (i=0; i<size; i++) {
298: starts[i] -= lsizes[i];
299: }
300: for (i=1; i<size; i++) {
301: sums[i] += sums[i-1];
302: starts[i] += sums[i-1];
303: }
305: /*
306: For each local index give it the new global number
307: */
308: PetscMalloc((n+1)*sizeof(int),&newi);
309: for (i=0; i<n; i++) {
310: newi[i] = starts[indices[i]]++;
311: }
312: PetscFree(lsizes);
314: ISRestoreIndices(part,&indices);
315: ISCreateGeneral(comm,n,newi,is);
316: PetscFree(newi);
317: ISSetPermutation(*is);
318: return(0);
319: }
321: /*@C
322: ISPartitioningCount - Takes a ISPartitioning and determines the number of
323: resulting elements on each processor
325: Collective on IS
327: Input Parameters:
328: . partitioning - a partitioning as generated by MatPartitioningApply()
330: Output Parameter:
331: . count - array of length size of communicator associated with IS, contains
332: the number of elements assigned to each processor
334: Level: advanced
336: .seealso: MatPartitioningCreate(), AOCreateBasic(), ISPartitioningToNumbering()
338: @*/
339: int ISPartitioningCount(IS part,int count[])
340: {
341: MPI_Comm comm;
342: int i,ierr,size,*indices,np,n,*lsizes;
345: PetscObjectGetComm((PetscObject)part,&comm);
346: MPI_Comm_size(comm,&size);
348: /* count the number of partitions,make sure <= size */
349: ISGetLocalSize(part,&n);
350: ISGetIndices(part,&indices);
351: np = 0;
352: for (i=0; i<n; i++) {
353: np = PetscMax(np,indices[i]);
354: }
355: if (np >= size) {
356: SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Number of partitions %d larger than number of processors %d",np,size);
357: }
359: /*
360: lsizes - number of elements of each partition on this particular processor
361: sums - total number of "previous" nodes for any particular partition
362: starts - global number of first element in each partition on this processor
363: */
364: PetscMalloc(size*sizeof(int),&lsizes);
365: ierr = PetscMemzero(lsizes,size*sizeof(int));
366: for (i=0; i<n; i++) {
367: lsizes[indices[i]]++;
368: }
369: ISRestoreIndices(part,&indices);
370: MPI_Allreduce(lsizes,count,size,MPI_INT,MPI_SUM,comm);
371: PetscFree(lsizes);
373: return(0);
374: }
376: /*@C
377: ISAllGather - Given an index set (IS) on each processor, generates a large
378: index set (same on each processor) by concatenating together each
379: processors index set.
381: Collective on IS
383: Input Parameter:
384: . is - the distributed index set
386: Output Parameter:
387: . isout - the concatenated index set (same on all processors)
389: Notes:
390: ISAllGather() is clearly not scalable for large index sets.
392: The IS created on each processor must be created with a common
393: communicator (e.g., PETSC_COMM_WORLD). If the index sets were created
394: with PETSC_COMM_SELF, this routine will not work as expected, since
395: each process will generate its own new IS that consists only of
396: itself.
398: Level: intermediate
400: Concepts: gather^index sets
401: Concepts: index sets^gathering to all processors
402: Concepts: IS^gathering to all processors
404: .seealso: ISCreateGeneral(), ISCreateStride(), ISCreateBlock()
405: @*/
406: int ISAllGather(IS is,IS *isout)
407: {
408: int *indices,*sizes,size,*offsets,n,*lindices,i,N,ierr;
409: MPI_Comm comm;
414: PetscObjectGetComm((PetscObject)is,&comm);
415: MPI_Comm_size(comm,&size);
416: PetscMalloc(2*size*sizeof(int),&sizes);
417: offsets = sizes + size;
418:
419: ISGetLocalSize(is,&n);
420: MPI_Allgather(&n,1,MPI_INT,sizes,1,MPI_INT,comm);
421: offsets[0] = 0;
422: for (i=1;i<size; i++) offsets[i] = offsets[i-1] + sizes[i-1];
423: N = offsets[size-1] + sizes[size-1];
425: PetscMalloc((N+1)*sizeof(int),&indices);
426: ISGetIndices(is,&lindices);
427: MPI_Allgatherv(lindices,n,MPI_INT,indices,sizes,offsets,MPI_INT,comm);
428: ISRestoreIndices(is,&lindices);
430: ISCreateGeneral(PETSC_COMM_SELF,N,indices,isout);
431: PetscFree(indices);
433: PetscFree(sizes);
434: return(0);
435: }