xref: /petsc/src/mat/matfd/fdmatrix.c (revision 0700a8246d308f50502909ba325e6169d3ee27eb)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
2be1d678aSKris Buschelman 
3bbf0e169SBarry Smith /*
4639f9d9dSBarry Smith    This is where the abstract matrix operations are defined that are
5639f9d9dSBarry Smith   used for finite difference computations of Jacobians using coloring.
6bbf0e169SBarry Smith */
7bbf0e169SBarry Smith 
87c4f633dSBarry Smith #include "private/matimpl.h"        /*I "petscmat.h" I*/
9bbf0e169SBarry Smith 
104a2ae208SSatish Balay #undef __FUNCT__
113a7fca6bSBarry Smith #define __FUNCT__ "MatFDColoringSetF"
12be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetF(MatFDColoring fd,Vec F)
133a7fca6bSBarry Smith {
143a7fca6bSBarry Smith   PetscFunctionBegin;
153a7fca6bSBarry Smith   fd->F = F;
163a7fca6bSBarry Smith   PetscFunctionReturn(0);
173a7fca6bSBarry Smith }
183a7fca6bSBarry Smith 
193a7fca6bSBarry Smith #undef __FUNCT__
204a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Draw_Zoom"
216849ba73SBarry Smith static PetscErrorCode MatFDColoringView_Draw_Zoom(PetscDraw draw,void *Aa)
229194eea9SBarry Smith {
239194eea9SBarry Smith   MatFDColoring  fd = (MatFDColoring)Aa;
24dfbe8321SBarry Smith   PetscErrorCode ierr;
2538baddfdSBarry Smith   PetscInt       i,j;
269194eea9SBarry Smith   PetscReal      x,y;
279194eea9SBarry Smith 
289194eea9SBarry Smith   PetscFunctionBegin;
299194eea9SBarry Smith 
309194eea9SBarry Smith   /* loop over colors  */
319194eea9SBarry Smith   for (i=0; i<fd->ncolors; i++) {
329194eea9SBarry Smith     for (j=0; j<fd->nrows[i]; j++) {
339194eea9SBarry Smith       y = fd->M - fd->rows[i][j] - fd->rstart;
349194eea9SBarry Smith       x = fd->columnsforrow[i][j];
35b0a32e0cSBarry Smith       ierr = PetscDrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);CHKERRQ(ierr);
369194eea9SBarry Smith     }
379194eea9SBarry Smith   }
389194eea9SBarry Smith   PetscFunctionReturn(0);
399194eea9SBarry Smith }
409194eea9SBarry Smith 
414a2ae208SSatish Balay #undef __FUNCT__
424a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Draw"
436849ba73SBarry Smith static PetscErrorCode MatFDColoringView_Draw(MatFDColoring fd,PetscViewer viewer)
44005c665bSBarry Smith {
45dfbe8321SBarry Smith   PetscErrorCode ierr;
46005c665bSBarry Smith   PetscTruth     isnull;
47b0a32e0cSBarry Smith   PetscDraw      draw;
4854d96fa3SBarry Smith   PetscReal      xr,yr,xl,yl,h,w;
49005c665bSBarry Smith 
503a40ed3dSBarry Smith   PetscFunctionBegin;
51b0a32e0cSBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
52b0a32e0cSBarry Smith   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
539194eea9SBarry Smith 
549194eea9SBarry Smith   ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
55005c665bSBarry Smith 
56005c665bSBarry Smith   xr  = fd->N; yr = fd->M; h = yr/10.0; w = xr/10.0;
57005c665bSBarry Smith   xr += w;     yr += h;    xl = -w;     yl = -h;
58b0a32e0cSBarry Smith   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
59b0a32e0cSBarry Smith   ierr = PetscDrawZoom(draw,MatFDColoringView_Draw_Zoom,fd);CHKERRQ(ierr);
609194eea9SBarry Smith   ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
613a40ed3dSBarry Smith   PetscFunctionReturn(0);
62005c665bSBarry Smith }
63005c665bSBarry Smith 
644a2ae208SSatish Balay #undef __FUNCT__
654a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView"
66bbf0e169SBarry Smith /*@C
67639f9d9dSBarry Smith    MatFDColoringView - Views a finite difference coloring context.
68bbf0e169SBarry Smith 
694c49b128SBarry Smith    Collective on MatFDColoring
70fee21e36SBarry Smith 
71ef5ee4d1SLois Curfman McInnes    Input  Parameters:
72ef5ee4d1SLois Curfman McInnes +  c - the coloring context
73ef5ee4d1SLois Curfman McInnes -  viewer - visualization context
74ef5ee4d1SLois Curfman McInnes 
7515091d37SBarry Smith    Level: intermediate
7615091d37SBarry Smith 
77b4fc646aSLois Curfman McInnes    Notes:
78b4fc646aSLois Curfman McInnes    The available visualization contexts include
79b0a32e0cSBarry Smith +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
80b0a32e0cSBarry Smith .     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
81ef5ee4d1SLois Curfman McInnes         output where only the first processor opens
82ef5ee4d1SLois Curfman McInnes         the file.  All other processors send their
83ef5ee4d1SLois Curfman McInnes         data to the first processor to print.
84b0a32e0cSBarry Smith -     PETSC_VIEWER_DRAW_WORLD - graphical display of nonzero structure
85bbf0e169SBarry Smith 
869194eea9SBarry Smith    Notes:
879194eea9SBarry Smith      Since PETSc uses only a small number of basic colors (currently 33), if the coloring
88b0a32e0cSBarry Smith    involves more than 33 then some seemingly identical colors are displayed making it look
899194eea9SBarry Smith    like an illegal coloring. This is just a graphical artifact.
909194eea9SBarry Smith 
91639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
92005c665bSBarry Smith 
93b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, view
94bbf0e169SBarry Smith @*/
95be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringView(MatFDColoring c,PetscViewer viewer)
96bbf0e169SBarry Smith {
976849ba73SBarry Smith   PetscErrorCode    ierr;
9838baddfdSBarry Smith   PetscInt          i,j;
9932077d6dSBarry Smith   PetscTruth        isdraw,iascii;
100f3ef73ceSBarry Smith   PetscViewerFormat format;
101bbf0e169SBarry Smith 
1023a40ed3dSBarry Smith   PetscFunctionBegin;
103*0700a824SBarry Smith   PetscValidHeaderSpecific(c,MAT_FDCOLORING_CLASSID,1);
1043050cee2SBarry Smith   if (!viewer) {
1057adad957SLisandro Dalcin     ierr = PetscViewerASCIIGetStdout(((PetscObject)c)->comm,&viewer);CHKERRQ(ierr);
1063050cee2SBarry Smith   }
107*0700a824SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
108c9780b6fSBarry Smith   PetscCheckSameComm(c,1,viewer,2);
109bbf0e169SBarry Smith 
110fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
11132077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
1120f5bd95cSBarry Smith   if (isdraw) {
113b4fc646aSLois Curfman McInnes     ierr = MatFDColoringView_Draw(c,viewer);CHKERRQ(ierr);
11432077d6dSBarry Smith   } else if (iascii) {
115b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"MatFDColoring Object:\n");CHKERRQ(ierr);
116a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Error tolerance=%G\n",c->error_rel);CHKERRQ(ierr);
117a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Umin=%G\n",c->umin);CHKERRQ(ierr);
11877431f27SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Number of colors=%D\n",c->ncolors);CHKERRQ(ierr);
119ae09f205SBarry Smith 
120b0a32e0cSBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
121fb9695e5SSatish Balay     if (format != PETSC_VIEWER_ASCII_INFO) {
122b4fc646aSLois Curfman McInnes       for (i=0; i<c->ncolors; i++) {
12377431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"  Information for color %D\n",i);CHKERRQ(ierr);
12477431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    Number of columns %D\n",c->ncolumns[i]);CHKERRQ(ierr);
125b4fc646aSLois Curfman McInnes         for (j=0; j<c->ncolumns[i]; j++) {
12677431f27SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"      %D\n",c->columns[i][j]);CHKERRQ(ierr);
127639f9d9dSBarry Smith         }
12877431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    Number of rows %D\n",c->nrows[i]);CHKERRQ(ierr);
129b4fc646aSLois Curfman McInnes         for (j=0; j<c->nrows[i]; j++) {
13077431f27SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"      %D %D \n",c->rows[i][j],c->columnsforrow[i][j]);CHKERRQ(ierr);
131b4fc646aSLois Curfman McInnes         }
132bbf0e169SBarry Smith       }
133bbf0e169SBarry Smith     }
134b0a32e0cSBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1355cd90555SBarry Smith   } else {
13679a5c55eSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for MatFDColoring",((PetscObject)viewer)->type_name);
137bbf0e169SBarry Smith   }
1383a40ed3dSBarry Smith   PetscFunctionReturn(0);
139639f9d9dSBarry Smith }
140639f9d9dSBarry Smith 
1414a2ae208SSatish Balay #undef __FUNCT__
1424a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetParameters"
143639f9d9dSBarry Smith /*@
144b4fc646aSLois Curfman McInnes    MatFDColoringSetParameters - Sets the parameters for the sparse approximation of
145b4fc646aSLois Curfman McInnes    a Jacobian matrix using finite differences.
146639f9d9dSBarry Smith 
147ef5ee4d1SLois Curfman McInnes    Collective on MatFDColoring
148ef5ee4d1SLois Curfman McInnes 
149ef5ee4d1SLois Curfman McInnes    The Jacobian is estimated with the differencing approximation
150ef5ee4d1SLois Curfman McInnes .vb
15165f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
152f23b5b22SLois Curfman McInnes        h = error_rel*u[i]                 if  abs(u[i]) > umin
153f23b5b22SLois Curfman McInnes          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
154ef5ee4d1SLois Curfman McInnes        dx_{i} = (0, ... 1, .... 0)
155ef5ee4d1SLois Curfman McInnes .ve
156639f9d9dSBarry Smith 
157639f9d9dSBarry Smith    Input Parameters:
158ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
159639f9d9dSBarry Smith .  error_rel - relative error
160f23b5b22SLois Curfman McInnes -  umin - minimum allowable u-value magnitude
161fee21e36SBarry Smith 
16215091d37SBarry Smith    Level: advanced
16315091d37SBarry Smith 
164b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, parameters
165b4fc646aSLois Curfman McInnes 
16695b89fc3SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringSetFromOptions()
16795b89fc3SBarry Smith 
168639f9d9dSBarry Smith @*/
169be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetParameters(MatFDColoring matfd,PetscReal error,PetscReal umin)
170639f9d9dSBarry Smith {
1713a40ed3dSBarry Smith   PetscFunctionBegin;
172*0700a824SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
173639f9d9dSBarry Smith 
174639f9d9dSBarry Smith   if (error != PETSC_DEFAULT) matfd->error_rel = error;
175639f9d9dSBarry Smith   if (umin != PETSC_DEFAULT)  matfd->umin      = umin;
1763a40ed3dSBarry Smith   PetscFunctionReturn(0);
177639f9d9dSBarry Smith }
178639f9d9dSBarry Smith 
179005c665bSBarry Smith 
180ff0cfa39SBarry Smith 
1814a2ae208SSatish Balay #undef __FUNCT__
1824a9d489dSBarry Smith #define __FUNCT__ "MatFDColoringGetFunction"
1834a9d489dSBarry Smith /*@C
1844a9d489dSBarry Smith    MatFDColoringGetFunction - Gets the function to use for computing the Jacobian.
1854a9d489dSBarry Smith 
1864a9d489dSBarry Smith    Collective on MatFDColoring
1874a9d489dSBarry Smith 
1884a9d489dSBarry Smith    Input Parameters:
1894a9d489dSBarry Smith .  coloring - the coloring context
1904a9d489dSBarry Smith 
1914a9d489dSBarry Smith    Output Parameters:
1924a9d489dSBarry Smith +  f - the function
1934a9d489dSBarry Smith -  fctx - the optional user-defined function context
1944a9d489dSBarry Smith 
1954a9d489dSBarry Smith    Level: intermediate
1964a9d489dSBarry Smith 
1974a9d489dSBarry Smith .keywords: Mat, Jacobian, finite differences, set, function
19895b89fc3SBarry Smith 
19995b89fc3SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringSetFunction(), MatFDColoringSetFromOptions()
20095b89fc3SBarry Smith 
2014a9d489dSBarry Smith @*/
2024a9d489dSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFunction(MatFDColoring matfd,PetscErrorCode (**f)(void),void **fctx)
2034a9d489dSBarry Smith {
2044a9d489dSBarry Smith   PetscFunctionBegin;
205*0700a824SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
2064a9d489dSBarry Smith   if (f) *f = matfd->f;
2074a9d489dSBarry Smith   if (fctx) *fctx = matfd->fctx;
2084a9d489dSBarry Smith   PetscFunctionReturn(0);
2094a9d489dSBarry Smith }
2104a9d489dSBarry Smith 
2114a9d489dSBarry Smith #undef __FUNCT__
2124a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFunction"
213d64ed03dSBarry Smith /*@C
214005c665bSBarry Smith    MatFDColoringSetFunction - Sets the function to use for computing the Jacobian.
215005c665bSBarry Smith 
216fee21e36SBarry Smith    Collective on MatFDColoring
217fee21e36SBarry Smith 
218ef5ee4d1SLois Curfman McInnes    Input Parameters:
219ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
220ef5ee4d1SLois Curfman McInnes .  f - the function
221ef5ee4d1SLois Curfman McInnes -  fctx - the optional user-defined function context
222ef5ee4d1SLois Curfman McInnes 
2237850c7c0SBarry Smith    Calling sequence of (*f) function:
2247850c7c0SBarry Smith     For SNES:    PetscErrorCode (*f)(SNES,Vec,Vec,void*)
2257850c7c0SBarry Smith     For TS:      PetscErrorCode (*f)(TS,PetscReal,Vec,Vec,void*)
2267850c7c0SBarry Smith     If not using SNES or TS: PetscErrorCode (*f)(void *dummy,Vec,Vec,void*) and dummy is ignored
22715091d37SBarry Smith 
2287850c7c0SBarry Smith    Level: advanced
2297850c7c0SBarry Smith 
2307850c7c0SBarry Smith    Notes: This function is usually used automatically by SNES or TS (when one uses SNESSetJacobian() with the argument
2317850c7c0SBarry Smith      SNESDefaultComputeJacobianColor() or TSSetRHSJacobian() with the argument TSDefaultComputeJacobianColor()) and only needs to be used
2327850c7c0SBarry Smith      by someone computing a matrix via coloring directly by calling MatFDColoringApply()
2337850c7c0SBarry Smith 
2347850c7c0SBarry Smith    Fortran Notes:
2357850c7c0SBarry Smith     In Fortran you must call MatFDColoringSetFunction() for a coloring object to
2367850c7c0SBarry Smith   be used without SNES or TS or within the SNES solvers and MatFDColoringSetFunctionTS() if it is to be used
2377850c7c0SBarry Smith   within the TS solvers.
238f881d145SBarry Smith 
239b4fc646aSLois Curfman McInnes .keywords: Mat, Jacobian, finite differences, set, function
24095b89fc3SBarry Smith 
24195b89fc3SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringGetFunction(), MatFDColoringSetFromOptions()
24295b89fc3SBarry Smith 
243005c665bSBarry Smith @*/
244be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFunction(MatFDColoring matfd,PetscErrorCode (*f)(void),void *fctx)
245005c665bSBarry Smith {
2463a40ed3dSBarry Smith   PetscFunctionBegin;
247*0700a824SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
248005c665bSBarry Smith   matfd->f    = f;
249005c665bSBarry Smith   matfd->fctx = fctx;
2503a40ed3dSBarry Smith   PetscFunctionReturn(0);
251005c665bSBarry Smith }
252005c665bSBarry Smith 
2534a2ae208SSatish Balay #undef __FUNCT__
2544a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFromOptions"
255639f9d9dSBarry Smith /*@
256b4fc646aSLois Curfman McInnes    MatFDColoringSetFromOptions - Sets coloring finite difference parameters from
257639f9d9dSBarry Smith    the options database.
258639f9d9dSBarry Smith 
259fee21e36SBarry Smith    Collective on MatFDColoring
260fee21e36SBarry Smith 
26165f2ba5bSLois Curfman McInnes    The Jacobian, F'(u), is estimated with the differencing approximation
262ef5ee4d1SLois Curfman McInnes .vb
26365f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
264f23b5b22SLois Curfman McInnes        h = error_rel*u[i]                 if  abs(u[i]) > umin
265f23b5b22SLois Curfman McInnes          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
266ef5ee4d1SLois Curfman McInnes        dx_{i} = (0, ... 1, .... 0)
267ef5ee4d1SLois Curfman McInnes .ve
268ef5ee4d1SLois Curfman McInnes 
269ef5ee4d1SLois Curfman McInnes    Input Parameter:
270ef5ee4d1SLois Curfman McInnes .  coloring - the coloring context
271ef5ee4d1SLois Curfman McInnes 
272b4fc646aSLois Curfman McInnes    Options Database Keys:
273d15ffeeaSSatish Balay +  -mat_fd_coloring_err <err> - Sets <err> (square root
274ef5ee4d1SLois Curfman McInnes            of relative error in the function)
275f23b5b22SLois Curfman McInnes .  -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude
2763ec795f1SBarry Smith .  -mat_fd_type - "wp" or "ds" (see MATMFFD_WP or MATMFFD_DS)
277ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_view - Activates basic viewing
278ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_view_info - Activates viewing info
279ef5ee4d1SLois Curfman McInnes -  -mat_fd_coloring_view_draw - Activates drawing
280639f9d9dSBarry Smith 
28115091d37SBarry Smith     Level: intermediate
28215091d37SBarry Smith 
283b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, parameters
284d1c39f55SBarry Smith 
285d1c39f55SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringView(), MatFDColoringSetParameters()
286d1c39f55SBarry Smith 
287639f9d9dSBarry Smith @*/
288be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFromOptions(MatFDColoring matfd)
289639f9d9dSBarry Smith {
290dfbe8321SBarry Smith   PetscErrorCode ierr;
291efb30889SBarry Smith   PetscTruth     flg;
292efb30889SBarry Smith   char           value[3];
2933a40ed3dSBarry Smith 
2943a40ed3dSBarry Smith   PetscFunctionBegin;
295*0700a824SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
296639f9d9dSBarry Smith 
2977adad957SLisandro Dalcin   ierr = PetscOptionsBegin(((PetscObject)matfd)->comm,((PetscObject)matfd)->prefix,"Jacobian computation via finite differences option","MatFD");CHKERRQ(ierr);
29887828ca2SBarry Smith     ierr = PetscOptionsReal("-mat_fd_coloring_err","Square root of relative error in function","MatFDColoringSetParameters",matfd->error_rel,&matfd->error_rel,0);CHKERRQ(ierr);
29987828ca2SBarry Smith     ierr = PetscOptionsReal("-mat_fd_coloring_umin","Minimum allowable u magnitude","MatFDColoringSetParameters",matfd->umin,&matfd->umin,0);CHKERRQ(ierr);
300efb30889SBarry Smith     ierr = PetscOptionsString("-mat_fd_type","Algorithm to compute h, wp or ds","MatFDColoringCreate",matfd->htype,value,2,&flg);CHKERRQ(ierr);
301efb30889SBarry Smith     if (flg) {
302efb30889SBarry Smith       if (value[0] == 'w' && value[1] == 'p') matfd->htype = "wp";
303efb30889SBarry Smith       else if (value[0] == 'd' && value[1] == 's') matfd->htype = "ds";
304efb30889SBarry Smith       else SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Unknown finite differencing type %s",value);
305efb30889SBarry Smith     }
306186905e3SBarry Smith     /* not used here; but so they are presented in the GUI */
307b0a32e0cSBarry Smith     ierr = PetscOptionsName("-mat_fd_coloring_view","Print entire datastructure used for Jacobian","None",0);CHKERRQ(ierr);
308b0a32e0cSBarry Smith     ierr = PetscOptionsName("-mat_fd_coloring_view_info","Print number of colors etc for Jacobian","None",0);CHKERRQ(ierr);
309b0a32e0cSBarry Smith     ierr = PetscOptionsName("-mat_fd_coloring_view_draw","Plot nonzero structure ofJacobian","None",0);CHKERRQ(ierr);
310b0a32e0cSBarry Smith   PetscOptionsEnd();CHKERRQ(ierr);
3113a40ed3dSBarry Smith   PetscFunctionReturn(0);
312005c665bSBarry Smith }
313005c665bSBarry Smith 
3144a2ae208SSatish Balay #undef __FUNCT__
3154a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Private"
316dfbe8321SBarry Smith PetscErrorCode MatFDColoringView_Private(MatFDColoring fd)
317005c665bSBarry Smith {
318dfbe8321SBarry Smith   PetscErrorCode ierr;
31990d69ab7SBarry Smith   PetscTruth     flg = PETSC_FALSE;
3203050cee2SBarry Smith   PetscViewer    viewer;
321005c665bSBarry Smith 
3223a40ed3dSBarry Smith   PetscFunctionBegin;
3237adad957SLisandro Dalcin   ierr = PetscViewerASCIIGetStdout(((PetscObject)fd)->comm,&viewer);CHKERRQ(ierr);
32490d69ab7SBarry Smith   ierr = PetscOptionsGetTruth(PETSC_NULL,"-mat_fd_coloring_view",&flg,PETSC_NULL);CHKERRQ(ierr);
325005c665bSBarry Smith   if (flg) {
3263050cee2SBarry Smith     ierr = MatFDColoringView(fd,viewer);CHKERRQ(ierr);
327005c665bSBarry Smith   }
32890d69ab7SBarry Smith   flg  = PETSC_FALSE;
32990d69ab7SBarry Smith   ierr = PetscOptionsGetTruth(PETSC_NULL,"-mat_fd_coloring_view_info",&flg,PETSC_NULL);CHKERRQ(ierr);
330ae09f205SBarry Smith   if (flg) {
3313050cee2SBarry Smith     ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
3323050cee2SBarry Smith     ierr = MatFDColoringView(fd,viewer);CHKERRQ(ierr);
3333050cee2SBarry Smith     ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
334ae09f205SBarry Smith   }
33590d69ab7SBarry Smith   flg  = PETSC_FALSE;
33690d69ab7SBarry Smith   ierr = PetscOptionsGetTruth(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg,PETSC_NULL);CHKERRQ(ierr);
337005c665bSBarry Smith   if (flg) {
3387adad957SLisandro Dalcin     ierr = MatFDColoringView(fd,PETSC_VIEWER_DRAW_(((PetscObject)fd)->comm));CHKERRQ(ierr);
3397adad957SLisandro Dalcin     ierr = PetscViewerFlush(PETSC_VIEWER_DRAW_(((PetscObject)fd)->comm));CHKERRQ(ierr);
340005c665bSBarry Smith   }
3413a40ed3dSBarry Smith   PetscFunctionReturn(0);
342bbf0e169SBarry Smith }
343bbf0e169SBarry Smith 
3444a2ae208SSatish Balay #undef __FUNCT__
3454a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringCreate"
34605869f15SSatish Balay /*@
347639f9d9dSBarry Smith    MatFDColoringCreate - Creates a matrix coloring context for finite difference
348639f9d9dSBarry Smith    computation of Jacobians.
349bbf0e169SBarry Smith 
350ef5ee4d1SLois Curfman McInnes    Collective on Mat
351ef5ee4d1SLois Curfman McInnes 
352639f9d9dSBarry Smith    Input Parameters:
353ef5ee4d1SLois Curfman McInnes +  mat - the matrix containing the nonzero structure of the Jacobian
354ef5ee4d1SLois Curfman McInnes -  iscoloring - the coloring of the matrix
355bbf0e169SBarry Smith 
356bbf0e169SBarry Smith     Output Parameter:
357639f9d9dSBarry Smith .   color - the new coloring context
358bbf0e169SBarry Smith 
35915091d37SBarry Smith     Level: intermediate
36015091d37SBarry Smith 
3616831982aSBarry Smith .seealso: MatFDColoringDestroy(),SNESDefaultComputeJacobianColor(), ISColoringCreate(),
362d1c39f55SBarry Smith           MatFDColoringSetFunction(), MatFDColoringSetFromOptions(), MatFDColoringApply(),
363ebd3b9afSBarry Smith           MatFDColoringView(), MatFDColoringSetParameters()
364bbf0e169SBarry Smith @*/
365be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
366bbf0e169SBarry Smith {
367639f9d9dSBarry Smith   MatFDColoring  c;
368639f9d9dSBarry Smith   MPI_Comm       comm;
369dfbe8321SBarry Smith   PetscErrorCode ierr;
37038baddfdSBarry Smith   PetscInt       M,N;
371b8f8c88eSHong Zhang   PetscMPIInt    size;
372639f9d9dSBarry Smith 
3733a40ed3dSBarry Smith   PetscFunctionBegin;
374d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
375639f9d9dSBarry Smith   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
37629bbc08cSBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,"Only for square matrices");
377639f9d9dSBarry Smith 
378f881d145SBarry Smith   ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
379*0700a824SBarry Smith   ierr = PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_CLASSID,0,"MatFDColoring",comm,MatFDColoringDestroy,MatFDColoringView);CHKERRQ(ierr);
380639f9d9dSBarry Smith 
381b8f8c88eSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
382b8f8c88eSHong Zhang   c->ctype = iscoloring->ctype;
383b8f8c88eSHong Zhang 
384f830108cSBarry Smith   if (mat->ops->fdcoloringcreate) {
385f830108cSBarry Smith     ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);CHKERRQ(ierr);
386639f9d9dSBarry Smith   } else {
38729bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Code not yet written for this matrix type");
388639f9d9dSBarry Smith   }
389639f9d9dSBarry Smith 
390b8f8c88eSHong Zhang   ierr = MatGetVecs(mat,PETSC_NULL,&c->w1);CHKERRQ(ierr);
391b8f8c88eSHong Zhang   ierr = PetscLogObjectParent(c,c->w1);CHKERRQ(ierr);
392b8f8c88eSHong Zhang   ierr = VecDuplicate(c->w1,&c->w2);CHKERRQ(ierr);
393b8f8c88eSHong Zhang   ierr = PetscLogObjectParent(c,c->w2);CHKERRQ(ierr);
394b8f8c88eSHong Zhang 
39577d8c4bbSBarry Smith   c->error_rel         = PETSC_SQRT_MACHINE_EPSILON;
39677d8c4bbSBarry Smith   c->umin              = 100.0*PETSC_SQRT_MACHINE_EPSILON;
39749b058dcSBarry Smith   c->currentcolor      = -1;
398efb30889SBarry Smith   c->htype             = "wp";
399639f9d9dSBarry Smith 
400639f9d9dSBarry Smith   *color = c;
401d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
4023a40ed3dSBarry Smith   PetscFunctionReturn(0);
403bbf0e169SBarry Smith }
404bbf0e169SBarry Smith 
4054a2ae208SSatish Balay #undef __FUNCT__
4064a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringDestroy"
40705869f15SSatish Balay /*@
408639f9d9dSBarry Smith     MatFDColoringDestroy - Destroys a matrix coloring context that was created
409639f9d9dSBarry Smith     via MatFDColoringCreate().
410bbf0e169SBarry Smith 
411ef5ee4d1SLois Curfman McInnes     Collective on MatFDColoring
412ef5ee4d1SLois Curfman McInnes 
413b4fc646aSLois Curfman McInnes     Input Parameter:
414639f9d9dSBarry Smith .   c - coloring context
415bbf0e169SBarry Smith 
41615091d37SBarry Smith     Level: intermediate
41715091d37SBarry Smith 
418639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
419bbf0e169SBarry Smith @*/
420be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringDestroy(MatFDColoring c)
421bbf0e169SBarry Smith {
4226849ba73SBarry Smith   PetscErrorCode ierr;
42338baddfdSBarry Smith   PetscInt       i;
424bbf0e169SBarry Smith 
4253a40ed3dSBarry Smith   PetscFunctionBegin;
4267adad957SLisandro Dalcin   if (--((PetscObject)c)->refct > 0) PetscFunctionReturn(0);
427d4bb536fSBarry Smith 
428639f9d9dSBarry Smith   for (i=0; i<c->ncolors; i++) {
42905b42c5fSBarry Smith     ierr = PetscFree(c->columns[i]);CHKERRQ(ierr);
43005b42c5fSBarry Smith     ierr = PetscFree(c->rows[i]);CHKERRQ(ierr);
43105b42c5fSBarry Smith     ierr = PetscFree(c->columnsforrow[i]);CHKERRQ(ierr);
43205b42c5fSBarry Smith     if (c->vscaleforrow) {ierr = PetscFree(c->vscaleforrow[i]);CHKERRQ(ierr);}
433bbf0e169SBarry Smith   }
434606d414cSSatish Balay   ierr = PetscFree(c->ncolumns);CHKERRQ(ierr);
435606d414cSSatish Balay   ierr = PetscFree(c->columns);CHKERRQ(ierr);
436606d414cSSatish Balay   ierr = PetscFree(c->nrows);CHKERRQ(ierr);
437606d414cSSatish Balay   ierr = PetscFree(c->rows);CHKERRQ(ierr);
438606d414cSSatish Balay   ierr = PetscFree(c->columnsforrow);CHKERRQ(ierr);
43905b42c5fSBarry Smith   ierr = PetscFree(c->vscaleforrow);CHKERRQ(ierr);
4404f113abfSBarry Smith   if (c->vscale)       {ierr = VecDestroy(c->vscale);CHKERRQ(ierr);}
441005c665bSBarry Smith   if (c->w1) {
442005c665bSBarry Smith     ierr = VecDestroy(c->w1);CHKERRQ(ierr);
443005c665bSBarry Smith     ierr = VecDestroy(c->w2);CHKERRQ(ierr);
444b8f8c88eSHong Zhang   }
445b8f8c88eSHong Zhang   if (c->w3){
446005c665bSBarry Smith     ierr = VecDestroy(c->w3);CHKERRQ(ierr);
447005c665bSBarry Smith   }
448d38fa0fbSBarry Smith   ierr = PetscHeaderDestroy(c);CHKERRQ(ierr);
4493a40ed3dSBarry Smith   PetscFunctionReturn(0);
450bbf0e169SBarry Smith }
45143a90d84SBarry Smith 
4524a2ae208SSatish Balay #undef __FUNCT__
45349b058dcSBarry Smith #define __FUNCT__ "MatFDColoringGetPerturbedColumns"
45449b058dcSBarry Smith /*@C
45549b058dcSBarry Smith     MatFDColoringGetPerturbedColumns - Returns the indices of the columns that
45649b058dcSBarry Smith       that are currently being perturbed.
45749b058dcSBarry Smith 
45849b058dcSBarry Smith     Not Collective
45949b058dcSBarry Smith 
46049b058dcSBarry Smith     Input Parameters:
46149b058dcSBarry Smith .   coloring - coloring context created with MatFDColoringCreate()
46249b058dcSBarry Smith 
46349b058dcSBarry Smith     Output Parameters:
46449b058dcSBarry Smith +   n - the number of local columns being perturbed
46549b058dcSBarry Smith -   cols - the column indices, in global numbering
46649b058dcSBarry Smith 
46749b058dcSBarry Smith    Level: intermediate
46849b058dcSBarry Smith 
46949b058dcSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringApply()
47049b058dcSBarry Smith 
47149b058dcSBarry Smith .keywords: coloring, Jacobian, finite differences
47249b058dcSBarry Smith @*/
473be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetPerturbedColumns(MatFDColoring coloring,PetscInt *n,PetscInt *cols[])
47449b058dcSBarry Smith {
47549b058dcSBarry Smith   PetscFunctionBegin;
47649b058dcSBarry Smith   if (coloring->currentcolor >= 0) {
47749b058dcSBarry Smith     *n    = coloring->ncolumns[coloring->currentcolor];
47849b058dcSBarry Smith     *cols = coloring->columns[coloring->currentcolor];
47949b058dcSBarry Smith   } else {
48049b058dcSBarry Smith     *n = 0;
48149b058dcSBarry Smith   }
48249b058dcSBarry Smith   PetscFunctionReturn(0);
48349b058dcSBarry Smith }
48449b058dcSBarry Smith 
48549b058dcSBarry Smith #undef __FUNCT__
4864a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringApply"
48743a90d84SBarry Smith /*@
488e0907662SLois Curfman McInnes     MatFDColoringApply - Given a matrix for which a MatFDColoring context
489e0907662SLois Curfman McInnes     has been created, computes the Jacobian for a function via finite differences.
49043a90d84SBarry Smith 
491fee21e36SBarry Smith     Collective on MatFDColoring
492fee21e36SBarry Smith 
493ef5ee4d1SLois Curfman McInnes     Input Parameters:
494ef5ee4d1SLois Curfman McInnes +   mat - location to store Jacobian
495ef5ee4d1SLois Curfman McInnes .   coloring - coloring context created with MatFDColoringCreate()
496ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
4977850c7c0SBarry Smith -   sctx - context required by function, if this is being used with the SNES solver then it is SNES object, otherwise it is null
498ef5ee4d1SLois Curfman McInnes 
4998bba8e72SBarry Smith     Options Database Keys:
500ebd3b9afSBarry Smith +    -mat_fd_type - "wp" or "ds"  (see MATMFFD_WP or MATMFFD_DS)
501b9722508SLois Curfman McInnes .    -mat_fd_coloring_view - Activates basic viewing or coloring
502b9722508SLois Curfman McInnes .    -mat_fd_coloring_view_draw - Activates drawing of coloring
503b9722508SLois Curfman McInnes -    -mat_fd_coloring_view_info - Activates viewing of coloring info
5048bba8e72SBarry Smith 
50515091d37SBarry Smith     Level: intermediate
50615091d37SBarry Smith 
5077850c7c0SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringSetFunction()
50843a90d84SBarry Smith 
50943a90d84SBarry Smith .keywords: coloring, Jacobian, finite differences
51043a90d84SBarry Smith @*/
511be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
51243a90d84SBarry Smith {
5133acb8795SBarry Smith   PetscErrorCode ierr;
5143acb8795SBarry Smith 
5153acb8795SBarry Smith   PetscFunctionBegin;
516*0700a824SBarry Smith   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
517*0700a824SBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_CLASSID,2);
518*0700a824SBarry Smith   PetscValidHeaderSpecific(x1,VEC_CLASSID,3);
5193acb8795SBarry Smith   if (!coloring->f) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetFunction()");
5203acb8795SBarry Smith   if (!J->ops->fdcoloringapply) SETERRQ1(PETSC_ERR_SUP,"Not supported for this matrix type %s",((PetscObject)J)->type_name);
5213acb8795SBarry Smith   ierr = (*J->ops->fdcoloringapply)(J,coloring,x1,flag,sctx);CHKERRQ(ierr);
5223acb8795SBarry Smith   PetscFunctionReturn(0);
5233acb8795SBarry Smith }
5243acb8795SBarry Smith 
5253acb8795SBarry Smith #undef __FUNCT__
5263acb8795SBarry Smith #define __FUNCT__ "MatFDColoringApply_AIJ"
5273acb8795SBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApply_AIJ(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
5283acb8795SBarry Smith {
5296849ba73SBarry Smith   PetscErrorCode (*f)(void*,Vec,Vec,void*) = (PetscErrorCode (*)(void*,Vec,Vec,void *))coloring->f;
5306849ba73SBarry Smith   PetscErrorCode ierr;
531b8f8c88eSHong Zhang   PetscInt       k,start,end,l,row,col,srow,**vscaleforrow,m1,m2;
532efb30889SBarry Smith   PetscScalar    dx,*y,*xx,*w3_array;
533ea709b57SSatish Balay   PetscScalar    *vscale_array;
534efb30889SBarry Smith   PetscReal      epsilon = coloring->error_rel,umin = coloring->umin,unorm;
535b8f8c88eSHong Zhang   Vec            w1=coloring->w1,w2=coloring->w2,w3;
536005c665bSBarry Smith   void           *fctx = coloring->fctx;
53790d69ab7SBarry Smith   PetscTruth     flg = PETSC_FALSE;
538705d48f7SSatish Balay   PetscInt       ctype=coloring->ctype,N,col_start=0,col_end=0;
539b8f8c88eSHong Zhang   Vec            x1_tmp;
540a2e34c3dSBarry Smith 
5413a40ed3dSBarry Smith   PetscFunctionBegin;
542*0700a824SBarry Smith   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
543*0700a824SBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_CLASSID,2);
544*0700a824SBarry Smith   PetscValidHeaderSpecific(x1,VEC_CLASSID,3);
545feba9168SBarry Smith   if (!f) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetFunction()");
546e0907662SLois Curfman McInnes 
547d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
5484c49b128SBarry Smith   ierr = MatSetUnfactored(J);CHKERRQ(ierr);
54990d69ab7SBarry Smith   ierr = PetscOptionsGetTruth(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg,PETSC_NULL);CHKERRQ(ierr);
550f1af5d2fSBarry Smith   if (flg) {
551ae15b995SBarry Smith     ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr);
552e0907662SLois Curfman McInnes   } else {
5530b9b6f31SBarry Smith     PetscTruth assembled;
5540b9b6f31SBarry Smith     ierr = MatAssembled(J,&assembled);CHKERRQ(ierr);
5550b9b6f31SBarry Smith     if (assembled) {
55643a90d84SBarry Smith       ierr = MatZeroEntries(J);CHKERRQ(ierr);
557e0907662SLois Curfman McInnes     }
5580b9b6f31SBarry Smith   }
55943a90d84SBarry Smith 
560b8f8c88eSHong Zhang   x1_tmp = x1;
561dac7f5fdSBarry Smith   if (!coloring->vscale){
562b8f8c88eSHong Zhang     ierr = VecDuplicate(x1_tmp,&coloring->vscale);CHKERRQ(ierr);
563b8f8c88eSHong Zhang   }
564be2fbe1fSBarry Smith 
5653a7fca6bSBarry Smith   /*
5663a7fca6bSBarry Smith     This is a horrible, horrible, hack. See DMMGComputeJacobian_Multigrid() it inproperly sets
5673a7fca6bSBarry Smith     coloring->F for the coarser grids from the finest
5683a7fca6bSBarry Smith   */
5693a7fca6bSBarry Smith   if (coloring->F) {
5703a7fca6bSBarry Smith     ierr = VecGetLocalSize(coloring->F,&m1);CHKERRQ(ierr);
5713a7fca6bSBarry Smith     ierr = VecGetLocalSize(w1,&m2);CHKERRQ(ierr);
5723a7fca6bSBarry Smith     if (m1 != m2) {
5733a7fca6bSBarry Smith       coloring->F = 0;
5743a7fca6bSBarry Smith       }
5753a7fca6bSBarry Smith     }
5763a7fca6bSBarry Smith 
577b8f8c88eSHong Zhang   if (coloring->htype[0] == 'w') { /* tacky test; need to make systematic if we add other approaches to computing h*/
578b8f8c88eSHong Zhang     ierr = VecNorm(x1_tmp,NORM_2,&unorm);CHKERRQ(ierr);
579b8f8c88eSHong Zhang   }
580b8f8c88eSHong Zhang   ierr = VecGetOwnershipRange(w1,&start,&end);CHKERRQ(ierr); /* OwnershipRange is used by ghosted x! */
581b8f8c88eSHong Zhang 
582b8f8c88eSHong Zhang   /* Set w1 = F(x1) */
583be2fbe1fSBarry Smith   if (coloring->F) {
584be2fbe1fSBarry Smith     w1          = coloring->F; /* use already computed value of function */
585be2fbe1fSBarry Smith     coloring->F = 0;
586be2fbe1fSBarry Smith   } else {
58766f9b7ceSBarry Smith     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
588b8f8c88eSHong Zhang     ierr = (*f)(sctx,x1_tmp,w1,fctx);CHKERRQ(ierr);
58966f9b7ceSBarry Smith     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
590be2fbe1fSBarry Smith   }
59143a90d84SBarry Smith 
592b8f8c88eSHong Zhang   if (!coloring->w3) {
593b8f8c88eSHong Zhang     ierr = VecDuplicate(x1_tmp,&coloring->w3);CHKERRQ(ierr);
594b8f8c88eSHong Zhang     ierr = PetscLogObjectParent(coloring,coloring->w3);CHKERRQ(ierr);
595efb30889SBarry Smith   }
596b8f8c88eSHong Zhang   w3 = coloring->w3;
597efb30889SBarry Smith 
598b8f8c88eSHong Zhang     /* Compute all the local scale factors, including ghost points */
599b8f8c88eSHong Zhang   ierr = VecGetLocalSize(x1_tmp,&N);CHKERRQ(ierr);
600b8f8c88eSHong Zhang   ierr = VecGetArray(x1_tmp,&xx);CHKERRQ(ierr);
601b8f8c88eSHong Zhang   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
602b8f8c88eSHong Zhang   if (ctype == IS_COLORING_GHOSTED){
603b8f8c88eSHong Zhang     col_start = 0; col_end = N;
6048ee2e534SBarry Smith   } else if (ctype == IS_COLORING_GLOBAL){
605b8f8c88eSHong Zhang     xx = xx - start;
606b8f8c88eSHong Zhang     vscale_array = vscale_array - start;
607b8f8c88eSHong Zhang     col_start = start; col_end = N + start;
608b8f8c88eSHong Zhang   }
609b8f8c88eSHong Zhang   for (col=col_start; col<col_end; col++){
610b8f8c88eSHong Zhang     /* Loop over each local column, vscale[col] = 1./(epsilon*dx[col]) */
611efb30889SBarry Smith     if (coloring->htype[0] == 'w') {
612efb30889SBarry Smith       dx = 1.0 + unorm;
613efb30889SBarry Smith     } else {
6149782cf98SBarry Smith       dx  = xx[col];
615efb30889SBarry Smith     }
6165b8514ebSBarry Smith     if (dx == 0.0) dx = 1.0;
6179782cf98SBarry Smith #if !defined(PETSC_USE_COMPLEX)
6189782cf98SBarry Smith     if (dx < umin && dx >= 0.0)      dx = umin;
6199782cf98SBarry Smith     else if (dx < 0.0 && dx > -umin) dx = -umin;
6209782cf98SBarry Smith #else
6219782cf98SBarry Smith     if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
6229782cf98SBarry Smith     else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
6239782cf98SBarry Smith #endif
6249782cf98SBarry Smith     dx               *= epsilon;
62530b34957SBarry Smith     vscale_array[col] = 1.0/dx;
6269782cf98SBarry Smith   }
6278ee2e534SBarry Smith   if (ctype == IS_COLORING_GLOBAL)  vscale_array = vscale_array + start;
628efb30889SBarry Smith   ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
6298ee2e534SBarry Smith   if (ctype == IS_COLORING_GLOBAL){
63030b34957SBarry Smith     ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
63130b34957SBarry Smith     ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
632b8f8c88eSHong Zhang   }
6339782cf98SBarry Smith 
634b8f8c88eSHong Zhang   if (coloring->vscaleforrow) {
635b8f8c88eSHong Zhang     vscaleforrow = coloring->vscaleforrow;
636b8f8c88eSHong Zhang   } else {
637b8f8c88eSHong Zhang     SETERRQ(PETSC_ERR_ARG_NULL,"Null Object: coloring->vscaleforrow");
638b8f8c88eSHong Zhang   }
639b0a32e0cSBarry Smith 
6409782cf98SBarry Smith   /*
64143a90d84SBarry Smith     Loop over each color
64243a90d84SBarry Smith   */
643b8f8c88eSHong Zhang   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
64443a90d84SBarry Smith   for (k=0; k<coloring->ncolors; k++) {
64549b058dcSBarry Smith     coloring->currentcolor = k;
646b8f8c88eSHong Zhang     ierr = VecCopy(x1_tmp,w3);CHKERRQ(ierr);
647b8f8c88eSHong Zhang     ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);
6488ee2e534SBarry Smith     if (ctype == IS_COLORING_GLOBAL) w3_array = w3_array - start;
64943a90d84SBarry Smith     /*
650b8f8c88eSHong Zhang       Loop over each column associated with color
651b8f8c88eSHong Zhang       adding the perturbation to the vector w3.
65243a90d84SBarry Smith     */
65343a90d84SBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
654b8f8c88eSHong Zhang       col = coloring->columns[k][l];    /* local column of the matrix we are probing for */
655efb30889SBarry Smith       if (coloring->htype[0] == 'w') {
656efb30889SBarry Smith         dx = 1.0 + unorm;
657efb30889SBarry Smith       } else {
65842460c72SBarry Smith         dx  = xx[col];
659efb30889SBarry Smith       }
6605b8514ebSBarry Smith       if (dx == 0.0) dx = 1.0;
661aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
662ae09f205SBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
663ae09f205SBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
66443a90d84SBarry Smith #else
665329f5518SBarry Smith       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
666329f5518SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
66743a90d84SBarry Smith #endif
66843a90d84SBarry Smith       dx            *= epsilon;
6691302d50aSBarry Smith       if (!PetscAbsScalar(dx)) SETERRQ(PETSC_ERR_PLIB,"Computed 0 differencing parameter");
67042460c72SBarry Smith       w3_array[col] += dx;
67143a90d84SBarry Smith     }
6728ee2e534SBarry Smith     if (ctype == IS_COLORING_GLOBAL) w3_array = w3_array + start;
673b8f8c88eSHong Zhang     ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
6743b28642cSBarry Smith 
67543a90d84SBarry Smith     /*
676b8f8c88eSHong Zhang       Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations)
677b8f8c88eSHong Zhang                            w2 = F(x1 + dx) - F(x1)
67843a90d84SBarry Smith     */
67966f9b7ceSBarry Smith     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
68043a90d84SBarry Smith     ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr);
68166f9b7ceSBarry Smith     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
682efb30889SBarry Smith     ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr);
6839782cf98SBarry Smith 
68443a90d84SBarry Smith     /*
685e0907662SLois Curfman McInnes       Loop over rows of vector, putting results into Jacobian matrix
68643a90d84SBarry Smith     */
6873b28642cSBarry Smith     ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
68843a90d84SBarry Smith     for (l=0; l<coloring->nrows[k]; l++) {
689b8f8c88eSHong Zhang       row    = coloring->rows[k][l];             /* local row index */
690b8f8c88eSHong Zhang       col    = coloring->columnsforrow[k][l];    /* global column index */
6915904e1b1SBarry Smith       y[row] *= vscale_array[vscaleforrow[k][l]];
69243a90d84SBarry Smith       srow   = row + start;
69343a90d84SBarry Smith       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
69443a90d84SBarry Smith     }
6953b28642cSBarry Smith     ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
696aeb7e98aSMatthew Knepley   } /* endof for each color */
6978ee2e534SBarry Smith   if (ctype == IS_COLORING_GLOBAL) xx = xx + start;
69830b34957SBarry Smith   ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
699b8f8c88eSHong Zhang   ierr = VecRestoreArray(x1_tmp,&xx);CHKERRQ(ierr);
700b8f8c88eSHong Zhang 
701b8f8c88eSHong Zhang   coloring->currentcolor = -1;
70243a90d84SBarry Smith   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
70343a90d84SBarry Smith   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
704d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
705a2e34c3dSBarry Smith 
70690d69ab7SBarry Smith   flg  = PETSC_FALSE;
70790d69ab7SBarry Smith   ierr = PetscOptionsGetTruth(PETSC_NULL,"-mat_null_space_test",&flg,PETSC_NULL);CHKERRQ(ierr);
708a2e34c3dSBarry Smith   if (flg) {
70995902228SMatthew Knepley     ierr = MatNullSpaceTest(J->nullsp,J,PETSC_NULL);CHKERRQ(ierr);
710a2e34c3dSBarry Smith   }
711b9722508SLois Curfman McInnes   ierr = MatFDColoringView_Private(coloring);CHKERRQ(ierr);
7123a40ed3dSBarry Smith   PetscFunctionReturn(0);
71343a90d84SBarry Smith }
714840b8ebdSBarry Smith 
7154a2ae208SSatish Balay #undef __FUNCT__
7164a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringApplyTS"
717840b8ebdSBarry Smith /*@
718840b8ebdSBarry Smith     MatFDColoringApplyTS - Given a matrix for which a MatFDColoring context
719840b8ebdSBarry Smith     has been created, computes the Jacobian for a function via finite differences.
720840b8ebdSBarry Smith 
721fee21e36SBarry Smith    Collective on Mat, MatFDColoring, and Vec
722fee21e36SBarry Smith 
723ef5ee4d1SLois Curfman McInnes     Input Parameters:
7243b28642cSBarry Smith +   mat - location to store Jacobian
725ef5ee4d1SLois Curfman McInnes .   coloring - coloring context created with MatFDColoringCreate()
726ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
7277850c7c0SBarry Smith -   sctx - context required by function, if this is being used with the TS solver then it is TS object, otherwise it is null
728ef5ee4d1SLois Curfman McInnes 
72915091d37SBarry Smith    Level: intermediate
73015091d37SBarry Smith 
7317850c7c0SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringSetFunction()
732840b8ebdSBarry Smith 
733840b8ebdSBarry Smith .keywords: coloring, Jacobian, finite differences
734840b8ebdSBarry Smith @*/
735be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApplyTS(Mat J,MatFDColoring coloring,PetscReal t,Vec x1,MatStructure *flag,void *sctx)
736840b8ebdSBarry Smith {
7376849ba73SBarry Smith   PetscErrorCode (*f)(void*,PetscReal,Vec,Vec,void*)=(PetscErrorCode (*)(void*,PetscReal,Vec,Vec,void *))coloring->f;
7386849ba73SBarry Smith   PetscErrorCode ierr;
73938baddfdSBarry Smith   PetscInt       k,N,start,end,l,row,col,srow,**vscaleforrow;
740efb30889SBarry Smith   PetscScalar    dx,*y,*xx,*w3_array;
741ea709b57SSatish Balay   PetscScalar    *vscale_array;
742329f5518SBarry Smith   PetscReal      epsilon = coloring->error_rel,umin = coloring->umin;
743ced766e8SHong Zhang   Vec            w1=coloring->w1,w2=coloring->w2,w3;
744840b8ebdSBarry Smith   void           *fctx = coloring->fctx;
745f1af5d2fSBarry Smith   PetscTruth     flg;
746840b8ebdSBarry Smith 
7473a40ed3dSBarry Smith   PetscFunctionBegin;
748*0700a824SBarry Smith   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
749*0700a824SBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_CLASSID,2);
750*0700a824SBarry Smith   PetscValidHeaderSpecific(x1,VEC_CLASSID,4);
751840b8ebdSBarry Smith 
752d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
753ced766e8SHong Zhang   if (!coloring->w3) {
754840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
75552e6d16bSBarry Smith     ierr = PetscLogObjectParent(coloring,coloring->w3);CHKERRQ(ierr);
756840b8ebdSBarry Smith   }
757ced766e8SHong Zhang   w3 = coloring->w3;
758840b8ebdSBarry Smith 
7595904e1b1SBarry Smith   ierr = MatSetUnfactored(J);CHKERRQ(ierr);
76090d69ab7SBarry Smith   flg  = PETSC_FALSE;
76190d69ab7SBarry Smith   ierr = PetscOptionsGetTruth(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg,PETSC_NULL);CHKERRQ(ierr);
762f1af5d2fSBarry Smith   if (flg) {
763ae15b995SBarry Smith     ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr);
764840b8ebdSBarry Smith   } else {
7650b9b6f31SBarry Smith     PetscTruth assembled;
7660b9b6f31SBarry Smith     ierr = MatAssembled(J,&assembled);CHKERRQ(ierr);
7670b9b6f31SBarry Smith     if (assembled) {
768840b8ebdSBarry Smith       ierr = MatZeroEntries(J);CHKERRQ(ierr);
769840b8ebdSBarry Smith     }
7700b9b6f31SBarry Smith   }
771840b8ebdSBarry Smith 
772840b8ebdSBarry Smith   ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr);
773840b8ebdSBarry Smith   ierr = VecGetSize(x1,&N);CHKERRQ(ierr);
77466f9b7ceSBarry Smith   ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
775840b8ebdSBarry Smith   ierr = (*f)(sctx,t,x1,w1,fctx);CHKERRQ(ierr);
77666f9b7ceSBarry Smith   ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
777840b8ebdSBarry Smith 
778840b8ebdSBarry Smith   /*
7795904e1b1SBarry Smith       Compute all the scale factors and share with other processors
780840b8ebdSBarry Smith   */
7815904e1b1SBarry Smith   ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);xx = xx - start;
7825904e1b1SBarry Smith   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);vscale_array = vscale_array - start;
783840b8ebdSBarry Smith   for (k=0; k<coloring->ncolors; k++) {
784840b8ebdSBarry Smith     /*
785840b8ebdSBarry Smith        Loop over each column associated with color adding the
786840b8ebdSBarry Smith        perturbation to the vector w3.
787840b8ebdSBarry Smith     */
788840b8ebdSBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
789840b8ebdSBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
7905904e1b1SBarry Smith       dx  = xx[col];
7915b8514ebSBarry Smith       if (dx == 0.0) dx = 1.0;
792aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
793840b8ebdSBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
794840b8ebdSBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
795840b8ebdSBarry Smith #else
796329f5518SBarry Smith       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
797329f5518SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
798840b8ebdSBarry Smith #endif
799840b8ebdSBarry Smith       dx                *= epsilon;
8005904e1b1SBarry Smith       vscale_array[col] = 1.0/dx;
801840b8ebdSBarry Smith     }
8025904e1b1SBarry Smith   }
8035904e1b1SBarry Smith   vscale_array = vscale_array - start;ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
8045904e1b1SBarry Smith   ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8055904e1b1SBarry Smith   ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8065904e1b1SBarry Smith 
8075904e1b1SBarry Smith   if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow;
8085904e1b1SBarry Smith   else                        vscaleforrow = coloring->columnsforrow;
8095904e1b1SBarry Smith 
8105904e1b1SBarry Smith   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
8115904e1b1SBarry Smith   /*
8125904e1b1SBarry Smith       Loop over each color
8135904e1b1SBarry Smith   */
8145904e1b1SBarry Smith   for (k=0; k<coloring->ncolors; k++) {
8155904e1b1SBarry Smith     ierr = VecCopy(x1,w3);CHKERRQ(ierr);
8165904e1b1SBarry Smith     ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);w3_array = w3_array - start;
8175904e1b1SBarry Smith     /*
8185904e1b1SBarry Smith        Loop over each column associated with color adding the
8195904e1b1SBarry Smith        perturbation to the vector w3.
8205904e1b1SBarry Smith     */
8215904e1b1SBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
8225904e1b1SBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
8235904e1b1SBarry Smith       dx  = xx[col];
8245b8514ebSBarry Smith       if (dx == 0.0) dx = 1.0;
8255904e1b1SBarry Smith #if !defined(PETSC_USE_COMPLEX)
8265904e1b1SBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
8275904e1b1SBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
8285904e1b1SBarry Smith #else
8295904e1b1SBarry Smith       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
8305904e1b1SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
8315904e1b1SBarry Smith #endif
8325904e1b1SBarry Smith       dx            *= epsilon;
8335904e1b1SBarry Smith       w3_array[col] += dx;
8345904e1b1SBarry Smith     }
8355904e1b1SBarry Smith     w3_array = w3_array + start; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
8365904e1b1SBarry Smith 
837840b8ebdSBarry Smith     /*
838840b8ebdSBarry Smith        Evaluate function at x1 + dx (here dx is a vector of perturbations)
839840b8ebdSBarry Smith     */
84066f9b7ceSBarry Smith     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
841840b8ebdSBarry Smith     ierr = (*f)(sctx,t,w3,w2,fctx);CHKERRQ(ierr);
84266f9b7ceSBarry Smith     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
843efb30889SBarry Smith     ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr);
8445904e1b1SBarry Smith 
845840b8ebdSBarry Smith     /*
846840b8ebdSBarry Smith        Loop over rows of vector, putting results into Jacobian matrix
847840b8ebdSBarry Smith     */
8483b28642cSBarry Smith     ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
849840b8ebdSBarry Smith     for (l=0; l<coloring->nrows[k]; l++) {
850840b8ebdSBarry Smith       row    = coloring->rows[k][l];
851840b8ebdSBarry Smith       col    = coloring->columnsforrow[k][l];
8525904e1b1SBarry Smith       y[row] *= vscale_array[vscaleforrow[k][l]];
853840b8ebdSBarry Smith       srow   = row + start;
854840b8ebdSBarry Smith       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
855840b8ebdSBarry Smith     }
8563b28642cSBarry Smith     ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
857840b8ebdSBarry Smith   }
8585904e1b1SBarry Smith   ierr  = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
8595904e1b1SBarry Smith   xx    = xx + start; ierr  = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
860840b8ebdSBarry Smith   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
861840b8ebdSBarry Smith   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
862d5ba7fb7SMatthew Knepley   ierr  = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
8633a40ed3dSBarry Smith   PetscFunctionReturn(0);
864840b8ebdSBarry Smith }
8653b28642cSBarry Smith 
8663b28642cSBarry Smith 
8673b28642cSBarry Smith 
868