Actual source code: greedy.c

petsc-3.10.5 2019-03-28
Report Typos and Errors
  1:  #include <petsc/private/matimpl.h>
  2:  #include <../src/mat/impls/aij/seq/aij.h>
  3:  #include <../src/mat/impls/aij/mpi/mpiaij.h>
  4:  #include <petscsf.h>

  6: typedef struct {
  7:   PetscBool symmetric;
  8: } MC_Greedy;

 10: static PetscErrorCode MatColoringDestroy_Greedy(MatColoring mc)
 11: {

 15:   PetscFree(mc->data);
 16:   return(0);
 17: }

 19: static PetscErrorCode GreedyColoringLocalDistanceOne_Private(MatColoring mc,PetscReal *wts,PetscInt *lperm,ISColoringValue *colors)
 20: {
 21:   PetscInt        i,j,k,s,e,n,no,nd,nd_global,n_global,idx,ncols,maxcolors,masksize,ccol,*mask;
 22:   PetscErrorCode  ierr;
 23:   Mat             m=mc->mat;
 24:   Mat_MPIAIJ      *aij = (Mat_MPIAIJ*)m->data;
 25:   Mat             md=NULL,mo=NULL;
 26:   const PetscInt  *md_i,*mo_i,*md_j,*mo_j;
 27:   PetscBool       isMPIAIJ,isSEQAIJ;
 28:   ISColoringValue pcol;
 29:   const PetscInt  *cidx;
 30:   PetscInt        *lcolors,*ocolors;
 31:   PetscReal       *owts=NULL;
 32:   PetscSF         sf;
 33:   PetscLayout     layout;

 36:   MatGetSize(m,&n_global,NULL);
 37:   MatGetOwnershipRange(m,&s,&e);
 38:   n=e-s;
 39:   masksize=20;
 40:   nd_global = 0;
 41:   /* get the matrix communication structures */
 42:   PetscObjectTypeCompare((PetscObject)m, MATMPIAIJ, &isMPIAIJ);
 43:   PetscObjectBaseTypeCompare((PetscObject)m, MATSEQAIJ, &isSEQAIJ);
 44:   if (isMPIAIJ) {
 45:     /* get the CSR data for on and off diagonal portions of m */
 46:     Mat_SeqAIJ *dseq;
 47:     Mat_SeqAIJ *oseq;
 48:     md=aij->A;
 49:     dseq = (Mat_SeqAIJ*)md->data;
 50:     mo=aij->B;
 51:     oseq = (Mat_SeqAIJ*)mo->data;
 52:     md_i = dseq->i;
 53:     md_j = dseq->j;
 54:     mo_i = oseq->i;
 55:     mo_j = oseq->j;
 56:   } else if (isSEQAIJ) {
 57:     /* get the CSR data for m */
 58:     Mat_SeqAIJ *dseq;
 59:     /* no off-processor nodes */
 60:     md=m;
 61:     dseq = (Mat_SeqAIJ*)md->data;
 62:     mo=NULL;
 63:     no=0;
 64:     md_i = dseq->i;
 65:     md_j = dseq->j;
 66:     mo_i = NULL;
 67:     mo_j = NULL;
 68:   } else SETERRQ(PetscObjectComm((PetscObject)mc),PETSC_ERR_ARG_WRONG,"Matrix must be AIJ for greedy coloring");
 69:   MatColoringGetMaxColors(mc,&maxcolors);
 70:   if (mo) {
 71:     VecGetSize(aij->lvec,&no);
 72:     PetscMalloc2(no,&ocolors,no,&owts);
 73:     for(i=0;i<no;i++) {
 74:       ocolors[i]=maxcolors;
 75:     }
 76:   }

 78:   PetscMalloc1(masksize,&mask);
 79:   PetscMalloc1(n,&lcolors);
 80:   for(i=0;i<n;i++) {
 81:     lcolors[i]=maxcolors;
 82:   }
 83:   for (i=0;i<masksize;i++) {
 84:     mask[i]=-1;
 85:   }
 86:   if (mo) {
 87:     /* transfer neighbor weights */
 88:     PetscSFCreate(PetscObjectComm((PetscObject)m),&sf);
 89:     MatGetLayouts(m,&layout,NULL);
 90:     PetscSFSetGraphLayout(sf,layout,no,NULL,PETSC_COPY_VALUES,aij->garray);
 91:     PetscSFBcastBegin(sf,MPIU_REAL,wts,owts);
 92:     PetscSFBcastEnd(sf,MPIU_REAL,wts,owts);
 93:   }
 94:   while (nd_global < n_global) {
 95:     nd=n;
 96:     /* assign lowest possible color to each local vertex */
 97:     PetscLogEventBegin(MATCOLORING_Local,mc,0,0,0);
 98:     for (i=0;i<n;i++) {
 99:       idx=lperm[i];
100:       if (lcolors[idx] == maxcolors) {
101:         ncols = md_i[idx+1]-md_i[idx];
102:         cidx = &(md_j[md_i[idx]]);
103:         for (j=0;j<ncols;j++) {
104:           if (lcolors[cidx[j]] != maxcolors) {
105:             ccol=lcolors[cidx[j]];
106:             if (ccol>=masksize) {
107:               PetscInt *newmask;
108:               PetscMalloc1(masksize*2,&newmask);
109:               for(k=0;k<2*masksize;k++) {
110:                 newmask[k]=-1;
111:               }
112:               for(k=0;k<masksize;k++) {
113:                 newmask[k]=mask[k];
114:               }
115:               PetscFree(mask);
116:               mask=newmask;
117:               masksize*=2;
118:             }
119:             mask[ccol]=idx;
120:           }
121:         }
122:         if (mo) {
123:           ncols = mo_i[idx+1]-mo_i[idx];
124:           cidx = &(mo_j[mo_i[idx]]);
125:           for (j=0;j<ncols;j++) {
126:             if (ocolors[cidx[j]] != maxcolors) {
127:               ccol=ocolors[cidx[j]];
128:               if (ccol>=masksize) {
129:                 PetscInt *newmask;
130:                 PetscMalloc1(masksize*2,&newmask);
131:                 for(k=0;k<2*masksize;k++) {
132:                   newmask[k]=-1;
133:                 }
134:                 for(k=0;k<masksize;k++) {
135:                   newmask[k]=mask[k];
136:                 }
137:                 PetscFree(mask);
138:                 mask=newmask;
139:                 masksize*=2;
140:               }
141:               mask[ccol]=idx;
142:             }
143:           }
144:         }
145:         for (j=0;j<masksize;j++) {
146:           if (mask[j]!=idx) {
147:             break;
148:           }
149:         }
150:         pcol=j;
151:         if (pcol>maxcolors)pcol=maxcolors;
152:         lcolors[idx]=pcol;
153:       }
154:     }
155:     PetscLogEventEnd(MATCOLORING_Local,mc,0,0,0);
156:     if (mo) {
157:       /* transfer neighbor colors */
158:       PetscLogEventBegin(MATCOLORING_Comm,mc,0,0,0);
159:       PetscSFBcastBegin(sf,MPIU_INT,lcolors,ocolors);
160:       PetscSFBcastEnd(sf,MPIU_INT,lcolors,ocolors);
161:       /* check for conflicts -- this is merely checking if any adjacent off-processor rows have the same color and marking the ones that are lower weight locally for changing */
162:       for (i=0;i<n;i++) {
163:         ncols = mo_i[i+1]-mo_i[i];
164:         cidx = &(mo_j[mo_i[i]]);
165:         for (j=0;j<ncols;j++) {
166:           /* in the case of conflicts, the highest weight one stays and the others go */
167:           if ((ocolors[cidx[j]] == lcolors[i]) && (owts[cidx[j]] > wts[i]) && lcolors[i] < maxcolors) {
168:             lcolors[i]=maxcolors;
169:             nd--;
170:           }
171:         }
172:       }
173:       nd_global=0;
174:     }
175:     MPIU_Allreduce(&nd,&nd_global,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)mc));
176:   }
177:   for (i=0;i<n;i++) {
178:     colors[i] = (ISColoringValue)lcolors[i];
179:   }
180:   PetscFree(mask);
181:   PetscFree(lcolors);
182:   if (mo) {
183:     PetscFree2(ocolors,owts);
184:     PetscSFDestroy(&sf);
185:   }
186:   return(0);
187: }

