Actual source code: fdmatrix.c

petsc-master 2019-08-22
Report Typos and Errors

  2: /*
  3:    This is where the abstract matrix operations are defined that are
  4:   used for finite difference computations of Jacobians using coloring.
  5: */

  7:  #include <petsc/private/matimpl.h>
  8:  #include <petsc/private/isimpl.h>

 10: PetscErrorCode  MatFDColoringSetF(MatFDColoring fd,Vec F)
 11: {

 15:   if (F) {
 16:     VecCopy(F,fd->w1);
 17:     fd->fset = PETSC_TRUE;
 18:   } else {
 19:     fd->fset = PETSC_FALSE;
 20:   }
 21:   return(0);
 22: }

 24:  #include <petscdraw.h>
 25: static PetscErrorCode MatFDColoringView_Draw_Zoom(PetscDraw draw,void *Aa)
 26: {
 27:   MatFDColoring  fd = (MatFDColoring)Aa;
 29:   PetscInt       i,j,nz,row;
 30:   PetscReal      x,y;
 31:   MatEntry       *Jentry=fd->matentry;

 34:   /* loop over colors  */
 35:   nz = 0;
 36:   for (i=0; i<fd->ncolors; i++) {
 37:     for (j=0; j<fd->nrows[i]; j++) {
 38:       row = Jentry[nz].row;
 39:       y   = fd->M - row - fd->rstart;
 40:       x   = (PetscReal)Jentry[nz++].col;
 41:       PetscDrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);
 42:     }
 43:   }
 44:   return(0);
 45: }

 47: static PetscErrorCode MatFDColoringView_Draw(MatFDColoring fd,PetscViewer viewer)
 48: {
 50:   PetscBool      isnull;
 51:   PetscDraw      draw;
 52:   PetscReal      xr,yr,xl,yl,h,w;

 55:   PetscViewerDrawGetDraw(viewer,0,&draw);
 56:   PetscDrawIsNull(draw,&isnull);
 57:   if (isnull) return(0);

 59:   xr   = fd->N; yr  = fd->M; h = yr/10.0; w = xr/10.0;
 60:   xr  += w;     yr += h;    xl = -w;     yl = -h;
 61:   PetscDrawSetCoordinates(draw,xl,yl,xr,yr);
 62:   PetscObjectCompose((PetscObject)fd,"Zoomviewer",(PetscObject)viewer);
 63:   PetscDrawZoom(draw,MatFDColoringView_Draw_Zoom,fd);
 64:   PetscObjectCompose((PetscObject)fd,"Zoomviewer",NULL);
 65:   PetscDrawSave(draw);
 66:   return(0);
 67: }

 69: /*@C
 70:    MatFDColoringView - Views a finite difference coloring context.

 72:    Collective on MatFDColoring

 74:    Input  Parameters:
 75: +  c - the coloring context
 76: -  viewer - visualization context

 78:    Level: intermediate

 80:    Notes:
 81:    The available visualization contexts include
 82: +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
 83: .     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
 84:         output where only the first processor opens
 85:         the file.  All other processors send their
 86:         data to the first processor to print.
 87: -     PETSC_VIEWER_DRAW_WORLD - graphical display of nonzero structure

 89:    Notes:
 90:      Since PETSc uses only a small number of basic colors (currently 33), if the coloring
 91:    involves more than 33 then some seemingly identical colors are displayed making it look
 92:    like an illegal coloring. This is just a graphical artifact.

 94: .seealso: MatFDColoringCreate()

 96: @*/
 97: PetscErrorCode  MatFDColoringView(MatFDColoring c,PetscViewer viewer)
 98: {
 99:   PetscErrorCode    ierr;
100:   PetscInt          i,j;
101:   PetscBool         isdraw,iascii;
102:   PetscViewerFormat format;

106:   if (!viewer) {
107:     PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)c),&viewer);
108:   }

112:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
113:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
114:   if (isdraw) {
115:     MatFDColoringView_Draw(c,viewer);
116:   } else if (iascii) {
117:     PetscObjectPrintClassNamePrefixType((PetscObject)c,viewer);
118:     PetscViewerASCIIPrintf(viewer,"  Error tolerance=%g\n",(double)c->error_rel);
119:     PetscViewerASCIIPrintf(viewer,"  Umin=%g\n",(double)c->umin);
120:     PetscViewerASCIIPrintf(viewer,"  Number of colors=%D\n",c->ncolors);

122:     PetscViewerGetFormat(viewer,&format);
123:     if (format != PETSC_VIEWER_ASCII_INFO) {
124:       PetscInt row,col,nz;
125:       nz = 0;
126:       for (i=0; i<c->ncolors; i++) {
127:         PetscViewerASCIIPrintf(viewer,"  Information for color %D\n",i);
128:         PetscViewerASCIIPrintf(viewer,"    Number of columns %D\n",c->ncolumns[i]);
129:         for (j=0; j<c->ncolumns[i]; j++) {
130:           PetscViewerASCIIPrintf(viewer,"      %D\n",c->columns[i][j]);
131:         }
132:         PetscViewerASCIIPrintf(viewer,"    Number of rows %D\n",c->nrows[i]);
133:         if (c->matentry) {
134:           for (j=0; j<c->nrows[i]; j++) {
135:             row  = c->matentry[nz].row;
136:             col  = c->matentry[nz++].col;
137:             PetscViewerASCIIPrintf(viewer,"      %D %D \n",row,col);
138:           }
139:         }
140:       }
141:     }
142:     PetscViewerFlush(viewer);
143:   }
144:   return(0);
145: }

147: /*@
148:    MatFDColoringSetParameters - Sets the parameters for the sparse approximation of
149:    a Jacobian matrix using finite differences.

151:    Logically Collective on MatFDColoring

153:    The Jacobian is estimated with the differencing approximation
154: .vb
155:        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
156:        htype = 'ds':
157:          h = error_rel*u[i]                 if  abs(u[i]) > umin
158:            = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
159:          dx_{i} = (0, ... 1, .... 0)

161:        htype = 'wp':
162:          h = error_rel * sqrt(1 + ||u||)
163: .ve

165:    Input Parameters:
166: +  coloring - the coloring context
167: .  error_rel - relative error
168: -  umin - minimum allowable u-value magnitude

170:    Level: advanced

172: .seealso: MatFDColoringCreate(), MatFDColoringSetFromOptions()

174: @*/
175: PetscErrorCode MatFDColoringSetParameters(MatFDColoring matfd,PetscReal error,PetscReal umin)
176: {
181:   if (error != PETSC_DEFAULT) matfd->error_rel = error;
182:   if (umin != PETSC_DEFAULT)  matfd->umin      = umin;
183:   return(0);
184: }

186: /*@
187:    MatFDColoringSetBlockSize - Sets block size for efficient inserting entries of Jacobian matrix.

189:    Logically Collective on MatFDColoring

191:    Input Parameters:
192: +  coloring - the coloring context
193: .  brows - number of rows in the block
194: -  bcols - number of columns in the block

196:    Level: intermediate

198: .seealso: MatFDColoringCreate(), MatFDColoringSetFromOptions()

200: @*/
201: PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring matfd,PetscInt brows,PetscInt bcols)
202: {
207:   if (brows != PETSC_DEFAULT) matfd->brows = brows;
208:   if (bcols != PETSC_DEFAULT) matfd->bcols = bcols;
209:   return(0);
210: }

212: /*@
213:    MatFDColoringSetUp - Sets up the internal data structures of matrix coloring context for the later use.

215:    Collective on Mat

217:    Input Parameters:
218: +  mat - the matrix containing the nonzero structure of the Jacobian
219: .  iscoloring - the coloring of the matrix; usually obtained with MatGetColoring() or DMCreateColoring()
220: -  color - the matrix coloring context

222:    Level: beginner

224: .seealso: MatFDColoringCreate(), MatFDColoringDestroy()
225: @*/
226: PetscErrorCode MatFDColoringSetUp(Mat mat,ISColoring iscoloring,MatFDColoring color)
227: {

233:   if (color->setupcalled) return(0);

235:   PetscLogEventBegin(MAT_FDColoringSetUp,mat,0,0,0);
236:   if (mat->ops->fdcoloringsetup) {
237:     (*mat->ops->fdcoloringsetup)(mat,iscoloring,color);
238:   } else SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Code not yet written for matrix type %s",((PetscObject)mat)->type_name);

240:   color->setupcalled = PETSC_TRUE;
241:    PetscLogEventEnd(MAT_FDColoringSetUp,mat,0,0,0);
242:   return(0);
243: }

245: /*@C
246:    MatFDColoringGetFunction - Gets the function to use for computing the Jacobian.

248:    Not Collective

250:    Input Parameters:
251: .  coloring - the coloring context

253:    Output Parameters:
254: +  f - the function
255: -  fctx - the optional user-defined function context

257:    Level: intermediate

259: .seealso: MatFDColoringCreate(), MatFDColoringSetFunction(), MatFDColoringSetFromOptions()

261: @*/
262: PetscErrorCode  MatFDColoringGetFunction(MatFDColoring matfd,PetscErrorCode (**f)(void),void **fctx)
263: {
266:   if (f) *f = matfd->f;
267:   if (fctx) *fctx = matfd->fctx;
268:   return(0);
269: }

271: /*@C
272:    MatFDColoringSetFunction - Sets the function to use for computing the Jacobian.

274:    Logically Collective on MatFDColoring

276:    Input Parameters:
277: +  coloring - the coloring context
278: .  f - the function
279: -  fctx - the optional user-defined function context

281:    Calling sequence of (*f) function:
282:     For SNES:    PetscErrorCode (*f)(SNES,Vec,Vec,void*)
283:     If not using SNES: PetscErrorCode (*f)(void *dummy,Vec,Vec,void*) and dummy is ignored

285:    Level: advanced

287:    Notes:
288:     This function is usually used automatically by SNES (when one uses SNESSetJacobian() with the argument
289:      SNESComputeJacobianDefaultColor()) and only needs to be used by someone computing a matrix via coloring directly by
290:      calling MatFDColoringApply()

292:    Fortran Notes:
293:     In Fortran you must call MatFDColoringSetFunction() for a coloring object to
294:   be used without SNES or within the SNES solvers.

296: .seealso: MatFDColoringCreate(), MatFDColoringGetFunction(), MatFDColoringSetFromOptions()

298: @*/
299: PetscErrorCode  MatFDColoringSetFunction(MatFDColoring matfd,PetscErrorCode (*f)(void),void *fctx)
300: {
303:   matfd->f    = f;
304:   matfd->fctx = fctx;
305:   return(0);
306: }

