Actual source code: fdmatrix.c
1: /*$Id: fdmatrix.c,v 1.92 2001/08/21 21:03:06 bsmith Exp $*/
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 src/mat/matimpl.h
10: /* Logging support */
11: int MAT_FDCOLORING_COOKIE;
15: int MatFDColoringSetF(MatFDColoring fd,Vec F)
16: {
18: fd->F = F;
19: return(0);
20: }
24: static int MatFDColoringView_Draw_Zoom(PetscDraw draw,void *Aa)
25: {
26: MatFDColoring fd = (MatFDColoring)Aa;
27: int ierr,i,j;
28: PetscReal x,y;
32: /* loop over colors */
33: for (i=0; i<fd->ncolors; i++) {
34: for (j=0; j<fd->nrows[i]; j++) {
35: y = fd->M - fd->rows[i][j] - fd->rstart;
36: x = fd->columnsforrow[i][j];
37: PetscDrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);
38: }
39: }
40: return(0);
41: }
45: static int MatFDColoringView_Draw(MatFDColoring fd,PetscViewer viewer)
46: {
47: int ierr;
48: PetscTruth isnull;
49: PetscDraw draw;
50: PetscReal xr,yr,xl,yl,h,w;
53: PetscViewerDrawGetDraw(viewer,0,&draw);
54: PetscDrawIsNull(draw,&isnull); if (isnull) return(0);
56: PetscObjectCompose((PetscObject)fd,"Zoomviewer",(PetscObject)viewer);
58: xr = fd->N; yr = fd->M; h = yr/10.0; w = xr/10.0;
59: xr += w; yr += h; xl = -w; yl = -h;
60: PetscDrawSetCoordinates(draw,xl,yl,xr,yr);
61: PetscDrawZoom(draw,MatFDColoringView_Draw_Zoom,fd);
62: PetscObjectCompose((PetscObject)fd,"Zoomviewer",PETSC_NULL);
63: return(0);
64: }
68: /*@C
69: MatFDColoringView - Views a finite difference coloring context.
71: Collective on MatFDColoring
73: Input Parameters:
74: + c - the coloring context
75: - viewer - visualization context
77: Level: intermediate
79: Notes:
80: The available visualization contexts include
81: + PETSC_VIEWER_STDOUT_SELF - standard output (default)
82: . PETSC_VIEWER_STDOUT_WORLD - synchronized standard
83: output where only the first processor opens
84: the file. All other processors send their
85: data to the first processor to print.
86: - PETSC_VIEWER_DRAW_WORLD - graphical display of nonzero structure
88: Notes:
89: Since PETSc uses only a small number of basic colors (currently 33), if the coloring
90: involves more than 33 then some seemingly identical colors are displayed making it look
91: like an illegal coloring. This is just a graphical artifact.
93: .seealso: MatFDColoringCreate()
95: .keywords: Mat, finite differences, coloring, view
96: @*/
97: int MatFDColoringView(MatFDColoring c,PetscViewer viewer)
98: {
99: int i,j,ierr;
100: PetscTruth isdraw,isascii;
101: PetscViewerFormat format;
105: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(c->comm);
109: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);
110: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
111: if (isdraw) {
112: MatFDColoringView_Draw(c,viewer);
113: } else if (isascii) {
114: PetscViewerASCIIPrintf(viewer,"MatFDColoring Object:\n");
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(1,"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: 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()
166: @*/
167: int MatFDColoringSetParameters(MatFDColoring matfd,PetscReal error,PetscReal umin)
168: {
172: if (error != PETSC_DEFAULT) matfd->error_rel = error;
173: if (umin != PETSC_DEFAULT) matfd->umin = umin;
174: return(0);
175: }
179: /*@
180: MatFDColoringSetFrequency - Sets the frequency for computing new Jacobian
181: matrices.
183: Collective on MatFDColoring
185: Input Parameters:
186: + coloring - the coloring context
187: - freq - frequency (default is 1)
189: Options Database Keys:
190: . -mat_fd_coloring_freq <freq> - Sets coloring frequency
192: Level: advanced
194: Notes:
195: Using a modified Newton strategy, where the Jacobian remains fixed for several
196: iterations, can be cost effective in terms of overall nonlinear solution
197: efficiency. This parameter indicates that a new Jacobian will be computed every
198: <freq> nonlinear iterations.
200: .keywords: Mat, finite differences, coloring, set, frequency
202: .seealso: MatFDColoringCreate(), MatFDColoringGetFrequency(), MatFDColoringSetRecompute()
203: @*/
204: int MatFDColoringSetFrequency(MatFDColoring matfd,int freq)
205: {
209: matfd->freq = freq;
210: return(0);
211: }
215: /*@
216: MatFDColoringGetFrequency - Gets the frequency for computing new Jacobian
217: matrices.
219: Not Collective
221: Input Parameters:
222: . coloring - the coloring context
224: Output Parameters:
225: . freq - frequency (default is 1)
227: Options Database Keys:
228: . -mat_fd_coloring_freq <freq> - Sets coloring frequency
230: Level: advanced
232: Notes:
233: Using a modified Newton strategy, where the Jacobian remains fixed for several
234: iterations, can be cost effective in terms of overall nonlinear solution
235: efficiency. This parameter indicates that a new Jacobian will be computed every
236: <freq> nonlinear iterations.
238: .keywords: Mat, finite differences, coloring, get, frequency
240: .seealso: MatFDColoringSetFrequency()
241: @*/
242: int MatFDColoringGetFrequency(MatFDColoring matfd,int *freq)
243: {
247: *freq = matfd->freq;
248: return(0);
249: }
253: /*@C
254: MatFDColoringSetFunction - Sets the function to use for computing the Jacobian.
256: Collective on MatFDColoring
258: Input Parameters:
259: + coloring - the coloring context
260: . f - the function
261: - fctx - the optional user-defined function context
263: Level: intermediate
265: Notes:
266: In Fortran you must call MatFDColoringSetFunctionSNES() for a coloring object to
267: be used with the SNES solvers and MatFDColoringSetFunctionTS() if it is to be used
268: with the TS solvers.
270: .keywords: Mat, Jacobian, finite differences, set, function
271: @*/
272: int MatFDColoringSetFunction(MatFDColoring matfd,int (*f)(void),void *fctx)
273: {
277: matfd->f = f;
278: matfd->fctx = fctx;
280: return(0);
281: }
285: /*@
286: MatFDColoringSetFromOptions - Sets coloring finite difference parameters from
287: the options database.
289: Collective on MatFDColoring
291: The Jacobian, F'(u), is estimated with the differencing approximation
292: .vb
293: F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
294: h = error_rel*u[i] if abs(u[i]) > umin
295: = +/- error_rel*umin otherwise, with +/- determined by the sign of u[i]
296: dx_{i} = (0, ... 1, .... 0)
297: .ve
299: Input Parameter:
300: . coloring - the coloring context
302: Options Database Keys:
303: + -mat_fd_coloring_err <err> - Sets <err> (square root
304: of relative error in the function)
305: . -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude
306: . -mat_fd_coloring_freq <freq> - Sets frequency of computing a new Jacobian
307: . -mat_fd_coloring_view - Activates basic viewing
308: . -mat_fd_coloring_view_info - Activates viewing info
309: - -mat_fd_coloring_view_draw - Activates drawing
311: Level: intermediate
313: .keywords: Mat, finite differences, parameters
315: .seealso: MatFDColoringCreate(), MatFDColoringView(), MatFDColoringSetParameters()
317: @*/
318: int MatFDColoringSetFromOptions(MatFDColoring matfd)
319: {
320: int ierr;
325: PetscOptionsBegin(matfd->comm,matfd->prefix,"Jacobian computation via finite differences option","MatFD");
326: PetscOptionsReal("-mat_fd_coloring_err","Square root of relative error in function","MatFDColoringSetParameters",matfd->error_rel,&matfd->error_rel,0);
327: PetscOptionsReal("-mat_fd_coloring_umin","Minimum allowable u magnitude","MatFDColoringSetParameters",matfd->umin,&matfd->umin,0);
328: PetscOptionsInt("-mat_fd_coloring_freq","How often Jacobian is recomputed","MatFDColoringSetFrequency",matfd->freq,&matfd->freq,0);
329: /* not used here; but so they are presented in the GUI */
330: PetscOptionsName("-mat_fd_coloring_view","Print entire datastructure used for Jacobian","None",0);
331: PetscOptionsName("-mat_fd_coloring_view_info","Print number of colors etc for Jacobian","None",0);
332: PetscOptionsName("-mat_fd_coloring_view_draw","Plot nonzero structure ofJacobian","None",0);
333: PetscOptionsEnd();
334: return(0);
335: }
339: int MatFDColoringView_Private(MatFDColoring fd)
340: {
341: int ierr;
342: PetscTruth flg;
345: PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view",&flg);
346: if (flg) {
347: MatFDColoringView(fd,PETSC_VIEWER_STDOUT_(fd->comm));
348: }
349: PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_info",&flg);
350: if (flg) {
351: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_(fd->comm),PETSC_VIEWER_ASCII_INFO);
352: MatFDColoringView(fd,PETSC_VIEWER_STDOUT_(fd->comm));
353: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_(fd->comm));
354: }
355: PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg);
356: if (flg) {
357: MatFDColoringView(fd,PETSC_VIEWER_DRAW_(fd->comm));
358: PetscViewerFlush(PETSC_VIEWER_DRAW_(fd->comm));
359: }
360: return(0);
361: }
365: /*@C
366: MatFDColoringCreate - Creates a matrix coloring context for finite difference
367: computation of Jacobians.
369: Collective on Mat
371: Input Parameters:
372: + mat - the matrix containing the nonzero structure of the Jacobian
373: - iscoloring - the coloring of the matrix
375: Output Parameter:
376: . color - the new coloring context
377:
378: Level: intermediate
380: .seealso: MatFDColoringDestroy(),SNESDefaultComputeJacobianColor(), ISColoringCreate(),
381: MatFDColoringSetFunction(), MatFDColoringSetFromOptions(), MatFDColoringApply(),
382: MatFDColoringSetFrequency(), MatFDColoringSetRecompute(), MatFDColoringView(),
383: MatFDColoringSetParameters()
384: @*/
385: int MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
386: {
387: MatFDColoring c;
388: MPI_Comm comm;
389: int ierr,M,N;
392: PetscLogEventBegin(MAT_FDColoringCreate,mat,0,0,0);
393: MatGetSize(mat,&M,&N);
394: if (M != N) SETERRQ(PETSC_ERR_SUP,"Only for square matrices");
396: PetscObjectGetComm((PetscObject)mat,&comm);
397: PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_COOKIE,0,"MatFDColoring",comm,MatFDColoringDestroy,MatFDColoringView);
398: PetscLogObjectCreate(c);
400: if (mat->ops->fdcoloringcreate) {
401: (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);
402: } else {
403: SETERRQ(PETSC_ERR_SUP,"Code not yet written for this matrix type");
404: }
406: c->error_rel = PETSC_SQRT_MACHINE_EPSILON;
407: c->umin = 100.0*PETSC_SQRT_MACHINE_EPSILON;
408: c->freq = 1;
409: c->usersetsrecompute = PETSC_FALSE;
410: c->recompute = PETSC_FALSE;
411: c->currentcolor = -1;
413: *color = c;
414: PetscLogEventEnd(MAT_FDColoringCreate,mat,0,0,0);
415: return(0);
416: }
420: /*@C
421: MatFDColoringDestroy - Destroys a matrix coloring context that was created
422: via MatFDColoringCreate().
424: Collective on MatFDColoring
426: Input Parameter:
427: . c - coloring context
429: Level: intermediate
431: .seealso: MatFDColoringCreate()
432: @*/
433: int MatFDColoringDestroy(MatFDColoring c)
434: {
435: int i,ierr;
438: if (--c->refct > 0) return(0);
440: for (i=0; i<c->ncolors; i++) {
441: if (c->columns[i]) {PetscFree(c->columns[i]);}
442: if (c->rows[i]) {PetscFree(c->rows[i]);}
443: if (c->columnsforrow[i]) {PetscFree(c->columnsforrow[i]);}
444: if (c->vscaleforrow && c->vscaleforrow[i]) {PetscFree(c->vscaleforrow[i]);}
445: }
446: PetscFree(c->ncolumns);
447: PetscFree(c->columns);
448: PetscFree(c->nrows);
449: PetscFree(c->rows);
450: PetscFree(c->columnsforrow);
451: if (c->vscaleforrow) {PetscFree(c->vscaleforrow);}
452: if (c->vscale) {VecDestroy(c->vscale);}
453: if (c->w1) {
454: VecDestroy(c->w1);
455: VecDestroy(c->w2);
456: VecDestroy(c->w3);
457: }
458: PetscLogObjectDestroy(c);
459: PetscHeaderDestroy(c);
460: return(0);
461: }
465: /*@C
466: MatFDColoringGetPerturbedColumns - Returns the indices of the columns that
467: that are currently being perturbed.
469: Not Collective
471: Input Parameters:
472: . coloring - coloring context created with MatFDColoringCreate()
474: Output Parameters:
475: + n - the number of local columns being perturbed
476: - cols - the column indices, in global numbering
478: Level: intermediate
480: .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringApply()
482: .keywords: coloring, Jacobian, finite differences
483: @*/
484: EXTERN int MatFDColoringGetPerturbedColumns(MatFDColoring coloring,int *n,int *cols[])
485: {
487: if (coloring->currentcolor >= 0) {
488: *n = coloring->ncolumns[coloring->currentcolor];
489: *cols = coloring->columns[coloring->currentcolor];
490: } else {
491: *n = 0;
492: }
493: return(0);
494: }
498: /*@
499: MatFDColoringApply - Given a matrix for which a MatFDColoring context
500: has been created, computes the Jacobian for a function via finite differences.
502: Collective on MatFDColoring
504: Input Parameters:
505: + mat - location to store Jacobian
506: . coloring - coloring context created with MatFDColoringCreate()
507: . x1 - location at which Jacobian is to be computed
508: - sctx - optional context required by function (actually a SNES context)
510: Options Database Keys:
511: + -mat_fd_coloring_freq <freq> - Sets coloring frequency
512: . -mat_fd_coloring_view - Activates basic viewing or coloring
513: . -mat_fd_coloring_view_draw - Activates drawing of coloring
514: - -mat_fd_coloring_view_info - Activates viewing of coloring info
516: Level: intermediate
518: .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
520: .keywords: coloring, Jacobian, finite differences
521: @*/
522: int MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
523: {
524: int (*f)(void *,Vec,Vec,void*) = (int (*)(void *,Vec,Vec,void *))coloring->f;
525: int k,ierr,N,start,end,l,row,col,srow,**vscaleforrow,m1,m2;
526: PetscScalar dx,mone = -1.0,*y,*xx,*w3_array;
527: PetscScalar *vscale_array;
528: PetscReal epsilon = coloring->error_rel,umin = coloring->umin;
529: Vec w1,w2,w3;
530: void *fctx = coloring->fctx;
531: PetscTruth flg;
539: if (coloring->usersetsrecompute) {
540: if (!coloring->recompute) {
541: *flag = SAME_PRECONDITIONER;
542: PetscLogInfo(J,"MatFDColoringApply:Skipping Jacobian, since user called MatFDColorSetRecompute()\n");
543: return(0);
544: } else {
545: coloring->recompute = PETSC_FALSE;
546: }
547: }
549: PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);
550: if (J->ops->fdcoloringapply) {
551: (*J->ops->fdcoloringapply)(J,coloring,x1,flag,sctx);
552: } else {
554: if (!coloring->w1) {
555: VecDuplicate(x1,&coloring->w1);
556: PetscLogObjectParent(coloring,coloring->w1);
557: VecDuplicate(x1,&coloring->w2);
558: PetscLogObjectParent(coloring,coloring->w2);
559: VecDuplicate(x1,&coloring->w3);
560: PetscLogObjectParent(coloring,coloring->w3);
561: }
562: w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
564: MatSetUnfactored(J);
565: PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);
566: if (flg) {
567: PetscLogInfo(coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n");
568: } else {
569: MatZeroEntries(J);
570: }
572: VecGetOwnershipRange(x1,&start,&end);
573: VecGetSize(x1,&N);
574:
575: /*
576: This is a horrible, horrible, hack. See DMMGComputeJacobian_Multigrid() it inproperly sets
577: coloring->F for the coarser grids from the finest
578: */
579: if (coloring->F) {
580: VecGetLocalSize(coloring->F,&m1);
581: VecGetLocalSize(w1,&m2);
582: if (m1 != m2) {
583: coloring->F = 0;
584: }
585: }
587: if (coloring->F) {
588: w1 = coloring->F; /* use already computed value of function */
589: coloring->F = 0;
590: } else {
591: PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);
592: (*f)(sctx,x1,w1,fctx);
593: PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);
594: }
596: /*
597: Compute all the scale factors and share with other processors
598: */
599: VecGetArray(x1,&xx);xx = xx - start;
600: VecGetArray(coloring->vscale,&vscale_array);vscale_array = vscale_array - start;
601: for (k=0; k<coloring->ncolors; k++) {
602: /*
603: Loop over each column associated with color adding the
604: perturbation to the vector w3.
605: */
606: for (l=0; l<coloring->ncolumns[k]; l++) {
607: col = coloring->columns[k][l]; /* column of the matrix we are probing for */
608: dx = xx[col];
609: if (dx == 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] = 1.0/dx;
619: }
620: }
621: vscale_array = vscale_array + start;VecRestoreArray(coloring->vscale,&vscale_array);
622: VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);
623: VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);
625: /* VecView(coloring->vscale,PETSC_VIEWER_STDOUT_WORLD);
626: VecView(x1,PETSC_VIEWER_STDOUT_WORLD);*/
628: if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow;
629: else vscaleforrow = coloring->columnsforrow;
631: VecGetArray(coloring->vscale,&vscale_array);
632: /*
633: Loop over each color
634: */
635: for (k=0; k<coloring->ncolors; k++) {
636: coloring->currentcolor = k;
637: VecCopy(x1,w3);
638: VecGetArray(w3,&w3_array);w3_array = w3_array - start;
639: /*
640: Loop over each column associated with color adding the
641: perturbation to the vector w3.
642: */
643: for (l=0; l<coloring->ncolumns[k]; l++) {
644: col = coloring->columns[k][l]; /* column of the matrix we are probing for */
645: dx = xx[col];
646: if (dx == 0.0) dx = 1.0;
647: #if !defined(PETSC_USE_COMPLEX)
648: if (dx < umin && dx >= 0.0) dx = umin;
649: else if (dx < 0.0 && dx > -umin) dx = -umin;
650: #else
651: if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin;
652: else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
653: #endif
654: dx *= epsilon;
655: if (!PetscAbsScalar(dx)) SETERRQ(1,"Computed 0 differencing parameter");
656: w3_array[col] += dx;
657: }
658: w3_array = w3_array + start; VecRestoreArray(w3,&w3_array);
660: /*
661: Evaluate function at x1 + dx (here dx is a vector of perturbations)
662: */
664: PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);
665: (*f)(sctx,w3,w2,fctx);
666: PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);
667: VecAXPY(&mone,w1,w2);
669: /*
670: Loop over rows of vector, putting results into Jacobian matrix
671: */
672: VecGetArray(w2,&y);
673: for (l=0; l<coloring->nrows[k]; l++) {
674: row = coloring->rows[k][l];
675: col = coloring->columnsforrow[k][l];
676: y[row] *= vscale_array[vscaleforrow[k][l]];
677: srow = row + start;
678: MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);
679: }
680: VecRestoreArray(w2,&y);
681: }
682: coloring->currentcolor = -1;
683: VecRestoreArray(coloring->vscale,&vscale_array);
684: xx = xx + start; VecRestoreArray(x1,&xx);
685: MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
686: MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
687: }
688: PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);
690: PetscOptionsHasName(PETSC_NULL,"-mat_null_space_test",&flg);
691: if (flg) {
692: MatNullSpaceTest(J->nullsp,J);
693: }
694: MatFDColoringView_Private(coloring);
696: return(0);
697: }
701: /*@
702: MatFDColoringApplyTS - Given a matrix for which a MatFDColoring context
703: has been created, computes the Jacobian for a function via finite differences.
705: Collective on Mat, MatFDColoring, and Vec
707: Input Parameters:
708: + mat - location to store Jacobian
709: . coloring - coloring context created with MatFDColoringCreate()
710: . x1 - location at which Jacobian is to be computed
711: - sctx - optional context required by function (actually a SNES context)
713: Options Database Keys:
714: . -mat_fd_coloring_freq <freq> - Sets coloring frequency
716: Level: intermediate
718: .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
720: .keywords: coloring, Jacobian, finite differences
721: @*/
722: int MatFDColoringApplyTS(Mat J,MatFDColoring coloring,PetscReal t,Vec x1,MatStructure *flag,void *sctx)
723: {
724: int (*f)(void *,PetscReal,Vec,Vec,void*)=(int (*)(void *,PetscReal,Vec,Vec,void *))coloring->f;
725: int k,ierr,N,start,end,l,row,col,srow,**vscaleforrow;
726: PetscScalar dx,mone = -1.0,*y,*xx,*w3_array;
727: PetscScalar *vscale_array;
728: PetscReal epsilon = coloring->error_rel,umin = coloring->umin;
729: Vec w1,w2,w3;
730: void *fctx = coloring->fctx;
731: PetscTruth flg;
738: PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);
739: if (!coloring->w1) {
740: VecDuplicate(x1,&coloring->w1);
741: PetscLogObjectParent(coloring,coloring->w1);
742: VecDuplicate(x1,&coloring->w2);
743: PetscLogObjectParent(coloring,coloring->w2);
744: VecDuplicate(x1,&coloring->w3);
745: PetscLogObjectParent(coloring,coloring->w3);
746: }
747: w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
749: MatSetUnfactored(J);
750: PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);
751: if (flg) {
752: PetscLogInfo(coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n");
753: } else {
754: MatZeroEntries(J);
755: }
757: VecGetOwnershipRange(x1,&start,&end);
758: VecGetSize(x1,&N);
759: PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);
760: (*f)(sctx,t,x1,w1,fctx);
761: PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);
763: /*
764: Compute all the scale factors and share with other processors
765: */
766: VecGetArray(x1,&xx);xx = xx - start;
767: VecGetArray(coloring->vscale,&vscale_array);vscale_array = vscale_array - start;
768: for (k=0; k<coloring->ncolors; k++) {
769: /*
770: Loop over each column associated with color adding the
771: perturbation to the vector w3.
772: */
773: for (l=0; l<coloring->ncolumns[k]; l++) {
774: col = coloring->columns[k][l]; /* column of the matrix we are probing for */
775: dx = xx[col];
776: if (dx == 0.0) dx = 1.0;
777: #if !defined(PETSC_USE_COMPLEX)
778: if (dx < umin && dx >= 0.0) dx = umin;
779: else if (dx < 0.0 && dx > -umin) dx = -umin;
780: #else
781: if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin;
782: else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
783: #endif
784: dx *= epsilon;
785: vscale_array[col] = 1.0/dx;
786: }
787: }
788: vscale_array = vscale_array - start;VecRestoreArray(coloring->vscale,&vscale_array);
789: VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);
790: VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);
792: if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow;
793: else vscaleforrow = coloring->columnsforrow;
795: VecGetArray(coloring->vscale,&vscale_array);
796: /*
797: Loop over each color
798: */
799: for (k=0; k<coloring->ncolors; k++) {
800: VecCopy(x1,w3);
801: VecGetArray(w3,&w3_array);w3_array = w3_array - start;
802: /*
803: Loop over each column associated with color adding the
804: perturbation to the vector w3.
805: */
806: for (l=0; l<coloring->ncolumns[k]; l++) {
807: col = coloring->columns[k][l]; /* column of the matrix we are probing for */
808: dx = xx[col];
809: if (dx == 0.0) dx = 1.0;
810: #if !defined(PETSC_USE_COMPLEX)
811: if (dx < umin && dx >= 0.0) dx = umin;
812: else if (dx < 0.0 && dx > -umin) dx = -umin;
813: #else
814: if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin;
815: else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
816: #endif
817: dx *= epsilon;
818: w3_array[col] += dx;
819: }
820: w3_array = w3_array + start; VecRestoreArray(w3,&w3_array);
822: /*
823: Evaluate function at x1 + dx (here dx is a vector of perturbations)
824: */
825: PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);
826: (*f)(sctx,t,w3,w2,fctx);
827: PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);
828: VecAXPY(&mone,w1,w2);
830: /*
831: Loop over rows of vector, putting results into Jacobian matrix
832: */
833: VecGetArray(w2,&y);
834: for (l=0; l<coloring->nrows[k]; l++) {
835: row = coloring->rows[k][l];
836: col = coloring->columnsforrow[k][l];
837: y[row] *= vscale_array[vscaleforrow[k][l]];
838: srow = row + start;
839: MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);
840: }
841: VecRestoreArray(w2,&y);
842: }
843: VecRestoreArray(coloring->vscale,&vscale_array);
844: xx = xx + start; VecRestoreArray(x1,&xx);
845: MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
846: MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
847: PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);
848: return(0);
849: }
854: /*@C
855: MatFDColoringSetRecompute - Indicates that the next time a Jacobian preconditioner
856: is needed it sholuld be recomputed. Once this is called and the new Jacobian is computed
857: no additional Jacobian's will be computed (the same one will be used) until this is
858: called again.
860: Collective on MatFDColoring
862: Input Parameters:
863: . c - the coloring context
865: Level: intermediate
867: Notes: The MatFDColoringSetFrequency() is ignored once this is called
869: .seealso: MatFDColoringCreate(), MatFDColoringSetFrequency()
871: .keywords: Mat, finite differences, coloring
872: @*/
873: int MatFDColoringSetRecompute(MatFDColoring c)
874: {
877: c->usersetsrecompute = PETSC_TRUE;
878: c->recompute = PETSC_TRUE;
879: return(0);
880: }