Actual source code: fdmatrix.c

petsc-3.3-p7 2013-05-11
  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>        /*I "petscmat.h" I*/

 11: PetscErrorCode  MatFDColoringSetF(MatFDColoring fd,Vec F)
 12: {
 14:   fd->F = F;
 15:   return(0);
 16: }

 20: static PetscErrorCode MatFDColoringView_Draw_Zoom(PetscDraw draw,void *Aa)
 21: {
 22:   MatFDColoring  fd = (MatFDColoring)Aa;
 24:   PetscInt       i,j;
 25:   PetscReal      x,y;


 29:   /* loop over colors  */
 30:   for (i=0; i<fd->ncolors; i++) {
 31:     for (j=0; j<fd->nrows[i]; j++) {
 32:       y = fd->M - fd->rows[i][j] - fd->rstart;
 33:       x = fd->columnsforrow[i][j];
 34:       PetscDrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);
 35:     }
 36:   }
 37:   return(0);
 38: }

 42: static PetscErrorCode MatFDColoringView_Draw(MatFDColoring fd,PetscViewer viewer)
 43: {
 45:   PetscBool      isnull;
 46:   PetscDraw      draw;
 47:   PetscReal      xr,yr,xl,yl,h,w;

 50:   PetscViewerDrawGetDraw(viewer,0,&draw);
 51:   PetscDrawIsNull(draw,&isnull); if (isnull) return(0);

 53:   PetscObjectCompose((PetscObject)fd,"Zoomviewer",(PetscObject)viewer);

 55:   xr  = fd->N; yr = fd->M; h = yr/10.0; w = xr/10.0;
 56:   xr += w;     yr += h;    xl = -w;     yl = -h;
 57:   PetscDrawSetCoordinates(draw,xl,yl,xr,yr);
 58:   PetscDrawZoom(draw,MatFDColoringView_Draw_Zoom,fd);
 59:   PetscObjectCompose((PetscObject)fd,"Zoomviewer",PETSC_NULL);
 60:   return(0);
 61: }

 65: /*@C
 66:    MatFDColoringView - Views a finite difference coloring context.

 68:    Collective on MatFDColoring

 70:    Input  Parameters:
 71: +  c - the coloring context
 72: -  viewer - visualization context

 74:    Level: intermediate

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

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

 90: .seealso: MatFDColoringCreate()

 92: .keywords: Mat, finite differences, coloring, view
 93: @*/
 94: PetscErrorCode  MatFDColoringView(MatFDColoring c,PetscViewer viewer)
 95: {
 96:   PetscErrorCode    ierr;
 97:   PetscInt          i,j;
 98:   PetscBool         isdraw,iascii;
 99:   PetscViewerFormat format;

103:   if (!viewer) {
104:     PetscViewerASCIIGetStdout(((PetscObject)c)->comm,&viewer);
105:   }

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

119:     PetscViewerGetFormat(viewer,&format);
120:     if (format != PETSC_VIEWER_ASCII_INFO) {
121:       for (i=0; i<c->ncolors; i++) {
122:         PetscViewerASCIIPrintf(viewer,"  Information for color %D\n",i);
123:         PetscViewerASCIIPrintf(viewer,"    Number of columns %D\n",c->ncolumns[i]);
124:         for (j=0; j<c->ncolumns[i]; j++) {
125:           PetscViewerASCIIPrintf(viewer,"      %D\n",c->columns[i][j]);
126:         }
127:         PetscViewerASCIIPrintf(viewer,"    Number of rows %D\n",c->nrows[i]);
128:         for (j=0; j<c->nrows[i]; j++) {
129:           PetscViewerASCIIPrintf(viewer,"      %D %D \n",c->rows[i][j],c->columnsforrow[i][j]);
130:         }
131:       }
132:     }
133:     PetscViewerFlush(viewer);
134:   } else {
135:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for MatFDColoring",((PetscObject)viewer)->type_name);
136:   }
137:   return(0);
138: }

