Actual source code: fdmatrix.c

  1: #define PETSCMAT_DLL

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

 8:  #include private/matimpl.h

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

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


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

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

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

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

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

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

 69:    Collective on MatFDColoring

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

 75:    Level: intermediate

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

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

 91: .seealso: MatFDColoringCreate()

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

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

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

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

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

147:    Collective on MatFDColoring

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

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

162:    Level: advanced

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

166: .seealso: MatFDColoringCreate()
167: @*/
168: PetscErrorCode  MatFDColoringSetParameters(MatFDColoring matfd,PetscReal error,PetscReal umin)
169: {

173:   if (error != PETSC_DEFAULT) matfd->error_rel = error;
174:   if (umin != PETSC_DEFAULT)  matfd->umin      = umin;
175:   return(0);
176: }



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

185:    Collective on MatFDColoring

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

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

194:    Level: intermediate

196: .keywords: Mat, Jacobian, finite differences, set, function
197: @*/
198: PetscErrorCode  MatFDColoringGetFunction(MatFDColoring matfd,PetscErrorCode (**f)(void),void **fctx)
199: {
202:   if (f) *f = matfd->f;
203:   if (fctx) *fctx = matfd->fctx;
204:   return(0);
205: }

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

212:    Collective on MatFDColoring

214:    Input Parameters:
215: +  coloring - the coloring context
216: .  f - the function
217: -  fctx - the optional user-defined function context

219:    Calling sequence of (*f) function:
220:     For SNES:    PetscErrorCode (*f)(SNES,Vec,Vec,void*)
221:     For TS:      PetscErrorCode (*f)(TS,PetscReal,Vec,Vec,void*)
222:     If not using SNES or TS: PetscErrorCode (*f)(void *dummy,Vec,Vec,void*) and dummy is ignored

224:    Level: advanced

226:    Notes: This function is usually used automatically by SNES or TS (when one uses SNESSetJacobian() with the argument 
227:      SNESDefaultComputeJacobianColor() or TSSetRHSJacobian() with the argument TSDefaultComputeJacobianColor()) and only needs to be used
228:      by someone computing a matrix via coloring directly by calling MatFDColoringApply()

230:    Fortran Notes:
231:     In Fortran you must call MatFDColoringSetFunction() for a coloring object to 
232:   be used without SNES or TS or within the SNES solvers and MatFDColoringSetFunctionTS() if it is to be used
233:   within the TS solvers.

235: .keywords: Mat, Jacobian, finite differences, set, function
236: @*/
237: PetscErrorCode  MatFDColoringSetFunction(MatFDColoring matfd,PetscErrorCode (*f)(void),void *fctx)
238: {
241:   matfd->f    = f;
242:   matfd->fctx = fctx;
243:   return(0);
244: }

248: /*@
249:    MatFDColoringSetFromOptions - Sets coloring finite difference parameters from 
250:    the options database.

252:    Collective on MatFDColoring

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

262:    Input Parameter:
263: .  coloring - the coloring context

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

274:     Level: intermediate

276: .keywords: Mat, finite differences, parameters

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

280: @*/
281: PetscErrorCode  MatFDColoringSetFromOptions(MatFDColoring matfd)
282: {
284:   PetscTruth     flg;
285:   char           value[3];


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

309: PetscErrorCode MatFDColoringView_Private(MatFDColoring fd)
310: {
312:   PetscTruth     flg;
313:   PetscViewer    viewer;

316:   PetscViewerASCIIGetStdout(((PetscObject)fd)->comm,&viewer);
317:   PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view",&flg);
318:   if (flg) {
319:     MatFDColoringView(fd,viewer);
320:   }
321:   PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_info",&flg);
322:   if (flg) {
323:     PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);
324:     MatFDColoringView(fd,viewer);
325:     PetscViewerPopFormat(viewer);
326:   }
327:   PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg);
328:   if (flg) {
329:     MatFDColoringView(fd,PETSC_VIEWER_DRAW_(((PetscObject)fd)->comm));
330:     PetscViewerFlush(PETSC_VIEWER_DRAW_(((PetscObject)fd)->comm));
331:   }
332:   return(0);
333: }