308: /*@
309:    MatFDColoringSetFromOptions - Sets coloring finite difference parameters from
310:    the options database.

312:    Collective on MatFDColoring

314:    The Jacobian, F'(u), is estimated with the differencing approximation
315: .vb
316:        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
317:        h = error_rel*u[i]                 if  abs(u[i]) > umin
318:          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
319:        dx_{i} = (0, ... 1, .... 0)
320: .ve

322:    Input Parameter:
323: .  coloring - the coloring context

325:    Options Database Keys:
326: +  -mat_fd_coloring_err <err> - Sets <err> (square root of relative error in the function)
327: .  -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude
328: .  -mat_fd_type - "wp" or "ds" (see MATMFFD_WP or MATMFFD_DS)
329: .  -mat_fd_coloring_view - Activates basic viewing
330: .  -mat_fd_coloring_view ::ascii_info - Activates viewing info
331: -  -mat_fd_coloring_view draw - Activates drawing

333:     Level: intermediate

335: .seealso: MatFDColoringCreate(), MatFDColoringView(), MatFDColoringSetParameters()

337: @*/
338: PetscErrorCode  MatFDColoringSetFromOptions(MatFDColoring matfd)
339: {
341:   PetscBool      flg;
342:   char           value[3];


347:   PetscObjectOptionsBegin((PetscObject)matfd);
348:   PetscOptionsReal("-mat_fd_coloring_err","Square root of relative error in function","MatFDColoringSetParameters",matfd->error_rel,&matfd->error_rel,0);
349:   PetscOptionsReal("-mat_fd_coloring_umin","Minimum allowable u magnitude","MatFDColoringSetParameters",matfd->umin,&matfd->umin,0);
350:   PetscOptionsString("-mat_fd_type","Algorithm to compute h, wp or ds","MatFDColoringCreate",matfd->htype,value,3,&flg);
351:   if (flg) {
352:     if (value[0] == 'w' && value[1] == 'p') matfd->htype = "wp";
353:     else if (value[0] == 'd' && value[1] == 's') matfd->htype = "ds";
354:     else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Unknown finite differencing type %s",value);
355:   }
356:   PetscOptionsInt("-mat_fd_coloring_brows","Number of block rows","MatFDColoringSetBlockSize",matfd->brows,&matfd->brows,NULL);
357:   PetscOptionsInt("-mat_fd_coloring_bcols","Number of block columns","MatFDColoringSetBlockSize",matfd->bcols,&matfd->bcols,&flg);
358:   if (flg && matfd->bcols > matfd->ncolors) {
359:     /* input bcols cannot be > matfd->ncolors, thus set it as ncolors */
360:     matfd->bcols = matfd->ncolors;
361:   }

363:   /* process any options handlers added with PetscObjectAddOptionsHandler() */
364:   PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)matfd);
365:   PetscOptionsEnd();
366:   return(0);
367: }

369: /*@C
370:    MatFDColoringSetType - Sets the approach for computing the finite difference parameter

372:    Collective on MatFDColoring

374:    Input Parameters:
375: +  coloring - the coloring context
376: -  type - either MATMFFD_WP or MATMFFD_DS

378:    Options Database Keys:
379: .  -mat_fd_type - "wp" or "ds"

381:    Note: It is goofy that the argument type is MatMFFDType since the MatFDColoring actually computes the matrix entries
382:          but the process of computing the entries is the same as as with the MatMFFD operation so we should reuse the names instead of
383:          introducing another one.

385:    Level: intermediate

387: .seealso: MatFDColoringCreate(), MatFDColoringView(), MatFDColoringSetParameters()

389: @*/
390: PetscErrorCode  MatFDColoringSetType(MatFDColoring matfd,MatMFFDType type)
391: {
394:   /*
395:      It is goofy to handle the strings this way but currently there is no code to free a dynamically created matfd->htype
396:      and this function is being provided as patch for a release so it shouldn't change the implementaton
397:   */
398:   if (type[0] == 'w' && type[1] == 'p') matfd->htype = "wp";
399:   else if (type[0] == 'd' && type[1] == 's') matfd->htype = "ds";
400:   else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Unknown finite differencing type %s",type);
401:   return(0);
402: }

404: PetscErrorCode MatFDColoringViewFromOptions(MatFDColoring fd,const char prefix[],const char optionname[])
405: {
406:   PetscErrorCode    ierr;
407:   PetscBool         flg;
408:   PetscViewer       viewer;
409:   PetscViewerFormat format;

412:   if (prefix) {
413:     PetscOptionsGetViewer(PetscObjectComm((PetscObject)fd),((PetscObject)fd)->options,prefix,optionname,&viewer,&format,&flg);
414:   } else {
415:     PetscOptionsGetViewer(PetscObjectComm((PetscObject)fd),((PetscObject)fd)->options,((PetscObject)fd)->prefix,optionname,&viewer,&format,&flg);
416:   }
417:   if (flg) {
418:     PetscViewerPushFormat(viewer,format);
419:     MatFDColoringView(fd,viewer);
420:     PetscViewerPopFormat(viewer);
421:     PetscViewerDestroy(&viewer);
422:   }
423:   return(0);
424: }

