Actual source code: color.c

  1: /*
  2:      Routines that call the kernel minpack coloring subroutines
  3: */

 5:  #include src/mat/matimpl.h
 6:  #include src/mat/color/color.h

  8: /*
  9:     MatFDColoringDegreeSequence_Minpack - Calls the MINPACK routine seqr() that
 10:       computes the degree sequence required by MINPACK coloring routines.
 11: */
 14: PetscErrorCode MatFDColoringDegreeSequence_Minpack(PetscInt m,PetscInt *cja, PetscInt *cia, PetscInt *rja, PetscInt *ria, PetscInt **seq)
 15: {
 16:   PetscInt       *work;

 20:   PetscMalloc(m*sizeof(PetscInt),&work);
 21:   PetscMalloc(m*sizeof(PetscInt),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: */
 36: PetscErrorCode MatFDColoringMinimumNumberofColors_Private(PetscInt m,PetscInt *ia,PetscInt *minc)
 37: {
 38:   PetscInt i,c = 0;

 41:   for (i=0; i<m; i++) {
 42:     c = PetscMax(c,ia[i+1]-ia[i]);
 43:   }
 44:   *minc = c;
 45:   return(0);
 46: }

 49: /* ----------------------------------------------------------------------------*/
 50: /*
 51:     MatFDColoringSL_Minpack - Uses the smallest-last (SL) coloring of minpack
 52: */
 55: PetscErrorCode MatFDColoringSL_Minpack(Mat mat,MatColoringType name,ISColoring *iscoloring)
 56: {
 58:   PetscInt        *list,*work,clique,*ria,*rja,*cia,*cja,*seq,*coloring,n;
 59:   PetscInt        ncolors,i;
 60:   PetscTruth      done;

 63:   MatGetRowIJ(mat,1,PETSC_FALSE,&n,&ria,&rja,&done);
 64:   MatGetColumnIJ(mat,1,PETSC_FALSE,&n,&cia,&cja,&done);
 65:   if (!done) SETERRQ(PETSC_ERR_SUP,"Ordering requires IJ");

 67:   MatFDColoringDegreeSequence_Minpack(n,cja,cia,rja,ria,&seq);

 69:   PetscMalloc(5*n*sizeof(PetscInt),&list);
 70:   work = list + n;

 72:   MINPACKslo(&n,cja,cia,rja,ria,seq,list,&clique,work,work+n,work+2*n,work+3*n);

 74:   PetscMalloc(n*sizeof(PetscInt),&coloring);
 75:   MINPACKseq(&n,cja,cia,rja,ria,list,coloring,&ncolors,work);

 77:   PetscFree(list);
 78:   PetscFree(seq);
 79:   MatRestoreRowIJ(mat,1,PETSC_FALSE,&n,&ria,&rja,&done);
 80:   MatRestoreColumnIJ(mat,1,PETSC_FALSE,&n,&cia,&cja,&done);

 82:   /* shift coloring numbers to start at zero and shorten */
 83:   if (ncolors > IS_COLORING_MAX-1) SETERRQ(PETSC_ERR_SUP,"Maximum color size exceeded");
 84:   {
 85:     ISColoringValue *s = (ISColoringValue*) coloring;
 86:     for (i=0; i<n; i++) {
 87:       s[i] = (ISColoringValue) (coloring[i]-1);
 88:     }
 89:     MatColoringPatch(mat,n,ncolors,s,iscoloring);
 90:   }
 91:   return(0);
 92: }

 96: /* ----------------------------------------------------------------------------*/
 97: /*
 98:     MatFDColoringLF_Minpack - 
 99: */
102: PetscErrorCode MatFDColoringLF_Minpack(Mat mat,MatColoringType name,ISColoring *iscoloring)
103: {
105:   PetscInt       *list,*work,*ria,*rja,*cia,*cja,*seq,*coloring,n;
106:   PetscInt       n1, none,ncolors,i;
107:   PetscTruth     done;

110:   MatGetRowIJ(mat,1,PETSC_FALSE,&n,&ria,&rja,&done);
111:   MatGetColumnIJ(mat,1,PETSC_FALSE,&n,&cia,&cja,&done);
112:   if (!done) SETERRQ(PETSC_ERR_SUP,"Ordering requires IJ");

114:   MatFDColoringDegreeSequence_Minpack(n,cja,cia,rja,ria,&seq);

116:   PetscMalloc(5*n*sizeof(PetscInt),&list);
117:   work = list + n;

119:   n1   = n - 1;
120:   none = -1;
121:   MINPACKnumsrt(&n,&n1,seq,&none,list,work+2*n,work+n);
122:   PetscMalloc(n*sizeof(PetscInt),&coloring);
123:   MINPACKseq(&n,cja,cia,rja,ria,list,coloring,&ncolors,work);

125:   PetscFree(list);
126:   PetscFree(seq);

128:   MatRestoreRowIJ(mat,1,PETSC_FALSE,&n,&ria,&rja,&done);
129:   MatRestoreColumnIJ(mat,1,PETSC_FALSE,&n,&cia,&cja,&done);

131:   /* shift coloring numbers to start at zero and shorten */
132:   if (ncolors > IS_COLORING_MAX-1) SETERRQ(PETSC_ERR_SUP,"Maximum color size exceeded");
133:   {
134:     ISColoringValue *s = (ISColoringValue*) coloring;
135:     for (i=0; i<n; i++) {
136:       s[i] = (ISColoringValue) (coloring[i]-1);
137:     }
138:     MatColoringPatch(mat,n,ncolors,s,iscoloring);
139:   }
140:   return(0);
141: }

145: /* ----------------------------------------------------------------------------*/
146: /*
147:     MatFDColoringID_Minpack - 
148: */
151: PetscErrorCode MatFDColoringID_Minpack(Mat mat,MatColoringType name,ISColoring *iscoloring)
152: {
154:   PetscInt       *list,*work,clique,*ria,*rja,*cia,*cja,*seq,*coloring,n;
155:   PetscInt       ncolors,i;
156:   PetscTruth     done;

159:   MatGetRowIJ(mat,1,PETSC_FALSE,&n,&ria,&rja,&done);
160:   MatGetColumnIJ(mat,1,PETSC_FALSE,&n,&cia,&cja,&done);
161:   if (!done) SETERRQ(PETSC_ERR_SUP,"Ordering requires IJ");

163:   MatFDColoringDegreeSequence_Minpack(n,cja,cia,rja,ria,&seq);

165:   PetscMalloc(5*n*sizeof(PetscInt),&list);
166:   work = list + n;

168:   MINPACKido(&n,&n,cja,cia,rja,ria,seq,list,&clique,work,work+n,work+2*n,work+3*n);

170:   PetscMalloc(n*sizeof(PetscInt),&coloring);
171:   MINPACKseq(&n,cja,cia,rja,ria,list,coloring,&ncolors,work);

173:   PetscFree(list);
174:   PetscFree(seq);

176:   MatRestoreRowIJ(mat,1,PETSC_FALSE,&n,&ria,&rja,&done);
177:   MatRestoreColumnIJ(mat,1,PETSC_FALSE,&n,&cia,&cja,&done);

179:   /* shift coloring numbers to start at zero and shorten */
180:   if (ncolors > IS_COLORING_MAX-1) SETERRQ(PETSC_ERR_SUP,"Maximum color size exceeded");
181:   {
182:     ISColoringValue *s = (ISColoringValue*) coloring;
183:     for (i=0; i<n; i++) {
184:       s[i] = (ISColoringValue) (coloring[i]-1);
185:     }
186:     MatColoringPatch(mat,n,ncolors,s,iscoloring);
187:   }
188:   return(0);
189: }

193: /*
194:    Simplest coloring, each column of the matrix gets its own unique color.
195: */
198: PetscErrorCode MatColoring_Natural(Mat mat,const MatColoringType color, ISColoring *iscoloring)
199: {
200:   PetscErrorCode  ierr;
201:   PetscInt        start,end,i;
202:   ISColoringValue *colors;
203:   MPI_Comm        comm;

206:   MatGetOwnershipRange(mat,&start,&end);
207:   PetscObjectGetComm((PetscObject)mat,&comm);
208:   PetscMalloc((end-start+1)*sizeof(PetscInt),&colors);
209:   for (i=start; i<end; i++) {
210:     colors[i-start] = i;
211:   }
212:   ISColoringCreate(comm,end-start,colors,iscoloring);

214:   return(0);
215: }
217: 
218: /* ===========================================================================================*/

220:  #include petscsys.h

222: PetscFList MatColoringList = 0;
223: PetscTruth MatColoringRegisterAllCalled = PETSC_FALSE;

227: PetscErrorCode MatColoringRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(Mat,const MatColoringType,ISColoring*))
228: {
230:   char           fullname[PETSC_MAX_PATH_LEN];

233:   PetscFListConcat(path,name,fullname);
234:   PetscFListAdd(&MatColoringList,sname,fullname,(void (*)(void))function);
235:   return(0);
236: }

240: /*@C
241:    MatColoringRegisterDestroy - Frees the list of coloringing routines.

243:    Not Collective

245:    Level: developer

247: .keywords: matrix, register, destroy

249: .seealso: MatColoringRegisterDynamic(), MatColoringRegisterAll()
250: @*/
251: PetscErrorCode MatColoringRegisterDestroy(void)
252: {

256:   if (MatColoringList) {
257:     PetscFListDestroy(&MatColoringList);
258:     MatColoringList = 0;
259:   }
260:   return(0);
261: }

263: EXTERN PetscErrorCode MatAdjustForInodes(Mat,IS *,IS *);

267: /*@C
268:    MatGetColoring - Gets a coloring for a matrix to reduce the number of function evaluations
269:    needed to compute a sparse Jacobian via differencing.

271:    Collective on Mat

273:    Input Parameters:
274: .  mat - the matrix
275: .  type - type of coloring, one of the following:
276: $      MATCOLORING_NATURAL - natural (one color for each column, very slow)
277: $      MATCOLORING_SL - smallest-last
278: $      MATCOLORING_LF - largest-first
279: $      MATCOLORING_ID - incidence-degree

281:    Output Parameters:
282: .   iscoloring - the coloring

284:    Options Database Keys:
285:    To specify the coloring through the options database, use one of
286:    the following 
287: $    -mat_coloring_type natural, -mat_coloring_type sl, -mat_coloring_type lf,
288: $    -mat_coloring_type id
289:    To see the coloring use
290: $    -mat_coloring_view

292:    Level: intermediate

294:    Notes:
295:      These compute the graph coloring of the graph of A^{T}A. The coloring used 
296:    for efficient (parallel or thread based) triangular solves etc is NOT yet 
297:    available. 

299:    The user can define additional colorings; see MatColoringRegisterDynamic().

301:    The sequential colorings SL, LF, and ID are obtained via the Minpack software that was
302:    converted to C using f2c.

304: .keywords: matrix, get, coloring

306: .seealso:  MatGetColoringTypeFromOptions(), MatColoringRegisterDynamic(), MatFDColoringCreate(),
307:            SNESDefaultComputeJacobianColor()
308: @*/
309: PetscErrorCode MatGetColoring(Mat mat,const MatColoringType type,ISColoring *iscoloring)
310: {
311:   PetscTruth     flag;
312:   PetscErrorCode ierr,(*r)(Mat,const MatColoringType,ISColoring *);
313:   char           tname[PETSC_MAX_PATH_LEN];

318:   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
319:   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
320:   if (!MatColoringRegisterAllCalled) {
321:     MatColoringRegisterAll(PETSC_NULL);
322:   }
323: 
324:   /* look for type on command line */
325:   PetscOptionsGetString(mat->prefix,"-mat_coloring_type",tname,256,&flag);
326:   if (flag) {
327:     type = tname;
328:   }

330:   PetscLogEventBegin(MAT_GetColoring,mat,0,0,0);
331:    PetscFListFind(mat->comm, MatColoringList, type,(void (**)(void)) &r);
332:   if (!r) {SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Unknown or unregistered type: %s",type);}
333:   (*r)(mat,type,iscoloring);
334:   PetscLogEventEnd(MAT_GetColoring,mat,0,0,0);

336:   PetscLogInfo((PetscObject)mat,"MatGetColoring:Number of colors %d\n",(*iscoloring)->n);
337:   PetscOptionsHasName(PETSC_NULL,"-mat_coloring_view",&flag);
338:   if (flag) {
339:     ISColoringView(*iscoloring,PETSC_VIEWER_STDOUT_((*iscoloring)->comm));
340:   }
341:   return(0);
342: }
343: