xref: /petsc/src/mat/matfd/fdmatrix.c (revision bcff55d9c22a11117c657f321125d3821e9fb4a8)
1 #ifdef PETSC_RCS_HEADER
2 static char vcid[] = "$Id: fdmatrix.c,v 1.10 1997/05/23 18:38:43 balay Exp balay $";
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" /* ADIC Ignore */
17 /*@C
18    MatFDColoringView - Views a finite difference coloring context.
19 
20    Input  Parameter:
21 .   color - the coloring context
22 
23 
24 .seealso: MatFDColoringCreate()
25 @*/
26 int MatFDColoringView(MatFDColoring color,Viewer viewer)
27 {
28   int i,j,format,ierr;
29 
30   if (!viewer) viewer = VIEWER_STDOUT_WORLD;
31 
32   ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr);
33 
34   if (format == VIEWER_FORMAT_ASCII_INFO) {
35     printf("MatFDColoring Object:\n");
36     printf("  Error tolerance %g\n",color->error_rel);
37     printf("  umin %g\n",color->umin);
38   } else {
39     printf("MatFDColoring Object:\n");
40     printf("  Error tolerance %g\n",color->error_rel);
41     printf("  umin %g\n",color->umin);
42     printf("Number of colors %d\n",color->ncolors);
43     for ( i=0; i<color->ncolors; i++ ) {
44       printf("Information for color %d\n",i);
45       printf("Number of columns %d\n",color->ncolumns[i]);
46       for ( j=0; j<color->ncolumns[i]; j++ ) {
47         printf("  %d\n",color->columns[i][j]);
48       }
49       printf("Number of rows %d\n",color->nrows[i]);
50       for ( j=0; j<color->nrows[i]; j++ ) {
51         printf("  %d %d \n",color->rows[i][j],color->columnsforrow[i][j]);
52       }
53     }
54   }
55   return 0;
56 }
57 
58 #undef __FUNC__
59 #define __FUNC__ "MatFDColoringSetParameters"
60 /*@
61    MatFDColoringSetParameters - Sets the parameters for the approximation of
62    Jacobian using finite differences.
63 
64 $       J(u)_{:,i} = [J(u+h*dx_{i}) - J(u)]/h where
65 $        h = error_rel*u[i]                    if  u[i] > umin
66 $          = error_rel*umin                    else
67 $
68 $   dx_{i} = (0, ... 1, .... 0)
69 
70    Input Parameters:
71 .  coloring - the jacobian coloring context
72 .  error_rel - relative error
73 .  umin - minimum allowable u-value
74 
75 .keywords: SNES, Jacobian, finite differences, parameters
76 @*/
77 int MatFDColoringSetParameters(MatFDColoring matfd,double error,double umin)
78 {
79   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
80 
81   if (error != PETSC_DEFAULT) matfd->error_rel = error;
82   if (umin != PETSC_DEFAULT)  matfd->umin      = umin;
83   return 0;
84 }
85 
86 #undef __FUNC__
87 #define __FUNC__ "MatFDColoringSetFromOptions" /* ADIC Ignore */
88 /*@
89    MatFDColoringSetFromOptions - Set coloring finite difference parameters from
90          the options database.
91 
92 $       J(u)_{:,i} = [J(u+h*dx_{i}) - J(u)]/h where
93 $        h = error_rel*u[i]                    if  u[i] > umin
94 $          = error_rel*.1                      else
95 $
96 $   dx_{i} = (0, ... 1, .... 0)
97 
98    Input Parameters:
99 .  coloring - the jacobian coloring context
100 
101    Options Database:
102 .  -mat_fd_coloring_error square root of relative error in function
103 .  -mat_fd_coloring_umin  see above
104 
105 .keywords: SNES, Jacobian, finite differences, parameters
106 @*/
107 int MatFDColoringSetFromOptions(MatFDColoring matfd)
108 {
109   int    ierr,flag;
110   double error = PETSC_DEFAULT,umin = PETSC_DEFAULT;
111   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
112 
113   ierr = OptionsGetDouble(matfd->prefix,"-mat_fd_coloring_err",&error,&flag);CHKERRQ(ierr);
114   ierr = OptionsGetDouble(matfd->prefix,"-mat_fd_coloring_umin",&umin,&flag);CHKERRQ(ierr);
115   ierr = OptionsHasName(PETSC_NULL,"-help",&flag); CHKERRQ(ierr);
116   ierr = MatFDColoringSetParameters(matfd,error,umin); CHKERRQ(ierr);
117   if (flag) {
118     ierr = MatFDColoringPrintHelp(matfd); CHKERRQ(ierr);
119   }
120   return 0;
121 }
122 
123 #undef __FUNC__
124 #define __FUNC__ "MatFDColoringPrintHelp" /* ADIC Ignore */
125 /*@
126     MatFDColoringPrintHelp - Prints help message for matrix finite difference calculations
127          using coloring.
128 
129    Input Parameter:
130 .   fdcoloring - the MatFDColoring context
131 
132 .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringSetFromOptions()
133 @*/
134 int MatFDColoringPrintHelp(MatFDColoring fd)
135 {
136   PetscValidHeaderSpecific(fd,MAT_FDCOLORING_COOKIE);
137 
138   PetscPrintf(fd->comm,"-mat_fd_coloring_err <err>: set sqrt rel tol in function, defaults to 1.e-8\n");
139   PetscPrintf(fd->comm,"-mat_fd_coloring_umin <umin>: see users manual, defaults to 1.e-8\n");
140   return 0;
141 }
142 
143 #undef __FUNC__
144 #define __FUNC__ "MatFDColoringCreate"
145 /*@C
146    MatFDColoringCreate - Creates a matrix coloring context for finite difference
147         computation of Jacobians.
148 
149    Input Parameters:
150 .    mat - the matrix containing the nonzero structure of the Jacobian
151 .    iscoloring - the coloring of the matrix
152 
153    Output Parameter:
154 .   color - the new coloring context
155 
156 
157 .seealso: MatFDColoringDestroy()
158 @*/
159 int MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
160 {
161   MatFDColoring c;
162   MPI_Comm      comm;
163   int           ierr,M,N;
164 
165   ierr = MatGetSize(mat,&M,&N); CHKERRQ(ierr);
166   if (M != N) SETERRQ(PETSC_ERR_SUP,0,"Only for square matrices");
167 
168   PetscObjectGetComm((PetscObject)mat,&comm);
169   PetscHeaderCreate(c,_p_MatFDColoring,MAT_FDCOLORING_COOKIE,0,comm);
170   PLogObjectCreate(c);
171 
172   if (mat->ops.fdcoloringcreate) {
173     ierr = (*mat->ops.fdcoloringcreate)(mat,iscoloring,c); CHKERRQ(ierr);
174   } else {
175     SETERRQ(PETSC_ERR_SUP,0,"Code not yet written for this matrix type");
176   }
177 
178   c->error_rel = 1.e-8;
179   c->umin      = 1.e-5;
180 
181   *color = c;
182 
183   return 0;
184 }
185 
186 #undef __FUNC__
187 #define __FUNC__ "MatFDColoringDestroy" /* ADIC Ignore */
188 /*@C
189     MatFDColoringDestroy - Destroys a matrix coloring context that was created
190          via MatFDColoringCreate().
191 
192     Input Paramter:
193 .   c - coloring context
194 
195 .seealso: MatFDColoringCreate()
196 @*/
197 int MatFDColoringDestroy(MatFDColoring c)
198 {
199   int i,ierr,flag;
200 
201   ierr = OptionsHasName(PETSC_NULL,"-matfdcoloring_view",&flag);
202   if (flag) {
203     Viewer viewer;
204     ierr = ViewerFileOpenASCII(c->comm,"stdout",&viewer);CHKERRQ(ierr);
205     ierr = MatFDColoringView(c,viewer);CHKERRQ(ierr);
206     ierr = ViewerDestroy(viewer); CHKERRQ(ierr);
207   }
208   ierr = OptionsHasName(PETSC_NULL,"-matfdcoloring_view_info",&flag);
209   if (flag) {
210     Viewer viewer;
211     ierr = ViewerFileOpenASCII(c->comm,"stdout",&viewer);CHKERRQ(ierr);
212     ierr = ViewerPushFormat(viewer,VIEWER_FORMAT_ASCII_INFO, PETSC_NULL); CHKERRQ(ierr);
213     ierr = MatFDColoringView(c,viewer);CHKERRQ(ierr);
214     ierr = ViewerPopFormat(viewer);
215     ierr = ViewerDestroy(viewer); CHKERRQ(ierr);
216   }
217 
218   for ( i=0; i<c->ncolors; i++ ) {
219     if (c->columns[i])       PetscFree(c->columns[i]);
220     if (c->rows[i])          PetscFree(c->rows[i]);
221     if (c->columnsforrow[i]) PetscFree(c->columnsforrow[i]);
222   }
223   PetscFree(c->ncolumns);
224   PetscFree(c->columns);
225   PetscFree(c->nrows);
226   PetscFree(c->rows);
227   PetscFree(c->columnsforrow);
228   PetscFree(c->scale);
229   PLogObjectDestroy(c);
230   PetscHeaderDestroy(c);
231   return 0;
232 }
233 
234 #undef __FUNC__
235 #define __FUNC__ "MatFDColoringApply"
236 /*@
237      MatFDColoringApply - Given a matrix for which a MatFDColoring has been created,
238          computes the Jacobian for a function via finite differences.
239 
240   Input Parameters:
241 .   mat - location to store Jacobian
242 .   coloring - coloring context created with MatFDColoringCreate()
243 .   x1 - location at which Jacobian is to be computed
244 .   w1,w2,w3 - three work vectors
245 .   f - function for which Jacobian is required
246 .   fctx - optional context required by function
247 
248 .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
249 
250 .keywords: coloring, Jacobian, finite differences
251 @*/
252 int MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,Vec w1,Vec w2,Vec w3,
253                        int (*f)(void *,Vec,Vec,void*),void *sctx,void *fctx)
254 {
255   int           k, ierr,N,start,end,l,row,col,srow;
256   Scalar        dx, mone = -1.0,*y,*scale = coloring->scale,*xx,*wscale = coloring->wscale;
257   double        epsilon = coloring->error_rel, umin = coloring->umin;
258   MPI_Comm      comm = coloring->comm;
259 
260   ierr = MatZeroEntries(J); CHKERRQ(ierr);
261 
262   ierr = VecGetOwnershipRange(x1,&start,&end); CHKERRQ(ierr);
263   ierr = VecGetSize(x1,&N); CHKERRQ(ierr);
264   ierr = VecGetArray(x1,&xx); CHKERRQ(ierr);
265   ierr = (*f)(sctx,x1,w1,fctx); CHKERRQ(ierr);
266 
267   PetscMemzero(wscale,N*sizeof(Scalar));
268   /*
269       Loop over each color
270   */
271 
272   for (k=0; k<coloring->ncolors; k++) {
273     ierr = VecCopy(x1,w3); CHKERRQ(ierr);
274     /*
275        Loop over each column associated with color adding the
276        perturbation to the vector w3.
277     */
278     for (l=0; l<coloring->ncolumns[k]; l++) {
279       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
280       dx  = xx[col-start];
281 #if !defined(PETSC_COMPLEX)
282       if (dx < umin && dx >= 0.0) dx = .1;
283       else if (dx < 0.0 && dx > -umin) dx = -.1;
284 #else
285       if (abs(dx) < umin && real(dx) >= 0.0) dx = .1;
286       else if (real(dx) < 0.0 && abs(dx) < umin) dx = -.1;
287 #endif
288       dx          *= epsilon;
289       wscale[col] = 1.0/dx;
290       VecSetValues(w3,1,&col,&dx,ADD_VALUES);
291     }
292     VecRestoreArray(x1,&xx);
293     /*
294        Evaluate function at x1 + dx (here dx is a vector, of perturbations)
295     */
296     ierr = (*f)(sctx,w3,w2,fctx); CHKERRQ(ierr);
297     ierr = VecAXPY(&mone,w1,w2); CHKERRQ(ierr);
298     /* Communicate scale to all processors */
299 #if !defined(PETSC_COMPLEX)
300     MPI_Allreduce(wscale,scale,N,MPI_DOUBLE,MPI_SUM,comm);
301 #else
302     MPI_Allreduce(wscale,scale,2*N,MPI_DOUBLE,MPI_SUM,comm);
303 #endif
304     /*
305        Loop over rows of vector putting results into Jacobian matrix
306     */
307     VecGetArray(w2,&y);
308     for (l=0; l<coloring->nrows[k]; l++) {
309       row    = coloring->rows[k][l];
310       col    = coloring->columnsforrow[k][l];
311       y[row] *= scale[col];
312       srow   = row + start;
313       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
314     }
315     VecRestoreArray(w2,&y);
316   }
317   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
318   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
319   return 0;
320 }
321