142: /*@
143:    MatFDColoringSetParameters - Sets the parameters for the sparse approximation of
144:    a Jacobian matrix using finite differences.

146:    Logically Collective on MatFDColoring

148:    The Jacobian is estimated with the differencing approximation
149: .vb
150:        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
151:        h = error_rel*u[i]                 if  abs(u[i]) > umin
152:          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
153:        dx_{i} = (0, ... 1, .... 0)
154: .ve

156:    Input Parameters:
157: +  coloring - the coloring context
158: .  error_rel - relative error
159: -  umin - minimum allowable u-value magnitude

161:    Level: advanced

163: .keywords: Mat, finite differences, coloring, set, parameters

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

167: @*/
168: PetscErrorCode  MatFDColoringSetParameters(MatFDColoring matfd,PetscReal error,PetscReal umin)
169: {
174:   if (error != PETSC_DEFAULT) matfd->error_rel = error;
175:   if (umin != PETSC_DEFAULT)  matfd->umin      = umin;
176:   return(0);
177: }



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

186:    Not Collective

188:    Input Parameters:
189: .  coloring - the coloring context

191:    Output Parameters:
192: +  f - the function
193: -  fctx - the optional user-defined function context

195:    Level: intermediate

197: .keywords: Mat, Jacobian, finite differences, set, function

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

201: @*/
202: PetscErrorCode  MatFDColoringGetFunction(MatFDColoring matfd,PetscErrorCode (**f)(void),void **fctx)
203: {
206:   if (f) *f = matfd->f;
207:   if (fctx) *fctx = matfd->fctx;
208:   return(0);
209: }

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

216:    Logically Collective on MatFDColoring

218:    Input Parameters:
219: +  coloring - the coloring context
220: .  f - the function
221: -  fctx - the optional user-defined function context

223:    Calling sequence of (*f) function:
224:     For SNES:    PetscErrorCode (*f)(SNES,Vec,Vec,void*)
225:     If not using SNES: PetscErrorCode (*f)(void *dummy,Vec,Vec,void*) and dummy is ignored

227:    Level: advanced

229:    Notes: This function is usually used automatically by SNES (when one uses SNESSetJacobian() with the argument
230:      SNESDefaultComputeJacobianColor()) and only needs to be used by someone computing a matrix via coloring directly by
231:      calling MatFDColoringApply()

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

237: .keywords: Mat, Jacobian, finite differences, set, function

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

241: @*/
242: PetscErrorCode  MatFDColoringSetFunction(MatFDColoring matfd,PetscErrorCode (*f)(void),void *fctx)
243: {
246:   matfd->f    = f;
247:   matfd->fctx = fctx;
248:   return(0);
249: }