426: /*@
427:    MatFDColoringCreate - Creates a matrix coloring context for finite difference
428:    computation of Jacobians.

430:    Collective on Mat

432:    Input Parameters:
433: +  mat - the matrix containing the nonzero structure of the Jacobian
434: -  iscoloring - the coloring of the matrix; usually obtained with MatColoringCreate() or DMCreateColoring()

436:     Output Parameter:
437: .   color - the new coloring context

439:     Level: intermediate

441: .seealso: MatFDColoringDestroy(),SNESComputeJacobianDefaultColor(), ISColoringCreate(),
442:           MatFDColoringSetFunction(), MatFDColoringSetFromOptions(), MatFDColoringApply(),
443:           MatFDColoringView(), MatFDColoringSetParameters(), MatColoringCreate(), DMCreateColoring()
444: @*/
445: PetscErrorCode  MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
446: {
447:   MatFDColoring  c;
448:   MPI_Comm       comm;
450:   PetscInt       M,N;

454:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Matrix must be assembled by calls to MatAssemblyBegin/End();");
455:   PetscLogEventBegin(MAT_FDColoringCreate,mat,0,0,0);
456:   MatGetSize(mat,&M,&N);
457:   if (M != N) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Only for square matrices");
458:   PetscObjectGetComm((PetscObject)mat,&comm);
459:   PetscHeaderCreate(c,MAT_FDCOLORING_CLASSID,"MatFDColoring","Jacobian computation via finite differences with coloring","Mat",comm,MatFDColoringDestroy,MatFDColoringView);

461:   c->ctype = iscoloring->ctype;

463:   if (mat->ops->fdcoloringcreate) {
464:     (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);
465:   } else SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Code not yet written for matrix type %s",((PetscObject)mat)->type_name);

467:   MatCreateVecs(mat,NULL,&c->w1);
468:   /* Vec is used instensively in particular piece of scalar CPU code; won't benifit from bouncing back and forth to the GPU */
469:   VecPinToCPU(c->w1,PETSC_TRUE);
470:   PetscLogObjectParent((PetscObject)c,(PetscObject)c->w1);
471:   VecDuplicate(c->w1,&c->w2);
472:   /* Vec is used instensively in particular piece of scalar CPU code; won't benifit from bouncing back and forth to the GPU */
473:   VecPinToCPU(c->w2,PETSC_TRUE);
474:   PetscLogObjectParent((PetscObject)c,(PetscObject)c->w2);

476:   c->error_rel    = PETSC_SQRT_MACHINE_EPSILON;
477:   c->umin         = 100.0*PETSC_SQRT_MACHINE_EPSILON;
478:   c->currentcolor = -1;
479:   c->htype        = "wp";
480:   c->fset         = PETSC_FALSE;
481:   c->setupcalled  = PETSC_FALSE;

483:   *color = c;
484:   PetscObjectCompose((PetscObject)mat,"SNESMatFDColoring",(PetscObject)c);
485:   PetscLogEventEnd(MAT_FDColoringCreate,mat,0,0,0);
486:   return(0);
487: }