189: static PetscErrorCode GreedyColoringLocalDistanceTwo_Private(MatColoring mc,PetscReal *wts,PetscInt *lperm,ISColoringValue *colors)
190: {
191:   MC_Greedy       *gr = (MC_Greedy *) mc->data;
192:   PetscInt        i,j,k,l,s,e,n,nd,nd_global,n_global,idx,ncols,maxcolors,mcol,mcol_global,nd1cols,*mask,masksize,*d1cols,*bad,*badnext,nbad,badsize,ccol,no,cbad;
193:   Mat             m = mc->mat, mt;
194:   Mat_MPIAIJ      *aij = (Mat_MPIAIJ*)m->data;
195:   Mat             md=NULL,mo=NULL;
196:   const PetscInt  *md_i,*mo_i,*md_j,*mo_j;
197:   const PetscInt  *rmd_i,*rmo_i,*rmd_j,*rmo_j;
198:   PetscBool       isMPIAIJ,isSEQAIJ;
199:   PetscInt        pcol,*dcolors,*ocolors;
200:   ISColoringValue *badidx;
201:   const PetscInt  *cidx;
202:   PetscReal       *owts,*colorweights;
203:   PetscInt        *oconf,*conf;
204:   PetscSF         sf;
205:   PetscLayout     layout;
206:   PetscErrorCode  ierr;

209:   MatGetSize(m,&n_global,NULL);
210:   MatGetOwnershipRange(m,&s,&e);
211:   n=e-s;
212:   nd_global = 0;
213:   /* get the matrix communication structures */
214:   PetscObjectBaseTypeCompare((PetscObject)m, MATMPIAIJ, &isMPIAIJ);
215:   PetscObjectBaseTypeCompare((PetscObject)m, MATSEQAIJ, &isSEQAIJ);
216:   if (isMPIAIJ) {
217:     Mat_SeqAIJ *dseq;
218:     Mat_SeqAIJ *oseq;
219:     md=aij->A;
220:     dseq = (Mat_SeqAIJ*)md->data;
221:     mo=aij->B;
222:     oseq = (Mat_SeqAIJ*)mo->data;
223:     md_i = dseq->i;
224:     md_j = dseq->j;
225:     mo_i = oseq->i;
226:     mo_j = oseq->j;
227:     rmd_i = dseq->i;
228:     rmd_j = dseq->j;
229:     rmo_i = oseq->i;
230:     rmo_j = oseq->j;
231:   } else if (isSEQAIJ) {
232:     Mat_SeqAIJ *dseq;
233:     /* no off-processor nodes */
234:     md=m;
235:     dseq = (Mat_SeqAIJ*)md->data;
236:     md_i = dseq->i;
237:     md_j = dseq->j;
238:     mo_i = NULL;
239:     mo_j = NULL;
240:     rmd_i = dseq->i;
241:     rmd_j = dseq->j;
242:     rmo_i = NULL;
243:     rmo_j = NULL;
244:   } else SETERRQ(PetscObjectComm((PetscObject)mc),PETSC_ERR_ARG_WRONG,"Matrix must be AIJ for greedy coloring");
245:   if (!gr->symmetric) {
246:     MatTranspose(m, MAT_INITIAL_MATRIX, &mt);
247:     if (isSEQAIJ) {
248:       Mat_SeqAIJ *dseq = (Mat_SeqAIJ*) mt->data;
249:       rmd_i = dseq->i;
250:       rmd_j = dseq->j;
251:       rmo_i = NULL;
252:       rmo_j = NULL;
253:     } else SETERRQ(PetscObjectComm((PetscObject) mc), PETSC_ERR_SUP, "Nonsymmetric greedy coloring only works in serial");
254:   }
255:   /* create the vectors and communication structures if necessary */
256:   no=0;
257:   if (mo) {
258:     VecGetLocalSize(aij->lvec,&no);
259:     PetscSFCreate(PetscObjectComm((PetscObject)m),&sf);
260:     MatGetLayouts(m,&layout,NULL);
261:     PetscSFSetGraphLayout(sf,layout,no,NULL,PETSC_COPY_VALUES,aij->garray);
262:   }
263:   MatColoringGetMaxColors(mc,&maxcolors);
264:   masksize=n;
265:   nbad=0;
266:   badsize=n;
267:   PetscMalloc1(masksize,&mask);
268:   PetscMalloc4(n,&d1cols,n,&dcolors,n,&conf,n,&bad);
269:   PetscMalloc2(badsize,&badidx,badsize,&badnext);
270:   for(i=0;i<masksize;i++) {
271:     mask[i]=-1;
272:   }
273:   for (i=0;i<n;i++) {
274:     dcolors[i]=maxcolors;
275:     bad[i]=-1;
276:   }
277:   for (i=0;i<badsize;i++) {
278:     badnext[i]=-1;
279:   }
280:   if (mo) {
281:     PetscMalloc3(no,&owts,no,&oconf,no,&ocolors);
282:     PetscSFBcastBegin(sf,MPIU_REAL,wts,owts);
283:     PetscSFBcastEnd(sf,MPIU_REAL,wts,owts);
284:     for (i=0;i<no;i++) {
285:       ocolors[i]=maxcolors;
286:     }
287:   } else {                      /* Appease overzealous -Wmaybe-initialized */
288:     owts = NULL;
289:     oconf = NULL;
290:     ocolors = NULL;
291:   }
292:   mcol=0;
293:   while (nd_global < n_global) {
294:     nd=n;
295:     /* assign lowest possible color to each local vertex */
296:     mcol_global=0;
297:     PetscLogEventBegin(MATCOLORING_Local,mc,0,0,0);
298:     for (i=0;i<n;i++) {
299:       idx=lperm[i];
300:       if (dcolors[idx] == maxcolors) {
301:         /* entries in bad */
302:         cbad=bad[idx];
303:         while (cbad>=0) {
304:           ccol=badidx[cbad];
305:           if (ccol>=masksize) {
306:             PetscInt *newmask;
307:             PetscMalloc1(masksize*2,&newmask);
308:             for(k=0;k<2*masksize;k++) {
309:               newmask[k]=-1;
310:             }
311:             for(k=0;k<masksize;k++) {
312:               newmask[k]=mask[k];
313:             }
314:             PetscFree(mask);
315:             mask=newmask;
316:             masksize*=2;
317:           }
318:           mask[ccol]=idx;
319:           cbad=badnext[cbad];
320:         }
321:         /* diagonal distance-one rows */
322:         nd1cols=0;
323:         ncols = rmd_i[idx+1]-rmd_i[idx];
324:         cidx = &(rmd_j[rmd_i[idx]]);
325:         for (j=0;j<ncols;j++) {
326:           d1cols[nd1cols] = cidx[j];
327:           nd1cols++;
328:           ccol=dcolors[cidx[j]];
329:           if (ccol != maxcolors) {
330:             if (ccol>=masksize) {
331:               PetscInt *newmask;
332:               PetscMalloc1(masksize*2,&newmask);
333:               for(k=0;k<2*masksize;k++) {
334:                 newmask[k]=-1;
335:               }
336:               for(k=0;k<masksize;k++) {
337:                 newmask[k]=mask[k];
338:               }
339:               PetscFree(mask);
340:               mask=newmask;
341:               masksize*=2;
342:             }
343:             mask[ccol]=idx;
344:           }
345:         }
346:         /* off-diagonal distance-one rows */
347:         if (mo) {
348:           ncols = rmo_i[idx+1]-rmo_i[idx];
349:           cidx = &(rmo_j[rmo_i[idx]]);
350:           for (j=0;j<ncols;j++) {
351:             ccol=ocolors[cidx[j]];
352:             if (ccol != maxcolors) {
353:               if (ccol>=masksize) {
354:                 PetscInt *newmask;
355:                 PetscMalloc1(masksize*2,&newmask);
356:                 for(k=0;k<2*masksize;k++) {
357:                   newmask[k]=-1;
358:                 }
359:                 for(k=0;k<masksize;k++) {
360:                   newmask[k]=mask[k];
361:                 }
362:                 PetscFree(mask);
363:                 mask=newmask;
364:                 masksize*=2;
365:               }
366:               mask[ccol]=idx;
367:             }
368:           }
369:         }
370:         /* diagonal distance-two rows */
371:         for (j=0;j<nd1cols;j++) {
372:           ncols = md_i[d1cols[j]+1]-md_i[d1cols[j]];
373:           cidx = &(md_j[md_i[d1cols[j]]]);
374:           for (l=0;l<ncols;l++) {
375:             ccol=dcolors[cidx[l]];
376:             if (ccol != maxcolors) {
377:               if (ccol>=masksize) {
378:                 PetscInt *newmask;
379:                 PetscMalloc1(masksize*2,&newmask);
380:                 for(k=0;k<2*masksize;k++) {
381:                   newmask[k]=-1;
382:                 }
383:                 for(k=0;k<masksize;k++) {
384:                   newmask[k]=mask[k];
385:                 }
386:                 PetscFree(mask);
387:                 mask=newmask;
388:                 masksize*=2;
389:               }
390:               mask[ccol]=idx;
391:             }
392:           }
393:         }
394:         /* off-diagonal distance-two rows */
395:         if (mo) {
396:           for (j=0;j<nd1cols;j++) {
397:             ncols = mo_i[d1cols[j]+1]-mo_i[d1cols[j]];
398:             cidx = &(mo_j[mo_i[d1cols[j]]]);
399:             for (l=0;l<ncols;l++) {
400:               ccol=ocolors[cidx[l]];
401:               if (ccol != maxcolors) {
402:                 if (ccol>=masksize) {
403:                   PetscInt *newmask;
404:                   PetscMalloc1(masksize*2,&newmask);
405:                   for(k=0;k<2*masksize;k++) {
406:                     newmask[k]=-1;
407:                   }
408:                   for(k=0;k<masksize;k++) {
409:                     newmask[k]=mask[k];
410:                   }
411:                   PetscFree(mask);
412:                   mask=newmask;
413:                   masksize*=2;
414:                 }
415:                 mask[ccol]=idx;
416:               }
417:             }
418:           }
419:         }
420:         /* assign this one the lowest color possible by seeing if there's a gap in the sequence of sorted neighbor colors */
421:         for (j=0;j<masksize;j++) {
422:           if (mask[j]!=idx) {
423:             break;
424:           }
425:         }
426:         pcol=j;
427:         if (pcol>maxcolors) pcol=maxcolors;
428:         dcolors[idx]=pcol;
429:         if (pcol>mcol) mcol=pcol;
430:       }
431:     }
432:     PetscLogEventEnd(MATCOLORING_Local,mc,0,0,0);
433:     if (mo) {
434:       /* transfer neighbor colors */
435:       PetscSFBcastBegin(sf,MPIU_INT,dcolors,ocolors);
436:       PetscSFBcastEnd(sf,MPIU_INT,dcolors,ocolors);
437:       /* find the maximum color assigned locally and allocate a mask */
438:       MPIU_Allreduce(&mcol,&mcol_global,1,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)mc));
439:       PetscMalloc1(mcol_global+1,&colorweights);
440:       /* check for conflicts */
441:       for (i=0;i<n;i++) {
442:         conf[i]=PETSC_FALSE;
443:       }
444:       for (i=0;i<no;i++) {
445:         oconf[i]=PETSC_FALSE;
446:       }
447:       for (i=0;i<n;i++) {
448:         ncols = mo_i[i+1]-mo_i[i];
449:         cidx = &(mo_j[mo_i[i]]);
450:         if (ncols > 0) {
451:           /* fill in the mask */
452:           for (j=0;j<mcol_global+1;j++) {
453:             colorweights[j]=0;
454:           }
455:           colorweights[dcolors[i]]=wts[i];
456:           /* fill in the off-diagonal part of the mask */
457:           for (j=0;j<ncols;j++) {
458:             ccol=ocolors[cidx[j]];
459:             if (ccol < maxcolors) {
460:               if (colorweights[ccol] < owts[cidx[j]]) {
461:                 colorweights[ccol] = owts[cidx[j]];
462:               }
463:             }
464:           }
465:           /* fill in the on-diagonal part of the mask */
466:           ncols = md_i[i+1]-md_i[i];
467:           cidx = &(md_j[md_i[i]]);
468:           for (j=0;j<ncols;j++) {
469:             ccol=dcolors[cidx[j]];
470:             if (ccol < maxcolors) {
471:               if (colorweights[ccol] < wts[cidx[j]]) {
472:                 colorweights[ccol] = wts[cidx[j]];
473:               }
474:             }
475:           }
476:           /* go back through and set up on and off-diagonal conflict vectors */
477:           ncols = md_i[i+1]-md_i[i];
478:           cidx = &(md_j[md_i[i]]);
479:           for (j=0;j<ncols;j++) {
480:             ccol=dcolors[cidx[j]];
481:             if (ccol < maxcolors) {
482:               if (colorweights[ccol] > wts[cidx[j]]) {
483:                 conf[cidx[j]]=PETSC_TRUE;
484:               }
485:             }
486:           }
487:           ncols = mo_i[i+1]-mo_i[i];
488:           cidx = &(mo_j[mo_i[i]]);
489:           for (j=0;j<ncols;j++) {
490:             ccol=ocolors[cidx[j]];
491:             if (ccol < maxcolors) {
492:               if (colorweights[ccol] > owts[cidx[j]]) {
493:                 oconf[cidx[j]]=PETSC_TRUE;
494:               }
495:             }
496:           }
497:         }
498:       }
499:       nd_global=0;
500:       PetscFree(colorweights);
501:       PetscLogEventBegin(MATCOLORING_Comm,mc,0,0,0);
502:       PetscSFReduceBegin(sf,MPIU_INT,oconf,conf,MPIU_SUM);
503:       PetscSFReduceEnd(sf,MPIU_INT,oconf,conf,MPIU_SUM);
504:       PetscLogEventEnd(MATCOLORING_Comm,mc,0,0,0);
505:       /* go through and unset local colors that have conflicts */
506:       for (i=0;i<n;i++) {
507:         if (conf[i]>0) {
508:           /* push this color onto the bad stack */
509:           badidx[nbad]=dcolors[i];
510:           badnext[nbad]=bad[i];
511:           bad[i]=nbad;
512:           nbad++;
513:           if (nbad>=badsize) {
514:             PetscInt *newbadnext;
515:             ISColoringValue *newbadidx;
516:             PetscMalloc2(badsize*2,&newbadidx,badsize*2,&newbadnext);
517:             for(k=0;k<2*badsize;k++) {
518:               newbadnext[k]=-1;
519:             }
520:             for(k=0;k<badsize;k++) {
521:               newbadidx[k]=badidx[k];
522:               newbadnext[k]=badnext[k];
523:             }
524:             PetscFree2(badidx,badnext);
525:             badidx=newbadidx;
526:             badnext=newbadnext;
527:             badsize*=2;
528:           }
529:           dcolors[i] = maxcolors;
530:           nd--;
531:         }
532:       }
533:     }
534:     MPIU_Allreduce(&nd,&nd_global,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)mc));
535:   }
536:   if (mo) {
537:     PetscSFDestroy(&sf);
538:     PetscFree3(owts,oconf,ocolors);
539:   }
540:   for (i=0;i<n;i++) {
541:     colors[i]=dcolors[i];
542:   }
543:   PetscFree(mask);
544:   PetscFree4(d1cols,dcolors,conf,bad);
545:   PetscFree2(badidx,badnext);
546:   if (!gr->symmetric) {MatDestroy(&mt);}
547:   return(0);
548: }

