Actual source code: matcoloring.c

petsc-3.10.5 2019-03-28
Report Typos and Errors
  1:  #include <petsc/private/matimpl.h>

  3: PetscFunctionList MatColoringList              = 0;
  4: PetscBool         MatColoringRegisterAllCalled = PETSC_FALSE;
  5: const char *const MatColoringWeightTypes[] = {"RANDOM","LEXICAL","LF","SL","MatColoringWeightType","MAT_COLORING_WEIGHT_",0};

  7: /*@C
  8:    MatColoringRegister - Adds a new sparse matrix coloring to the  matrix package.

 10:    Not Collective

 12:    Input Parameters:
 13: +  sname - name of Coloring (for example MATCOLORINGSL)
 14: -  function - function pointer that creates the coloring

 16:    Level: developer

 18:    Sample usage:
 19: .vb
 20:    MatColoringRegister("my_color",MyColor);
 21: .ve

 23:    Then, your partitioner can be chosen with the procedural interface via
 24: $     MatColoringSetType(part,"my_color")
 25:    or at runtime via the option
 26: $     -mat_coloring_type my_color

 28: .keywords: matrix, Coloring, register

 30: .seealso: MatColoringRegisterDestroy(), MatColoringRegisterAll()
 31: @*/
 32: PetscErrorCode  MatColoringRegister(const char sname[],PetscErrorCode (*function)(MatColoring))
 33: {

 37:   MatInitializePackage();
 38:   PetscFunctionListAdd(&MatColoringList,sname,function);
 39:   return(0);
 40: }

 42: /*@
 43:    MatColoringCreate - Creates a matrix coloring context.

 45:    Collective on MatColoring

 47:    Input Parameters:
 48: .  comm - MPI communicator

 50:    Output Parameter:
 51: .  mcptr - the new MatColoring context

 53:    Options Database Keys:
 54: +   -mat_coloring_type - the type of coloring algorithm used
 55: .   -mat_coloring_maxcolors - the maximum number of relevant colors, all nodes not in a color are in maxcolors+1
 56: .   -mat_coloring_distance - compute a distance 1,2,... coloring.
 57: .   -mat_coloring_view - print information about the coloring and the produced index sets
 58: .   -mat_coloring_test - debugging option that prints all coloring incompatibilities
 59: -   -mat_is_coloring_test - debugging option that throws an error if MatColoringApply() generates an incorrect iscoloring

 61:    Level: beginner

 63:    Notes:
 64:     A distance one coloring is useful, for example, multi-color SOR. A distance two coloring is for the finite difference computation of Jacobians 
 65:           (see MatFDColoringCreate()).

 67:           Some coloring types only support distance two colorings

 69: .keywords: Coloring, Matrix

 71: .seealso: MatColoring, MatColoringApply(), MatFDColoringCreate()
 72: @*/
 73: PetscErrorCode MatColoringCreate(Mat m,MatColoring *mcptr)
 74: {
 75:   MatColoring    mc;

 81:   *mcptr = NULL;

 83: #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
 84:   MatInitializePackage();
 85: #endif
 86:   PetscHeaderCreate(mc, MAT_COLORING_CLASSID,"MatColoring","Matrix coloring", "MatColoring",PetscObjectComm((PetscObject)m),MatColoringDestroy, MatColoringView);
 87:   PetscObjectReference((PetscObject)m);
 88:   mc->mat       = m;
 89:   mc->dist      = 2; /* default to Jacobian computation case */
 90:   mc->maxcolors = IS_COLORING_MAX;
 91:   *mcptr        = mc;
 92:   mc->valid     = PETSC_FALSE;
 93:   mc->weight_type = MAT_COLORING_WEIGHT_RANDOM;
 94:   mc->user_weights = NULL;
 95:   mc->user_lperm = NULL;
 96:   return(0);
 97: }


