xref: /petsc/src/mat/matfd/fdmatrix.c (revision 3acb87959449f26430bfb44c593b059e1dfb39fd)
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;
1034482741eSBarry Smith   PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE,1);
1043050cee2SBarry Smith   if (!viewer) {
1057adad957SLisandro Dalcin     ierr = PetscViewerASCIIGetStdout(((PetscObject)c)->comm,&viewer);CHKERRQ(ierr);
1063050cee2SBarry Smith   }
1074482741eSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE,2);
108c9780b6fSBarry Smith   PetscCheckSameComm(c,1,viewer,2);
109bbf0e169SBarry Smith 
110fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
11132077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
1120f5bd95cSBarry Smith   if (isdraw) {
113b4fc646aSLois Curfman McInnes     ierr = MatFDColoringView_Draw(c,viewer);CHKERRQ(ierr);
11432077d6dSBarry Smith   } else if (iascii) {
115b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"MatFDColoring Object:\n");CHKERRQ(ierr);
116a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Error tolerance=%G\n",c->error_rel);CHKERRQ(ierr);
117a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Umin=%G\n",c->umin);CHKERRQ(ierr);
11877431f27SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Number of colors=%D\n",c->ncolors);CHKERRQ(ierr);
119ae09f205SBarry Smith 
120b0a32e0cSBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
121fb9695e5SSatish Balay     if (format != PETSC_VIEWER_ASCII_INFO) {
122b4fc646aSLois Curfman McInnes       for (i=0; i<c->ncolors; i++) {
12377431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"  Information for color %D\n",i);CHKERRQ(ierr);
12477431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    Number of columns %D\n",c->ncolumns[i]);CHKERRQ(ierr);
125b4fc646aSLois Curfman McInnes         for (j=0; j<c->ncolumns[i]; j++) {
12677431f27SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"      %D\n",c->columns[i][j]);CHKERRQ(ierr);
127639f9d9dSBarry Smith         }
12877431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    Number of rows %D\n",c->nrows[i]);CHKERRQ(ierr);
129b4fc646aSLois Curfman McInnes         for (j=0; j<c->nrows[i]; j++) {
13077431f27SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"      %D %D \n",c->rows[i][j],c->columnsforrow[i][j]);CHKERRQ(ierr);
131b4fc646aSLois Curfman McInnes         }
132bbf0e169SBarry Smith       }
133bbf0e169SBarry Smith     }
134b0a32e0cSBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1355cd90555SBarry Smith   } else {
13679a5c55eSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for MatFDColoring",((PetscObject)viewer)->type_name);
137bbf0e169SBarry Smith   }
1383a40ed3dSBarry Smith   PetscFunctionReturn(0);
139639f9d9dSBarry Smith }
140639f9d9dSBarry Smith 
1414a2ae208SSatish Balay #undef __FUNCT__
1424a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetParameters"
143639f9d9dSBarry Smith /*@
144b4fc646aSLois Curfman McInnes    MatFDColoringSetParameters - Sets the parameters for the sparse approximation of
145b4fc646aSLois Curfman McInnes    a Jacobian matrix using finite differences.
146639f9d9dSBarry Smith 
147ef5ee4d1SLois Curfman McInnes    Collective on MatFDColoring
148ef5ee4d1SLois Curfman McInnes 
149ef5ee4d1SLois Curfman McInnes    The Jacobian is estimated with the differencing approximation
150ef5ee4d1SLois Curfman McInnes .vb
15165f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
152f23b5b22SLois Curfman McInnes        h = error_rel*u[i]                 if  abs(u[i]) > umin
153f23b5b22SLois Curfman McInnes          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
154ef5ee4d1SLois Curfman McInnes        dx_{i} = (0, ... 1, .... 0)
155ef5ee4d1SLois Curfman McInnes .ve
156639f9d9dSBarry Smith 
157639f9d9dSBarry Smith    Input Parameters:
158ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
159639f9d9dSBarry Smith .  error_rel - relative error
160f23b5b22SLois Curfman McInnes -  umin - minimum allowable u-value magnitude
161fee21e36SBarry Smith 
16215091d37SBarry Smith    Level: advanced
16315091d37SBarry Smith 
164b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, parameters
165b4fc646aSLois Curfman McInnes 
166b4fc646aSLois Curfman McInnes .seealso: MatFDColoringCreate()
167639f9d9dSBarry Smith @*/
168be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetParameters(MatFDColoring matfd,PetscReal error,PetscReal umin)
169639f9d9dSBarry Smith {
1703a40ed3dSBarry Smith   PetscFunctionBegin;
1714482741eSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1);
172639f9d9dSBarry Smith 
173639f9d9dSBarry Smith   if (error != PETSC_DEFAULT) matfd->error_rel = error;
174639f9d9dSBarry Smith   if (umin != PETSC_DEFAULT)  matfd->umin      = umin;
1753a40ed3dSBarry Smith   PetscFunctionReturn(0);
176639f9d9dSBarry Smith }
177639f9d9dSBarry Smith 
178005c665bSBarry Smith 
179ff0cfa39SBarry Smith 
1804a2ae208SSatish Balay #undef __FUNCT__
1814a9d489dSBarry Smith #define __FUNCT__ "MatFDColoringGetFunction"
1824a9d489dSBarry Smith /*@C
1834a9d489dSBarry Smith    MatFDColoringGetFunction - Gets the function to use for computing the Jacobian.
1844a9d489dSBarry Smith 
1854a9d489dSBarry Smith    Collective on MatFDColoring
1864a9d489dSBarry Smith 
1874a9d489dSBarry Smith    Input Parameters:
1884a9d489dSBarry Smith .  coloring - the coloring context
1894a9d489dSBarry Smith 
1904a9d489dSBarry Smith    Output Parameters:
1914a9d489dSBarry Smith +  f - the function
1924a9d489dSBarry Smith -  fctx - the optional user-defined function context
1934a9d489dSBarry Smith 
1944a9d489dSBarry Smith    Level: intermediate
1954a9d489dSBarry Smith 
1964a9d489dSBarry Smith .keywords: Mat, Jacobian, finite differences, set, function
1974a9d489dSBarry Smith @*/
1984a9d489dSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFunction(MatFDColoring matfd,PetscErrorCode (**f)(void),void **fctx)
1994a9d489dSBarry Smith {
2004a9d489dSBarry Smith   PetscFunctionBegin;
2014a9d489dSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1);
2024a9d489dSBarry Smith   if (f) *f = matfd->f;
2034a9d489dSBarry Smith   if (fctx) *fctx = matfd->fctx;
2044a9d489dSBarry Smith   PetscFunctionReturn(0);
2054a9d489dSBarry Smith }
2064a9d489dSBarry Smith 
2074a9d489dSBarry Smith #undef __FUNCT__
2084a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFunction"
209d64ed03dSBarry Smith /*@C
210005c665bSBarry Smith    MatFDColoringSetFunction - Sets the function to use for computing the Jacobian.
211005c665bSBarry Smith 
212fee21e36SBarry Smith    Collective on MatFDColoring
213fee21e36SBarry Smith 
214ef5ee4d1SLois Curfman McInnes    Input Parameters:
215ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
216ef5ee4d1SLois Curfman McInnes .  f - the function
217ef5ee4d1SLois Curfman McInnes -  fctx - the optional user-defined function context
218ef5ee4d1SLois Curfman McInnes 
2197850c7c0SBarry Smith    Calling sequence of (*f) function:
2207850c7c0SBarry Smith     For SNES:    PetscErrorCode (*f)(SNES,Vec,Vec,void*)
2217850c7c0SBarry Smith     For TS:      PetscErrorCode (*f)(TS,PetscReal,Vec,Vec,void*)
2227850c7c0SBarry Smith     If not using SNES or TS: PetscErrorCode (*f)(void *dummy,Vec,Vec,void*) and dummy is ignored
22315091d37SBarry Smith 
2247850c7c0SBarry Smith    Level: advanced
2257850c7c0SBarry Smith 
2267850c7c0SBarry Smith    Notes: This function is usually used automatically by SNES or TS (when one uses SNESSetJacobian() with the argument
2277850c7c0SBarry Smith      SNESDefaultComputeJacobianColor() or TSSetRHSJacobian() with the argument TSDefaultComputeJacobianColor()) and only needs to be used
2287850c7c0SBarry Smith      by someone computing a matrix via coloring directly by calling MatFDColoringApply()
2297850c7c0SBarry Smith 
2307850c7c0SBarry Smith    Fortran Notes:
2317850c7c0SBarry Smith     In Fortran you must call MatFDColoringSetFunction() for a coloring object to
2327850c7c0SBarry Smith   be used without SNES or TS or within the SNES solvers and MatFDColoringSetFunctionTS() if it is to be used
2337850c7c0SBarry Smith   within the TS solvers.
234f881d145SBarry Smith 
235b4fc646aSLois Curfman McInnes .keywords: Mat, Jacobian, finite differences, set, function
236005c665bSBarry Smith @*/
237be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFunction(MatFDColoring matfd,PetscErrorCode (*f)(void),void *fctx)
238005c665bSBarry Smith {
2393a40ed3dSBarry Smith   PetscFunctionBegin;
2404482741eSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1);
241005c665bSBarry Smith   matfd->f    = f;
242005c665bSBarry Smith   matfd->fctx = fctx;
2433a40ed3dSBarry Smith   PetscFunctionReturn(0);
244005c665bSBarry Smith }
245005c665bSBarry Smith 
2464a2ae208SSatish Balay #undef __FUNCT__
2474a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFromOptions"
248639f9d9dSBarry Smith /*@
249b4fc646aSLois Curfman McInnes    MatFDColoringSetFromOptions - Sets coloring finite difference parameters from
250639f9d9dSBarry Smith    the options database.
251639f9d9dSBarry Smith 
252fee21e36SBarry Smith    Collective on MatFDColoring
253fee21e36SBarry Smith 
25465f2ba5bSLois Curfman McInnes    The Jacobian, F'(u), is estimated with the differencing approximation
255ef5ee4d1SLois Curfman McInnes .vb
25665f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
257f23b5b22SLois Curfman McInnes        h = error_rel*u[i]                 if  abs(u[i]) > umin
258f23b5b22SLois Curfman McInnes          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
259ef5ee4d1SLois Curfman McInnes        dx_{i} = (0, ... 1, .... 0)
260ef5ee4d1SLois Curfman McInnes .ve
261ef5ee4d1SLois Curfman McInnes 
262ef5ee4d1SLois Curfman McInnes    Input Parameter:
263ef5ee4d1SLois Curfman McInnes .  coloring - the coloring context
264ef5ee4d1SLois Curfman McInnes 
265b4fc646aSLois Curfman McInnes    Options Database Keys:
266d15ffeeaSSatish Balay +  -mat_fd_coloring_err <err> - Sets <err> (square root
267ef5ee4d1SLois Curfman McInnes            of relative error in the function)
268f23b5b22SLois Curfman McInnes .  -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude
2693ec795f1SBarry Smith .  -mat_fd_type - "wp" or "ds" (see MATMFFD_WP or MATMFFD_DS)
270ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_view - Activates basic viewing
271ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_view_info - Activates viewing info
272ef5ee4d1SLois Curfman McInnes -  -mat_fd_coloring_view_draw - Activates drawing
273639f9d9dSBarry Smith 
27415091d37SBarry Smith     Level: intermediate
27515091d37SBarry Smith 
276b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, parameters
277d1c39f55SBarry Smith 
278d1c39f55SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringView(), MatFDColoringSetParameters()
279d1c39f55SBarry Smith 
280639f9d9dSBarry Smith @*/
281be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFromOptions(MatFDColoring matfd)
282639f9d9dSBarry Smith {
283dfbe8321SBarry Smith   PetscErrorCode ierr;
284efb30889SBarry Smith   PetscTruth     flg;
285efb30889SBarry Smith   char           value[3];
2863a40ed3dSBarry Smith 
2873a40ed3dSBarry Smith   PetscFunctionBegin;
2884482741eSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1);
289639f9d9dSBarry Smith 
2907adad957SLisandro Dalcin   ierr = PetscOptionsBegin(((PetscObject)matfd)->comm,((PetscObject)matfd)->prefix,"Jacobian computation via finite differences option","MatFD");CHKERRQ(ierr);
29187828ca2SBarry Smith     ierr = PetscOptionsReal("-mat_fd_coloring_err","Square root of relative error in function","MatFDColoringSetParameters",matfd->error_rel,&matfd->error_rel,0);CHKERRQ(ierr);
29287828ca2SBarry Smith     ierr = PetscOptionsReal("-mat_fd_coloring_umin","Minimum allowable u magnitude","MatFDColoringSetParameters",matfd->umin,&matfd->umin,0);CHKERRQ(ierr);
293efb30889SBarry Smith     ierr = PetscOptionsString("-mat_fd_type","Algorithm to compute h, wp or ds","MatFDColoringCreate",matfd->htype,value,2,&flg);CHKERRQ(ierr);
294efb30889SBarry Smith     if (flg) {
295efb30889SBarry Smith       if (value[0] == 'w' && value[1] == 'p') matfd->htype = "wp";
296efb30889SBarry Smith       else if (value[0] == 'd' && value[1] == 's') matfd->htype = "ds";
297efb30889SBarry Smith       else SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Unknown finite differencing type %s",value);
298efb30889SBarry Smith     }
299186905e3SBarry Smith     /* not used here; but so they are presented in the GUI */
300b0a32e0cSBarry Smith     ierr = PetscOptionsName("-mat_fd_coloring_view","Print entire datastructure used for Jacobian","None",0);CHKERRQ(ierr);
301b0a32e0cSBarry Smith     ierr = PetscOptionsName("-mat_fd_coloring_view_info","Print number of colors etc for Jacobian","None",0);CHKERRQ(ierr);
302b0a32e0cSBarry Smith     ierr = PetscOptionsName("-mat_fd_coloring_view_draw","Plot nonzero structure ofJacobian","None",0);CHKERRQ(ierr);
303b0a32e0cSBarry Smith   PetscOptionsEnd();CHKERRQ(ierr);
3043a40ed3dSBarry Smith   PetscFunctionReturn(0);
305005c665bSBarry Smith }
306005c665bSBarry Smith 
3074a2ae208SSatish Balay #undef __FUNCT__
3084a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Private"
309dfbe8321SBarry Smith PetscErrorCode MatFDColoringView_Private(MatFDColoring fd)
310005c665bSBarry Smith {
311dfbe8321SBarry Smith   PetscErrorCode ierr;
31290d69ab7SBarry Smith   PetscTruth     flg = PETSC_FALSE;
3133050cee2SBarry Smith   PetscViewer    viewer;
314005c665bSBarry Smith 
3153a40ed3dSBarry Smith   PetscFunctionBegin;
3167adad957SLisandro Dalcin   ierr = PetscViewerASCIIGetStdout(((PetscObject)fd)->comm,&viewer);CHKERRQ(ierr);
31790d69ab7SBarry Smith   ierr = PetscOptionsGetTruth(PETSC_NULL,"-mat_fd_coloring_view",&flg,PETSC_NULL);CHKERRQ(ierr);
318005c665bSBarry Smith   if (flg) {
3193050cee2SBarry Smith     ierr = MatFDColoringView(fd,viewer);CHKERRQ(ierr);
320005c665bSBarry Smith   }
32190d69ab7SBarry Smith   flg  = PETSC_FALSE;
32290d69ab7SBarry Smith   ierr = PetscOptionsGetTruth(PETSC_NULL,"-mat_fd_coloring_view_info",&flg,PETSC_NULL);CHKERRQ(ierr);
323ae09f205SBarry Smith   if (flg) {
3243050cee2SBarry Smith     ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
3253050cee2SBarry Smith     ierr = MatFDColoringView(fd,viewer);CHKERRQ(ierr);
3263050cee2SBarry Smith     ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
327ae09f205SBarry Smith   }
32890d69ab7SBarry Smith   flg  = PETSC_FALSE;
32990d69ab7SBarry Smith   ierr = PetscOptionsGetTruth(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg,PETSC_NULL);CHKERRQ(ierr);
330005c665bSBarry Smith   if (flg) {
3317adad957SLisandro Dalcin     ierr = MatFDColoringView(fd,PETSC_VIEWER_DRAW_(((PetscObject)fd)->comm));CHKERRQ(ierr);
3327adad957SLisandro Dalcin     ierr = PetscViewerFlush(PETSC_VIEWER_DRAW_(((PetscObject)fd)->comm));CHKERRQ(ierr);
333005c665bSBarry Smith   }
3343a40ed3dSBarry Smith   PetscFunctionReturn(0);
335bbf0e169SBarry Smith }
336bbf0e169SBarry Smith 
3374a2ae208SSatish Balay #undef __FUNCT__
3384a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringCreate"
33905869f15SSatish Balay /*@
340639f9d9dSBarry Smith    MatFDColoringCreate - Creates a matrix coloring context for finite difference
341639f9d9dSBarry Smith    computation of Jacobians.
342bbf0e169SBarry Smith 
343ef5ee4d1SLois Curfman McInnes    Collective on Mat
344ef5ee4d1SLois Curfman McInnes 
345639f9d9dSBarry Smith    Input Parameters:
346ef5ee4d1SLois Curfman McInnes +  mat - the matrix containing the nonzero structure of the Jacobian
347ef5ee4d1SLois Curfman McInnes -  iscoloring - the coloring of the matrix
348bbf0e169SBarry Smith 
349bbf0e169SBarry Smith     Output Parameter:
350639f9d9dSBarry Smith .   color - the new coloring context
351bbf0e169SBarry Smith 
35215091d37SBarry Smith     Level: intermediate
35315091d37SBarry Smith 
3546831982aSBarry Smith .seealso: MatFDColoringDestroy(),SNESDefaultComputeJacobianColor(), ISColoringCreate(),
355d1c39f55SBarry Smith           MatFDColoringSetFunction(), MatFDColoringSetFromOptions(), MatFDColoringApply(),
356ebd3b9afSBarry Smith           MatFDColoringView(), MatFDColoringSetParameters()
357bbf0e169SBarry Smith @*/
358be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
359bbf0e169SBarry Smith {
360639f9d9dSBarry Smith   MatFDColoring  c;
361639f9d9dSBarry Smith   MPI_Comm       comm;
362dfbe8321SBarry Smith   PetscErrorCode ierr;
36338baddfdSBarry Smith   PetscInt       M,N;
364b8f8c88eSHong Zhang   PetscMPIInt    size;
365639f9d9dSBarry Smith 
3663a40ed3dSBarry Smith   PetscFunctionBegin;
367d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
368639f9d9dSBarry Smith   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
36929bbc08cSBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,"Only for square matrices");
370639f9d9dSBarry Smith 
371f881d145SBarry Smith   ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
37252e6d16bSBarry Smith   ierr = PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_COOKIE,0,"MatFDColoring",comm,MatFDColoringDestroy,MatFDColoringView);CHKERRQ(ierr);
373639f9d9dSBarry Smith 
374b8f8c88eSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
375b8f8c88eSHong Zhang   c->ctype = iscoloring->ctype;
376b8f8c88eSHong Zhang 
377f830108cSBarry Smith   if (mat->ops->fdcoloringcreate) {
378f830108cSBarry Smith     ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);CHKERRQ(ierr);
379639f9d9dSBarry Smith   } else {
38029bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Code not yet written for this matrix type");
381639f9d9dSBarry Smith   }
382639f9d9dSBarry Smith 
383b8f8c88eSHong Zhang   ierr = MatGetVecs(mat,PETSC_NULL,&c->w1);CHKERRQ(ierr);
384b8f8c88eSHong Zhang   ierr = PetscLogObjectParent(c,c->w1);CHKERRQ(ierr);
385b8f8c88eSHong Zhang   ierr = VecDuplicate(c->w1,&c->w2);CHKERRQ(ierr);
386b8f8c88eSHong Zhang   ierr = PetscLogObjectParent(c,c->w2);CHKERRQ(ierr);
387b8f8c88eSHong Zhang 
38877d8c4bbSBarry Smith   c->error_rel         = PETSC_SQRT_MACHINE_EPSILON;
38977d8c4bbSBarry Smith   c->umin              = 100.0*PETSC_SQRT_MACHINE_EPSILON;
39049b058dcSBarry Smith   c->currentcolor      = -1;
391efb30889SBarry Smith   c->htype             = "wp";
392639f9d9dSBarry Smith 
393639f9d9dSBarry Smith   *color = c;
394d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
3953a40ed3dSBarry Smith   PetscFunctionReturn(0);
396bbf0e169SBarry Smith }
397bbf0e169SBarry Smith 
3984a2ae208SSatish Balay #undef __FUNCT__
3994a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringDestroy"
40005869f15SSatish Balay /*@
401639f9d9dSBarry Smith     MatFDColoringDestroy - Destroys a matrix coloring context that was created
402639f9d9dSBarry Smith     via MatFDColoringCreate().
403bbf0e169SBarry Smith 
404ef5ee4d1SLois Curfman McInnes     Collective on MatFDColoring
405ef5ee4d1SLois Curfman McInnes 
406b4fc646aSLois Curfman McInnes     Input Parameter:
407639f9d9dSBarry Smith .   c - coloring context
408bbf0e169SBarry Smith 
40915091d37SBarry Smith     Level: intermediate
41015091d37SBarry Smith 
411639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
412bbf0e169SBarry Smith @*/
413be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringDestroy(MatFDColoring c)
414bbf0e169SBarry Smith {
4156849ba73SBarry Smith   PetscErrorCode ierr;
41638baddfdSBarry Smith   PetscInt       i;
417bbf0e169SBarry Smith 
4183a40ed3dSBarry Smith   PetscFunctionBegin;
4197adad957SLisandro Dalcin   if (--((PetscObject)c)->refct > 0) PetscFunctionReturn(0);
420d4bb536fSBarry Smith 
421639f9d9dSBarry Smith   for (i=0; i<c->ncolors; i++) {
42205b42c5fSBarry Smith     ierr = PetscFree(c->columns[i]);CHKERRQ(ierr);
42305b42c5fSBarry Smith     ierr = PetscFree(c->rows[i]);CHKERRQ(ierr);
42405b42c5fSBarry Smith     ierr = PetscFree(c->columnsforrow[i]);CHKERRQ(ierr);
42505b42c5fSBarry Smith     if (c->vscaleforrow) {ierr = PetscFree(c->vscaleforrow[i]);CHKERRQ(ierr);}
426bbf0e169SBarry Smith   }
427606d414cSSatish Balay   ierr = PetscFree(c->ncolumns);CHKERRQ(ierr);
428606d414cSSatish Balay   ierr = PetscFree(c->columns);CHKERRQ(ierr);
429606d414cSSatish Balay   ierr = PetscFree(c->nrows);CHKERRQ(ierr);
430606d414cSSatish Balay   ierr = PetscFree(c->rows);CHKERRQ(ierr);
431606d414cSSatish Balay   ierr = PetscFree(c->columnsforrow);CHKERRQ(ierr);
43205b42c5fSBarry Smith   ierr = PetscFree(c->vscaleforrow);CHKERRQ(ierr);
4334f113abfSBarry Smith   if (c->vscale)       {ierr = VecDestroy(c->vscale);CHKERRQ(ierr);}
434005c665bSBarry Smith   if (c->w1) {
435005c665bSBarry Smith     ierr = VecDestroy(c->w1);CHKERRQ(ierr);
436005c665bSBarry Smith     ierr = VecDestroy(c->w2);CHKERRQ(ierr);
437b8f8c88eSHong Zhang   }
438b8f8c88eSHong Zhang   if (c->w3){
439005c665bSBarry Smith     ierr = VecDestroy(c->w3);CHKERRQ(ierr);
440005c665bSBarry Smith   }
441d38fa0fbSBarry Smith   ierr = PetscHeaderDestroy(c);CHKERRQ(ierr);
4423a40ed3dSBarry Smith   PetscFunctionReturn(0);
443bbf0e169SBarry Smith }
44443a90d84SBarry Smith 
4454a2ae208SSatish Balay #undef __FUNCT__
44649b058dcSBarry Smith #define __FUNCT__ "MatFDColoringGetPerturbedColumns"
44749b058dcSBarry Smith /*@C
44849b058dcSBarry Smith     MatFDColoringGetPerturbedColumns - Returns the indices of the columns that
44949b058dcSBarry Smith       that are currently being perturbed.
45049b058dcSBarry Smith 
45149b058dcSBarry Smith     Not Collective
45249b058dcSBarry Smith 
45349b058dcSBarry Smith     Input Parameters:
45449b058dcSBarry Smith .   coloring - coloring context created with MatFDColoringCreate()
45549b058dcSBarry Smith 
45649b058dcSBarry Smith     Output Parameters:
45749b058dcSBarry Smith +   n - the number of local columns being perturbed
45849b058dcSBarry Smith -   cols - the column indices, in global numbering
45949b058dcSBarry Smith 
46049b058dcSBarry Smith    Level: intermediate
46149b058dcSBarry Smith 
46249b058dcSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringApply()
46349b058dcSBarry Smith 
46449b058dcSBarry Smith .keywords: coloring, Jacobian, finite differences
46549b058dcSBarry Smith @*/
466be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetPerturbedColumns(MatFDColoring coloring,PetscInt *n,PetscInt *cols[])
46749b058dcSBarry Smith {
46849b058dcSBarry Smith   PetscFunctionBegin;
46949b058dcSBarry Smith   if (coloring->currentcolor >= 0) {
47049b058dcSBarry Smith     *n    = coloring->ncolumns[coloring->currentcolor];
47149b058dcSBarry Smith     *cols = coloring->columns[coloring->currentcolor];
47249b058dcSBarry Smith   } else {
47349b058dcSBarry Smith     *n = 0;
47449b058dcSBarry Smith   }
47549b058dcSBarry Smith   PetscFunctionReturn(0);
47649b058dcSBarry Smith }
47749b058dcSBarry Smith 
478b8f8c88eSHong Zhang #include "petscda.h"      /*I      "petscda.h"    I*/
479b8f8c88eSHong Zhang 
48049b058dcSBarry Smith #undef __FUNCT__
4814a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringApply"
48243a90d84SBarry Smith /*@
483e0907662SLois Curfman McInnes     MatFDColoringApply - Given a matrix for which a MatFDColoring context
484e0907662SLois Curfman McInnes     has been created, computes the Jacobian for a function via finite differences.
48543a90d84SBarry Smith 
486fee21e36SBarry Smith     Collective on MatFDColoring
487fee21e36SBarry Smith 
488ef5ee4d1SLois Curfman McInnes     Input Parameters:
489ef5ee4d1SLois Curfman McInnes +   mat - location to store Jacobian
490ef5ee4d1SLois Curfman McInnes .   coloring - coloring context created with MatFDColoringCreate()
491ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
4927850c7c0SBarry Smith -   sctx - context required by function, if this is being used with the SNES solver then it is SNES object, otherwise it is null
493ef5ee4d1SLois Curfman McInnes 
4948bba8e72SBarry Smith     Options Database Keys:
495ebd3b9afSBarry Smith +    -mat_fd_type - "wp" or "ds"  (see MATMFFD_WP or MATMFFD_DS)
496b9722508SLois Curfman McInnes .    -mat_fd_coloring_view - Activates basic viewing or coloring
497b9722508SLois Curfman McInnes .    -mat_fd_coloring_view_draw - Activates drawing of coloring
498b9722508SLois Curfman McInnes -    -mat_fd_coloring_view_info - Activates viewing of coloring info
4998bba8e72SBarry Smith 
50015091d37SBarry Smith     Level: intermediate
50115091d37SBarry Smith 
5027850c7c0SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringSetFunction()
50343a90d84SBarry Smith 
50443a90d84SBarry Smith .keywords: coloring, Jacobian, finite differences
50543a90d84SBarry Smith @*/
506be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
50743a90d84SBarry Smith {
508*3acb8795SBarry Smith   PetscErrorCode ierr;
509*3acb8795SBarry Smith 
510*3acb8795SBarry Smith   PetscFunctionBegin;
511*3acb8795SBarry Smith   PetscValidHeaderSpecific(J,MAT_COOKIE,1);
512*3acb8795SBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE,2);
513*3acb8795SBarry Smith   PetscValidHeaderSpecific(x1,VEC_COOKIE,3);
514*3acb8795SBarry Smith   if (!coloring->f) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetFunction()");
515*3acb8795SBarry Smith   if (!J->ops->fdcoloringapply) SETERRQ1(PETSC_ERR_SUP,"Not supported for this matrix type %s",((PetscObject)J)->type_name);
516*3acb8795SBarry Smith   ierr = (*J->ops->fdcoloringapply)(J,coloring,x1,flag,sctx);CHKERRQ(ierr);
517*3acb8795SBarry Smith   PetscFunctionReturn(0);
518*3acb8795SBarry Smith }
519*3acb8795SBarry Smith 
520*3acb8795SBarry Smith #undef __FUNCT__
521*3acb8795SBarry Smith #define __FUNCT__ "MatFDColoringApply_AIJ"
522*3acb8795SBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApply_AIJ(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
523*3acb8795SBarry Smith {
5246849ba73SBarry Smith   PetscErrorCode (*f)(void*,Vec,Vec,void*) = (PetscErrorCode (*)(void*,Vec,Vec,void *))coloring->f;
5256849ba73SBarry Smith   PetscErrorCode ierr;
526b8f8c88eSHong Zhang   PetscInt       k,start,end,l,row,col,srow,**vscaleforrow,m1,m2;
527efb30889SBarry Smith   PetscScalar    dx,*y,*xx,*w3_array;
528ea709b57SSatish Balay   PetscScalar    *vscale_array;
529efb30889SBarry Smith   PetscReal      epsilon = coloring->error_rel,umin = coloring->umin,unorm;
530b8f8c88eSHong Zhang   Vec            w1=coloring->w1,w2=coloring->w2,w3;
531005c665bSBarry Smith   void           *fctx = coloring->fctx;
53290d69ab7SBarry Smith   PetscTruth     flg = PETSC_FALSE;
533705d48f7SSatish Balay   PetscInt       ctype=coloring->ctype,N,col_start=0,col_end=0;
534b8f8c88eSHong Zhang   Vec            x1_tmp;
535a2e34c3dSBarry Smith 
5363a40ed3dSBarry Smith   PetscFunctionBegin;
5374482741eSBarry Smith   PetscValidHeaderSpecific(J,MAT_COOKIE,1);
5384482741eSBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE,2);
5394482741eSBarry Smith   PetscValidHeaderSpecific(x1,VEC_COOKIE,3);
540feba9168SBarry Smith   if (!f) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetFunction()");
541e0907662SLois Curfman McInnes 
542d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
5434c49b128SBarry Smith   ierr = MatSetUnfactored(J);CHKERRQ(ierr);
54490d69ab7SBarry Smith   ierr = PetscOptionsGetTruth(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg,PETSC_NULL);CHKERRQ(ierr);
545f1af5d2fSBarry Smith   if (flg) {
546ae15b995SBarry Smith     ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr);
547e0907662SLois Curfman McInnes   } else {
5480b9b6f31SBarry Smith     PetscTruth assembled;
5490b9b6f31SBarry Smith     ierr = MatAssembled(J,&assembled);CHKERRQ(ierr);
5500b9b6f31SBarry Smith     if (assembled) {
55143a90d84SBarry Smith       ierr = MatZeroEntries(J);CHKERRQ(ierr);
552e0907662SLois Curfman McInnes     }
5530b9b6f31SBarry Smith   }
55443a90d84SBarry Smith 
555b8f8c88eSHong Zhang   x1_tmp = x1;
556dac7f5fdSBarry Smith   if (!coloring->vscale){
557b8f8c88eSHong Zhang     ierr = VecDuplicate(x1_tmp,&coloring->vscale);CHKERRQ(ierr);
558b8f8c88eSHong Zhang   }
559be2fbe1fSBarry Smith 
5603a7fca6bSBarry Smith   /*
5613a7fca6bSBarry Smith     This is a horrible, horrible, hack. See DMMGComputeJacobian_Multigrid() it inproperly sets
5623a7fca6bSBarry Smith     coloring->F for the coarser grids from the finest
5633a7fca6bSBarry Smith   */
5643a7fca6bSBarry Smith   if (coloring->F) {
5653a7fca6bSBarry Smith     ierr = VecGetLocalSize(coloring->F,&m1);CHKERRQ(ierr);
5663a7fca6bSBarry Smith     ierr = VecGetLocalSize(w1,&m2);CHKERRQ(ierr);
5673a7fca6bSBarry Smith     if (m1 != m2) {
5683a7fca6bSBarry Smith       coloring->F = 0;
5693a7fca6bSBarry Smith       }
5703a7fca6bSBarry Smith     }
5713a7fca6bSBarry Smith 
572b8f8c88eSHong Zhang   if (coloring->htype[0] == 'w') { /* tacky test; need to make systematic if we add other approaches to computing h*/
573b8f8c88eSHong Zhang     ierr = VecNorm(x1_tmp,NORM_2,&unorm);CHKERRQ(ierr);
574b8f8c88eSHong Zhang   }
575b8f8c88eSHong Zhang   ierr = VecGetOwnershipRange(w1,&start,&end);CHKERRQ(ierr); /* OwnershipRange is used by ghosted x! */
576b8f8c88eSHong Zhang 
577b8f8c88eSHong Zhang   /* Set w1 = F(x1) */
578be2fbe1fSBarry Smith   if (coloring->F) {
579be2fbe1fSBarry Smith     w1          = coloring->F; /* use already computed value of function */
580be2fbe1fSBarry Smith     coloring->F = 0;
581be2fbe1fSBarry Smith   } else {
58266f9b7ceSBarry Smith     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
583b8f8c88eSHong Zhang     ierr = (*f)(sctx,x1_tmp,w1,fctx);CHKERRQ(ierr);
58466f9b7ceSBarry Smith     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
585be2fbe1fSBarry Smith   }
58643a90d84SBarry Smith 
587b8f8c88eSHong Zhang   if (!coloring->w3) {
588b8f8c88eSHong Zhang     ierr = VecDuplicate(x1_tmp,&coloring->w3);CHKERRQ(ierr);
589b8f8c88eSHong Zhang     ierr = PetscLogObjectParent(coloring,coloring->w3);CHKERRQ(ierr);
590efb30889SBarry Smith   }
591b8f8c88eSHong Zhang   w3 = coloring->w3;
592efb30889SBarry Smith 
593b8f8c88eSHong Zhang     /* Compute all the local scale factors, including ghost points */
594b8f8c88eSHong Zhang   ierr = VecGetLocalSize(x1_tmp,&N);CHKERRQ(ierr);
595b8f8c88eSHong Zhang   ierr = VecGetArray(x1_tmp,&xx);CHKERRQ(ierr);
596b8f8c88eSHong Zhang   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
597b8f8c88eSHong Zhang   if (ctype == IS_COLORING_GHOSTED){
598b8f8c88eSHong Zhang     col_start = 0; col_end = N;
5998ee2e534SBarry Smith   } else if (ctype == IS_COLORING_GLOBAL){
600b8f8c88eSHong Zhang     xx = xx - start;
601b8f8c88eSHong Zhang     vscale_array = vscale_array - start;
602b8f8c88eSHong Zhang     col_start = start; col_end = N + start;
603b8f8c88eSHong Zhang   }
604b8f8c88eSHong Zhang   for (col=col_start; col<col_end; col++){
605b8f8c88eSHong Zhang     /* Loop over each local column, vscale[col] = 1./(epsilon*dx[col]) */
606efb30889SBarry Smith     if (coloring->htype[0] == 'w') {
607efb30889SBarry Smith       dx = 1.0 + unorm;
608efb30889SBarry Smith     } else {
6099782cf98SBarry Smith       dx  = xx[col];
610efb30889SBarry Smith     }
6115b8514ebSBarry Smith     if (dx == 0.0) dx = 1.0;
6129782cf98SBarry Smith #if !defined(PETSC_USE_COMPLEX)
6139782cf98SBarry Smith     if (dx < umin && dx >= 0.0)      dx = umin;
6149782cf98SBarry Smith     else if (dx < 0.0 && dx > -umin) dx = -umin;
6159782cf98SBarry Smith #else
6169782cf98SBarry Smith     if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
6179782cf98SBarry Smith     else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
6189782cf98SBarry Smith #endif
6199782cf98SBarry Smith     dx               *= epsilon;
62030b34957SBarry Smith     vscale_array[col] = 1.0/dx;
6219782cf98SBarry Smith   }
6228ee2e534SBarry Smith   if (ctype == IS_COLORING_GLOBAL)  vscale_array = vscale_array + start;
623efb30889SBarry Smith   ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
6248ee2e534SBarry Smith   if (ctype == IS_COLORING_GLOBAL){
62530b34957SBarry Smith     ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
62630b34957SBarry Smith     ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
627b8f8c88eSHong Zhang   }
6289782cf98SBarry Smith 
629b8f8c88eSHong Zhang   if (coloring->vscaleforrow) {
630b8f8c88eSHong Zhang     vscaleforrow = coloring->vscaleforrow;
631b8f8c88eSHong Zhang   } else {
632b8f8c88eSHong Zhang     SETERRQ(PETSC_ERR_ARG_NULL,"Null Object: coloring->vscaleforrow");
633b8f8c88eSHong Zhang   }
634b0a32e0cSBarry Smith 
6359782cf98SBarry Smith   /*
63643a90d84SBarry Smith     Loop over each color
63743a90d84SBarry Smith   */
638b8f8c88eSHong Zhang   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
63943a90d84SBarry Smith   for (k=0; k<coloring->ncolors; k++) {
64049b058dcSBarry Smith     coloring->currentcolor = k;
641b8f8c88eSHong Zhang     ierr = VecCopy(x1_tmp,w3);CHKERRQ(ierr);
642b8f8c88eSHong Zhang     ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);
6438ee2e534SBarry Smith     if (ctype == IS_COLORING_GLOBAL) w3_array = w3_array - start;
64443a90d84SBarry Smith     /*
645b8f8c88eSHong Zhang       Loop over each column associated with color
646b8f8c88eSHong Zhang       adding the perturbation to the vector w3.
64743a90d84SBarry Smith     */
64843a90d84SBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
649b8f8c88eSHong Zhang       col = coloring->columns[k][l];    /* local column of the matrix we are probing for */
650efb30889SBarry Smith       if (coloring->htype[0] == 'w') {
651efb30889SBarry Smith         dx = 1.0 + unorm;
652efb30889SBarry Smith       } else {
65342460c72SBarry Smith         dx  = xx[col];
654efb30889SBarry Smith       }
6555b8514ebSBarry Smith       if (dx == 0.0) dx = 1.0;
656aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
657ae09f205SBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
658ae09f205SBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
65943a90d84SBarry Smith #else
660329f5518SBarry Smith       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
661329f5518SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
66243a90d84SBarry Smith #endif
66343a90d84SBarry Smith       dx            *= epsilon;
6641302d50aSBarry Smith       if (!PetscAbsScalar(dx)) SETERRQ(PETSC_ERR_PLIB,"Computed 0 differencing parameter");
66542460c72SBarry Smith       w3_array[col] += dx;
66643a90d84SBarry Smith     }
6678ee2e534SBarry Smith     if (ctype == IS_COLORING_GLOBAL) w3_array = w3_array + start;
668b8f8c88eSHong Zhang     ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
6693b28642cSBarry Smith 
67043a90d84SBarry Smith     /*
671b8f8c88eSHong Zhang       Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations)
672b8f8c88eSHong Zhang                            w2 = F(x1 + dx) - F(x1)
67343a90d84SBarry Smith     */
67466f9b7ceSBarry Smith     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
67543a90d84SBarry Smith     ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr);
67666f9b7ceSBarry Smith     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
677efb30889SBarry Smith     ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr);
6789782cf98SBarry Smith 
67943a90d84SBarry Smith     /*
680e0907662SLois Curfman McInnes       Loop over rows of vector, putting results into Jacobian matrix
68143a90d84SBarry Smith     */
6823b28642cSBarry Smith     ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
68343a90d84SBarry Smith     for (l=0; l<coloring->nrows[k]; l++) {
684b8f8c88eSHong Zhang       row    = coloring->rows[k][l];             /* local row index */
685b8f8c88eSHong Zhang       col    = coloring->columnsforrow[k][l];    /* global column index */
6865904e1b1SBarry Smith       y[row] *= vscale_array[vscaleforrow[k][l]];
68743a90d84SBarry Smith       srow   = row + start;
68843a90d84SBarry Smith       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
68943a90d84SBarry Smith     }
6903b28642cSBarry Smith     ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
691aeb7e98aSMatthew Knepley   } /* endof for each color */
6928ee2e534SBarry Smith   if (ctype == IS_COLORING_GLOBAL) xx = xx + start;
69330b34957SBarry Smith   ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
694b8f8c88eSHong Zhang   ierr = VecRestoreArray(x1_tmp,&xx);CHKERRQ(ierr);
695b8f8c88eSHong Zhang 
696b8f8c88eSHong Zhang   coloring->currentcolor = -1;
69743a90d84SBarry Smith   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
69843a90d84SBarry Smith   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
699d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
700a2e34c3dSBarry Smith 
70190d69ab7SBarry Smith   flg  = PETSC_FALSE;
70290d69ab7SBarry Smith   ierr = PetscOptionsGetTruth(PETSC_NULL,"-mat_null_space_test",&flg,PETSC_NULL);CHKERRQ(ierr);
703a2e34c3dSBarry Smith   if (flg) {
70495902228SMatthew Knepley     ierr = MatNullSpaceTest(J->nullsp,J,PETSC_NULL);CHKERRQ(ierr);
705a2e34c3dSBarry Smith   }
706b9722508SLois Curfman McInnes   ierr = MatFDColoringView_Private(coloring);CHKERRQ(ierr);
7073a40ed3dSBarry Smith   PetscFunctionReturn(0);
70843a90d84SBarry Smith }
709840b8ebdSBarry Smith 
7104a2ae208SSatish Balay #undef __FUNCT__
7114a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringApplyTS"
712840b8ebdSBarry Smith /*@
713840b8ebdSBarry Smith     MatFDColoringApplyTS - Given a matrix for which a MatFDColoring context
714840b8ebdSBarry Smith     has been created, computes the Jacobian for a function via finite differences.
715840b8ebdSBarry Smith 
716fee21e36SBarry Smith    Collective on Mat, MatFDColoring, and Vec
717fee21e36SBarry Smith 
718ef5ee4d1SLois Curfman McInnes     Input Parameters:
7193b28642cSBarry Smith +   mat - location to store Jacobian
720ef5ee4d1SLois Curfman McInnes .   coloring - coloring context created with MatFDColoringCreate()
721ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
7227850c7c0SBarry Smith -   sctx - context required by function, if this is being used with the TS solver then it is TS object, otherwise it is null
723ef5ee4d1SLois Curfman McInnes 
72415091d37SBarry Smith    Level: intermediate
72515091d37SBarry Smith 
7267850c7c0SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringSetFunction()
727840b8ebdSBarry Smith 
728840b8ebdSBarry Smith .keywords: coloring, Jacobian, finite differences
729840b8ebdSBarry Smith @*/
730be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApplyTS(Mat J,MatFDColoring coloring,PetscReal t,Vec x1,MatStructure *flag,void *sctx)
731840b8ebdSBarry Smith {
7326849ba73SBarry Smith   PetscErrorCode (*f)(void*,PetscReal,Vec,Vec,void*)=(PetscErrorCode (*)(void*,PetscReal,Vec,Vec,void *))coloring->f;
7336849ba73SBarry Smith   PetscErrorCode ierr;
73438baddfdSBarry Smith   PetscInt       k,N,start,end,l,row,col,srow,**vscaleforrow;
735efb30889SBarry Smith   PetscScalar    dx,*y,*xx,*w3_array;
736ea709b57SSatish Balay   PetscScalar    *vscale_array;
737329f5518SBarry Smith   PetscReal      epsilon = coloring->error_rel,umin = coloring->umin;
738ced766e8SHong Zhang   Vec            w1=coloring->w1,w2=coloring->w2,w3;
739840b8ebdSBarry Smith   void           *fctx = coloring->fctx;
740f1af5d2fSBarry Smith   PetscTruth     flg;
741840b8ebdSBarry Smith 
7423a40ed3dSBarry Smith   PetscFunctionBegin;
7434482741eSBarry Smith   PetscValidHeaderSpecific(J,MAT_COOKIE,1);
7444482741eSBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE,2);
7454482741eSBarry Smith   PetscValidHeaderSpecific(x1,VEC_COOKIE,4);
746840b8ebdSBarry Smith 
747d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
748ced766e8SHong Zhang   if (!coloring->w3) {
749840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
75052e6d16bSBarry Smith     ierr = PetscLogObjectParent(coloring,coloring->w3);CHKERRQ(ierr);
751840b8ebdSBarry Smith   }
752ced766e8SHong Zhang   w3 = coloring->w3;
753840b8ebdSBarry Smith 
7545904e1b1SBarry Smith   ierr = MatSetUnfactored(J);CHKERRQ(ierr);
75590d69ab7SBarry Smith   flg  = PETSC_FALSE;
75690d69ab7SBarry Smith   ierr = PetscOptionsGetTruth(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg,PETSC_NULL);CHKERRQ(ierr);
757f1af5d2fSBarry Smith   if (flg) {
758ae15b995SBarry Smith     ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr);
759840b8ebdSBarry Smith   } else {
7600b9b6f31SBarry Smith     PetscTruth assembled;
7610b9b6f31SBarry Smith     ierr = MatAssembled(J,&assembled);CHKERRQ(ierr);
7620b9b6f31SBarry Smith     if (assembled) {
763840b8ebdSBarry Smith       ierr = MatZeroEntries(J);CHKERRQ(ierr);
764840b8ebdSBarry Smith     }
7650b9b6f31SBarry Smith   }
766840b8ebdSBarry Smith 
767840b8ebdSBarry Smith   ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr);
768840b8ebdSBarry Smith   ierr = VecGetSize(x1,&N);CHKERRQ(ierr);
76966f9b7ceSBarry Smith   ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
770840b8ebdSBarry Smith   ierr = (*f)(sctx,t,x1,w1,fctx);CHKERRQ(ierr);
77166f9b7ceSBarry Smith   ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
772840b8ebdSBarry Smith 
773840b8ebdSBarry Smith   /*
7745904e1b1SBarry Smith       Compute all the scale factors and share with other processors
775840b8ebdSBarry Smith   */
7765904e1b1SBarry Smith   ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);xx = xx - start;
7775904e1b1SBarry Smith   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);vscale_array = vscale_array - start;
778840b8ebdSBarry Smith   for (k=0; k<coloring->ncolors; k++) {
779840b8ebdSBarry Smith     /*
780840b8ebdSBarry Smith        Loop over each column associated with color adding the
781840b8ebdSBarry Smith        perturbation to the vector w3.
782840b8ebdSBarry Smith     */
783840b8ebdSBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
784840b8ebdSBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
7855904e1b1SBarry Smith       dx  = xx[col];
7865b8514ebSBarry Smith       if (dx == 0.0) dx = 1.0;
787aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
788840b8ebdSBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
789840b8ebdSBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
790840b8ebdSBarry Smith #else
791329f5518SBarry Smith       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
792329f5518SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
793840b8ebdSBarry Smith #endif
794840b8ebdSBarry Smith       dx                *= epsilon;
7955904e1b1SBarry Smith       vscale_array[col] = 1.0/dx;
796840b8ebdSBarry Smith     }
7975904e1b1SBarry Smith   }
7985904e1b1SBarry Smith   vscale_array = vscale_array - start;ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
7995904e1b1SBarry Smith   ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8005904e1b1SBarry Smith   ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8015904e1b1SBarry Smith 
8025904e1b1SBarry Smith   if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow;
8035904e1b1SBarry Smith   else                        vscaleforrow = coloring->columnsforrow;
8045904e1b1SBarry Smith 
8055904e1b1SBarry Smith   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
8065904e1b1SBarry Smith   /*
8075904e1b1SBarry Smith       Loop over each color
8085904e1b1SBarry Smith   */
8095904e1b1SBarry Smith   for (k=0; k<coloring->ncolors; k++) {
8105904e1b1SBarry Smith     ierr = VecCopy(x1,w3);CHKERRQ(ierr);
8115904e1b1SBarry Smith     ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);w3_array = w3_array - start;
8125904e1b1SBarry Smith     /*
8135904e1b1SBarry Smith        Loop over each column associated with color adding the
8145904e1b1SBarry Smith        perturbation to the vector w3.
8155904e1b1SBarry Smith     */
8165904e1b1SBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
8175904e1b1SBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
8185904e1b1SBarry Smith       dx  = xx[col];
8195b8514ebSBarry Smith       if (dx == 0.0) dx = 1.0;
8205904e1b1SBarry Smith #if !defined(PETSC_USE_COMPLEX)
8215904e1b1SBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
8225904e1b1SBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
8235904e1b1SBarry Smith #else
8245904e1b1SBarry Smith       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
8255904e1b1SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
8265904e1b1SBarry Smith #endif
8275904e1b1SBarry Smith       dx            *= epsilon;
8285904e1b1SBarry Smith       w3_array[col] += dx;
8295904e1b1SBarry Smith     }
8305904e1b1SBarry Smith     w3_array = w3_array + start; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
8315904e1b1SBarry Smith 
832840b8ebdSBarry Smith     /*
833840b8ebdSBarry Smith        Evaluate function at x1 + dx (here dx is a vector of perturbations)
834840b8ebdSBarry Smith     */
83566f9b7ceSBarry Smith     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
836840b8ebdSBarry Smith     ierr = (*f)(sctx,t,w3,w2,fctx);CHKERRQ(ierr);
83766f9b7ceSBarry Smith     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
838efb30889SBarry Smith     ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr);
8395904e1b1SBarry Smith 
840840b8ebdSBarry Smith     /*
841840b8ebdSBarry Smith        Loop over rows of vector, putting results into Jacobian matrix
842840b8ebdSBarry Smith     */
8433b28642cSBarry Smith     ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
844840b8ebdSBarry Smith     for (l=0; l<coloring->nrows[k]; l++) {
845840b8ebdSBarry Smith       row    = coloring->rows[k][l];
846840b8ebdSBarry Smith       col    = coloring->columnsforrow[k][l];
8475904e1b1SBarry Smith       y[row] *= vscale_array[vscaleforrow[k][l]];
848840b8ebdSBarry Smith       srow   = row + start;
849840b8ebdSBarry Smith       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
850840b8ebdSBarry Smith     }
8513b28642cSBarry Smith     ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
852840b8ebdSBarry Smith   }
8535904e1b1SBarry Smith   ierr  = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
8545904e1b1SBarry Smith   xx    = xx + start; ierr  = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
855840b8ebdSBarry Smith   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
856840b8ebdSBarry Smith   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
857d5ba7fb7SMatthew Knepley   ierr  = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
8583a40ed3dSBarry Smith   PetscFunctionReturn(0);
859840b8ebdSBarry Smith }
8603b28642cSBarry Smith 
8613b28642cSBarry Smith 
8623b28642cSBarry Smith 
863