Actual source code: color.c
petsc-3.7.1 2016-05-15
2: /*
3: Routines that call the kernel minpack coloring subroutines
4: */
6: #include <petsc/private/matimpl.h>
7: #include <petsc/private/isimpl.h>
8: #include <../src/mat/color/impls/minpack/color.h>
10: /*
11: MatFDColoringDegreeSequence_Minpack - Calls the MINPACK routine seqr() that
12: computes the degree sequence required by MINPACK coloring routines.
13: */
16: PETSC_INTERN PetscErrorCode MatFDColoringDegreeSequence_Minpack(PetscInt m,const PetscInt *cja,const PetscInt *cia,const PetscInt *rja,const PetscInt *ria,PetscInt **seq)
17: {
18: PetscInt *work;
22: PetscMalloc1(m,&work);
23: PetscMalloc1(m,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++) c = PetscMax(c,ia[i+1]-ia[i]);
44: *minc = c;
45: return(0);
46: }
50: static PetscErrorCode MatColoringApply_SL(MatColoring mc,ISColoring *iscoloring)
51: {
52: PetscErrorCode ierr;
53: PetscInt *list,*work,clique,*seq,*coloring,n;
54: const PetscInt *ria,*rja,*cia,*cja;
55: PetscInt ncolors,i;
56: PetscBool done;
57: Mat mat = mc->mat;
58: Mat mat_seq = mat;
59: PetscMPIInt size;
60: MPI_Comm comm;
61: ISColoring iscoloring_seq;
62: PetscInt bs = 1,rstart,rend,N_loc,nc;
63: ISColoringValue *colors_loc;
64: PetscBool flg1,flg2;
67: if (mc->dist != 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"SL may only do distance 2 coloring");
68: /* this is ugly way to get blocksize but cannot call MatGetBlockSize() because AIJ can have bs > 1 */
69: PetscObjectTypeCompare((PetscObject)mat,MATSEQBAIJ,&flg1);
70: PetscObjectTypeCompare((PetscObject)mat,MATMPIBAIJ,&flg2);
71: if (flg1 || flg2) {
72: MatGetBlockSize(mat,&bs);
73: }
75: PetscObjectGetComm((PetscObject)mat,&comm);
76: MPI_Comm_size(comm,&size);
77: if (size > 1) {
78: /* create a sequential iscoloring on all processors */
79: MatGetSeqNonzeroStructure(mat,&mat_seq);
80: }
82: MatGetRowIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,&n,&ria,&rja,&done);
83: MatGetColumnIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,&n,&cia,&cja,&done);
84: if (!done) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Ordering requires IJ");
86: MatFDColoringDegreeSequence_Minpack(n,cja,cia,rja,ria,&seq);
88: PetscMalloc2(n,&list,4*n,&work);
90: MINPACKslo(&n,cja,cia,rja,ria,seq,list,&clique,work,work+n,work+2*n,work+3*n);
92: PetscMalloc1(n,&coloring);
93: MINPACKseq(&n,cja,cia,rja,ria,list,coloring,&ncolors,work);
95: PetscFree2(list,work);
96: PetscFree(seq);
97: MatRestoreRowIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,NULL,&ria,&rja,&done);
98: MatRestoreColumnIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,NULL,&cia,&cja,&done);
100: /* shift coloring numbers to start at zero and shorten */
101: if (ncolors > IS_COLORING_MAX-1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Maximum color size exceeded");
102: {
103: ISColoringValue *s = (ISColoringValue*) coloring;
104: for (i=0; i<n; i++) {
105: s[i] = (ISColoringValue) (coloring[i]-1);
106: }
107: MatColoringPatch(mat_seq,ncolors,n,s,iscoloring);
108: }
110: if (size > 1) {
111: MatDestroySeqNonzeroStructure(&mat_seq);
113: /* convert iscoloring_seq to a parallel iscoloring */
114: iscoloring_seq = *iscoloring;
115: rstart = mat->rmap->rstart/bs;
116: rend = mat->rmap->rend/bs;
117: N_loc = rend - rstart; /* number of local nodes */
119: /* get local colors for each local node */
120: PetscMalloc1(N_loc+1,&colors_loc);
121: for (i=rstart; i<rend; i++) {
122: colors_loc[i-rstart] = iscoloring_seq->colors[i];
123: }
124: /* create a parallel iscoloring */
125: nc = iscoloring_seq->n;
126: ISColoringCreate(comm,nc,N_loc,colors_loc,PETSC_OWN_POINTER,iscoloring);
127: ISColoringDestroy(&iscoloring_seq);
128: }
129: return(0);
130: }
134: PETSC_EXTERN PetscErrorCode MatColoringCreate_SL(MatColoring mc)
135: {
137: mc->dist = 2;
138: mc->data = NULL;
139: mc->ops->apply = MatColoringApply_SL;
140: mc->ops->view = NULL;
141: mc->ops->destroy = NULL;
142: mc->ops->setfromoptions = NULL;
143: return(0);
144: }
148: static PetscErrorCode MatColoringApply_LF(MatColoring mc,ISColoring *iscoloring)
149: {
150: PetscErrorCode ierr;
151: PetscInt *list,*work,*seq,*coloring,n;
152: const PetscInt *ria,*rja,*cia,*cja;
153: PetscInt n1, none,ncolors,i;
154: PetscBool done;
155: Mat mat = mc->mat;
156: Mat mat_seq = mat;
157: PetscMPIInt size;
158: MPI_Comm comm;
159: ISColoring iscoloring_seq;
160: PetscInt bs = 1,rstart,rend,N_loc,nc;
161: ISColoringValue *colors_loc;
162: PetscBool flg1,flg2;
165: if (mc->dist != 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"LF may only do distance 2 coloring");
166: /* this is ugly way to get blocksize but cannot call MatGetBlockSize() because AIJ can have bs > 1 */
167: PetscObjectTypeCompare((PetscObject)mat,MATSEQBAIJ,&flg1);
168: PetscObjectTypeCompare((PetscObject)mat,MATMPIBAIJ,&flg2);
169: if (flg1 || flg2) {
170: MatGetBlockSize(mat,&bs);
171: }
173: PetscObjectGetComm((PetscObject)mat,&comm);
174: MPI_Comm_size(comm,&size);
175: if (size > 1) {
176: /* create a sequential iscoloring on all processors */
177: MatGetSeqNonzeroStructure(mat,&mat_seq);
178: }
180: MatGetRowIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,&n,&ria,&rja,&done);
181: MatGetColumnIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,&n,&cia,&cja,&done);
182: if (!done) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Ordering requires IJ");
184: MatFDColoringDegreeSequence_Minpack(n,cja,cia,rja,ria,&seq);
186: PetscMalloc2(n,&list,4*n,&work);
188: n1 = n - 1;
189: none = -1;
190: MINPACKnumsrt(&n,&n1,seq,&none,list,work+2*n,work+n);
191: PetscMalloc1(n,&coloring);
192: MINPACKseq(&n,cja,cia,rja,ria,list,coloring,&ncolors,work);
194: PetscFree2(list,work);
195: PetscFree(seq);
197: MatRestoreRowIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,NULL,&ria,&rja,&done);
198: MatRestoreColumnIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,NULL,&cia,&cja,&done);
200: /* shift coloring numbers to start at zero and shorten */
201: if (ncolors > IS_COLORING_MAX-1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Maximum color size exceeded");
202: {
203: ISColoringValue *s = (ISColoringValue*) coloring;
204: for (i=0; i<n; i++) s[i] = (ISColoringValue) (coloring[i]-1);
205: MatColoringPatch(mat_seq,ncolors,n,s,iscoloring);
206: }
208: if (size > 1) {
209: MatDestroySeqNonzeroStructure(&mat_seq);
211: /* convert iscoloring_seq to a parallel iscoloring */
212: iscoloring_seq = *iscoloring;
213: rstart = mat->rmap->rstart/bs;
214: rend = mat->rmap->rend/bs;
215: N_loc = rend - rstart; /* number of local nodes */
217: /* get local colors for each local node */
218: PetscMalloc1(N_loc+1,&colors_loc);
219: for (i=rstart; i<rend; i++) colors_loc[i-rstart] = iscoloring_seq->colors[i];
221: /* create a parallel iscoloring */
222: nc = iscoloring_seq->n;
223: ISColoringCreate(comm,nc,N_loc,colors_loc,PETSC_OWN_POINTER,iscoloring);
224: ISColoringDestroy(&iscoloring_seq);
225: }
226: return(0);
227: }
231: PETSC_EXTERN PetscErrorCode MatColoringCreate_LF(MatColoring mc)
232: {
234: mc->dist = 2;
235: mc->data = NULL;
236: mc->ops->apply = MatColoringApply_LF;
237: mc->ops->view = NULL;
238: mc->ops->destroy = NULL;
239: mc->ops->setfromoptions = NULL;
240: return(0);
241: }
245: static PetscErrorCode MatColoringApply_ID(MatColoring mc,ISColoring *iscoloring)
246: {
247: PetscErrorCode ierr;
248: PetscInt *list,*work,clique,*seq,*coloring,n;
249: const PetscInt *ria,*rja,*cia,*cja;
250: PetscInt ncolors,i;
251: PetscBool done;
252: Mat mat = mc->mat;
253: Mat mat_seq = mat;
254: PetscMPIInt size;
255: MPI_Comm comm;
256: ISColoring iscoloring_seq;
257: PetscInt bs = 1,rstart,rend,N_loc,nc;
258: ISColoringValue *colors_loc;
259: PetscBool flg1,flg2;
262: if (mc->dist != 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"IDO may only do distance 2 coloring");
263: /* this is ugly way to get blocksize but cannot call MatGetBlockSize() because AIJ can have bs > 1 */
264: PetscObjectTypeCompare((PetscObject)mat,MATSEQBAIJ,&flg1);
265: PetscObjectTypeCompare((PetscObject)mat,MATMPIBAIJ,&flg2);
266: if (flg1 || flg2) {
267: MatGetBlockSize(mat,&bs);
268: }
270: PetscObjectGetComm((PetscObject)mat,&comm);
271: MPI_Comm_size(comm,&size);
272: if (size > 1) {
273: /* create a sequential iscoloring on all processors */
274: MatGetSeqNonzeroStructure(mat,&mat_seq);
275: }
277: MatGetRowIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,&n,&ria,&rja,&done);
278: MatGetColumnIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,&n,&cia,&cja,&done);
279: if (!done) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Ordering requires IJ");
281: MatFDColoringDegreeSequence_Minpack(n,cja,cia,rja,ria,&seq);
283: PetscMalloc2(n,&list,4*n,&work);
285: MINPACKido(&n,&n,cja,cia,rja,ria,seq,list,&clique,work,work+n,work+2*n,work+3*n);
287: PetscMalloc1(n,&coloring);
288: MINPACKseq(&n,cja,cia,rja,ria,list,coloring,&ncolors,work);
290: PetscFree2(list,work);
291: PetscFree(seq);
293: MatRestoreRowIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,NULL,&ria,&rja,&done);
294: MatRestoreColumnIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,NULL,&cia,&cja,&done);
296: /* shift coloring numbers to start at zero and shorten */
297: if (ncolors > IS_COLORING_MAX-1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Maximum color size exceeded");
298: {
299: ISColoringValue *s = (ISColoringValue*) coloring;
300: for (i=0; i<n; i++) {
301: s[i] = (ISColoringValue) (coloring[i]-1);
302: }
303: MatColoringPatch(mat_seq,ncolors,n,s,iscoloring);
304: }
306: if (size > 1) {
307: MatDestroySeqNonzeroStructure(&mat_seq);
309: /* convert iscoloring_seq to a parallel iscoloring */
310: iscoloring_seq = *iscoloring;
311: rstart = mat->rmap->rstart/bs;
312: rend = mat->rmap->rend/bs;
313: N_loc = rend - rstart; /* number of local nodes */
315: /* get local colors for each local node */
316: PetscMalloc1(N_loc+1,&colors_loc);
317: for (i=rstart; i<rend; i++) {
318: colors_loc[i-rstart] = iscoloring_seq->colors[i];
319: }
320: /* create a parallel iscoloring */
321: nc = iscoloring_seq->n;
322: ISColoringCreate(comm,nc,N_loc,colors_loc,PETSC_OWN_POINTER,iscoloring);
323: ISColoringDestroy(&iscoloring_seq);
324: }
325: return(0);
326: }
331: PETSC_EXTERN PetscErrorCode MatColoringCreate_ID(MatColoring mc)
332: {
334: mc->dist = 2;
335: mc->data = NULL;
336: mc->ops->apply = MatColoringApply_ID;
337: mc->ops->view = NULL;
338: mc->ops->destroy = NULL;
339: mc->ops->setfromoptions = NULL;
340: return(0);
341: }