Actual source code: greedy.c

petsc-3.7.7 2017-09-25
Report Typos and Errors
  1: #include <petsc/private/matimpl.h>      /*I "petscmat.h"  I*/
  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;

 12: static PetscErrorCode MatColoringDestroy_Greedy(MatColoring mc)
 13: {

 17:   PetscFree(mc->data);
 18:   return(0);
 19: }

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

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

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

195: static PetscErrorCode GreedyColoringLocalDistanceTwo_Private(MatColoring mc,PetscReal *wts,PetscInt *lperm,ISColoringValue *colors)
196: {
197:   MC_Greedy       *gr = (MC_Greedy *) mc->data;
198:   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;
199:   Mat             m = mc->mat, mt;
200:   Mat_MPIAIJ      *aij = (Mat_MPIAIJ*)m->data;
201:   Mat             md=NULL,mo=NULL;
202:   const PetscInt  *md_i,*mo_i,*md_j,*mo_j;
203:   const PetscInt  *rmd_i,*rmo_i,*rmd_j,*rmo_j;
204:   PetscBool       isMPIAIJ,isSEQAIJ;
205:   PetscInt        pcol,*dcolors,*ocolors;
206:   ISColoringValue *badidx;
207:   const PetscInt  *cidx;
208:   PetscReal       *owts,*colorweights;
209:   PetscInt        *oconf,*conf;
210:   PetscSF         sf;
211:   PetscLayout     layout;
212:   PetscErrorCode  ierr;

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

558: static PetscErrorCode MatColoringApply_Greedy(MatColoring mc,ISColoring *iscoloring)
559: {
560:   PetscErrorCode  ierr;
561:   PetscInt        finalcolor,finalcolor_global;
562:   ISColoringValue *colors;
563:   PetscInt        ncolstotal,ncols;
564:   PetscReal       *wts;
565:   PetscInt        i,*lperm;

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

600: static PetscErrorCode MatColoringSetFromOptions_Greedy(PetscOptionItems *PetscOptionsObject, MatColoring mc)
601: {
602:   MC_Greedy     *gr = (MC_Greedy *) mc->data;

606:   PetscOptionsHead(PetscOptionsObject, "Greedy options");
607:   PetscOptionsBool("-mat_coloring_greedy_symmetric", "Flag for assuming a symmetric matrix", "", gr->symmetric, &gr->symmetric, NULL);
608:   PetscOptionsTail();
609:   return(0);
610: }

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

617:    Level: beginner

619:    Notes:

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

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

634: .seealso: MatColoringCreate(), MatColoring, MatColoringSetType()
635: M*/
636: PETSC_EXTERN PetscErrorCode MatColoringCreate_Greedy(MatColoring mc)
637: {
638:   MC_Greedy      *gr;

642:   PetscNewLog(mc,&gr);
643:   mc->data                = gr;
644:   mc->ops->apply          = MatColoringApply_Greedy;
645:   mc->ops->view           = NULL;
646:   mc->ops->destroy        = MatColoringDestroy_Greedy;
647:   mc->ops->setfromoptions = MatColoringSetFromOptions_Greedy;

649:   gr->symmetric = PETSC_TRUE;
650:   return(0);
651: }