xref: /petsc/src/mat/matfd/fdmatrix.c (revision 4e269d778a4980f937dbf53e8c867b2c3b45d02b)
1be1d678aSKris Buschelman 
2bbf0e169SBarry Smith /*
3639f9d9dSBarry Smith    This is where the abstract matrix operations are defined that are
4639f9d9dSBarry Smith   used for finite difference computations of Jacobians using coloring.
5bbf0e169SBarry Smith */
6bbf0e169SBarry Smith 
7b45d2f2cSJed Brown #include <petsc-private/matimpl.h>        /*I "petscmat.h" I*/
8bbf0e169SBarry Smith 
94a2ae208SSatish Balay #undef __FUNCT__
103a7fca6bSBarry Smith #define __FUNCT__ "MatFDColoringSetF"
117087cfbeSBarry Smith PetscErrorCode  MatFDColoringSetF(MatFDColoring fd,Vec F)
123a7fca6bSBarry Smith {
13*4e269d77SPeter Brune   PetscErrorCode ierr;
143a7fca6bSBarry Smith   PetscFunctionBegin;
15*4e269d77SPeter Brune   if (F) {
16*4e269d77SPeter Brune     ierr = VecCopy(F,fd->w1);CHKERRQ(ierr);
17*4e269d77SPeter Brune     fd->fset = PETSC_TRUE;
18*4e269d77SPeter Brune   } else {
19*4e269d77SPeter Brune     fd->fset = PETSC_FALSE;
20*4e269d77SPeter Brune   }
213a7fca6bSBarry Smith   PetscFunctionReturn(0);
223a7fca6bSBarry Smith }
233a7fca6bSBarry Smith 
243a7fca6bSBarry Smith #undef __FUNCT__
254a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Draw_Zoom"
266849ba73SBarry Smith static PetscErrorCode MatFDColoringView_Draw_Zoom(PetscDraw draw,void *Aa)
279194eea9SBarry Smith {
289194eea9SBarry Smith   MatFDColoring  fd = (MatFDColoring)Aa;
29dfbe8321SBarry Smith   PetscErrorCode ierr;
3038baddfdSBarry Smith   PetscInt       i,j;
319194eea9SBarry Smith   PetscReal      x,y;
329194eea9SBarry Smith 
339194eea9SBarry Smith   PetscFunctionBegin;
349194eea9SBarry Smith 
359194eea9SBarry Smith   /* loop over colors  */
369194eea9SBarry Smith   for (i=0; i<fd->ncolors; i++) {
379194eea9SBarry Smith     for (j=0; j<fd->nrows[i]; j++) {
389194eea9SBarry Smith       y = fd->M - fd->rows[i][j] - fd->rstart;
399194eea9SBarry Smith       x = fd->columnsforrow[i][j];
40b0a32e0cSBarry Smith       ierr = PetscDrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);CHKERRQ(ierr);
419194eea9SBarry Smith     }
429194eea9SBarry Smith   }
439194eea9SBarry Smith   PetscFunctionReturn(0);
449194eea9SBarry Smith }
459194eea9SBarry Smith 
464a2ae208SSatish Balay #undef __FUNCT__
474a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Draw"
486849ba73SBarry Smith static PetscErrorCode MatFDColoringView_Draw(MatFDColoring fd,PetscViewer viewer)
49005c665bSBarry Smith {
50dfbe8321SBarry Smith   PetscErrorCode ierr;
51ace3abfcSBarry Smith   PetscBool      isnull;
52b0a32e0cSBarry Smith   PetscDraw      draw;
5354d96fa3SBarry Smith   PetscReal      xr,yr,xl,yl,h,w;
54005c665bSBarry Smith 
553a40ed3dSBarry Smith   PetscFunctionBegin;
56b0a32e0cSBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
57b0a32e0cSBarry Smith   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
589194eea9SBarry Smith 
599194eea9SBarry Smith   ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
60005c665bSBarry Smith 
61005c665bSBarry Smith   xr  = fd->N; yr = fd->M; h = yr/10.0; w = xr/10.0;
62005c665bSBarry Smith   xr += w;     yr += h;    xl = -w;     yl = -h;
63b0a32e0cSBarry Smith   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
64b0a32e0cSBarry Smith   ierr = PetscDrawZoom(draw,MatFDColoringView_Draw_Zoom,fd);CHKERRQ(ierr);
659194eea9SBarry Smith   ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
663a40ed3dSBarry Smith   PetscFunctionReturn(0);
67005c665bSBarry Smith }
68005c665bSBarry Smith 
694a2ae208SSatish Balay #undef __FUNCT__
704a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView"
71bbf0e169SBarry Smith /*@C
72639f9d9dSBarry Smith    MatFDColoringView - Views a finite difference coloring context.
73bbf0e169SBarry Smith 
744c49b128SBarry Smith    Collective on MatFDColoring
75fee21e36SBarry Smith 
76ef5ee4d1SLois Curfman McInnes    Input  Parameters:
77ef5ee4d1SLois Curfman McInnes +  c - the coloring context
78ef5ee4d1SLois Curfman McInnes -  viewer - visualization context
79ef5ee4d1SLois Curfman McInnes 
8015091d37SBarry Smith    Level: intermediate
8115091d37SBarry Smith 
82b4fc646aSLois Curfman McInnes    Notes:
83b4fc646aSLois Curfman McInnes    The available visualization contexts include
84b0a32e0cSBarry Smith +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
85b0a32e0cSBarry Smith .     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
86ef5ee4d1SLois Curfman McInnes         output where only the first processor opens
87ef5ee4d1SLois Curfman McInnes         the file.  All other processors send their
88ef5ee4d1SLois Curfman McInnes         data to the first processor to print.
89b0a32e0cSBarry Smith -     PETSC_VIEWER_DRAW_WORLD - graphical display of nonzero structure
90bbf0e169SBarry Smith 
919194eea9SBarry Smith    Notes:
929194eea9SBarry Smith      Since PETSc uses only a small number of basic colors (currently 33), if the coloring
93b0a32e0cSBarry Smith    involves more than 33 then some seemingly identical colors are displayed making it look
949194eea9SBarry Smith    like an illegal coloring. This is just a graphical artifact.
959194eea9SBarry Smith 
96639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
97005c665bSBarry Smith 
98b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, view
99bbf0e169SBarry Smith @*/
1007087cfbeSBarry Smith PetscErrorCode  MatFDColoringView(MatFDColoring c,PetscViewer viewer)
101bbf0e169SBarry Smith {
1026849ba73SBarry Smith   PetscErrorCode    ierr;
10338baddfdSBarry Smith   PetscInt          i,j;
104ace3abfcSBarry Smith   PetscBool         isdraw,iascii;
105f3ef73ceSBarry Smith   PetscViewerFormat format;
106bbf0e169SBarry Smith 
1073a40ed3dSBarry Smith   PetscFunctionBegin;
1080700a824SBarry Smith   PetscValidHeaderSpecific(c,MAT_FDCOLORING_CLASSID,1);
1093050cee2SBarry Smith   if (!viewer) {
1107adad957SLisandro Dalcin     ierr = PetscViewerASCIIGetStdout(((PetscObject)c)->comm,&viewer);CHKERRQ(ierr);
1113050cee2SBarry Smith   }
1120700a824SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
113c9780b6fSBarry Smith   PetscCheckSameComm(c,1,viewer,2);
114bbf0e169SBarry Smith 
115251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
116251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1170f5bd95cSBarry Smith   if (isdraw) {
118b4fc646aSLois Curfman McInnes     ierr = MatFDColoringView_Draw(c,viewer);CHKERRQ(ierr);
11932077d6dSBarry Smith   } else if (iascii) {
120317d6ea6SBarry Smith     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)c,viewer,"MatFDColoring Object");CHKERRQ(ierr);
121a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Error tolerance=%G\n",c->error_rel);CHKERRQ(ierr);
122a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Umin=%G\n",c->umin);CHKERRQ(ierr);
12377431f27SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Number of colors=%D\n",c->ncolors);CHKERRQ(ierr);
124ae09f205SBarry Smith 
125b0a32e0cSBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
126fb9695e5SSatish Balay     if (format != PETSC_VIEWER_ASCII_INFO) {
127b4fc646aSLois Curfman McInnes       for (i=0; i<c->ncolors; i++) {
12877431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"  Information for color %D\n",i);CHKERRQ(ierr);
12977431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    Number of columns %D\n",c->ncolumns[i]);CHKERRQ(ierr);
130b4fc646aSLois Curfman McInnes         for (j=0; j<c->ncolumns[i]; j++) {
13177431f27SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"      %D\n",c->columns[i][j]);CHKERRQ(ierr);
132639f9d9dSBarry Smith         }
13377431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    Number of rows %D\n",c->nrows[i]);CHKERRQ(ierr);
134b4fc646aSLois Curfman McInnes         for (j=0; j<c->nrows[i]; j++) {
13577431f27SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"      %D %D \n",c->rows[i][j],c->columnsforrow[i][j]);CHKERRQ(ierr);
136b4fc646aSLois Curfman McInnes         }
137bbf0e169SBarry Smith       }
138bbf0e169SBarry Smith     }
139b0a32e0cSBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1405cd90555SBarry Smith   } else {
141e32f2f54SBarry Smith     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for MatFDColoring",((PetscObject)viewer)->type_name);
142bbf0e169SBarry Smith   }
1433a40ed3dSBarry Smith   PetscFunctionReturn(0);
144639f9d9dSBarry Smith }
145639f9d9dSBarry Smith 
1464a2ae208SSatish Balay #undef __FUNCT__
1474a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetParameters"
148639f9d9dSBarry Smith /*@
149b4fc646aSLois Curfman McInnes    MatFDColoringSetParameters - Sets the parameters for the sparse approximation of
150b4fc646aSLois Curfman McInnes    a Jacobian matrix using finite differences.
151639f9d9dSBarry Smith 
1523f9fe445SBarry Smith    Logically Collective on MatFDColoring
153ef5ee4d1SLois Curfman McInnes 
154ef5ee4d1SLois Curfman McInnes    The Jacobian is estimated with the differencing approximation
155ef5ee4d1SLois Curfman McInnes .vb
15665f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
157f23b5b22SLois Curfman McInnes        h = error_rel*u[i]                 if  abs(u[i]) > umin
158f23b5b22SLois Curfman McInnes          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
159ef5ee4d1SLois Curfman McInnes        dx_{i} = (0, ... 1, .... 0)
160ef5ee4d1SLois Curfman McInnes .ve
161639f9d9dSBarry Smith 
162639f9d9dSBarry Smith    Input Parameters:
163ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
164639f9d9dSBarry Smith .  error_rel - relative error
165f23b5b22SLois Curfman McInnes -  umin - minimum allowable u-value magnitude
166fee21e36SBarry Smith 
16715091d37SBarry Smith    Level: advanced
16815091d37SBarry Smith 
169b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, parameters
170b4fc646aSLois Curfman McInnes 
17195b89fc3SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringSetFromOptions()
17295b89fc3SBarry Smith 
173639f9d9dSBarry Smith @*/
1747087cfbeSBarry Smith PetscErrorCode  MatFDColoringSetParameters(MatFDColoring matfd,PetscReal error,PetscReal umin)
175639f9d9dSBarry Smith {
1763a40ed3dSBarry Smith   PetscFunctionBegin;
1770700a824SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
178c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(matfd,error,2);
179c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(matfd,umin,3);
180639f9d9dSBarry Smith   if (error != PETSC_DEFAULT) matfd->error_rel = error;
181639f9d9dSBarry Smith   if (umin != PETSC_DEFAULT)  matfd->umin      = umin;
1823a40ed3dSBarry Smith   PetscFunctionReturn(0);
183639f9d9dSBarry Smith }
184639f9d9dSBarry Smith 
185005c665bSBarry Smith 
186ff0cfa39SBarry Smith 
1874a2ae208SSatish Balay #undef __FUNCT__
1884a9d489dSBarry Smith #define __FUNCT__ "MatFDColoringGetFunction"
1894a9d489dSBarry Smith /*@C
1904a9d489dSBarry Smith    MatFDColoringGetFunction - Gets the function to use for computing the Jacobian.
1914a9d489dSBarry Smith 
1923f9fe445SBarry Smith    Not Collective
1934a9d489dSBarry Smith 
1944a9d489dSBarry Smith    Input Parameters:
1954a9d489dSBarry Smith .  coloring - the coloring context
1964a9d489dSBarry Smith 
1974a9d489dSBarry Smith    Output Parameters:
1984a9d489dSBarry Smith +  f - the function
1994a9d489dSBarry Smith -  fctx - the optional user-defined function context
2004a9d489dSBarry Smith 
2014a9d489dSBarry Smith    Level: intermediate
2024a9d489dSBarry Smith 
2034a9d489dSBarry Smith .keywords: Mat, Jacobian, finite differences, set, function
20495b89fc3SBarry Smith 
20595b89fc3SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringSetFunction(), MatFDColoringSetFromOptions()
20695b89fc3SBarry Smith 
2074a9d489dSBarry Smith @*/
2087087cfbeSBarry Smith PetscErrorCode  MatFDColoringGetFunction(MatFDColoring matfd,PetscErrorCode (**f)(void),void **fctx)
2094a9d489dSBarry Smith {
2104a9d489dSBarry Smith   PetscFunctionBegin;
2110700a824SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
2124a9d489dSBarry Smith   if (f) *f = matfd->f;
2134a9d489dSBarry Smith   if (fctx) *fctx = matfd->fctx;
2144a9d489dSBarry Smith   PetscFunctionReturn(0);
2154a9d489dSBarry Smith }
2164a9d489dSBarry Smith 
2174a9d489dSBarry Smith #undef __FUNCT__
2184a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFunction"
219d64ed03dSBarry Smith /*@C
220005c665bSBarry Smith    MatFDColoringSetFunction - Sets the function to use for computing the Jacobian.
221005c665bSBarry Smith 
2223f9fe445SBarry Smith    Logically Collective on MatFDColoring
223fee21e36SBarry Smith 
224ef5ee4d1SLois Curfman McInnes    Input Parameters:
225ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
226ef5ee4d1SLois Curfman McInnes .  f - the function
227ef5ee4d1SLois Curfman McInnes -  fctx - the optional user-defined function context
228ef5ee4d1SLois Curfman McInnes 
2297850c7c0SBarry Smith    Calling sequence of (*f) function:
2307850c7c0SBarry Smith     For SNES:    PetscErrorCode (*f)(SNES,Vec,Vec,void*)
231ab637aeaSJed Brown     If not using SNES: PetscErrorCode (*f)(void *dummy,Vec,Vec,void*) and dummy is ignored
23215091d37SBarry Smith 
2337850c7c0SBarry Smith    Level: advanced
2347850c7c0SBarry Smith 
235ab637aeaSJed Brown    Notes: This function is usually used automatically by SNES (when one uses SNESSetJacobian() with the argument
236ab637aeaSJed Brown      SNESDefaultComputeJacobianColor()) and only needs to be used by someone computing a matrix via coloring directly by
237ab637aeaSJed Brown      calling MatFDColoringApply()
2387850c7c0SBarry Smith 
2397850c7c0SBarry Smith    Fortran Notes:
2407850c7c0SBarry Smith     In Fortran you must call MatFDColoringSetFunction() for a coloring object to
241ab637aeaSJed Brown   be used without SNES or within the SNES solvers.
242f881d145SBarry Smith 
243b4fc646aSLois Curfman McInnes .keywords: Mat, Jacobian, finite differences, set, function
24495b89fc3SBarry Smith 
24595b89fc3SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringGetFunction(), MatFDColoringSetFromOptions()
24695b89fc3SBarry Smith 
247005c665bSBarry Smith @*/
2487087cfbeSBarry Smith PetscErrorCode  MatFDColoringSetFunction(MatFDColoring matfd,PetscErrorCode (*f)(void),void *fctx)
249005c665bSBarry Smith {
2503a40ed3dSBarry Smith   PetscFunctionBegin;
2510700a824SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
252005c665bSBarry Smith   matfd->f    = f;
253005c665bSBarry Smith   matfd->fctx = fctx;
2543a40ed3dSBarry Smith   PetscFunctionReturn(0);
255005c665bSBarry Smith }
256005c665bSBarry Smith 
2574a2ae208SSatish Balay #undef __FUNCT__
2584a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFromOptions"
259639f9d9dSBarry Smith /*@
260b4fc646aSLois Curfman McInnes    MatFDColoringSetFromOptions - Sets coloring finite difference parameters from
261639f9d9dSBarry Smith    the options database.
262639f9d9dSBarry Smith 
263fee21e36SBarry Smith    Collective on MatFDColoring
264fee21e36SBarry Smith 
26565f2ba5bSLois Curfman McInnes    The Jacobian, F'(u), is estimated with the differencing approximation
266ef5ee4d1SLois Curfman McInnes .vb
26765f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
268f23b5b22SLois Curfman McInnes        h = error_rel*u[i]                 if  abs(u[i]) > umin
269f23b5b22SLois Curfman McInnes          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
270ef5ee4d1SLois Curfman McInnes        dx_{i} = (0, ... 1, .... 0)
271ef5ee4d1SLois Curfman McInnes .ve
272ef5ee4d1SLois Curfman McInnes 
273ef5ee4d1SLois Curfman McInnes    Input Parameter:
274ef5ee4d1SLois Curfman McInnes .  coloring - the coloring context
275ef5ee4d1SLois Curfman McInnes 
276b4fc646aSLois Curfman McInnes    Options Database Keys:
277d15ffeeaSSatish Balay +  -mat_fd_coloring_err <err> - Sets <err> (square root
278ef5ee4d1SLois Curfman McInnes            of relative error in the function)
279f23b5b22SLois Curfman McInnes .  -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude
2803ec795f1SBarry Smith .  -mat_fd_type - "wp" or "ds" (see MATMFFD_WP or MATMFFD_DS)
281ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_view - Activates basic viewing
282ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_view_info - Activates viewing info
283ef5ee4d1SLois Curfman McInnes -  -mat_fd_coloring_view_draw - Activates drawing
284639f9d9dSBarry Smith 
28515091d37SBarry Smith     Level: intermediate
28615091d37SBarry Smith 
287b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, parameters
288d1c39f55SBarry Smith 
289d1c39f55SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringView(), MatFDColoringSetParameters()
290d1c39f55SBarry Smith 
291639f9d9dSBarry Smith @*/
2927087cfbeSBarry Smith PetscErrorCode  MatFDColoringSetFromOptions(MatFDColoring matfd)
293639f9d9dSBarry Smith {
294dfbe8321SBarry Smith   PetscErrorCode ierr;
295ace3abfcSBarry Smith   PetscBool      flg;
296efb30889SBarry Smith   char           value[3];
2973a40ed3dSBarry Smith 
2983a40ed3dSBarry Smith   PetscFunctionBegin;
2990700a824SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
300639f9d9dSBarry Smith 
3013194b578SJed Brown   ierr = PetscObjectOptionsBegin((PetscObject)matfd);CHKERRQ(ierr);
30287828ca2SBarry Smith     ierr = PetscOptionsReal("-mat_fd_coloring_err","Square root of relative error in function","MatFDColoringSetParameters",matfd->error_rel,&matfd->error_rel,0);CHKERRQ(ierr);
30387828ca2SBarry Smith     ierr = PetscOptionsReal("-mat_fd_coloring_umin","Minimum allowable u magnitude","MatFDColoringSetParameters",matfd->umin,&matfd->umin,0);CHKERRQ(ierr);
3041bc75a8dSBarry Smith     ierr = PetscOptionsString("-mat_fd_type","Algorithm to compute h, wp or ds","MatFDColoringCreate",matfd->htype,value,3,&flg);CHKERRQ(ierr);
305efb30889SBarry Smith     if (flg) {
306efb30889SBarry Smith       if (value[0] == 'w' && value[1] == 'p') matfd->htype = "wp";
307efb30889SBarry Smith       else if (value[0] == 'd' && value[1] == 's') matfd->htype = "ds";
308e32f2f54SBarry Smith       else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Unknown finite differencing type %s",value);
309efb30889SBarry Smith     }
310186905e3SBarry Smith     /* not used here; but so they are presented in the GUI */
311b0a32e0cSBarry Smith     ierr = PetscOptionsName("-mat_fd_coloring_view","Print entire datastructure used for Jacobian","None",0);CHKERRQ(ierr);
312b0a32e0cSBarry Smith     ierr = PetscOptionsName("-mat_fd_coloring_view_info","Print number of colors etc for Jacobian","None",0);CHKERRQ(ierr);
313b0a32e0cSBarry Smith     ierr = PetscOptionsName("-mat_fd_coloring_view_draw","Plot nonzero structure ofJacobian","None",0);CHKERRQ(ierr);
3145d973c19SBarry Smith 
3155d973c19SBarry Smith     /* process any options handlers added with PetscObjectAddOptionsHandler() */
3165d973c19SBarry Smith     ierr = PetscObjectProcessOptionsHandlers((PetscObject)matfd);CHKERRQ(ierr);
317b0a32e0cSBarry Smith   PetscOptionsEnd();CHKERRQ(ierr);
3183a40ed3dSBarry Smith   PetscFunctionReturn(0);
319005c665bSBarry Smith }
320005c665bSBarry Smith 
3214a2ae208SSatish Balay #undef __FUNCT__
3224a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Private"
323dfbe8321SBarry Smith PetscErrorCode MatFDColoringView_Private(MatFDColoring fd)
324005c665bSBarry Smith {
325dfbe8321SBarry Smith   PetscErrorCode ierr;
326ace3abfcSBarry Smith   PetscBool      flg = PETSC_FALSE;
3273050cee2SBarry Smith   PetscViewer    viewer;
328005c665bSBarry Smith 
3293a40ed3dSBarry Smith   PetscFunctionBegin;
3307adad957SLisandro Dalcin   ierr = PetscViewerASCIIGetStdout(((PetscObject)fd)->comm,&viewer);CHKERRQ(ierr);
331acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-mat_fd_coloring_view",&flg,PETSC_NULL);CHKERRQ(ierr);
332005c665bSBarry Smith   if (flg) {
3333050cee2SBarry Smith     ierr = MatFDColoringView(fd,viewer);CHKERRQ(ierr);
334005c665bSBarry Smith   }
33590d69ab7SBarry Smith   flg  = PETSC_FALSE;
336acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-mat_fd_coloring_view_info",&flg,PETSC_NULL);CHKERRQ(ierr);
337ae09f205SBarry Smith   if (flg) {
3383050cee2SBarry Smith     ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
3393050cee2SBarry Smith     ierr = MatFDColoringView(fd,viewer);CHKERRQ(ierr);
3403050cee2SBarry Smith     ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
341ae09f205SBarry Smith   }
34290d69ab7SBarry Smith   flg  = PETSC_FALSE;
343acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg,PETSC_NULL);CHKERRQ(ierr);
344005c665bSBarry Smith   if (flg) {
3457adad957SLisandro Dalcin     ierr = MatFDColoringView(fd,PETSC_VIEWER_DRAW_(((PetscObject)fd)->comm));CHKERRQ(ierr);
3467adad957SLisandro Dalcin     ierr = PetscViewerFlush(PETSC_VIEWER_DRAW_(((PetscObject)fd)->comm));CHKERRQ(ierr);
347005c665bSBarry Smith   }
3483a40ed3dSBarry Smith   PetscFunctionReturn(0);
349bbf0e169SBarry Smith }
350bbf0e169SBarry Smith 
3514a2ae208SSatish Balay #undef __FUNCT__
3524a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringCreate"
35305869f15SSatish Balay /*@
354639f9d9dSBarry Smith    MatFDColoringCreate - Creates a matrix coloring context for finite difference
355639f9d9dSBarry Smith    computation of Jacobians.
356bbf0e169SBarry Smith 
357ef5ee4d1SLois Curfman McInnes    Collective on Mat
358ef5ee4d1SLois Curfman McInnes 
359639f9d9dSBarry Smith    Input Parameters:
360ef5ee4d1SLois Curfman McInnes +  mat - the matrix containing the nonzero structure of the Jacobian
361e727c939SJed Brown -  iscoloring - the coloring of the matrix; usually obtained with MatGetColoring() or DMCreateColoring()
362bbf0e169SBarry Smith 
363bbf0e169SBarry Smith     Output Parameter:
364639f9d9dSBarry Smith .   color - the new coloring context
365bbf0e169SBarry Smith 
36615091d37SBarry Smith     Level: intermediate
36715091d37SBarry Smith 
3686831982aSBarry Smith .seealso: MatFDColoringDestroy(),SNESDefaultComputeJacobianColor(), ISColoringCreate(),
369d1c39f55SBarry Smith           MatFDColoringSetFunction(), MatFDColoringSetFromOptions(), MatFDColoringApply(),
370e727c939SJed Brown           MatFDColoringView(), MatFDColoringSetParameters(), MatGetColoring(), DMCreateColoring()
371bbf0e169SBarry Smith @*/
3727087cfbeSBarry Smith PetscErrorCode  MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
373bbf0e169SBarry Smith {
374639f9d9dSBarry Smith   MatFDColoring  c;
375639f9d9dSBarry Smith   MPI_Comm       comm;
376dfbe8321SBarry Smith   PetscErrorCode ierr;
37738baddfdSBarry Smith   PetscInt       M,N;
378639f9d9dSBarry Smith 
3793a40ed3dSBarry Smith   PetscFunctionBegin;
380d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
381639f9d9dSBarry Smith   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
38217186662SBarry Smith   if (M != N) SETERRQ(((PetscObject)mat)->comm,PETSC_ERR_SUP,"Only for square matrices");
383f881d145SBarry Smith   ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
3843194b578SJed Brown   ierr = PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_CLASSID,0,"MatFDColoring","Jacobian computation via finite differences with coloring","Mat",comm,MatFDColoringDestroy,MatFDColoringView);CHKERRQ(ierr);
385639f9d9dSBarry Smith 
386b8f8c88eSHong Zhang   c->ctype = iscoloring->ctype;
387b8f8c88eSHong Zhang 
388f830108cSBarry Smith   if (mat->ops->fdcoloringcreate) {
389f830108cSBarry Smith     ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);CHKERRQ(ierr);
390577b14d0SJed Brown   } else SETERRQ1(((PetscObject)mat)->comm,PETSC_ERR_SUP,"Code not yet written for matrix type %s",((PetscObject)mat)->type_name);
391639f9d9dSBarry Smith 
392b8f8c88eSHong Zhang   ierr = MatGetVecs(mat,PETSC_NULL,&c->w1);CHKERRQ(ierr);
393b8f8c88eSHong Zhang   ierr = PetscLogObjectParent(c,c->w1);CHKERRQ(ierr);
394b8f8c88eSHong Zhang   ierr = VecDuplicate(c->w1,&c->w2);CHKERRQ(ierr);
395b8f8c88eSHong Zhang   ierr = PetscLogObjectParent(c,c->w2);CHKERRQ(ierr);
396b8f8c88eSHong Zhang 
39777d8c4bbSBarry Smith   c->error_rel         = PETSC_SQRT_MACHINE_EPSILON;
39877d8c4bbSBarry Smith   c->umin              = 100.0*PETSC_SQRT_MACHINE_EPSILON;
39949b058dcSBarry Smith   c->currentcolor      = -1;
400efb30889SBarry Smith   c->htype             = "wp";
401*4e269d77SPeter Brune   c->fset              = PETSC_FALSE;
402639f9d9dSBarry Smith 
403639f9d9dSBarry Smith   *color = c;
404*4e269d77SPeter Brune   ierr = PetscObjectCompose((PetscObject)mat,"SNESMatFDColoring",(PetscObject)c);CHKERRQ(ierr);
405d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
4063a40ed3dSBarry Smith   PetscFunctionReturn(0);
407bbf0e169SBarry Smith }
408bbf0e169SBarry Smith 
4094a2ae208SSatish Balay #undef __FUNCT__
4104a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringDestroy"
41105869f15SSatish Balay /*@
412639f9d9dSBarry Smith     MatFDColoringDestroy - Destroys a matrix coloring context that was created
413639f9d9dSBarry Smith     via MatFDColoringCreate().
414bbf0e169SBarry Smith 
415ef5ee4d1SLois Curfman McInnes     Collective on MatFDColoring
416ef5ee4d1SLois Curfman McInnes 
417b4fc646aSLois Curfman McInnes     Input Parameter:
418639f9d9dSBarry Smith .   c - coloring context
419bbf0e169SBarry Smith 
42015091d37SBarry Smith     Level: intermediate
42115091d37SBarry Smith 
422639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
423bbf0e169SBarry Smith @*/
4246bf464f9SBarry Smith PetscErrorCode  MatFDColoringDestroy(MatFDColoring *c)
425bbf0e169SBarry Smith {
4266849ba73SBarry Smith   PetscErrorCode ierr;
42738baddfdSBarry Smith   PetscInt       i;
428bbf0e169SBarry Smith 
4293a40ed3dSBarry Smith   PetscFunctionBegin;
4306bf464f9SBarry Smith   if (!*c) PetscFunctionReturn(0);
4316bf464f9SBarry Smith   if (--((PetscObject)(*c))->refct > 0) {*c = 0; PetscFunctionReturn(0);}
432d4bb536fSBarry Smith 
4336bf464f9SBarry Smith   for (i=0; i<(*c)->ncolors; i++) {
4346bf464f9SBarry Smith     ierr = PetscFree((*c)->columns[i]);CHKERRQ(ierr);
4356bf464f9SBarry Smith     ierr = PetscFree((*c)->rows[i]);CHKERRQ(ierr);
4366bf464f9SBarry Smith     ierr = PetscFree((*c)->columnsforrow[i]);CHKERRQ(ierr);
4376bf464f9SBarry Smith     if ((*c)->vscaleforrow) {ierr = PetscFree((*c)->vscaleforrow[i]);CHKERRQ(ierr);}
438bbf0e169SBarry Smith   }
4396bf464f9SBarry Smith   ierr = PetscFree((*c)->ncolumns);CHKERRQ(ierr);
4406bf464f9SBarry Smith   ierr = PetscFree((*c)->columns);CHKERRQ(ierr);
4416bf464f9SBarry Smith   ierr = PetscFree((*c)->nrows);CHKERRQ(ierr);
4426bf464f9SBarry Smith   ierr = PetscFree((*c)->rows);CHKERRQ(ierr);
4436bf464f9SBarry Smith   ierr = PetscFree((*c)->columnsforrow);CHKERRQ(ierr);
4446bf464f9SBarry Smith   ierr = PetscFree((*c)->vscaleforrow);CHKERRQ(ierr);
4456bf464f9SBarry Smith   ierr = VecDestroy(&(*c)->vscale);CHKERRQ(ierr);
4466bf464f9SBarry Smith   ierr = VecDestroy(&(*c)->w1);CHKERRQ(ierr);
4476bf464f9SBarry Smith   ierr = VecDestroy(&(*c)->w2);CHKERRQ(ierr);
4486bf464f9SBarry Smith   ierr = VecDestroy(&(*c)->w3);CHKERRQ(ierr);
449d38fa0fbSBarry Smith   ierr = PetscHeaderDestroy(c);CHKERRQ(ierr);
4503a40ed3dSBarry Smith   PetscFunctionReturn(0);
451bbf0e169SBarry Smith }
45243a90d84SBarry Smith 
4534a2ae208SSatish Balay #undef __FUNCT__
45449b058dcSBarry Smith #define __FUNCT__ "MatFDColoringGetPerturbedColumns"
45549b058dcSBarry Smith /*@C
45649b058dcSBarry Smith     MatFDColoringGetPerturbedColumns - Returns the indices of the columns that
45749b058dcSBarry Smith       that are currently being perturbed.
45849b058dcSBarry Smith 
45949b058dcSBarry Smith     Not Collective
46049b058dcSBarry Smith 
46149b058dcSBarry Smith     Input Parameters:
46249b058dcSBarry Smith .   coloring - coloring context created with MatFDColoringCreate()
46349b058dcSBarry Smith 
46449b058dcSBarry Smith     Output Parameters:
46549b058dcSBarry Smith +   n - the number of local columns being perturbed
46649b058dcSBarry Smith -   cols - the column indices, in global numbering
46749b058dcSBarry Smith 
46849b058dcSBarry Smith    Level: intermediate
46949b058dcSBarry Smith 
47049b058dcSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringApply()
47149b058dcSBarry Smith 
47249b058dcSBarry Smith .keywords: coloring, Jacobian, finite differences
47349b058dcSBarry Smith @*/
4747087cfbeSBarry Smith PetscErrorCode  MatFDColoringGetPerturbedColumns(MatFDColoring coloring,PetscInt *n,PetscInt *cols[])
47549b058dcSBarry Smith {
47649b058dcSBarry Smith   PetscFunctionBegin;
47749b058dcSBarry Smith   if (coloring->currentcolor >= 0) {
47849b058dcSBarry Smith     *n    = coloring->ncolumns[coloring->currentcolor];
47949b058dcSBarry Smith     *cols = coloring->columns[coloring->currentcolor];
48049b058dcSBarry Smith   } else {
48149b058dcSBarry Smith     *n = 0;
48249b058dcSBarry Smith   }
48349b058dcSBarry Smith   PetscFunctionReturn(0);
48449b058dcSBarry Smith }
48549b058dcSBarry Smith 
48649b058dcSBarry Smith #undef __FUNCT__
4874a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringApply"
48843a90d84SBarry Smith /*@
489e0907662SLois Curfman McInnes     MatFDColoringApply - Given a matrix for which a MatFDColoring context
490e0907662SLois Curfman McInnes     has been created, computes the Jacobian for a function via finite differences.
49143a90d84SBarry Smith 
492fee21e36SBarry Smith     Collective on MatFDColoring
493fee21e36SBarry Smith 
494ef5ee4d1SLois Curfman McInnes     Input Parameters:
495ef5ee4d1SLois Curfman McInnes +   mat - location to store Jacobian
496ef5ee4d1SLois Curfman McInnes .   coloring - coloring context created with MatFDColoringCreate()
497ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
4987850c7c0SBarry Smith -   sctx - context required by function, if this is being used with the SNES solver then it is SNES object, otherwise it is null
499ef5ee4d1SLois Curfman McInnes 
5008bba8e72SBarry Smith     Options Database Keys:
501ebd3b9afSBarry Smith +    -mat_fd_type - "wp" or "ds"  (see MATMFFD_WP or MATMFFD_DS)
502b9722508SLois Curfman McInnes .    -mat_fd_coloring_view - Activates basic viewing or coloring
503b9722508SLois Curfman McInnes .    -mat_fd_coloring_view_draw - Activates drawing of coloring
504b9722508SLois Curfman McInnes -    -mat_fd_coloring_view_info - Activates viewing of coloring info
5058bba8e72SBarry Smith 
50615091d37SBarry Smith     Level: intermediate
50798d79826SSatish Balay 
5087850c7c0SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringSetFunction()
50943a90d84SBarry Smith 
51043a90d84SBarry Smith .keywords: coloring, Jacobian, finite differences
51143a90d84SBarry Smith @*/
5127087cfbeSBarry Smith PetscErrorCode  MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
51343a90d84SBarry Smith {
5143acb8795SBarry Smith   PetscErrorCode ierr;
5153acb8795SBarry Smith 
5163acb8795SBarry Smith   PetscFunctionBegin;
5170700a824SBarry Smith   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
5180700a824SBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_CLASSID,2);
5190700a824SBarry Smith   PetscValidHeaderSpecific(x1,VEC_CLASSID,3);
52017186662SBarry Smith   if (!coloring->f) SETERRQ(((PetscObject)J)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetFunction()");
521e32f2f54SBarry Smith   if (!J->ops->fdcoloringapply) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not supported for this matrix type %s",((PetscObject)J)->type_name);
5223acb8795SBarry Smith   ierr = (*J->ops->fdcoloringapply)(J,coloring,x1,flag,sctx);CHKERRQ(ierr);
5233acb8795SBarry Smith   PetscFunctionReturn(0);
5243acb8795SBarry Smith }
5253acb8795SBarry Smith 
5263acb8795SBarry Smith #undef __FUNCT__
5273acb8795SBarry Smith #define __FUNCT__ "MatFDColoringApply_AIJ"
5287087cfbeSBarry Smith PetscErrorCode  MatFDColoringApply_AIJ(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
5293acb8795SBarry Smith {
5306849ba73SBarry Smith   PetscErrorCode (*f)(void*,Vec,Vec,void*) = (PetscErrorCode (*)(void*,Vec,Vec,void *))coloring->f;
5316849ba73SBarry Smith   PetscErrorCode ierr;
532*4e269d77SPeter Brune   PetscInt       k,start,end,l,row,col,srow,**vscaleforrow;
533efb30889SBarry Smith   PetscScalar    dx,*y,*xx,*w3_array;
534ea709b57SSatish Balay   PetscScalar    *vscale_array;
535efb30889SBarry Smith   PetscReal      epsilon = coloring->error_rel,umin = coloring->umin,unorm;
536b8f8c88eSHong Zhang   Vec            w1=coloring->w1,w2=coloring->w2,w3;
537005c665bSBarry Smith   void           *fctx = coloring->fctx;
538ace3abfcSBarry Smith   PetscBool      flg = PETSC_FALSE;
539705d48f7SSatish Balay   PetscInt       ctype=coloring->ctype,N,col_start=0,col_end=0;
540b8f8c88eSHong Zhang   Vec            x1_tmp;
541a2e34c3dSBarry Smith 
5423a40ed3dSBarry Smith   PetscFunctionBegin;
5430700a824SBarry Smith   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
5440700a824SBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_CLASSID,2);
5450700a824SBarry Smith   PetscValidHeaderSpecific(x1,VEC_CLASSID,3);
54617186662SBarry Smith   if (!f) SETERRQ(((PetscObject)J)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetFunction()");
547e0907662SLois Curfman McInnes 
548d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
5494c49b128SBarry Smith   ierr = MatSetUnfactored(J);CHKERRQ(ierr);
550acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg,PETSC_NULL);CHKERRQ(ierr);
551f1af5d2fSBarry Smith   if (flg) {
552ae15b995SBarry Smith     ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr);
553e0907662SLois Curfman McInnes   } else {
554ace3abfcSBarry Smith     PetscBool  assembled;
5550b9b6f31SBarry Smith     ierr = MatAssembled(J,&assembled);CHKERRQ(ierr);
5560b9b6f31SBarry Smith     if (assembled) {
55743a90d84SBarry Smith       ierr = MatZeroEntries(J);CHKERRQ(ierr);
558e0907662SLois Curfman McInnes     }
5590b9b6f31SBarry Smith   }
56043a90d84SBarry Smith 
561b8f8c88eSHong Zhang   x1_tmp = x1;
562dac7f5fdSBarry Smith   if (!coloring->vscale){
563b8f8c88eSHong Zhang     ierr = VecDuplicate(x1_tmp,&coloring->vscale);CHKERRQ(ierr);
564b8f8c88eSHong Zhang   }
565be2fbe1fSBarry Smith 
566b8f8c88eSHong Zhang   if (coloring->htype[0] == 'w') { /* tacky test; need to make systematic if we add other approaches to computing h*/
567b8f8c88eSHong Zhang     ierr = VecNorm(x1_tmp,NORM_2,&unorm);CHKERRQ(ierr);
568b8f8c88eSHong Zhang   }
569b8f8c88eSHong Zhang   ierr = VecGetOwnershipRange(w1,&start,&end);CHKERRQ(ierr); /* OwnershipRange is used by ghosted x! */
570b8f8c88eSHong Zhang 
571b8f8c88eSHong Zhang   /* Set w1 = F(x1) */
572*4e269d77SPeter Brune   if (!coloring->fset) {
57366f9b7ceSBarry Smith     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
574b8f8c88eSHong Zhang     ierr = (*f)(sctx,x1_tmp,w1,fctx);CHKERRQ(ierr);
57566f9b7ceSBarry Smith     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
576*4e269d77SPeter Brune   } else {
577*4e269d77SPeter Brune     coloring->fset = PETSC_FALSE;
578be2fbe1fSBarry Smith }
57943a90d84SBarry Smith 
580b8f8c88eSHong Zhang   if (!coloring->w3) {
581b8f8c88eSHong Zhang     ierr = VecDuplicate(x1_tmp,&coloring->w3);CHKERRQ(ierr);
582b8f8c88eSHong Zhang     ierr = PetscLogObjectParent(coloring,coloring->w3);CHKERRQ(ierr);
583efb30889SBarry Smith   }
584b8f8c88eSHong Zhang   w3 = coloring->w3;
585efb30889SBarry Smith 
586b8f8c88eSHong Zhang     /* Compute all the local scale factors, including ghost points */
587b8f8c88eSHong Zhang   ierr = VecGetLocalSize(x1_tmp,&N);CHKERRQ(ierr);
588b8f8c88eSHong Zhang   ierr = VecGetArray(x1_tmp,&xx);CHKERRQ(ierr);
589b8f8c88eSHong Zhang   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
590b8f8c88eSHong Zhang   if (ctype == IS_COLORING_GHOSTED){
591b8f8c88eSHong Zhang     col_start = 0; col_end = N;
5928ee2e534SBarry Smith   } else if (ctype == IS_COLORING_GLOBAL){
593b8f8c88eSHong Zhang     xx = xx - start;
594b8f8c88eSHong Zhang     vscale_array = vscale_array - start;
595b8f8c88eSHong Zhang     col_start = start; col_end = N + start;
596b8f8c88eSHong Zhang   }
597b8f8c88eSHong Zhang   for (col=col_start; col<col_end; col++){
598b8f8c88eSHong Zhang     /* Loop over each local column, vscale[col] = 1./(epsilon*dx[col]) */
599efb30889SBarry Smith     if (coloring->htype[0] == 'w') {
600efb30889SBarry Smith       dx = 1.0 + unorm;
601efb30889SBarry Smith     } else {
6029782cf98SBarry Smith       dx  = xx[col];
603efb30889SBarry Smith     }
604d4a378daSJed Brown     if (dx == (PetscScalar)0.0) dx = 1.0;
6059782cf98SBarry Smith #if !defined(PETSC_USE_COMPLEX)
6069782cf98SBarry Smith     if (dx < umin && dx >= 0.0)      dx = umin;
6079782cf98SBarry Smith     else if (dx < 0.0 && dx > -umin) dx = -umin;
6089782cf98SBarry Smith #else
6099782cf98SBarry Smith     if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
6109782cf98SBarry Smith     else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
6119782cf98SBarry Smith #endif
6129782cf98SBarry Smith     dx               *= epsilon;
613d4a378daSJed Brown     vscale_array[col] = (PetscScalar)1.0/dx;
6149782cf98SBarry Smith   }
6158ee2e534SBarry Smith   if (ctype == IS_COLORING_GLOBAL)  vscale_array = vscale_array + start;
616efb30889SBarry Smith   ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
6178ee2e534SBarry Smith   if (ctype == IS_COLORING_GLOBAL){
61830b34957SBarry Smith     ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
61930b34957SBarry Smith     ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
620b8f8c88eSHong Zhang   }
6219782cf98SBarry Smith 
622b8f8c88eSHong Zhang   if (coloring->vscaleforrow) {
623b8f8c88eSHong Zhang     vscaleforrow = coloring->vscaleforrow;
62417186662SBarry Smith   } else SETERRQ(((PetscObject)J)->comm,PETSC_ERR_ARG_NULL,"Null Object: coloring->vscaleforrow");
625b0a32e0cSBarry Smith 
6269782cf98SBarry Smith   /*
62743a90d84SBarry Smith     Loop over each color
62843a90d84SBarry Smith   */
629b8f8c88eSHong Zhang   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
63043a90d84SBarry Smith   for (k=0; k<coloring->ncolors; k++) {
63149b058dcSBarry Smith     coloring->currentcolor = k;
632b8f8c88eSHong Zhang     ierr = VecCopy(x1_tmp,w3);CHKERRQ(ierr);
633b8f8c88eSHong Zhang     ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);
6348ee2e534SBarry Smith     if (ctype == IS_COLORING_GLOBAL) w3_array = w3_array - start;
63543a90d84SBarry Smith     /*
636b8f8c88eSHong Zhang       Loop over each column associated with color
637b8f8c88eSHong Zhang       adding the perturbation to the vector w3.
63843a90d84SBarry Smith     */
63943a90d84SBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
640b8f8c88eSHong Zhang       col = coloring->columns[k][l];    /* local column of the matrix we are probing for */
641efb30889SBarry Smith       if (coloring->htype[0] == 'w') {
642efb30889SBarry Smith         dx = 1.0 + unorm;
643efb30889SBarry Smith       } else {
64442460c72SBarry Smith         dx  = xx[col];
645efb30889SBarry Smith       }
646d4a378daSJed Brown       if (dx == (PetscScalar)0.0) dx = 1.0;
647aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
648ae09f205SBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
649ae09f205SBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
65043a90d84SBarry Smith #else
651329f5518SBarry Smith       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
652329f5518SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
65343a90d84SBarry Smith #endif
65443a90d84SBarry Smith       dx            *= epsilon;
655e32f2f54SBarry Smith       if (!PetscAbsScalar(dx)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Computed 0 differencing parameter");
65642460c72SBarry Smith       w3_array[col] += dx;
65743a90d84SBarry Smith     }
65898d79826SSatish Balay     if (ctype == IS_COLORING_GLOBAL) w3_array = w3_array + start;
659b8f8c88eSHong Zhang     ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
6603b28642cSBarry Smith 
66143a90d84SBarry Smith     /*
662b8f8c88eSHong Zhang       Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations)
663b8f8c88eSHong Zhang                            w2 = F(x1 + dx) - F(x1)
66443a90d84SBarry Smith     */
66566f9b7ceSBarry Smith     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
66643a90d84SBarry Smith     ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr);
66766f9b7ceSBarry Smith     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
668efb30889SBarry Smith     ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr);
6699782cf98SBarry Smith 
67043a90d84SBarry Smith     /*
671e0907662SLois Curfman McInnes       Loop over rows of vector, putting results into Jacobian matrix
67243a90d84SBarry Smith     */
6733b28642cSBarry Smith     ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
67443a90d84SBarry Smith     for (l=0; l<coloring->nrows[k]; l++) {
675b8f8c88eSHong Zhang       row    = coloring->rows[k][l];             /* local row index */
676b8f8c88eSHong Zhang       col    = coloring->columnsforrow[k][l];    /* global column index */
6775904e1b1SBarry Smith       y[row] *= vscale_array[vscaleforrow[k][l]];
67843a90d84SBarry Smith       srow   = row + start;
67943a90d84SBarry Smith       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
68043a90d84SBarry Smith     }
6813b28642cSBarry Smith     ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
682aeb7e98aSMatthew Knepley   } /* endof for each color */
6838ee2e534SBarry Smith   if (ctype == IS_COLORING_GLOBAL) xx = xx + start;
68430b34957SBarry Smith   ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
685b8f8c88eSHong Zhang   ierr = VecRestoreArray(x1_tmp,&xx);CHKERRQ(ierr);
686b8f8c88eSHong Zhang 
687b8f8c88eSHong Zhang   coloring->currentcolor = -1;
68843a90d84SBarry Smith   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
68943a90d84SBarry Smith   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
690d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
691a2e34c3dSBarry Smith 
69290d69ab7SBarry Smith   flg  = PETSC_FALSE;
693acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-mat_null_space_test",&flg,PETSC_NULL);CHKERRQ(ierr);
694a2e34c3dSBarry Smith   if (flg) {
69595902228SMatthew Knepley     ierr = MatNullSpaceTest(J->nullsp,J,PETSC_NULL);CHKERRQ(ierr);
696a2e34c3dSBarry Smith   }
697b9722508SLois Curfman McInnes   ierr = MatFDColoringView_Private(coloring);CHKERRQ(ierr);
6983a40ed3dSBarry Smith   PetscFunctionReturn(0);
69943a90d84SBarry Smith }
700