Actual source code: color.c
1: #define PETSCMAT_DLL
3: /*
4: Routines that call the kernel minpack coloring subroutines
5: */
7: #include private/matimpl.h
8: #include ../src/mat/color/color.h
10: /*
11: MatFDColoringDegreeSequence_Minpack - Calls the MINPACK routine seqr() that
12: computes the degree sequence required by MINPACK coloring routines.
13: */
16: PetscErrorCode MatFDColoringDegreeSequence_Minpack(PetscInt m,PetscInt *cja, PetscInt *cia, PetscInt *rja, PetscInt *ria, PetscInt **seq)
17: {
18: PetscInt *work;
22: PetscMalloc(m*sizeof(PetscInt),&work);
23: PetscMalloc(m*sizeof(PetscInt),seq);
25: MINPACKdegr(&m,cja,cia,rja,ria,*seq,work);
27: PetscFree(work);
28: return(0);
29: }
31: /*
32: MatFDColoringMinimumNumberofColors_Private - For a given sparse
33: matrix computes the minimum number of colors needed.
35: */
38: PetscErrorCode MatFDColoringMinimumNumberofColors_Private(PetscInt m,PetscInt *ia,PetscInt *minc)
39: {
40: PetscInt i,c = 0;
43: for (i=0; i<m; i++) {
44: c = PetscMax(c,ia[i+1]-ia[i]);
45: }
46: *minc = c;
47: return(0);
48: }
51: /* ----------------------------------------------------------------------------*/
52: /*
53: MatGetColoring_SL_Minpack - Uses the smallest-last (SL) coloring of minpack
54: */
57: PetscErrorCode MatGetColoring_SL_Minpack(Mat mat,MatColoringType name,ISColoring *iscoloring)
58: {
59: PetscErrorCode ierr;
60: PetscInt *list,*work,clique,*ria,*rja,*cia,*cja,*seq,*coloring,n;
61: PetscInt ncolors,i;
62: PetscTruth done;
63: Mat mat_seq = mat;
64: PetscMPIInt size;
65: MPI_Comm comm;
66: ISColoring iscoloring_seq;
67: PetscInt bs = 1,rstart,rend,N_loc,nc;
68: ISColoringValue *colors_loc;
69: PetscTruth flg1,flg2;
72: /* this is ugly way to get blocksize but cannot call MatGetBlockSize() because AIJ can have bs > 1 */
73: PetscTypeCompare((PetscObject)mat,MATSEQBAIJ,&flg1);
74: PetscTypeCompare((PetscObject)mat,MATMPIBAIJ,&flg2);
75: if (flg1 || flg2) {
76: MatGetBlockSize(mat,&bs);
77: }
79: PetscObjectGetComm((PetscObject)mat,&comm);
80: MPI_Comm_size(comm,&size);
81: if (size > 1){
82: /* create a sequential iscoloring on all processors */
83: MatGetSeqNonzeroStructure(mat,&mat_seq);
84: }
86: MatGetRowIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,&n,&ria,&rja,&done);
87: MatGetColumnIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,&n,&cia,&cja,&done);
88: if (!done) SETERRQ(PETSC_ERR_SUP,"Ordering requires IJ");
90: MatFDColoringDegreeSequence_Minpack(n,cja,cia,rja,ria,&seq);
92: PetscMalloc2(n,PetscInt,&list,4*n,PetscInt,&work);
94: MINPACKslo(&n,cja,cia,rja,ria,seq,list,&clique,work,work+n,work+2*n,work+3*n);
96: PetscMalloc(n*sizeof(PetscInt),&coloring);
97: MINPACKseq(&n,cja,cia,rja,ria,list,coloring,&ncolors,work);
99: PetscFree2(list,work);
100: PetscFree(seq);
101: MatRestoreRowIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,&n,&ria,&rja,&done);
102: MatRestoreColumnIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,&n,&cia,&cja,&done);
104: /* shift coloring numbers to start at zero and shorten */
105: if (ncolors > IS_COLORING_MAX-1) SETERRQ(PETSC_ERR_SUP,"Maximum color size exceeded");
106: {
107: ISColoringValue *s = (ISColoringValue*) coloring;
108: for (i=0; i<n; i++) {
109: s[i] = (ISColoringValue) (coloring[i]-1);
110: }
111: MatColoringPatch(mat_seq,ncolors,n,s,iscoloring);
112: }
114: if (size > 1) {
115: MatDestroySeqNonzeroStructure(&mat_seq);
117: /* convert iscoloring_seq to a parallel iscoloring */
118: iscoloring_seq = *iscoloring;
119: rstart = mat->rmap->rstart/bs;
120: rend = mat->rmap->rend/bs;
121: N_loc = rend - rstart; /* number of local nodes */
123: /* get local colors for each local node */
124: PetscMalloc((N_loc+1)*sizeof(ISColoringValue),&colors_loc);
125: for (i=rstart; i<rend; i++){
126: colors_loc[i-rstart] = iscoloring_seq->colors[i];
127: }
128: /* create a parallel iscoloring */
129: nc=iscoloring_seq->n;
130: ISColoringCreate(comm,nc,N_loc,colors_loc,iscoloring);
131: ISColoringDestroy(iscoloring_seq);
132: }
133: return(0);
134: }
138: /* ----------------------------------------------------------------------------*/
139: /*
140: MatGetColoring_LF_Minpack -
141: */
144: PetscErrorCode MatGetColoring_LF_Minpack(Mat mat,MatColoringType name,ISColoring *iscoloring)
145: {
146: PetscErrorCode ierr;
147: PetscInt *list,*work,*ria,*rja,*cia,*cja,*seq,*coloring,n;
148: PetscInt n1, none,ncolors,i;
149: PetscTruth done;
150: Mat mat_seq = mat;
151: PetscMPIInt size;
152: MPI_Comm comm;
153: ISColoring iscoloring_seq;
154: PetscInt bs = 1,rstart,rend,N_loc,nc;
155: ISColoringValue *colors_loc;
156: PetscTruth flg1,flg2;
159: /* this is ugly way to get blocksize but cannot call MatGetBlockSize() because AIJ can have bs > 1 */
160: PetscTypeCompare((PetscObject)mat,MATSEQBAIJ,&flg1);
161: PetscTypeCompare((PetscObject)mat,MATMPIBAIJ,&flg2);
162: if (flg1 || flg2) {
163: MatGetBlockSize(mat,&bs);
164: }
166: PetscObjectGetComm((PetscObject)mat,&comm);
167: MPI_Comm_size(comm,&size);
168: if (size > 1){
169: /* create a sequential iscoloring on all processors */
170: MatGetSeqNonzeroStructure(mat,&mat_seq);
171: }
173: MatGetRowIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,&n,&ria,&rja,&done);
174: MatGetColumnIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,&n,&cia,&cja,&done);
175: if (!done) SETERRQ(PETSC_ERR_SUP,"Ordering requires IJ");
177: MatFDColoringDegreeSequence_Minpack(n,cja,cia,rja,ria,&seq);
179: PetscMalloc2(n,PetscInt,&list,4*n,PetscInt,&work);
181: n1 = n - 1;
182: none = -1;
183: MINPACKnumsrt(&n,&n1,seq,&none,list,work+2*n,work+n);
184: PetscMalloc(n*sizeof(PetscInt),&coloring);
185: MINPACKseq(&n,cja,cia,rja,ria,list,coloring,&ncolors,work);
187: PetscFree2(list,work);
188: PetscFree(seq);
190: MatRestoreRowIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,&n,&ria,&rja,&done);
191: MatRestoreColumnIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,&n,&cia,&cja,&done);
193: /* shift coloring numbers to start at zero and shorten */
194: if (ncolors > IS_COLORING_MAX-1) SETERRQ(PETSC_ERR_SUP,"Maximum color size exceeded");
195: {
196: ISColoringValue *s = (ISColoringValue*) coloring;
197: for (i=0; i<n; i++) {
198: s[i] = (ISColoringValue) (coloring[i]-1);
199: }
200: MatColoringPatch(mat_seq,ncolors,n,s,iscoloring);
201: }
203: if (size > 1) {
204: MatDestroySeqNonzeroStructure(&mat_seq);
206: /* convert iscoloring_seq to a parallel iscoloring */
207: iscoloring_seq = *iscoloring;
208: rstart = mat->rmap->rstart/bs;
209: rend = mat->rmap->rend/bs;
210: N_loc = rend - rstart; /* number of local nodes */
212: /* get local colors for each local node */
213: PetscMalloc((N_loc+1)*sizeof(ISColoringValue),&colors_loc);
214: for (i=rstart; i<rend; i++){
215: colors_loc[i-rstart] = iscoloring_seq->colors[i];
216: }
217: /* create a parallel iscoloring */
218: nc=iscoloring_seq->n;
219: ISColoringCreate(comm,nc,N_loc,colors_loc,iscoloring);
220: ISColoringDestroy(iscoloring_seq);
221: }
222: return(0);
223: }
227: /* ----------------------------------------------------------------------------*/
228: /*
229: MatGetColoring_ID_Minpack -
230: */
233: PetscErrorCode MatGetColoring_ID_Minpack(Mat mat,MatColoringType name,ISColoring *iscoloring)
234: {
235: PetscErrorCode ierr;
236: PetscInt *list,*work,clique,*ria,*rja,*cia,*cja,*seq,*coloring,n;
237: PetscInt ncolors,i;
238: PetscTruth done;
239: Mat mat_seq = mat;
240: PetscMPIInt size;
241: MPI_Comm comm;
242: ISColoring iscoloring_seq;
243: PetscInt bs = 1,rstart,rend,N_loc,nc;
244: ISColoringValue *colors_loc;
245: PetscTruth flg1,flg2;
248: /* this is ugly way to get blocksize but cannot call MatGetBlockSize() because AIJ can have bs > 1 */
249: PetscTypeCompare((PetscObject)mat,MATSEQBAIJ,&flg1);
250: PetscTypeCompare((PetscObject)mat,MATMPIBAIJ,&flg2);
251: if (flg1 || flg2) {
252: MatGetBlockSize(mat,&bs);
253: }
255: PetscObjectGetComm((PetscObject)mat,&comm);
256: MPI_Comm_size(comm,&size);
257: if (size > 1){
258: /* create a sequential iscoloring on all processors */
259: MatGetSeqNonzeroStructure(mat,&mat_seq);
260: }
262: MatGetRowIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,&n,&ria,&rja,&done);
263: MatGetColumnIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,&n,&cia,&cja,&done);
264: if (!done) SETERRQ(PETSC_ERR_SUP,"Ordering requires IJ");
266: MatFDColoringDegreeSequence_Minpack(n,cja,cia,rja,ria,&seq);
268: PetscMalloc2(n,PetscInt,&list,4*n,PetscInt,&work);
270: MINPACKido(&n,&n,cja,cia,rja,ria,seq,list,&clique,work,work+n,work+2*n,work+3*n);
272: PetscMalloc(n*sizeof(PetscInt),&coloring);
273: MINPACKseq(&n,cja,cia,rja,ria,list,coloring,&ncolors,work);
275: PetscFree2(list,work);
276: PetscFree(seq);
278: MatRestoreRowIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,&n,&ria,&rja,&done);
279: MatRestoreColumnIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,&n,&cia,&cja,&done);
281: /* shift coloring numbers to start at zero and shorten */
282: if (ncolors > IS_COLORING_MAX-1) SETERRQ(PETSC_ERR_SUP,"Maximum color size exceeded");
283: {
284: ISColoringValue *s = (ISColoringValue*) coloring;
285: for (i=0; i<n; i++) {
286: s[i] = (ISColoringValue) (coloring[i]-1);
287: }
288: MatColoringPatch(mat_seq,ncolors,n,s,iscoloring);
289: }
291: if (size > 1) {
292: MatDestroySeqNonzeroStructure(&mat_seq);
294: /* convert iscoloring_seq to a parallel iscoloring */
295: iscoloring_seq = *iscoloring;
296: rstart = mat->rmap->rstart/bs;
297: rend = mat->rmap->rend/bs;
298: N_loc = rend - rstart; /* number of local nodes */
300: /* get local colors for each local node */
301: PetscMalloc((N_loc+1)*sizeof(ISColoringValue),&colors_loc);
302: for (i=rstart; i<rend; i++){
303: colors_loc[i-rstart] = iscoloring_seq->colors[i];
304: }
305: /* create a parallel iscoloring */
306: nc=iscoloring_seq->n;
307: ISColoringCreate(comm,nc,N_loc,colors_loc,iscoloring);
308: ISColoringDestroy(iscoloring_seq);
309: }
310: return(0);
311: }
315: /*
316: Simplest coloring, each column of the matrix gets its own unique color.
317: */
320: PetscErrorCode MatGetColoring_Natural(Mat mat,MatColoringType color, ISColoring *iscoloring)
321: {
322: PetscErrorCode ierr;
323: PetscInt start,end,i,bs = 1,n;
324: ISColoringValue *colors;
325: MPI_Comm comm;
326: PetscTruth flg1,flg2;
327: Mat mat_seq = mat;
328: PetscMPIInt size;
329: ISColoring iscoloring_seq;
330: ISColoringValue *colors_loc;
331: PetscInt rstart,rend,N_loc,nc;
334: /* this is ugly way to get blocksize but cannot call MatGetBlockSize() because AIJ can have bs > 1 */
335: PetscTypeCompare((PetscObject)mat,MATSEQBAIJ,&flg1);
336: PetscTypeCompare((PetscObject)mat,MATMPIBAIJ,&flg2);
337: if (flg1 || flg2) {
338: MatGetBlockSize(mat,&bs);
339: }
341: PetscObjectGetComm((PetscObject)mat,&comm);
342: MPI_Comm_size(comm,&size);
343: if (size > 1){
344: /* create a sequential iscoloring on all processors */
345: MatGetSeqNonzeroStructure(mat,&mat_seq);
346: }
348: MatGetSize(mat_seq,PETSC_NULL,&n);
349: MatGetOwnershipRange(mat_seq,&start,&end);
350: n = n/bs;
351: if (n > IS_COLORING_MAX-1) SETERRQ(PETSC_ERR_SUP,"Maximum color size exceeded");
353: start = start/bs;
354: end = end/bs;
355: PetscMalloc((end-start+1)*sizeof(PetscInt),&colors);
356: for (i=start; i<end; i++) {
357: colors[i-start] = (ISColoringValue)i;
358: }
359: ISColoringCreate(comm,n,end-start,colors,iscoloring);
361: if (size > 1) {
362: MatDestroySeqNonzeroStructure(&mat_seq);
364: /* convert iscoloring_seq to a parallel iscoloring */
365: iscoloring_seq = *iscoloring;
366: rstart = mat->rmap->rstart/bs;
367: rend = mat->rmap->rend/bs;
368: N_loc = rend - rstart; /* number of local nodes */
370: /* get local colors for each local node */
371: PetscMalloc((N_loc+1)*sizeof(ISColoringValue),&colors_loc);
372: for (i=rstart; i<rend; i++){
373: colors_loc[i-rstart] = iscoloring_seq->colors[i];
374: }
375: /* create a parallel iscoloring */
376: nc=iscoloring_seq->n;
377: ISColoringCreate(comm,nc,N_loc,colors_loc,iscoloring);
378: ISColoringDestroy(iscoloring_seq);
379: }
380: return(0);
381: }
383:
384: /* ===========================================================================================*/
386: PetscFList MatColoringList = 0;
387: PetscTruth MatColoringRegisterAllCalled = PETSC_FALSE;
391: PetscErrorCode MatColoringRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(Mat,MatColoringType,ISColoring*))
392: {
394: char fullname[PETSC_MAX_PATH_LEN];
397: PetscFListConcat(path,name,fullname);
398: PetscFListAdd(&MatColoringList,sname,fullname,(void (*)(void))function);
399: return(0);
400: }
404: /*@C
405: MatColoringRegisterDestroy - Frees the list of coloringing routines.
407: Not Collective
409: Level: developer
411: .keywords: matrix, register, destroy
413: .seealso: MatColoringRegisterDynamic(), MatColoringRegisterAll()
414: @*/
415: PetscErrorCode MatColoringRegisterDestroy(void)
416: {
420: PetscFListDestroy(&MatColoringList);
421: MatColoringRegisterAllCalled = PETSC_FALSE;
422: return(0);
423: }
427: /*@C
428: MatGetColoring - Gets a coloring for a matrix, from its sparsity structure,
429: to reduce the number of function evaluations needed to compute a sparse Jacobian via differencing.
431: Collective on Mat
433: Input Parameters:
434: . mat - the matrix
435: . type - type of coloring, one of the following:
436: $ MATCOLORING_NATURAL - natural (one color for each column, very slow)
437: $ MATCOLORING_SL - smallest-last
438: $ MATCOLORING_LF - largest-first
439: $ MATCOLORING_ID - incidence-degree
441: Output Parameters:
442: . iscoloring - the coloring
444: Options Database Keys:
445: To specify the coloring through the options database, use one of
446: the following
447: $ -mat_coloring_type natural, -mat_coloring_type sl, -mat_coloring_type lf,
448: $ -mat_coloring_type id
449: To see the coloring use
450: $ -mat_coloring_view
452: Level: intermediate
454: Notes:
455: These compute the graph coloring of the graph of A^{T}A. The coloring used
456: for efficient (parallel or thread based) triangular solves etc is NOT yet
457: available.
459: The user can define additional colorings; see MatColoringRegisterDynamic().
461: For parallel matrices currently converts to sequential matrix and uses the sequential coloring
462: on that.
464: The colorings SL, LF, and ID are obtained via the Minpack software that was
465: converted to C using f2c.
467: For BAIJ matrices this colors the blocks. The true number of colors would be block size times the number of colors
468: returned here.
470: References:
471: $ Thomas F. Coleman and Jorge J. More, Estimation of Sparse {J}acobian Matrices and Graph Coloring Problems,
472: $ SIAM Journal on Numerical Analysis, 1983, pages 187-209, volume 20
473: $ Jorge J. Mor\'{e} and Danny C. Sorenson and Burton S. Garbow and Kenneth E. Hillstrom, The {MINPACK} Project,
474: $ Sources and Development of Mathematical Software, Wayne R. Cowell editor, 1984, pages 88-111
476: .keywords: matrix, get, coloring
478: .seealso: MatGetColoringTypeFromOptions(), MatColoringRegisterDynamic(), MatFDColoringCreate(),
479: SNESDefaultComputeJacobianColor()
480: @*/
481: PetscErrorCode MatGetColoring(Mat mat,const MatColoringType type,ISColoring *iscoloring)
482: {
483: PetscTruth flag;
484: PetscErrorCode ierr,(*r)(Mat,const MatColoringType,ISColoring *);
485: char tname[PETSC_MAX_PATH_LEN];
486: MPI_Comm comm;
491: if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
492: if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
494: /* look for type on command line */
495: if (!MatColoringRegisterAllCalled) {MatColoringRegisterAll(PETSC_NULL);}
496: PetscOptionsGetString(((PetscObject)mat)->prefix,"-mat_coloring_type",tname,256,&flag);
497: if (flag) { type = tname; }
499: PetscObjectGetComm((PetscObject)mat,&comm);
500: PetscFListFind(MatColoringList,comm, type,(void (**)(void)) &r);
501: if (!r) {SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Unknown or unregistered type: %s",type);}
503: PetscLogEventBegin(MAT_GetColoring,mat,0,0,0);
504: (*r)(mat,type,iscoloring);
505: PetscLogEventEnd(MAT_GetColoring,mat,0,0,0);
507: PetscInfo1(mat,"Number of colors %d\n",(*iscoloring)->n);
508: flag = PETSC_FALSE;
509: PetscOptionsGetTruth(PETSC_NULL,"-mat_coloring_view",&flag,PETSC_NULL);
510: if (flag) {
511: PetscViewer viewer;
512: PetscViewerASCIIGetStdout((*iscoloring)->comm,&viewer);
513: ISColoringView(*iscoloring,viewer);
514: }
515: return(0);
516: }
517: