Actual source code: fdmatrix.c
1: /*$Id: fdmatrix.c,v 1.85 2001/04/10 19:35:59 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 petsc.h
9: #include src/mat/matimpl.h
11: static int MatFDColoringView_Draw_Zoom(PetscDraw draw,void *Aa)
12: {
13: MatFDColoring fd = (MatFDColoring)Aa;
14: int ierr,i,j;
15: PetscReal x,y;
19: /* loop over colors */
20: for (i=0; i<fd->ncolors; i++) {
21: for (j=0; j<fd->nrows[i]; j++) {
22: y = fd->M - fd->rows[i][j] - fd->rstart;
23: x = fd->columnsforrow[i][j];
24: PetscDrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);
25: }
26: }
27: return(0);
28: }
30: static int MatFDColoringView_Draw(MatFDColoring fd,PetscViewer viewer)
31: {
32: int ierr;
33: PetscTruth isnull;
34: PetscDraw draw;
35: PetscReal xr,yr,xl,yl,h,w;
38: PetscViewerDrawGetDraw(viewer,0,&draw);
39: PetscDrawIsNull(draw,&isnull); if (isnull) return(0);
41: PetscObjectCompose((PetscObject)fd,"Zoomviewer",(PetscObject)viewer);
43: xr = fd->N; yr = fd->M; h = yr/10.0; w = xr/10.0;
44: xr += w; yr += h; xl = -w; yl = -h;
45: PetscDrawSetCoordinates(draw,xl,yl,xr,yr);
46: PetscDrawZoom(draw,MatFDColoringView_Draw_Zoom,fd);
47: PetscObjectCompose((PetscObject)fd,"Zoomviewer",PETSC_NULL);
48: return(0);
49: }
51: /*@C
52: MatFDColoringView - Views a finite difference coloring context.
54: Collective on MatFDColoring
56: Input Parameters:
57: + c - the coloring context
58: - viewer - visualization context
60: Level: intermediate
62: Notes:
63: The available visualization contexts include
64: + PETSC_VIEWER_STDOUT_SELF - standard output (default)
65: . PETSC_VIEWER_STDOUT_WORLD - synchronized standard
66: output where only the first processor opens
67: the file. All other processors send their
68: data to the first processor to print.
69: - PETSC_VIEWER_DRAW_WORLD - graphical display of nonzero structure
71: Notes:
72: Since PETSc uses only a small number of basic colors (currently 33), if the coloring
73: involves more than 33 then some seemingly identical colors are displayed making it look
74: like an illegal coloring. This is just a graphical artifact.
76: .seealso: MatFDColoringCreate()
78: .keywords: Mat, finite differences, coloring, view
79: @*/
80: int MatFDColoringView(MatFDColoring c,PetscViewer viewer)
81: {
82: int i,j,ierr;
83: PetscTruth isdraw,isascii;
84: PetscViewerFormat format;
88: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(c->comm);
92: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);
93: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
94: if (isdraw) {
95: MatFDColoringView_Draw(c,viewer);
96: } else if (isascii) {
97: PetscViewerASCIIPrintf(viewer,"MatFDColoring Object:n");
98: PetscViewerASCIIPrintf(viewer," Error tolerance=%gn",c->error_rel);
99: PetscViewerASCIIPrintf(viewer," Umin=%gn",c->umin);
100: PetscViewerASCIIPrintf(viewer," Number of colors=%dn",c->ncolors);
102: PetscViewerGetFormat(viewer,&format);
103: if (format != PETSC_VIEWER_ASCII_INFO) {
104: for (i=0; i<c->ncolors; i++) {
105: PetscViewerASCIIPrintf(viewer," Information for color %dn",i);
106: PetscViewerASCIIPrintf(viewer," Number of columns %dn",c->ncolumns[i]);
107: for (j=0; j<c->ncolumns[i]; j++) {
108: PetscViewerASCIIPrintf(viewer," %dn",c->columns[i][j]);
109: }
110: PetscViewerASCIIPrintf(viewer," Number of rows %dn",c->nrows[i]);
111: for (j=0; j<c->nrows[i]; j++) {
112: PetscViewerASCIIPrintf(viewer," %d %d n",c->rows[i][j],c->columnsforrow[i][j]);
113: }
114: }
115: }
116: PetscViewerFlush(viewer);
117: } else {
118: SETERRQ1(1,"Viewer type %s not supported for MatFDColoring",((PetscObject)viewer)->type_name);
119: }
120: return(0);
121: }
123: /*@
124: MatFDColoringSetParameters - Sets the parameters for the sparse approximation of
125: a Jacobian matrix using finite differences.
127: Collective on MatFDColoring
129: The Jacobian is estimated with the differencing approximation
130: .vb
131: F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
132: h = error_rel*u[i] if abs(u[i]) > umin
133: = +/- error_rel*umin otherwise, with +/- determined by the sign of u[i]
134: dx_{i} = (0, ... 1, .... 0)
135: .ve
137: Input Parameters:
138: + coloring - the coloring context
139: . error_rel - relative error
140: - umin - minimum allowable u-value magnitude
142: Level: advanced
144: .keywords: Mat, finite differences, coloring, set, parameters
146: .seealso: MatFDColoringCreate()
147: @*/
148: int MatFDColoringSetParameters(MatFDColoring matfd,PetscReal error,PetscReal umin)
149: {
153: if (error != PETSC_DEFAULT) matfd->error_rel = error;
154: if (umin != PETSC_DEFAULT) matfd->umin = umin;
155: return(0);
156: }
158: /*@
159: MatFDColoringSetFrequency - Sets the frequency for computing new Jacobian
160: matrices.
162: Collective on MatFDColoring
164: Input Parameters:
165: + coloring - the coloring context
166: - freq - frequency (default is 1)
168: Options Database Keys:
169: . -mat_fd_coloring_freq <freq> - Sets coloring frequency
171: Level: advanced
173: Notes:
174: Using a modified Newton strategy, where the Jacobian remains fixed for several
175: iterations, can be cost effective in terms of overall nonlinear solution
176: efficiency. This parameter indicates that a new Jacobian will be computed every
177: <freq> nonlinear iterations.
179: .keywords: Mat, finite differences, coloring, set, frequency
181: .seealso: MatFDColoringCreate(), MatFDColoringGetFrequency(), MatFDColoringSetRecompute()
182: @*/
183: int MatFDColoringSetFrequency(MatFDColoring matfd,int freq)
184: {
188: matfd->freq = freq;
189: return(0);
190: }
192: /*@
193: MatFDColoringGetFrequency - Gets the frequency for computing new Jacobian
194: matrices.
196: Not Collective
198: Input Parameters:
199: . coloring - the coloring context
201: Output Parameters:
202: . freq - frequency (default is 1)
204: Options Database Keys:
205: . -mat_fd_coloring_freq <freq> - Sets coloring frequency
207: Level: advanced
209: Notes:
210: Using a modified Newton strategy, where the Jacobian remains fixed for several
211: iterations, can be cost effective in terms of overall nonlinear solution
212: efficiency. This parameter indicates that a new Jacobian will be computed every
213: <freq> nonlinear iterations.
215: .keywords: Mat, finite differences, coloring, get, frequency
217: .seealso: MatFDColoringSetFrequency()
218: @*/
219: int MatFDColoringGetFrequency(MatFDColoring matfd,int *freq)
220: {
224: *freq = matfd->freq;
225: return(0);
226: }
228: /*@C
229: MatFDColoringSetFunction - Sets the function to use for computing the Jacobian.
231: Collective on MatFDColoring
233: Input Parameters:
234: + coloring - the coloring context
235: . f - the function
236: - fctx - the optional user-defined function context
238: Level: intermediate
240: Notes:
241: In Fortran you must call MatFDColoringSetFunctionSNES() for a coloring object to
242: be used with the SNES solvers and MatFDColoringSetFunctionTS() if it is to be used
243: with the TS solvers.
245: .keywords: Mat, Jacobian, finite differences, set, function
246: @*/
247: int MatFDColoringSetFunction(MatFDColoring matfd,int (*f)(void),void *fctx)
248: {
252: matfd->f = f;
253: matfd->fctx = fctx;
255: return(0);
256: }
258: /*@
259: MatFDColoringSetFromOptions - Sets coloring finite difference parameters from
260: the options database.
262: Collective on MatFDColoring
264: The Jacobian, F'(u), is estimated with the differencing approximation
265: .vb
266: F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
267: h = error_rel*u[i] if abs(u[i]) > umin
268: = +/- error_rel*umin otherwise, with +/- determined by the sign of u[i]
269: dx_{i} = (0, ... 1, .... 0)
270: .ve
272: Input Parameter:
273: . coloring - the coloring context
275: Options Database Keys:
276: + -mat_fd_coloring_err <err> - Sets <err> (square root
277: of relative error in the function)
278: . -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude
279: . -mat_fd_coloring_freq <freq> - Sets frequency of computing a new Jacobian
280: . -mat_fd_coloring_view - Activates basic viewing
281: . -mat_fd_coloring_view_info - Activates viewing info
282: - -mat_fd_coloring_view_draw - Activates drawing
284: Level: intermediate
286: .keywords: Mat, finite differences, parameters
287: @*/
288: int MatFDColoringSetFromOptions(MatFDColoring matfd)
289: {
290: int ierr;
295: PetscOptionsBegin(matfd->comm,matfd->prefix,"Jacobian computation via finite differences option","MatFD");
296: PetscOptionsDouble("-mat_fd_coloring_err","Square root of relative error in function","MatFDColoringSetParameters",matfd->error_rel,&matfd->error_rel,0);
297: PetscOptionsDouble("-mat_fd_coloring_umin","Minimum allowable u magnitude","MatFDColoringSetParameters",matfd->umin,&matfd->umin,0);
298: PetscOptionsInt("-mat_fd_coloring_freq","How often Jacobian is recomputed","MatFDColoringSetFrequency",matfd->freq,&matfd->freq,0);
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: }
307: int MatFDColoringView_Private(MatFDColoring fd)
308: {
309: int ierr;
310: PetscTruth flg;
313: PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view",&flg);
314: if (flg) {
315: MatFDColoringView(fd,PETSC_VIEWER_STDOUT_(fd->comm));
316: }
317: PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_info",&flg);
318: if (flg) {
319: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_(fd->comm),PETSC_VIEWER_ASCII_INFO);
320: MatFDColoringView(fd,PETSC_VIEWER_STDOUT_(fd->comm));
321: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_(fd->comm));
322: }
323: PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg);
324: if (flg) {
325: MatFDColoringView(fd,PETSC_VIEWER_DRAW_(fd->comm));
326: PetscViewerFlush(PETSC_VIEWER_DRAW_(fd->comm));
327: }
328: return(0);
329: }
331: /*@C
332: MatFDColoringCreate - Creates a matrix coloring context for finite difference
333: computation of Jacobians.
335: Collective on Mat
337: Input Parameters:
338: + mat - the matrix containing the nonzero structure of the Jacobian
339: - iscoloring - the coloring of the matrix
341: Output Parameter:
342: . color - the new coloring context
343:
344: Options Database Keys:
345: + -mat_fd_coloring_view - Activates basic viewing or coloring
346: . -mat_fd_coloring_view_draw - Activates drawing of coloring
347: - -mat_fd_coloring_view_info - Activates viewing of coloring info
349: Level: intermediate
351: .seealso: MatFDColoringDestroy(),SNESDefaultComputeJacobianColor(), ISColoringCreate(),
352: MatFDColoringSetFunction(), MatFDColoringSetFromOptions()
353: @*/
354: int MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
355: {
356: MatFDColoring c;
357: MPI_Comm comm;
358: int ierr,M,N;
361: PetscLogEventBegin(MAT_FDColoringCreate,mat,0,0,0);
362: MatGetSize(mat,&M,&N);
363: if (M != N) SETERRQ(PETSC_ERR_SUP,"Only for square matrices");
365: PetscObjectGetComm((PetscObject)mat,&comm);
366: PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_COOKIE,0,"MatFDColoring",comm,MatFDColoringDestroy,MatFDColoringView);
367: PetscLogObjectCreate(c);
369: if (mat->ops->fdcoloringcreate) {
370: (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);
371: } else {
372: SETERRQ(PETSC_ERR_SUP,"Code not yet written for this matrix type");
373: }
375: c->error_rel = 1.e-8;
376: c->umin = 1.e-6;
377: c->freq = 1;
378: c->usersetsrecompute = PETSC_FALSE;
379: c->recompute = PETSC_FALSE;
381: MatFDColoringView_Private(c);
383: *color = c;
384: PetscLogEventEnd(MAT_FDColoringCreate,mat,0,0,0);
385: return(0);
386: }
388: /*@C
389: MatFDColoringDestroy - Destroys a matrix coloring context that was created
390: via MatFDColoringCreate().
392: Collective on MatFDColoring
394: Input Parameter:
395: . c - coloring context
397: Level: intermediate
399: .seealso: MatFDColoringCreate()
400: @*/
401: int MatFDColoringDestroy(MatFDColoring c)
402: {
403: int i,ierr;
406: if (--c->refct > 0) return(0);
408: for (i=0; i<c->ncolors; i++) {
409: if (c->columns[i]) {PetscFree(c->columns[i]);}
410: if (c->rows[i]) {PetscFree(c->rows[i]);}
411: if (c->columnsforrow[i]) {PetscFree(c->columnsforrow[i]);}
412: if (c->vscaleforrow && c->vscaleforrow[i]) {PetscFree(c->vscaleforrow[i]);}
413: }
414: PetscFree(c->ncolumns);
415: PetscFree(c->columns);
416: PetscFree(c->nrows);
417: PetscFree(c->rows);
418: PetscFree(c->columnsforrow);
419: if (c->vscaleforrow) {PetscFree(c->vscaleforrow);}
420: if (c->vscale) {VecDestroy(c->vscale);}
421: if (c->w1) {
422: VecDestroy(c->w1);
423: VecDestroy(c->w2);
424: VecDestroy(c->w3);
425: }
426: PetscLogObjectDestroy(c);
427: PetscHeaderDestroy(c);
428: return(0);
429: }
431: /*@
432: MatFDColoringApply - Given a matrix for which a MatFDColoring context
433: has been created, computes the Jacobian for a function via finite differences.
435: Collective on MatFDColoring
437: Input Parameters:
438: + mat - location to store Jacobian
439: . coloring - coloring context created with MatFDColoringCreate()
440: . x1 - location at which Jacobian is to be computed
441: - sctx - optional context required by function (actually a SNES context)
443: Options Database Keys:
444: . -mat_fd_coloring_freq <freq> - Sets coloring frequency
446: Level: intermediate
448: .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
450: .keywords: coloring, Jacobian, finite differences
451: @*/
452: int MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
453: {
454: int (*f)(void *,Vec,Vec,void*) = (int (*)(void *,Vec,Vec,void *))coloring->f;
455: int k,ierr,N,start,end,l,row,col,srow,**vscaleforrow;
456: Scalar dx,mone = -1.0,*y,*xx,*w3_array;
457: Scalar *vscale_array;
458: PetscReal epsilon = coloring->error_rel,umin = coloring->umin;
459: Vec w1,w2,w3;
460: void *fctx = coloring->fctx;
461: PetscTruth flg;
469: if (coloring->usersetsrecompute) {
470: if (!coloring->recompute) {
471: *flag = SAME_PRECONDITIONER;
472: return(0);
473: } else {
474: coloring->recompute = PETSC_FALSE;
475: }
476: }
478: PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);
479: if (!coloring->w1) {
480: VecDuplicate(x1,&coloring->w1);
481: PetscLogObjectParent(coloring,coloring->w1);
482: VecDuplicate(x1,&coloring->w2);
483: PetscLogObjectParent(coloring,coloring->w2);
484: VecDuplicate(x1,&coloring->w3);
485: PetscLogObjectParent(coloring,coloring->w3);
486: }
487: w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
489: MatSetUnfactored(J);
490: PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);
491: if (flg) {
492: PetscLogInfo(coloring,"MatFDColoringApply: Not calling MatZeroEntries()n");
493: } else {
494: MatZeroEntries(J);
495: }
497: VecGetOwnershipRange(x1,&start,&end);
498: VecGetSize(x1,&N);
499: (*f)(sctx,x1,w1,fctx);
501: /*
502: Compute all the scale factors and share with other processors
503: */
504: VecGetArray(x1,&xx);xx = xx - start;
505: VecGetArray(coloring->vscale,&vscale_array);vscale_array = vscale_array - start;
506: for (k=0; k<coloring->ncolors; k++) {
507: /*
508: Loop over each column associated with color adding the
509: perturbation to the vector w3.
510: */
511: for (l=0; l<coloring->ncolumns[k]; l++) {
512: col = coloring->columns[k][l]; /* column of the matrix we are probing for */
513: dx = xx[col];
514: if (dx == 0.0) dx = 1.0;
515: #if !defined(PETSC_USE_COMPLEX)
516: if (dx < umin && dx >= 0.0) dx = umin;
517: else if (dx < 0.0 && dx > -umin) dx = -umin;
518: #else
519: if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin;
520: else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
521: #endif
522: dx *= epsilon;
523: vscale_array[col] = 1.0/dx;
524: }
525: }
526: vscale_array = vscale_array + start;VecRestoreArray(coloring->vscale,&vscale_array);
527: VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);
528: VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);
530: /* VecView(coloring->vscale,PETSC_VIEWER_STDOUT_WORLD);
531: VecView(x1,PETSC_VIEWER_STDOUT_WORLD);*/
533: if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow;
534: else vscaleforrow = coloring->columnsforrow;
536: VecGetArray(coloring->vscale,&vscale_array);
537: /*
538: Loop over each color
539: */
540: for (k=0; k<coloring->ncolors; k++) {
541: VecCopy(x1,w3);
542: VecGetArray(w3,&w3_array);w3_array = w3_array - start;
543: /*
544: Loop over each column associated with color adding the
545: perturbation to the vector w3.
546: */
547: for (l=0; l<coloring->ncolumns[k]; l++) {
548: col = coloring->columns[k][l]; /* column of the matrix we are probing for */
549: dx = xx[col];
550: if (dx == 0.0) dx = 1.0;
551: #if !defined(PETSC_USE_COMPLEX)
552: if (dx < umin && dx >= 0.0) dx = umin;
553: else if (dx < 0.0 && dx > -umin) dx = -umin;
554: #else
555: if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin;
556: else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
557: #endif
558: dx *= epsilon;
559: if (!PetscAbsScalar(dx)) SETERRQ(1,"Computed 0 differencing parameter");
560: w3_array[col] += dx;
561: }
562: w3_array = w3_array + start; VecRestoreArray(w3,&w3_array);
564: /*
565: Evaluate function at x1 + dx (here dx is a vector of perturbations)
566: */
568: (*f)(sctx,w3,w2,fctx);
569: VecAXPY(&mone,w1,w2);
571: /*
572: Loop over rows of vector, putting results into Jacobian matrix
573: */
574: VecGetArray(w2,&y);
575: for (l=0; l<coloring->nrows[k]; l++) {
576: row = coloring->rows[k][l];
577: col = coloring->columnsforrow[k][l];
578: y[row] *= vscale_array[vscaleforrow[k][l]];
579: srow = row + start;
580: ierr = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);
581: }
582: VecRestoreArray(w2,&y);
583: }
584: VecRestoreArray(coloring->vscale,&vscale_array);
585: xx = xx + start; ierr = VecRestoreArray(x1,&xx);
586: ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
587: ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
588: PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);
590: PetscOptionsHasName(PETSC_NULL,"-mat_null_space_test",&flg);
591: if (flg) {
592: MatNullSpaceTest(J->nullsp,J);
593: }
594: return(0);
595: }
597: /*@
598: MatFDColoringApplyTS - Given a matrix for which a MatFDColoring context
599: has been created, computes the Jacobian for a function via finite differences.
601: Collective on Mat, MatFDColoring, and Vec
603: Input Parameters:
604: + mat - location to store Jacobian
605: . coloring - coloring context created with MatFDColoringCreate()
606: . x1 - location at which Jacobian is to be computed
607: - sctx - optional context required by function (actually a SNES context)
609: Options Database Keys:
610: . -mat_fd_coloring_freq <freq> - Sets coloring frequency
612: Level: intermediate
614: .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
616: .keywords: coloring, Jacobian, finite differences
617: @*/
618: int MatFDColoringApplyTS(Mat J,MatFDColoring coloring,PetscReal t,Vec x1,MatStructure *flag,void *sctx)
619: {
620: int (*f)(void *,PetscReal,Vec,Vec,void*)=(int (*)(void *,PetscReal,Vec,Vec,void *))coloring->f;
621: int k,ierr,N,start,end,l,row,col,srow,**vscaleforrow;
622: Scalar dx,mone = -1.0,*y,*xx,*w3_array;
623: Scalar *vscale_array;
624: PetscReal epsilon = coloring->error_rel,umin = coloring->umin;
625: Vec w1,w2,w3;
626: void *fctx = coloring->fctx;
627: PetscTruth flg;
634: PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);
635: if (!coloring->w1) {
636: VecDuplicate(x1,&coloring->w1);
637: PetscLogObjectParent(coloring,coloring->w1);
638: VecDuplicate(x1,&coloring->w2);
639: PetscLogObjectParent(coloring,coloring->w2);
640: VecDuplicate(x1,&coloring->w3);
641: PetscLogObjectParent(coloring,coloring->w3);
642: }
643: w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
645: MatSetUnfactored(J);
646: PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);
647: if (flg) {
648: PetscLogInfo(coloring,"MatFDColoringApply: Not calling MatZeroEntries()n");
649: } else {
650: MatZeroEntries(J);
651: }
653: VecGetOwnershipRange(x1,&start,&end);
654: VecGetSize(x1,&N);
655: (*f)(sctx,t,x1,w1,fctx);
657: /*
658: Compute all the scale factors and share with other processors
659: */
660: VecGetArray(x1,&xx);xx = xx - start;
661: VecGetArray(coloring->vscale,&vscale_array);vscale_array = vscale_array - start;
662: for (k=0; k<coloring->ncolors; k++) {
663: /*
664: Loop over each column associated with color adding the
665: perturbation to the vector w3.
666: */
667: for (l=0; l<coloring->ncolumns[k]; l++) {
668: col = coloring->columns[k][l]; /* column of the matrix we are probing for */
669: dx = xx[col];
670: if (dx == 0.0) dx = 1.0;
671: #if !defined(PETSC_USE_COMPLEX)
672: if (dx < umin && dx >= 0.0) dx = umin;
673: else if (dx < 0.0 && dx > -umin) dx = -umin;
674: #else
675: if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin;
676: else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
677: #endif
678: dx *= epsilon;
679: vscale_array[col] = 1.0/dx;
680: }
681: }
682: vscale_array = vscale_array - start;VecRestoreArray(coloring->vscale,&vscale_array);
683: VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);
684: VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);
686: if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow;
687: else vscaleforrow = coloring->columnsforrow;
689: VecGetArray(coloring->vscale,&vscale_array);
690: /*
691: Loop over each color
692: */
693: for (k=0; k<coloring->ncolors; k++) {
694: VecCopy(x1,w3);
695: VecGetArray(w3,&w3_array);w3_array = w3_array - start;
696: /*
697: Loop over each column associated with color adding the
698: perturbation to the vector w3.
699: */
700: for (l=0; l<coloring->ncolumns[k]; l++) {
701: col = coloring->columns[k][l]; /* column of the matrix we are probing for */
702: dx = xx[col];
703: if (dx == 0.0) dx = 1.0;
704: #if !defined(PETSC_USE_COMPLEX)
705: if (dx < umin && dx >= 0.0) dx = umin;
706: else if (dx < 0.0 && dx > -umin) dx = -umin;
707: #else
708: if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin;
709: else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
710: #endif
711: dx *= epsilon;
712: w3_array[col] += dx;
713: }
714: w3_array = w3_array + start; VecRestoreArray(w3,&w3_array);
716: /*
717: Evaluate function at x1 + dx (here dx is a vector of perturbations)
718: */
719: (*f)(sctx,t,w3,w2,fctx);
720: VecAXPY(&mone,w1,w2);
722: /*
723: Loop over rows of vector, putting results into Jacobian matrix
724: */
725: VecGetArray(w2,&y);
726: for (l=0; l<coloring->nrows[k]; l++) {
727: row = coloring->rows[k][l];
728: col = coloring->columnsforrow[k][l];
729: y[row] *= vscale_array[vscaleforrow[k][l]];
730: srow = row + start;
731: ierr = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);
732: }
733: VecRestoreArray(w2,&y);
734: }
735: ierr = VecRestoreArray(coloring->vscale,&vscale_array);
736: xx = xx + start; ierr = VecRestoreArray(x1,&xx);
737: ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
738: ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
739: ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);
740: return(0);
741: }
744: /*@C
745: MatFDColoringSetRecompute - Indicates that the next time a Jacobian preconditioner
746: is needed it sholuld be recomputed. Once this is called and the new Jacobian is computed
747: no additional Jacobian's will be computed (the same one will be used) until this is
748: called again.
750: Collective on MatFDColoring
752: Input Parameters:
753: . c - the coloring context
755: Level: intermediate
757: Notes: The MatFDColoringSetFrequency() is ignored once this is called
759: .seealso: MatFDColoringCreate(), MatFDColoringSetFrequency()
761: .keywords: Mat, finite differences, coloring
762: @*/
763: int MatFDColoringSetRecompute(MatFDColoring c)
764: {
767: c->usersetsrecompute = PETSC_TRUE;
768: c->recompute = PETSC_TRUE;
769: return(0);
770: }