xref: /petsc/src/mat/matfd/fdmatrix.c (revision ae09f205e0a50ea5bbede23f787c5c49f4a90e16)
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*ae09f205SBarry Smith static char vcid[] = "$Id: fdmatrix.c,v 1.14 1997/08/24 23:27:24 curfman Exp bsmith $";
3bbf0e169SBarry Smith #endif
4bbf0e169SBarry Smith 
5bbf0e169SBarry Smith /*
6639f9d9dSBarry Smith    This is where the abstract matrix operations are defined that are
7639f9d9dSBarry Smith   used for finite difference computations of Jacobians using coloring.
8bbf0e169SBarry Smith */
9bbf0e169SBarry Smith 
10bbf0e169SBarry Smith #include "petsc.h"
11bbf0e169SBarry Smith #include "src/mat/matimpl.h"        /*I "mat.h" I*/
12bbf0e169SBarry Smith #include "src/vec/vecimpl.h"
13bbf0e169SBarry Smith #include "pinclude/pviewer.h"
14bbf0e169SBarry Smith 
155615d1e5SSatish Balay #undef __FUNC__
16d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringView_Draw"
17005c665bSBarry Smith static int MatFDColoringView_Draw(MatFDColoring fd,Viewer viewer)
18005c665bSBarry Smith {
19005c665bSBarry Smith   int         ierr,i,j,pause;
20005c665bSBarry Smith   PetscTruth  isnull;
21005c665bSBarry Smith   Draw        draw;
22d4bb536fSBarry Smith   double      xr,yr,xl,yl,h,w,x,y,xc,yc,scale = 0.0;
23005c665bSBarry Smith   DrawButton  button;
24005c665bSBarry Smith 
25005c665bSBarry Smith   ierr = ViewerDrawGetDraw(viewer,&draw); CHKERRQ(ierr);
26005c665bSBarry Smith   ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); if (isnull) return 0;
27005c665bSBarry Smith   ierr = DrawSyncClear(draw); CHKERRQ(ierr);
28005c665bSBarry Smith 
29005c665bSBarry Smith   xr  = fd->N; yr = fd->M; h = yr/10.0; w = xr/10.0;
30005c665bSBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
31005c665bSBarry Smith   ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
32005c665bSBarry Smith 
33005c665bSBarry Smith 
34005c665bSBarry Smith   /* loop over colors  */
35005c665bSBarry Smith   for (i=0; i<fd->ncolors; i++ ) {
36005c665bSBarry Smith     for ( j=0; j<fd->nrows[i]; j++ ) {
37005c665bSBarry Smith       y = fd->M - fd->rows[i][j] - fd->rstart;
38005c665bSBarry Smith       x = fd->columnsforrow[i][j];
39005c665bSBarry Smith       DrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);
40005c665bSBarry Smith     }
41005c665bSBarry Smith   }
42005c665bSBarry Smith   ierr = DrawSyncFlush(draw); CHKERRQ(ierr);
43005c665bSBarry Smith   ierr = DrawGetPause(draw,&pause); CHKERRQ(ierr);
44005c665bSBarry Smith   if (pause >= 0) { PetscSleep(pause); return 0;}
45005c665bSBarry Smith   ierr = DrawCheckResizedWindow(draw);
46005c665bSBarry Smith   ierr = DrawSyncGetMouseButton(draw,&button,&xc,&yc,0,0);
47005c665bSBarry Smith   while (button != BUTTON_RIGHT) {
48005c665bSBarry Smith     ierr = DrawSyncClear(draw); CHKERRQ(ierr);
49005c665bSBarry Smith     if (button == BUTTON_LEFT) scale = .5;
50005c665bSBarry Smith     else if (button == BUTTON_CENTER) scale = 2.;
51005c665bSBarry Smith     xl = scale*(xl + w - xc) + xc - w*scale;
52005c665bSBarry Smith     xr = scale*(xr - w - xc) + xc + w*scale;
53005c665bSBarry Smith     yl = scale*(yl + h - yc) + yc - h*scale;
54005c665bSBarry Smith     yr = scale*(yr - h - yc) + yc + h*scale;
55005c665bSBarry Smith     w *= scale; h *= scale;
56005c665bSBarry Smith     ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
57005c665bSBarry Smith     /* loop over colors  */
58005c665bSBarry Smith     for (i=0; i<fd->ncolors; i++ ) {
59005c665bSBarry Smith       for ( j=0; j<fd->nrows[i]; j++ ) {
60005c665bSBarry Smith         y = fd->M - fd->rows[i][j] - fd->rstart;
61005c665bSBarry Smith         x = fd->columnsforrow[i][j];
62005c665bSBarry Smith         DrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);
63005c665bSBarry Smith       }
64005c665bSBarry Smith     }
65005c665bSBarry Smith     ierr = DrawCheckResizedWindow(draw);
66005c665bSBarry Smith     ierr = DrawSyncGetMouseButton(draw,&button,&xc,&yc,0,0);
67005c665bSBarry Smith   }
68005c665bSBarry Smith 
69005c665bSBarry Smith   return 0;
70005c665bSBarry Smith }
71005c665bSBarry Smith 
72005c665bSBarry Smith 
73005c665bSBarry Smith #undef __FUNC__
74d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringView"
75bbf0e169SBarry Smith /*@C
76639f9d9dSBarry Smith    MatFDColoringView - Views a finite difference coloring context.
77bbf0e169SBarry Smith 
78639f9d9dSBarry Smith    Input  Parameter:
79639f9d9dSBarry Smith .   color - the coloring context
80bbf0e169SBarry Smith 
81639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
82005c665bSBarry Smith 
83bbf0e169SBarry Smith @*/
84639f9d9dSBarry Smith int MatFDColoringView(MatFDColoring color,Viewer viewer)
85bbf0e169SBarry Smith {
86005c665bSBarry Smith   ViewerType  vtype;
87639f9d9dSBarry Smith   int         i,j,format,ierr;
88bbf0e169SBarry Smith 
89639f9d9dSBarry Smith   if (!viewer) viewer = VIEWER_STDOUT_WORLD;
90bbf0e169SBarry Smith 
91005c665bSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
92005c665bSBarry Smith   if (vtype == DRAW_VIEWER) {
93005c665bSBarry Smith     ierr = MatFDColoringView_Draw(color,viewer); CHKERRQ(ierr);
94005c665bSBarry Smith     return 0;
95005c665bSBarry Smith   }
96005c665bSBarry Smith 
97639f9d9dSBarry Smith   printf("MatFDColoring Object:\n");
98639f9d9dSBarry Smith   printf("  Error tolerance %g\n",color->error_rel);
99*ae09f205SBarry Smith   printf("  Umin %g\n",color->umin);
100639f9d9dSBarry Smith   printf("  Number of colors %d\n",color->ncolors);
101*ae09f205SBarry Smith 
102*ae09f205SBarry Smith   ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr);
103*ae09f205SBarry Smith   if (format != VIEWER_FORMAT_ASCII_INFO) {
104639f9d9dSBarry Smith     for ( i=0; i<color->ncolors; i++ ) {
105639f9d9dSBarry Smith       printf("Information for color %d\n",i);
106639f9d9dSBarry Smith       printf("Number of columns %d\n",color->ncolumns[i]);
107639f9d9dSBarry Smith       for ( j=0; j<color->ncolumns[i]; j++ ) {
108639f9d9dSBarry Smith         printf("  %d\n",color->columns[i][j]);
109639f9d9dSBarry Smith       }
110639f9d9dSBarry Smith       printf("Number of rows %d\n",color->nrows[i]);
111639f9d9dSBarry Smith       for ( j=0; j<color->nrows[i]; j++ ) {
112639f9d9dSBarry Smith         printf("  %d %d \n",color->rows[i][j],color->columnsforrow[i][j]);
113bbf0e169SBarry Smith       }
114bbf0e169SBarry Smith     }
115bbf0e169SBarry Smith   }
116639f9d9dSBarry Smith   return 0;
117639f9d9dSBarry Smith }
118639f9d9dSBarry Smith 
1195615d1e5SSatish Balay #undef __FUNC__
1205615d1e5SSatish Balay #define __FUNC__ "MatFDColoringSetParameters"
121639f9d9dSBarry Smith /*@
122639f9d9dSBarry Smith    MatFDColoringSetParameters - Sets the parameters for the approximation of
123639f9d9dSBarry Smith    Jacobian using finite differences.
124639f9d9dSBarry Smith 
125639f9d9dSBarry Smith $       J(u)_{:,i} = [J(u+h*dx_{i}) - J(u)]/h where
126639f9d9dSBarry Smith $        h = error_rel*u[i]                    if  u[i] > umin
127639f9d9dSBarry Smith $          = error_rel*umin                    else
128639f9d9dSBarry Smith $
129639f9d9dSBarry Smith $   dx_{i} = (0, ... 1, .... 0)
130639f9d9dSBarry Smith 
131639f9d9dSBarry Smith    Input Parameters:
132639f9d9dSBarry Smith .  coloring - the jacobian coloring context
133639f9d9dSBarry Smith .  error_rel - relative error
134639f9d9dSBarry Smith .  umin - minimum allowable u-value
135639f9d9dSBarry Smith 
136639f9d9dSBarry Smith .keywords: SNES, Jacobian, finite differences, parameters
137639f9d9dSBarry Smith @*/
138639f9d9dSBarry Smith int MatFDColoringSetParameters(MatFDColoring matfd,double error,double umin)
139639f9d9dSBarry Smith {
140639f9d9dSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
141639f9d9dSBarry Smith 
142639f9d9dSBarry Smith   if (error != PETSC_DEFAULT) matfd->error_rel = error;
143639f9d9dSBarry Smith   if (umin != PETSC_DEFAULT)  matfd->umin      = umin;
144639f9d9dSBarry Smith   return 0;
145639f9d9dSBarry Smith }
146639f9d9dSBarry Smith 
1475615d1e5SSatish Balay #undef __FUNC__
148005c665bSBarry Smith #define __FUNC__ "MatFDColoringSetFrequency"
149005c665bSBarry Smith /*@
150005c665bSBarry Smith    MatFDColoringSetFrequency - Sets the frequency at which new Jacobians
151005c665bSBarry Smith      are computed.
152005c665bSBarry Smith 
153005c665bSBarry Smith    Input Parameters:
154005c665bSBarry Smith .  coloring - the jacobian coloring context
155005c665bSBarry Smith .  freq - frequency
156005c665bSBarry Smith 
157005c665bSBarry Smith .keywords: SNES, Jacobian, finite differences, parameters
158005c665bSBarry Smith @*/
159005c665bSBarry Smith int MatFDColoringSetFrequency(MatFDColoring matfd,int freq)
160005c665bSBarry Smith {
161005c665bSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
162005c665bSBarry Smith 
163005c665bSBarry Smith   matfd->freq = freq;
164005c665bSBarry Smith   return 0;
165005c665bSBarry Smith }
166005c665bSBarry Smith 
167005c665bSBarry Smith #undef __FUNC__
168005c665bSBarry Smith #define __FUNC__ "MatFDColoringSetFunction"
169005c665bSBarry Smith /*@
170005c665bSBarry Smith    MatFDColoringSetFunction - Sets the function to use for computing the Jacobian.
171005c665bSBarry Smith 
172005c665bSBarry Smith    Input Parameters:
173005c665bSBarry Smith .  coloring - the jacobian coloring context
174005c665bSBarry Smith .  f - the function
175005c665bSBarry Smith .  fctx - the function context
176005c665bSBarry Smith 
177005c665bSBarry Smith 
178005c665bSBarry Smith .keywords: SNES, Jacobian, finite differences, parameters, function
179005c665bSBarry Smith @*/
180005c665bSBarry Smith int MatFDColoringSetFunction(MatFDColoring matfd,int (*f)(void *,Vec,Vec,void *),void *fctx)
181005c665bSBarry Smith {
182005c665bSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
183005c665bSBarry Smith 
184005c665bSBarry Smith   matfd->f    = f;
185005c665bSBarry Smith   matfd->fctx = fctx;
186005c665bSBarry Smith 
187005c665bSBarry Smith   return 0;
188005c665bSBarry Smith }
189005c665bSBarry Smith 
190005c665bSBarry Smith #undef __FUNC__
191d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringSetFromOptions"
192639f9d9dSBarry Smith /*@
193639f9d9dSBarry Smith    MatFDColoringSetFromOptions - Set coloring finite difference parameters from
194639f9d9dSBarry Smith          the options database.
195639f9d9dSBarry Smith 
196639f9d9dSBarry Smith $       J(u)_{:,i} = [J(u+h*dx_{i}) - J(u)]/h where
197639f9d9dSBarry Smith $        h = error_rel*u[i]                    if  u[i] > umin
198*ae09f205SBarry Smith $          = error_rel*umin                      else
199639f9d9dSBarry Smith $
200639f9d9dSBarry Smith $   dx_{i} = (0, ... 1, .... 0)
201639f9d9dSBarry Smith 
202639f9d9dSBarry Smith    Input Parameters:
203639f9d9dSBarry Smith .  coloring - the jacobian coloring context
204639f9d9dSBarry Smith 
205639f9d9dSBarry Smith    Options Database:
206639f9d9dSBarry Smith .  -mat_fd_coloring_error square root of relative error in function
207639f9d9dSBarry Smith .  -mat_fd_coloring_umin  see above
208005c665bSBarry Smith .  -mat_fd_coloring_freq frequency at which Jacobian is computed
209639f9d9dSBarry Smith 
210639f9d9dSBarry Smith .keywords: SNES, Jacobian, finite differences, parameters
211639f9d9dSBarry Smith @*/
212639f9d9dSBarry Smith int MatFDColoringSetFromOptions(MatFDColoring matfd)
213639f9d9dSBarry Smith {
214005c665bSBarry Smith   int    ierr,flag,freq = 1;
215639f9d9dSBarry Smith   double error = PETSC_DEFAULT,umin = PETSC_DEFAULT;
216639f9d9dSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
217639f9d9dSBarry Smith 
218639f9d9dSBarry Smith   ierr = OptionsGetDouble(matfd->prefix,"-mat_fd_coloring_err",&error,&flag);CHKERRQ(ierr);
219639f9d9dSBarry Smith   ierr = OptionsGetDouble(matfd->prefix,"-mat_fd_coloring_umin",&umin,&flag);CHKERRQ(ierr);
220639f9d9dSBarry Smith   ierr = MatFDColoringSetParameters(matfd,error,umin); CHKERRQ(ierr);
221005c665bSBarry Smith   ierr = OptionsGetInt(matfd->prefix,"-mat_fd_coloring_freq",&freq,&flag);CHKERRQ(ierr);
222005c665bSBarry Smith   ierr = MatFDColoringSetFrequency(matfd,freq);CHKERRQ(ierr);
223005c665bSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help",&flag); CHKERRQ(ierr);
224639f9d9dSBarry Smith   if (flag) {
225639f9d9dSBarry Smith     ierr = MatFDColoringPrintHelp(matfd); CHKERRQ(ierr);
226639f9d9dSBarry Smith   }
227639f9d9dSBarry Smith   return 0;
228639f9d9dSBarry Smith }
229639f9d9dSBarry Smith 
2305615d1e5SSatish Balay #undef __FUNC__
231d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringPrintHelp"
232639f9d9dSBarry Smith /*@
233639f9d9dSBarry Smith     MatFDColoringPrintHelp - Prints help message for matrix finite difference calculations
234639f9d9dSBarry Smith          using coloring.
235639f9d9dSBarry Smith 
236639f9d9dSBarry Smith    Input Parameter:
237639f9d9dSBarry Smith .   fdcoloring - the MatFDColoring context
238639f9d9dSBarry Smith 
239639f9d9dSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringSetFromOptions()
240639f9d9dSBarry Smith @*/
241639f9d9dSBarry Smith int MatFDColoringPrintHelp(MatFDColoring fd)
242639f9d9dSBarry Smith {
243639f9d9dSBarry Smith   PetscValidHeaderSpecific(fd,MAT_FDCOLORING_COOKIE);
244639f9d9dSBarry Smith 
2450f665d81SLois Curfman McInnes   PetscPrintf(fd->comm,"-mat_fd_coloring_err <err>: set sqrt rel tol in function, defaults to 1.e-8\n");
2460f665d81SLois Curfman McInnes   PetscPrintf(fd->comm,"-mat_fd_coloring_umin <umin>: see users manual, defaults to 1.e-8\n");
247005c665bSBarry Smith   PetscPrintf(fd->comm,"-mat_fd_coloring_freq <freq>: frequency that Jacobian is recomputed, defaults 1\n");
248005c665bSBarry Smith   PetscPrintf(fd->comm,"-mat_fd_coloring_view\n");
249005c665bSBarry Smith   PetscPrintf(fd->comm,"-mat_fd_coloring_view_draw\n");
250*ae09f205SBarry Smith   PetscPrintf(fd->comm,"-mat_fd_coloring_view_info\n");
251005c665bSBarry Smith   return 0;
252005c665bSBarry Smith }
253005c665bSBarry Smith 
254005c665bSBarry Smith int MatFDColoringView_Private(MatFDColoring fd)
255005c665bSBarry Smith {
256005c665bSBarry Smith   int ierr,flg;
257005c665bSBarry Smith 
258005c665bSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view",&flg); CHKERRQ(ierr);
259005c665bSBarry Smith   if (flg) {
260005c665bSBarry Smith     Viewer viewer;
261005c665bSBarry Smith     ierr = ViewerFileOpenASCII(fd->comm,"stdout",&viewer);CHKERRQ(ierr);
262005c665bSBarry Smith     ierr = MatFDColoringView(fd,viewer); CHKERRQ(ierr);
263005c665bSBarry Smith     ierr = ViewerDestroy(viewer); CHKERRQ(ierr);
264005c665bSBarry Smith   }
265*ae09f205SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_info",&flg); CHKERRQ(ierr);
266*ae09f205SBarry Smith   if (flg) {
267*ae09f205SBarry Smith     Viewer viewer;
268*ae09f205SBarry Smith     ierr = ViewerFileOpenASCII(fd->comm,"stdout",&viewer);CHKERRQ(ierr);
269*ae09f205SBarry Smith     ierr = ViewerPushFormat(viewer,VIEWER_FORMAT_ASCII_INFO,PETSC_NULL);CHKERRQ(ierr);
270*ae09f205SBarry Smith     ierr = MatFDColoringView(fd,viewer); CHKERRQ(ierr);
271*ae09f205SBarry Smith     ierr = ViewerPopFormat(viewer);CHKERRQ(ierr);
272*ae09f205SBarry Smith     ierr = ViewerDestroy(viewer); CHKERRQ(ierr);
273*ae09f205SBarry Smith   }
274005c665bSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg); CHKERRQ(ierr);
275005c665bSBarry Smith   if (flg) {
276005c665bSBarry Smith     ierr = MatFDColoringView(fd,VIEWER_DRAWX_(fd->comm)); CHKERRQ(ierr);
277005c665bSBarry Smith     ierr = ViewerFlush(VIEWER_DRAWX_(fd->comm)); CHKERRQ(ierr);
278005c665bSBarry Smith   }
279bbf0e169SBarry Smith   return 0;
280bbf0e169SBarry Smith }
281bbf0e169SBarry Smith 
2825615d1e5SSatish Balay #undef __FUNC__
2835615d1e5SSatish Balay #define __FUNC__ "MatFDColoringCreate"
284bbf0e169SBarry Smith /*@C
285639f9d9dSBarry Smith    MatFDColoringCreate - Creates a matrix coloring context for finite difference
286639f9d9dSBarry Smith         computation of Jacobians.
287bbf0e169SBarry Smith 
288639f9d9dSBarry Smith    Input Parameters:
289639f9d9dSBarry Smith .    mat - the matrix containing the nonzero structure of the Jacobian
290639f9d9dSBarry Smith .    iscoloring - the coloring of the matrix
291bbf0e169SBarry Smith 
292bbf0e169SBarry Smith    Output Parameter:
293639f9d9dSBarry Smith .   color - the new coloring context
294bbf0e169SBarry Smith 
295005c665bSBarry Smith    Options Database:
296005c665bSBarry Smith .  -mat_fd_coloring_error square root of relative error in function
297005c665bSBarry Smith .  -mat_fd_coloring_umin  see above
298005c665bSBarry Smith .  -mat_fd_coloring_view
299005c665bSBarry Smith .  -mat_fd_coloring_view_draw
300639f9d9dSBarry Smith 
301639f9d9dSBarry Smith .seealso: MatFDColoringDestroy()
302bbf0e169SBarry Smith @*/
303639f9d9dSBarry Smith int MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
304bbf0e169SBarry Smith {
305639f9d9dSBarry Smith   MatFDColoring c;
306639f9d9dSBarry Smith   MPI_Comm      comm;
307639f9d9dSBarry Smith   int           ierr,M,N;
308639f9d9dSBarry Smith 
309639f9d9dSBarry Smith   ierr = MatGetSize(mat,&M,&N); CHKERRQ(ierr);
310e3372554SBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,0,"Only for square matrices");
311639f9d9dSBarry Smith 
312639f9d9dSBarry Smith   PetscObjectGetComm((PetscObject)mat,&comm);
313d4bb536fSBarry Smith   PetscHeaderCreate(c,_p_MatFDColoring,MAT_FDCOLORING_COOKIE,0,comm,MatFDColoringDestroy,MatFDColoringView);
314639f9d9dSBarry Smith   PLogObjectCreate(c);
315639f9d9dSBarry Smith 
316639f9d9dSBarry Smith   if (mat->ops.fdcoloringcreate) {
317639f9d9dSBarry Smith     ierr = (*mat->ops.fdcoloringcreate)(mat,iscoloring,c); CHKERRQ(ierr);
318639f9d9dSBarry Smith   } else {
319e3372554SBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"Code not yet written for this matrix type");
320639f9d9dSBarry Smith   }
321639f9d9dSBarry Smith 
322639f9d9dSBarry Smith   c->error_rel = 1.e-8;
323*ae09f205SBarry Smith   c->umin      = 1.e-6;
324005c665bSBarry Smith   c->freq      = 1;
325005c665bSBarry Smith 
326005c665bSBarry Smith   ierr = MatFDColoringView_Private(c); CHKERRQ(ierr);
327639f9d9dSBarry Smith 
328639f9d9dSBarry Smith   *color = c;
329639f9d9dSBarry Smith 
330bbf0e169SBarry Smith   return 0;
331bbf0e169SBarry Smith }
332bbf0e169SBarry Smith 
3335615d1e5SSatish Balay #undef __FUNC__
334d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringDestroy"
335bbf0e169SBarry Smith /*@C
336639f9d9dSBarry Smith     MatFDColoringDestroy - Destroys a matrix coloring context that was created
337639f9d9dSBarry Smith     via MatFDColoringCreate().
338bbf0e169SBarry Smith 
339639f9d9dSBarry Smith     Input Paramter:
340639f9d9dSBarry Smith .   c - coloring context
341bbf0e169SBarry Smith 
342639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
343bbf0e169SBarry Smith @*/
344639f9d9dSBarry Smith int MatFDColoringDestroy(MatFDColoring c)
345bbf0e169SBarry Smith {
346639f9d9dSBarry Smith   int i,ierr,flag;
347bbf0e169SBarry Smith 
348d4bb536fSBarry Smith   if (--c->refct > 0) return 0;
349d4bb536fSBarry Smith 
350639f9d9dSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-matfdcoloring_view",&flag);
351639f9d9dSBarry Smith   if (flag) {
352bbf0e169SBarry Smith     Viewer viewer;
353639f9d9dSBarry Smith     ierr = ViewerFileOpenASCII(c->comm,"stdout",&viewer);CHKERRQ(ierr);
354639f9d9dSBarry Smith     ierr = MatFDColoringView(c,viewer);CHKERRQ(ierr);
355bbf0e169SBarry Smith     ierr = ViewerDestroy(viewer); CHKERRQ(ierr);
356bbf0e169SBarry Smith   }
357639f9d9dSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-matfdcoloring_view_info",&flag);
358639f9d9dSBarry Smith   if (flag) {
359bbf0e169SBarry Smith     Viewer viewer;
360639f9d9dSBarry Smith     ierr = ViewerFileOpenASCII(c->comm,"stdout",&viewer);CHKERRQ(ierr);
361639f9d9dSBarry Smith     ierr = ViewerPushFormat(viewer,VIEWER_FORMAT_ASCII_INFO, PETSC_NULL); CHKERRQ(ierr);
362639f9d9dSBarry Smith     ierr = MatFDColoringView(c,viewer);CHKERRQ(ierr);
363639f9d9dSBarry Smith     ierr = ViewerPopFormat(viewer);
364bbf0e169SBarry Smith     ierr = ViewerDestroy(viewer); CHKERRQ(ierr);
365bbf0e169SBarry Smith   }
366639f9d9dSBarry Smith 
367639f9d9dSBarry Smith   for ( i=0; i<c->ncolors; i++ ) {
368639f9d9dSBarry Smith     if (c->columns[i])       PetscFree(c->columns[i]);
369639f9d9dSBarry Smith     if (c->rows[i])          PetscFree(c->rows[i]);
370639f9d9dSBarry Smith     if (c->columnsforrow[i]) PetscFree(c->columnsforrow[i]);
371bbf0e169SBarry Smith   }
372639f9d9dSBarry Smith   PetscFree(c->ncolumns);
373639f9d9dSBarry Smith   PetscFree(c->columns);
374639f9d9dSBarry Smith   PetscFree(c->nrows);
375639f9d9dSBarry Smith   PetscFree(c->rows);
376639f9d9dSBarry Smith   PetscFree(c->columnsforrow);
377639f9d9dSBarry Smith   PetscFree(c->scale);
378005c665bSBarry Smith   if (c->w1) {
379005c665bSBarry Smith     ierr = VecDestroy(c->w1); CHKERRQ(ierr);
380005c665bSBarry Smith     ierr = VecDestroy(c->w2); CHKERRQ(ierr);
381005c665bSBarry Smith     ierr = VecDestroy(c->w3); CHKERRQ(ierr);
382005c665bSBarry Smith   }
383639f9d9dSBarry Smith   PLogObjectDestroy(c);
384639f9d9dSBarry Smith   PetscHeaderDestroy(c);
385bbf0e169SBarry Smith   return 0;
386bbf0e169SBarry Smith }
38743a90d84SBarry Smith 
388005c665bSBarry Smith #include "snes.h"
389005c665bSBarry Smith 
3905615d1e5SSatish Balay #undef __FUNC__
3915615d1e5SSatish Balay #define __FUNC__ "MatFDColoringApply"
39243a90d84SBarry Smith /*@
39343a90d84SBarry Smith     MatFDColoringApply - Given a matrix for which a MatFDColoring has been created,
39443a90d84SBarry Smith     computes the Jacobian for a function via finite differences.
39543a90d84SBarry Smith 
39643a90d84SBarry Smith     Input Parameters:
39743a90d84SBarry Smith .   mat - location to store Jacobian
39843a90d84SBarry Smith .   coloring - coloring context created with MatFDColoringCreate()
39943a90d84SBarry Smith .   x1 - location at which Jacobian is to be computed
400005c665bSBarry Smith .   sctx - optional context required by function (actually a SNES context)
40143a90d84SBarry Smith 
40243a90d84SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
40343a90d84SBarry Smith 
40443a90d84SBarry Smith .keywords: coloring, Jacobian, finite differences
40543a90d84SBarry Smith @*/
406005c665bSBarry Smith int MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
40743a90d84SBarry Smith {
408dd9fa9a5SLois Curfman McInnes   int           k, ierr,N,start,end,l,row,col,srow;
40943a90d84SBarry Smith   Scalar        dx, mone = -1.0,*y,*scale = coloring->scale,*xx,*wscale = coloring->wscale;
41043a90d84SBarry Smith   double        epsilon = coloring->error_rel, umin = coloring->umin;
41143a90d84SBarry Smith   MPI_Comm      comm = coloring->comm;
412005c665bSBarry Smith   Vec           w1,w2,w3;
413dd9fa9a5SLois Curfman McInnes   int           (*f)(void *,Vec,Vec,void *) = coloring->f;
414005c665bSBarry Smith   void          *fctx = coloring->fctx;
415005c665bSBarry Smith 
416d4bb536fSBarry Smith   /*
417005c665bSBarry Smith   ierr = SNESGetIterationNumber((SNES)sctx,&it); CHKERRQ(ierr);
418005c665bSBarry Smith   if ((freq > 1) && ((it % freq) != 1)) {
419005c665bSBarry Smith     PLogInfo(coloring,"MatFDColoringApply:Skipping Jacobian, iteration %d, freq %d\n",it,freq);
420005c665bSBarry Smith     *flag = SAME_PRECONDITIONER;
421005c665bSBarry Smith     return 0;
422005c665bSBarry Smith   } else {
423005c665bSBarry Smith     PLogInfo(coloring,"MatFDColoringApply:Computing Jacobian, iteration %d, freq %d\n",it,freq);
424005c665bSBarry Smith     *flag = SAME_NONZERO_PATTERN;
425d4bb536fSBarry Smith   }*/
426005c665bSBarry Smith 
427005c665bSBarry Smith   if (!coloring->w1) {
428005c665bSBarry Smith     ierr = VecDuplicate(x1,&coloring->w1); CHKERRQ(ierr);
429005c665bSBarry Smith     PLogObjectParent(coloring,coloring->w1);
430005c665bSBarry Smith     ierr = VecDuplicate(x1,&coloring->w2); CHKERRQ(ierr);
431005c665bSBarry Smith     PLogObjectParent(coloring,coloring->w2);
432005c665bSBarry Smith     ierr = VecDuplicate(x1,&coloring->w3); CHKERRQ(ierr);
433005c665bSBarry Smith     PLogObjectParent(coloring,coloring->w3);
434005c665bSBarry Smith   }
435005c665bSBarry Smith   w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
43643a90d84SBarry Smith 
43743a90d84SBarry Smith   ierr = MatZeroEntries(J); CHKERRQ(ierr);
43843a90d84SBarry Smith 
43943a90d84SBarry Smith   ierr = VecGetOwnershipRange(x1,&start,&end); CHKERRQ(ierr);
44043a90d84SBarry Smith   ierr = VecGetSize(x1,&N); CHKERRQ(ierr);
44143a90d84SBarry Smith   ierr = VecGetArray(x1,&xx); CHKERRQ(ierr);
44243a90d84SBarry Smith   ierr = (*f)(sctx,x1,w1,fctx); CHKERRQ(ierr);
44343a90d84SBarry Smith 
44443a90d84SBarry Smith   PetscMemzero(wscale,N*sizeof(Scalar));
44543a90d84SBarry Smith   /*
44643a90d84SBarry Smith       Loop over each color
44743a90d84SBarry Smith   */
44843a90d84SBarry Smith 
44943a90d84SBarry Smith   for (k=0; k<coloring->ncolors; k++) {
45043a90d84SBarry Smith     ierr = VecCopy(x1,w3); CHKERRQ(ierr);
45143a90d84SBarry Smith     /*
45243a90d84SBarry Smith        Loop over each column associated with color adding the
45343a90d84SBarry Smith        perturbation to the vector w3.
45443a90d84SBarry Smith     */
45543a90d84SBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
45643a90d84SBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
45743a90d84SBarry Smith       dx  = xx[col-start];
458*ae09f205SBarry Smith       if (dx == 0.0) dx = 1.0;
45943a90d84SBarry Smith #if !defined(PETSC_COMPLEX)
460*ae09f205SBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
461*ae09f205SBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
46243a90d84SBarry Smith #else
463*ae09f205SBarry Smith       if (abs(dx) < umin && real(dx) >= 0.0)     dx = umin;
464*ae09f205SBarry Smith       else if (real(dx) < 0.0 && abs(dx) < umin) dx = -umin;
46543a90d84SBarry Smith #endif
46643a90d84SBarry Smith       dx          *= epsilon;
46743a90d84SBarry Smith       wscale[col] = 1.0/dx;
46843a90d84SBarry Smith       VecSetValues(w3,1,&col,&dx,ADD_VALUES);
46943a90d84SBarry Smith     }
47043a90d84SBarry Smith     VecRestoreArray(x1,&xx);
47143a90d84SBarry Smith     /*
47243a90d84SBarry Smith        Evaluate function at x1 + dx (here dx is a vector, of perturbations)
47343a90d84SBarry Smith     */
47443a90d84SBarry Smith     ierr = (*f)(sctx,w3,w2,fctx); CHKERRQ(ierr);
47543a90d84SBarry Smith     ierr = VecAXPY(&mone,w1,w2); CHKERRQ(ierr);
47643a90d84SBarry Smith     /* Communicate scale to all processors */
47743a90d84SBarry Smith #if !defined(PETSC_COMPLEX)
47843a90d84SBarry Smith     MPI_Allreduce(wscale,scale,N,MPI_DOUBLE,MPI_SUM,comm);
47943a90d84SBarry Smith #else
48043a90d84SBarry Smith     MPI_Allreduce(wscale,scale,2*N,MPI_DOUBLE,MPI_SUM,comm);
48143a90d84SBarry Smith #endif
48243a90d84SBarry Smith     /*
48343a90d84SBarry Smith        Loop over rows of vector putting results into Jacobian matrix
48443a90d84SBarry Smith     */
48543a90d84SBarry Smith     VecGetArray(w2,&y);
48643a90d84SBarry Smith     for (l=0; l<coloring->nrows[k]; l++) {
48743a90d84SBarry Smith       row    = coloring->rows[k][l];
48843a90d84SBarry Smith       col    = coloring->columnsforrow[k][l];
48943a90d84SBarry Smith       y[row] *= scale[col];
49043a90d84SBarry Smith       srow   = row + start;
49143a90d84SBarry Smith       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
49243a90d84SBarry Smith     }
49343a90d84SBarry Smith     VecRestoreArray(w2,&y);
49443a90d84SBarry Smith   }
49543a90d84SBarry Smith   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
49643a90d84SBarry Smith   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
49743a90d84SBarry Smith   return 0;
49843a90d84SBarry Smith }
499