xref: /petsc/src/mat/matfd/fdmatrix.c (revision ce94432eddcd14845bc7e8083b7f8ea723b9bf7d)
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 {
134e269d77SPeter Brune   PetscErrorCode ierr;
146e111a19SKarl Rupp 
153a7fca6bSBarry Smith   PetscFunctionBegin;
164e269d77SPeter Brune   if (F) {
174e269d77SPeter Brune     ierr     = VecCopy(F,fd->w1);CHKERRQ(ierr);
184e269d77SPeter Brune     fd->fset = PETSC_TRUE;
194e269d77SPeter Brune   } else {
204e269d77SPeter Brune     fd->fset = PETSC_FALSE;
214e269d77SPeter Brune   }
223a7fca6bSBarry Smith   PetscFunctionReturn(0);
233a7fca6bSBarry Smith }
243a7fca6bSBarry Smith 
253a7fca6bSBarry Smith #undef __FUNCT__
264a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Draw_Zoom"
276849ba73SBarry Smith static PetscErrorCode MatFDColoringView_Draw_Zoom(PetscDraw draw,void *Aa)
289194eea9SBarry Smith {
299194eea9SBarry Smith   MatFDColoring  fd = (MatFDColoring)Aa;
30dfbe8321SBarry Smith   PetscErrorCode ierr;
3138baddfdSBarry Smith   PetscInt       i,j;
329194eea9SBarry Smith   PetscReal      x,y;
339194eea9SBarry Smith 
349194eea9SBarry Smith   PetscFunctionBegin;
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);
65f77a5242SKarl Rupp   ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",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) {
110*ce94432eSBarry Smith     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)c),&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);
140bbf0e169SBarry Smith   }
1413a40ed3dSBarry Smith   PetscFunctionReturn(0);
142639f9d9dSBarry Smith }
143639f9d9dSBarry Smith 
1444a2ae208SSatish Balay #undef __FUNCT__
1454a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetParameters"
146639f9d9dSBarry Smith /*@
147b4fc646aSLois Curfman McInnes    MatFDColoringSetParameters - Sets the parameters for the sparse approximation of
148b4fc646aSLois Curfman McInnes    a Jacobian matrix using finite differences.
149639f9d9dSBarry Smith 
1503f9fe445SBarry Smith    Logically Collective on MatFDColoring
151ef5ee4d1SLois Curfman McInnes 
152ef5ee4d1SLois Curfman McInnes    The Jacobian is estimated with the differencing approximation
153ef5ee4d1SLois Curfman McInnes .vb
15465f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
155f23b5b22SLois Curfman McInnes        h = error_rel*u[i]                 if  abs(u[i]) > umin
156f23b5b22SLois Curfman McInnes          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
157ef5ee4d1SLois Curfman McInnes        dx_{i} = (0, ... 1, .... 0)
158ef5ee4d1SLois Curfman McInnes .ve
159639f9d9dSBarry Smith 
160639f9d9dSBarry Smith    Input Parameters:
161ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
162639f9d9dSBarry Smith .  error_rel - relative error
163f23b5b22SLois Curfman McInnes -  umin - minimum allowable u-value magnitude
164fee21e36SBarry Smith 
16515091d37SBarry Smith    Level: advanced
16615091d37SBarry Smith 
167b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, parameters
168b4fc646aSLois Curfman McInnes 
16995b89fc3SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringSetFromOptions()
17095b89fc3SBarry Smith 
171639f9d9dSBarry Smith @*/
1727087cfbeSBarry Smith PetscErrorCode  MatFDColoringSetParameters(MatFDColoring matfd,PetscReal error,PetscReal umin)
173639f9d9dSBarry Smith {
1743a40ed3dSBarry Smith   PetscFunctionBegin;
1750700a824SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
176c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(matfd,error,2);
177c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(matfd,umin,3);
178639f9d9dSBarry Smith   if (error != PETSC_DEFAULT) matfd->error_rel = error;
179639f9d9dSBarry Smith   if (umin != PETSC_DEFAULT)  matfd->umin      = umin;
1803a40ed3dSBarry Smith   PetscFunctionReturn(0);
181639f9d9dSBarry Smith }
182639f9d9dSBarry Smith 
183005c665bSBarry Smith 
184ff0cfa39SBarry Smith 
1854a2ae208SSatish Balay #undef __FUNCT__
1864a9d489dSBarry Smith #define __FUNCT__ "MatFDColoringGetFunction"
1874a9d489dSBarry Smith /*@C
1884a9d489dSBarry Smith    MatFDColoringGetFunction - Gets the function to use for computing the Jacobian.
1894a9d489dSBarry Smith 
1903f9fe445SBarry Smith    Not Collective
1914a9d489dSBarry Smith 
1924a9d489dSBarry Smith    Input Parameters:
1934a9d489dSBarry Smith .  coloring - the coloring context
1944a9d489dSBarry Smith 
1954a9d489dSBarry Smith    Output Parameters:
1964a9d489dSBarry Smith +  f - the function
1974a9d489dSBarry Smith -  fctx - the optional user-defined function context
1984a9d489dSBarry Smith 
1994a9d489dSBarry Smith    Level: intermediate
2004a9d489dSBarry Smith 
2014a9d489dSBarry Smith .keywords: Mat, Jacobian, finite differences, set, function
20295b89fc3SBarry Smith 
20395b89fc3SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringSetFunction(), MatFDColoringSetFromOptions()
20495b89fc3SBarry Smith 
2054a9d489dSBarry Smith @*/
2067087cfbeSBarry Smith PetscErrorCode  MatFDColoringGetFunction(MatFDColoring matfd,PetscErrorCode (**f)(void),void **fctx)
2074a9d489dSBarry Smith {
2084a9d489dSBarry Smith   PetscFunctionBegin;
2090700a824SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
2104a9d489dSBarry Smith   if (f) *f = matfd->f;
2114a9d489dSBarry Smith   if (fctx) *fctx = matfd->fctx;
2124a9d489dSBarry Smith   PetscFunctionReturn(0);
2134a9d489dSBarry Smith }
2144a9d489dSBarry Smith 
2154a9d489dSBarry Smith #undef __FUNCT__
2164a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFunction"
217d64ed03dSBarry Smith /*@C
218005c665bSBarry Smith    MatFDColoringSetFunction - Sets the function to use for computing the Jacobian.
219005c665bSBarry Smith 
2203f9fe445SBarry Smith    Logically Collective on MatFDColoring
221fee21e36SBarry Smith 
222ef5ee4d1SLois Curfman McInnes    Input Parameters:
223ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
224ef5ee4d1SLois Curfman McInnes .  f - the function
225ef5ee4d1SLois Curfman McInnes -  fctx - the optional user-defined function context
226ef5ee4d1SLois Curfman McInnes 
2277850c7c0SBarry Smith    Calling sequence of (*f) function:
2287850c7c0SBarry Smith     For SNES:    PetscErrorCode (*f)(SNES,Vec,Vec,void*)
229ab637aeaSJed Brown     If not using SNES: PetscErrorCode (*f)(void *dummy,Vec,Vec,void*) and dummy is ignored
23015091d37SBarry Smith 
2317850c7c0SBarry Smith    Level: advanced
2327850c7c0SBarry Smith 
233ab637aeaSJed Brown    Notes: This function is usually used automatically by SNES (when one uses SNESSetJacobian() with the argument
234ab637aeaSJed Brown      SNESDefaultComputeJacobianColor()) and only needs to be used by someone computing a matrix via coloring directly by
235ab637aeaSJed Brown      calling MatFDColoringApply()
2367850c7c0SBarry Smith 
2377850c7c0SBarry Smith    Fortran Notes:
2387850c7c0SBarry Smith     In Fortran you must call MatFDColoringSetFunction() for a coloring object to
239ab637aeaSJed Brown   be used without SNES or within the SNES solvers.
240f881d145SBarry Smith 
241b4fc646aSLois Curfman McInnes .keywords: Mat, Jacobian, finite differences, set, function
24295b89fc3SBarry Smith 
24395b89fc3SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringGetFunction(), MatFDColoringSetFromOptions()
24495b89fc3SBarry Smith 
245005c665bSBarry Smith @*/
2467087cfbeSBarry Smith PetscErrorCode  MatFDColoringSetFunction(MatFDColoring matfd,PetscErrorCode (*f)(void),void *fctx)
247005c665bSBarry Smith {
2483a40ed3dSBarry Smith   PetscFunctionBegin;
2490700a824SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
250005c665bSBarry Smith   matfd->f    = f;
251005c665bSBarry Smith   matfd->fctx = fctx;
2523a40ed3dSBarry Smith   PetscFunctionReturn(0);
253005c665bSBarry Smith }
254005c665bSBarry Smith 
2554a2ae208SSatish Balay #undef __FUNCT__
2564a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFromOptions"
257639f9d9dSBarry Smith /*@
258b4fc646aSLois Curfman McInnes    MatFDColoringSetFromOptions - Sets coloring finite difference parameters from
259639f9d9dSBarry Smith    the options database.
260639f9d9dSBarry Smith 
261fee21e36SBarry Smith    Collective on MatFDColoring
262fee21e36SBarry Smith 
26365f2ba5bSLois Curfman McInnes    The Jacobian, F'(u), is estimated with the differencing approximation
264ef5ee4d1SLois Curfman McInnes .vb
26565f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
266f23b5b22SLois Curfman McInnes        h = error_rel*u[i]                 if  abs(u[i]) > umin
267f23b5b22SLois Curfman McInnes          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
268ef5ee4d1SLois Curfman McInnes        dx_{i} = (0, ... 1, .... 0)
269ef5ee4d1SLois Curfman McInnes .ve
270ef5ee4d1SLois Curfman McInnes 
271ef5ee4d1SLois Curfman McInnes    Input Parameter:
272ef5ee4d1SLois Curfman McInnes .  coloring - the coloring context
273ef5ee4d1SLois Curfman McInnes 
274b4fc646aSLois Curfman McInnes    Options Database Keys:
275d15ffeeaSSatish Balay +  -mat_fd_coloring_err <err> - Sets <err> (square root
276ef5ee4d1SLois Curfman McInnes            of relative error in the function)
277f23b5b22SLois Curfman McInnes .  -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude
2783ec795f1SBarry Smith .  -mat_fd_type - "wp" or "ds" (see MATMFFD_WP or MATMFFD_DS)
279ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_view - Activates basic viewing
280e350db43SBarry Smith .  -mat_fd_coloring_view ::ascii_info - Activates viewing info
281e350db43SBarry Smith -  -mat_fd_coloring_view draw - Activates drawing
282639f9d9dSBarry Smith 
28315091d37SBarry Smith     Level: intermediate
28415091d37SBarry Smith 
285b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, parameters
286d1c39f55SBarry Smith 
287d1c39f55SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringView(), MatFDColoringSetParameters()
288d1c39f55SBarry Smith 
289639f9d9dSBarry Smith @*/
2907087cfbeSBarry Smith PetscErrorCode  MatFDColoringSetFromOptions(MatFDColoring matfd)
291639f9d9dSBarry Smith {
292dfbe8321SBarry Smith   PetscErrorCode ierr;
293ace3abfcSBarry Smith   PetscBool      flg;
294efb30889SBarry Smith   char           value[3];
2953a40ed3dSBarry Smith 
2963a40ed3dSBarry Smith   PetscFunctionBegin;
2970700a824SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
298639f9d9dSBarry Smith 
2993194b578SJed Brown   ierr = PetscObjectOptionsBegin((PetscObject)matfd);CHKERRQ(ierr);
30087828ca2SBarry Smith   ierr = PetscOptionsReal("-mat_fd_coloring_err","Square root of relative error in function","MatFDColoringSetParameters",matfd->error_rel,&matfd->error_rel,0);CHKERRQ(ierr);
30187828ca2SBarry Smith   ierr = PetscOptionsReal("-mat_fd_coloring_umin","Minimum allowable u magnitude","MatFDColoringSetParameters",matfd->umin,&matfd->umin,0);CHKERRQ(ierr);
3021bc75a8dSBarry Smith   ierr = PetscOptionsString("-mat_fd_type","Algorithm to compute h, wp or ds","MatFDColoringCreate",matfd->htype,value,3,&flg);CHKERRQ(ierr);
303efb30889SBarry Smith   if (flg) {
304efb30889SBarry Smith     if (value[0] == 'w' && value[1] == 'p') matfd->htype = "wp";
305efb30889SBarry Smith     else if (value[0] == 'd' && value[1] == 's') matfd->htype = "ds";
306e32f2f54SBarry Smith     else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Unknown finite differencing type %s",value);
307efb30889SBarry Smith   }
3085d973c19SBarry Smith   /* process any options handlers added with PetscObjectAddOptionsHandler() */
3095d973c19SBarry Smith   ierr = PetscObjectProcessOptionsHandlers((PetscObject)matfd);CHKERRQ(ierr);
310b0a32e0cSBarry Smith   PetscOptionsEnd();CHKERRQ(ierr);
3113a40ed3dSBarry Smith   PetscFunctionReturn(0);
312005c665bSBarry Smith }
313005c665bSBarry Smith 
3144a2ae208SSatish Balay #undef __FUNCT__
315e350db43SBarry Smith #define __FUNCT__ "MatFDColoringViewFromOptions"
316e350db43SBarry Smith PetscErrorCode MatFDColoringViewFromOptions(MatFDColoring fd,const char optionname[])
317005c665bSBarry Smith {
318dfbe8321SBarry Smith   PetscErrorCode    ierr;
319e350db43SBarry Smith   PetscBool         flg;
3203050cee2SBarry Smith   PetscViewer       viewer;
321cffb1e40SBarry Smith   PetscViewerFormat format;
322005c665bSBarry Smith 
3233a40ed3dSBarry Smith   PetscFunctionBegin;
324*ce94432eSBarry Smith   ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)fd),((PetscObject)fd)->prefix,optionname,&viewer,&format,&flg);CHKERRQ(ierr);
325005c665bSBarry Smith   if (flg) {
326cffb1e40SBarry Smith     ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr);
3273050cee2SBarry Smith     ierr = MatFDColoringView(fd,viewer);CHKERRQ(ierr);
328cffb1e40SBarry Smith     ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
329cffb1e40SBarry Smith     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
330005c665bSBarry Smith   }
3313a40ed3dSBarry Smith   PetscFunctionReturn(0);
332bbf0e169SBarry Smith }
333bbf0e169SBarry Smith 
3344a2ae208SSatish Balay #undef __FUNCT__
3354a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringCreate"
33605869f15SSatish Balay /*@
337639f9d9dSBarry Smith    MatFDColoringCreate - Creates a matrix coloring context for finite difference
338639f9d9dSBarry Smith    computation of Jacobians.
339bbf0e169SBarry Smith 
340ef5ee4d1SLois Curfman McInnes    Collective on Mat
341ef5ee4d1SLois Curfman McInnes 
342639f9d9dSBarry Smith    Input Parameters:
343ef5ee4d1SLois Curfman McInnes +  mat - the matrix containing the nonzero structure of the Jacobian
344e727c939SJed Brown -  iscoloring - the coloring of the matrix; usually obtained with MatGetColoring() or DMCreateColoring()
345bbf0e169SBarry Smith 
346bbf0e169SBarry Smith     Output Parameter:
347639f9d9dSBarry Smith .   color - the new coloring context
348bbf0e169SBarry Smith 
34915091d37SBarry Smith     Level: intermediate
35015091d37SBarry Smith 
3516831982aSBarry Smith .seealso: MatFDColoringDestroy(),SNESDefaultComputeJacobianColor(), ISColoringCreate(),
352d1c39f55SBarry Smith           MatFDColoringSetFunction(), MatFDColoringSetFromOptions(), MatFDColoringApply(),
353e727c939SJed Brown           MatFDColoringView(), MatFDColoringSetParameters(), MatGetColoring(), DMCreateColoring()
354bbf0e169SBarry Smith @*/
3557087cfbeSBarry Smith PetscErrorCode  MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
356bbf0e169SBarry Smith {
357639f9d9dSBarry Smith   MatFDColoring  c;
358639f9d9dSBarry Smith   MPI_Comm       comm;
359dfbe8321SBarry Smith   PetscErrorCode ierr;
36038baddfdSBarry Smith   PetscInt       M,N;
361639f9d9dSBarry Smith 
3623a40ed3dSBarry Smith   PetscFunctionBegin;
363d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
364639f9d9dSBarry Smith   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
365*ce94432eSBarry Smith   if (M != N) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Only for square matrices");
366f881d145SBarry Smith   ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
36767c2884eSBarry Smith   ierr = PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_CLASSID,"MatFDColoring","Jacobian computation via finite differences with coloring","Mat",comm,MatFDColoringDestroy,MatFDColoringView);CHKERRQ(ierr);
368639f9d9dSBarry Smith 
369b8f8c88eSHong Zhang   c->ctype = iscoloring->ctype;
370b8f8c88eSHong Zhang 
371f830108cSBarry Smith   if (mat->ops->fdcoloringcreate) {
372f830108cSBarry Smith     ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);CHKERRQ(ierr);
373*ce94432eSBarry Smith   } else SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Code not yet written for matrix type %s",((PetscObject)mat)->type_name);
374639f9d9dSBarry Smith 
375f77a5242SKarl Rupp   ierr = MatGetVecs(mat,NULL,&c->w1);CHKERRQ(ierr);
376b8f8c88eSHong Zhang   ierr = PetscLogObjectParent(c,c->w1);CHKERRQ(ierr);
377b8f8c88eSHong Zhang   ierr = VecDuplicate(c->w1,&c->w2);CHKERRQ(ierr);
378b8f8c88eSHong Zhang   ierr = PetscLogObjectParent(c,c->w2);CHKERRQ(ierr);
379b8f8c88eSHong Zhang 
38077d8c4bbSBarry Smith   c->error_rel    = PETSC_SQRT_MACHINE_EPSILON;
38177d8c4bbSBarry Smith   c->umin         = 100.0*PETSC_SQRT_MACHINE_EPSILON;
38249b058dcSBarry Smith   c->currentcolor = -1;
383efb30889SBarry Smith   c->htype        = "wp";
3844e269d77SPeter Brune   c->fset         = PETSC_FALSE;
385639f9d9dSBarry Smith 
386639f9d9dSBarry Smith   *color = c;
3874e269d77SPeter Brune   ierr   = PetscObjectCompose((PetscObject)mat,"SNESMatFDColoring",(PetscObject)c);CHKERRQ(ierr);
388d5ba7fb7SMatthew Knepley   ierr   = PetscLogEventEnd(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
3893a40ed3dSBarry Smith   PetscFunctionReturn(0);
390bbf0e169SBarry Smith }
391bbf0e169SBarry Smith 
3924a2ae208SSatish Balay #undef __FUNCT__
3934a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringDestroy"
39405869f15SSatish Balay /*@
395639f9d9dSBarry Smith     MatFDColoringDestroy - Destroys a matrix coloring context that was created
396639f9d9dSBarry Smith     via MatFDColoringCreate().
397bbf0e169SBarry Smith 
398ef5ee4d1SLois Curfman McInnes     Collective on MatFDColoring
399ef5ee4d1SLois Curfman McInnes 
400b4fc646aSLois Curfman McInnes     Input Parameter:
401639f9d9dSBarry Smith .   c - coloring context
402bbf0e169SBarry Smith 
40315091d37SBarry Smith     Level: intermediate
40415091d37SBarry Smith 
405639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
406bbf0e169SBarry Smith @*/
4076bf464f9SBarry Smith PetscErrorCode  MatFDColoringDestroy(MatFDColoring *c)
408bbf0e169SBarry Smith {
4096849ba73SBarry Smith   PetscErrorCode ierr;
41038baddfdSBarry Smith   PetscInt       i;
411bbf0e169SBarry Smith 
4123a40ed3dSBarry Smith   PetscFunctionBegin;
4136bf464f9SBarry Smith   if (!*c) PetscFunctionReturn(0);
4146bf464f9SBarry Smith   if (--((PetscObject)(*c))->refct > 0) {*c = 0; PetscFunctionReturn(0);}
415d4bb536fSBarry Smith 
4166bf464f9SBarry Smith   for (i=0; i<(*c)->ncolors; i++) {
4176bf464f9SBarry Smith     ierr = PetscFree((*c)->columns[i]);CHKERRQ(ierr);
4186bf464f9SBarry Smith     ierr = PetscFree((*c)->rows[i]);CHKERRQ(ierr);
4196bf464f9SBarry Smith     ierr = PetscFree((*c)->columnsforrow[i]);CHKERRQ(ierr);
4206bf464f9SBarry Smith     if ((*c)->vscaleforrow) {ierr = PetscFree((*c)->vscaleforrow[i]);CHKERRQ(ierr);}
421bbf0e169SBarry Smith   }
4226bf464f9SBarry Smith   ierr = PetscFree((*c)->ncolumns);CHKERRQ(ierr);
4236bf464f9SBarry Smith   ierr = PetscFree((*c)->columns);CHKERRQ(ierr);
4246bf464f9SBarry Smith   ierr = PetscFree((*c)->nrows);CHKERRQ(ierr);
4256bf464f9SBarry Smith   ierr = PetscFree((*c)->rows);CHKERRQ(ierr);
4266bf464f9SBarry Smith   ierr = PetscFree((*c)->columnsforrow);CHKERRQ(ierr);
4276bf464f9SBarry Smith   ierr = PetscFree((*c)->vscaleforrow);CHKERRQ(ierr);
4286bf464f9SBarry Smith   ierr = VecDestroy(&(*c)->vscale);CHKERRQ(ierr);
4296bf464f9SBarry Smith   ierr = VecDestroy(&(*c)->w1);CHKERRQ(ierr);
4306bf464f9SBarry Smith   ierr = VecDestroy(&(*c)->w2);CHKERRQ(ierr);
4316bf464f9SBarry Smith   ierr = VecDestroy(&(*c)->w3);CHKERRQ(ierr);
432d38fa0fbSBarry Smith   ierr = PetscHeaderDestroy(c);CHKERRQ(ierr);
4333a40ed3dSBarry Smith   PetscFunctionReturn(0);
434bbf0e169SBarry Smith }
43543a90d84SBarry Smith 
4364a2ae208SSatish Balay #undef __FUNCT__
43749b058dcSBarry Smith #define __FUNCT__ "MatFDColoringGetPerturbedColumns"
43849b058dcSBarry Smith /*@C
43949b058dcSBarry Smith     MatFDColoringGetPerturbedColumns - Returns the indices of the columns that
44049b058dcSBarry Smith       that are currently being perturbed.
44149b058dcSBarry Smith 
44249b058dcSBarry Smith     Not Collective
44349b058dcSBarry Smith 
44449b058dcSBarry Smith     Input Parameters:
44549b058dcSBarry Smith .   coloring - coloring context created with MatFDColoringCreate()
44649b058dcSBarry Smith 
44749b058dcSBarry Smith     Output Parameters:
44849b058dcSBarry Smith +   n - the number of local columns being perturbed
44949b058dcSBarry Smith -   cols - the column indices, in global numbering
45049b058dcSBarry Smith 
45149b058dcSBarry Smith    Level: intermediate
45249b058dcSBarry Smith 
45349b058dcSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringApply()
45449b058dcSBarry Smith 
45549b058dcSBarry Smith .keywords: coloring, Jacobian, finite differences
45649b058dcSBarry Smith @*/
4577087cfbeSBarry Smith PetscErrorCode  MatFDColoringGetPerturbedColumns(MatFDColoring coloring,PetscInt *n,PetscInt *cols[])
45849b058dcSBarry Smith {
45949b058dcSBarry Smith   PetscFunctionBegin;
46049b058dcSBarry Smith   if (coloring->currentcolor >= 0) {
46149b058dcSBarry Smith     *n    = coloring->ncolumns[coloring->currentcolor];
46249b058dcSBarry Smith     *cols = coloring->columns[coloring->currentcolor];
46349b058dcSBarry Smith   } else {
46449b058dcSBarry Smith     *n = 0;
46549b058dcSBarry Smith   }
46649b058dcSBarry Smith   PetscFunctionReturn(0);
46749b058dcSBarry Smith }
46849b058dcSBarry Smith 
46949b058dcSBarry Smith #undef __FUNCT__
4704a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringApply"
47143a90d84SBarry Smith /*@
472e0907662SLois Curfman McInnes     MatFDColoringApply - Given a matrix for which a MatFDColoring context
473e0907662SLois Curfman McInnes     has been created, computes the Jacobian for a function via finite differences.
47443a90d84SBarry Smith 
475fee21e36SBarry Smith     Collective on MatFDColoring
476fee21e36SBarry Smith 
477ef5ee4d1SLois Curfman McInnes     Input Parameters:
478ef5ee4d1SLois Curfman McInnes +   mat - location to store Jacobian
479ef5ee4d1SLois Curfman McInnes .   coloring - coloring context created with MatFDColoringCreate()
480ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
4817850c7c0SBarry Smith -   sctx - context required by function, if this is being used with the SNES solver then it is SNES object, otherwise it is null
482ef5ee4d1SLois Curfman McInnes 
4838bba8e72SBarry Smith     Options Database Keys:
484ebd3b9afSBarry Smith +    -mat_fd_type - "wp" or "ds"  (see MATMFFD_WP or MATMFFD_DS)
485b9722508SLois Curfman McInnes .    -mat_fd_coloring_view - Activates basic viewing or coloring
486e350db43SBarry Smith .    -mat_fd_coloring_view draw - Activates drawing of coloring
487e350db43SBarry Smith -    -mat_fd_coloring_view ::ascii_info - Activates viewing of coloring info
4888bba8e72SBarry Smith 
48915091d37SBarry Smith     Level: intermediate
49098d79826SSatish Balay 
4917850c7c0SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringSetFunction()
49243a90d84SBarry Smith 
49343a90d84SBarry Smith .keywords: coloring, Jacobian, finite differences
49443a90d84SBarry Smith @*/
4957087cfbeSBarry Smith PetscErrorCode  MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
49643a90d84SBarry Smith {
4973acb8795SBarry Smith   PetscErrorCode ierr;
4983acb8795SBarry Smith 
4993acb8795SBarry Smith   PetscFunctionBegin;
5000700a824SBarry Smith   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
5010700a824SBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_CLASSID,2);
5020700a824SBarry Smith   PetscValidHeaderSpecific(x1,VEC_CLASSID,3);
503*ce94432eSBarry Smith   if (!coloring->f) SETERRQ(PetscObjectComm((PetscObject)J),PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetFunction()");
504e32f2f54SBarry Smith   if (!J->ops->fdcoloringapply) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not supported for this matrix type %s",((PetscObject)J)->type_name);
5053acb8795SBarry Smith   ierr = (*J->ops->fdcoloringapply)(J,coloring,x1,flag,sctx);CHKERRQ(ierr);
5063acb8795SBarry Smith   PetscFunctionReturn(0);
5073acb8795SBarry Smith }
5083acb8795SBarry Smith 
5093acb8795SBarry Smith #undef __FUNCT__
5103acb8795SBarry Smith #define __FUNCT__ "MatFDColoringApply_AIJ"
5117087cfbeSBarry Smith PetscErrorCode  MatFDColoringApply_AIJ(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
5123acb8795SBarry Smith {
5136849ba73SBarry Smith   PetscErrorCode (*f)(void*,Vec,Vec,void*) = (PetscErrorCode (*)(void*,Vec,Vec,void*))coloring->f;
5146849ba73SBarry Smith   PetscErrorCode ierr;
5154e269d77SPeter Brune   PetscInt       k,start,end,l,row,col,srow,**vscaleforrow;
516efb30889SBarry Smith   PetscScalar    dx,*y,*xx,*w3_array;
517ea709b57SSatish Balay   PetscScalar    *vscale_array;
518efb30889SBarry Smith   PetscReal      epsilon = coloring->error_rel,umin = coloring->umin,unorm;
519b8f8c88eSHong Zhang   Vec            w1      = coloring->w1,w2=coloring->w2,w3;
520005c665bSBarry Smith   void           *fctx   = coloring->fctx;
521ace3abfcSBarry Smith   PetscBool      flg     = PETSC_FALSE;
522705d48f7SSatish Balay   PetscInt       ctype   = coloring->ctype,N,col_start=0,col_end=0;
523b8f8c88eSHong Zhang   Vec            x1_tmp;
524a2e34c3dSBarry Smith 
5253a40ed3dSBarry Smith   PetscFunctionBegin;
5260700a824SBarry Smith   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
5270700a824SBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_CLASSID,2);
5280700a824SBarry Smith   PetscValidHeaderSpecific(x1,VEC_CLASSID,3);
529*ce94432eSBarry Smith   if (!f) SETERRQ(PetscObjectComm((PetscObject)J),PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetFunction()");
530e0907662SLois Curfman McInnes 
531d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
5324c49b128SBarry Smith   ierr = MatSetUnfactored(J);CHKERRQ(ierr);
533f77a5242SKarl Rupp   ierr = PetscOptionsGetBool(NULL,"-mat_fd_coloring_dont_rezero",&flg,NULL);CHKERRQ(ierr);
534f1af5d2fSBarry Smith   if (flg) {
535ae15b995SBarry Smith     ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr);
536e0907662SLois Curfman McInnes   } else {
537ace3abfcSBarry Smith     PetscBool assembled;
5380b9b6f31SBarry Smith     ierr = MatAssembled(J,&assembled);CHKERRQ(ierr);
5390b9b6f31SBarry Smith     if (assembled) {
54043a90d84SBarry Smith       ierr = MatZeroEntries(J);CHKERRQ(ierr);
541e0907662SLois Curfman McInnes     }
5420b9b6f31SBarry Smith   }
54343a90d84SBarry Smith 
544b8f8c88eSHong Zhang   x1_tmp = x1;
545dac7f5fdSBarry Smith   if (!coloring->vscale) {
546b8f8c88eSHong Zhang     ierr = VecDuplicate(x1_tmp,&coloring->vscale);CHKERRQ(ierr);
547b8f8c88eSHong Zhang   }
548be2fbe1fSBarry Smith 
549b8f8c88eSHong Zhang   if (coloring->htype[0] == 'w') { /* tacky test; need to make systematic if we add other approaches to computing h*/
550b8f8c88eSHong Zhang     ierr = VecNorm(x1_tmp,NORM_2,&unorm);CHKERRQ(ierr);
551b8f8c88eSHong Zhang   }
552b8f8c88eSHong Zhang   ierr = VecGetOwnershipRange(w1,&start,&end);CHKERRQ(ierr); /* OwnershipRange is used by ghosted x! */
553b8f8c88eSHong Zhang 
554b8f8c88eSHong Zhang   /* Set w1 = F(x1) */
5554e269d77SPeter Brune   if (!coloring->fset) {
55666f9b7ceSBarry Smith     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
557b8f8c88eSHong Zhang     ierr = (*f)(sctx,x1_tmp,w1,fctx);CHKERRQ(ierr);
55866f9b7ceSBarry Smith     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
5594e269d77SPeter Brune   } else {
5604e269d77SPeter Brune     coloring->fset = PETSC_FALSE;
561be2fbe1fSBarry Smith   }
56243a90d84SBarry Smith 
563b8f8c88eSHong Zhang   if (!coloring->w3) {
564b8f8c88eSHong Zhang     ierr = VecDuplicate(x1_tmp,&coloring->w3);CHKERRQ(ierr);
565b8f8c88eSHong Zhang     ierr = PetscLogObjectParent(coloring,coloring->w3);CHKERRQ(ierr);
566efb30889SBarry Smith   }
567b8f8c88eSHong Zhang   w3 = coloring->w3;
568efb30889SBarry Smith 
569b8f8c88eSHong Zhang   /* Compute all the local scale factors, including ghost points */
570b8f8c88eSHong Zhang   ierr = VecGetLocalSize(x1_tmp,&N);CHKERRQ(ierr);
571b8f8c88eSHong Zhang   ierr = VecGetArray(x1_tmp,&xx);CHKERRQ(ierr);
572b8f8c88eSHong Zhang   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
573b8f8c88eSHong Zhang   if (ctype == IS_COLORING_GHOSTED) {
574b8f8c88eSHong Zhang     col_start = 0; col_end = N;
5758ee2e534SBarry Smith   } else if (ctype == IS_COLORING_GLOBAL) {
576b8f8c88eSHong Zhang     xx           = xx - start;
577b8f8c88eSHong Zhang     vscale_array = vscale_array - start;
578b8f8c88eSHong Zhang     col_start    = start; col_end = N + start;
579b8f8c88eSHong Zhang   }
580b8f8c88eSHong Zhang   for (col=col_start; col<col_end; col++) {
581b8f8c88eSHong Zhang     /* Loop over each local column, vscale[col] = 1./(epsilon*dx[col]) */
582efb30889SBarry Smith     if (coloring->htype[0] == 'w') {
583efb30889SBarry Smith       dx = 1.0 + unorm;
584efb30889SBarry Smith     } else {
5859782cf98SBarry Smith       dx = xx[col];
586efb30889SBarry Smith     }
587d4a378daSJed Brown     if (dx == (PetscScalar)0.0) dx = 1.0;
5889782cf98SBarry Smith     if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
5899782cf98SBarry Smith     else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
5909782cf98SBarry Smith     dx               *= epsilon;
591d4a378daSJed Brown     vscale_array[col] = (PetscScalar)1.0/dx;
5929782cf98SBarry Smith   }
5938ee2e534SBarry Smith   if (ctype == IS_COLORING_GLOBAL) vscale_array = vscale_array + start;
594efb30889SBarry Smith   ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
5958ee2e534SBarry Smith   if (ctype == IS_COLORING_GLOBAL) {
59630b34957SBarry Smith     ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
59730b34957SBarry Smith     ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
598b8f8c88eSHong Zhang   }
5999782cf98SBarry Smith 
600b8f8c88eSHong Zhang   if (coloring->vscaleforrow) {
601b8f8c88eSHong Zhang     vscaleforrow = coloring->vscaleforrow;
602*ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)J),PETSC_ERR_ARG_NULL,"Null Object: coloring->vscaleforrow");
603b0a32e0cSBarry Smith 
6049782cf98SBarry Smith   /*
60543a90d84SBarry Smith     Loop over each color
60643a90d84SBarry Smith   */
607b8f8c88eSHong Zhang   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
60843a90d84SBarry Smith   for (k=0; k<coloring->ncolors; k++) {
60949b058dcSBarry Smith     coloring->currentcolor = k;
61026fbe8dcSKarl Rupp 
611b8f8c88eSHong Zhang     ierr = VecCopy(x1_tmp,w3);CHKERRQ(ierr);
612b8f8c88eSHong Zhang     ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);
6138ee2e534SBarry Smith     if (ctype == IS_COLORING_GLOBAL) w3_array = w3_array - start;
61443a90d84SBarry Smith     /*
615b8f8c88eSHong Zhang       Loop over each column associated with color
616b8f8c88eSHong Zhang       adding the perturbation to the vector w3.
61743a90d84SBarry Smith     */
61843a90d84SBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
619b8f8c88eSHong Zhang       col = coloring->columns[k][l];    /* local column of the matrix we are probing for */
620efb30889SBarry Smith       if (coloring->htype[0] == 'w') {
621efb30889SBarry Smith         dx = 1.0 + unorm;
622efb30889SBarry Smith       } else {
62342460c72SBarry Smith         dx = xx[col];
624efb30889SBarry Smith       }
625d4a378daSJed Brown       if (dx == (PetscScalar)0.0) dx = 1.0;
626329f5518SBarry Smith       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
627329f5518SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
62843a90d84SBarry Smith       dx            *= epsilon;
629e32f2f54SBarry Smith       if (!PetscAbsScalar(dx)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Computed 0 differencing parameter");
63042460c72SBarry Smith       w3_array[col] += dx;
63143a90d84SBarry Smith     }
63298d79826SSatish Balay     if (ctype == IS_COLORING_GLOBAL) w3_array = w3_array + start;
633b8f8c88eSHong Zhang     ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
6343b28642cSBarry Smith 
63543a90d84SBarry Smith     /*
636b8f8c88eSHong Zhang       Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations)
637b8f8c88eSHong Zhang                            w2 = F(x1 + dx) - F(x1)
63843a90d84SBarry Smith     */
63966f9b7ceSBarry Smith     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
64043a90d84SBarry Smith     ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr);
64166f9b7ceSBarry Smith     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
642efb30889SBarry Smith     ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr);
6439782cf98SBarry Smith 
64443a90d84SBarry Smith     /*
645e0907662SLois Curfman McInnes       Loop over rows of vector, putting results into Jacobian matrix
64643a90d84SBarry Smith     */
6473b28642cSBarry Smith     ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
64843a90d84SBarry Smith     for (l=0; l<coloring->nrows[k]; l++) {
649b8f8c88eSHong Zhang       row     = coloring->rows[k][l];            /* local row index */
650b8f8c88eSHong Zhang       col     = coloring->columnsforrow[k][l];   /* global column index */
6515904e1b1SBarry Smith       y[row] *= vscale_array[vscaleforrow[k][l]];
65243a90d84SBarry Smith       srow    = row + start;
65343a90d84SBarry Smith       ierr    = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
65443a90d84SBarry Smith     }
6553b28642cSBarry Smith     ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
656aeb7e98aSMatthew Knepley   } /* endof for each color */
6578ee2e534SBarry Smith   if (ctype == IS_COLORING_GLOBAL) xx = xx + start;
65830b34957SBarry Smith   ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
659b8f8c88eSHong Zhang   ierr = VecRestoreArray(x1_tmp,&xx);CHKERRQ(ierr);
660b8f8c88eSHong Zhang 
661b8f8c88eSHong Zhang   coloring->currentcolor = -1;
6628865f1eaSKarl Rupp 
66343a90d84SBarry Smith   ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
66443a90d84SBarry Smith   ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
665d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
666a2e34c3dSBarry Smith 
667e350db43SBarry Smith   ierr = MatFDColoringViewFromOptions(coloring,"-mat_fd_coloring_view");CHKERRQ(ierr);
6683a40ed3dSBarry Smith   PetscFunctionReturn(0);
66943a90d84SBarry Smith }
670