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 src/mat/matimpl.h
10: /* Logging support */
11: PetscCookie MAT_FDCOLORING_COOKIE = 0;
15: PetscErrorCode MatFDColoringSetF(MatFDColoring fd,Vec F)
16: {
18: fd->F = F;
19: return(0);
20: }
24: static PetscErrorCode MatFDColoringView_Draw_Zoom(PetscDraw draw,void *Aa)
25: {
26: MatFDColoring fd = (MatFDColoring)Aa;
28: PetscInt i,j;
29: PetscReal x,y;
33: /* loop over colors */
34: for (i=0; i<fd->ncolors; i++) {
35: for (j=0; j<fd->nrows[i]; j++) {
36: y = fd->M - fd->rows[i][j] - fd->rstart;
37: x = fd->columnsforrow[i][j];
38: PetscDrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);
39: }
40: }
41: return(0);
42: }
46: static PetscErrorCode MatFDColoringView_Draw(MatFDColoring fd,PetscViewer viewer)
47: {
49: PetscTruth isnull;
50: PetscDraw draw;
51: PetscReal xr,yr,xl,yl,h,w;
54: PetscViewerDrawGetDraw(viewer,0,&draw);
55: PetscDrawIsNull(draw,&isnull); if (isnull) return(0);
57: PetscObjectCompose((PetscObject)fd,"Zoomviewer",(PetscObject)viewer);
59: xr = fd->N; yr = fd->M; h = yr/10.0; w = xr/10.0;
60: xr += w; yr += h; xl = -w; yl = -h;
61: PetscDrawSetCoordinates(draw,xl,yl,xr,yr);
62: PetscDrawZoom(draw,MatFDColoringView_Draw_Zoom,fd);
63: PetscObjectCompose((PetscObject)fd,"Zoomviewer",PETSC_NULL);
64: return(0);
65: }
69: /*@C
70: MatFDColoringView - Views a finite difference coloring context.
72: Collective on MatFDColoring
74: Input Parameters:
75: + c - the coloring context
76: - viewer - visualization context
78: Level: intermediate
80: Notes:
81: The available visualization contexts include
82: + PETSC_VIEWER_STDOUT_SELF - standard output (default)
83: . PETSC_VIEWER_STDOUT_WORLD - synchronized standard
84: output where only the first processor opens
85: the file. All other processors send their
86: data to the first processor to print.
87: - PETSC_VIEWER_DRAW_WORLD - graphical display of nonzero structure
89: Notes:
90: Since PETSc uses only a small number of basic colors (currently 33), if the coloring
91: involves more than 33 then some seemingly identical colors are displayed making it look
92: like an illegal coloring. This is just a graphical artifact.
94: .seealso: MatFDColoringCreate()
96: .keywords: Mat, finite differences, coloring, view
97: @*/
98: PetscErrorCode MatFDColoringView(MatFDColoring c,PetscViewer viewer)
99: {
100: PetscErrorCode ierr;
101: PetscInt i,j;
102: PetscTruth isdraw,iascii;
103: PetscViewerFormat format;
107: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(c->comm);
111: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);
112: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
113: if (isdraw) {
114: MatFDColoringView_Draw(c,viewer);
115: } else if (iascii) {
116: PetscViewerASCIIPrintf(viewer,"MatFDColoring Object:\n");
117: PetscViewerASCIIPrintf(viewer," Error tolerance=%G\n",c->error_rel);
118: PetscViewerASCIIPrintf(viewer," Umin=%G\n",c->umin);
119: PetscViewerASCIIPrintf(viewer," Number of colors=%D\n",c->ncolors);
121: PetscViewerGetFormat(viewer,&format);
122: if (format != PETSC_VIEWER_ASCII_INFO) {
123: for (i=0; i<c->ncolors; i++) {
124: PetscViewerASCIIPrintf(viewer," Information for color %D\n",i);
125: PetscViewerASCIIPrintf(viewer," Number of columns %D\n",c->ncolumns[i]);
126: for (j=0; j<c->ncolumns[i]; j++) {
127: PetscViewerASCIIPrintf(viewer," %D\n",c->columns[i][j]);
128: }
129: PetscViewerASCIIPrintf(viewer," Number of rows %D\n",c->nrows[i]);
130: for (j=0; j<c->nrows[i]; j++) {
131: PetscViewerASCIIPrintf(viewer," %D %D \n",c->rows[i][j],c->columnsforrow[i][j]);
132: }
133: }
134: }
135: PetscViewerFlush(viewer);
136: } else {
137: SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for MatFDColoring",((PetscObject)viewer)->type_name);
138: }
139: return(0);
140: }
144: /*@
145: MatFDColoringSetParameters - Sets the parameters for the sparse approximation of
146: a Jacobian matrix using finite differences.
148: Collective on MatFDColoring
150: The Jacobian is estimated with the differencing approximation
151: .vb
152: F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
153: h = error_rel*u[i] if abs(u[i]) > umin
154: = +/- error_rel*umin otherwise, with +/- determined by the sign of u[i]
155: dx_{i} = (0, ... 1, .... 0)
156: .ve
158: Input Parameters:
159: + coloring - the coloring context
160: . error_rel - relative error
161: - umin - minimum allowable u-value magnitude
163: Level: advanced
165: .keywords: Mat, finite differences, coloring, set, parameters
167: .seealso: MatFDColoringCreate()
168: @*/
169: PetscErrorCode MatFDColoringSetParameters(MatFDColoring matfd,PetscReal error,PetscReal umin)
170: {
174: if (error != PETSC_DEFAULT) matfd->error_rel = error;
175: if (umin != PETSC_DEFAULT) matfd->umin = umin;
176: return(0);
177: }
181: /*@
182: MatFDColoringSetFrequency - Sets the frequency for computing new Jacobian
183: matrices.
185: Collective on MatFDColoring
187: Input Parameters:
188: + coloring - the coloring context
189: - freq - frequency (default is 1)
191: Options Database Keys:
192: . -mat_fd_coloring_freq <freq> - Sets coloring frequency
194: Level: advanced
196: Notes:
197: Using a modified Newton strategy, where the Jacobian remains fixed for several
198: iterations, can be cost effective in terms of overall nonlinear solution
199: efficiency. This parameter indicates that a new Jacobian will be computed every
200: <freq> nonlinear iterations.
202: .keywords: Mat, finite differences, coloring, set, frequency
204: .seealso: MatFDColoringCreate(), MatFDColoringGetFrequency(), MatFDColoringSetRecompute()
205: @*/
206: PetscErrorCode MatFDColoringSetFrequency(MatFDColoring matfd,PetscInt freq)
207: {
211: matfd->freq = freq;
212: return(0);
213: }
217: /*@
218: MatFDColoringGetFrequency - Gets the frequency for computing new Jacobian
219: matrices.
221: Not Collective
223: Input Parameters:
224: . coloring - the coloring context
226: Output Parameters:
227: . freq - frequency (default is 1)
229: Options Database Keys:
230: . -mat_fd_coloring_freq <freq> - Sets coloring frequency
232: Level: advanced
234: Notes:
235: Using a modified Newton strategy, where the Jacobian remains fixed for several
236: iterations, can be cost effective in terms of overall nonlinear solution
237: efficiency. This parameter indicates that a new Jacobian will be computed every
238: <freq> nonlinear iterations.
240: .keywords: Mat, finite differences, coloring, get, frequency
242: .seealso: MatFDColoringSetFrequency()
243: @*/
244: PetscErrorCode MatFDColoringGetFrequency(MatFDColoring matfd,PetscInt *freq)
245: {
248: *freq = matfd->freq;
249: return(0);
250: }
254: /*@C
255: MatFDColoringGetFunction - Gets the function to use for computing the Jacobian.
257: Collective on MatFDColoring
259: Input Parameters:
260: . coloring - the coloring context
262: Output Parameters:
263: + f - the function
264: - fctx - the optional user-defined function context
266: Level: intermediate
268: .keywords: Mat, Jacobian, finite differences, set, function
269: @*/
270: PetscErrorCode MatFDColoringGetFunction(MatFDColoring matfd,PetscErrorCode (**f)(void),void **fctx)
271: {
274: if (f) *f = matfd->f;
275: if (fctx) *fctx = matfd->fctx;
276: return(0);
277: }
281: /*@C
282: MatFDColoringSetFunction - Sets the function to use for computing the Jacobian.
284: Collective on MatFDColoring
286: Input Parameters:
287: + coloring - the coloring context
288: . f - the function
289: - fctx - the optional user-defined function context
291: Level: intermediate
293: Notes:
294: In Fortran you must call MatFDColoringSetFunctionSNES() for a coloring object to
295: be used with the SNES solvers and MatFDColoringSetFunctionTS() if it is to be used
296: with the TS solvers.
298: .keywords: Mat, Jacobian, finite differences, set, function
299: @*/
300: PetscErrorCode MatFDColoringSetFunction(MatFDColoring matfd,PetscErrorCode (*f)(void),void *fctx)
301: {
304: matfd->f = f;
305: matfd->fctx = fctx;
306: return(0);
307: }
311: /*@
312: MatFDColoringSetFromOptions - Sets coloring finite difference parameters from
313: the options database.
315: Collective on MatFDColoring
317: The Jacobian, F'(u), is estimated with the differencing approximation
318: .vb
319: F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
320: h = error_rel*u[i] if abs(u[i]) > umin
321: = +/- error_rel*umin otherwise, with +/- determined by the sign of u[i]
322: dx_{i} = (0, ... 1, .... 0)
323: .ve
325: Input Parameter:
326: . coloring - the coloring context
328: Options Database Keys:
329: + -mat_fd_coloring_err <err> - Sets <err> (square root
330: of relative error in the function)
331: . -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude
332: . -mat_fd_coloring_freq <freq> - Sets frequency of computing a new Jacobian
333: . -mat_fd_type - "wp" or "ds" (see MATSNESMF_WP or MATSNESMF_DS)
334: . -mat_fd_coloring_view - Activates basic viewing
335: . -mat_fd_coloring_view_info - Activates viewing info
336: - -mat_fd_coloring_view_draw - Activates drawing
338: Level: intermediate
340: .keywords: Mat, finite differences, parameters
342: .seealso: MatFDColoringCreate(), MatFDColoringView(), MatFDColoringSetParameters()
344: @*/
345: PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring matfd)
346: {
348: PetscTruth flg;
349: char value[3];
354: PetscOptionsBegin(matfd->comm,matfd->prefix,"Jacobian computation via finite differences option","MatFD");
355: PetscOptionsReal("-mat_fd_coloring_err","Square root of relative error in function","MatFDColoringSetParameters",matfd->error_rel,&matfd->error_rel,0);
356: PetscOptionsReal("-mat_fd_coloring_umin","Minimum allowable u magnitude","MatFDColoringSetParameters",matfd->umin,&matfd->umin,0);
357: PetscOptionsInt("-mat_fd_coloring_freq","How often Jacobian is recomputed","MatFDColoringSetFrequency",matfd->freq,&matfd->freq,0);
358: PetscOptionsString("-mat_fd_type","Algorithm to compute h, wp or ds","MatFDColoringCreate",matfd->htype,value,2,&flg);
359: if (flg) {
360: if (value[0] == 'w' && value[1] == 'p') matfd->htype = "wp";
361: else if (value[0] == 'd' && value[1] == 's') matfd->htype = "ds";
362: else SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Unknown finite differencing type %s",value);
363: }
364: /* not used here; but so they are presented in the GUI */
365: PetscOptionsName("-mat_fd_coloring_view","Print entire datastructure used for Jacobian","None",0);
366: PetscOptionsName("-mat_fd_coloring_view_info","Print number of colors etc for Jacobian","None",0);
367: PetscOptionsName("-mat_fd_coloring_view_draw","Plot nonzero structure ofJacobian","None",0);
368: PetscOptionsEnd();
369: return(0);
370: }
374: PetscErrorCode MatFDColoringView_Private(MatFDColoring fd)
375: {
377: PetscTruth flg;
380: PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view",&flg);
381: if (flg) {
382: MatFDColoringView(fd,PETSC_VIEWER_STDOUT_(fd->comm));
383: }
384: PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_info",&flg);
385: if (flg) {
386: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_(fd->comm),PETSC_VIEWER_ASCII_INFO);
387: MatFDColoringView(fd,PETSC_VIEWER_STDOUT_(fd->comm));
388: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_(fd->comm));
389: }
390: PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg);
391: if (flg) {
392: MatFDColoringView(fd,PETSC_VIEWER_DRAW_(fd->comm));
393: PetscViewerFlush(PETSC_VIEWER_DRAW_(fd->comm));
394: }
395: return(0);
396: }
400: /*@
401: MatFDColoringCreate - Creates a matrix coloring context for finite difference
402: computation of Jacobians.
404: Collective on Mat
406: Input Parameters:
407: + mat - the matrix containing the nonzero structure of the Jacobian
408: - iscoloring - the coloring of the matrix
410: Output Parameter:
411: . color - the new coloring context
412:
413: Level: intermediate
415: .seealso: MatFDColoringDestroy(),SNESDefaultComputeJacobianColor(), ISColoringCreate(),
416: MatFDColoringSetFunction(), MatFDColoringSetFromOptions(), MatFDColoringApply(),
417: MatFDColoringSetFrequency(), MatFDColoringSetRecompute(), MatFDColoringView(),
418: MatFDColoringSetParameters()
419: @*/
420: PetscErrorCode MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
421: {
422: MatFDColoring c;
423: MPI_Comm comm;
425: PetscInt M,N;
426: PetscMPIInt size;
430: MatGetSize(mat,&M,&N);
431: if (M != N) SETERRQ(PETSC_ERR_SUP,"Only for square matrices");
433: PetscObjectGetComm((PetscObject)mat,&comm);
434: PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_COOKIE,0,"MatFDColoring",comm,MatFDColoringDestroy,MatFDColoringView);
436: MPI_Comm_size(comm,&size);
437: c->ctype = iscoloring->ctype;
439: if (mat->ops->fdcoloringcreate) {
440: (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);
441: } else {
442: SETERRQ(PETSC_ERR_SUP,"Code not yet written for this matrix type");
443: }
445: MatGetVecs(mat,PETSC_NULL,&c->w1);
446: PetscLogObjectParent(c,c->w1);
447: VecDuplicate(c->w1,&c->w2);
448: PetscLogObjectParent(c,c->w2);
450: c->error_rel = PETSC_SQRT_MACHINE_EPSILON;
451: c->umin = 100.0*PETSC_SQRT_MACHINE_EPSILON;
452: c->freq = 1;
453: c->usersetsrecompute = PETSC_FALSE;
454: c->recompute = PETSC_FALSE;
455: c->currentcolor = -1;
456: c->htype = "wp";
458: *color = c;
460: return(0);
461: }
465: /*@
466: MatFDColoringDestroy - Destroys a matrix coloring context that was created
467: via MatFDColoringCreate().
469: Collective on MatFDColoring
471: Input Parameter:
472: . c - coloring context
474: Level: intermediate
476: .seealso: MatFDColoringCreate()
477: @*/
478: PetscErrorCode MatFDColoringDestroy(MatFDColoring c)
479: {
481: PetscInt i;
484: if (--c->refct > 0) return(0);
486: for (i=0; i<c->ncolors; i++) {
487: PetscFree(c->columns[i]);
488: PetscFree(c->rows[i]);
489: PetscFree(c->columnsforrow[i]);
490: if (c->vscaleforrow) {PetscFree(c->vscaleforrow[i]);}
491: }
492: PetscFree(c->ncolumns);
493: PetscFree(c->columns);
494: PetscFree(c->nrows);
495: PetscFree(c->rows);
496: PetscFree(c->columnsforrow);
497: PetscFree(c->vscaleforrow);
498: if (c->vscale) {VecDestroy(c->vscale);}
499: if (c->w1) {
500: VecDestroy(c->w1);
501: VecDestroy(c->w2);
502: }
503: if (c->w3){
504: VecDestroy(c->w3);
505: }
506: PetscHeaderDestroy(c);
507: return(0);
508: }
512: /*@C
513: MatFDColoringGetPerturbedColumns - Returns the indices of the columns that
514: that are currently being perturbed.
516: Not Collective
518: Input Parameters:
519: . coloring - coloring context created with MatFDColoringCreate()
521: Output Parameters:
522: + n - the number of local columns being perturbed
523: - cols - the column indices, in global numbering
525: Level: intermediate
527: .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringApply()
529: .keywords: coloring, Jacobian, finite differences
530: @*/
531: PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring coloring,PetscInt *n,PetscInt *cols[])
532: {
534: if (coloring->currentcolor >= 0) {
535: *n = coloring->ncolumns[coloring->currentcolor];
536: *cols = coloring->columns[coloring->currentcolor];
537: } else {
538: *n = 0;
539: }
540: return(0);
541: }
543: #include "petscda.h" /*I "petscda.h" I*/
547: /*@
548: MatFDColoringApply - Given a matrix for which a MatFDColoring context
549: has been created, computes the Jacobian for a function via finite differences.
551: Collective on MatFDColoring
553: Input Parameters:
554: + mat - location to store Jacobian
555: . coloring - coloring context created with MatFDColoringCreate()
556: . x1 - location at which Jacobian is to be computed
557: - sctx - optional context required by function (actually a SNES context)
559: Options Database Keys:
560: + -mat_fd_coloring_freq <freq> - Sets coloring frequency
561: . -mat_fd_type - "wp" or "ds" (see MATSNESMF_WP or MATSNESMF_DS)
562: . -mat_fd_coloring_view - Activates basic viewing or coloring
563: . -mat_fd_coloring_view_draw - Activates drawing of coloring
564: - -mat_fd_coloring_view_info - Activates viewing of coloring info
566: Level: intermediate
568: .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
570: .keywords: coloring, Jacobian, finite differences
571: @*/
572: PetscErrorCode MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
573: {
574: PetscErrorCode (*f)(void*,Vec,Vec,void*) = (PetscErrorCode (*)(void*,Vec,Vec,void *))coloring->f;
576: PetscInt k,start,end,l,row,col,srow,**vscaleforrow,m1,m2;
577: PetscScalar dx,*y,*xx,*w3_array;
578: PetscScalar *vscale_array;
579: PetscReal epsilon = coloring->error_rel,umin = coloring->umin,unorm;
580: Vec w1=coloring->w1,w2=coloring->w2,w3;
581: void *fctx = coloring->fctx;
582: PetscTruth flg;
583: PetscInt ctype=coloring->ctype,N,col_start=0,col_end=0;
584: Vec x1_tmp;
590: if (!f) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetFunction()");
592: if (coloring->usersetsrecompute) {
593: if (!coloring->recompute) {
594: *flag = SAME_PRECONDITIONER;
595: PetscInfo(J,"Skipping Jacobian, since user called MatFDColorSetRecompute()\n");
596: return(0);
597: } else {
598: coloring->recompute = PETSC_FALSE;
599: }
600: }
603: MatSetUnfactored(J);
604: PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);
605: if (flg) {
606: PetscInfo(coloring,"Not calling MatZeroEntries()\n");
607: } else {
608: PetscTruth assembled;
609: MatAssembled(J,&assembled);
610: if (assembled) {
611: MatZeroEntries(J);
612: }
613: }
615: x1_tmp = x1;
616: if (!coloring->vscale){
617: VecDuplicate(x1_tmp,&coloring->vscale);
618: }
619:
620: /*
621: This is a horrible, horrible, hack. See DMMGComputeJacobian_Multigrid() it inproperly sets
622: coloring->F for the coarser grids from the finest
623: */
624: if (coloring->F) {
625: VecGetLocalSize(coloring->F,&m1);
626: VecGetLocalSize(w1,&m2);
627: if (m1 != m2) {
628: coloring->F = 0;
629: }
630: }
632: if (coloring->htype[0] == 'w') { /* tacky test; need to make systematic if we add other approaches to computing h*/
633: VecNorm(x1_tmp,NORM_2,&unorm);
634: }
635: VecGetOwnershipRange(w1,&start,&end); /* OwnershipRange is used by ghosted x! */
636:
637: /* Set w1 = F(x1) */
638: if (coloring->F) {
639: w1 = coloring->F; /* use already computed value of function */
640: coloring->F = 0;
641: } else {
643: (*f)(sctx,x1_tmp,w1,fctx);
645: }
646:
647: if (!coloring->w3) {
648: VecDuplicate(x1_tmp,&coloring->w3);
649: PetscLogObjectParent(coloring,coloring->w3);
650: }
651: w3 = coloring->w3;
653: /* Compute all the local scale factors, including ghost points */
654: VecGetLocalSize(x1_tmp,&N);
655: VecGetArray(x1_tmp,&xx);
656: VecGetArray(coloring->vscale,&vscale_array);
657: if (ctype == IS_COLORING_GHOSTED){
658: col_start = 0; col_end = N;
659: } else if (ctype == IS_COLORING_LOCAL){
660: xx = xx - start;
661: vscale_array = vscale_array - start;
662: col_start = start; col_end = N + start;
663: }
664: for (col=col_start; col<col_end; col++){
665: /* Loop over each local column, vscale[col] = 1./(epsilon*dx[col]) */
666: if (coloring->htype[0] == 'w') {
667: dx = 1.0 + unorm;
668: } else {
669: dx = xx[col];
670: }
671: if (dx == 0.0) dx = 1.0;
672: #if !defined(PETSC_USE_COMPLEX)
673: if (dx < umin && dx >= 0.0) dx = umin;
674: else if (dx < 0.0 && dx > -umin) dx = -umin;
675: #else
676: if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin;
677: else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
678: #endif
679: dx *= epsilon;
680: vscale_array[col] = 1.0/dx;
681: }
682: if (ctype == IS_COLORING_LOCAL) vscale_array = vscale_array + start;
683: VecRestoreArray(coloring->vscale,&vscale_array);
684: if (ctype == IS_COLORING_LOCAL){
685: VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);
686: VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);
687: }
688:
689: if (coloring->vscaleforrow) {
690: vscaleforrow = coloring->vscaleforrow;
691: } else {
692: SETERRQ(PETSC_ERR_ARG_NULL,"Null Object: coloring->vscaleforrow");
693: }
695: /*
696: Loop over each color
697: */
698: VecGetArray(coloring->vscale,&vscale_array);
699: for (k=0; k<coloring->ncolors; k++) {
700: coloring->currentcolor = k;
701: VecCopy(x1_tmp,w3);
702: VecGetArray(w3,&w3_array);
703: if (ctype == IS_COLORING_LOCAL) w3_array = w3_array - start;
704: /*
705: Loop over each column associated with color
706: adding the perturbation to the vector w3.
707: */
708: for (l=0; l<coloring->ncolumns[k]; l++) {
709: col = coloring->columns[k][l]; /* local column of the matrix we are probing for */
710: if (coloring->htype[0] == 'w') {
711: dx = 1.0 + unorm;
712: } else {
713: dx = xx[col];
714: }
715: if (dx == 0.0) dx = 1.0;
716: #if !defined(PETSC_USE_COMPLEX)
717: if (dx < umin && dx >= 0.0) dx = umin;
718: else if (dx < 0.0 && dx > -umin) dx = -umin;
719: #else
720: if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin;
721: else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
722: #endif
723: dx *= epsilon;
724: if (!PetscAbsScalar(dx)) SETERRQ(PETSC_ERR_PLIB,"Computed 0 differencing parameter");
725: w3_array[col] += dx;
726: }
727: if (ctype == IS_COLORING_LOCAL) w3_array = w3_array + start;
728: VecRestoreArray(w3,&w3_array);
730: /*
731: Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations)
732: w2 = F(x1 + dx) - F(x1)
733: */
735: (*f)(sctx,w3,w2,fctx);
737: VecAXPY(w2,-1.0,w1);
738:
739: /*
740: Loop over rows of vector, putting results into Jacobian matrix
741: */
742: VecGetArray(w2,&y);
743: for (l=0; l<coloring->nrows[k]; l++) {
744: row = coloring->rows[k][l]; /* local row index */
745: col = coloring->columnsforrow[k][l]; /* global column index */
746: y[row] *= vscale_array[vscaleforrow[k][l]];
747: srow = row + start;
748: MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);
749: }
750: VecRestoreArray(w2,&y);
751: } /* endof for each color */
752: if (ctype == IS_COLORING_LOCAL) xx = xx + start;
753: VecRestoreArray(coloring->vscale,&vscale_array);
754: VecRestoreArray(x1_tmp,&xx);
755:
756: coloring->currentcolor = -1;
757: MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
758: MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
761: PetscOptionsHasName(PETSC_NULL,"-mat_null_space_test",&flg);
762: if (flg) {
763: MatNullSpaceTest(J->nullsp,J);
764: }
765: MatFDColoringView_Private(coloring);
766: return(0);
767: }
771: /*@
772: MatFDColoringApplyTS - Given a matrix for which a MatFDColoring context
773: has been created, computes the Jacobian for a function via finite differences.
775: Collective on Mat, MatFDColoring, and Vec
777: Input Parameters:
778: + mat - location to store Jacobian
779: . coloring - coloring context created with MatFDColoringCreate()
780: . x1 - location at which Jacobian is to be computed
781: - sctx - optional context required by function (actually a SNES context)
783: Options Database Keys:
784: . -mat_fd_coloring_freq <freq> - Sets coloring frequency
786: Level: intermediate
788: .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
790: .keywords: coloring, Jacobian, finite differences
791: @*/
792: PetscErrorCode MatFDColoringApplyTS(Mat J,MatFDColoring coloring,PetscReal t,Vec x1,MatStructure *flag,void *sctx)
793: {
794: PetscErrorCode (*f)(void*,PetscReal,Vec,Vec,void*)=(PetscErrorCode (*)(void*,PetscReal,Vec,Vec,void *))coloring->f;
796: PetscInt k,N,start,end,l,row,col,srow,**vscaleforrow;
797: PetscScalar dx,*y,*xx,*w3_array;
798: PetscScalar *vscale_array;
799: PetscReal epsilon = coloring->error_rel,umin = coloring->umin;
800: Vec w1=coloring->w1,w2=coloring->w2,w3;
801: void *fctx = coloring->fctx;
802: PetscTruth flg;
810: if (!coloring->w3) {
811: VecDuplicate(x1,&coloring->w3);
812: PetscLogObjectParent(coloring,coloring->w3);
813: }
814: w3 = coloring->w3;
816: MatSetUnfactored(J);
817: PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);
818: if (flg) {
819: PetscInfo(coloring,"Not calling MatZeroEntries()\n");
820: } else {
821: PetscTruth assembled;
822: MatAssembled(J,&assembled);
823: if (assembled) {
824: MatZeroEntries(J);
825: }
826: }
828: VecGetOwnershipRange(x1,&start,&end);
829: VecGetSize(x1,&N);
831: (*f)(sctx,t,x1,w1,fctx);
834: /*
835: Compute all the scale factors and share with other processors
836: */
837: VecGetArray(x1,&xx);xx = xx - start;
838: VecGetArray(coloring->vscale,&vscale_array);vscale_array = vscale_array - start;
839: for (k=0; k<coloring->ncolors; k++) {
840: /*
841: Loop over each column associated with color adding the
842: perturbation to the vector w3.
843: */
844: for (l=0; l<coloring->ncolumns[k]; l++) {
845: col = coloring->columns[k][l]; /* column of the matrix we are probing for */
846: dx = xx[col];
847: if (dx == 0.0) dx = 1.0;
848: #if !defined(PETSC_USE_COMPLEX)
849: if (dx < umin && dx >= 0.0) dx = umin;
850: else if (dx < 0.0 && dx > -umin) dx = -umin;
851: #else
852: if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin;
853: else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
854: #endif
855: dx *= epsilon;
856: vscale_array[col] = 1.0/dx;
857: }
858: }
859: vscale_array = vscale_array - start;VecRestoreArray(coloring->vscale,&vscale_array);
860: VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);
861: VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);
863: if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow;
864: else vscaleforrow = coloring->columnsforrow;
866: VecGetArray(coloring->vscale,&vscale_array);
867: /*
868: Loop over each color
869: */
870: for (k=0; k<coloring->ncolors; k++) {
871: VecCopy(x1,w3);
872: VecGetArray(w3,&w3_array);w3_array = w3_array - start;
873: /*
874: Loop over each column associated with color adding the
875: perturbation to the vector w3.
876: */
877: for (l=0; l<coloring->ncolumns[k]; l++) {
878: col = coloring->columns[k][l]; /* column of the matrix we are probing for */
879: dx = xx[col];
880: if (dx == 0.0) dx = 1.0;
881: #if !defined(PETSC_USE_COMPLEX)
882: if (dx < umin && dx >= 0.0) dx = umin;
883: else if (dx < 0.0 && dx > -umin) dx = -umin;
884: #else
885: if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin;
886: else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
887: #endif
888: dx *= epsilon;
889: w3_array[col] += dx;
890: }
891: w3_array = w3_array + start; VecRestoreArray(w3,&w3_array);
893: /*
894: Evaluate function at x1 + dx (here dx is a vector of perturbations)
895: */
897: (*f)(sctx,t,w3,w2,fctx);
899: VecAXPY(w2,-1.0,w1);
901: /*
902: Loop over rows of vector, putting results into Jacobian matrix
903: */
904: VecGetArray(w2,&y);
905: for (l=0; l<coloring->nrows[k]; l++) {
906: row = coloring->rows[k][l];
907: col = coloring->columnsforrow[k][l];
908: y[row] *= vscale_array[vscaleforrow[k][l]];
909: srow = row + start;
910: MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);
911: }
912: VecRestoreArray(w2,&y);
913: }
914: VecRestoreArray(coloring->vscale,&vscale_array);
915: xx = xx + start; VecRestoreArray(x1,&xx);
916: MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
917: MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
919: return(0);
920: }
925: /*@C
926: MatFDColoringSetRecompute - Indicates that the next time a Jacobian preconditioner
927: is needed it sholuld be recomputed. Once this is called and the new Jacobian is computed
928: no additional Jacobian's will be computed (the same one will be used) until this is
929: called again.
931: Collective on MatFDColoring
933: Input Parameters:
934: . c - the coloring context
936: Level: intermediate
938: Notes: The MatFDColoringSetFrequency() is ignored once this is called
940: .seealso: MatFDColoringCreate(), MatFDColoringSetFrequency()
942: .keywords: Mat, finite differences, coloring
943: @*/
944: PetscErrorCode MatFDColoringSetRecompute(MatFDColoring c)
945: {
948: c->usersetsrecompute = PETSC_TRUE;
949: c->recompute = PETSC_TRUE;
950: return(0);
951: }