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