xref: /petsc/src/mat/matfd/fdmatrix.c (revision 5d973c197c8b7d308ce7b868381cbe50f94f9f77)
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;
1030700a824SBarry Smith   PetscValidHeaderSpecific(c,MAT_FDCOLORING_CLASSID,1);
1043050cee2SBarry Smith   if (!viewer) {
1057adad957SLisandro Dalcin     ierr = PetscViewerASCIIGetStdout(((PetscObject)c)->comm,&viewer);CHKERRQ(ierr);
1063050cee2SBarry Smith   }
1070700a824SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
108c9780b6fSBarry Smith   PetscCheckSameComm(c,1,viewer,2);
109bbf0e169SBarry Smith 
1102692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
1112692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&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 {
136e32f2f54SBarry Smith     SETERRQ1(PETSC_COMM_SELF,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;
1720700a824SBarry 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;
2050700a824SBarry 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;
2470700a824SBarry 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;
2950700a824SBarry 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);
3001bc75a8dSBarry Smith     ierr = PetscOptionsString("-mat_fd_type","Algorithm to compute h, wp or ds","MatFDColoringCreate",matfd->htype,value,3,&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";
304e32f2f54SBarry Smith       else SETERRQ1(PETSC_COMM_SELF,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);
310*5d973c19SBarry Smith 
311*5d973c19SBarry Smith     /* process any options handlers added with PetscObjectAddOptionsHandler() */
312*5d973c19SBarry Smith     ierr = PetscObjectProcessOptionsHandlers((PetscObject)matfd);CHKERRQ(ierr);
313b0a32e0cSBarry Smith   PetscOptionsEnd();CHKERRQ(ierr);
3143a40ed3dSBarry Smith   PetscFunctionReturn(0);
315005c665bSBarry Smith }
316005c665bSBarry Smith 
3174a2ae208SSatish Balay #undef __FUNCT__
3184a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Private"
319dfbe8321SBarry Smith PetscErrorCode MatFDColoringView_Private(MatFDColoring fd)
320005c665bSBarry Smith {
321dfbe8321SBarry Smith   PetscErrorCode ierr;
32290d69ab7SBarry Smith   PetscTruth     flg = PETSC_FALSE;
3233050cee2SBarry Smith   PetscViewer    viewer;
324005c665bSBarry Smith 
3253a40ed3dSBarry Smith   PetscFunctionBegin;
3267adad957SLisandro Dalcin   ierr = PetscViewerASCIIGetStdout(((PetscObject)fd)->comm,&viewer);CHKERRQ(ierr);
32790d69ab7SBarry Smith   ierr = PetscOptionsGetTruth(PETSC_NULL,"-mat_fd_coloring_view",&flg,PETSC_NULL);CHKERRQ(ierr);
328005c665bSBarry Smith   if (flg) {
3293050cee2SBarry Smith     ierr = MatFDColoringView(fd,viewer);CHKERRQ(ierr);
330005c665bSBarry Smith   }
33190d69ab7SBarry Smith   flg  = PETSC_FALSE;
33290d69ab7SBarry Smith   ierr = PetscOptionsGetTruth(PETSC_NULL,"-mat_fd_coloring_view_info",&flg,PETSC_NULL);CHKERRQ(ierr);
333ae09f205SBarry Smith   if (flg) {
3343050cee2SBarry Smith     ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
3353050cee2SBarry Smith     ierr = MatFDColoringView(fd,viewer);CHKERRQ(ierr);
3363050cee2SBarry Smith     ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
337ae09f205SBarry Smith   }
33890d69ab7SBarry Smith   flg  = PETSC_FALSE;
33990d69ab7SBarry Smith   ierr = PetscOptionsGetTruth(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg,PETSC_NULL);CHKERRQ(ierr);
340005c665bSBarry Smith   if (flg) {
3417adad957SLisandro Dalcin     ierr = MatFDColoringView(fd,PETSC_VIEWER_DRAW_(((PetscObject)fd)->comm));CHKERRQ(ierr);
3427adad957SLisandro Dalcin     ierr = PetscViewerFlush(PETSC_VIEWER_DRAW_(((PetscObject)fd)->comm));CHKERRQ(ierr);
343005c665bSBarry Smith   }
3443a40ed3dSBarry Smith   PetscFunctionReturn(0);
345bbf0e169SBarry Smith }
346bbf0e169SBarry Smith 
3474a2ae208SSatish Balay #undef __FUNCT__
3484a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringCreate"
34905869f15SSatish Balay /*@
350639f9d9dSBarry Smith    MatFDColoringCreate - Creates a matrix coloring context for finite difference
351639f9d9dSBarry Smith    computation of Jacobians.
352bbf0e169SBarry Smith 
353ef5ee4d1SLois Curfman McInnes    Collective on Mat
354ef5ee4d1SLois Curfman McInnes 
355639f9d9dSBarry Smith    Input Parameters:
356ef5ee4d1SLois Curfman McInnes +  mat - the matrix containing the nonzero structure of the Jacobian
3571d0fab5eSBarry Smith -  iscoloring - the coloring of the matrix; usually obtained with MatGetColoring() or DAGetColoring()
358bbf0e169SBarry Smith 
359bbf0e169SBarry Smith     Output Parameter:
360639f9d9dSBarry Smith .   color - the new coloring context
361bbf0e169SBarry Smith 
36215091d37SBarry Smith     Level: intermediate
36315091d37SBarry Smith 
3646831982aSBarry Smith .seealso: MatFDColoringDestroy(),SNESDefaultComputeJacobianColor(), ISColoringCreate(),
365d1c39f55SBarry Smith           MatFDColoringSetFunction(), MatFDColoringSetFromOptions(), MatFDColoringApply(),
3661d0fab5eSBarry Smith           MatFDColoringView(), MatFDColoringSetParameters(), MatGetColoring(), DAGetColoring()
367bbf0e169SBarry Smith @*/
368be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
369bbf0e169SBarry Smith {
370639f9d9dSBarry Smith   MatFDColoring  c;
371639f9d9dSBarry Smith   MPI_Comm       comm;
372dfbe8321SBarry Smith   PetscErrorCode ierr;
37338baddfdSBarry Smith   PetscInt       M,N;
374639f9d9dSBarry Smith 
3753a40ed3dSBarry Smith   PetscFunctionBegin;
376d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
377639f9d9dSBarry Smith   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
37817186662SBarry Smith   if (M != N) SETERRQ(((PetscObject)mat)->comm,PETSC_ERR_SUP,"Only for square matrices");
379639f9d9dSBarry Smith 
380f881d145SBarry Smith   ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
3810700a824SBarry Smith   ierr = PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_CLASSID,0,"MatFDColoring",comm,MatFDColoringDestroy,MatFDColoringView);CHKERRQ(ierr);
382639f9d9dSBarry Smith 
383b8f8c88eSHong Zhang   c->ctype = iscoloring->ctype;
384b8f8c88eSHong Zhang 
385f830108cSBarry Smith   if (mat->ops->fdcoloringcreate) {
386f830108cSBarry Smith     ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);CHKERRQ(ierr);
38717186662SBarry Smith   } else SETERRQ(((PetscObject)mat)->comm,PETSC_ERR_SUP,"Code not yet written for this matrix type");
388639f9d9dSBarry Smith 
389b8f8c88eSHong Zhang   ierr = MatGetVecs(mat,PETSC_NULL,&c->w1);CHKERRQ(ierr);
390b8f8c88eSHong Zhang   ierr = PetscLogObjectParent(c,c->w1);CHKERRQ(ierr);
391b8f8c88eSHong Zhang   ierr = VecDuplicate(c->w1,&c->w2);CHKERRQ(ierr);
392b8f8c88eSHong Zhang   ierr = PetscLogObjectParent(c,c->w2);CHKERRQ(ierr);
393b8f8c88eSHong Zhang 
39477d8c4bbSBarry Smith   c->error_rel         = PETSC_SQRT_MACHINE_EPSILON;
39577d8c4bbSBarry Smith   c->umin              = 100.0*PETSC_SQRT_MACHINE_EPSILON;
39649b058dcSBarry Smith   c->currentcolor      = -1;
397efb30889SBarry Smith   c->htype             = "wp";
398639f9d9dSBarry Smith 
399639f9d9dSBarry Smith   *color = c;
400d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
4013a40ed3dSBarry Smith   PetscFunctionReturn(0);
402bbf0e169SBarry Smith }
403bbf0e169SBarry Smith 
4044a2ae208SSatish Balay #undef __FUNCT__
4054a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringDestroy"
40605869f15SSatish Balay /*@
407639f9d9dSBarry Smith     MatFDColoringDestroy - Destroys a matrix coloring context that was created
408639f9d9dSBarry Smith     via MatFDColoringCreate().
409bbf0e169SBarry Smith 
410ef5ee4d1SLois Curfman McInnes     Collective on MatFDColoring
411ef5ee4d1SLois Curfman McInnes 
412b4fc646aSLois Curfman McInnes     Input Parameter:
413639f9d9dSBarry Smith .   c - coloring context
414bbf0e169SBarry Smith 
41515091d37SBarry Smith     Level: intermediate
41615091d37SBarry Smith 
417639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
418bbf0e169SBarry Smith @*/
419be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringDestroy(MatFDColoring c)
420bbf0e169SBarry Smith {
4216849ba73SBarry Smith   PetscErrorCode ierr;
42238baddfdSBarry Smith   PetscInt       i;
423bbf0e169SBarry Smith 
4243a40ed3dSBarry Smith   PetscFunctionBegin;
4257adad957SLisandro Dalcin   if (--((PetscObject)c)->refct > 0) PetscFunctionReturn(0);
426d4bb536fSBarry Smith 
427639f9d9dSBarry Smith   for (i=0; i<c->ncolors; i++) {
42805b42c5fSBarry Smith     ierr = PetscFree(c->columns[i]);CHKERRQ(ierr);
42905b42c5fSBarry Smith     ierr = PetscFree(c->rows[i]);CHKERRQ(ierr);
43005b42c5fSBarry Smith     ierr = PetscFree(c->columnsforrow[i]);CHKERRQ(ierr);
43105b42c5fSBarry Smith     if (c->vscaleforrow) {ierr = PetscFree(c->vscaleforrow[i]);CHKERRQ(ierr);}
432bbf0e169SBarry Smith   }
433606d414cSSatish Balay   ierr = PetscFree(c->ncolumns);CHKERRQ(ierr);
434606d414cSSatish Balay   ierr = PetscFree(c->columns);CHKERRQ(ierr);
435606d414cSSatish Balay   ierr = PetscFree(c->nrows);CHKERRQ(ierr);
436606d414cSSatish Balay   ierr = PetscFree(c->rows);CHKERRQ(ierr);
437606d414cSSatish Balay   ierr = PetscFree(c->columnsforrow);CHKERRQ(ierr);
43805b42c5fSBarry Smith   ierr = PetscFree(c->vscaleforrow);CHKERRQ(ierr);
4394f113abfSBarry Smith   if (c->vscale)       {ierr = VecDestroy(c->vscale);CHKERRQ(ierr);}
440005c665bSBarry Smith   if (c->w1) {
441005c665bSBarry Smith     ierr = VecDestroy(c->w1);CHKERRQ(ierr);
442005c665bSBarry Smith     ierr = VecDestroy(c->w2);CHKERRQ(ierr);
443b8f8c88eSHong Zhang   }
444b8f8c88eSHong Zhang   if (c->w3){
445005c665bSBarry Smith     ierr = VecDestroy(c->w3);CHKERRQ(ierr);
446005c665bSBarry Smith   }
447d38fa0fbSBarry Smith   ierr = PetscHeaderDestroy(c);CHKERRQ(ierr);
4483a40ed3dSBarry Smith   PetscFunctionReturn(0);
449bbf0e169SBarry Smith }
45043a90d84SBarry Smith 
4514a2ae208SSatish Balay #undef __FUNCT__
45249b058dcSBarry Smith #define __FUNCT__ "MatFDColoringGetPerturbedColumns"
45349b058dcSBarry Smith /*@C
45449b058dcSBarry Smith     MatFDColoringGetPerturbedColumns - Returns the indices of the columns that
45549b058dcSBarry Smith       that are currently being perturbed.
45649b058dcSBarry Smith 
45749b058dcSBarry Smith     Not Collective
45849b058dcSBarry Smith 
45949b058dcSBarry Smith     Input Parameters:
46049b058dcSBarry Smith .   coloring - coloring context created with MatFDColoringCreate()
46149b058dcSBarry Smith 
46249b058dcSBarry Smith     Output Parameters:
46349b058dcSBarry Smith +   n - the number of local columns being perturbed
46449b058dcSBarry Smith -   cols - the column indices, in global numbering
46549b058dcSBarry Smith 
46649b058dcSBarry Smith    Level: intermediate
46749b058dcSBarry Smith 
46849b058dcSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringApply()
46949b058dcSBarry Smith 
47049b058dcSBarry Smith .keywords: coloring, Jacobian, finite differences
47149b058dcSBarry Smith @*/
472be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetPerturbedColumns(MatFDColoring coloring,PetscInt *n,PetscInt *cols[])
47349b058dcSBarry Smith {
47449b058dcSBarry Smith   PetscFunctionBegin;
47549b058dcSBarry Smith   if (coloring->currentcolor >= 0) {
47649b058dcSBarry Smith     *n    = coloring->ncolumns[coloring->currentcolor];
47749b058dcSBarry Smith     *cols = coloring->columns[coloring->currentcolor];
47849b058dcSBarry Smith   } else {
47949b058dcSBarry Smith     *n = 0;
48049b058dcSBarry Smith   }
48149b058dcSBarry Smith   PetscFunctionReturn(0);
48249b058dcSBarry Smith }
48349b058dcSBarry Smith 
48449b058dcSBarry Smith #undef __FUNCT__
4854a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringApply"
48643a90d84SBarry Smith /*@
487e0907662SLois Curfman McInnes     MatFDColoringApply - Given a matrix for which a MatFDColoring context
488e0907662SLois Curfman McInnes     has been created, computes the Jacobian for a function via finite differences.
48943a90d84SBarry Smith 
490fee21e36SBarry Smith     Collective on MatFDColoring
491fee21e36SBarry Smith 
492ef5ee4d1SLois Curfman McInnes     Input Parameters:
493ef5ee4d1SLois Curfman McInnes +   mat - location to store Jacobian
494ef5ee4d1SLois Curfman McInnes .   coloring - coloring context created with MatFDColoringCreate()
495ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
4967850c7c0SBarry Smith -   sctx - context required by function, if this is being used with the SNES solver then it is SNES object, otherwise it is null
497ef5ee4d1SLois Curfman McInnes 
4988bba8e72SBarry Smith     Options Database Keys:
499ebd3b9afSBarry Smith +    -mat_fd_type - "wp" or "ds"  (see MATMFFD_WP or MATMFFD_DS)
500b9722508SLois Curfman McInnes .    -mat_fd_coloring_view - Activates basic viewing or coloring
501b9722508SLois Curfman McInnes .    -mat_fd_coloring_view_draw - Activates drawing of coloring
502b9722508SLois Curfman McInnes -    -mat_fd_coloring_view_info - Activates viewing of coloring info
5038bba8e72SBarry Smith 
50415091d37SBarry Smith     Level: intermediate
50515091d37SBarry Smith 
5067850c7c0SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringSetFunction()
50743a90d84SBarry Smith 
50843a90d84SBarry Smith .keywords: coloring, Jacobian, finite differences
50943a90d84SBarry Smith @*/
510be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
51143a90d84SBarry Smith {
5123acb8795SBarry Smith   PetscErrorCode ierr;
5133acb8795SBarry Smith 
5143acb8795SBarry Smith   PetscFunctionBegin;
5150700a824SBarry Smith   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
5160700a824SBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_CLASSID,2);
5170700a824SBarry Smith   PetscValidHeaderSpecific(x1,VEC_CLASSID,3);
51817186662SBarry Smith   if (!coloring->f) SETERRQ(((PetscObject)J)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetFunction()");
519e32f2f54SBarry Smith   if (!J->ops->fdcoloringapply) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not supported for this matrix type %s",((PetscObject)J)->type_name);
5203acb8795SBarry Smith   ierr = (*J->ops->fdcoloringapply)(J,coloring,x1,flag,sctx);CHKERRQ(ierr);
5213acb8795SBarry Smith   PetscFunctionReturn(0);
5223acb8795SBarry Smith }
5233acb8795SBarry Smith 
5243acb8795SBarry Smith #undef __FUNCT__
5253acb8795SBarry Smith #define __FUNCT__ "MatFDColoringApply_AIJ"
5263acb8795SBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApply_AIJ(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
5273acb8795SBarry Smith {
5286849ba73SBarry Smith   PetscErrorCode (*f)(void*,Vec,Vec,void*) = (PetscErrorCode (*)(void*,Vec,Vec,void *))coloring->f;
5296849ba73SBarry Smith   PetscErrorCode ierr;
530b8f8c88eSHong Zhang   PetscInt       k,start,end,l,row,col,srow,**vscaleforrow,m1,m2;
531efb30889SBarry Smith   PetscScalar    dx,*y,*xx,*w3_array;
532ea709b57SSatish Balay   PetscScalar    *vscale_array;
533efb30889SBarry Smith   PetscReal      epsilon = coloring->error_rel,umin = coloring->umin,unorm;
534b8f8c88eSHong Zhang   Vec            w1=coloring->w1,w2=coloring->w2,w3;
535005c665bSBarry Smith   void           *fctx = coloring->fctx;
53690d69ab7SBarry Smith   PetscTruth     flg = PETSC_FALSE;
537705d48f7SSatish Balay   PetscInt       ctype=coloring->ctype,N,col_start=0,col_end=0;
538b8f8c88eSHong Zhang   Vec            x1_tmp;
539a2e34c3dSBarry Smith 
5403a40ed3dSBarry Smith   PetscFunctionBegin;
5410700a824SBarry Smith   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
5420700a824SBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_CLASSID,2);
5430700a824SBarry Smith   PetscValidHeaderSpecific(x1,VEC_CLASSID,3);
54417186662SBarry Smith   if (!f) SETERRQ(((PetscObject)J)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetFunction()");
545e0907662SLois Curfman McInnes 
546d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
5474c49b128SBarry Smith   ierr = MatSetUnfactored(J);CHKERRQ(ierr);
54890d69ab7SBarry Smith   ierr = PetscOptionsGetTruth(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg,PETSC_NULL);CHKERRQ(ierr);
549f1af5d2fSBarry Smith   if (flg) {
550ae15b995SBarry Smith     ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr);
551e0907662SLois Curfman McInnes   } else {
5520b9b6f31SBarry Smith     PetscTruth assembled;
5530b9b6f31SBarry Smith     ierr = MatAssembled(J,&assembled);CHKERRQ(ierr);
5540b9b6f31SBarry Smith     if (assembled) {
55543a90d84SBarry Smith       ierr = MatZeroEntries(J);CHKERRQ(ierr);
556e0907662SLois Curfman McInnes     }
5570b9b6f31SBarry Smith   }
55843a90d84SBarry Smith 
559b8f8c88eSHong Zhang   x1_tmp = x1;
560dac7f5fdSBarry Smith   if (!coloring->vscale){
561b8f8c88eSHong Zhang     ierr = VecDuplicate(x1_tmp,&coloring->vscale);CHKERRQ(ierr);
562b8f8c88eSHong Zhang   }
563be2fbe1fSBarry Smith 
5643a7fca6bSBarry Smith   /*
5653a7fca6bSBarry Smith     This is a horrible, horrible, hack. See DMMGComputeJacobian_Multigrid() it inproperly sets
5663a7fca6bSBarry Smith     coloring->F for the coarser grids from the finest
5673a7fca6bSBarry Smith   */
5683a7fca6bSBarry Smith   if (coloring->F) {
5693a7fca6bSBarry Smith     ierr = VecGetLocalSize(coloring->F,&m1);CHKERRQ(ierr);
5703a7fca6bSBarry Smith     ierr = VecGetLocalSize(w1,&m2);CHKERRQ(ierr);
5713a7fca6bSBarry Smith     if (m1 != m2) {
5723a7fca6bSBarry Smith       coloring->F = 0;
5733a7fca6bSBarry Smith       }
5743a7fca6bSBarry Smith     }
5753a7fca6bSBarry Smith 
576b8f8c88eSHong Zhang   if (coloring->htype[0] == 'w') { /* tacky test; need to make systematic if we add other approaches to computing h*/
577b8f8c88eSHong Zhang     ierr = VecNorm(x1_tmp,NORM_2,&unorm);CHKERRQ(ierr);
578b8f8c88eSHong Zhang   }
579b8f8c88eSHong Zhang   ierr = VecGetOwnershipRange(w1,&start,&end);CHKERRQ(ierr); /* OwnershipRange is used by ghosted x! */
580b8f8c88eSHong Zhang 
581b8f8c88eSHong Zhang   /* Set w1 = F(x1) */
582be2fbe1fSBarry Smith   if (coloring->F) {
583be2fbe1fSBarry Smith     w1          = coloring->F; /* use already computed value of function */
584be2fbe1fSBarry Smith     coloring->F = 0;
585be2fbe1fSBarry Smith   } else {
58666f9b7ceSBarry Smith     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
587b8f8c88eSHong Zhang     ierr = (*f)(sctx,x1_tmp,w1,fctx);CHKERRQ(ierr);
58866f9b7ceSBarry Smith     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
589be2fbe1fSBarry Smith   }
59043a90d84SBarry Smith 
591b8f8c88eSHong Zhang   if (!coloring->w3) {
592b8f8c88eSHong Zhang     ierr = VecDuplicate(x1_tmp,&coloring->w3);CHKERRQ(ierr);
593b8f8c88eSHong Zhang     ierr = PetscLogObjectParent(coloring,coloring->w3);CHKERRQ(ierr);
594efb30889SBarry Smith   }
595b8f8c88eSHong Zhang   w3 = coloring->w3;
596efb30889SBarry Smith 
597b8f8c88eSHong Zhang     /* Compute all the local scale factors, including ghost points */
598b8f8c88eSHong Zhang   ierr = VecGetLocalSize(x1_tmp,&N);CHKERRQ(ierr);
599b8f8c88eSHong Zhang   ierr = VecGetArray(x1_tmp,&xx);CHKERRQ(ierr);
600b8f8c88eSHong Zhang   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
601b8f8c88eSHong Zhang   if (ctype == IS_COLORING_GHOSTED){
602b8f8c88eSHong Zhang     col_start = 0; col_end = N;
6038ee2e534SBarry Smith   } else if (ctype == IS_COLORING_GLOBAL){
604b8f8c88eSHong Zhang     xx = xx - start;
605b8f8c88eSHong Zhang     vscale_array = vscale_array - start;
606b8f8c88eSHong Zhang     col_start = start; col_end = N + start;
607b8f8c88eSHong Zhang   }
608b8f8c88eSHong Zhang   for (col=col_start; col<col_end; col++){
609b8f8c88eSHong Zhang     /* Loop over each local column, vscale[col] = 1./(epsilon*dx[col]) */
610efb30889SBarry Smith     if (coloring->htype[0] == 'w') {
611efb30889SBarry Smith       dx = 1.0 + unorm;
612efb30889SBarry Smith     } else {
6139782cf98SBarry Smith       dx  = xx[col];
614efb30889SBarry Smith     }
6155b8514ebSBarry Smith     if (dx == 0.0) dx = 1.0;
6169782cf98SBarry Smith #if !defined(PETSC_USE_COMPLEX)
6179782cf98SBarry Smith     if (dx < umin && dx >= 0.0)      dx = umin;
6189782cf98SBarry Smith     else if (dx < 0.0 && dx > -umin) dx = -umin;
6199782cf98SBarry Smith #else
6209782cf98SBarry Smith     if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
6219782cf98SBarry Smith     else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
6229782cf98SBarry Smith #endif
6239782cf98SBarry Smith     dx               *= epsilon;
62430b34957SBarry Smith     vscale_array[col] = 1.0/dx;
6259782cf98SBarry Smith   }
6268ee2e534SBarry Smith   if (ctype == IS_COLORING_GLOBAL)  vscale_array = vscale_array + start;
627efb30889SBarry Smith   ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
6288ee2e534SBarry Smith   if (ctype == IS_COLORING_GLOBAL){
62930b34957SBarry Smith     ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
63030b34957SBarry Smith     ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
631b8f8c88eSHong Zhang   }
6329782cf98SBarry Smith 
633b8f8c88eSHong Zhang   if (coloring->vscaleforrow) {
634b8f8c88eSHong Zhang     vscaleforrow = coloring->vscaleforrow;
63517186662SBarry Smith   } else SETERRQ(((PetscObject)J)->comm,PETSC_ERR_ARG_NULL,"Null Object: coloring->vscaleforrow");
636b0a32e0cSBarry Smith 
6379782cf98SBarry Smith   /*
63843a90d84SBarry Smith     Loop over each color
63943a90d84SBarry Smith   */
640b8f8c88eSHong Zhang   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
64143a90d84SBarry Smith   for (k=0; k<coloring->ncolors; k++) {
64249b058dcSBarry Smith     coloring->currentcolor = k;
643b8f8c88eSHong Zhang     ierr = VecCopy(x1_tmp,w3);CHKERRQ(ierr);
644b8f8c88eSHong Zhang     ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);
6458ee2e534SBarry Smith     if (ctype == IS_COLORING_GLOBAL) w3_array = w3_array - start;
64643a90d84SBarry Smith     /*
647b8f8c88eSHong Zhang       Loop over each column associated with color
648b8f8c88eSHong Zhang       adding the perturbation to the vector w3.
64943a90d84SBarry Smith     */
65043a90d84SBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
651b8f8c88eSHong Zhang       col = coloring->columns[k][l];    /* local column of the matrix we are probing for */
652efb30889SBarry Smith       if (coloring->htype[0] == 'w') {
653efb30889SBarry Smith         dx = 1.0 + unorm;
654efb30889SBarry Smith       } else {
65542460c72SBarry Smith         dx  = xx[col];
656efb30889SBarry Smith       }
6575b8514ebSBarry Smith       if (dx == 0.0) dx = 1.0;
658aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
659ae09f205SBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
660ae09f205SBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
66143a90d84SBarry Smith #else
662329f5518SBarry Smith       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
663329f5518SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
66443a90d84SBarry Smith #endif
66543a90d84SBarry Smith       dx            *= epsilon;
666e32f2f54SBarry Smith       if (!PetscAbsScalar(dx)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Computed 0 differencing parameter");
66742460c72SBarry Smith       w3_array[col] += dx;
66843a90d84SBarry Smith     }
6698ee2e534SBarry Smith     if (ctype == IS_COLORING_GLOBAL) w3_array = w3_array + start;
670b8f8c88eSHong Zhang     ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
6713b28642cSBarry Smith 
67243a90d84SBarry Smith     /*
673b8f8c88eSHong Zhang       Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations)
674b8f8c88eSHong Zhang                            w2 = F(x1 + dx) - F(x1)
67543a90d84SBarry Smith     */
67666f9b7ceSBarry Smith     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
67743a90d84SBarry Smith     ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr);
67866f9b7ceSBarry Smith     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
679efb30889SBarry Smith     ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr);
6809782cf98SBarry Smith 
68143a90d84SBarry Smith     /*
682e0907662SLois Curfman McInnes       Loop over rows of vector, putting results into Jacobian matrix
68343a90d84SBarry Smith     */
6843b28642cSBarry Smith     ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
68543a90d84SBarry Smith     for (l=0; l<coloring->nrows[k]; l++) {
686b8f8c88eSHong Zhang       row    = coloring->rows[k][l];             /* local row index */
687b8f8c88eSHong Zhang       col    = coloring->columnsforrow[k][l];    /* global column index */
6885904e1b1SBarry Smith       y[row] *= vscale_array[vscaleforrow[k][l]];
68943a90d84SBarry Smith       srow   = row + start;
69043a90d84SBarry Smith       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
69143a90d84SBarry Smith     }
6923b28642cSBarry Smith     ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
693aeb7e98aSMatthew Knepley   } /* endof for each color */
6948ee2e534SBarry Smith   if (ctype == IS_COLORING_GLOBAL) xx = xx + start;
69530b34957SBarry Smith   ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
696b8f8c88eSHong Zhang   ierr = VecRestoreArray(x1_tmp,&xx);CHKERRQ(ierr);
697b8f8c88eSHong Zhang 
698b8f8c88eSHong Zhang   coloring->currentcolor = -1;
69943a90d84SBarry Smith   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
70043a90d84SBarry Smith   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
701d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
702a2e34c3dSBarry Smith 
70390d69ab7SBarry Smith   flg  = PETSC_FALSE;
70490d69ab7SBarry Smith   ierr = PetscOptionsGetTruth(PETSC_NULL,"-mat_null_space_test",&flg,PETSC_NULL);CHKERRQ(ierr);
705a2e34c3dSBarry Smith   if (flg) {
70695902228SMatthew Knepley     ierr = MatNullSpaceTest(J->nullsp,J,PETSC_NULL);CHKERRQ(ierr);
707a2e34c3dSBarry Smith   }
708b9722508SLois Curfman McInnes   ierr = MatFDColoringView_Private(coloring);CHKERRQ(ierr);
7093a40ed3dSBarry Smith   PetscFunctionReturn(0);
71043a90d84SBarry Smith }
711840b8ebdSBarry Smith 
7124a2ae208SSatish Balay #undef __FUNCT__
7134a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringApplyTS"
714840b8ebdSBarry Smith /*@
715840b8ebdSBarry Smith     MatFDColoringApplyTS - Given a matrix for which a MatFDColoring context
716840b8ebdSBarry Smith     has been created, computes the Jacobian for a function via finite differences.
717840b8ebdSBarry Smith 
718fee21e36SBarry Smith    Collective on Mat, MatFDColoring, and Vec
719fee21e36SBarry Smith 
720ef5ee4d1SLois Curfman McInnes     Input Parameters:
7213b28642cSBarry Smith +   mat - location to store Jacobian
722ef5ee4d1SLois Curfman McInnes .   coloring - coloring context created with MatFDColoringCreate()
723ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
7247850c7c0SBarry Smith -   sctx - context required by function, if this is being used with the TS solver then it is TS object, otherwise it is null
725ef5ee4d1SLois Curfman McInnes 
72615091d37SBarry Smith    Level: intermediate
72715091d37SBarry Smith 
7287850c7c0SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringSetFunction()
729840b8ebdSBarry Smith 
730840b8ebdSBarry Smith .keywords: coloring, Jacobian, finite differences
731840b8ebdSBarry Smith @*/
732be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApplyTS(Mat J,MatFDColoring coloring,PetscReal t,Vec x1,MatStructure *flag,void *sctx)
733840b8ebdSBarry Smith {
7346849ba73SBarry Smith   PetscErrorCode (*f)(void*,PetscReal,Vec,Vec,void*)=(PetscErrorCode (*)(void*,PetscReal,Vec,Vec,void *))coloring->f;
7356849ba73SBarry Smith   PetscErrorCode ierr;
73638baddfdSBarry Smith   PetscInt       k,N,start,end,l,row,col,srow,**vscaleforrow;
737efb30889SBarry Smith   PetscScalar    dx,*y,*xx,*w3_array;
738ea709b57SSatish Balay   PetscScalar    *vscale_array;
739329f5518SBarry Smith   PetscReal      epsilon = coloring->error_rel,umin = coloring->umin;
740ced766e8SHong Zhang   Vec            w1=coloring->w1,w2=coloring->w2,w3;
741840b8ebdSBarry Smith   void           *fctx = coloring->fctx;
742f1af5d2fSBarry Smith   PetscTruth     flg;
743840b8ebdSBarry Smith 
7443a40ed3dSBarry Smith   PetscFunctionBegin;
7450700a824SBarry Smith   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
7460700a824SBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_CLASSID,2);
7470700a824SBarry Smith   PetscValidHeaderSpecific(x1,VEC_CLASSID,4);
748840b8ebdSBarry Smith 
749d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
750ced766e8SHong Zhang   if (!coloring->w3) {
751840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
75252e6d16bSBarry Smith     ierr = PetscLogObjectParent(coloring,coloring->w3);CHKERRQ(ierr);
753840b8ebdSBarry Smith   }
754ced766e8SHong Zhang   w3 = coloring->w3;
755840b8ebdSBarry Smith 
7565904e1b1SBarry Smith   ierr = MatSetUnfactored(J);CHKERRQ(ierr);
75790d69ab7SBarry Smith   flg  = PETSC_FALSE;
75890d69ab7SBarry Smith   ierr = PetscOptionsGetTruth(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg,PETSC_NULL);CHKERRQ(ierr);
759f1af5d2fSBarry Smith   if (flg) {
760ae15b995SBarry Smith     ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr);
761840b8ebdSBarry Smith   } else {
7620b9b6f31SBarry Smith     PetscTruth assembled;
7630b9b6f31SBarry Smith     ierr = MatAssembled(J,&assembled);CHKERRQ(ierr);
7640b9b6f31SBarry Smith     if (assembled) {
765840b8ebdSBarry Smith       ierr = MatZeroEntries(J);CHKERRQ(ierr);
766840b8ebdSBarry Smith     }
7670b9b6f31SBarry Smith   }
768840b8ebdSBarry Smith 
769840b8ebdSBarry Smith   ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr);
770840b8ebdSBarry Smith   ierr = VecGetSize(x1,&N);CHKERRQ(ierr);
77166f9b7ceSBarry Smith   ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
772840b8ebdSBarry Smith   ierr = (*f)(sctx,t,x1,w1,fctx);CHKERRQ(ierr);
77366f9b7ceSBarry Smith   ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
774840b8ebdSBarry Smith 
775840b8ebdSBarry Smith   /*
7765904e1b1SBarry Smith       Compute all the scale factors and share with other processors
777840b8ebdSBarry Smith   */
7785904e1b1SBarry Smith   ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);xx = xx - start;
7795904e1b1SBarry Smith   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);vscale_array = vscale_array - start;
780840b8ebdSBarry Smith   for (k=0; k<coloring->ncolors; k++) {
781840b8ebdSBarry Smith     /*
782840b8ebdSBarry Smith        Loop over each column associated with color adding the
783840b8ebdSBarry Smith        perturbation to the vector w3.
784840b8ebdSBarry Smith     */
785840b8ebdSBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
786840b8ebdSBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
7875904e1b1SBarry Smith       dx  = xx[col];
7885b8514ebSBarry Smith       if (dx == 0.0) dx = 1.0;
789aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
790840b8ebdSBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
791840b8ebdSBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
792840b8ebdSBarry Smith #else
793329f5518SBarry Smith       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
794329f5518SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
795840b8ebdSBarry Smith #endif
796840b8ebdSBarry Smith       dx                *= epsilon;
7975904e1b1SBarry Smith       vscale_array[col] = 1.0/dx;
798840b8ebdSBarry Smith     }
7995904e1b1SBarry Smith   }
8005904e1b1SBarry Smith   vscale_array = vscale_array - start;ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
8015904e1b1SBarry Smith   ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8025904e1b1SBarry Smith   ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8035904e1b1SBarry Smith 
8045904e1b1SBarry Smith   if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow;
8055904e1b1SBarry Smith   else                        vscaleforrow = coloring->columnsforrow;
8065904e1b1SBarry Smith 
8075904e1b1SBarry Smith   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
8085904e1b1SBarry Smith   /*
8095904e1b1SBarry Smith       Loop over each color
8105904e1b1SBarry Smith   */
8115904e1b1SBarry Smith   for (k=0; k<coloring->ncolors; k++) {
8125904e1b1SBarry Smith     ierr = VecCopy(x1,w3);CHKERRQ(ierr);
8135904e1b1SBarry Smith     ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);w3_array = w3_array - start;
8145904e1b1SBarry Smith     /*
8155904e1b1SBarry Smith        Loop over each column associated with color adding the
8165904e1b1SBarry Smith        perturbation to the vector w3.
8175904e1b1SBarry Smith     */
8185904e1b1SBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
8195904e1b1SBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
8205904e1b1SBarry Smith       dx  = xx[col];
8215b8514ebSBarry Smith       if (dx == 0.0) dx = 1.0;
8225904e1b1SBarry Smith #if !defined(PETSC_USE_COMPLEX)
8235904e1b1SBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
8245904e1b1SBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
8255904e1b1SBarry Smith #else
8265904e1b1SBarry Smith       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
8275904e1b1SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
8285904e1b1SBarry Smith #endif
8295904e1b1SBarry Smith       dx            *= epsilon;
8305904e1b1SBarry Smith       w3_array[col] += dx;
8315904e1b1SBarry Smith     }
8325904e1b1SBarry Smith     w3_array = w3_array + start; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
8335904e1b1SBarry Smith 
834840b8ebdSBarry Smith     /*
835840b8ebdSBarry Smith        Evaluate function at x1 + dx (here dx is a vector of perturbations)
836840b8ebdSBarry Smith     */
83766f9b7ceSBarry Smith     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
838840b8ebdSBarry Smith     ierr = (*f)(sctx,t,w3,w2,fctx);CHKERRQ(ierr);
83966f9b7ceSBarry Smith     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
840efb30889SBarry Smith     ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr);
8415904e1b1SBarry Smith 
842840b8ebdSBarry Smith     /*
843840b8ebdSBarry Smith        Loop over rows of vector, putting results into Jacobian matrix
844840b8ebdSBarry Smith     */
8453b28642cSBarry Smith     ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
846840b8ebdSBarry Smith     for (l=0; l<coloring->nrows[k]; l++) {
847840b8ebdSBarry Smith       row    = coloring->rows[k][l];
848840b8ebdSBarry Smith       col    = coloring->columnsforrow[k][l];
8495904e1b1SBarry Smith       y[row] *= vscale_array[vscaleforrow[k][l]];
850840b8ebdSBarry Smith       srow   = row + start;
851840b8ebdSBarry Smith       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
852840b8ebdSBarry Smith     }
8533b28642cSBarry Smith     ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
854840b8ebdSBarry Smith   }
8555904e1b1SBarry Smith   ierr  = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
8565904e1b1SBarry Smith   xx    = xx + start; ierr  = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
857840b8ebdSBarry Smith   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
858840b8ebdSBarry Smith   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
859d5ba7fb7SMatthew Knepley   ierr  = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
8603a40ed3dSBarry Smith   PetscFunctionReturn(0);
861840b8ebdSBarry Smith }
8623b28642cSBarry Smith 
8633b28642cSBarry Smith 
8643b28642cSBarry Smith 
865