100: /*@
101:    MatColoringDestroy - Destroys the matrix coloring context

103:    Collective on MatColoring

105:    Input Parameter:
106: .  mc - the MatColoring context

108:    Level: beginner

110: .keywords: Coloring, destroy

112: .seealso: MatColoringCreate(), MatColoringApply()
113: @*/
114: PetscErrorCode MatColoringDestroy(MatColoring *mc)
115: {

119:   if (--((PetscObject)(*mc))->refct > 0) {*mc = 0; return(0);}
120:   MatDestroy(&(*mc)->mat);
121:   if ((*mc)->ops->destroy) {(*((*mc)->ops->destroy))(*mc);}
122:   if ((*mc)->user_weights) {PetscFree((*mc)->user_weights);}
123:   if ((*mc)->user_lperm) {PetscFree((*mc)->user_lperm);}
124:   PetscHeaderDestroy(mc);
125:   return(0);
126: }

128: /*@C
129:    MatColoringSetType - Sets the type of coloring algorithm used

131:    Collective on MatColoring

133:    Input Parameter:
134: +  mc - the MatColoring context
135: -  type - the type of coloring

137:    Level: beginner

139:    Notes:
140:     Possible types include the sequential types MATCOLORINGLF,
141:    MATCOLORINGSL, and MATCOLORINGID from the MINPACK package as well
142:    as a parallel MATCOLORINGMIS algorithm.

144: .keywords: Coloring, type

146: .seealso: MatColoringCreate(), MatColoringApply()
147: @*/
148: PetscErrorCode MatColoringSetType(MatColoring mc,MatColoringType type)
149: {
150:   PetscBool      match;
151:   PetscErrorCode ierr,(*r)(MatColoring);

156:   PetscObjectTypeCompare((PetscObject)mc,type,&match);
157:   if (match) return(0);
158:    PetscFunctionListFind(MatColoringList,type,&r);
159:   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested MatColoring type %s",type);
160:   if (mc->ops->destroy) {
161:     (*(mc)->ops->destroy)(mc);
162:     mc->ops->destroy = NULL;
163:   }
164:   mc->ops->apply            = 0;
165:   mc->ops->view             = 0;
166:   mc->ops->setfromoptions   = 0;
167:   mc->ops->destroy          = 0;

169:   PetscObjectChangeTypeName((PetscObject)mc,type);
170:   (*r)(mc);
171:   return(0);
172: }

174: /*@
175:    MatColoringSetFromOptions - Sets MatColoring options from user parameters

177:    Collective on MatColoring

179:    Input Parameters:
180: .  mc - MatColoring context

182:    Options Database Keys:
183: +   -mat_coloring_type - the type of coloring algorithm used
184: .   -mat_coloring_maxcolors - the maximum number of relevant colors, all nodes not in a color are in maxcolors+1
185: .   -mat_coloring_distance - compute a distance 1,2,... coloring.
186: .   -mat_coloring_view - print information about the coloring and the produced index sets

188:    Level: beginner

190: .keywords: Coloring, Matrix

192: .seealso: MatColoring, MatColoringApply()
193: @*/
194: PetscErrorCode MatColoringSetFromOptions(MatColoring mc)
195: {
196:   PetscBool      flg;
197:   MatColoringType deft = MATCOLORINGSL;
198:   char           type[256];
200:   PetscInt       dist,maxcolors;

204:   MatColoringGetDistance(mc,&dist);
205:   if (dist == 2) deft = MATCOLORINGSL;
206:   else           deft = MATCOLORINGGREEDY;
207:   MatColoringGetMaxColors(mc,&maxcolors);
208:   MatColoringRegisterAll();
209:   PetscObjectOptionsBegin((PetscObject)mc);
210:   if (((PetscObject)mc)->type_name) deft = ((PetscObject)mc)->type_name;
211:   PetscOptionsFList("-mat_coloring_type","The coloring method used","MatColoringSetType",MatColoringList,deft,type,256,&flg);
212:   if (flg) {
213:     MatColoringSetType(mc,type);
214:   } else if (!((PetscObject)mc)->type_name) {
215:     MatColoringSetType(mc,deft);
216:   }
217:   PetscOptionsInt("-mat_coloring_distance","Distance of the coloring","MatColoringSetDistance",dist,&dist,&flg);
218:   if (flg) {MatColoringSetDistance(mc,dist);}
219:   PetscOptionsInt("-mat_coloring_maxcolors","Maximum colors returned at the end. 1 returns an independent set","MatColoringSetMaxColors",maxcolors,&maxcolors,&flg);
220:   if (flg) {MatColoringSetMaxColors(mc,maxcolors);}
221:   if (mc->ops->setfromoptions) {
222:     (*mc->ops->setfromoptions)(PetscOptionsObject,mc);
223:   }
224:   PetscOptionsBool("-mat_coloring_test","Check that a valid coloring has been produced","",mc->valid,&mc->valid,NULL);
225:   PetscOptionsBool("-mat_is_coloring_test","Check that a valid iscoloring has been produced","",mc->valid_iscoloring,&mc->valid_iscoloring,NULL);
226:   PetscOptionsEnum("-mat_coloring_weight_type","Sets the type of vertex weighting used","MatColoringSetWeightType",MatColoringWeightTypes,(PetscEnum)mc->weight_type,(PetscEnum*)&mc->weight_type,NULL);
227:   PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)mc);
228:   PetscOptionsEnd();
229:   return(0);
230: }