253: /*@
254:    MatFDColoringSetFromOptions - Sets coloring finite difference parameters from 
255:    the options database.

257:    Collective on MatFDColoring

259:    The Jacobian, F'(u), is estimated with the differencing approximation
260: .vb
261:        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
262:        h = error_rel*u[i]                 if  abs(u[i]) > umin
263:          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
264:        dx_{i} = (0, ... 1, .... 0)
265: .ve

267:    Input Parameter:
268: .  coloring - the coloring context

270:    Options Database Keys:
271: +  -mat_fd_coloring_err <err> - Sets <err> (square root
272:            of relative error in the function)
273: .  -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude
274: .  -mat_fd_type - "wp" or "ds" (see MATMFFD_WP or MATMFFD_DS)
275: .  -mat_fd_coloring_view - Activates basic viewing
276: .  -mat_fd_coloring_view_info - Activates viewing info
277: -  -mat_fd_coloring_view_draw - Activates drawing

279:     Level: intermediate

281: .keywords: Mat, finite differences, parameters

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

285: @*/
286: PetscErrorCode  MatFDColoringSetFromOptions(MatFDColoring matfd)
287: {
289:   PetscBool      flg;
290:   char           value[3];


295:   PetscObjectOptionsBegin((PetscObject)matfd);
296:     PetscOptionsReal("-mat_fd_coloring_err","Square root of relative error in function","MatFDColoringSetParameters",matfd->error_rel,&matfd->error_rel,0);
297:     PetscOptionsReal("-mat_fd_coloring_umin","Minimum allowable u magnitude","MatFDColoringSetParameters",matfd->umin,&matfd->umin,0);
298:     PetscOptionsString("-mat_fd_type","Algorithm to compute h, wp or ds","MatFDColoringCreate",matfd->htype,value,3,&flg);
299:     if (flg) {
300:       if (value[0] == 'w' && value[1] == 'p') matfd->htype = "wp";
301:       else if (value[0] == 'd' && value[1] == 's') matfd->htype = "ds";
302:       else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Unknown finite differencing type %s",value);
303:     }
304:     /* not used here; but so they are presented in the GUI */
305:     PetscOptionsName("-mat_fd_coloring_view","Print entire datastructure used for Jacobian","None",0);
306:     PetscOptionsName("-mat_fd_coloring_view_info","Print number of colors etc for Jacobian","None",0);
307:     PetscOptionsName("-mat_fd_coloring_view_draw","Plot nonzero structure ofJacobian","None",0);

309:     /* process any options handlers added with PetscObjectAddOptionsHandler() */
310:     PetscObjectProcessOptionsHandlers((PetscObject)matfd);
311:   PetscOptionsEnd();
312:   return(0);
313: }

317: PetscErrorCode MatFDColoringView_Private(MatFDColoring fd)
318: {
320:   PetscBool      flg = PETSC_FALSE;
321:   PetscViewer    viewer;

324:   PetscViewerASCIIGetStdout(((PetscObject)fd)->comm,&viewer);
325:   PetscOptionsGetBool(PETSC_NULL,"-mat_fd_coloring_view",&flg,PETSC_NULL);
326:   if (flg) {
327:     MatFDColoringView(fd,viewer);
328:   }
329:   flg  = PETSC_FALSE;
330:   PetscOptionsGetBool(PETSC_NULL,"-mat_fd_coloring_view_info",&flg,PETSC_NULL);
331:   if (flg) {
332:     PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);
333:     MatFDColoringView(fd,viewer);
334:     PetscViewerPopFormat(viewer);
335:   }
336:   flg  = PETSC_FALSE;
337:   PetscOptionsGetBool(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg,PETSC_NULL);
338:   if (flg) {
339:     MatFDColoringView(fd,PETSC_VIEWER_DRAW_(((PetscObject)fd)->comm));
340:     PetscViewerFlush(PETSC_VIEWER_DRAW_(((PetscObject)fd)->comm));
341:   }
342:   return(0);
343: }

347: /*@
348:    MatFDColoringCreate - Creates a matrix coloring context for finite difference 
349:    computation of Jacobians.

351:    Collective on Mat

353:    Input Parameters:
354: +  mat - the matrix containing the nonzero structure of the Jacobian
355: -  iscoloring - the coloring of the matrix; usually obtained with MatGetColoring() or DMCreateColoring()

357:     Output Parameter:
358: .   color - the new coloring context
359:    
360:     Level: intermediate

362: .seealso: MatFDColoringDestroy(),SNESDefaultComputeJacobianColor(), ISColoringCreate(),
363:           MatFDColoringSetFunction(), MatFDColoringSetFromOptions(), MatFDColoringApply(),
364:           MatFDColoringView(), MatFDColoringSetParameters(), MatGetColoring(), DMCreateColoring()
365: @*/
366: PetscErrorCode  MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
367: {
368:   MatFDColoring  c;
369:   MPI_Comm       comm;
371:   PetscInt       M,N;

374:   PetscLogEventBegin(MAT_FDColoringCreate,mat,0,0,0);
375:   MatGetSize(mat,&M,&N);
376:   if (M != N) SETERRQ(((PetscObject)mat)->comm,PETSC_ERR_SUP,"Only for square matrices");

378:   PetscObjectGetComm((PetscObject)mat,&comm);
379:   PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_CLASSID,0,"MatFDColoring","Jacobian computation via finite differences with coloring","Mat",comm,MatFDColoringDestroy,MatFDColoringView);

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