550: static PetscErrorCode MatColoringApply_Greedy(MatColoring mc,ISColoring *iscoloring)
551: {
552:   PetscErrorCode  ierr;
553:   PetscInt        finalcolor,finalcolor_global;
554:   ISColoringValue *colors;
555:   PetscInt        ncolstotal,ncols;
556:   PetscReal       *wts;
557:   PetscInt        i,*lperm;

560:   MatGetSize(mc->mat,NULL,&ncolstotal);
561:   MatGetLocalSize(mc->mat,NULL,&ncols);
562:   if (!mc->user_weights) {
563:     MatColoringCreateWeights(mc,&wts,&lperm);
564:   } else {
565:     wts = mc->user_weights;
566:     lperm = mc->user_lperm;
567:   }
568:   PetscMalloc1(ncols,&colors);
569:   if (mc->dist == 1) {
570:     GreedyColoringLocalDistanceOne_Private(mc,wts,lperm,colors);
571:   } else if (mc->dist == 2) {
572:     GreedyColoringLocalDistanceTwo_Private(mc,wts,lperm,colors);
573:   } else SETERRQ(PetscObjectComm((PetscObject)mc),PETSC_ERR_ARG_OUTOFRANGE,"Only distance 1 and distance 2 supported by MatColoringGreedy");
574:   finalcolor=0;
575:   for (i=0;i<ncols;i++) {
576:     if (colors[i] > finalcolor) finalcolor=colors[i];
577:   }
578:   finalcolor_global=0;
579:   MPIU_Allreduce(&finalcolor,&finalcolor_global,1,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)mc));
580:   PetscLogEventBegin(MATCOLORING_ISCreate,mc,0,0,0);
581:   ISColoringCreate(PetscObjectComm((PetscObject)mc),finalcolor_global+1,ncols,colors,PETSC_OWN_POINTER,iscoloring);
582:   PetscLogEventEnd(MATCOLORING_ISCreate,mc,0,0,0);
583:   if (!mc->user_weights) {
584:     PetscFree(wts);
585:     PetscFree(lperm);
586:   }
587:   return(0);
588: }