337: /*@
338:    MatFDColoringCreate - Creates a matrix coloring context for finite difference 
339:    computation of Jacobians.

341:    Collective on Mat

343:    Input Parameters:
344: +  mat - the matrix containing the nonzero structure of the Jacobian
345: -  iscoloring - the coloring of the matrix

347:     Output Parameter:
348: .   color - the new coloring context
349:    
350:     Level: intermediate

352: .seealso: MatFDColoringDestroy(),SNESDefaultComputeJacobianColor(), ISColoringCreate(),
353:           MatFDColoringSetFunction(), MatFDColoringSetFromOptions(), MatFDColoringApply(),
354:           MatFDColoringView(), MatFDColoringSetParameters()
355: @*/
356: PetscErrorCode  MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
357: {
358:   MatFDColoring  c;
359:   MPI_Comm       comm;
361:   PetscInt       M,N;
362:   PetscMPIInt    size;

365:   PetscLogEventBegin(MAT_FDColoringCreate,mat,0,0,0);
366:   MatGetSize(mat,&M,&N);
367:   if (M != N) SETERRQ(PETSC_ERR_SUP,"Only for square matrices");

369:   PetscObjectGetComm((PetscObject)mat,&comm);
370:   PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_COOKIE,0,"MatFDColoring",comm,MatFDColoringDestroy,MatFDColoringView);

372:   MPI_Comm_size(comm,&size);
373:   c->ctype = iscoloring->ctype;

375:   if (mat->ops->fdcoloringcreate) {
376:     (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);
377:   } else {
378:     SETERRQ(PETSC_ERR_SUP,"Code not yet written for this matrix type");
379:   }

381:   MatGetVecs(mat,PETSC_NULL,&c->w1);
382:   PetscLogObjectParent(c,c->w1);
383:   VecDuplicate(c->w1,&c->w2);
384:   PetscLogObjectParent(c,c->w2);

386:   c->error_rel         = PETSC_SQRT_MACHINE_EPSILON;
387:   c->umin              = 100.0*PETSC_SQRT_MACHINE_EPSILON;
388:   c->currentcolor      = -1;
389:   c->htype             = "wp";

391:   *color = c;
392:   PetscLogEventEnd(MAT_FDColoringCreate,mat,0,0,0);
393:   return(0);
394: }

398: /*@
399:     MatFDColoringDestroy - Destroys a matrix coloring context that was created
400:     via MatFDColoringCreate().

402:     Collective on MatFDColoring

404:     Input Parameter:
405: .   c - coloring context

407:     Level: intermediate

409: .seealso: MatFDColoringCreate()
410: @*/
411: PetscErrorCode  MatFDColoringDestroy(MatFDColoring c)
412: {
414:   PetscInt       i;

417:   if (--((PetscObject)c)->refct > 0) return(0);

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

445: /*@C
446:     MatFDColoringGetPerturbedColumns - Returns the indices of the columns that
447:       that are currently being perturbed.

449:     Not Collective

451:     Input Parameters:
452: .   coloring - coloring context created with MatFDColoringCreate()

454:     Output Parameters:
455: +   n - the number of local columns being perturbed
456: -   cols - the column indices, in global numbering

458:    Level: intermediate

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

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

476: #include "petscda.h"      /*I      "petscda.h"    I*/ 

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

484:     Collective on MatFDColoring

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

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

498:     Level: intermediate

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

502: .keywords: coloring, Jacobian, finite differences
503: @*/
504: PetscErrorCode  MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
505: {
506:   PetscErrorCode (*f)(void*,Vec,Vec,void*) = (PetscErrorCode (*)(void*,Vec,Vec,void *))coloring->f;
508:   PetscInt       k,start,end,l,row,col,srow,**vscaleforrow,m1,m2;
509:   PetscScalar    dx,*y,*xx,*w3_array;
510:   PetscScalar    *vscale_array;
511:   PetscReal      epsilon = coloring->error_rel,umin = coloring->umin,unorm;
512:   Vec            w1=coloring->w1,w2=coloring->w2,w3;
513:   void           *fctx = coloring->fctx;
514:   PetscTruth     flg;
515:   PetscInt       ctype=coloring->ctype,N,col_start=0,col_end=0;
516:   Vec            x1_tmp;

522:   if (!f) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetFunction()");

524:   PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);
525:   MatSetUnfactored(J);
526:   PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);
527:   if (flg) {
528:     PetscInfo(coloring,"Not calling MatZeroEntries()\n");
529:   } else {
530:     PetscTruth assembled;
531:     MatAssembled(J,&assembled);
532:     if (assembled) {
533:       MatZeroEntries(J);
534:     }
535:   }

537:   x1_tmp = x1;
538:   if (!coloring->vscale){
539:     VecDuplicate(x1_tmp,&coloring->vscale);
540:   }
541: 
542:   /*
543:     This is a horrible, horrible, hack. See DMMGComputeJacobian_Multigrid() it inproperly sets
544:     coloring->F for the coarser grids from the finest
545:   */
546:   if (coloring->F) {
547:     VecGetLocalSize(coloring->F,&m1);
548:     VecGetLocalSize(w1,&m2);
549:     if (m1 != m2) {
550:       coloring->F = 0;
551:       }
552:     }

554:   if (coloring->htype[0] == 'w') { /* tacky test; need to make systematic if we add other approaches to computing h*/
555:     VecNorm(x1_tmp,NORM_2,&unorm);
556:   }
557:   VecGetOwnershipRange(w1,&start,&end); /* OwnershipRange is used by ghosted x! */
558: 
559:   /* Set w1 = F(x1) */
560:   if (coloring->F) {
561:     w1          = coloring->F; /* use already computed value of function */
562:     coloring->F = 0;
563:   } else {
564:     PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);
565:     (*f)(sctx,x1_tmp,w1,fctx);
566:     PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);
567:   }
568: 
569:   if (!coloring->w3) {
570:     VecDuplicate(x1_tmp,&coloring->w3);
571:     PetscLogObjectParent(coloring,coloring->w3);
572:   }
573:   w3 = coloring->w3;

575:     /* Compute all the local scale factors, including ghost points */
576:   VecGetLocalSize(x1_tmp,&N);
577:   VecGetArray(x1_tmp,&xx);
578:   VecGetArray(coloring->vscale,&vscale_array);
579:   if (ctype == IS_COLORING_GHOSTED){
580:     col_start = 0; col_end = N;
581:   } else if (ctype == IS_COLORING_GLOBAL){
582:     xx = xx - start;
583:     vscale_array = vscale_array - start;
584:     col_start = start; col_end = N + start;
585:   }
586:   for (col=col_start; col<col_end; col++){
587:     /* Loop over each local column, vscale[col] = 1./(epsilon*dx[col]) */
588:     if (coloring->htype[0] == 'w') {
589:       dx = 1.0 + unorm;
590:     } else {
591:       dx  = xx[col];
592:     }
593:     if (dx == 0.0) dx = 1.0;
594: #if !defined(PETSC_USE_COMPLEX)
595:     if (dx < umin && dx >= 0.0)      dx = umin;
596:     else if (dx < 0.0 && dx > -umin) dx = -umin;
597: #else
598:     if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
599:     else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
600: #endif
601:     dx               *= epsilon;
602:     vscale_array[col] = 1.0/dx;
603:   }
604:   if (ctype == IS_COLORING_GLOBAL)  vscale_array = vscale_array + start;
605:   VecRestoreArray(coloring->vscale,&vscale_array);
606:   if (ctype == IS_COLORING_GLOBAL){
607:     VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);
608:     VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);
609:   }
610: 
611:   if (coloring->vscaleforrow) {
612:     vscaleforrow = coloring->vscaleforrow;
613:   } else {
614:     SETERRQ(PETSC_ERR_ARG_NULL,"Null Object: coloring->vscaleforrow");
615:   }