383:   if (mat->ops->fdcoloringcreate) {
384:     (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);
385:   } else SETERRQ1(((PetscObject)mat)->comm,PETSC_ERR_SUP,"Code not yet written for matrix type %s",((PetscObject)mat)->type_name);

387:   MatGetVecs(mat,PETSC_NULL,&c->w1);
388:   PetscLogObjectParent(c,c->w1);
389:   VecDuplicate(c->w1,&c->w2);
390:   PetscLogObjectParent(c,c->w2);

392:   c->error_rel         = PETSC_SQRT_MACHINE_EPSILON;
393:   c->umin              = 100.0*PETSC_SQRT_MACHINE_EPSILON;
394:   c->currentcolor      = -1;
395:   c->htype             = "wp";

397:   *color = c;
398:   PetscLogEventEnd(MAT_FDColoringCreate,mat,0,0,0);
399:   return(0);
400: }

404: /*@
405:     MatFDColoringDestroy - Destroys a matrix coloring context that was created
406:     via MatFDColoringCreate().

408:     Collective on MatFDColoring

410:     Input Parameter:
411: .   c - coloring context

413:     Level: intermediate

415: .seealso: MatFDColoringCreate()
416: @*/
417: PetscErrorCode  MatFDColoringDestroy(MatFDColoring *c)
418: {
420:   PetscInt       i;

423:   if (!*c) return(0);
424:   if (--((PetscObject)(*c))->refct > 0) {*c = 0; return(0);}

426:   for (i=0; i<(*c)->ncolors; i++) {
427:     PetscFree((*c)->columns[i]);
428:     PetscFree((*c)->rows[i]);
429:     PetscFree((*c)->columnsforrow[i]);
430:     if ((*c)->vscaleforrow) {PetscFree((*c)->vscaleforrow[i]);}
431:   }
432:   PetscFree((*c)->ncolumns);
433:   PetscFree((*c)->columns);
434:   PetscFree((*c)->nrows);
435:   PetscFree((*c)->rows);
436:   PetscFree((*c)->columnsforrow);
437:   PetscFree((*c)->vscaleforrow);
438:   VecDestroy(&(*c)->vscale);
439:   VecDestroy(&(*c)->w1);
440:   VecDestroy(&(*c)->w2);
441:   VecDestroy(&(*c)->w3);
442:   PetscHeaderDestroy(c);
443:   return(0);
444: }

448: /*@C
449:     MatFDColoringGetPerturbedColumns - Returns the indices of the columns that
450:       that are currently being perturbed.

452:     Not Collective

454:     Input Parameters:
455: .   coloring - coloring context created with MatFDColoringCreate()

457:     Output Parameters:
458: +   n - the number of local columns being perturbed
459: -   cols - the column indices, in global numbering

461:    Level: intermediate

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

465: .keywords: coloring, Jacobian, finite differences
466: @*/
467: PetscErrorCode  MatFDColoringGetPerturbedColumns(MatFDColoring coloring,PetscInt *n,PetscInt *cols[])
468: {
470:   if (coloring->currentcolor >= 0) {
471:     *n    = coloring->ncolumns[coloring->currentcolor];
472:     *cols = coloring->columns[coloring->currentcolor];
473:   } else {
474:     *n = 0;
475:   }
476:   return(0);
477: }

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

485:     Collective on MatFDColoring

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

493:     Options Database Keys:
494: +    -mat_fd_type - "wp" or "ds"  (see MATMFFD_WP or MATMFFD_DS)
495: .    -mat_fd_coloring_view - Activates basic viewing or coloring
496: .    -mat_fd_coloring_view_draw - Activates drawing of coloring
497: -    -mat_fd_coloring_view_info - Activates viewing of coloring info

499:     Level: intermediate

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

