Actual source code: color.c
petsc-3.8.3 2017-12-09
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: */
14: PETSC_INTERN PetscErrorCode MatFDColoringDegreeSequence_Minpack(PetscInt m,const PetscInt *cja,const PetscInt *cia,const PetscInt *rja,const PetscInt *ria,PetscInt **seq)
15: {
16: PetscInt *work;
20: PetscMalloc1(m,&work);
21: PetscMalloc1(m,seq);
23: MINPACKdegr(&m,cja,cia,rja,ria,*seq,work);
25: PetscFree(work);
26: return(0);
27: }
29: /*
30: MatFDColoringMinimumNumberofColors_Private - For a given sparse
31: matrix computes the minimum number of colors needed.
33: */
34: PetscErrorCode MatFDColoringMinimumNumberofColors_Private(PetscInt m,PetscInt *ia,PetscInt *minc)
35: {
36: PetscInt i,c = 0;
39: for (i=0; i<m; i++) c = PetscMax(c,ia[i+1]-ia[i]);
40: *minc = c;
41: return(0);
42: }
44: static PetscErrorCode MatColoringApply_SL(MatColoring mc,ISColoring *iscoloring)
45: {
46: PetscErrorCode ierr;
47: PetscInt *list,*work,clique,*seq,*coloring,n;
48: const PetscInt *ria,*rja,*cia,*cja;
49: PetscInt ncolors,i;
50: PetscBool done;
51: Mat mat = mc->mat;
52: Mat mat_seq = mat;
53: PetscMPIInt size;
54: MPI_Comm comm;
55: ISColoring iscoloring_seq;
56: PetscInt bs = 1,rstart,rend,N_loc,nc;
57: ISColoringValue *colors_loc;
58: PetscBool flg1,flg2;
61: if (mc->dist != 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"SL may only do distance 2 coloring");
62: /* this is ugly way to get blocksize but cannot call MatGetBlockSize() because AIJ can have bs > 1 */
63: PetscObjectTypeCompare((PetscObject)mat,MATSEQBAIJ,&flg1);
64: PetscObjectTypeCompare((PetscObject)mat,MATMPIBAIJ,&flg2);
65: if (flg1 || flg2) {
66: MatGetBlockSize(mat,&bs);
67: }
69: PetscObjectGetComm((PetscObject)mat,&comm);
70: MPI_Comm_size(comm,&size);
71: if (size > 1) {
72: /* create a sequential iscoloring on all processors */
73: MatGetSeqNonzeroStructure(mat,&mat_seq);
74: }
76: MatGetRowIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,&n,&ria,&rja,&done);
77: MatGetColumnIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,&n,&cia,&cja,&done);
78: if (!done) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Ordering requires IJ");
80: MatFDColoringDegreeSequence_Minpack(n,cja,cia,rja,ria,&seq);
82: PetscMalloc2(n,&list,4*n,&work);
84: MINPACKslo(&n,cja,cia,rja,ria,seq,list,&clique,work,work+n,work+2*n,work+3*n);
86: PetscMalloc1(n,&coloring);
87: MINPACKseq(&n,cja,cia,rja,ria,list,coloring,&ncolors,work);
89: PetscFree2(list,work);
90: PetscFree(seq);
91: MatRestoreRowIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,NULL,&ria,&rja,&done);
92: MatRestoreColumnIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,NULL,&cia,&cja,&done);
94: /* shift coloring numbers to start at zero and shorten */
95: if (ncolors > IS_COLORING_MAX-1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Maximum color size exceeded");
96: {
97: ISColoringValue *s = (ISColoringValue*) coloring;
98: for (i=0; i<n; i++) {
99: s[i] = (ISColoringValue) (coloring[i]-1);
100: }
101: MatColoringPatch(mat_seq,ncolors,n,s,iscoloring);
102: }
104: if (size > 1) {
105: MatDestroySeqNonzeroStructure(&mat_seq);
107: /* convert iscoloring_seq to a parallel iscoloring */
108: iscoloring_seq = *iscoloring;
109: rstart = mat->rmap->rstart/bs;
110: rend = mat->rmap->rend/bs;
111: N_loc = rend - rstart; /* number of local nodes */
113: /* get local colors for each local node */
114: PetscMalloc1(N_loc+1,&colors_loc);
115: for (i=rstart; i<rend; i++) {
116: colors_loc[i-rstart] = iscoloring_seq->colors[i];
117: }
118: /* create a parallel iscoloring */
119: nc = iscoloring_seq->n;
120: ISColoringCreate(comm,nc,N_loc,colors_loc,PETSC_OWN_POINTER,iscoloring);
121: ISColoringDestroy(&iscoloring_seq);
122: }
123: return(0);
124: }
126: PETSC_EXTERN PetscErrorCode MatColoringCreate_SL(MatColoring mc)
127: {
129: mc->dist = 2;
130: mc->data = NULL;
131: mc->ops->apply = MatColoringApply_SL;
132: mc->ops->view = NULL;
133: mc->ops->destroy = NULL;
134: mc->ops->setfromoptions = NULL;
135: return(0);
136: }
138: static PetscErrorCode MatColoringApply_LF(MatColoring mc,ISColoring *iscoloring)
139: {
140: PetscErrorCode ierr;
141: PetscInt *list,*work,*seq,*coloring,n;
142: const PetscInt *ria,*rja,*cia,*cja;
143: PetscInt n1, none,ncolors,i;
144: PetscBool done;
145: Mat mat = mc->mat;
146: Mat mat_seq = mat;
147: PetscMPIInt size;
148: MPI_Comm comm;
149: ISColoring iscoloring_seq;
150: PetscInt bs = 1,rstart,rend,N_loc,nc;
151: ISColoringValue *colors_loc;
152: PetscBool flg1,flg2;
155: if (mc->dist != 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"LF may only do distance 2 coloring");
156: /* this is ugly way to get blocksize but cannot call MatGetBlockSize() because AIJ can have bs > 1 */
157: PetscObjectTypeCompare((PetscObject)mat,MATSEQBAIJ,&flg1);
158: PetscObjectTypeCompare((PetscObject)mat,MATMPIBAIJ,&flg2);
159: if (flg1 || flg2) {
160: MatGetBlockSize(mat,&bs);
161: }
163: PetscObjectGetComm((PetscObject)mat,&comm);
164: MPI_Comm_size(comm,&size);
165: if (size > 1) {
166: /* create a sequential iscoloring on all processors */
167: MatGetSeqNonzeroStructure(mat,&mat_seq);
168: }
170: MatGetRowIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,&n,&ria,&rja,&done);
171: MatGetColumnIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,&n,&cia,&cja,&done);
172: if (!done) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Ordering requires IJ");
174: MatFDColoringDegreeSequence_Minpack(n,cja,cia,rja,ria,&seq);
176: PetscMalloc2(n,&list,4*n,&work);
178: n1 = n - 1;
179: none = -1;
180: MINPACKnumsrt(&n,&n1,seq,&none,list,work+2*n,work+n);
181: PetscMalloc1(n,&coloring);
182: MINPACKseq(&n,cja,cia,rja,ria,list,coloring,&ncolors,work);
184: PetscFree2(list,work);
185: PetscFree(seq);
187: MatRestoreRowIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,NULL,&ria,&rja,&done);
188: MatRestoreColumnIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,NULL,&cia,&cja,&done);
190: /* shift coloring numbers to start at zero and shorten */
191: if (ncolors > IS_COLORING_MAX-1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Maximum color size exceeded");
192: {
193: ISColoringValue *s = (ISColoringValue*) coloring;
194: for (i=0; i<n; i++) s[i] = (ISColoringValue) (coloring[i]-1);
195: MatColoringPatch(mat_seq,ncolors,n,s,iscoloring);
196: }
198: if (size > 1) {
199: MatDestroySeqNonzeroStructure(&mat_seq);
201: /* convert iscoloring_seq to a parallel iscoloring */
202: iscoloring_seq = *iscoloring;
203: rstart = mat->rmap->rstart/bs;
204: rend = mat->rmap->rend/bs;
205: N_loc = rend - rstart; /* number of local nodes */
207: /* get local colors for each local node */
208: PetscMalloc1(N_loc+1,&colors_loc);
209: for (i=rstart; i<rend; i++) colors_loc[i-rstart] = iscoloring_seq->colors[i];
211: /* create a parallel iscoloring */
212: nc = iscoloring_seq->n;
213: ISColoringCreate(comm,nc,N_loc,colors_loc,PETSC_OWN_POINTER,iscoloring);
214: ISColoringDestroy(&iscoloring_seq);
215: }
216: return(0);
217: }
219: PETSC_EXTERN PetscErrorCode MatColoringCreate_LF(MatColoring mc)
220: {
222: mc->dist = 2;
223: mc->data = NULL;
224: mc->ops->apply = MatColoringApply_LF;
225: mc->ops->view = NULL;
226: mc->ops->destroy = NULL;
227: mc->ops->setfromoptions = NULL;
228: return(0);
229: }
231: static PetscErrorCode MatColoringApply_ID(MatColoring mc,ISColoring *iscoloring)
232: {
233: PetscErrorCode ierr;
234: PetscInt *list,*work,clique,*seq,*coloring,n;
235: const PetscInt *ria,*rja,*cia,*cja;
236: PetscInt ncolors,i;
237: PetscBool done;
238: Mat mat = mc->mat;
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: PetscBool flg1,flg2;
248: if (mc->dist != 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"IDO may only do distance 2 coloring");
249: /* this is ugly way to get blocksize but cannot call MatGetBlockSize() because AIJ can have bs > 1 */
250: PetscObjectTypeCompare((PetscObject)mat,MATSEQBAIJ,&flg1);
251: PetscObjectTypeCompare((PetscObject)mat,MATMPIBAIJ,&flg2);
252: if (flg1 || flg2) {
253: MatGetBlockSize(mat,&bs);
254: }
256: PetscObjectGetComm((PetscObject)mat,&comm);
257: MPI_Comm_size(comm,&size);
258: if (size > 1) {
259: /* create a sequential iscoloring on all processors */
260: MatGetSeqNonzeroStructure(mat,&mat_seq);
261: }
263: MatGetRowIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,&n,&ria,&rja,&done);
264: MatGetColumnIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,&n,&cia,&cja,&done);
265: if (!done) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Ordering requires IJ");
267: MatFDColoringDegreeSequence_Minpack(n,cja,cia,rja,ria,&seq);
269: PetscMalloc2(n,&list,4*n,&work);
271: MINPACKido(&n,&n,cja,cia,rja,ria,seq,list,&clique,work,work+n,work+2*n,work+3*n);
273: PetscMalloc1(n,&coloring);
274: MINPACKseq(&n,cja,cia,rja,ria,list,coloring,&ncolors,work);
276: PetscFree2(list,work);
277: PetscFree(seq);
279: MatRestoreRowIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,NULL,&ria,&rja,&done);
280: MatRestoreColumnIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,NULL,&cia,&cja,&done);
282: /* shift coloring numbers to start at zero and shorten */
283: if (ncolors > IS_COLORING_MAX-1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Maximum color size exceeded");
284: {
285: ISColoringValue *s = (ISColoringValue*) coloring;
286: for (i=0; i<n; i++) {
287: s[i] = (ISColoringValue) (coloring[i]-1);
288: }
289: MatColoringPatch(mat_seq,ncolors,n,s,iscoloring);
290: }
292: if (size > 1) {
293: MatDestroySeqNonzeroStructure(&mat_seq);
295: /* convert iscoloring_seq to a parallel iscoloring */
296: iscoloring_seq = *iscoloring;
297: rstart = mat->rmap->rstart/bs;
298: rend = mat->rmap->rend/bs;
299: N_loc = rend - rstart; /* number of local nodes */
301: /* get local colors for each local node */
302: PetscMalloc1(N_loc+1,&colors_loc);
303: for (i=rstart; i<rend; i++) {
304: colors_loc[i-rstart] = iscoloring_seq->colors[i];
305: }
306: /* create a parallel iscoloring */
307: nc = iscoloring_seq->n;
308: ISColoringCreate(comm,nc,N_loc,colors_loc,PETSC_OWN_POINTER,iscoloring);
309: ISColoringDestroy(&iscoloring_seq);
310: }
311: return(0);
312: }
315: PETSC_EXTERN PetscErrorCode MatColoringCreate_ID(MatColoring mc)
316: {
318: mc->dist = 2;
319: mc->data = NULL;
320: mc->ops->apply = MatColoringApply_ID;
321: mc->ops->view = NULL;
322: mc->ops->destroy = NULL;
323: mc->ops->setfromoptions = NULL;
324: return(0);
325: }