Actual source code: color.c
petsc-3.5.4 2015-05-23
2: /*
3: Routines that call the kernel minpack coloring subroutines
4: */
6: #include <petsc-private/matimpl.h>
7: #include <../src/mat/color/impls/minpack/color.h>
9: /*
10: MatFDColoringDegreeSequence_Minpack - Calls the MINPACK routine seqr() that
11: computes the degree sequence required by MINPACK coloring routines.
12: */
15: PetscErrorCode MatFDColoringDegreeSequence_Minpack(PetscInt m,const PetscInt *cja,const PetscInt *cia,const PetscInt *rja,const PetscInt *ria,PetscInt **seq)
16: {
17: PetscInt *work;
21: PetscMalloc1(m,&work);
22: PetscMalloc1(m,seq);
24: MINPACKdegr(&m,cja,cia,rja,ria,*seq,work);
26: PetscFree(work);
27: return(0);
28: }
30: /*
31: MatFDColoringMinimumNumberofColors_Private - For a given sparse
32: matrix computes the minimum number of colors needed.
34: */
37: PetscErrorCode MatFDColoringMinimumNumberofColors_Private(PetscInt m,PetscInt *ia,PetscInt *minc)
38: {
39: PetscInt i,c = 0;
42: for (i=0; i<m; i++) c = PetscMax(c,ia[i+1]-ia[i]);
43: *minc = c;
44: return(0);
45: }
49: PETSC_EXTERN PetscErrorCode MatColoringApply_SL(MatColoring mc,ISColoring *iscoloring)
50: {
51: PetscErrorCode ierr;
52: PetscInt *list,*work,clique,*seq,*coloring,n;
53: const PetscInt *ria,*rja,*cia,*cja;
54: PetscInt ncolors,i;
55: PetscBool done;
56: Mat mat = mc->mat;
57: Mat mat_seq = mat;
58: PetscMPIInt size;
59: MPI_Comm comm;
60: ISColoring iscoloring_seq;
61: PetscInt bs = 1,rstart,rend,N_loc,nc;
62: ISColoringValue *colors_loc;
63: PetscBool flg1,flg2;
66: if (mc->dist != 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"SL may only do distance 2 coloring");
67: /* this is ugly way to get blocksize but cannot call MatGetBlockSize() because AIJ can have bs > 1 */
68: PetscObjectTypeCompare((PetscObject)mat,MATSEQBAIJ,&flg1);
69: PetscObjectTypeCompare((PetscObject)mat,MATMPIBAIJ,&flg2);
70: if (flg1 || flg2) {
71: MatGetBlockSize(mat,&bs);
72: }
74: PetscObjectGetComm((PetscObject)mat,&comm);
75: MPI_Comm_size(comm,&size);
76: if (size > 1) {
77: /* create a sequential iscoloring on all processors */
78: MatGetSeqNonzeroStructure(mat,&mat_seq);
79: }
81: MatGetRowIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,&n,&ria,&rja,&done);
82: MatGetColumnIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,&n,&cia,&cja,&done);
83: if (!done) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Ordering requires IJ");
85: MatFDColoringDegreeSequence_Minpack(n,cja,cia,rja,ria,&seq);
87: PetscMalloc2(n,&list,4*n,&work);
89: MINPACKslo(&n,cja,cia,rja,ria,seq,list,&clique,work,work+n,work+2*n,work+3*n);
91: PetscMalloc1(n,&coloring);
92: MINPACKseq(&n,cja,cia,rja,ria,list,coloring,&ncolors,work);
94: PetscFree2(list,work);
95: PetscFree(seq);
96: MatRestoreRowIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,NULL,&ria,&rja,&done);
97: MatRestoreColumnIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,NULL,&cia,&cja,&done);
99: /* shift coloring numbers to start at zero and shorten */
100: if (ncolors > IS_COLORING_MAX-1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Maximum color size exceeded");
101: {
102: ISColoringValue *s = (ISColoringValue*) coloring;
103: for (i=0; i<n; i++) {
104: s[i] = (ISColoringValue) (coloring[i]-1);
105: }
106: MatColoringPatch(mat_seq,ncolors,n,s,iscoloring);
107: }
109: if (size > 1) {
110: MatDestroySeqNonzeroStructure(&mat_seq);
112: /* convert iscoloring_seq to a parallel iscoloring */
113: iscoloring_seq = *iscoloring;
114: rstart = mat->rmap->rstart/bs;
115: rend = mat->rmap->rend/bs;
116: N_loc = rend - rstart; /* number of local nodes */
118: /* get local colors for each local node */
119: PetscMalloc1((N_loc+1),&colors_loc);
120: for (i=rstart; i<rend; i++) {
121: colors_loc[i-rstart] = iscoloring_seq->colors[i];
122: }
123: /* create a parallel iscoloring */
124: nc = iscoloring_seq->n;
125: ISColoringCreate(comm,nc,N_loc,colors_loc,iscoloring);
126: ISColoringDestroy(&iscoloring_seq);
127: }
128: return(0);
129: }
133: PETSC_EXTERN PetscErrorCode MatColoringCreate_SL(MatColoring mc)
134: {
136: mc->dist = 2;
137: mc->data = NULL;
138: mc->ops->apply = MatColoringApply_SL;
139: mc->ops->view = NULL;
140: mc->ops->destroy = NULL;
141: mc->ops->setfromoptions = NULL;
142: return(0);
143: }
147: PETSC_EXTERN PetscErrorCode MatColoringApply_LF(MatColoring mc,ISColoring *iscoloring)
148: {
149: PetscErrorCode ierr;
150: PetscInt *list,*work,*seq,*coloring,n;
151: const PetscInt *ria,*rja,*cia,*cja;
152: PetscInt n1, none,ncolors,i;
153: PetscBool done;
154: Mat mat = mc->mat;
155: Mat mat_seq = mat;
156: PetscMPIInt size;
157: MPI_Comm comm;
158: ISColoring iscoloring_seq;
159: PetscInt bs = 1,rstart,rend,N_loc,nc;
160: ISColoringValue *colors_loc;
161: PetscBool flg1,flg2;
164: if (mc->dist != 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"LF may only do distance 2 coloring");
165: /* this is ugly way to get blocksize but cannot call MatGetBlockSize() because AIJ can have bs > 1 */
166: PetscObjectTypeCompare((PetscObject)mat,MATSEQBAIJ,&flg1);
167: PetscObjectTypeCompare((PetscObject)mat,MATMPIBAIJ,&flg2);
168: if (flg1 || flg2) {
169: MatGetBlockSize(mat,&bs);
170: }
172: PetscObjectGetComm((PetscObject)mat,&comm);
173: MPI_Comm_size(comm,&size);
174: if (size > 1) {
175: /* create a sequential iscoloring on all processors */
176: MatGetSeqNonzeroStructure(mat,&mat_seq);
177: }
179: MatGetRowIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,&n,&ria,&rja,&done);
180: MatGetColumnIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,&n,&cia,&cja,&done);
181: if (!done) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Ordering requires IJ");
183: MatFDColoringDegreeSequence_Minpack(n,cja,cia,rja,ria,&seq);
185: PetscMalloc2(n,&list,4*n,&work);
187: n1 = n - 1;
188: none = -1;
189: MINPACKnumsrt(&n,&n1,seq,&none,list,work+2*n,work+n);
190: PetscMalloc1(n,&coloring);
191: MINPACKseq(&n,cja,cia,rja,ria,list,coloring,&ncolors,work);
193: PetscFree2(list,work);
194: PetscFree(seq);
196: MatRestoreRowIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,NULL,&ria,&rja,&done);
197: MatRestoreColumnIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,NULL,&cia,&cja,&done);
199: /* shift coloring numbers to start at zero and shorten */
200: if (ncolors > IS_COLORING_MAX-1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Maximum color size exceeded");
201: {
202: ISColoringValue *s = (ISColoringValue*) coloring;
203: for (i=0; i<n; i++) s[i] = (ISColoringValue) (coloring[i]-1);
204: MatColoringPatch(mat_seq,ncolors,n,s,iscoloring);
205: }
207: if (size > 1) {
208: MatDestroySeqNonzeroStructure(&mat_seq);
210: /* convert iscoloring_seq to a parallel iscoloring */
211: iscoloring_seq = *iscoloring;
212: rstart = mat->rmap->rstart/bs;
213: rend = mat->rmap->rend/bs;
214: N_loc = rend - rstart; /* number of local nodes */
216: /* get local colors for each local node */
217: PetscMalloc1((N_loc+1),&colors_loc);
218: for (i=rstart; i<rend; i++) colors_loc[i-rstart] = iscoloring_seq->colors[i];
220: /* create a parallel iscoloring */
221: nc = iscoloring_seq->n;
222: ISColoringCreate(comm,nc,N_loc,colors_loc,iscoloring);
223: ISColoringDestroy(&iscoloring_seq);
224: }
225: return(0);
226: }
230: PETSC_EXTERN PetscErrorCode MatColoringCreate_LF(MatColoring mc)
231: {
233: mc->dist = 2;
234: mc->data = NULL;
235: mc->ops->apply = MatColoringApply_LF;
236: mc->ops->view = NULL;
237: mc->ops->destroy = NULL;
238: mc->ops->setfromoptions = NULL;
239: return(0);
240: }
244: PETSC_EXTERN PetscErrorCode MatColoringApply_ID(MatColoring mc,ISColoring *iscoloring)
245: {
246: PetscErrorCode ierr;
247: PetscInt *list,*work,clique,*seq,*coloring,n;
248: const PetscInt *ria,*rja,*cia,*cja;
249: PetscInt ncolors,i;
250: PetscBool done;
251: Mat mat = mc->mat;
252: Mat mat_seq = mat;
253: PetscMPIInt size;
254: MPI_Comm comm;
255: ISColoring iscoloring_seq;
256: PetscInt bs = 1,rstart,rend,N_loc,nc;
257: ISColoringValue *colors_loc;
258: PetscBool flg1,flg2;
261: if (mc->dist != 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"IDO may only do distance 2 coloring");
262: /* this is ugly way to get blocksize but cannot call MatGetBlockSize() because AIJ can have bs > 1 */
263: PetscObjectTypeCompare((PetscObject)mat,MATSEQBAIJ,&flg1);
264: PetscObjectTypeCompare((PetscObject)mat,MATMPIBAIJ,&flg2);
265: if (flg1 || flg2) {
266: MatGetBlockSize(mat,&bs);
267: }
269: PetscObjectGetComm((PetscObject)mat,&comm);
270: MPI_Comm_size(comm,&size);
271: if (size > 1) {
272: /* create a sequential iscoloring on all processors */
273: MatGetSeqNonzeroStructure(mat,&mat_seq);
274: }
276: MatGetRowIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,&n,&ria,&rja,&done);
277: MatGetColumnIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,&n,&cia,&cja,&done);
278: if (!done) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Ordering requires IJ");
280: MatFDColoringDegreeSequence_Minpack(n,cja,cia,rja,ria,&seq);
282: PetscMalloc2(n,&list,4*n,&work);
284: MINPACKido(&n,&n,cja,cia,rja,ria,seq,list,&clique,work,work+n,work+2*n,work+3*n);
286: PetscMalloc1(n,&coloring);
287: MINPACKseq(&n,cja,cia,rja,ria,list,coloring,&ncolors,work);
289: PetscFree2(list,work);
290: PetscFree(seq);
292: MatRestoreRowIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,NULL,&ria,&rja,&done);
293: MatRestoreColumnIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,NULL,&cia,&cja,&done);
295: /* shift coloring numbers to start at zero and shorten */
296: if (ncolors > IS_COLORING_MAX-1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Maximum color size exceeded");
297: {
298: ISColoringValue *s = (ISColoringValue*) coloring;
299: for (i=0; i<n; i++) {
300: s[i] = (ISColoringValue) (coloring[i]-1);
301: }
302: MatColoringPatch(mat_seq,ncolors,n,s,iscoloring);
303: }
305: if (size > 1) {
306: MatDestroySeqNonzeroStructure(&mat_seq);
308: /* convert iscoloring_seq to a parallel iscoloring */
309: iscoloring_seq = *iscoloring;
310: rstart = mat->rmap->rstart/bs;
311: rend = mat->rmap->rend/bs;
312: N_loc = rend - rstart; /* number of local nodes */
314: /* get local colors for each local node */
315: PetscMalloc1((N_loc+1),&colors_loc);
316: for (i=rstart; i<rend; i++) {
317: colors_loc[i-rstart] = iscoloring_seq->colors[i];
318: }
319: /* create a parallel iscoloring */
320: nc = iscoloring_seq->n;
321: ISColoringCreate(comm,nc,N_loc,colors_loc,iscoloring);
322: ISColoringDestroy(&iscoloring_seq);
323: }
324: return(0);
325: }
330: PETSC_EXTERN PetscErrorCode MatColoringCreate_ID(MatColoring mc)
331: {
333: mc->dist = 2;
334: mc->data = NULL;
335: mc->ops->apply = MatColoringApply_ID;
336: mc->ops->view = NULL;
337: mc->ops->destroy = NULL;
338: mc->ops->setfromoptions = NULL;
339: return(0);
340: }