503: .keywords: coloring, Jacobian, finite differences
504: @*/
505: PetscErrorCode  MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
506: {

513:   if (!coloring->f) SETERRQ(((PetscObject)J)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetFunction()");
514:   if (!J->ops->fdcoloringapply) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not supported for this matrix type %s",((PetscObject)J)->type_name);
515:   (*J->ops->fdcoloringapply)(J,coloring,x1,flag,sctx);
516:   return(0);
517: }

521: PetscErrorCode  MatFDColoringApply_AIJ(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
522: {
523:   PetscErrorCode (*f)(void*,Vec,Vec,void*) = (PetscErrorCode (*)(void*,Vec,Vec,void *))coloring->f;
525:   PetscInt       k,start,end,l,row,col,srow,**vscaleforrow,m1,m2;
526:   PetscScalar    dx,*y,*xx,*w3_array;
527:   PetscScalar    *vscale_array;
528:   PetscReal      epsilon = coloring->error_rel,umin = coloring->umin,unorm;
529:   Vec            w1=coloring->w1,w2=coloring->w2,w3;
530:   void           *fctx = coloring->fctx;
531:   PetscBool      flg = PETSC_FALSE;
532:   PetscInt       ctype=coloring->ctype,N,col_start=0,col_end=0;
533:   Vec            x1_tmp;

539:   if (!f) SETERRQ(((PetscObject)J)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetFunction()");

541:   PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);
542:   MatSetUnfactored(J);
543:   PetscOptionsGetBool(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg,PETSC_NULL);
544:   if (flg) {
545:     PetscInfo(coloring,"Not calling MatZeroEntries()\n");
546:   } else {
547:     PetscBool  assembled;
548:     MatAssembled(J,&assembled);
549:     if (assembled) {
550:       MatZeroEntries(J);
551:     }
552:   }

554:   x1_tmp = x1;
555:   if (!coloring->vscale){
556:     VecDuplicate(x1_tmp,&coloring->vscale);
557:   }
558: 
559:   /*
560:     This is a horrible, horrible, hack. 
561:   */
562:   if (coloring->F) {
563:     VecGetLocalSize(coloring->F,&m1);
564:     VecGetLocalSize(w1,&m2);
565:     if (m1 != m2) {
566:       coloring->F = 0;
567:       }
568:     }

570:   if (coloring->htype[0] == 'w') { /* tacky test; need to make systematic if we add other approaches to computing h*/
571:     VecNorm(x1_tmp,NORM_2,&unorm);
572:   }
573:   VecGetOwnershipRange(w1,&start,&end); /* OwnershipRange is used by ghosted x! */
574: 
575:   /* Set w1 = F(x1) */
576:   if (coloring->F) {
577:     w1          = coloring->F; /* use already computed value of function */
578:     coloring->F = 0;
579:   } else {
580:     PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);
581:     (*f)(sctx,x1_tmp,w1,fctx);
582:     PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);
583:   }
584: 
585:   if (!coloring->w3) {
586:     VecDuplicate(x1_tmp,&coloring->w3);
587:     PetscLogObjectParent(coloring,coloring->w3);
588:   }
589:   w3 = coloring->w3;

591:     /* Compute all the local scale factors, including ghost points */
592:   VecGetLocalSize(x1_tmp,&N);
593:   VecGetArray(x1_tmp,&xx);
594:   VecGetArray(coloring->vscale,&vscale_array);
595:   if (ctype == IS_COLORING_GHOSTED){
596:     col_start = 0; col_end = N;
597:   } else if (ctype == IS_COLORING_GLOBAL){
598:     xx = xx - start;
599:     vscale_array = vscale_array - start;
600:     col_start = start; col_end = N + start;
601:   }
602:   for (col=col_start; col<col_end; col++){
603:     /* Loop over each local column, vscale[col] = 1./(epsilon*dx[col]) */
604:     if (coloring->htype[0] == 'w') {
605:       dx = 1.0 + unorm;
606:     } else {
607:       dx  = xx[col];
608:     }
609:     if (dx == (PetscScalar)0.0) dx = 1.0;
610: #if !defined(PETSC_USE_COMPLEX)
611:     if (dx < umin && dx >= 0.0)      dx = umin;
612:     else if (dx < 0.0 && dx > -umin) dx = -umin;
613: #else
614:     if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
615:     else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
616: #endif
617:     dx               *= epsilon;
618:     vscale_array[col] = (PetscScalar)1.0/dx;
619:   }
620:   if (ctype == IS_COLORING_GLOBAL)  vscale_array = vscale_array + start;
621:   VecRestoreArray(coloring->vscale,&vscale_array);
622:   if (ctype == IS_COLORING_GLOBAL){
623:     VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);
624:     VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);
625:   }
626: 
627:   if (coloring->vscaleforrow) {
628:     vscaleforrow = coloring->vscaleforrow;
629:   } else SETERRQ(((PetscObject)J)->comm,PETSC_ERR_ARG_NULL,"Null Object: coloring->vscaleforrow");

