xref: /petsc/src/mat/matfd/fdmatrix.c (revision 43a90d8472dff9fb76015431d971d6cf7e7ec88a)
1bbf0e169SBarry Smith 
2bbf0e169SBarry Smith #ifndef lint
3*43a90d84SBarry Smith static char vcid[] = "$Id: fdmatrix.c,v 1.2 1996/11/07 15:09:08 bsmith Exp bsmith $";
4bbf0e169SBarry Smith #endif
5bbf0e169SBarry Smith 
6bbf0e169SBarry Smith /*
7639f9d9dSBarry Smith    This is where the abstract matrix operations are defined that are
8639f9d9dSBarry Smith   used for finite difference computations of Jacobians using coloring.
9bbf0e169SBarry Smith */
10bbf0e169SBarry Smith 
11bbf0e169SBarry Smith #include "petsc.h"
12bbf0e169SBarry Smith #include "src/mat/matimpl.h"        /*I "mat.h" I*/
13bbf0e169SBarry Smith #include "src/vec/vecimpl.h"
14bbf0e169SBarry Smith #include "pinclude/pviewer.h"
15bbf0e169SBarry Smith 
16bbf0e169SBarry Smith /*@C
17639f9d9dSBarry Smith    MatFDColoringView - Views a finite difference coloring context.
18bbf0e169SBarry Smith 
19639f9d9dSBarry Smith    Input  Parameter:
20639f9d9dSBarry Smith .   color - the coloring context
21bbf0e169SBarry Smith 
22bbf0e169SBarry Smith 
23639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
24bbf0e169SBarry Smith @*/
25639f9d9dSBarry Smith int MatFDColoringView(MatFDColoring color,Viewer viewer)
26bbf0e169SBarry Smith {
27639f9d9dSBarry Smith   int i,j,format,ierr;
28bbf0e169SBarry Smith 
29639f9d9dSBarry Smith   if (!viewer) viewer = VIEWER_STDOUT_WORLD;
30bbf0e169SBarry Smith 
31bbf0e169SBarry Smith   ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr);
32639f9d9dSBarry Smith 
33639f9d9dSBarry Smith   if (format == VIEWER_FORMAT_ASCII_INFO) {
34639f9d9dSBarry Smith     printf("MatFDColoring Object:\n");
35639f9d9dSBarry Smith     printf("  Error tolerance %g\n",color->error_rel);
36639f9d9dSBarry Smith     printf("  umin %g\n",color->umin);
37639f9d9dSBarry Smith   } else {
38639f9d9dSBarry Smith     printf("MatFDColoring Object:\n");
39639f9d9dSBarry Smith     printf("  Error tolerance %g\n",color->error_rel);
40639f9d9dSBarry Smith     printf("  umin %g\n",color->umin);
41639f9d9dSBarry Smith     printf("Number of colors %d\n",color->ncolors);
42639f9d9dSBarry Smith     for ( i=0; i<color->ncolors; i++ ) {
43639f9d9dSBarry Smith       printf("Information for color %d\n",i);
44639f9d9dSBarry Smith       printf("Number of columns %d\n",color->ncolumns[i]);
45639f9d9dSBarry Smith       for ( j=0; j<color->ncolumns[i]; j++ ) {
46639f9d9dSBarry Smith         printf("  %d\n",color->columns[i][j]);
47639f9d9dSBarry Smith       }
48639f9d9dSBarry Smith       printf("Number of rows %d\n",color->nrows[i]);
49639f9d9dSBarry Smith       for ( j=0; j<color->nrows[i]; j++ ) {
50639f9d9dSBarry Smith         printf("  %d %d \n",color->rows[i][j],color->columnsforrow[i][j]);
51bbf0e169SBarry Smith       }
52bbf0e169SBarry Smith     }
53bbf0e169SBarry Smith   }
54639f9d9dSBarry Smith   return 0;
55639f9d9dSBarry Smith }
56639f9d9dSBarry Smith 
57639f9d9dSBarry Smith /*@
58639f9d9dSBarry Smith    MatFDColoringSetParameters - Sets the parameters for the approximation of
59639f9d9dSBarry Smith    Jacobian using finite differences.
60639f9d9dSBarry Smith 
61639f9d9dSBarry Smith $       J(u)_{:,i} = [J(u+h*dx_{i}) - J(u)]/h where
62639f9d9dSBarry Smith $        h = error_rel*u[i]                    if  u[i] > umin
63639f9d9dSBarry Smith $          = error_rel*umin                    else
64639f9d9dSBarry Smith $
65639f9d9dSBarry Smith $   dx_{i} = (0, ... 1, .... 0)
66639f9d9dSBarry Smith 
67639f9d9dSBarry Smith    Input Parameters:
68639f9d9dSBarry Smith .  coloring - the jacobian coloring context
69639f9d9dSBarry Smith .  error_rel - relative error
70639f9d9dSBarry Smith .  umin - minimum allowable u-value
71639f9d9dSBarry Smith 
72639f9d9dSBarry Smith .keywords: SNES, Jacobian, finite differences, parameters
73639f9d9dSBarry Smith @*/
74639f9d9dSBarry Smith int MatFDColoringSetParameters(MatFDColoring matfd,double error,double umin)
75639f9d9dSBarry Smith {
76639f9d9dSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
77639f9d9dSBarry Smith 
78639f9d9dSBarry Smith   if (error != PETSC_DEFAULT) matfd->error_rel = error;
79639f9d9dSBarry Smith   if (umin != PETSC_DEFAULT)  matfd->umin      = umin;
80639f9d9dSBarry Smith   return 0;
81639f9d9dSBarry Smith }
82639f9d9dSBarry Smith 
83639f9d9dSBarry Smith /*@
84639f9d9dSBarry Smith    MatFDColoringSetFromOptions - Set coloring finite difference parameters from
85639f9d9dSBarry Smith          the options database.
86639f9d9dSBarry Smith 
87639f9d9dSBarry Smith $       J(u)_{:,i} = [J(u+h*dx_{i}) - J(u)]/h where
88639f9d9dSBarry Smith $        h = error_rel*u[i]                    if  u[i] > umin
89639f9d9dSBarry Smith $          = error_rel*.1                      else
90639f9d9dSBarry Smith $
91639f9d9dSBarry Smith $   dx_{i} = (0, ... 1, .... 0)
92639f9d9dSBarry Smith 
93639f9d9dSBarry Smith    Input Parameters:
94639f9d9dSBarry Smith .  coloring - the jacobian coloring context
95639f9d9dSBarry Smith 
96639f9d9dSBarry Smith    Options Database:
97639f9d9dSBarry Smith .  -mat_fd_coloring_error square root of relative error in function
98639f9d9dSBarry Smith .  -mat_fd_coloring_umin  see above
99639f9d9dSBarry Smith 
100639f9d9dSBarry Smith .keywords: SNES, Jacobian, finite differences, parameters
101639f9d9dSBarry Smith @*/
102639f9d9dSBarry Smith int MatFDColoringSetFromOptions(MatFDColoring matfd)
103639f9d9dSBarry Smith {
104639f9d9dSBarry Smith   int    ierr,flag;
105639f9d9dSBarry Smith   double error = PETSC_DEFAULT,umin = PETSC_DEFAULT;
106639f9d9dSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
107639f9d9dSBarry Smith 
108639f9d9dSBarry Smith   ierr = OptionsGetDouble(matfd->prefix,"-mat_fd_coloring_err",&error,&flag);CHKERRQ(ierr);
109639f9d9dSBarry Smith   ierr = OptionsGetDouble(matfd->prefix,"-mat_fd_coloring_umin",&umin,&flag);CHKERRQ(ierr);
110639f9d9dSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help",&flag); CHKERRQ(ierr);
111639f9d9dSBarry Smith   ierr = MatFDColoringSetParameters(matfd,error,umin); CHKERRQ(ierr);
112639f9d9dSBarry Smith   if (flag) {
113639f9d9dSBarry Smith     ierr = MatFDColoringPrintHelp(matfd); CHKERRQ(ierr);
114639f9d9dSBarry Smith   }
115639f9d9dSBarry Smith   return 0;
116639f9d9dSBarry Smith }
117639f9d9dSBarry Smith 
118639f9d9dSBarry Smith /*@
119639f9d9dSBarry Smith     MatFDColoringPrintHelp - Prints help message for matrix finite difference calculations
120639f9d9dSBarry Smith          using coloring.
121639f9d9dSBarry Smith 
122639f9d9dSBarry Smith    Input Parameter:
123639f9d9dSBarry Smith .   fdcoloring - the MatFDColoring context
124639f9d9dSBarry Smith 
125639f9d9dSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringSetFromOptions()
126639f9d9dSBarry Smith @*/
127639f9d9dSBarry Smith int MatFDColoringPrintHelp(MatFDColoring fd)
128639f9d9dSBarry Smith {
129639f9d9dSBarry Smith   PetscValidHeaderSpecific(fd,MAT_FDCOLORING_COOKIE);
130639f9d9dSBarry Smith 
131639f9d9dSBarry Smith   PetscPrintf(fd->comm,"-mat_fd_coloring_err <err> set sqrt rel tol in function. Default 1.e-8\n");
132639f9d9dSBarry Smith   PetscPrintf(fd->comm,"-mat_fd_coloring_umin <umin> see users manual. Default 1.e-8\n");
133bbf0e169SBarry Smith   return 0;
134bbf0e169SBarry Smith }
135bbf0e169SBarry Smith 
136bbf0e169SBarry Smith /*@C
137639f9d9dSBarry Smith    MatFDColoringCreate - Creates a matrix coloring context for finite difference
138639f9d9dSBarry Smith         computation of Jacobians.
139bbf0e169SBarry Smith 
140639f9d9dSBarry Smith    Input Parameters:
141639f9d9dSBarry Smith .    mat - the matrix containing the nonzero structure of the Jacobian
142639f9d9dSBarry Smith .    iscoloring - the coloring of the matrix
143bbf0e169SBarry Smith 
144bbf0e169SBarry Smith    Output Parameter:
145639f9d9dSBarry Smith .   color - the new coloring context
146bbf0e169SBarry Smith 
147639f9d9dSBarry Smith 
148639f9d9dSBarry Smith .seealso: MatFDColoringDestroy()
149bbf0e169SBarry Smith @*/
150639f9d9dSBarry Smith int MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
151bbf0e169SBarry Smith {
152639f9d9dSBarry Smith   MatFDColoring c;
153639f9d9dSBarry Smith   MPI_Comm      comm;
154639f9d9dSBarry Smith   int           ierr,M,N;
155639f9d9dSBarry Smith 
156639f9d9dSBarry Smith   ierr = MatGetSize(mat,&M,&N); CHKERRQ(ierr);
157639f9d9dSBarry Smith   if (M != N) SETERRQ(1,"MatFDColoringCreate:Only for square matrices");
158639f9d9dSBarry Smith 
159639f9d9dSBarry Smith   PetscObjectGetComm((PetscObject)mat,&comm);
160639f9d9dSBarry Smith   PetscHeaderCreate(c,_MatFDColoring,MAT_FDCOLORING_COOKIE,0,comm);
161639f9d9dSBarry Smith   PLogObjectCreate(c);
162639f9d9dSBarry Smith 
163639f9d9dSBarry Smith   if (mat->ops.fdcoloringcreate) {
164639f9d9dSBarry Smith     ierr = (*mat->ops.fdcoloringcreate)(mat,iscoloring,c); CHKERRQ(ierr);
165639f9d9dSBarry Smith   } else {
166639f9d9dSBarry Smith     SETERRQ(1,"MatFDColoringCreate:Code not yet written for this matrix type");
167639f9d9dSBarry Smith   }
168639f9d9dSBarry Smith 
169639f9d9dSBarry Smith   c->error_rel = 1.e-8;
170639f9d9dSBarry Smith   c->umin      = 1.e-5;
171639f9d9dSBarry Smith 
172639f9d9dSBarry Smith   *color = c;
173639f9d9dSBarry Smith 
174bbf0e169SBarry Smith   return 0;
175bbf0e169SBarry Smith }
176bbf0e169SBarry Smith 
177bbf0e169SBarry Smith /*@C
178639f9d9dSBarry Smith     MatFDColoringDestroy - Destroys a matrix coloring context that was created
179639f9d9dSBarry Smith          via MatFDColoringCreate().
180bbf0e169SBarry Smith 
181639f9d9dSBarry Smith     Input Paramter:
182639f9d9dSBarry Smith .   c - coloring context
183bbf0e169SBarry Smith 
184639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
185bbf0e169SBarry Smith @*/
186639f9d9dSBarry Smith int MatFDColoringDestroy(MatFDColoring c)
187bbf0e169SBarry Smith {
188639f9d9dSBarry Smith   int i,ierr,flag;
189bbf0e169SBarry Smith 
190639f9d9dSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-matfdcoloring_view",&flag);
191639f9d9dSBarry Smith   if (flag) {
192bbf0e169SBarry Smith     Viewer viewer;
193639f9d9dSBarry Smith     ierr = ViewerFileOpenASCII(c->comm,"stdout",&viewer);CHKERRQ(ierr);
194639f9d9dSBarry Smith     ierr = MatFDColoringView(c,viewer);CHKERRQ(ierr);
195bbf0e169SBarry Smith     ierr = ViewerDestroy(viewer); CHKERRQ(ierr);
196bbf0e169SBarry Smith   }
197639f9d9dSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-matfdcoloring_view_info",&flag);
198639f9d9dSBarry Smith   if (flag) {
199bbf0e169SBarry Smith     Viewer viewer;
200639f9d9dSBarry Smith     ierr = ViewerFileOpenASCII(c->comm,"stdout",&viewer);CHKERRQ(ierr);
201639f9d9dSBarry Smith     ierr = ViewerPushFormat(viewer,VIEWER_FORMAT_ASCII_INFO, PETSC_NULL); CHKERRQ(ierr);
202639f9d9dSBarry Smith     ierr = MatFDColoringView(c,viewer);CHKERRQ(ierr);
203639f9d9dSBarry Smith     ierr = ViewerPopFormat(viewer);
204bbf0e169SBarry Smith     ierr = ViewerDestroy(viewer); CHKERRQ(ierr);
205bbf0e169SBarry Smith   }
206639f9d9dSBarry Smith 
207639f9d9dSBarry Smith   for ( i=0; i<c->ncolors; i++ ) {
208639f9d9dSBarry Smith     if (c->columns[i])       PetscFree(c->columns[i]);
209639f9d9dSBarry Smith     if (c->rows[i])          PetscFree(c->rows[i]);
210639f9d9dSBarry Smith     if (c->columnsforrow[i]) PetscFree(c->columnsforrow[i]);
211bbf0e169SBarry Smith   }
212639f9d9dSBarry Smith   PetscFree(c->ncolumns);
213639f9d9dSBarry Smith   PetscFree(c->columns);
214639f9d9dSBarry Smith   PetscFree(c->nrows);
215639f9d9dSBarry Smith   PetscFree(c->rows);
216639f9d9dSBarry Smith   PetscFree(c->columnsforrow);
217639f9d9dSBarry Smith   PetscFree(c->scale);
218639f9d9dSBarry Smith   PLogObjectDestroy(c);
219639f9d9dSBarry Smith   PetscHeaderDestroy(c);
220bbf0e169SBarry Smith   return 0;
221bbf0e169SBarry Smith }
222*43a90d84SBarry Smith 
223*43a90d84SBarry Smith /*@
224*43a90d84SBarry Smith      MatFDColoringApply - Given a matrix for which a MatFDColoring has been created,
225*43a90d84SBarry Smith          computes the Jacobian for a function via finite differences.
226*43a90d84SBarry Smith 
227*43a90d84SBarry Smith   Input Parameters:
228*43a90d84SBarry Smith .   mat - location to store Jacobian
229*43a90d84SBarry Smith .   coloring - coloring context created with MatFDColoringCreate()
230*43a90d84SBarry Smith .   x1 - location at which Jacobian is to be computed
231*43a90d84SBarry Smith .   w1,w2,w3 - three work vectors
232*43a90d84SBarry Smith .   f - function for which Jacobian is required
233*43a90d84SBarry Smith .   fctx - optional context required by function
234*43a90d84SBarry Smith 
235*43a90d84SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
236*43a90d84SBarry Smith 
237*43a90d84SBarry Smith .keywords: coloring, Jacobian, finite differences
238*43a90d84SBarry Smith @*/
239*43a90d84SBarry Smith int MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,Vec w1,Vec w2,Vec w3,
240*43a90d84SBarry Smith                        int (*f)(void *,Vec,Vec,void*),void *sctx,void *fctx)
241*43a90d84SBarry Smith {
242*43a90d84SBarry Smith   int           k, ierr,N,start,end,l,row,col,srow;
243*43a90d84SBarry Smith   Scalar        dx, mone = -1.0,*y,*scale = coloring->scale,*xx,*wscale = coloring->wscale;
244*43a90d84SBarry Smith   double        epsilon = coloring->error_rel, umin = coloring->umin;
245*43a90d84SBarry Smith   MPI_Comm      comm = coloring->comm;
246*43a90d84SBarry Smith 
247*43a90d84SBarry Smith   ierr = MatZeroEntries(J); CHKERRQ(ierr);
248*43a90d84SBarry Smith 
249*43a90d84SBarry Smith   ierr = VecGetOwnershipRange(x1,&start,&end); CHKERRQ(ierr);
250*43a90d84SBarry Smith   ierr = VecGetSize(x1,&N); CHKERRQ(ierr);
251*43a90d84SBarry Smith   ierr = VecGetArray(x1,&xx); CHKERRQ(ierr);
252*43a90d84SBarry Smith   ierr = (*f)(sctx,x1,w1,fctx); CHKERRQ(ierr);
253*43a90d84SBarry Smith 
254*43a90d84SBarry Smith   PetscMemzero(wscale,N*sizeof(Scalar));
255*43a90d84SBarry Smith   /*
256*43a90d84SBarry Smith       Loop over each color
257*43a90d84SBarry Smith   */
258*43a90d84SBarry Smith 
259*43a90d84SBarry Smith   for (k=0; k<coloring->ncolors; k++) {
260*43a90d84SBarry Smith     ierr = VecCopy(x1,w3); CHKERRQ(ierr);
261*43a90d84SBarry Smith     /*
262*43a90d84SBarry Smith        Loop over each column associated with color adding the
263*43a90d84SBarry Smith        perturbation to the vector w3.
264*43a90d84SBarry Smith     */
265*43a90d84SBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
266*43a90d84SBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
267*43a90d84SBarry Smith       dx  = xx[col-start];
268*43a90d84SBarry Smith #if !defined(PETSC_COMPLEX)
269*43a90d84SBarry Smith       if (dx < umin && dx >= 0.0) dx = .1;
270*43a90d84SBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -.1;
271*43a90d84SBarry Smith #else
272*43a90d84SBarry Smith       if (abs(dx) < umin && real(dx) >= 0.0) dx = .1;
273*43a90d84SBarry Smith       else if (real(dx) < 0.0 && abs(dx) < umin) dx = -.1;
274*43a90d84SBarry Smith #endif
275*43a90d84SBarry Smith       dx          *= epsilon;
276*43a90d84SBarry Smith       wscale[col] = 1.0/dx;
277*43a90d84SBarry Smith       VecSetValues(w3,1,&col,&dx,ADD_VALUES);
278*43a90d84SBarry Smith     }
279*43a90d84SBarry Smith     VecRestoreArray(x1,&xx);
280*43a90d84SBarry Smith     /*
281*43a90d84SBarry Smith        Evaluate function at x1 + dx (here dx is a vector, of perturbations)
282*43a90d84SBarry Smith     */
283*43a90d84SBarry Smith     ierr = (*f)(sctx,w3,w2,fctx); CHKERRQ(ierr);
284*43a90d84SBarry Smith     ierr = VecAXPY(&mone,w1,w2); CHKERRQ(ierr);
285*43a90d84SBarry Smith     /* Communicate scale to all processors */
286*43a90d84SBarry Smith #if !defined(PETSC_COMPLEX)
287*43a90d84SBarry Smith     MPI_Allreduce(wscale,scale,N,MPI_DOUBLE,MPI_SUM,comm);
288*43a90d84SBarry Smith #else
289*43a90d84SBarry Smith     MPI_Allreduce(wscale,scale,2*N,MPI_DOUBLE,MPI_SUM,comm);
290*43a90d84SBarry Smith #endif
291*43a90d84SBarry Smith     /*
292*43a90d84SBarry Smith        Loop over rows of vector putting results into Jacobian matrix
293*43a90d84SBarry Smith     */
294*43a90d84SBarry Smith     VecGetArray(w2,&y);
295*43a90d84SBarry Smith     for (l=0; l<coloring->nrows[k]; l++) {
296*43a90d84SBarry Smith       row    = coloring->rows[k][l];
297*43a90d84SBarry Smith       col    = coloring->columnsforrow[k][l];
298*43a90d84SBarry Smith       y[row] *= scale[col];
299*43a90d84SBarry Smith       srow   = row + start;
300*43a90d84SBarry Smith       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
301*43a90d84SBarry Smith     }
302*43a90d84SBarry Smith     VecRestoreArray(w2,&y);
303*43a90d84SBarry Smith   }
304*43a90d84SBarry Smith   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
305*43a90d84SBarry Smith   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
306*43a90d84SBarry Smith   return 0;
307*43a90d84SBarry Smith }
308