xref: /petsc/src/mat/matfd/fdmatrix.c (revision d4a378da529a33c9fd3485de6780ea57674d47bd)
1be1d678aSKris Buschelman 
2bbf0e169SBarry Smith /*
3639f9d9dSBarry Smith    This is where the abstract matrix operations are defined that are
4639f9d9dSBarry Smith   used for finite difference computations of Jacobians using coloring.
5bbf0e169SBarry Smith */
6bbf0e169SBarry Smith 
7c6db04a5SJed Brown #include <private/matimpl.h>        /*I "petscmat.h" I*/
8bbf0e169SBarry Smith 
94a2ae208SSatish Balay #undef __FUNCT__
103a7fca6bSBarry Smith #define __FUNCT__ "MatFDColoringSetF"
117087cfbeSBarry Smith PetscErrorCode  MatFDColoringSetF(MatFDColoring fd,Vec F)
123a7fca6bSBarry Smith {
133a7fca6bSBarry Smith   PetscFunctionBegin;
143a7fca6bSBarry Smith   fd->F = F;
153a7fca6bSBarry Smith   PetscFunctionReturn(0);
163a7fca6bSBarry Smith }
173a7fca6bSBarry Smith 
183a7fca6bSBarry Smith #undef __FUNCT__
194a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Draw_Zoom"
206849ba73SBarry Smith static PetscErrorCode MatFDColoringView_Draw_Zoom(PetscDraw draw,void *Aa)
219194eea9SBarry Smith {
229194eea9SBarry Smith   MatFDColoring  fd = (MatFDColoring)Aa;
23dfbe8321SBarry Smith   PetscErrorCode ierr;
2438baddfdSBarry Smith   PetscInt       i,j;
259194eea9SBarry Smith   PetscReal      x,y;
269194eea9SBarry Smith 
279194eea9SBarry Smith   PetscFunctionBegin;
289194eea9SBarry Smith 
299194eea9SBarry Smith   /* loop over colors  */
309194eea9SBarry Smith   for (i=0; i<fd->ncolors; i++) {
319194eea9SBarry Smith     for (j=0; j<fd->nrows[i]; j++) {
329194eea9SBarry Smith       y = fd->M - fd->rows[i][j] - fd->rstart;
339194eea9SBarry Smith       x = fd->columnsforrow[i][j];
34b0a32e0cSBarry Smith       ierr = PetscDrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);CHKERRQ(ierr);
359194eea9SBarry Smith     }
369194eea9SBarry Smith   }
379194eea9SBarry Smith   PetscFunctionReturn(0);
389194eea9SBarry Smith }
399194eea9SBarry Smith 
404a2ae208SSatish Balay #undef __FUNCT__
414a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Draw"
426849ba73SBarry Smith static PetscErrorCode MatFDColoringView_Draw(MatFDColoring fd,PetscViewer viewer)
43005c665bSBarry Smith {
44dfbe8321SBarry Smith   PetscErrorCode ierr;
45ace3abfcSBarry Smith   PetscBool      isnull;
46b0a32e0cSBarry Smith   PetscDraw      draw;
4754d96fa3SBarry Smith   PetscReal      xr,yr,xl,yl,h,w;
48005c665bSBarry Smith 
493a40ed3dSBarry Smith   PetscFunctionBegin;
50b0a32e0cSBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
51b0a32e0cSBarry Smith   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
529194eea9SBarry Smith 
539194eea9SBarry Smith   ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
54005c665bSBarry Smith 
55005c665bSBarry Smith   xr  = fd->N; yr = fd->M; h = yr/10.0; w = xr/10.0;
56005c665bSBarry Smith   xr += w;     yr += h;    xl = -w;     yl = -h;
57b0a32e0cSBarry Smith   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
58b0a32e0cSBarry Smith   ierr = PetscDrawZoom(draw,MatFDColoringView_Draw_Zoom,fd);CHKERRQ(ierr);
599194eea9SBarry Smith   ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
603a40ed3dSBarry Smith   PetscFunctionReturn(0);
61005c665bSBarry Smith }
62005c665bSBarry Smith 
634a2ae208SSatish Balay #undef __FUNCT__
644a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView"
65bbf0e169SBarry Smith /*@C
66639f9d9dSBarry Smith    MatFDColoringView - Views a finite difference coloring context.
67bbf0e169SBarry Smith 
684c49b128SBarry Smith    Collective on MatFDColoring
69fee21e36SBarry Smith 
70ef5ee4d1SLois Curfman McInnes    Input  Parameters:
71ef5ee4d1SLois Curfman McInnes +  c - the coloring context
72ef5ee4d1SLois Curfman McInnes -  viewer - visualization context
73ef5ee4d1SLois Curfman McInnes 
7415091d37SBarry Smith    Level: intermediate
7515091d37SBarry Smith 
76b4fc646aSLois Curfman McInnes    Notes:
77b4fc646aSLois Curfman McInnes    The available visualization contexts include
78b0a32e0cSBarry Smith +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
79b0a32e0cSBarry Smith .     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
80ef5ee4d1SLois Curfman McInnes         output where only the first processor opens
81ef5ee4d1SLois Curfman McInnes         the file.  All other processors send their
82ef5ee4d1SLois Curfman McInnes         data to the first processor to print.
83b0a32e0cSBarry Smith -     PETSC_VIEWER_DRAW_WORLD - graphical display of nonzero structure
84bbf0e169SBarry Smith 
859194eea9SBarry Smith    Notes:
869194eea9SBarry Smith      Since PETSc uses only a small number of basic colors (currently 33), if the coloring
87b0a32e0cSBarry Smith    involves more than 33 then some seemingly identical colors are displayed making it look
889194eea9SBarry Smith    like an illegal coloring. This is just a graphical artifact.
899194eea9SBarry Smith 
90639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
91005c665bSBarry Smith 
92b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, view
93bbf0e169SBarry Smith @*/
947087cfbeSBarry Smith PetscErrorCode  MatFDColoringView(MatFDColoring c,PetscViewer viewer)
95bbf0e169SBarry Smith {
966849ba73SBarry Smith   PetscErrorCode    ierr;
9738baddfdSBarry Smith   PetscInt          i,j;
98ace3abfcSBarry Smith   PetscBool         isdraw,iascii;
99f3ef73ceSBarry Smith   PetscViewerFormat format;
100bbf0e169SBarry Smith 
1013a40ed3dSBarry Smith   PetscFunctionBegin;
1020700a824SBarry Smith   PetscValidHeaderSpecific(c,MAT_FDCOLORING_CLASSID,1);
1033050cee2SBarry Smith   if (!viewer) {
1047adad957SLisandro Dalcin     ierr = PetscViewerASCIIGetStdout(((PetscObject)c)->comm,&viewer);CHKERRQ(ierr);
1053050cee2SBarry Smith   }
1060700a824SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
107c9780b6fSBarry Smith   PetscCheckSameComm(c,1,viewer,2);
108bbf0e169SBarry Smith 
1092692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
1102692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1110f5bd95cSBarry Smith   if (isdraw) {
112b4fc646aSLois Curfman McInnes     ierr = MatFDColoringView_Draw(c,viewer);CHKERRQ(ierr);
11332077d6dSBarry Smith   } else if (iascii) {
114317d6ea6SBarry Smith     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)c,viewer,"MatFDColoring Object");CHKERRQ(ierr);
115a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Error tolerance=%G\n",c->error_rel);CHKERRQ(ierr);
116a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Umin=%G\n",c->umin);CHKERRQ(ierr);
11777431f27SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Number of colors=%D\n",c->ncolors);CHKERRQ(ierr);
118ae09f205SBarry Smith 
119b0a32e0cSBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
120fb9695e5SSatish Balay     if (format != PETSC_VIEWER_ASCII_INFO) {
121b4fc646aSLois Curfman McInnes       for (i=0; i<c->ncolors; i++) {
12277431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"  Information for color %D\n",i);CHKERRQ(ierr);
12377431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    Number of columns %D\n",c->ncolumns[i]);CHKERRQ(ierr);
124b4fc646aSLois Curfman McInnes         for (j=0; j<c->ncolumns[i]; j++) {
12577431f27SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"      %D\n",c->columns[i][j]);CHKERRQ(ierr);
126639f9d9dSBarry Smith         }
12777431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    Number of rows %D\n",c->nrows[i]);CHKERRQ(ierr);
128b4fc646aSLois Curfman McInnes         for (j=0; j<c->nrows[i]; j++) {
12977431f27SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"      %D %D \n",c->rows[i][j],c->columnsforrow[i][j]);CHKERRQ(ierr);
130b4fc646aSLois Curfman McInnes         }
131bbf0e169SBarry Smith       }
132bbf0e169SBarry Smith     }
133b0a32e0cSBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1345cd90555SBarry Smith   } else {
135e32f2f54SBarry Smith     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for MatFDColoring",((PetscObject)viewer)->type_name);
136bbf0e169SBarry Smith   }
1373a40ed3dSBarry Smith   PetscFunctionReturn(0);
138639f9d9dSBarry Smith }
139639f9d9dSBarry Smith 
1404a2ae208SSatish Balay #undef __FUNCT__
1414a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetParameters"
142639f9d9dSBarry Smith /*@
143b4fc646aSLois Curfman McInnes    MatFDColoringSetParameters - Sets the parameters for the sparse approximation of
144b4fc646aSLois Curfman McInnes    a Jacobian matrix using finite differences.
145639f9d9dSBarry Smith 
1463f9fe445SBarry Smith    Logically Collective on MatFDColoring
147ef5ee4d1SLois Curfman McInnes 
148ef5ee4d1SLois Curfman McInnes    The Jacobian is estimated with the differencing approximation
149ef5ee4d1SLois Curfman McInnes .vb
15065f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
151f23b5b22SLois Curfman McInnes        h = error_rel*u[i]                 if  abs(u[i]) > umin
152f23b5b22SLois Curfman McInnes          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
153ef5ee4d1SLois Curfman McInnes        dx_{i} = (0, ... 1, .... 0)
154ef5ee4d1SLois Curfman McInnes .ve
155639f9d9dSBarry Smith 
156639f9d9dSBarry Smith    Input Parameters:
157ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
158639f9d9dSBarry Smith .  error_rel - relative error
159f23b5b22SLois Curfman McInnes -  umin - minimum allowable u-value magnitude
160fee21e36SBarry Smith 
16115091d37SBarry Smith    Level: advanced
16215091d37SBarry Smith 
163b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, parameters
164b4fc646aSLois Curfman McInnes 
16595b89fc3SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringSetFromOptions()
16695b89fc3SBarry Smith 
167639f9d9dSBarry Smith @*/
1687087cfbeSBarry Smith PetscErrorCode  MatFDColoringSetParameters(MatFDColoring matfd,PetscReal error,PetscReal umin)
169639f9d9dSBarry Smith {
1703a40ed3dSBarry Smith   PetscFunctionBegin;
1710700a824SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
172c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(matfd,error,2);
173c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(matfd,umin,3);
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 
1863f9fe445SBarry Smith    Not Collective
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 @*/
2027087cfbeSBarry Smith PetscErrorCode  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 
2163f9fe445SBarry Smith    Logically 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 @*/
2447087cfbeSBarry Smith PetscErrorCode  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 @*/
2887087cfbeSBarry Smith PetscErrorCode  MatFDColoringSetFromOptions(MatFDColoring matfd)
289639f9d9dSBarry Smith {
290dfbe8321SBarry Smith   PetscErrorCode ierr;
291ace3abfcSBarry Smith   PetscBool      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);
3105d973c19SBarry Smith 
3115d973c19SBarry Smith     /* process any options handlers added with PetscObjectAddOptionsHandler() */
3125d973c19SBarry 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;
322ace3abfcSBarry Smith   PetscBool      flg = PETSC_FALSE;
3233050cee2SBarry Smith   PetscViewer    viewer;
324005c665bSBarry Smith 
3253a40ed3dSBarry Smith   PetscFunctionBegin;
3267adad957SLisandro Dalcin   ierr = PetscViewerASCIIGetStdout(((PetscObject)fd)->comm,&viewer);CHKERRQ(ierr);
327acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(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;
332acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(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;
339acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(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
35794013140SBarry Smith -  iscoloring - the coloring of the matrix; usually obtained with MatGetColoring() or DMGetColoring()
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(),
36694013140SBarry Smith           MatFDColoringView(), MatFDColoringSetParameters(), MatGetColoring(), DMGetColoring()
367bbf0e169SBarry Smith @*/
3687087cfbeSBarry Smith PetscErrorCode  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 @*/
4196bf464f9SBarry Smith PetscErrorCode  MatFDColoringDestroy(MatFDColoring *c)
420bbf0e169SBarry Smith {
4216849ba73SBarry Smith   PetscErrorCode ierr;
42238baddfdSBarry Smith   PetscInt       i;
423bbf0e169SBarry Smith 
4243a40ed3dSBarry Smith   PetscFunctionBegin;
4256bf464f9SBarry Smith   if (!*c) PetscFunctionReturn(0);
4266bf464f9SBarry Smith   if (--((PetscObject)(*c))->refct > 0) {*c = 0; PetscFunctionReturn(0);}
427d4bb536fSBarry Smith 
4286bf464f9SBarry Smith   for (i=0; i<(*c)->ncolors; i++) {
4296bf464f9SBarry Smith     ierr = PetscFree((*c)->columns[i]);CHKERRQ(ierr);
4306bf464f9SBarry Smith     ierr = PetscFree((*c)->rows[i]);CHKERRQ(ierr);
4316bf464f9SBarry Smith     ierr = PetscFree((*c)->columnsforrow[i]);CHKERRQ(ierr);
4326bf464f9SBarry Smith     if ((*c)->vscaleforrow) {ierr = PetscFree((*c)->vscaleforrow[i]);CHKERRQ(ierr);}
433bbf0e169SBarry Smith   }
4346bf464f9SBarry Smith   ierr = PetscFree((*c)->ncolumns);CHKERRQ(ierr);
4356bf464f9SBarry Smith   ierr = PetscFree((*c)->columns);CHKERRQ(ierr);
4366bf464f9SBarry Smith   ierr = PetscFree((*c)->nrows);CHKERRQ(ierr);
4376bf464f9SBarry Smith   ierr = PetscFree((*c)->rows);CHKERRQ(ierr);
4386bf464f9SBarry Smith   ierr = PetscFree((*c)->columnsforrow);CHKERRQ(ierr);
4396bf464f9SBarry Smith   ierr = PetscFree((*c)->vscaleforrow);CHKERRQ(ierr);
4406bf464f9SBarry Smith   ierr = VecDestroy(&(*c)->vscale);CHKERRQ(ierr);
4416bf464f9SBarry Smith   ierr = VecDestroy(&(*c)->w1);CHKERRQ(ierr);
4426bf464f9SBarry Smith   ierr = VecDestroy(&(*c)->w2);CHKERRQ(ierr);
4436bf464f9SBarry Smith   ierr = VecDestroy(&(*c)->w3);CHKERRQ(ierr);
444d38fa0fbSBarry Smith   ierr = PetscHeaderDestroy(c);CHKERRQ(ierr);
4453a40ed3dSBarry Smith   PetscFunctionReturn(0);
446bbf0e169SBarry Smith }
44743a90d84SBarry Smith 
4484a2ae208SSatish Balay #undef __FUNCT__
44949b058dcSBarry Smith #define __FUNCT__ "MatFDColoringGetPerturbedColumns"
45049b058dcSBarry Smith /*@C
45149b058dcSBarry Smith     MatFDColoringGetPerturbedColumns - Returns the indices of the columns that
45249b058dcSBarry Smith       that are currently being perturbed.
45349b058dcSBarry Smith 
45449b058dcSBarry Smith     Not Collective
45549b058dcSBarry Smith 
45649b058dcSBarry Smith     Input Parameters:
45749b058dcSBarry Smith .   coloring - coloring context created with MatFDColoringCreate()
45849b058dcSBarry Smith 
45949b058dcSBarry Smith     Output Parameters:
46049b058dcSBarry Smith +   n - the number of local columns being perturbed
46149b058dcSBarry Smith -   cols - the column indices, in global numbering
46249b058dcSBarry Smith 
46349b058dcSBarry Smith    Level: intermediate
46449b058dcSBarry Smith 
46549b058dcSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringApply()
46649b058dcSBarry Smith 
46749b058dcSBarry Smith .keywords: coloring, Jacobian, finite differences
46849b058dcSBarry Smith @*/
4697087cfbeSBarry Smith PetscErrorCode  MatFDColoringGetPerturbedColumns(MatFDColoring coloring,PetscInt *n,PetscInt *cols[])
47049b058dcSBarry Smith {
47149b058dcSBarry Smith   PetscFunctionBegin;
47249b058dcSBarry Smith   if (coloring->currentcolor >= 0) {
47349b058dcSBarry Smith     *n    = coloring->ncolumns[coloring->currentcolor];
47449b058dcSBarry Smith     *cols = coloring->columns[coloring->currentcolor];
47549b058dcSBarry Smith   } else {
47649b058dcSBarry Smith     *n = 0;
47749b058dcSBarry Smith   }
47849b058dcSBarry Smith   PetscFunctionReturn(0);
47949b058dcSBarry Smith }
48049b058dcSBarry Smith 
48149b058dcSBarry Smith #undef __FUNCT__
4824a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringApply"
48343a90d84SBarry Smith /*@
484e0907662SLois Curfman McInnes     MatFDColoringApply - Given a matrix for which a MatFDColoring context
485e0907662SLois Curfman McInnes     has been created, computes the Jacobian for a function via finite differences.
48643a90d84SBarry Smith 
487fee21e36SBarry Smith     Collective on MatFDColoring
488fee21e36SBarry Smith 
489ef5ee4d1SLois Curfman McInnes     Input Parameters:
490ef5ee4d1SLois Curfman McInnes +   mat - location to store Jacobian
491ef5ee4d1SLois Curfman McInnes .   coloring - coloring context created with MatFDColoringCreate()
492ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
4937850c7c0SBarry Smith -   sctx - context required by function, if this is being used with the SNES solver then it is SNES object, otherwise it is null
494ef5ee4d1SLois Curfman McInnes 
4958bba8e72SBarry Smith     Options Database Keys:
496ebd3b9afSBarry Smith +    -mat_fd_type - "wp" or "ds"  (see MATMFFD_WP or MATMFFD_DS)
497b9722508SLois Curfman McInnes .    -mat_fd_coloring_view - Activates basic viewing or coloring
498b9722508SLois Curfman McInnes .    -mat_fd_coloring_view_draw - Activates drawing of coloring
499b9722508SLois Curfman McInnes -    -mat_fd_coloring_view_info - Activates viewing of coloring info
5008bba8e72SBarry Smith 
50115091d37SBarry Smith     Level: intermediate
50215091d37SBarry Smith 
5037850c7c0SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringSetFunction()
50443a90d84SBarry Smith 
50543a90d84SBarry Smith .keywords: coloring, Jacobian, finite differences
50643a90d84SBarry Smith @*/
5077087cfbeSBarry Smith PetscErrorCode  MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
50843a90d84SBarry Smith {
5093acb8795SBarry Smith   PetscErrorCode ierr;
5103acb8795SBarry Smith 
5113acb8795SBarry Smith   PetscFunctionBegin;
5120700a824SBarry Smith   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
5130700a824SBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_CLASSID,2);
5140700a824SBarry Smith   PetscValidHeaderSpecific(x1,VEC_CLASSID,3);
51517186662SBarry Smith   if (!coloring->f) SETERRQ(((PetscObject)J)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetFunction()");
516e32f2f54SBarry Smith   if (!J->ops->fdcoloringapply) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not supported for this matrix type %s",((PetscObject)J)->type_name);
5173acb8795SBarry Smith   ierr = (*J->ops->fdcoloringapply)(J,coloring,x1,flag,sctx);CHKERRQ(ierr);
5183acb8795SBarry Smith   PetscFunctionReturn(0);
5193acb8795SBarry Smith }
5203acb8795SBarry Smith 
5213acb8795SBarry Smith #undef __FUNCT__
5223acb8795SBarry Smith #define __FUNCT__ "MatFDColoringApply_AIJ"
5237087cfbeSBarry Smith PetscErrorCode  MatFDColoringApply_AIJ(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
5243acb8795SBarry Smith {
5256849ba73SBarry Smith   PetscErrorCode (*f)(void*,Vec,Vec,void*) = (PetscErrorCode (*)(void*,Vec,Vec,void *))coloring->f;
5266849ba73SBarry Smith   PetscErrorCode ierr;
527b8f8c88eSHong Zhang   PetscInt       k,start,end,l,row,col,srow,**vscaleforrow,m1,m2;
528efb30889SBarry Smith   PetscScalar    dx,*y,*xx,*w3_array;
529ea709b57SSatish Balay   PetscScalar    *vscale_array;
530efb30889SBarry Smith   PetscReal      epsilon = coloring->error_rel,umin = coloring->umin,unorm;
531b8f8c88eSHong Zhang   Vec            w1=coloring->w1,w2=coloring->w2,w3;
532005c665bSBarry Smith   void           *fctx = coloring->fctx;
533ace3abfcSBarry Smith   PetscBool      flg = PETSC_FALSE;
534705d48f7SSatish Balay   PetscInt       ctype=coloring->ctype,N,col_start=0,col_end=0;
535b8f8c88eSHong Zhang   Vec            x1_tmp;
536a2e34c3dSBarry Smith 
5373a40ed3dSBarry Smith   PetscFunctionBegin;
5380700a824SBarry Smith   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
5390700a824SBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_CLASSID,2);
5400700a824SBarry Smith   PetscValidHeaderSpecific(x1,VEC_CLASSID,3);
54117186662SBarry Smith   if (!f) SETERRQ(((PetscObject)J)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetFunction()");
542e0907662SLois Curfman McInnes 
543d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
5444c49b128SBarry Smith   ierr = MatSetUnfactored(J);CHKERRQ(ierr);
545acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg,PETSC_NULL);CHKERRQ(ierr);
546f1af5d2fSBarry Smith   if (flg) {
547ae15b995SBarry Smith     ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr);
548e0907662SLois Curfman McInnes   } else {
549ace3abfcSBarry Smith     PetscBool  assembled;
5500b9b6f31SBarry Smith     ierr = MatAssembled(J,&assembled);CHKERRQ(ierr);
5510b9b6f31SBarry Smith     if (assembled) {
55243a90d84SBarry Smith       ierr = MatZeroEntries(J);CHKERRQ(ierr);
553e0907662SLois Curfman McInnes     }
5540b9b6f31SBarry Smith   }
55543a90d84SBarry Smith 
556b8f8c88eSHong Zhang   x1_tmp = x1;
557dac7f5fdSBarry Smith   if (!coloring->vscale){
558b8f8c88eSHong Zhang     ierr = VecDuplicate(x1_tmp,&coloring->vscale);CHKERRQ(ierr);
559b8f8c88eSHong Zhang   }
560be2fbe1fSBarry Smith 
5613a7fca6bSBarry Smith   /*
5623a7fca6bSBarry Smith     This is a horrible, horrible, hack. See DMMGComputeJacobian_Multigrid() it inproperly sets
5633a7fca6bSBarry Smith     coloring->F for the coarser grids from the finest
5643a7fca6bSBarry Smith   */
5653a7fca6bSBarry Smith   if (coloring->F) {
5663a7fca6bSBarry Smith     ierr = VecGetLocalSize(coloring->F,&m1);CHKERRQ(ierr);
5673a7fca6bSBarry Smith     ierr = VecGetLocalSize(w1,&m2);CHKERRQ(ierr);
5683a7fca6bSBarry Smith     if (m1 != m2) {
5693a7fca6bSBarry Smith       coloring->F = 0;
5703a7fca6bSBarry Smith       }
5713a7fca6bSBarry Smith     }
5723a7fca6bSBarry Smith 
573b8f8c88eSHong Zhang   if (coloring->htype[0] == 'w') { /* tacky test; need to make systematic if we add other approaches to computing h*/
574b8f8c88eSHong Zhang     ierr = VecNorm(x1_tmp,NORM_2,&unorm);CHKERRQ(ierr);
575b8f8c88eSHong Zhang   }
576b8f8c88eSHong Zhang   ierr = VecGetOwnershipRange(w1,&start,&end);CHKERRQ(ierr); /* OwnershipRange is used by ghosted x! */
577b8f8c88eSHong Zhang 
578b8f8c88eSHong Zhang   /* Set w1 = F(x1) */
579be2fbe1fSBarry Smith   if (coloring->F) {
580be2fbe1fSBarry Smith     w1          = coloring->F; /* use already computed value of function */
581be2fbe1fSBarry Smith     coloring->F = 0;
582be2fbe1fSBarry Smith   } else {
58366f9b7ceSBarry Smith     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
584b8f8c88eSHong Zhang     ierr = (*f)(sctx,x1_tmp,w1,fctx);CHKERRQ(ierr);
58566f9b7ceSBarry Smith     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
586be2fbe1fSBarry Smith   }
58743a90d84SBarry Smith 
588b8f8c88eSHong Zhang   if (!coloring->w3) {
589b8f8c88eSHong Zhang     ierr = VecDuplicate(x1_tmp,&coloring->w3);CHKERRQ(ierr);
590b8f8c88eSHong Zhang     ierr = PetscLogObjectParent(coloring,coloring->w3);CHKERRQ(ierr);
591efb30889SBarry Smith   }
592b8f8c88eSHong Zhang   w3 = coloring->w3;
593efb30889SBarry Smith 
594b8f8c88eSHong Zhang     /* Compute all the local scale factors, including ghost points */
595b8f8c88eSHong Zhang   ierr = VecGetLocalSize(x1_tmp,&N);CHKERRQ(ierr);
596b8f8c88eSHong Zhang   ierr = VecGetArray(x1_tmp,&xx);CHKERRQ(ierr);
597b8f8c88eSHong Zhang   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
598b8f8c88eSHong Zhang   if (ctype == IS_COLORING_GHOSTED){
599b8f8c88eSHong Zhang     col_start = 0; col_end = N;
6008ee2e534SBarry Smith   } else if (ctype == IS_COLORING_GLOBAL){
601b8f8c88eSHong Zhang     xx = xx - start;
602b8f8c88eSHong Zhang     vscale_array = vscale_array - start;
603b8f8c88eSHong Zhang     col_start = start; col_end = N + start;
604b8f8c88eSHong Zhang   }
605b8f8c88eSHong Zhang   for (col=col_start; col<col_end; col++){
606b8f8c88eSHong Zhang     /* Loop over each local column, vscale[col] = 1./(epsilon*dx[col]) */
607efb30889SBarry Smith     if (coloring->htype[0] == 'w') {
608efb30889SBarry Smith       dx = 1.0 + unorm;
609efb30889SBarry Smith     } else {
6109782cf98SBarry Smith       dx  = xx[col];
611efb30889SBarry Smith     }
612*d4a378daSJed Brown     if (dx == (PetscScalar)0.0) dx = 1.0;
6139782cf98SBarry Smith #if !defined(PETSC_USE_COMPLEX)
6149782cf98SBarry Smith     if (dx < umin && dx >= 0.0)      dx = umin;
6159782cf98SBarry Smith     else if (dx < 0.0 && dx > -umin) dx = -umin;
6169782cf98SBarry Smith #else
6179782cf98SBarry Smith     if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
6189782cf98SBarry Smith     else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
6199782cf98SBarry Smith #endif
6209782cf98SBarry Smith     dx               *= epsilon;
621*d4a378daSJed Brown     vscale_array[col] = (PetscScalar)1.0/dx;
6229782cf98SBarry Smith   }
6238ee2e534SBarry Smith   if (ctype == IS_COLORING_GLOBAL)  vscale_array = vscale_array + start;
624efb30889SBarry Smith   ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
6258ee2e534SBarry Smith   if (ctype == IS_COLORING_GLOBAL){
62630b34957SBarry Smith     ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
62730b34957SBarry Smith     ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
628b8f8c88eSHong Zhang   }
6299782cf98SBarry Smith 
630b8f8c88eSHong Zhang   if (coloring->vscaleforrow) {
631b8f8c88eSHong Zhang     vscaleforrow = coloring->vscaleforrow;
63217186662SBarry Smith   } else SETERRQ(((PetscObject)J)->comm,PETSC_ERR_ARG_NULL,"Null Object: coloring->vscaleforrow");
633b0a32e0cSBarry Smith 
6349782cf98SBarry Smith   /*
63543a90d84SBarry Smith     Loop over each color
63643a90d84SBarry Smith   */
637b8f8c88eSHong Zhang   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
63843a90d84SBarry Smith   for (k=0; k<coloring->ncolors; k++) {
63949b058dcSBarry Smith     coloring->currentcolor = k;
640b8f8c88eSHong Zhang     ierr = VecCopy(x1_tmp,w3);CHKERRQ(ierr);
641b8f8c88eSHong Zhang     ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);
6428ee2e534SBarry Smith     if (ctype == IS_COLORING_GLOBAL) w3_array = w3_array - start;
64343a90d84SBarry Smith     /*
644b8f8c88eSHong Zhang       Loop over each column associated with color
645b8f8c88eSHong Zhang       adding the perturbation to the vector w3.
64643a90d84SBarry Smith     */
64743a90d84SBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
648b8f8c88eSHong Zhang       col = coloring->columns[k][l];    /* local column of the matrix we are probing for */
649efb30889SBarry Smith       if (coloring->htype[0] == 'w') {
650efb30889SBarry Smith         dx = 1.0 + unorm;
651efb30889SBarry Smith       } else {
65242460c72SBarry Smith         dx  = xx[col];
653efb30889SBarry Smith       }
654*d4a378daSJed Brown       if (dx == (PetscScalar)0.0) dx = 1.0;
655aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
656ae09f205SBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
657ae09f205SBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
65843a90d84SBarry Smith #else
659329f5518SBarry Smith       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
660329f5518SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
66143a90d84SBarry Smith #endif
66243a90d84SBarry Smith       dx            *= epsilon;
663e32f2f54SBarry Smith       if (!PetscAbsScalar(dx)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Computed 0 differencing parameter");
66442460c72SBarry Smith       w3_array[col] += dx;
66543a90d84SBarry Smith     }
6668ee2e534SBarry Smith     if (ctype == IS_COLORING_GLOBAL) w3_array = w3_array + start;
667b8f8c88eSHong Zhang     ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
6683b28642cSBarry Smith 
66943a90d84SBarry Smith     /*
670b8f8c88eSHong Zhang       Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations)
671b8f8c88eSHong Zhang                            w2 = F(x1 + dx) - F(x1)
67243a90d84SBarry Smith     */
67366f9b7ceSBarry Smith     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
67443a90d84SBarry Smith     ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr);
67566f9b7ceSBarry Smith     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
676efb30889SBarry Smith     ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr);
6779782cf98SBarry Smith 
67843a90d84SBarry Smith     /*
679e0907662SLois Curfman McInnes       Loop over rows of vector, putting results into Jacobian matrix
68043a90d84SBarry Smith     */
6813b28642cSBarry Smith     ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
68243a90d84SBarry Smith     for (l=0; l<coloring->nrows[k]; l++) {
683b8f8c88eSHong Zhang       row    = coloring->rows[k][l];             /* local row index */
684b8f8c88eSHong Zhang       col    = coloring->columnsforrow[k][l];    /* global column index */
6855904e1b1SBarry Smith       y[row] *= vscale_array[vscaleforrow[k][l]];
68643a90d84SBarry Smith       srow   = row + start;
68743a90d84SBarry Smith       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
68843a90d84SBarry Smith     }
6893b28642cSBarry Smith     ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
690aeb7e98aSMatthew Knepley   } /* endof for each color */
6918ee2e534SBarry Smith   if (ctype == IS_COLORING_GLOBAL) xx = xx + start;
69230b34957SBarry Smith   ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
693b8f8c88eSHong Zhang   ierr = VecRestoreArray(x1_tmp,&xx);CHKERRQ(ierr);
694b8f8c88eSHong Zhang 
695b8f8c88eSHong Zhang   coloring->currentcolor = -1;
69643a90d84SBarry Smith   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
69743a90d84SBarry Smith   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
698d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
699a2e34c3dSBarry Smith 
70090d69ab7SBarry Smith   flg  = PETSC_FALSE;
701acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-mat_null_space_test",&flg,PETSC_NULL);CHKERRQ(ierr);
702a2e34c3dSBarry Smith   if (flg) {
70395902228SMatthew Knepley     ierr = MatNullSpaceTest(J->nullsp,J,PETSC_NULL);CHKERRQ(ierr);
704a2e34c3dSBarry Smith   }
705b9722508SLois Curfman McInnes   ierr = MatFDColoringView_Private(coloring);CHKERRQ(ierr);
7063a40ed3dSBarry Smith   PetscFunctionReturn(0);
70743a90d84SBarry Smith }
708840b8ebdSBarry Smith 
7094a2ae208SSatish Balay #undef __FUNCT__
7104a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringApplyTS"
711840b8ebdSBarry Smith /*@
712840b8ebdSBarry Smith     MatFDColoringApplyTS - Given a matrix for which a MatFDColoring context
713840b8ebdSBarry Smith     has been created, computes the Jacobian for a function via finite differences.
714840b8ebdSBarry Smith 
715fee21e36SBarry Smith    Collective on Mat, MatFDColoring, and Vec
716fee21e36SBarry Smith 
717ef5ee4d1SLois Curfman McInnes     Input Parameters:
7183b28642cSBarry Smith +   mat - location to store Jacobian
719ef5ee4d1SLois Curfman McInnes .   coloring - coloring context created with MatFDColoringCreate()
720ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
7217850c7c0SBarry Smith -   sctx - context required by function, if this is being used with the TS solver then it is TS object, otherwise it is null
722ef5ee4d1SLois Curfman McInnes 
72315091d37SBarry Smith    Level: intermediate
72415091d37SBarry Smith 
7257850c7c0SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringSetFunction()
726840b8ebdSBarry Smith 
727840b8ebdSBarry Smith .keywords: coloring, Jacobian, finite differences
728840b8ebdSBarry Smith @*/
7297087cfbeSBarry Smith PetscErrorCode  MatFDColoringApplyTS(Mat J,MatFDColoring coloring,PetscReal t,Vec x1,MatStructure *flag,void *sctx)
730840b8ebdSBarry Smith {
7316849ba73SBarry Smith   PetscErrorCode (*f)(void*,PetscReal,Vec,Vec,void*)=(PetscErrorCode (*)(void*,PetscReal,Vec,Vec,void *))coloring->f;
7326849ba73SBarry Smith   PetscErrorCode ierr;
73338baddfdSBarry Smith   PetscInt       k,N,start,end,l,row,col,srow,**vscaleforrow;
734efb30889SBarry Smith   PetscScalar    dx,*y,*xx,*w3_array;
735ea709b57SSatish Balay   PetscScalar    *vscale_array;
736329f5518SBarry Smith   PetscReal      epsilon = coloring->error_rel,umin = coloring->umin;
737ced766e8SHong Zhang   Vec            w1=coloring->w1,w2=coloring->w2,w3;
738840b8ebdSBarry Smith   void           *fctx = coloring->fctx;
739ace3abfcSBarry Smith   PetscBool      flg;
740840b8ebdSBarry Smith 
7413a40ed3dSBarry Smith   PetscFunctionBegin;
7420700a824SBarry Smith   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
7430700a824SBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_CLASSID,2);
7440700a824SBarry Smith   PetscValidHeaderSpecific(x1,VEC_CLASSID,4);
745840b8ebdSBarry Smith 
746d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
747ced766e8SHong Zhang   if (!coloring->w3) {
748840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
74952e6d16bSBarry Smith     ierr = PetscLogObjectParent(coloring,coloring->w3);CHKERRQ(ierr);
750840b8ebdSBarry Smith   }
751ced766e8SHong Zhang   w3 = coloring->w3;
752840b8ebdSBarry Smith 
7535904e1b1SBarry Smith   ierr = MatSetUnfactored(J);CHKERRQ(ierr);
75490d69ab7SBarry Smith   flg  = PETSC_FALSE;
755acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg,PETSC_NULL);CHKERRQ(ierr);
756f1af5d2fSBarry Smith   if (flg) {
757ae15b995SBarry Smith     ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr);
758840b8ebdSBarry Smith   } else {
759ace3abfcSBarry Smith     PetscBool  assembled;
7600b9b6f31SBarry Smith     ierr = MatAssembled(J,&assembled);CHKERRQ(ierr);
7610b9b6f31SBarry Smith     if (assembled) {
762840b8ebdSBarry Smith       ierr = MatZeroEntries(J);CHKERRQ(ierr);
763840b8ebdSBarry Smith     }
7640b9b6f31SBarry Smith   }
765840b8ebdSBarry Smith 
766840b8ebdSBarry Smith   ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr);
767840b8ebdSBarry Smith   ierr = VecGetSize(x1,&N);CHKERRQ(ierr);
76866f9b7ceSBarry Smith   ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
769840b8ebdSBarry Smith   ierr = (*f)(sctx,t,x1,w1,fctx);CHKERRQ(ierr);
77066f9b7ceSBarry Smith   ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
771840b8ebdSBarry Smith 
772840b8ebdSBarry Smith   /*
7735904e1b1SBarry Smith       Compute all the scale factors and share with other processors
774840b8ebdSBarry Smith   */
7755904e1b1SBarry Smith   ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);xx = xx - start;
7765904e1b1SBarry Smith   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);vscale_array = vscale_array - start;
777840b8ebdSBarry Smith   for (k=0; k<coloring->ncolors; k++) {
778840b8ebdSBarry Smith     /*
779840b8ebdSBarry Smith        Loop over each column associated with color adding the
780840b8ebdSBarry Smith        perturbation to the vector w3.
781840b8ebdSBarry Smith     */
782840b8ebdSBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
783840b8ebdSBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
7845904e1b1SBarry Smith       dx  = xx[col];
785*d4a378daSJed Brown       if (dx == (PetscScalar)0.0) dx = 1.0;
786aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
787840b8ebdSBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
788840b8ebdSBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
789840b8ebdSBarry Smith #else
790329f5518SBarry Smith       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
791329f5518SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
792840b8ebdSBarry Smith #endif
793840b8ebdSBarry Smith       dx                *= epsilon;
794*d4a378daSJed Brown       vscale_array[col] = (PetscScalar)1.0/dx;
795840b8ebdSBarry Smith     }
7965904e1b1SBarry Smith   }
7975904e1b1SBarry Smith   vscale_array = vscale_array - start;ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
7985904e1b1SBarry Smith   ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7995904e1b1SBarry Smith   ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8005904e1b1SBarry Smith 
8015904e1b1SBarry Smith   if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow;
8025904e1b1SBarry Smith   else                        vscaleforrow = coloring->columnsforrow;
8035904e1b1SBarry Smith 
8045904e1b1SBarry Smith   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
8055904e1b1SBarry Smith   /*
8065904e1b1SBarry Smith       Loop over each color
8075904e1b1SBarry Smith   */
8085904e1b1SBarry Smith   for (k=0; k<coloring->ncolors; k++) {
8095904e1b1SBarry Smith     ierr = VecCopy(x1,w3);CHKERRQ(ierr);
8105904e1b1SBarry Smith     ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);w3_array = w3_array - start;
8115904e1b1SBarry Smith     /*
8125904e1b1SBarry Smith        Loop over each column associated with color adding the
8135904e1b1SBarry Smith        perturbation to the vector w3.
8145904e1b1SBarry Smith     */
8155904e1b1SBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
8165904e1b1SBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
8175904e1b1SBarry Smith       dx  = xx[col];
818*d4a378daSJed Brown       if (dx == (PetscScalar)0.0) dx = 1.0;
8195904e1b1SBarry Smith #if !defined(PETSC_USE_COMPLEX)
8205904e1b1SBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
8215904e1b1SBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
8225904e1b1SBarry Smith #else
8235904e1b1SBarry Smith       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
8245904e1b1SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
8255904e1b1SBarry Smith #endif
8265904e1b1SBarry Smith       dx            *= epsilon;
8275904e1b1SBarry Smith       w3_array[col] += dx;
8285904e1b1SBarry Smith     }
8295904e1b1SBarry Smith     w3_array = w3_array + start; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
8305904e1b1SBarry Smith 
831840b8ebdSBarry Smith     /*
832840b8ebdSBarry Smith        Evaluate function at x1 + dx (here dx is a vector of perturbations)
833840b8ebdSBarry Smith     */
83466f9b7ceSBarry Smith     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
835840b8ebdSBarry Smith     ierr = (*f)(sctx,t,w3,w2,fctx);CHKERRQ(ierr);
83666f9b7ceSBarry Smith     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
837efb30889SBarry Smith     ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr);
8385904e1b1SBarry Smith 
839840b8ebdSBarry Smith     /*
840840b8ebdSBarry Smith        Loop over rows of vector, putting results into Jacobian matrix
841840b8ebdSBarry Smith     */
8423b28642cSBarry Smith     ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
843840b8ebdSBarry Smith     for (l=0; l<coloring->nrows[k]; l++) {
844840b8ebdSBarry Smith       row    = coloring->rows[k][l];
845840b8ebdSBarry Smith       col    = coloring->columnsforrow[k][l];
8465904e1b1SBarry Smith       y[row] *= vscale_array[vscaleforrow[k][l]];
847840b8ebdSBarry Smith       srow   = row + start;
848840b8ebdSBarry Smith       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
849840b8ebdSBarry Smith     }
8503b28642cSBarry Smith     ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
851840b8ebdSBarry Smith   }
8525904e1b1SBarry Smith   ierr  = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
8535904e1b1SBarry Smith   xx    = xx + start; ierr  = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
854840b8ebdSBarry Smith   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
855840b8ebdSBarry Smith   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
856d5ba7fb7SMatthew Knepley   ierr  = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
8573a40ed3dSBarry Smith   PetscFunctionReturn(0);
858840b8ebdSBarry Smith }
8593b28642cSBarry Smith 
8603b28642cSBarry Smith 
8613b28642cSBarry Smith 
862