xref: /petsc/src/mat/matfd/fdmatrix.c (revision f141eedf9782c96c0d2810308bf7683d831c0843)
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 
259804daf3SBarry Smith #include <petscdraw.h>
263a7fca6bSBarry Smith #undef __FUNCT__
274a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Draw_Zoom"
286849ba73SBarry Smith static PetscErrorCode MatFDColoringView_Draw_Zoom(PetscDraw draw,void *Aa)
299194eea9SBarry Smith {
309194eea9SBarry Smith   MatFDColoring  fd = (MatFDColoring)Aa;
31dfbe8321SBarry Smith   PetscErrorCode ierr;
3238baddfdSBarry Smith   PetscInt       i,j;
339194eea9SBarry Smith   PetscReal      x,y;
349194eea9SBarry Smith 
359194eea9SBarry Smith   PetscFunctionBegin;
369194eea9SBarry Smith   /* loop over colors  */
379194eea9SBarry Smith   for (i=0; i<fd->ncolors; i++) {
389194eea9SBarry Smith     for (j=0; j<fd->nrows[i]; j++) {
399194eea9SBarry Smith       y    = fd->M - fd->rows[i][j] - fd->rstart;
409194eea9SBarry Smith       x    = fd->columnsforrow[i][j];
41b0a32e0cSBarry Smith       ierr = PetscDrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);CHKERRQ(ierr);
429194eea9SBarry Smith     }
439194eea9SBarry Smith   }
449194eea9SBarry Smith   PetscFunctionReturn(0);
459194eea9SBarry Smith }
469194eea9SBarry Smith 
474a2ae208SSatish Balay #undef __FUNCT__
484a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Draw"
496849ba73SBarry Smith static PetscErrorCode MatFDColoringView_Draw(MatFDColoring fd,PetscViewer viewer)
50005c665bSBarry Smith {
51dfbe8321SBarry Smith   PetscErrorCode ierr;
52ace3abfcSBarry Smith   PetscBool      isnull;
53b0a32e0cSBarry Smith   PetscDraw      draw;
5454d96fa3SBarry Smith   PetscReal      xr,yr,xl,yl,h,w;
55005c665bSBarry Smith 
563a40ed3dSBarry Smith   PetscFunctionBegin;
57b0a32e0cSBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
58b0a32e0cSBarry Smith   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
599194eea9SBarry Smith 
609194eea9SBarry Smith   ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
61005c665bSBarry Smith 
62005c665bSBarry Smith   xr   = fd->N; yr  = fd->M; h = yr/10.0; w = xr/10.0;
63005c665bSBarry Smith   xr  += w;     yr += h;    xl = -w;     yl = -h;
64b0a32e0cSBarry Smith   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
65b0a32e0cSBarry Smith   ierr = PetscDrawZoom(draw,MatFDColoringView_Draw_Zoom,fd);CHKERRQ(ierr);
66f77a5242SKarl Rupp   ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",NULL);CHKERRQ(ierr);
673a40ed3dSBarry Smith   PetscFunctionReturn(0);
68005c665bSBarry Smith }
69005c665bSBarry Smith 
704a2ae208SSatish Balay #undef __FUNCT__
714a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView"
72bbf0e169SBarry Smith /*@C
73639f9d9dSBarry Smith    MatFDColoringView - Views a finite difference coloring context.
74bbf0e169SBarry Smith 
754c49b128SBarry Smith    Collective on MatFDColoring
76fee21e36SBarry Smith 
77ef5ee4d1SLois Curfman McInnes    Input  Parameters:
78ef5ee4d1SLois Curfman McInnes +  c - the coloring context
79ef5ee4d1SLois Curfman McInnes -  viewer - visualization context
80ef5ee4d1SLois Curfman McInnes 
8115091d37SBarry Smith    Level: intermediate
8215091d37SBarry Smith 
83b4fc646aSLois Curfman McInnes    Notes:
84b4fc646aSLois Curfman McInnes    The available visualization contexts include
85b0a32e0cSBarry Smith +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
86b0a32e0cSBarry Smith .     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
87ef5ee4d1SLois Curfman McInnes         output where only the first processor opens
88ef5ee4d1SLois Curfman McInnes         the file.  All other processors send their
89ef5ee4d1SLois Curfman McInnes         data to the first processor to print.
90b0a32e0cSBarry Smith -     PETSC_VIEWER_DRAW_WORLD - graphical display of nonzero structure
91bbf0e169SBarry Smith 
929194eea9SBarry Smith    Notes:
939194eea9SBarry Smith      Since PETSc uses only a small number of basic colors (currently 33), if the coloring
94b0a32e0cSBarry Smith    involves more than 33 then some seemingly identical colors are displayed making it look
959194eea9SBarry Smith    like an illegal coloring. This is just a graphical artifact.
969194eea9SBarry Smith 
97639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
98005c665bSBarry Smith 
99b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, view
100bbf0e169SBarry Smith @*/
1017087cfbeSBarry Smith PetscErrorCode  MatFDColoringView(MatFDColoring c,PetscViewer viewer)
102bbf0e169SBarry Smith {
1036849ba73SBarry Smith   PetscErrorCode    ierr;
10438baddfdSBarry Smith   PetscInt          i,j;
105ace3abfcSBarry Smith   PetscBool         isdraw,iascii;
106f3ef73ceSBarry Smith   PetscViewerFormat format;
107bbf0e169SBarry Smith 
1083a40ed3dSBarry Smith   PetscFunctionBegin;
1090700a824SBarry Smith   PetscValidHeaderSpecific(c,MAT_FDCOLORING_CLASSID,1);
1103050cee2SBarry Smith   if (!viewer) {
111ce94432eSBarry Smith     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)c),&viewer);CHKERRQ(ierr);
1123050cee2SBarry Smith   }
1130700a824SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
114c9780b6fSBarry Smith   PetscCheckSameComm(c,1,viewer,2);
115bbf0e169SBarry Smith 
116251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
117251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1180f5bd95cSBarry Smith   if (isdraw) {
119b4fc646aSLois Curfman McInnes     ierr = MatFDColoringView_Draw(c,viewer);CHKERRQ(ierr);
12032077d6dSBarry Smith   } else if (iascii) {
121dae58748SBarry Smith     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)c,viewer);CHKERRQ(ierr);
122a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Error tolerance=%G\n",c->error_rel);CHKERRQ(ierr);
123a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Umin=%G\n",c->umin);CHKERRQ(ierr);
12477431f27SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Number of colors=%D\n",c->ncolors);CHKERRQ(ierr);
125ae09f205SBarry Smith 
126b0a32e0cSBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
127fb9695e5SSatish Balay     if (format != PETSC_VIEWER_ASCII_INFO) {
128b4fc646aSLois Curfman McInnes       for (i=0; i<c->ncolors; i++) {
12977431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"  Information for color %D\n",i);CHKERRQ(ierr);
13077431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    Number of columns %D\n",c->ncolumns[i]);CHKERRQ(ierr);
131b4fc646aSLois Curfman McInnes         for (j=0; j<c->ncolumns[i]; j++) {
13277431f27SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"      %D\n",c->columns[i][j]);CHKERRQ(ierr);
133639f9d9dSBarry Smith         }
13477431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    Number of rows %D\n",c->nrows[i]);CHKERRQ(ierr);
135b4fc646aSLois Curfman McInnes         for (j=0; j<c->nrows[i]; j++) {
13677431f27SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"      %D %D \n",c->rows[i][j],c->columnsforrow[i][j]);CHKERRQ(ierr);
137b4fc646aSLois Curfman McInnes         }
138bbf0e169SBarry Smith       }
139bbf0e169SBarry Smith     }
140b0a32e0cSBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
141bbf0e169SBarry Smith   }
1423a40ed3dSBarry Smith   PetscFunctionReturn(0);
143639f9d9dSBarry Smith }
144639f9d9dSBarry Smith 
1454a2ae208SSatish Balay #undef __FUNCT__
1464a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetParameters"
147639f9d9dSBarry Smith /*@
148b4fc646aSLois Curfman McInnes    MatFDColoringSetParameters - Sets the parameters for the sparse approximation of
149b4fc646aSLois Curfman McInnes    a Jacobian matrix using finite differences.
150639f9d9dSBarry Smith 
1513f9fe445SBarry Smith    Logically Collective on MatFDColoring
152ef5ee4d1SLois Curfman McInnes 
153ef5ee4d1SLois Curfman McInnes    The Jacobian is estimated with the differencing approximation
154ef5ee4d1SLois Curfman McInnes .vb
15565f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
156f23b5b22SLois Curfman McInnes        h = error_rel*u[i]                 if  abs(u[i]) > umin
157f23b5b22SLois Curfman McInnes          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
158ef5ee4d1SLois Curfman McInnes        dx_{i} = (0, ... 1, .... 0)
159ef5ee4d1SLois Curfman McInnes .ve
160639f9d9dSBarry Smith 
161639f9d9dSBarry Smith    Input Parameters:
162ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
163639f9d9dSBarry Smith .  error_rel - relative error
164f23b5b22SLois Curfman McInnes -  umin - minimum allowable u-value magnitude
165fee21e36SBarry Smith 
16615091d37SBarry Smith    Level: advanced
16715091d37SBarry Smith 
168b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, parameters
169b4fc646aSLois Curfman McInnes 
17095b89fc3SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringSetFromOptions()
17195b89fc3SBarry Smith 
172639f9d9dSBarry Smith @*/
1737087cfbeSBarry Smith PetscErrorCode  MatFDColoringSetParameters(MatFDColoring matfd,PetscReal error,PetscReal umin)
174639f9d9dSBarry Smith {
1753a40ed3dSBarry Smith   PetscFunctionBegin;
1760700a824SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
177c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(matfd,error,2);
178c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(matfd,umin,3);
179639f9d9dSBarry Smith   if (error != PETSC_DEFAULT) matfd->error_rel = error;
180639f9d9dSBarry Smith   if (umin != PETSC_DEFAULT)  matfd->umin      = umin;
1813a40ed3dSBarry Smith   PetscFunctionReturn(0);
182639f9d9dSBarry Smith }
183639f9d9dSBarry Smith 
184005c665bSBarry Smith 
185ff0cfa39SBarry Smith 
1864a2ae208SSatish Balay #undef __FUNCT__
1874a9d489dSBarry Smith #define __FUNCT__ "MatFDColoringGetFunction"
1884a9d489dSBarry Smith /*@C
1894a9d489dSBarry Smith    MatFDColoringGetFunction - Gets the function to use for computing the Jacobian.
1904a9d489dSBarry Smith 
1913f9fe445SBarry Smith    Not Collective
1924a9d489dSBarry Smith 
1934a9d489dSBarry Smith    Input Parameters:
1944a9d489dSBarry Smith .  coloring - the coloring context
1954a9d489dSBarry Smith 
1964a9d489dSBarry Smith    Output Parameters:
1974a9d489dSBarry Smith +  f - the function
1984a9d489dSBarry Smith -  fctx - the optional user-defined function context
1994a9d489dSBarry Smith 
2004a9d489dSBarry Smith    Level: intermediate
2014a9d489dSBarry Smith 
2024a9d489dSBarry Smith .keywords: Mat, Jacobian, finite differences, set, function
20395b89fc3SBarry Smith 
20495b89fc3SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringSetFunction(), MatFDColoringSetFromOptions()
20595b89fc3SBarry Smith 
2064a9d489dSBarry Smith @*/
2077087cfbeSBarry Smith PetscErrorCode  MatFDColoringGetFunction(MatFDColoring matfd,PetscErrorCode (**f)(void),void **fctx)
2084a9d489dSBarry Smith {
2094a9d489dSBarry Smith   PetscFunctionBegin;
2100700a824SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
2114a9d489dSBarry Smith   if (f) *f = matfd->f;
2124a9d489dSBarry Smith   if (fctx) *fctx = matfd->fctx;
2134a9d489dSBarry Smith   PetscFunctionReturn(0);
2144a9d489dSBarry Smith }
2154a9d489dSBarry Smith 
2164a9d489dSBarry Smith #undef __FUNCT__
2174a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFunction"
218d64ed03dSBarry Smith /*@C
219005c665bSBarry Smith    MatFDColoringSetFunction - Sets the function to use for computing the Jacobian.
220005c665bSBarry Smith 
2213f9fe445SBarry Smith    Logically Collective on MatFDColoring
222fee21e36SBarry Smith 
223ef5ee4d1SLois Curfman McInnes    Input Parameters:
224ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
225ef5ee4d1SLois Curfman McInnes .  f - the function
226ef5ee4d1SLois Curfman McInnes -  fctx - the optional user-defined function context
227ef5ee4d1SLois Curfman McInnes 
2287850c7c0SBarry Smith    Calling sequence of (*f) function:
2297850c7c0SBarry Smith     For SNES:    PetscErrorCode (*f)(SNES,Vec,Vec,void*)
230ab637aeaSJed Brown     If not using SNES: PetscErrorCode (*f)(void *dummy,Vec,Vec,void*) and dummy is ignored
23115091d37SBarry Smith 
2327850c7c0SBarry Smith    Level: advanced
2337850c7c0SBarry Smith 
234ab637aeaSJed Brown    Notes: This function is usually used automatically by SNES (when one uses SNESSetJacobian() with the argument
2358d359177SBarry Smith      SNESComputeJacobianDefaultColor()) and only needs to be used by someone computing a matrix via coloring directly by
236ab637aeaSJed Brown      calling MatFDColoringApply()
2377850c7c0SBarry Smith 
2387850c7c0SBarry Smith    Fortran Notes:
2397850c7c0SBarry Smith     In Fortran you must call MatFDColoringSetFunction() for a coloring object to
240ab637aeaSJed Brown   be used without SNES or within the SNES solvers.
241f881d145SBarry Smith 
242b4fc646aSLois Curfman McInnes .keywords: Mat, Jacobian, finite differences, set, function
24395b89fc3SBarry Smith 
24495b89fc3SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringGetFunction(), MatFDColoringSetFromOptions()
24595b89fc3SBarry Smith 
246005c665bSBarry Smith @*/
2477087cfbeSBarry Smith PetscErrorCode  MatFDColoringSetFunction(MatFDColoring matfd,PetscErrorCode (*f)(void),void *fctx)
248005c665bSBarry Smith {
2493a40ed3dSBarry Smith   PetscFunctionBegin;
2500700a824SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
251005c665bSBarry Smith   matfd->f    = f;
252005c665bSBarry Smith   matfd->fctx = fctx;
2533a40ed3dSBarry Smith   PetscFunctionReturn(0);
254005c665bSBarry Smith }
255005c665bSBarry Smith 
2564a2ae208SSatish Balay #undef __FUNCT__
2574a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFromOptions"
258639f9d9dSBarry Smith /*@
259b4fc646aSLois Curfman McInnes    MatFDColoringSetFromOptions - Sets coloring finite difference parameters from
260639f9d9dSBarry Smith    the options database.
261639f9d9dSBarry Smith 
262fee21e36SBarry Smith    Collective on MatFDColoring
263fee21e36SBarry Smith 
26465f2ba5bSLois Curfman McInnes    The Jacobian, F'(u), is estimated with the differencing approximation
265ef5ee4d1SLois Curfman McInnes .vb
26665f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
267f23b5b22SLois Curfman McInnes        h = error_rel*u[i]                 if  abs(u[i]) > umin
268f23b5b22SLois Curfman McInnes          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
269ef5ee4d1SLois Curfman McInnes        dx_{i} = (0, ... 1, .... 0)
270ef5ee4d1SLois Curfman McInnes .ve
271ef5ee4d1SLois Curfman McInnes 
272ef5ee4d1SLois Curfman McInnes    Input Parameter:
273ef5ee4d1SLois Curfman McInnes .  coloring - the coloring context
274ef5ee4d1SLois Curfman McInnes 
275b4fc646aSLois Curfman McInnes    Options Database Keys:
276d15ffeeaSSatish Balay +  -mat_fd_coloring_err <err> - Sets <err> (square root
277ef5ee4d1SLois Curfman McInnes            of relative error in the function)
278f23b5b22SLois Curfman McInnes .  -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude
2793ec795f1SBarry Smith .  -mat_fd_type - "wp" or "ds" (see MATMFFD_WP or MATMFFD_DS)
280ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_view - Activates basic viewing
281e350db43SBarry Smith .  -mat_fd_coloring_view ::ascii_info - Activates viewing info
282e350db43SBarry Smith -  -mat_fd_coloring_view draw - Activates drawing
283639f9d9dSBarry Smith 
28415091d37SBarry Smith     Level: intermediate
28515091d37SBarry Smith 
286b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, parameters
287d1c39f55SBarry Smith 
288d1c39f55SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringView(), MatFDColoringSetParameters()
289d1c39f55SBarry Smith 
290639f9d9dSBarry Smith @*/
2917087cfbeSBarry Smith PetscErrorCode  MatFDColoringSetFromOptions(MatFDColoring matfd)
292639f9d9dSBarry Smith {
293dfbe8321SBarry Smith   PetscErrorCode ierr;
294ace3abfcSBarry Smith   PetscBool      flg;
295efb30889SBarry Smith   char           value[3];
2963a40ed3dSBarry Smith 
2973a40ed3dSBarry Smith   PetscFunctionBegin;
2980700a824SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
299639f9d9dSBarry Smith 
3003194b578SJed Brown   ierr = PetscObjectOptionsBegin((PetscObject)matfd);CHKERRQ(ierr);
30187828ca2SBarry Smith   ierr = PetscOptionsReal("-mat_fd_coloring_err","Square root of relative error in function","MatFDColoringSetParameters",matfd->error_rel,&matfd->error_rel,0);CHKERRQ(ierr);
30287828ca2SBarry Smith   ierr = PetscOptionsReal("-mat_fd_coloring_umin","Minimum allowable u magnitude","MatFDColoringSetParameters",matfd->umin,&matfd->umin,0);CHKERRQ(ierr);
3031bc75a8dSBarry Smith   ierr = PetscOptionsString("-mat_fd_type","Algorithm to compute h, wp or ds","MatFDColoringCreate",matfd->htype,value,3,&flg);CHKERRQ(ierr);
304efb30889SBarry Smith   if (flg) {
305efb30889SBarry Smith     if (value[0] == 'w' && value[1] == 'p') matfd->htype = "wp";
306efb30889SBarry Smith     else if (value[0] == 'd' && value[1] == 's') matfd->htype = "ds";
307e32f2f54SBarry Smith     else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Unknown finite differencing type %s",value);
308efb30889SBarry Smith   }
3095d973c19SBarry Smith   /* process any options handlers added with PetscObjectAddOptionsHandler() */
3105d973c19SBarry Smith   ierr = PetscObjectProcessOptionsHandlers((PetscObject)matfd);CHKERRQ(ierr);
311b0a32e0cSBarry Smith   PetscOptionsEnd();CHKERRQ(ierr);
3123a40ed3dSBarry Smith   PetscFunctionReturn(0);
313005c665bSBarry Smith }
314005c665bSBarry Smith 
3154a2ae208SSatish Balay #undef __FUNCT__
316e350db43SBarry Smith #define __FUNCT__ "MatFDColoringViewFromOptions"
317146574abSBarry Smith PetscErrorCode MatFDColoringViewFromOptions(MatFDColoring fd,const char prefix[],const char optionname[])
318005c665bSBarry Smith {
319dfbe8321SBarry Smith   PetscErrorCode    ierr;
320e350db43SBarry Smith   PetscBool         flg;
3213050cee2SBarry Smith   PetscViewer       viewer;
322cffb1e40SBarry Smith   PetscViewerFormat format;
323005c665bSBarry Smith 
3243a40ed3dSBarry Smith   PetscFunctionBegin;
325146574abSBarry Smith   if (prefix) {
326146574abSBarry Smith     ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)fd),prefix,optionname,&viewer,&format,&flg);CHKERRQ(ierr);
327146574abSBarry Smith   } else {
328ce94432eSBarry Smith     ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)fd),((PetscObject)fd)->prefix,optionname,&viewer,&format,&flg);CHKERRQ(ierr);
329146574abSBarry Smith   }
330005c665bSBarry Smith   if (flg) {
331cffb1e40SBarry Smith     ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr);
3323050cee2SBarry Smith     ierr = MatFDColoringView(fd,viewer);CHKERRQ(ierr);
333cffb1e40SBarry Smith     ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
334cffb1e40SBarry Smith     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
335005c665bSBarry Smith   }
3363a40ed3dSBarry Smith   PetscFunctionReturn(0);
337bbf0e169SBarry Smith }
338bbf0e169SBarry Smith 
3394a2ae208SSatish Balay #undef __FUNCT__
3404a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringCreate"
34105869f15SSatish Balay /*@
342639f9d9dSBarry Smith    MatFDColoringCreate - Creates a matrix coloring context for finite difference
343639f9d9dSBarry Smith    computation of Jacobians.
344bbf0e169SBarry Smith 
345ef5ee4d1SLois Curfman McInnes    Collective on Mat
346ef5ee4d1SLois Curfman McInnes 
347639f9d9dSBarry Smith    Input Parameters:
348ef5ee4d1SLois Curfman McInnes +  mat - the matrix containing the nonzero structure of the Jacobian
349e727c939SJed Brown -  iscoloring - the coloring of the matrix; usually obtained with MatGetColoring() or DMCreateColoring()
350bbf0e169SBarry Smith 
351bbf0e169SBarry Smith     Output Parameter:
352639f9d9dSBarry Smith .   color - the new coloring context
353bbf0e169SBarry Smith 
35415091d37SBarry Smith     Level: intermediate
35515091d37SBarry Smith 
3568d359177SBarry Smith .seealso: MatFDColoringDestroy(),SNESComputeJacobianDefaultColor(), ISColoringCreate(),
357d1c39f55SBarry Smith           MatFDColoringSetFunction(), MatFDColoringSetFromOptions(), MatFDColoringApply(),
358e727c939SJed Brown           MatFDColoringView(), MatFDColoringSetParameters(), MatGetColoring(), DMCreateColoring()
359bbf0e169SBarry Smith @*/
3607087cfbeSBarry Smith PetscErrorCode  MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
361bbf0e169SBarry Smith {
362639f9d9dSBarry Smith   MatFDColoring  c;
363639f9d9dSBarry Smith   MPI_Comm       comm;
364dfbe8321SBarry Smith   PetscErrorCode ierr;
36538baddfdSBarry Smith   PetscInt       M,N;
366639f9d9dSBarry Smith 
3673a40ed3dSBarry Smith   PetscFunctionBegin;
368*f141eedfSHong Zhang   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Matrix must be assembled by calls to MatAssemblyBegin/End();");
369d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
370639f9d9dSBarry Smith   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
371ce94432eSBarry Smith   if (M != N) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Only for square matrices");
372f881d145SBarry Smith   ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
37367c2884eSBarry Smith   ierr = PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_CLASSID,"MatFDColoring","Jacobian computation via finite differences with coloring","Mat",comm,MatFDColoringDestroy,MatFDColoringView);CHKERRQ(ierr);
374639f9d9dSBarry Smith 
375b8f8c88eSHong Zhang   c->ctype = iscoloring->ctype;
376b8f8c88eSHong Zhang 
377f830108cSBarry Smith   if (mat->ops->fdcoloringcreate) {
378f830108cSBarry Smith     ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);CHKERRQ(ierr);
379ce94432eSBarry Smith   } else SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Code not yet written for matrix type %s",((PetscObject)mat)->type_name);
380639f9d9dSBarry Smith 
381f77a5242SKarl Rupp   ierr = MatGetVecs(mat,NULL,&c->w1);CHKERRQ(ierr);
3823bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)c,(PetscObject)c->w1);CHKERRQ(ierr);
383b8f8c88eSHong Zhang   ierr = VecDuplicate(c->w1,&c->w2);CHKERRQ(ierr);
3843bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)c,(PetscObject)c->w2);CHKERRQ(ierr);
385b8f8c88eSHong Zhang 
38677d8c4bbSBarry Smith   c->error_rel    = PETSC_SQRT_MACHINE_EPSILON;
38777d8c4bbSBarry Smith   c->umin         = 100.0*PETSC_SQRT_MACHINE_EPSILON;
38849b058dcSBarry Smith   c->currentcolor = -1;
389efb30889SBarry Smith   c->htype        = "wp";
3904e269d77SPeter Brune   c->fset         = PETSC_FALSE;
391639f9d9dSBarry Smith 
392639f9d9dSBarry Smith   *color = c;
3934e269d77SPeter Brune   ierr   = PetscObjectCompose((PetscObject)mat,"SNESMatFDColoring",(PetscObject)c);CHKERRQ(ierr);
394d5ba7fb7SMatthew Knepley   ierr   = PetscLogEventEnd(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
3953a40ed3dSBarry Smith   PetscFunctionReturn(0);
396bbf0e169SBarry Smith }
397bbf0e169SBarry Smith 
3984a2ae208SSatish Balay #undef __FUNCT__
3994a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringDestroy"
40005869f15SSatish Balay /*@
401639f9d9dSBarry Smith     MatFDColoringDestroy - Destroys a matrix coloring context that was created
402639f9d9dSBarry Smith     via MatFDColoringCreate().
403bbf0e169SBarry Smith 
404ef5ee4d1SLois Curfman McInnes     Collective on MatFDColoring
405ef5ee4d1SLois Curfman McInnes 
406b4fc646aSLois Curfman McInnes     Input Parameter:
407639f9d9dSBarry Smith .   c - coloring context
408bbf0e169SBarry Smith 
40915091d37SBarry Smith     Level: intermediate
41015091d37SBarry Smith 
411639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
412bbf0e169SBarry Smith @*/
4136bf464f9SBarry Smith PetscErrorCode  MatFDColoringDestroy(MatFDColoring *c)
414bbf0e169SBarry Smith {
4156849ba73SBarry Smith   PetscErrorCode ierr;
41638baddfdSBarry Smith   PetscInt       i;
417bbf0e169SBarry Smith 
4183a40ed3dSBarry Smith   PetscFunctionBegin;
4196bf464f9SBarry Smith   if (!*c) PetscFunctionReturn(0);
4206bf464f9SBarry Smith   if (--((PetscObject)(*c))->refct > 0) {*c = 0; PetscFunctionReturn(0);}
421d4bb536fSBarry Smith 
4226bf464f9SBarry Smith   for (i=0; i<(*c)->ncolors; i++) {
4236bf464f9SBarry Smith     ierr = PetscFree((*c)->columns[i]);CHKERRQ(ierr);
4246bf464f9SBarry Smith     ierr = PetscFree((*c)->rows[i]);CHKERRQ(ierr);
4256bf464f9SBarry Smith     ierr = PetscFree((*c)->columnsforrow[i]);CHKERRQ(ierr);
4266bf464f9SBarry Smith     if ((*c)->vscaleforrow) {ierr = PetscFree((*c)->vscaleforrow[i]);CHKERRQ(ierr);}
427bbf0e169SBarry Smith   }
4286bf464f9SBarry Smith   ierr = PetscFree((*c)->ncolumns);CHKERRQ(ierr);
4296bf464f9SBarry Smith   ierr = PetscFree((*c)->columns);CHKERRQ(ierr);
4306bf464f9SBarry Smith   ierr = PetscFree((*c)->nrows);CHKERRQ(ierr);
4316bf464f9SBarry Smith   ierr = PetscFree((*c)->rows);CHKERRQ(ierr);
4326bf464f9SBarry Smith   ierr = PetscFree((*c)->columnsforrow);CHKERRQ(ierr);
4336bf464f9SBarry Smith   ierr = PetscFree((*c)->vscaleforrow);CHKERRQ(ierr);
43479c1e64dSHong Zhang   ierr = PetscFree((*c)->den2sp);CHKERRQ(ierr);
4356bf464f9SBarry Smith   ierr = VecDestroy(&(*c)->vscale);CHKERRQ(ierr);
4366bf464f9SBarry Smith   ierr = VecDestroy(&(*c)->w1);CHKERRQ(ierr);
4376bf464f9SBarry Smith   ierr = VecDestroy(&(*c)->w2);CHKERRQ(ierr);
4386bf464f9SBarry Smith   ierr = VecDestroy(&(*c)->w3);CHKERRQ(ierr);
439d38fa0fbSBarry Smith   ierr = PetscHeaderDestroy(c);CHKERRQ(ierr);
4403a40ed3dSBarry Smith   PetscFunctionReturn(0);
441bbf0e169SBarry Smith }
44243a90d84SBarry Smith 
4434a2ae208SSatish Balay #undef __FUNCT__
44449b058dcSBarry Smith #define __FUNCT__ "MatFDColoringGetPerturbedColumns"
44549b058dcSBarry Smith /*@C
44649b058dcSBarry Smith     MatFDColoringGetPerturbedColumns - Returns the indices of the columns that
44749b058dcSBarry Smith       that are currently being perturbed.
44849b058dcSBarry Smith 
44949b058dcSBarry Smith     Not Collective
45049b058dcSBarry Smith 
45149b058dcSBarry Smith     Input Parameters:
45249b058dcSBarry Smith .   coloring - coloring context created with MatFDColoringCreate()
45349b058dcSBarry Smith 
45449b058dcSBarry Smith     Output Parameters:
45549b058dcSBarry Smith +   n - the number of local columns being perturbed
45649b058dcSBarry Smith -   cols - the column indices, in global numbering
45749b058dcSBarry Smith 
45849b058dcSBarry Smith    Level: intermediate
45949b058dcSBarry Smith 
46049b058dcSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringApply()
46149b058dcSBarry Smith 
46249b058dcSBarry Smith .keywords: coloring, Jacobian, finite differences
46349b058dcSBarry Smith @*/
4647087cfbeSBarry Smith PetscErrorCode  MatFDColoringGetPerturbedColumns(MatFDColoring coloring,PetscInt *n,PetscInt *cols[])
46549b058dcSBarry Smith {
46649b058dcSBarry Smith   PetscFunctionBegin;
46749b058dcSBarry Smith   if (coloring->currentcolor >= 0) {
46849b058dcSBarry Smith     *n    = coloring->ncolumns[coloring->currentcolor];
46949b058dcSBarry Smith     *cols = coloring->columns[coloring->currentcolor];
47049b058dcSBarry Smith   } else {
47149b058dcSBarry Smith     *n = 0;
47249b058dcSBarry Smith   }
47349b058dcSBarry Smith   PetscFunctionReturn(0);
47449b058dcSBarry Smith }
47549b058dcSBarry Smith 
47649b058dcSBarry Smith #undef __FUNCT__
4774a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringApply"
47843a90d84SBarry Smith /*@
479e0907662SLois Curfman McInnes     MatFDColoringApply - Given a matrix for which a MatFDColoring context
480e0907662SLois Curfman McInnes     has been created, computes the Jacobian for a function via finite differences.
48143a90d84SBarry Smith 
482fee21e36SBarry Smith     Collective on MatFDColoring
483fee21e36SBarry Smith 
484ef5ee4d1SLois Curfman McInnes     Input Parameters:
485ef5ee4d1SLois Curfman McInnes +   mat - location to store Jacobian
486ef5ee4d1SLois Curfman McInnes .   coloring - coloring context created with MatFDColoringCreate()
487ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
4887850c7c0SBarry Smith -   sctx - context required by function, if this is being used with the SNES solver then it is SNES object, otherwise it is null
489ef5ee4d1SLois Curfman McInnes 
4908bba8e72SBarry Smith     Options Database Keys:
491ebd3b9afSBarry Smith +    -mat_fd_type - "wp" or "ds"  (see MATMFFD_WP or MATMFFD_DS)
492b9722508SLois Curfman McInnes .    -mat_fd_coloring_view - Activates basic viewing or coloring
493e350db43SBarry Smith .    -mat_fd_coloring_view draw - Activates drawing of coloring
494e350db43SBarry Smith -    -mat_fd_coloring_view ::ascii_info - Activates viewing of coloring info
4958bba8e72SBarry Smith 
49615091d37SBarry Smith     Level: intermediate
49798d79826SSatish Balay 
4987850c7c0SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringSetFunction()
49943a90d84SBarry Smith 
50043a90d84SBarry Smith .keywords: coloring, Jacobian, finite differences
50143a90d84SBarry Smith @*/
5027087cfbeSBarry Smith PetscErrorCode  MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
50343a90d84SBarry Smith {
5043acb8795SBarry Smith   PetscErrorCode ierr;
5053acb8795SBarry Smith 
5063acb8795SBarry Smith   PetscFunctionBegin;
5070700a824SBarry Smith   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
5080700a824SBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_CLASSID,2);
5090700a824SBarry Smith   PetscValidHeaderSpecific(x1,VEC_CLASSID,3);
510ce94432eSBarry Smith   if (!coloring->f) SETERRQ(PetscObjectComm((PetscObject)J),PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetFunction()");
511e32f2f54SBarry Smith   if (!J->ops->fdcoloringapply) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not supported for this matrix type %s",((PetscObject)J)->type_name);
5123acb8795SBarry Smith   ierr = (*J->ops->fdcoloringapply)(J,coloring,x1,flag,sctx);CHKERRQ(ierr);
5133acb8795SBarry Smith   PetscFunctionReturn(0);
5143acb8795SBarry Smith }
5153acb8795SBarry Smith 
5169a3c1acfSHong Zhang /* #define JACOBIANCOLOROPT */
5179a3c1acfSHong Zhang #if defined(JACOBIANCOLOROPT)
5189a3c1acfSHong Zhang #include <petsctime.h>
5199a3c1acfSHong Zhang #endif
5203acb8795SBarry Smith #undef __FUNCT__
5213acb8795SBarry Smith #define __FUNCT__ "MatFDColoringApply_AIJ"
5227087cfbeSBarry Smith PetscErrorCode  MatFDColoringApply_AIJ(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
5233acb8795SBarry Smith {
5246849ba73SBarry Smith   PetscErrorCode (*f)(void*,Vec,Vec,void*) = (PetscErrorCode (*)(void*,Vec,Vec,void*))coloring->f;
5256849ba73SBarry Smith   PetscErrorCode ierr;
5264e269d77SPeter Brune   PetscInt       k,start,end,l,row,col,srow,**vscaleforrow;
527efb30889SBarry Smith   PetscScalar    dx,*y,*xx,*w3_array;
528ea709b57SSatish Balay   PetscScalar    *vscale_array;
529efb30889SBarry Smith   PetscReal      epsilon = coloring->error_rel,umin = coloring->umin,unorm;
530b8f8c88eSHong Zhang   Vec            w1      = coloring->w1,w2=coloring->w2,w3;
531005c665bSBarry Smith   void           *fctx   = coloring->fctx;
532ace3abfcSBarry Smith   PetscBool      flg     = PETSC_FALSE;
533705d48f7SSatish Balay   PetscInt       ctype   = coloring->ctype,N,col_start=0,col_end=0;
534b8f8c88eSHong Zhang   Vec            x1_tmp;
5359a3c1acfSHong Zhang #if defined(JACOBIANCOLOROPT)
5369a3c1acfSHong Zhang   PetscLogDouble t0,t1,time_setvalues=0.0;
5379a3c1acfSHong Zhang #endif
538a2e34c3dSBarry Smith 
5393a40ed3dSBarry Smith   PetscFunctionBegin;
5400700a824SBarry Smith   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
5410700a824SBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_CLASSID,2);
5420700a824SBarry Smith   PetscValidHeaderSpecific(x1,VEC_CLASSID,3);
543ce94432eSBarry Smith   if (!f) SETERRQ(PetscObjectComm((PetscObject)J),PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetFunction()");
544e0907662SLois Curfman McInnes 
545d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
5464c49b128SBarry Smith   ierr = MatSetUnfactored(J);CHKERRQ(ierr);
547f77a5242SKarl Rupp   ierr = PetscOptionsGetBool(NULL,"-mat_fd_coloring_dont_rezero",&flg,NULL);CHKERRQ(ierr);
548f1af5d2fSBarry Smith   if (flg) {
549ae15b995SBarry Smith     ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr);
550e0907662SLois Curfman McInnes   } else {
551ace3abfcSBarry Smith     PetscBool assembled;
5520b9b6f31SBarry Smith     ierr = MatAssembled(J,&assembled);CHKERRQ(ierr);
5530b9b6f31SBarry Smith     if (assembled) {
55443a90d84SBarry Smith       ierr = MatZeroEntries(J);CHKERRQ(ierr);
555e0907662SLois Curfman McInnes     }
5560b9b6f31SBarry Smith   }
55743a90d84SBarry Smith 
558b8f8c88eSHong Zhang   x1_tmp = x1;
559dac7f5fdSBarry Smith   if (!coloring->vscale) {
560b8f8c88eSHong Zhang     ierr = VecDuplicate(x1_tmp,&coloring->vscale);CHKERRQ(ierr);
561b8f8c88eSHong Zhang   }
562be2fbe1fSBarry Smith 
563b8f8c88eSHong Zhang   if (coloring->htype[0] == 'w') { /* tacky test; need to make systematic if we add other approaches to computing h*/
564b8f8c88eSHong Zhang     ierr = VecNorm(x1_tmp,NORM_2,&unorm);CHKERRQ(ierr);
565b8f8c88eSHong Zhang   }
566b8f8c88eSHong Zhang   ierr = VecGetOwnershipRange(w1,&start,&end);CHKERRQ(ierr); /* OwnershipRange is used by ghosted x! */
567b8f8c88eSHong Zhang 
568b8f8c88eSHong Zhang   /* Set w1 = F(x1) */
5694e269d77SPeter Brune   if (!coloring->fset) {
57066f9b7ceSBarry Smith     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
571b8f8c88eSHong Zhang     ierr = (*f)(sctx,x1_tmp,w1,fctx);CHKERRQ(ierr);
57266f9b7ceSBarry Smith     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
5734e269d77SPeter Brune   } else {
5744e269d77SPeter Brune     coloring->fset = PETSC_FALSE;
575be2fbe1fSBarry Smith   }
57643a90d84SBarry Smith 
577b8f8c88eSHong Zhang   if (!coloring->w3) {
578b8f8c88eSHong Zhang     ierr = VecDuplicate(x1_tmp,&coloring->w3);CHKERRQ(ierr);
5793bb1ff40SBarry Smith     ierr = PetscLogObjectParent((PetscObject)coloring,(PetscObject)coloring->w3);CHKERRQ(ierr);
580efb30889SBarry Smith   }
581b8f8c88eSHong Zhang   w3 = coloring->w3;
582efb30889SBarry Smith 
583b8f8c88eSHong Zhang   /* Compute all the local scale factors, including ghost points */
584b8f8c88eSHong Zhang   ierr = VecGetLocalSize(x1_tmp,&N);CHKERRQ(ierr);
585b8f8c88eSHong Zhang   ierr = VecGetArray(x1_tmp,&xx);CHKERRQ(ierr);
586b8f8c88eSHong Zhang   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
587b8f8c88eSHong Zhang   if (ctype == IS_COLORING_GHOSTED) {
588b8f8c88eSHong Zhang     col_start = 0; col_end = N;
5898ee2e534SBarry Smith   } else if (ctype == IS_COLORING_GLOBAL) {
590b8f8c88eSHong Zhang     xx           = xx - start;
591b8f8c88eSHong Zhang     vscale_array = vscale_array - start;
592b8f8c88eSHong Zhang     col_start    = start; col_end = N + start;
593b8f8c88eSHong Zhang   }
594b8f8c88eSHong Zhang   for (col=col_start; col<col_end; col++) {
595b8f8c88eSHong Zhang     /* Loop over each local column, vscale[col] = 1./(epsilon*dx[col]) */
596efb30889SBarry Smith     if (coloring->htype[0] == 'w') {
597efb30889SBarry Smith       dx = 1.0 + unorm;
598efb30889SBarry Smith     } else {
5999782cf98SBarry Smith       dx = xx[col];
600efb30889SBarry Smith     }
601d4a378daSJed Brown     if (dx == (PetscScalar)0.0) dx = 1.0;
6029782cf98SBarry Smith     if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
6039782cf98SBarry Smith     else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
6049782cf98SBarry Smith     dx               *= epsilon;
605d4a378daSJed Brown     vscale_array[col] = (PetscScalar)1.0/dx;
6069782cf98SBarry Smith   }
6078ee2e534SBarry Smith   if (ctype == IS_COLORING_GLOBAL) vscale_array = vscale_array + start;
608efb30889SBarry Smith   ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
6098ee2e534SBarry Smith   if (ctype == IS_COLORING_GLOBAL) {
61030b34957SBarry Smith     ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
61130b34957SBarry Smith     ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
612b8f8c88eSHong Zhang   }
6139782cf98SBarry Smith 
614b8f8c88eSHong Zhang   if (coloring->vscaleforrow) {
615b8f8c88eSHong Zhang     vscaleforrow = coloring->vscaleforrow;
616ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)J),PETSC_ERR_ARG_NULL,"Null Object: coloring->vscaleforrow");
617b0a32e0cSBarry Smith 
6189782cf98SBarry Smith   /*
61943a90d84SBarry Smith     Loop over each color
62043a90d84SBarry Smith   */
621b8f8c88eSHong Zhang   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
62243a90d84SBarry Smith   for (k=0; k<coloring->ncolors; k++) {
62349b058dcSBarry Smith     coloring->currentcolor = k;
62426fbe8dcSKarl Rupp 
625b8f8c88eSHong Zhang     ierr = VecCopy(x1_tmp,w3);CHKERRQ(ierr);
626b8f8c88eSHong Zhang     ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);
6278ee2e534SBarry Smith     if (ctype == IS_COLORING_GLOBAL) w3_array = w3_array - start;
62843a90d84SBarry Smith     /*
629b8f8c88eSHong Zhang       Loop over each column associated with color
630b8f8c88eSHong Zhang       adding the perturbation to the vector w3.
63143a90d84SBarry Smith     */
63243a90d84SBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
633b8f8c88eSHong Zhang       col = coloring->columns[k][l];    /* local column of the matrix we are probing for */
634efb30889SBarry Smith       if (coloring->htype[0] == 'w') {
635efb30889SBarry Smith         dx = 1.0 + unorm;
636efb30889SBarry Smith       } else {
63742460c72SBarry Smith         dx = xx[col];
638efb30889SBarry Smith       }
639d4a378daSJed Brown       if (dx == (PetscScalar)0.0) dx = 1.0;
640329f5518SBarry Smith       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
641329f5518SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
64243a90d84SBarry Smith       dx            *= epsilon;
643e32f2f54SBarry Smith       if (!PetscAbsScalar(dx)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Computed 0 differencing parameter");
64442460c72SBarry Smith       w3_array[col] += dx;
64543a90d84SBarry Smith     }
64698d79826SSatish Balay     if (ctype == IS_COLORING_GLOBAL) w3_array = w3_array + start;
647b8f8c88eSHong Zhang     ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
6483b28642cSBarry Smith 
64943a90d84SBarry Smith     /*
650b8f8c88eSHong Zhang       Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations)
651b8f8c88eSHong Zhang                            w2 = F(x1 + dx) - F(x1)
65243a90d84SBarry Smith     */
65366f9b7ceSBarry Smith     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
65443a90d84SBarry Smith     ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr);
65566f9b7ceSBarry Smith     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
656efb30889SBarry Smith     ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr);
6579782cf98SBarry Smith 
65843a90d84SBarry Smith     /*
659e0907662SLois Curfman McInnes       Loop over rows of vector, putting results into Jacobian matrix
66043a90d84SBarry Smith     */
6619a3c1acfSHong Zhang #if defined(JACOBIANCOLOROPT)
6629a3c1acfSHong Zhang     ierr = PetscTime(&t0);CHKERRQ(ierr);
6639a3c1acfSHong Zhang #endif
6643b28642cSBarry Smith     ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
66543a90d84SBarry Smith     for (l=0; l<coloring->nrows[k]; l++) {
666b8f8c88eSHong Zhang       row     = coloring->rows[k][l];            /* local row index */
667b8f8c88eSHong Zhang       col     = coloring->columnsforrow[k][l];   /* global column index */
6685904e1b1SBarry Smith       y[row] *= vscale_array[vscaleforrow[k][l]];
66943a90d84SBarry Smith       srow    = row + start;
67043a90d84SBarry Smith       ierr    = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
67143a90d84SBarry Smith     }
6723b28642cSBarry Smith     ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
6739a3c1acfSHong Zhang #if defined(JACOBIANCOLOROPT)
6749a3c1acfSHong Zhang     ierr = PetscTime(&t1);CHKERRQ(ierr);
6759a3c1acfSHong Zhang     time_setvalues += t1-t0;
6769a3c1acfSHong Zhang #endif
677aeb7e98aSMatthew Knepley   } /* endof for each color */
6789a3c1acfSHong Zhang #if defined(JACOBIANCOLOROPT)
6799a3c1acfSHong Zhang   printf("     MatFDColoringApply_AIJ: time_setvalues %g\n",time_setvalues);
6809a3c1acfSHong Zhang #endif
6818ee2e534SBarry Smith   if (ctype == IS_COLORING_GLOBAL) xx = xx + start;
68230b34957SBarry Smith   ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
683b8f8c88eSHong Zhang   ierr = VecRestoreArray(x1_tmp,&xx);CHKERRQ(ierr);
684b8f8c88eSHong Zhang 
685b8f8c88eSHong Zhang   coloring->currentcolor = -1;
6868865f1eaSKarl Rupp 
68743a90d84SBarry Smith   ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
68843a90d84SBarry Smith   ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
689d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
690a2e34c3dSBarry Smith 
691146574abSBarry Smith   ierr = MatFDColoringViewFromOptions(coloring,NULL,"-mat_fd_coloring_view");CHKERRQ(ierr);
6923a40ed3dSBarry Smith   PetscFunctionReturn(0);
69343a90d84SBarry Smith }
694