617:   /*
618:     Loop over each color
619:   */
620:   VecGetArray(coloring->vscale,&vscale_array);
621:   for (k=0; k<coloring->ncolors; k++) {
622:     coloring->currentcolor = k;
623:     VecCopy(x1_tmp,w3);
624:     VecGetArray(w3,&w3_array);
625:     if (ctype == IS_COLORING_GLOBAL) w3_array = w3_array - start;
626:     /*
627:       Loop over each column associated with color 
628:       adding the perturbation to the vector w3.
629:     */
630:     for (l=0; l<coloring->ncolumns[k]; l++) {
631:       col = coloring->columns[k][l];    /* local column of the matrix we are probing for */
632:       if (coloring->htype[0] == 'w') {
633:         dx = 1.0 + unorm;
634:       } else {
635:         dx  = xx[col];
636:       }
637:       if (dx == 0.0) dx = 1.0;
638: #if !defined(PETSC_USE_COMPLEX)
639:       if (dx < umin && dx >= 0.0)      dx = umin;
640:       else if (dx < 0.0 && dx > -umin) dx = -umin;
641: #else
642:       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
643:       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
644: #endif
645:       dx            *= epsilon;
646:       if (!PetscAbsScalar(dx)) SETERRQ(PETSC_ERR_PLIB,"Computed 0 differencing parameter");
647:       w3_array[col] += dx;
648:     }
649:     if (ctype == IS_COLORING_GLOBAL) w3_array = w3_array + start;
650:     VecRestoreArray(w3,&w3_array);

652:     /*
653:       Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations)
654:                            w2 = F(x1 + dx) - F(x1)
655:     */
656:     PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);
657:     (*f)(sctx,w3,w2,fctx);
658:     PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);
659:     VecAXPY(w2,-1.0,w1);
660: 
661:     /*
662:       Loop over rows of vector, putting results into Jacobian matrix
663:     */
664:     VecGetArray(w2,&y);
665:     for (l=0; l<coloring->nrows[k]; l++) {
666:       row    = coloring->rows[k][l];             /* local row index */
667:       col    = coloring->columnsforrow[k][l];    /* global column index */
668:       y[row] *= vscale_array[vscaleforrow[k][l]];
669:       srow   = row + start;
670:       MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);
671:     }
672:     VecRestoreArray(w2,&y);
673:   } /* endof for each color */
674:   if (ctype == IS_COLORING_GLOBAL) xx = xx + start;
675:   VecRestoreArray(coloring->vscale,&vscale_array);
676:   VecRestoreArray(x1_tmp,&xx);
677: 
678:   coloring->currentcolor = -1;
679:   MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
680:   MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
681:   PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);

683:   PetscOptionsHasName(PETSC_NULL,"-mat_null_space_test",&flg);
684:   if (flg) {
685:     MatNullSpaceTest(J->nullsp,J,PETSC_NULL);
686:   }
687:   MatFDColoringView_Private(coloring);
688:   return(0);
689: }

693: /*@
694:     MatFDColoringApplyTS - Given a matrix for which a MatFDColoring context 
695:     has been created, computes the Jacobian for a function via finite differences.

697:    Collective on Mat, MatFDColoring, and Vec

699:     Input Parameters:
700: +   mat - location to store Jacobian
701: .   coloring - coloring context created with MatFDColoringCreate()
702: .   x1 - location at which Jacobian is to be computed
703: -   sctx - context required by function, if this is being used with the TS solver then it is TS object, otherwise it is null

705:    Level: intermediate

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

709: .keywords: coloring, Jacobian, finite differences
710: @*/
711: PetscErrorCode  MatFDColoringApplyTS(Mat J,MatFDColoring coloring,PetscReal t,Vec x1,MatStructure *flag,void *sctx)
712: {
713:   PetscErrorCode (*f)(void*,PetscReal,Vec,Vec,void*)=(PetscErrorCode (*)(void*,PetscReal,Vec,Vec,void *))coloring->f;
715:   PetscInt       k,N,start,end,l,row,col,srow,**vscaleforrow;
716:   PetscScalar    dx,*y,*xx,*w3_array;
717:   PetscScalar    *vscale_array;
718:   PetscReal      epsilon = coloring->error_rel,umin = coloring->umin;
719:   Vec            w1=coloring->w1,w2=coloring->w2,w3;
720:   void           *fctx = coloring->fctx;
721:   PetscTruth     flg;


728:   PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);
729:   if (!coloring->w3) {
730:     VecDuplicate(x1,&coloring->w3);
731:     PetscLogObjectParent(coloring,coloring->w3);
732:   }
733:   w3 = coloring->w3;