232: /*@
233:    MatColoringSetDistance - Sets the distance of the coloring

235:    Logically Collective on MatColoring

237:    Input Parameter:
238: .  mc - the MatColoring context
239: .  dist - the distance the coloring should compute

241:    Level: beginner

243:    Notes:
244:     The distance of the coloring denotes the minimum number
245:    of edges in the graph induced by the matrix any two vertices
246:    of the same color are.  Distance-1 colorings are the classical
247:    coloring, where no two vertices of the same color are adjacent.
248:    distance-2 colorings are useful for the computation of Jacobians.

250: .keywords: Coloring, distance, Jacobian

252: .seealso: MatColoringGetDistance(), MatColoringApply()
253: @*/
254: PetscErrorCode MatColoringSetDistance(MatColoring mc,PetscInt dist)
255: {
258:   mc->dist = dist;
259:   return(0);
260: }

262: /*@
263:    MatColoringGetDistance - Gets the distance of the coloring

265:    Logically Collective on MatColoring

267:    Input Parameter:
268: .  mc - the MatColoring context

270:    Output Paramter:
271: .  dist - the current distance being used for the coloring.

273:    Level: beginner

275: .keywords: Coloring, distance

277: .seealso: MatColoringSetDistance(), MatColoringApply()
278: @*/
279: PetscErrorCode MatColoringGetDistance(MatColoring mc,PetscInt *dist)
280: {
283:   if (dist) *dist = mc->dist;
284:   return(0);
285: }

287: /*@
288:    MatColoringSetMaxColors - Sets the maximum number of colors

290:    Logically Collective on MatColoring

292:    Input Parameter:
293: +  mc - the MatColoring context
294: -  maxcolors - the maximum number of colors to produce

296:    Level: beginner

298:    Notes:
299:     This may be used to compute a certain number of
300:    independent sets from the graph.  For instance, while using
301:    MATCOLORINGMIS and maxcolors = 1, one gets out an MIS.  Vertices
302:    not in a color are set to have color maxcolors+1, which is not
303:    a valid color as they may be adjacent.

305: .keywords: Coloring

307: .seealso: MatColoringGetMaxColors(), MatColoringApply()
308: @*/
309: PetscErrorCode MatColoringSetMaxColors(MatColoring mc,PetscInt maxcolors)
310: {
313:   mc->maxcolors = maxcolors;
314:   return(0);
315: }

317: /*@
318:    MatColoringGetMaxColors - Gets the maximum number of colors

320:    Logically Collective on MatColoring

322:    Input Parameter:
323: .  mc - the MatColoring context

325:    Output Paramter:
326: .  maxcolors - the current maximum number of colors to produce

328:    Level: beginner

330: .keywords: Coloring

332: .seealso: MatColoringSetMaxColors(), MatColoringApply()
333: @*/
334: PetscErrorCode MatColoringGetMaxColors(MatColoring mc,PetscInt *maxcolors)
335: {
338:   if (maxcolors) *maxcolors = mc->maxcolors;
339:   return(0);
340: }