590: static PetscErrorCode MatColoringSetFromOptions_Greedy(PetscOptionItems *PetscOptionsObject, MatColoring mc)
591: {
592:   MC_Greedy     *gr = (MC_Greedy *) mc->data;

596:   PetscOptionsHead(PetscOptionsObject, "Greedy options");
597:   PetscOptionsBool("-mat_coloring_greedy_symmetric", "Flag for assuming a symmetric matrix", "", gr->symmetric, &gr->symmetric, NULL);
598:   PetscOptionsTail();
599:   return(0);
600: }

602: /*MC
603:   MATCOLORINGGREEDY - Greedy-with-conflict correction based Matrix Coloring for distance 1 and 2.

605:    Level: beginner

607:    Notes:

609:    These algorithms proceed in two phases -- local coloring and conflict resolution.  The local coloring
610:    tentatively colors all vertices at the distance required given what's known of the global coloring.  Then,
611:    the updated colors are transferred to different processors at distance one.  In the distance one case, each
612:    vertex with nonlocal neighbors is then checked to see if it conforms, with the vertex being
613:    marked for recoloring if its lower weight than its same colored neighbor.  In the distance two case,
614:    each boundary vertex's immediate star is checked for validity of the coloring.  Lower-weight conflict
615:    vertices are marked, and then the conflicts are gathered back on owning processors.  In both cases
616:    this is done until each column has received a valid color.

618:    References:
619: .  1. - Bozdag et al. "A Parallel Distance 2 Graph Coloring Algorithm for Distributed Memory Computers"
620:    HPCC'05 Proceedings of the First international conference on High Performance Computing and Communications

622: .seealso: MatColoringCreate(), MatColoring, MatColoringSetType(), MatColoringType
623: M*/
624: PETSC_EXTERN PetscErrorCode MatColoringCreate_Greedy(MatColoring mc)
625: {
626:   MC_Greedy      *gr;

630:   PetscNewLog(mc,&gr);
631:   mc->data                = gr;
632:   mc->ops->apply          = MatColoringApply_Greedy;
633:   mc->ops->view           = NULL;
634:   mc->ops->destroy        = MatColoringDestroy_Greedy;
635:   mc->ops->setfromoptions = MatColoringSetFromOptions_Greedy;

637:   gr->symmetric = PETSC_TRUE;
638:   return(0);
639: }