489: /*@
490:     MatFDColoringDestroy - Destroys a matrix coloring context that was created
491:     via MatFDColoringCreate().

493:     Collective on MatFDColoring

495:     Input Parameter:
496: .   c - coloring context

498:     Level: intermediate

500: .seealso: MatFDColoringCreate()
501: @*/
502: PetscErrorCode  MatFDColoringDestroy(MatFDColoring *c)
503: {
505:   PetscInt       i;
506:   MatFDColoring  color = *c;

509:   if (!*c) return(0);
510:   if (--((PetscObject)color)->refct > 0) {*c = 0; return(0);}

512:   /* we do not free the column arrays since their entries are owned by the ISs in color->isa */
513:   for (i=0; i<color->ncolors; i++) {
514:     ISDestroy(&color->isa[i]);
515:   }
516:   PetscFree(color->isa);
517:   PetscFree2(color->ncolumns,color->columns);
518:   PetscFree(color->nrows);
519:   if (color->htype[0] == 'w') {
520:     PetscFree(color->matentry2);
521:   } else {
522:     PetscFree(color->matentry);
523:   }
524:   PetscFree(color->dy);
525:   if (color->vscale) {VecDestroy(&color->vscale);}
526:   VecDestroy(&color->w1);
527:   VecDestroy(&color->w2);
528:   VecDestroy(&color->w3);
529:   PetscHeaderDestroy(c);
530:   return(0);
531: }

533: /*@C
534:     MatFDColoringGetPerturbedColumns - Returns the indices of the columns that
535:       that are currently being perturbed.

537:     Not Collective

539:     Input Parameters:
540: .   coloring - coloring context created with MatFDColoringCreate()

542:     Output Parameters:
543: +   n - the number of local columns being perturbed
544: -   cols - the column indices, in global numbering

546:    Level: advanced

548:    Note: IF the matrix type is BAIJ, then the block column indices are returned

550:    Fortran Note:
551:    This routine has a different interface for Fortran
552: $     #include <petsc/finclude/petscmat.h>
553: $          use petscmat
554: $          PetscInt, pointer :: array(:)
555: $          PetscErrorCode  ierr
556: $          MatFDColoring   i
557: $          call MatFDColoringGetPerturbedColumnsF90(i,array,ierr)
558: $      use the entries of array ...
559: $          call MatFDColoringRestorePerturbedColumnsF90(i,array,ierr)

561: .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringApply()

563: @*/
564: PetscErrorCode  MatFDColoringGetPerturbedColumns(MatFDColoring coloring,PetscInt *n,const PetscInt *cols[])
565: {
567:   if (coloring->currentcolor >= 0) {
568:     *n    = coloring->ncolumns[coloring->currentcolor];
569:     *cols = coloring->columns[coloring->currentcolor];
570:   } else {
571:     *n = 0;
572:   }
573:   return(0);
574: }

576: /*@
577:     MatFDColoringApply - Given a matrix for which a MatFDColoring context
578:     has been created, computes the Jacobian for a function via finite differences.

580:     Collective on MatFDColoring

582:     Input Parameters:
583: +   mat - location to store Jacobian
584: .   coloring - coloring context created with MatFDColoringCreate()
585: .   x1 - location at which Jacobian is to be computed
586: -   sctx - context required by function, if this is being used with the SNES solver then it is SNES object, otherwise it is null

588:     Options Database Keys:
589: +    -mat_fd_type - "wp" or "ds"  (see MATMFFD_WP or MATMFFD_DS)
590: .    -mat_fd_coloring_view - Activates basic viewing or coloring
591: .    -mat_fd_coloring_view draw - Activates drawing of coloring
592: -    -mat_fd_coloring_view ::ascii_info - Activates viewing of coloring info

594:     Level: intermediate

596: .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringSetFunction()

598: @*/
599: PetscErrorCode  MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,void *sctx)
600: {

607:   if (!coloring->f) SETERRQ(PetscObjectComm((PetscObject)J),PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetFunction()");
608:   if (!J->ops->fdcoloringapply) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not supported for this matrix type %s",((PetscObject)J)->type_name);
609:   if (!coloring->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetUp()");

611:   MatSetUnfactored(J);
612:   PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);
613:   (*J->ops->fdcoloringapply)(J,coloring,x1,sctx);
614:   PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);
615:   if (!coloring->viewed) {
616:     MatFDColoringViewFromOptions(coloring,NULL,"-mat_fd_coloring_view");
617:     coloring->viewed = PETSC_TRUE;
618:   }
619:   return(0);
620: }