xref: /petsc/src/mat/matfd/fdmatrix.c (revision 005c665b4cb72e0f49f0521e223449c941bc47c1)
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*005c665bSBarry Smith static char vcid[] = "$Id: fdmatrix.c,v 1.11 1997/07/09 20:53:23 balay 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__
16*005c665bSBarry Smith #define __FUNC__ "MatFDColoringView_Draw" /* ADIC Ignore */
17*005c665bSBarry Smith static int MatFDColoringView_Draw(MatFDColoring fd,Viewer viewer)
18*005c665bSBarry Smith {
19*005c665bSBarry Smith   int         ierr,i,j,pause;
20*005c665bSBarry Smith   PetscTruth  isnull;
21*005c665bSBarry Smith   Draw        draw;
22*005c665bSBarry Smith   double      xr,yr,xl,yl,h,w,x,y,xc,yc,scale;
23*005c665bSBarry Smith   DrawButton  button;
24*005c665bSBarry Smith 
25*005c665bSBarry Smith   ierr = ViewerDrawGetDraw(viewer,&draw); CHKERRQ(ierr);
26*005c665bSBarry Smith   ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); if (isnull) return 0;
27*005c665bSBarry Smith   ierr = DrawSyncClear(draw); CHKERRQ(ierr);
28*005c665bSBarry Smith 
29*005c665bSBarry Smith   xr  = fd->N; yr = fd->M; h = yr/10.0; w = xr/10.0;
30*005c665bSBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
31*005c665bSBarry Smith   ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
32*005c665bSBarry Smith 
33*005c665bSBarry Smith 
34*005c665bSBarry Smith   /* loop over colors  */
35*005c665bSBarry Smith   for (i=0; i<fd->ncolors; i++ ) {
36*005c665bSBarry Smith     for ( j=0; j<fd->nrows[i]; j++ ) {
37*005c665bSBarry Smith       y = fd->M - fd->rows[i][j] - fd->rstart;
38*005c665bSBarry Smith       x = fd->columnsforrow[i][j];
39*005c665bSBarry Smith       DrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);
40*005c665bSBarry Smith     }
41*005c665bSBarry Smith   }
42*005c665bSBarry Smith   ierr = DrawSyncFlush(draw); CHKERRQ(ierr);
43*005c665bSBarry Smith   ierr = DrawGetPause(draw,&pause); CHKERRQ(ierr);
44*005c665bSBarry Smith   if (pause >= 0) { PetscSleep(pause); return 0;}
45*005c665bSBarry Smith   ierr = DrawCheckResizedWindow(draw);
46*005c665bSBarry Smith   ierr = DrawSyncGetMouseButton(draw,&button,&xc,&yc,0,0);
47*005c665bSBarry Smith   while (button != BUTTON_RIGHT) {
48*005c665bSBarry Smith     ierr = DrawSyncClear(draw); CHKERRQ(ierr);
49*005c665bSBarry Smith     if (button == BUTTON_LEFT) scale = .5;
50*005c665bSBarry Smith     else if (button == BUTTON_CENTER) scale = 2.;
51*005c665bSBarry Smith     xl = scale*(xl + w - xc) + xc - w*scale;
52*005c665bSBarry Smith     xr = scale*(xr - w - xc) + xc + w*scale;
53*005c665bSBarry Smith     yl = scale*(yl + h - yc) + yc - h*scale;
54*005c665bSBarry Smith     yr = scale*(yr - h - yc) + yc + h*scale;
55*005c665bSBarry Smith     w *= scale; h *= scale;
56*005c665bSBarry Smith     ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
57*005c665bSBarry Smith     /* loop over colors  */
58*005c665bSBarry Smith     for (i=0; i<fd->ncolors; i++ ) {
59*005c665bSBarry Smith       for ( j=0; j<fd->nrows[i]; j++ ) {
60*005c665bSBarry Smith         y = fd->M - fd->rows[i][j] - fd->rstart;
61*005c665bSBarry Smith         x = fd->columnsforrow[i][j];
62*005c665bSBarry Smith         DrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);
63*005c665bSBarry Smith       }
64*005c665bSBarry Smith     }
65*005c665bSBarry Smith     ierr = DrawCheckResizedWindow(draw);
66*005c665bSBarry Smith     ierr = DrawSyncGetMouseButton(draw,&button,&xc,&yc,0,0);
67*005c665bSBarry Smith   }
68*005c665bSBarry Smith 
69*005c665bSBarry Smith   return 0;
70*005c665bSBarry Smith }
71*005c665bSBarry Smith 
72*005c665bSBarry Smith 
73*005c665bSBarry Smith #undef __FUNC__
745eea60f9SBarry Smith #define __FUNC__ "MatFDColoringView" /* ADIC Ignore */
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()
82*005c665bSBarry Smith 
83bbf0e169SBarry Smith @*/
84639f9d9dSBarry Smith int MatFDColoringView(MatFDColoring color,Viewer viewer)
85bbf0e169SBarry Smith {
86*005c665bSBarry Smith   ViewerType  vtype;
87639f9d9dSBarry Smith   int         i,j,format,ierr;
88bbf0e169SBarry Smith 
89639f9d9dSBarry Smith   if (!viewer) viewer = VIEWER_STDOUT_WORLD;
90bbf0e169SBarry Smith 
91*005c665bSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
92*005c665bSBarry Smith   if (vtype == DRAW_VIEWER) {
93*005c665bSBarry Smith     ierr = MatFDColoringView_Draw(color,viewer); CHKERRQ(ierr);
94*005c665bSBarry Smith     return 0;
95*005c665bSBarry Smith   }
96*005c665bSBarry Smith 
97bbf0e169SBarry Smith   ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr);
98639f9d9dSBarry Smith 
99639f9d9dSBarry Smith   if (format == VIEWER_FORMAT_ASCII_INFO) {
100639f9d9dSBarry Smith     printf("MatFDColoring Object:\n");
101639f9d9dSBarry Smith     printf("  Error tolerance %g\n",color->error_rel);
102639f9d9dSBarry Smith     printf("  umin %g\n",color->umin);
103639f9d9dSBarry Smith   } else {
104639f9d9dSBarry Smith     printf("MatFDColoring Object:\n");
105639f9d9dSBarry Smith     printf("  Error tolerance %g\n",color->error_rel);
106639f9d9dSBarry Smith     printf("  umin %g\n",color->umin);
107639f9d9dSBarry Smith     printf("Number of colors %d\n",color->ncolors);
108639f9d9dSBarry Smith     for ( i=0; i<color->ncolors; i++ ) {
109639f9d9dSBarry Smith       printf("Information for color %d\n",i);
110639f9d9dSBarry Smith       printf("Number of columns %d\n",color->ncolumns[i]);
111639f9d9dSBarry Smith       for ( j=0; j<color->ncolumns[i]; j++ ) {
112639f9d9dSBarry Smith         printf("  %d\n",color->columns[i][j]);
113639f9d9dSBarry Smith       }
114639f9d9dSBarry Smith       printf("Number of rows %d\n",color->nrows[i]);
115639f9d9dSBarry Smith       for ( j=0; j<color->nrows[i]; j++ ) {
116639f9d9dSBarry Smith         printf("  %d %d \n",color->rows[i][j],color->columnsforrow[i][j]);
117bbf0e169SBarry Smith       }
118bbf0e169SBarry Smith     }
119bbf0e169SBarry Smith   }
120639f9d9dSBarry Smith   return 0;
121639f9d9dSBarry Smith }
122639f9d9dSBarry Smith 
1235615d1e5SSatish Balay #undef __FUNC__
1245615d1e5SSatish Balay #define __FUNC__ "MatFDColoringSetParameters"
125639f9d9dSBarry Smith /*@
126639f9d9dSBarry Smith    MatFDColoringSetParameters - Sets the parameters for the approximation of
127639f9d9dSBarry Smith    Jacobian using finite differences.
128639f9d9dSBarry Smith 
129639f9d9dSBarry Smith $       J(u)_{:,i} = [J(u+h*dx_{i}) - J(u)]/h where
130639f9d9dSBarry Smith $        h = error_rel*u[i]                    if  u[i] > umin
131639f9d9dSBarry Smith $          = error_rel*umin                    else
132639f9d9dSBarry Smith $
133639f9d9dSBarry Smith $   dx_{i} = (0, ... 1, .... 0)
134639f9d9dSBarry Smith 
135639f9d9dSBarry Smith    Input Parameters:
136639f9d9dSBarry Smith .  coloring - the jacobian coloring context
137639f9d9dSBarry Smith .  error_rel - relative error
138639f9d9dSBarry Smith .  umin - minimum allowable u-value
139639f9d9dSBarry Smith 
140639f9d9dSBarry Smith .keywords: SNES, Jacobian, finite differences, parameters
141639f9d9dSBarry Smith @*/
142639f9d9dSBarry Smith int MatFDColoringSetParameters(MatFDColoring matfd,double error,double umin)
143639f9d9dSBarry Smith {
144639f9d9dSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
145639f9d9dSBarry Smith 
146639f9d9dSBarry Smith   if (error != PETSC_DEFAULT) matfd->error_rel = error;
147639f9d9dSBarry Smith   if (umin != PETSC_DEFAULT)  matfd->umin      = umin;
148639f9d9dSBarry Smith   return 0;
149639f9d9dSBarry Smith }
150639f9d9dSBarry Smith 
1515615d1e5SSatish Balay #undef __FUNC__
152*005c665bSBarry Smith #define __FUNC__ "MatFDColoringSetFrequency"
153*005c665bSBarry Smith /*@
154*005c665bSBarry Smith    MatFDColoringSetFrequency - Sets the frequency at which new Jacobians
155*005c665bSBarry Smith      are computed.
156*005c665bSBarry Smith 
157*005c665bSBarry Smith    Input Parameters:
158*005c665bSBarry Smith .  coloring - the jacobian coloring context
159*005c665bSBarry Smith .  freq - frequency
160*005c665bSBarry Smith 
161*005c665bSBarry Smith .keywords: SNES, Jacobian, finite differences, parameters
162*005c665bSBarry Smith @*/
163*005c665bSBarry Smith int MatFDColoringSetFrequency(MatFDColoring matfd,int freq)
164*005c665bSBarry Smith {
165*005c665bSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
166*005c665bSBarry Smith 
167*005c665bSBarry Smith   matfd->freq = freq;
168*005c665bSBarry Smith   return 0;
169*005c665bSBarry Smith }
170*005c665bSBarry Smith 
171*005c665bSBarry Smith #undef __FUNC__
172*005c665bSBarry Smith #define __FUNC__ "MatFDColoringSetFunction"
173*005c665bSBarry Smith /*@
174*005c665bSBarry Smith    MatFDColoringSetFunction - Sets the function to use for computing the Jacobian.
175*005c665bSBarry Smith 
176*005c665bSBarry Smith    Input Parameters:
177*005c665bSBarry Smith .  coloring - the jacobian coloring context
178*005c665bSBarry Smith .  f - the function
179*005c665bSBarry Smith .  fctx - the function context
180*005c665bSBarry Smith 
181*005c665bSBarry Smith 
182*005c665bSBarry Smith .keywords: SNES, Jacobian, finite differences, parameters, function
183*005c665bSBarry Smith @*/
184*005c665bSBarry Smith int MatFDColoringSetFunction(MatFDColoring matfd,int (*f)(void *,Vec,Vec,void *),void *fctx)
185*005c665bSBarry Smith {
186*005c665bSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
187*005c665bSBarry Smith 
188*005c665bSBarry Smith   matfd->f    = f;
189*005c665bSBarry Smith   matfd->fctx = fctx;
190*005c665bSBarry Smith 
191*005c665bSBarry Smith   return 0;
192*005c665bSBarry Smith }
193*005c665bSBarry Smith 
194*005c665bSBarry Smith #undef __FUNC__
1955eea60f9SBarry Smith #define __FUNC__ "MatFDColoringSetFromOptions" /* ADIC Ignore */
196639f9d9dSBarry Smith /*@
197639f9d9dSBarry Smith    MatFDColoringSetFromOptions - Set coloring finite difference parameters from
198639f9d9dSBarry Smith          the options database.
199639f9d9dSBarry Smith 
200639f9d9dSBarry Smith $       J(u)_{:,i} = [J(u+h*dx_{i}) - J(u)]/h where
201639f9d9dSBarry Smith $        h = error_rel*u[i]                    if  u[i] > umin
202639f9d9dSBarry Smith $          = error_rel*.1                      else
203639f9d9dSBarry Smith $
204639f9d9dSBarry Smith $   dx_{i} = (0, ... 1, .... 0)
205639f9d9dSBarry Smith 
206639f9d9dSBarry Smith    Input Parameters:
207639f9d9dSBarry Smith .  coloring - the jacobian coloring context
208639f9d9dSBarry Smith 
209639f9d9dSBarry Smith    Options Database:
210639f9d9dSBarry Smith .  -mat_fd_coloring_error square root of relative error in function
211639f9d9dSBarry Smith .  -mat_fd_coloring_umin  see above
212*005c665bSBarry Smith .  -mat_fd_coloring_freq frequency at which Jacobian is computed
213639f9d9dSBarry Smith 
214639f9d9dSBarry Smith .keywords: SNES, Jacobian, finite differences, parameters
215639f9d9dSBarry Smith @*/
216639f9d9dSBarry Smith int MatFDColoringSetFromOptions(MatFDColoring matfd)
217639f9d9dSBarry Smith {
218*005c665bSBarry Smith   int    ierr,flag,freq = 1;
219639f9d9dSBarry Smith   double error = PETSC_DEFAULT,umin = PETSC_DEFAULT;
220639f9d9dSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
221639f9d9dSBarry Smith 
222639f9d9dSBarry Smith   ierr = OptionsGetDouble(matfd->prefix,"-mat_fd_coloring_err",&error,&flag);CHKERRQ(ierr);
223639f9d9dSBarry Smith   ierr = OptionsGetDouble(matfd->prefix,"-mat_fd_coloring_umin",&umin,&flag);CHKERRQ(ierr);
224639f9d9dSBarry Smith   ierr = MatFDColoringSetParameters(matfd,error,umin); CHKERRQ(ierr);
225*005c665bSBarry Smith   ierr = OptionsGetInt(matfd->prefix,"-mat_fd_coloring_freq",&freq,&flag);CHKERRQ(ierr);
226*005c665bSBarry Smith   ierr = MatFDColoringSetFrequency(matfd,freq);CHKERRQ(ierr);
227*005c665bSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help",&flag); CHKERRQ(ierr);
228639f9d9dSBarry Smith   if (flag) {
229639f9d9dSBarry Smith     ierr = MatFDColoringPrintHelp(matfd); CHKERRQ(ierr);
230639f9d9dSBarry Smith   }
231639f9d9dSBarry Smith   return 0;
232639f9d9dSBarry Smith }
233639f9d9dSBarry Smith 
2345615d1e5SSatish Balay #undef __FUNC__
2355eea60f9SBarry Smith #define __FUNC__ "MatFDColoringPrintHelp" /* ADIC Ignore */
236639f9d9dSBarry Smith /*@
237639f9d9dSBarry Smith     MatFDColoringPrintHelp - Prints help message for matrix finite difference calculations
238639f9d9dSBarry Smith          using coloring.
239639f9d9dSBarry Smith 
240639f9d9dSBarry Smith    Input Parameter:
241639f9d9dSBarry Smith .   fdcoloring - the MatFDColoring context
242639f9d9dSBarry Smith 
243639f9d9dSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringSetFromOptions()
244639f9d9dSBarry Smith @*/
245639f9d9dSBarry Smith int MatFDColoringPrintHelp(MatFDColoring fd)
246639f9d9dSBarry Smith {
247639f9d9dSBarry Smith   PetscValidHeaderSpecific(fd,MAT_FDCOLORING_COOKIE);
248639f9d9dSBarry Smith 
2490f665d81SLois Curfman McInnes   PetscPrintf(fd->comm,"-mat_fd_coloring_err <err>: set sqrt rel tol in function, defaults to 1.e-8\n");
2500f665d81SLois Curfman McInnes   PetscPrintf(fd->comm,"-mat_fd_coloring_umin <umin>: see users manual, defaults to 1.e-8\n");
251*005c665bSBarry Smith   PetscPrintf(fd->comm,"-mat_fd_coloring_freq <freq>: frequency that Jacobian is recomputed, defaults 1\n");
252*005c665bSBarry Smith   PetscPrintf(fd->comm,"-mat_fd_coloring_view\n");
253*005c665bSBarry Smith   PetscPrintf(fd->comm,"-mat_fd_coloring_view_draw\n");
254*005c665bSBarry Smith   return 0;
255*005c665bSBarry Smith }
256*005c665bSBarry Smith 
257*005c665bSBarry Smith int MatFDColoringView_Private(MatFDColoring fd)
258*005c665bSBarry Smith {
259*005c665bSBarry Smith   int ierr,flg;
260*005c665bSBarry Smith 
261*005c665bSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view",&flg); CHKERRQ(ierr);
262*005c665bSBarry Smith   if (flg) {
263*005c665bSBarry Smith     Viewer viewer;
264*005c665bSBarry Smith     ierr = ViewerFileOpenASCII(fd->comm,"stdout",&viewer);CHKERRQ(ierr);
265*005c665bSBarry Smith     ierr = MatFDColoringView(fd,viewer); CHKERRQ(ierr);
266*005c665bSBarry Smith     ierr = ViewerDestroy(viewer); CHKERRQ(ierr);
267*005c665bSBarry Smith   }
268*005c665bSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg); CHKERRQ(ierr);
269*005c665bSBarry Smith   if (flg) {
270*005c665bSBarry Smith     ierr = MatFDColoringView(fd,VIEWER_DRAWX_(fd->comm)); CHKERRQ(ierr);
271*005c665bSBarry Smith     ierr = ViewerFlush(VIEWER_DRAWX_(fd->comm)); CHKERRQ(ierr);
272*005c665bSBarry Smith   }
273bbf0e169SBarry Smith   return 0;
274bbf0e169SBarry Smith }
275bbf0e169SBarry Smith 
2765615d1e5SSatish Balay #undef __FUNC__
2775615d1e5SSatish Balay #define __FUNC__ "MatFDColoringCreate"
278bbf0e169SBarry Smith /*@C
279639f9d9dSBarry Smith    MatFDColoringCreate - Creates a matrix coloring context for finite difference
280639f9d9dSBarry Smith         computation of Jacobians.
281bbf0e169SBarry Smith 
282639f9d9dSBarry Smith    Input Parameters:
283639f9d9dSBarry Smith .    mat - the matrix containing the nonzero structure of the Jacobian
284639f9d9dSBarry Smith .    iscoloring - the coloring of the matrix
285bbf0e169SBarry Smith 
286bbf0e169SBarry Smith    Output Parameter:
287639f9d9dSBarry Smith .   color - the new coloring context
288bbf0e169SBarry Smith 
289*005c665bSBarry Smith    Options Database:
290*005c665bSBarry Smith .  -mat_fd_coloring_error square root of relative error in function
291*005c665bSBarry Smith .  -mat_fd_coloring_umin  see above
292*005c665bSBarry Smith .  -mat_fd_coloring_view
293*005c665bSBarry Smith .  -mat_fd_coloring_view_draw
294639f9d9dSBarry Smith 
295639f9d9dSBarry Smith .seealso: MatFDColoringDestroy()
296bbf0e169SBarry Smith @*/
297639f9d9dSBarry Smith int MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
298bbf0e169SBarry Smith {
299639f9d9dSBarry Smith   MatFDColoring c;
300639f9d9dSBarry Smith   MPI_Comm      comm;
301639f9d9dSBarry Smith   int           ierr,M,N;
302639f9d9dSBarry Smith 
303639f9d9dSBarry Smith   ierr = MatGetSize(mat,&M,&N); CHKERRQ(ierr);
304e3372554SBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,0,"Only for square matrices");
305639f9d9dSBarry Smith 
306639f9d9dSBarry Smith   PetscObjectGetComm((PetscObject)mat,&comm);
307f09e8eb9SSatish Balay   PetscHeaderCreate(c,_p_MatFDColoring,MAT_FDCOLORING_COOKIE,0,comm);
308639f9d9dSBarry Smith   PLogObjectCreate(c);
309639f9d9dSBarry Smith 
310639f9d9dSBarry Smith   if (mat->ops.fdcoloringcreate) {
311639f9d9dSBarry Smith     ierr = (*mat->ops.fdcoloringcreate)(mat,iscoloring,c); CHKERRQ(ierr);
312639f9d9dSBarry Smith   } else {
313e3372554SBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"Code not yet written for this matrix type");
314639f9d9dSBarry Smith   }
315639f9d9dSBarry Smith 
316639f9d9dSBarry Smith   c->error_rel = 1.e-8;
317639f9d9dSBarry Smith   c->umin      = 1.e-5;
318*005c665bSBarry Smith   c->freq      = 1;
319*005c665bSBarry Smith 
320*005c665bSBarry Smith   ierr = MatFDColoringView_Private(c); CHKERRQ(ierr);
321639f9d9dSBarry Smith 
322639f9d9dSBarry Smith   *color = c;
323639f9d9dSBarry Smith 
324bbf0e169SBarry Smith   return 0;
325bbf0e169SBarry Smith }
326bbf0e169SBarry Smith 
3275615d1e5SSatish Balay #undef __FUNC__
3285eea60f9SBarry Smith #define __FUNC__ "MatFDColoringDestroy" /* ADIC Ignore */
329bbf0e169SBarry Smith /*@C
330639f9d9dSBarry Smith     MatFDColoringDestroy - Destroys a matrix coloring context that was created
331639f9d9dSBarry Smith          via MatFDColoringCreate().
332bbf0e169SBarry Smith 
333639f9d9dSBarry Smith     Input Paramter:
334639f9d9dSBarry Smith .   c - coloring context
335bbf0e169SBarry Smith 
336639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
337bbf0e169SBarry Smith @*/
338639f9d9dSBarry Smith int MatFDColoringDestroy(MatFDColoring c)
339bbf0e169SBarry Smith {
340639f9d9dSBarry Smith   int i,ierr,flag;
341bbf0e169SBarry Smith 
342639f9d9dSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-matfdcoloring_view",&flag);
343639f9d9dSBarry Smith   if (flag) {
344bbf0e169SBarry Smith     Viewer viewer;
345639f9d9dSBarry Smith     ierr = ViewerFileOpenASCII(c->comm,"stdout",&viewer);CHKERRQ(ierr);
346639f9d9dSBarry Smith     ierr = MatFDColoringView(c,viewer);CHKERRQ(ierr);
347bbf0e169SBarry Smith     ierr = ViewerDestroy(viewer); CHKERRQ(ierr);
348bbf0e169SBarry Smith   }
349639f9d9dSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-matfdcoloring_view_info",&flag);
350639f9d9dSBarry Smith   if (flag) {
351bbf0e169SBarry Smith     Viewer viewer;
352639f9d9dSBarry Smith     ierr = ViewerFileOpenASCII(c->comm,"stdout",&viewer);CHKERRQ(ierr);
353639f9d9dSBarry Smith     ierr = ViewerPushFormat(viewer,VIEWER_FORMAT_ASCII_INFO, PETSC_NULL); CHKERRQ(ierr);
354639f9d9dSBarry Smith     ierr = MatFDColoringView(c,viewer);CHKERRQ(ierr);
355639f9d9dSBarry Smith     ierr = ViewerPopFormat(viewer);
356bbf0e169SBarry Smith     ierr = ViewerDestroy(viewer); CHKERRQ(ierr);
357bbf0e169SBarry Smith   }
358639f9d9dSBarry Smith 
359639f9d9dSBarry Smith   for ( i=0; i<c->ncolors; i++ ) {
360639f9d9dSBarry Smith     if (c->columns[i])       PetscFree(c->columns[i]);
361639f9d9dSBarry Smith     if (c->rows[i])          PetscFree(c->rows[i]);
362639f9d9dSBarry Smith     if (c->columnsforrow[i]) PetscFree(c->columnsforrow[i]);
363bbf0e169SBarry Smith   }
364639f9d9dSBarry Smith   PetscFree(c->ncolumns);
365639f9d9dSBarry Smith   PetscFree(c->columns);
366639f9d9dSBarry Smith   PetscFree(c->nrows);
367639f9d9dSBarry Smith   PetscFree(c->rows);
368639f9d9dSBarry Smith   PetscFree(c->columnsforrow);
369639f9d9dSBarry Smith   PetscFree(c->scale);
370*005c665bSBarry Smith   if (c->w1) {
371*005c665bSBarry Smith     ierr = VecDestroy(c->w1); CHKERRQ(ierr);
372*005c665bSBarry Smith     ierr = VecDestroy(c->w2); CHKERRQ(ierr);
373*005c665bSBarry Smith     ierr = VecDestroy(c->w3); CHKERRQ(ierr);
374*005c665bSBarry Smith   }
375639f9d9dSBarry Smith   PLogObjectDestroy(c);
376639f9d9dSBarry Smith   PetscHeaderDestroy(c);
377bbf0e169SBarry Smith   return 0;
378bbf0e169SBarry Smith }
37943a90d84SBarry Smith 
380*005c665bSBarry Smith #include "snes.h"
381*005c665bSBarry Smith 
3825615d1e5SSatish Balay #undef __FUNC__
3835615d1e5SSatish Balay #define __FUNC__ "MatFDColoringApply"
38443a90d84SBarry Smith /*@
38543a90d84SBarry Smith      MatFDColoringApply - Given a matrix for which a MatFDColoring has been created,
38643a90d84SBarry Smith          computes the Jacobian for a function via finite differences.
38743a90d84SBarry Smith 
38843a90d84SBarry Smith   Input Parameters:
38943a90d84SBarry Smith .   mat - location to store Jacobian
39043a90d84SBarry Smith .   coloring - coloring context created with MatFDColoringCreate()
39143a90d84SBarry Smith .   x1 - location at which Jacobian is to be computed
392*005c665bSBarry Smith .   sctx - optional context required by function (actually a SNES context)
39343a90d84SBarry Smith 
39443a90d84SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
39543a90d84SBarry Smith 
39643a90d84SBarry Smith .keywords: coloring, Jacobian, finite differences
39743a90d84SBarry Smith @*/
398*005c665bSBarry Smith int MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
39943a90d84SBarry Smith {
400*005c665bSBarry Smith   int           k, ierr,N,start,end,l,row,col,srow, it;
40143a90d84SBarry Smith   Scalar        dx, mone = -1.0,*y,*scale = coloring->scale,*xx,*wscale = coloring->wscale;
40243a90d84SBarry Smith   double        epsilon = coloring->error_rel, umin = coloring->umin;
40343a90d84SBarry Smith   MPI_Comm      comm = coloring->comm;
404*005c665bSBarry Smith   Vec           w1,w2,w3;
405*005c665bSBarry Smith   int           (*f)(void *,Vec,Vec,void *) = coloring->f,freq = coloring->freq;
406*005c665bSBarry Smith   void          *fctx = coloring->fctx;
407*005c665bSBarry Smith 
408*005c665bSBarry Smith   ierr = SNESGetIterationNumber((SNES)sctx,&it); CHKERRQ(ierr);
409*005c665bSBarry Smith   if ((freq > 1) && ((it % freq) != 1)) {
410*005c665bSBarry Smith     PLogInfo(coloring,"MatFDColoringApply:Skipping Jacobian, iteration %d, freq %d\n",it,freq);
411*005c665bSBarry Smith     *flag = SAME_PRECONDITIONER;
412*005c665bSBarry Smith     return 0;
413*005c665bSBarry Smith   } else {
414*005c665bSBarry Smith     PLogInfo(coloring,"MatFDColoringApply:Computing Jacobian, iteration %d, freq %d\n",it,freq);
415*005c665bSBarry Smith     *flag = SAME_NONZERO_PATTERN;
416*005c665bSBarry Smith   }
417*005c665bSBarry Smith 
418*005c665bSBarry Smith 
419*005c665bSBarry Smith 
420*005c665bSBarry Smith   if (!coloring->w1) {
421*005c665bSBarry Smith     ierr = VecDuplicate(x1,&coloring->w1); CHKERRQ(ierr);
422*005c665bSBarry Smith     PLogObjectParent(coloring,coloring->w1);
423*005c665bSBarry Smith     ierr = VecDuplicate(x1,&coloring->w2); CHKERRQ(ierr);
424*005c665bSBarry Smith     PLogObjectParent(coloring,coloring->w2);
425*005c665bSBarry Smith     ierr = VecDuplicate(x1,&coloring->w3); CHKERRQ(ierr);
426*005c665bSBarry Smith     PLogObjectParent(coloring,coloring->w3);
427*005c665bSBarry Smith   }
428*005c665bSBarry Smith   w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
42943a90d84SBarry Smith 
43043a90d84SBarry Smith   ierr = MatZeroEntries(J); CHKERRQ(ierr);
43143a90d84SBarry Smith 
43243a90d84SBarry Smith   ierr = VecGetOwnershipRange(x1,&start,&end); CHKERRQ(ierr);
43343a90d84SBarry Smith   ierr = VecGetSize(x1,&N); CHKERRQ(ierr);
43443a90d84SBarry Smith   ierr = VecGetArray(x1,&xx); CHKERRQ(ierr);
43543a90d84SBarry Smith   ierr = (*f)(sctx,x1,w1,fctx); CHKERRQ(ierr);
43643a90d84SBarry Smith 
43743a90d84SBarry Smith   PetscMemzero(wscale,N*sizeof(Scalar));
43843a90d84SBarry Smith   /*
43943a90d84SBarry Smith       Loop over each color
44043a90d84SBarry Smith   */
44143a90d84SBarry Smith 
44243a90d84SBarry Smith   for (k=0; k<coloring->ncolors; k++) {
44343a90d84SBarry Smith     ierr = VecCopy(x1,w3); CHKERRQ(ierr);
44443a90d84SBarry Smith     /*
44543a90d84SBarry Smith        Loop over each column associated with color adding the
44643a90d84SBarry Smith        perturbation to the vector w3.
44743a90d84SBarry Smith     */
44843a90d84SBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
44943a90d84SBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
45043a90d84SBarry Smith       dx  = xx[col-start];
45143a90d84SBarry Smith #if !defined(PETSC_COMPLEX)
45243a90d84SBarry Smith       if (dx < umin && dx >= 0.0) dx = .1;
45343a90d84SBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -.1;
45443a90d84SBarry Smith #else
45543a90d84SBarry Smith       if (abs(dx) < umin && real(dx) >= 0.0) dx = .1;
45643a90d84SBarry Smith       else if (real(dx) < 0.0 && abs(dx) < umin) dx = -.1;
45743a90d84SBarry Smith #endif
45843a90d84SBarry Smith       dx          *= epsilon;
45943a90d84SBarry Smith       wscale[col] = 1.0/dx;
46043a90d84SBarry Smith       VecSetValues(w3,1,&col,&dx,ADD_VALUES);
46143a90d84SBarry Smith     }
46243a90d84SBarry Smith     VecRestoreArray(x1,&xx);
46343a90d84SBarry Smith     /*
46443a90d84SBarry Smith        Evaluate function at x1 + dx (here dx is a vector, of perturbations)
46543a90d84SBarry Smith     */
46643a90d84SBarry Smith     ierr = (*f)(sctx,w3,w2,fctx); CHKERRQ(ierr);
46743a90d84SBarry Smith     ierr = VecAXPY(&mone,w1,w2); CHKERRQ(ierr);
46843a90d84SBarry Smith     /* Communicate scale to all processors */
46943a90d84SBarry Smith #if !defined(PETSC_COMPLEX)
47043a90d84SBarry Smith     MPI_Allreduce(wscale,scale,N,MPI_DOUBLE,MPI_SUM,comm);
47143a90d84SBarry Smith #else
47243a90d84SBarry Smith     MPI_Allreduce(wscale,scale,2*N,MPI_DOUBLE,MPI_SUM,comm);
47343a90d84SBarry Smith #endif
47443a90d84SBarry Smith     /*
47543a90d84SBarry Smith        Loop over rows of vector putting results into Jacobian matrix
47643a90d84SBarry Smith     */
47743a90d84SBarry Smith     VecGetArray(w2,&y);
47843a90d84SBarry Smith     for (l=0; l<coloring->nrows[k]; l++) {
47943a90d84SBarry Smith       row    = coloring->rows[k][l];
48043a90d84SBarry Smith       col    = coloring->columnsforrow[k][l];
48143a90d84SBarry Smith       y[row] *= scale[col];
48243a90d84SBarry Smith       srow   = row + start;
48343a90d84SBarry Smith       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
48443a90d84SBarry Smith     }
48543a90d84SBarry Smith     VecRestoreArray(w2,&y);
48643a90d84SBarry Smith   }
48743a90d84SBarry Smith   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
48843a90d84SBarry Smith   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
48943a90d84SBarry Smith   return 0;
49043a90d84SBarry Smith }
491