735:   MatSetUnfactored(J);
736:   PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);
737:   if (flg) {
738:     PetscInfo(coloring,"Not calling MatZeroEntries()\n");
739:   } else {
740:     PetscTruth assembled;
741:     MatAssembled(J,&assembled);
742:     if (assembled) {
743:       MatZeroEntries(J);
744:     }
745:   }

747:   VecGetOwnershipRange(x1,&start,&end);
748:   VecGetSize(x1,&N);
749:   PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);
750:   (*f)(sctx,t,x1,w1,fctx);
751:   PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);

753:   /* 
754:       Compute all the scale factors and share with other processors
755:   */
756:   VecGetArray(x1,&xx);xx = xx - start;
757:   VecGetArray(coloring->vscale,&vscale_array);vscale_array = vscale_array - start;
758:   for (k=0; k<coloring->ncolors; k++) {
759:     /*
760:        Loop over each column associated with color adding the 
761:        perturbation to the vector w3.
762:     */
763:     for (l=0; l<coloring->ncolumns[k]; l++) {
764:       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
765:       dx  = xx[col];
766:       if (dx == 0.0) dx = 1.0;
767: #if !defined(PETSC_USE_COMPLEX)
768:       if (dx < umin && dx >= 0.0)      dx = umin;
769:       else if (dx < 0.0 && dx > -umin) dx = -umin;
770: #else
771:       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
772:       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
773: #endif
774:       dx                *= epsilon;
775:       vscale_array[col] = 1.0/dx;
776:     }
777:   }
778:   vscale_array = vscale_array - start;VecRestoreArray(coloring->vscale,&vscale_array);
779:   VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);
780:   VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);

782:   if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow;
783:   else                        vscaleforrow = coloring->columnsforrow;

785:   VecGetArray(coloring->vscale,&vscale_array);
786:   /*
787:       Loop over each color
788:   */
789:   for (k=0; k<coloring->ncolors; k++) {
790:     VecCopy(x1,w3);
791:     VecGetArray(w3,&w3_array);w3_array = w3_array - start;
792:     /*
793:        Loop over each column associated with color adding the 
794:        perturbation to the vector w3.
795:     */
796:     for (l=0; l<coloring->ncolumns[k]; l++) {
797:       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
798:       dx  = xx[col];
799:       if (dx == 0.0) dx = 1.0;
800: #if !defined(PETSC_USE_COMPLEX)
801:       if (dx < umin && dx >= 0.0)      dx = umin;
802:       else if (dx < 0.0 && dx > -umin) dx = -umin;
803: #else
804:       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
805:       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
806: #endif
807:       dx            *= epsilon;
808:       w3_array[col] += dx;
809:     }
810:     w3_array = w3_array + start; VecRestoreArray(w3,&w3_array);

812:     /*
813:        Evaluate function at x1 + dx (here dx is a vector of perturbations)
814:     */
815:     PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);
816:     (*f)(sctx,t,w3,w2,fctx);
817:     PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);
818:     VecAXPY(w2,-1.0,w1);

820:     /*
821:        Loop over rows of vector, putting results into Jacobian matrix
822:     */
823:     VecGetArray(w2,&y);
824:     for (l=0; l<coloring->nrows[k]; l++) {
825:       row    = coloring->rows[k][l];
826:       col    = coloring->columnsforrow[k][l];
827:       y[row] *= vscale_array[vscaleforrow[k][l]];
828:       srow   = row + start;
829:       MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);
830:     }
831:     VecRestoreArray(w2,&y);
832:   }
833:   VecRestoreArray(coloring->vscale,&vscale_array);
834:   xx    = xx + start; VecRestoreArray(x1,&xx);
835:   MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
836:   MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
837:   PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);
838:   return(0);
839: }