xref: /petsc/src/mat/matfd/fdmatrix.c (revision 2bdab257d34f9d407eac123c854c8b30dbb6719c)
1 #ifdef PETSC_RCS_HEADER
2 static char vcid[] = "$Id: fdmatrix.c,v 1.19 1997/10/01 22:45:08 bsmith Exp bsmith $";
3 #endif
4 
5 /*
6    This is where the abstract matrix operations are defined that are
7   used for finite difference computations of Jacobians using coloring.
8 */
9 
10 #include "petsc.h"
11 #include "src/mat/matimpl.h"        /*I "mat.h" I*/
12 #include "src/vec/vecimpl.h"
13 #include "pinclude/pviewer.h"
14 
15 #undef __FUNC__
16 #define __FUNC__ "MatFDColoringView_Draw"
17 static int MatFDColoringView_Draw(MatFDColoring fd,Viewer viewer)
18 {
19   int         ierr,i,j,pause;
20   PetscTruth  isnull;
21   Draw        draw;
22   double      xr,yr,xl,yl,h,w,x,y,xc,yc,scale = 0.0;
23   DrawButton  button;
24 
25   ierr = ViewerDrawGetDraw(viewer,&draw); CHKERRQ(ierr);
26   ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); if (isnull) return 0;
27   ierr = DrawSynchronizedClear(draw); CHKERRQ(ierr);
28 
29   xr  = fd->N; yr = fd->M; h = yr/10.0; w = xr/10.0;
30   xr += w;    yr += h;  xl = -w;     yl = -h;
31   ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
32 
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       DrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);
39     }
40   }
41   ierr = DrawSynchronizedFlush(draw); CHKERRQ(ierr);
42   ierr = DrawGetPause(draw,&pause); CHKERRQ(ierr);
43   if (pause >= 0) { PetscSleep(pause); return 0;}
44   ierr = DrawCheckResizedWindow(draw);
45   ierr = DrawSynchronizedGetMouseButton(draw,&button,&xc,&yc,0,0);
46   while (button != BUTTON_RIGHT) {
47     ierr = DrawSynchronizedClear(draw); CHKERRQ(ierr);
48     if (button == BUTTON_LEFT) scale = .5;
49     else if (button == BUTTON_CENTER) scale = 2.;
50     xl = scale*(xl + w - xc) + xc - w*scale;
51     xr = scale*(xr - w - xc) + xc + w*scale;
52     yl = scale*(yl + h - yc) + yc - h*scale;
53     yr = scale*(yr - h - yc) + yc + h*scale;
54     w *= scale; h *= scale;
55     ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
56     /* loop over colors  */
57     for (i=0; i<fd->ncolors; i++ ) {
58       for ( j=0; j<fd->nrows[i]; j++ ) {
59         y = fd->M - fd->rows[i][j] - fd->rstart;
60         x = fd->columnsforrow[i][j];
61         DrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);
62       }
63     }
64     ierr = DrawCheckResizedWindow(draw);
65     ierr = DrawSynchronizedGetMouseButton(draw,&button,&xc,&yc,0,0);
66   }
67 
68   return 0;
69 }
70 
71 #undef __FUNC__
72 #define __FUNC__ "MatFDColoringView"
73 /*@C
74    MatFDColoringView - Views a finite difference coloring context.
75 
76    Input  Parameters:
77 .  c - the coloring context
78 .  viewer - visualization context
79 
80    Notes:
81    The available visualization contexts include
82 $     VIEWER_STDOUT_SELF - standard output (default)
83 $     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 $     VIEWER_DRAWX_WORLD - graphical display of nonzero structure
88 
89 .seealso: MatFDColoringCreate()
90 
91 .keywords: Mat, finite differences, coloring, view
92 @*/
93 int MatFDColoringView(MatFDColoring c,Viewer viewer)
94 {
95   ViewerType vtype;
96   int        i,j,format,ierr;
97   FILE       *fd;
98 
99   PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE);
100   if (viewer) {PetscValidHeader(viewer);}
101   else {viewer = VIEWER_STDOUT_SELF;}
102 
103   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
104   if (vtype == DRAW_VIEWER) {
105     ierr = MatFDColoringView_Draw(c,viewer); CHKERRQ(ierr);
106     return 0;
107   }
108   else if (vtype  == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) {
109     ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
110     PetscFPrintf(c->comm,fd,"MatFDColoring Object:\n");
111     PetscFPrintf(c->comm,fd,"  Error tolerance=%g\n",c->error_rel);
112     PetscFPrintf(c->comm,fd,"  Umin=%g\n",c->umin);
113     PetscFPrintf(c->comm,fd,"  Number of colors=%d\n",c->ncolors);
114 
115     ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr);
116     if (format != VIEWER_FORMAT_ASCII_INFO) {
117       for ( i=0; i<c->ncolors; i++ ) {
118         PetscFPrintf(c->comm,fd,"  Information for color %d\n",i);
119         PetscFPrintf(c->comm,fd,"    Number of columns %d\n",c->ncolumns[i]);
120         for ( j=0; j<c->ncolumns[i]; j++ ) {
121           PetscFPrintf(c->comm,fd,"      %d\n",c->columns[i][j]);
122         }
123         PetscFPrintf(c->comm,fd,"    Number of rows %d\n",c->nrows[i]);
124         for ( j=0; j<c->nrows[i]; j++ ) {
125           PetscFPrintf(c->comm,fd,"      %d %d \n",c->rows[i][j],c->columnsforrow[i][j]);
126         }
127       }
128     }
129   }
130   return 0;
131 }
132 
133 #undef __FUNC__
134 #define __FUNC__ "MatFDColoringSetParameters"
135 /*@
136    MatFDColoringSetParameters - Sets the parameters for the sparse approximation of
137    a Jacobian matrix using finite differences.
138 
139 $       J(u)_{:,i} = [J(u+h*dx_{i}) - J(u)]/h where
140 $        h = error_rel*u[i]                    if  u[i] > umin
141 $          = error_rel*umin                    else
142 $
143 $   dx_{i} = (0, ... 1, .... 0)
144 
145    Input Parameters:
146 .  coloring - the coloring context
147 .  error_rel - relative error
148 .  umin - minimum allowable u-value
149 
150 .keywords: Mat, finite differences, coloring, set, parameters
151 
152 .seealso: MatFDColoringCreate()
153 @*/
154 int MatFDColoringSetParameters(MatFDColoring matfd,double error,double umin)
155 {
156   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
157 
158   if (error != PETSC_DEFAULT) matfd->error_rel = error;
159   if (umin != PETSC_DEFAULT)  matfd->umin      = umin;
160   return 0;
161 }
162 
163 #undef __FUNC__
164 #define __FUNC__ "MatFDColoringSetFrequency"
165 /*@
166    MatFDColoringSetFrequency - Sets the frequency for computing new Jacobian
167    matrices.
168 
169    Input Parameters:
170 .  coloring - the coloring context
171 .  freq - frequency (default is 1)
172 
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.
178 
179    Options Database Keys:
180 $  -mat_fd_coloring_freq <freq>
181 
182 .keywords: Mat, finite differences, coloring, set, frequency
183 @*/
184 int MatFDColoringSetFrequency(MatFDColoring matfd,int freq)
185 {
186   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
187 
188   matfd->freq = freq;
189   return 0;
190 }
191 
192 #undef __FUNC__
193 #define __FUNC__ "MatFDColoringSetFunction"
194 /*@
195    MatFDColoringSetFunction - Sets the function to use for computing the Jacobian.
196 
197    Input Parameters:
198 .  coloring - the coloring context
199 .  f - the function
200 .  fctx - the function context
201 
202 .keywords: Mat, Jacobian, finite differences, set, function
203 @*/
204 int MatFDColoringSetFunction(MatFDColoring matfd,int (*f)(void *,Vec,Vec,void *),void *fctx)
205 {
206   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
207 
208   matfd->f    = f;
209   matfd->fctx = fctx;
210 
211   return 0;
212 }
213 
214 #undef __FUNC__
215 #define __FUNC__ "MatFDColoringSetFromOptions"
216 /*@
217    MatFDColoringSetFromOptions - Sets coloring finite difference parameters from
218    the options database.
219 
220    The Jacobian is estimated with the differencing approximation
221 $       J(u)_{:,i} = [J(u+h*dx_{i}) - J(u)]/h where
222 $        h = error_rel*u[i]                    if  u[i] > umin
223 $          = error_rel*umin                      else
224 $
225 $   dx_{i} = (0, ... 1, .... 0)
226 
227    Input Parameters:
228 .  coloring - the coloring context
229 
230    Options Database Keys:
231 $  -mat_fd_coloring_error <err>, where <err> is the square root
232 $           of relative error in the function
233 $  -mat_fd_coloring_umin  <umin>, where umin is described above
234 $  -mat_fd_coloring_freq <freq> where <freq> is the frequency of
235 $           computing a new Jacobian
236 $  -mat_fd_coloring_view
237 $  -mat_fd_coloring_view_info
238 $  -mat_fd_coloring_view_draw
239 
240 .keywords: Mat, finite differences, parameters
241 @*/
242 int MatFDColoringSetFromOptions(MatFDColoring matfd)
243 {
244   int    ierr,flag,freq = 1;
245   double error = PETSC_DEFAULT,umin = PETSC_DEFAULT;
246   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
247 
248   ierr = OptionsGetDouble(matfd->prefix,"-mat_fd_coloring_err",&error,&flag);CHKERRQ(ierr);
249   ierr = OptionsGetDouble(matfd->prefix,"-mat_fd_coloring_umin",&umin,&flag);CHKERRQ(ierr);
250   ierr = MatFDColoringSetParameters(matfd,error,umin); CHKERRQ(ierr);
251   ierr = OptionsGetInt(matfd->prefix,"-mat_fd_coloring_freq",&freq,&flag);CHKERRQ(ierr);
252   ierr = MatFDColoringSetFrequency(matfd,freq);CHKERRQ(ierr);
253   ierr = OptionsHasName(PETSC_NULL,"-help",&flag); CHKERRQ(ierr);
254   if (flag) {
255     ierr = MatFDColoringPrintHelp(matfd); CHKERRQ(ierr);
256   }
257   return 0;
258 }
259 
260 #undef __FUNC__
261 #define __FUNC__ "MatFDColoringPrintHelp"
262 /*@
263     MatFDColoringPrintHelp - Prints help message for matrix finite difference calculations
264     using coloring.
265 
266     Input Parameter:
267 .   fdcoloring - the MatFDColoring context
268 
269 .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringSetFromOptions()
270 @*/
271 int MatFDColoringPrintHelp(MatFDColoring fd)
272 {
273   PetscValidHeaderSpecific(fd,MAT_FDCOLORING_COOKIE);
274 
275   PetscPrintf(fd->comm,"-mat_fd_coloring_err <err>: set sqrt rel tol in function, defaults to %g\n",fd->error_rel);
276   PetscPrintf(fd->comm,"-mat_fd_coloring_umin <umin>: see users manual, defaults to %d\n",fd->umin);
277   PetscPrintf(fd->comm,"-mat_fd_coloring_freq <freq>: frequency that Jacobian is recomputed, defaults to %d\n",fd->freq);
278   PetscPrintf(fd->comm,"-mat_fd_coloring_view\n");
279   PetscPrintf(fd->comm,"-mat_fd_coloring_view_draw\n");
280   PetscPrintf(fd->comm,"-mat_fd_coloring_view_info\n");
281   return 0;
282 }
283 
284 int MatFDColoringView_Private(MatFDColoring fd)
285 {
286   int ierr,flg;
287 
288   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view",&flg); CHKERRQ(ierr);
289   if (flg) {
290     ierr = MatFDColoringView(fd,VIEWER_STDOUT_(fd->comm)); CHKERRQ(ierr);
291   }
292   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_info",&flg); CHKERRQ(ierr);
293   if (flg) {
294     ierr = ViewerPushFormat(VIEWER_STDOUT_(fd->comm),VIEWER_FORMAT_ASCII_INFO,PETSC_NULL);CHKERRQ(ierr);
295     ierr = MatFDColoringView(fd,VIEWER_STDOUT_(fd->comm)); CHKERRQ(ierr);
296     ierr = ViewerPopFormat(VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr);
297   }
298   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg); CHKERRQ(ierr);
299   if (flg) {
300     ierr = MatFDColoringView(fd,VIEWER_DRAWX_(fd->comm)); CHKERRQ(ierr);
301     ierr = ViewerFlush(VIEWER_DRAWX_(fd->comm)); CHKERRQ(ierr);
302   }
303   return 0;
304 }
305 
306 #undef __FUNC__
307 #define __FUNC__ "MatFDColoringCreate"
308 /*@C
309    MatFDColoringCreate - Creates a matrix coloring context for finite difference
310    computation of Jacobians.
311 
312    Input Parameters:
313 .  mat - the matrix containing the nonzero structure of the Jacobian
314 .  iscoloring - the coloring of the matrix
315 
316     Output Parameter:
317 .   color - the new coloring context
318 
319     Options Database Keys:
320 $    -mat_fd_coloring_view
321 $    -mat_fd_coloring_view_draw
322 $    -mat_fd_coloring_view_info
323 
324 .seealso: MatFDColoringDestroy()
325 @*/
326 int MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
327 {
328   MatFDColoring c;
329   MPI_Comm      comm;
330   int           ierr,M,N;
331 
332   ierr = MatGetSize(mat,&M,&N); CHKERRQ(ierr);
333   if (M != N) SETERRQ(PETSC_ERR_SUP,0,"Only for square matrices");
334 
335   PetscObjectGetComm((PetscObject)mat,&comm);
336   PetscHeaderCreate(c,_p_MatFDColoring,MAT_FDCOLORING_COOKIE,0,comm,MatFDColoringDestroy,MatFDColoringView);
337   PLogObjectCreate(c);
338 
339   if (mat->ops.fdcoloringcreate) {
340     ierr = (*mat->ops.fdcoloringcreate)(mat,iscoloring,c); CHKERRQ(ierr);
341   } else {
342     SETERRQ(PETSC_ERR_SUP,0,"Code not yet written for this matrix type");
343   }
344 
345   c->error_rel = 1.e-8;
346   c->umin      = 1.e-6;
347   c->freq      = 1;
348 
349   ierr = MatFDColoringView_Private(c); CHKERRQ(ierr);
350 
351   *color = c;
352 
353   return 0;
354 }
355 
356 #undef __FUNC__
357 #define __FUNC__ "MatFDColoringDestroy"
358 /*@C
359     MatFDColoringDestroy - Destroys a matrix coloring context that was created
360     via MatFDColoringCreate().
361 
362     Input Parameter:
363 .   c - coloring context
364 
365 .seealso: MatFDColoringCreate()
366 @*/
367 int MatFDColoringDestroy(MatFDColoring c)
368 {
369   int i,ierr,flag;
370 
371   if (--c->refct > 0) return 0;
372 
373   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view",&flag);
374   if (flag) {
375     ierr = MatFDColoringView(c,VIEWER_STDOUT_(c->comm));CHKERRQ(ierr);
376   }
377   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_info",&flag);
378   if (flag) {
379     ierr = ViewerPushFormat(VIEWER_STDOUT_(c->comm),VIEWER_FORMAT_ASCII_INFO,PETSC_NULL);CHKERRQ(ierr);
380     ierr = MatFDColoringView(c,VIEWER_STDOUT_(c->comm));CHKERRQ(ierr);
381     ierr = ViewerPopFormat(VIEWER_STDOUT_(c->comm));
382   }
383 
384   for ( i=0; i<c->ncolors; i++ ) {
385     if (c->columns[i])       PetscFree(c->columns[i]);
386     if (c->rows[i])          PetscFree(c->rows[i]);
387     if (c->columnsforrow[i]) PetscFree(c->columnsforrow[i]);
388   }
389   PetscFree(c->ncolumns);
390   PetscFree(c->columns);
391   PetscFree(c->nrows);
392   PetscFree(c->rows);
393   PetscFree(c->columnsforrow);
394   PetscFree(c->scale);
395   if (c->w1) {
396     ierr = VecDestroy(c->w1); CHKERRQ(ierr);
397     ierr = VecDestroy(c->w2); CHKERRQ(ierr);
398     ierr = VecDestroy(c->w3); CHKERRQ(ierr);
399   }
400   PLogObjectDestroy(c);
401   PetscHeaderDestroy(c);
402   return 0;
403 }
404 
405 #include "snes.h"
406 
407 #undef __FUNC__
408 #define __FUNC__ "MatFDColoringApply"
409 /*@
410     MatFDColoringApply - Given a matrix for which a MatFDColoring context
411     has been created, computes the Jacobian for a function via finite differences.
412 
413     Input Parameters:
414 .   mat - location to store Jacobian
415 .   coloring - coloring context created with MatFDColoringCreate()
416 .   x1 - location at which Jacobian is to be computed
417 .   sctx - optional context required by function (actually a SNES context)
418 
419 .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
420 
421 .keywords: coloring, Jacobian, finite differences
422 @*/
423 int MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
424 {
425   int           k,fg,ierr,N,start,end,l,row,col,srow;
426   Scalar        dx, mone = -1.0,*y,*scale = coloring->scale,*xx,*wscale = coloring->wscale;
427   double        epsilon = coloring->error_rel, umin = coloring->umin;
428   MPI_Comm      comm = coloring->comm;
429   Vec           w1,w2,w3;
430   int           (*f)(void *,Vec,Vec,void *) = coloring->f;
431   void          *fctx = coloring->fctx;
432 
433   PetscValidHeaderSpecific(J,MAT_COOKIE);
434   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE);
435   PetscValidHeaderSpecific(x1,VEC_COOKIE);
436 
437   /*
438   ierr = SNESGetIterationNumber((SNES)sctx,&it); CHKERRQ(ierr);
439   if ((freq > 1) && ((it % freq) != 1)) {
440     PLogInfo(coloring,"MatFDColoringApply:Skipping Jacobian, iteration %d, freq %d\n",it,freq);
441     *flag = SAME_PRECONDITIONER;
442     return 0;
443   } else {
444     PLogInfo(coloring,"MatFDColoringApply:Computing Jacobian, iteration %d, freq %d\n",it,freq);
445     *flag = SAME_NONZERO_PATTERN;
446   }*/
447 
448   if (!coloring->w1) {
449     ierr = VecDuplicate(x1,&coloring->w1); CHKERRQ(ierr);
450     PLogObjectParent(coloring,coloring->w1);
451     ierr = VecDuplicate(x1,&coloring->w2); CHKERRQ(ierr);
452     PLogObjectParent(coloring,coloring->w2);
453     ierr = VecDuplicate(x1,&coloring->w3); CHKERRQ(ierr);
454     PLogObjectParent(coloring,coloring->w3);
455   }
456   w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
457 
458   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&fg); CHKERRQ(ierr);
459   if (fg) {
460     PLogInfo(coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n");
461   } else {
462     ierr = MatZeroEntries(J); CHKERRQ(ierr);
463   }
464 
465   ierr = VecGetOwnershipRange(x1,&start,&end); CHKERRQ(ierr);
466   ierr = VecGetSize(x1,&N); CHKERRQ(ierr);
467   ierr = VecGetArray(x1,&xx); CHKERRQ(ierr);
468   ierr = (*f)(sctx,x1,w1,fctx); CHKERRQ(ierr);
469 
470   PetscMemzero(wscale,N*sizeof(Scalar));
471   /*
472       Loop over each color
473   */
474 
475   for (k=0; k<coloring->ncolors; k++) {
476     ierr = VecCopy(x1,w3); CHKERRQ(ierr);
477     /*
478        Loop over each column associated with color adding the
479        perturbation to the vector w3.
480     */
481     for (l=0; l<coloring->ncolumns[k]; l++) {
482       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
483       dx  = xx[col-start];
484       if (dx == 0.0) dx = 1.0;
485 #if !defined(PETSC_COMPLEX)
486       if (dx < umin && dx >= 0.0)      dx = umin;
487       else if (dx < 0.0 && dx > -umin) dx = -umin;
488 #else
489       if (abs(dx) < umin && real(dx) >= 0.0)     dx = umin;
490       else if (real(dx) < 0.0 && abs(dx) < umin) dx = -umin;
491 #endif
492       dx          *= epsilon;
493       wscale[col] = 1.0/dx;
494       VecSetValues(w3,1,&col,&dx,ADD_VALUES);
495     }
496     VecRestoreArray(x1,&xx);
497     /*
498        Evaluate function at x1 + dx (here dx is a vector of perturbations)
499     */
500     ierr = (*f)(sctx,w3,w2,fctx); CHKERRQ(ierr);
501     ierr = VecAXPY(&mone,w1,w2); CHKERRQ(ierr);
502     /* Communicate scale to all processors */
503 #if !defined(PETSC_COMPLEX)
504     MPI_Allreduce(wscale,scale,N,MPI_DOUBLE,MPI_SUM,comm);
505 #else
506     MPI_Allreduce(wscale,scale,2*N,MPI_DOUBLE,MPI_SUM,comm);
507 #endif
508     /*
509        Loop over rows of vector, putting results into Jacobian matrix
510     */
511     VecGetArray(w2,&y);
512     for (l=0; l<coloring->nrows[k]; l++) {
513       row    = coloring->rows[k][l];
514       col    = coloring->columnsforrow[k][l];
515       y[row] *= scale[col];
516       srow   = row + start;
517       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
518     }
519     VecRestoreArray(w2,&y);
520   }
521   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
522   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
523   return 0;
524 }
525