631:   /*
632:     Loop over each color
633:   */
634:   VecGetArray(coloring->vscale,&vscale_array);
635:   for (k=0; k<coloring->ncolors; k++) {
636:     coloring->currentcolor = k;
637:     VecCopy(x1_tmp,w3);
638:     VecGetArray(w3,&w3_array);
639:     if (ctype == IS_COLORING_GLOBAL) w3_array = w3_array - start;
640:     /*
641:       Loop over each column associated with color 
642:       adding the perturbation to the vector w3.
643:     */
644:     for (l=0; l<coloring->ncolumns[k]; l++) {
645:       col = coloring->columns[k][l];    /* local column of the matrix we are probing for */
646:       if (coloring->htype[0] == 'w') {
647:         dx = 1.0 + unorm;
648:       } else {
649:         dx  = xx[col];
650:       }
651:       if (dx == (PetscScalar)0.0) dx = 1.0;
652: #if !defined(PETSC_USE_COMPLEX)
653:       if (dx < umin && dx >= 0.0)      dx = umin;
654:       else if (dx < 0.0 && dx > -umin) dx = -umin;
655: #else
656:       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
657:       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
658: #endif
659:       dx            *= epsilon;
660:       if (!PetscAbsScalar(dx)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Computed 0 differencing parameter");
661:       w3_array[col] += dx;
662:     }
663:     if (ctype == IS_COLORING_GLOBAL) w3_array = w3_array + start;
664:     VecRestoreArray(w3,&w3_array);

666:     /*
667:       Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations)
668:                            w2 = F(x1 + dx) - F(x1)
669:     */
670:     PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);
671:     (*f)(sctx,w3,w2,fctx);
672:     PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);
673:     VecAXPY(w2,-1.0,w1);
674: 
675:     /*
676:       Loop over rows of vector, putting results into Jacobian matrix
677:     */
678:     VecGetArray(w2,&y);
679:     for (l=0; l<coloring->nrows[k]; l++) {
680:       row    = coloring->rows[k][l];             /* local row index */
681:       col    = coloring->columnsforrow[k][l];    /* global column index */
682:       y[row] *= vscale_array[vscaleforrow[k][l]];
683:       srow   = row + start;
684:       MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);
685:     }
686:     VecRestoreArray(w2,&y);
687:   } /* endof for each color */
688:   if (ctype == IS_COLORING_GLOBAL) xx = xx + start;
689:   VecRestoreArray(coloring->vscale,&vscale_array);
690:   VecRestoreArray(x1_tmp,&xx);
691: 
692:   coloring->currentcolor = -1;
693:   MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
694:   MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
695:   PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);

697:   flg  = PETSC_FALSE;
698:   PetscOptionsGetBool(PETSC_NULL,"-mat_null_space_test",&flg,PETSC_NULL);
699:   if (flg) {
700:     MatNullSpaceTest(J->nullsp,J,PETSC_NULL);
701:   }
702:   MatFDColoringView_Private(coloring);
703:   return(0);
704: }