342: /*@
343:    MatColoringApply - Apply the coloring to the matrix, producing index
344:    sets corresponding to a number of independent sets in the induced
345:    graph.

347:    Collective on MatColoring

349:    Input Parameters:
350: .  mc - the MatColoring context

352:    Output Parameter:
353: .  coloring - the ISColoring instance containing the coloring

355:    Level: beginner

357: .keywords: Coloring, Apply

359: .seealso: MatColoring, MatColoringCreate()
360: @*/
361: PetscErrorCode MatColoringApply(MatColoring mc,ISColoring *coloring)
362: {
363:   PetscErrorCode    ierr;
364:   PetscBool         flg;
365:   PetscViewerFormat format;
366:   PetscViewer       viewer;
367:   PetscInt          nc,ncolors;

371:   PetscLogEventBegin(MATCOLORING_Apply,mc,0,0,0);
372:   (*mc->ops->apply)(mc,coloring);
373:   PetscLogEventEnd(MATCOLORING_Apply,mc,0,0,0);

375:   /* valid */
376:   if (mc->valid) {
377:     MatColoringTest(mc,*coloring);
378:   }
379:   if (mc->valid_iscoloring) {
380:     MatISColoringTest(mc->mat,*coloring);
381:   }

383:   /* view */
384:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)mc),((PetscObject)mc)->prefix,"-mat_coloring_view",&viewer,&format,&flg);
385:   if (flg && !PetscPreLoadingOn) {
386:     PetscViewerPushFormat(viewer,format);
387:     MatColoringView(mc,viewer);
388:     MatGetSize(mc->mat,NULL,&nc);
389:     ISColoringGetIS(*coloring,&ncolors,NULL);
390:     PetscViewerASCIIPrintf(viewer,"  Number of colors %d\n",ncolors);
391:     PetscViewerASCIIPrintf(viewer,"  Number of total columns %d\n",nc);
392:     if (nc <= 1000) {ISColoringView(*coloring,viewer);}
393:     PetscViewerPopFormat(viewer);
394:     PetscViewerDestroy(&viewer);
395:   }
396:   return(0);
397: }

399: /*@
400:    MatColoringView - Output details about the MatColoring.

402:    Collective on MatColoring

404:    Input Parameters:
405: -  mc - the MatColoring context
406: +  viewer - the Viewer context

408:    Level: beginner

410: .keywords: Coloring, view

412: .seealso: MatColoring, MatColoringApply()
413: @*/
414: PetscErrorCode MatColoringView(MatColoring mc,PetscViewer viewer)
415: {
417:   PetscBool      iascii;

421:   if (!viewer) {
422:     PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)mc),&viewer);
423:   }

427:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
428:   if (iascii) {
429:     PetscObjectPrintClassNamePrefixType((PetscObject)mc,viewer);
430:     PetscViewerASCIIPrintf(viewer,"  Weight type: %s\n",MatColoringWeightTypes[mc->weight_type]);
431:     if (mc->maxcolors > 0) {
432:       PetscViewerASCIIPrintf(viewer,"  Distance %D, Max. Colors %D\n",mc->dist,mc->maxcolors);
433:     } else {
434:       PetscViewerASCIIPrintf(viewer,"  Distance %d\n",mc->dist);
435:     }
436:   }
437:   return(0);
438: }

440: /*@
441:    MatColoringSetWeightType - Set the type of weight computation used.

443:    Logically collective on MatColoring

445:    Input Parameters:
446: -  mc - the MatColoring context
447: +  wt - the weight type

449:    Level: beginner

451: .keywords: Coloring, view

453: .seealso: MatColoring, MatColoringWeightType
454: @*/
455: PetscErrorCode MatColoringSetWeightType(MatColoring mc,MatColoringWeightType wt)
456: {
458:   mc->weight_type = wt;
459:   return(0);

461: }