xref: /petsc/src/mat/matfd/fdmatrix.c (revision 684f2004e99e45d6318bc6d6b71aeefd2c18360e)
1be1d678aSKris Buschelman 
2bbf0e169SBarry Smith /*
3639f9d9dSBarry Smith    This is where the abstract matrix operations are defined that are
4639f9d9dSBarry Smith   used for finite difference computations of Jacobians using coloring.
5bbf0e169SBarry Smith */
6bbf0e169SBarry Smith 
7b45d2f2cSJed Brown #include <petsc-private/matimpl.h>        /*I "petscmat.h" I*/
8bbf0e169SBarry Smith 
94a2ae208SSatish Balay #undef __FUNCT__
103a7fca6bSBarry Smith #define __FUNCT__ "MatFDColoringSetF"
117087cfbeSBarry Smith PetscErrorCode  MatFDColoringSetF(MatFDColoring fd,Vec F)
123a7fca6bSBarry Smith {
134e269d77SPeter Brune   PetscErrorCode ierr;
146e111a19SKarl Rupp 
153a7fca6bSBarry Smith   PetscFunctionBegin;
164e269d77SPeter Brune   if (F) {
174e269d77SPeter Brune     ierr     = VecCopy(F,fd->w1);CHKERRQ(ierr);
184e269d77SPeter Brune     fd->fset = PETSC_TRUE;
194e269d77SPeter Brune   } else {
204e269d77SPeter Brune     fd->fset = PETSC_FALSE;
214e269d77SPeter Brune   }
223a7fca6bSBarry Smith   PetscFunctionReturn(0);
233a7fca6bSBarry Smith }
243a7fca6bSBarry Smith 
259804daf3SBarry Smith #include <petscdraw.h>
263a7fca6bSBarry Smith #undef __FUNCT__
274a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Draw_Zoom"
286849ba73SBarry Smith static PetscErrorCode MatFDColoringView_Draw_Zoom(PetscDraw draw,void *Aa)
299194eea9SBarry Smith {
309194eea9SBarry Smith   MatFDColoring  fd = (MatFDColoring)Aa;
31dfbe8321SBarry Smith   PetscErrorCode ierr;
32b312708bSHong Zhang   PetscInt       i,j,nz,row;
339194eea9SBarry Smith   PetscReal      x,y;
34b312708bSHong Zhang   MatEntry       *Jentry=fd->matentry;
359194eea9SBarry Smith 
369194eea9SBarry Smith   PetscFunctionBegin;
379194eea9SBarry Smith   /* loop over colors  */
38b312708bSHong Zhang   nz = 0;
399194eea9SBarry Smith   for (i=0; i<fd->ncolors; i++) {
409194eea9SBarry Smith     for (j=0; j<fd->nrows[i]; j++) {
41b312708bSHong Zhang       row = Jentry[nz].row;
42b312708bSHong Zhang       y   = fd->M - row - fd->rstart;
43b312708bSHong Zhang       x   = (PetscReal)Jentry[nz++].col;
44b0a32e0cSBarry Smith       ierr = PetscDrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);CHKERRQ(ierr);
459194eea9SBarry Smith     }
469194eea9SBarry Smith   }
479194eea9SBarry Smith   PetscFunctionReturn(0);
489194eea9SBarry Smith }
499194eea9SBarry Smith 
504a2ae208SSatish Balay #undef __FUNCT__
514a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Draw"
526849ba73SBarry Smith static PetscErrorCode MatFDColoringView_Draw(MatFDColoring fd,PetscViewer viewer)
53005c665bSBarry Smith {
54dfbe8321SBarry Smith   PetscErrorCode ierr;
55ace3abfcSBarry Smith   PetscBool      isnull;
56b0a32e0cSBarry Smith   PetscDraw      draw;
5754d96fa3SBarry Smith   PetscReal      xr,yr,xl,yl,h,w;
58005c665bSBarry Smith 
593a40ed3dSBarry Smith   PetscFunctionBegin;
60b0a32e0cSBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
61b0a32e0cSBarry Smith   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
629194eea9SBarry Smith 
639194eea9SBarry Smith   ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
64005c665bSBarry Smith 
65005c665bSBarry Smith   xr   = fd->N; yr  = fd->M; h = yr/10.0; w = xr/10.0;
66005c665bSBarry Smith   xr  += w;     yr += h;    xl = -w;     yl = -h;
67b0a32e0cSBarry Smith   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
68b0a32e0cSBarry Smith   ierr = PetscDrawZoom(draw,MatFDColoringView_Draw_Zoom,fd);CHKERRQ(ierr);
69f77a5242SKarl Rupp   ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",NULL);CHKERRQ(ierr);
703a40ed3dSBarry Smith   PetscFunctionReturn(0);
71005c665bSBarry Smith }
72005c665bSBarry Smith 
734a2ae208SSatish Balay #undef __FUNCT__
744a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView"
75bbf0e169SBarry Smith /*@C
76639f9d9dSBarry Smith    MatFDColoringView - Views a finite difference coloring context.
77bbf0e169SBarry Smith 
784c49b128SBarry Smith    Collective on MatFDColoring
79fee21e36SBarry Smith 
80ef5ee4d1SLois Curfman McInnes    Input  Parameters:
81ef5ee4d1SLois Curfman McInnes +  c - the coloring context
82ef5ee4d1SLois Curfman McInnes -  viewer - visualization context
83ef5ee4d1SLois Curfman McInnes 
8415091d37SBarry Smith    Level: intermediate
8515091d37SBarry Smith 
86b4fc646aSLois Curfman McInnes    Notes:
87b4fc646aSLois Curfman McInnes    The available visualization contexts include
88b0a32e0cSBarry Smith +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
89b0a32e0cSBarry Smith .     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
90ef5ee4d1SLois Curfman McInnes         output where only the first processor opens
91ef5ee4d1SLois Curfman McInnes         the file.  All other processors send their
92ef5ee4d1SLois Curfman McInnes         data to the first processor to print.
93b0a32e0cSBarry Smith -     PETSC_VIEWER_DRAW_WORLD - graphical display of nonzero structure
94bbf0e169SBarry Smith 
959194eea9SBarry Smith    Notes:
969194eea9SBarry Smith      Since PETSc uses only a small number of basic colors (currently 33), if the coloring
97b0a32e0cSBarry Smith    involves more than 33 then some seemingly identical colors are displayed making it look
989194eea9SBarry Smith    like an illegal coloring. This is just a graphical artifact.
999194eea9SBarry Smith 
100639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
101005c665bSBarry Smith 
102b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, view
103bbf0e169SBarry Smith @*/
1047087cfbeSBarry Smith PetscErrorCode  MatFDColoringView(MatFDColoring c,PetscViewer viewer)
105bbf0e169SBarry Smith {
1066849ba73SBarry Smith   PetscErrorCode    ierr;
10738baddfdSBarry Smith   PetscInt          i,j;
108ace3abfcSBarry Smith   PetscBool         isdraw,iascii;
109f3ef73ceSBarry Smith   PetscViewerFormat format;
110bbf0e169SBarry Smith 
1113a40ed3dSBarry Smith   PetscFunctionBegin;
1120700a824SBarry Smith   PetscValidHeaderSpecific(c,MAT_FDCOLORING_CLASSID,1);
1133050cee2SBarry Smith   if (!viewer) {
114ce94432eSBarry Smith     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)c),&viewer);CHKERRQ(ierr);
1153050cee2SBarry Smith   }
1160700a824SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
117c9780b6fSBarry Smith   PetscCheckSameComm(c,1,viewer,2);
118bbf0e169SBarry Smith 
119251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
120251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1210f5bd95cSBarry Smith   if (isdraw) {
122b4fc646aSLois Curfman McInnes     ierr = MatFDColoringView_Draw(c,viewer);CHKERRQ(ierr);
12332077d6dSBarry Smith   } else if (iascii) {
124dae58748SBarry Smith     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)c,viewer);CHKERRQ(ierr);
125a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Error tolerance=%G\n",c->error_rel);CHKERRQ(ierr);
126a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Umin=%G\n",c->umin);CHKERRQ(ierr);
12777431f27SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Number of colors=%D\n",c->ncolors);CHKERRQ(ierr);
128ae09f205SBarry Smith 
129b0a32e0cSBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
130fb9695e5SSatish Balay     if (format != PETSC_VIEWER_ASCII_INFO) {
131b312708bSHong Zhang       PetscInt row,col,nz;
132b312708bSHong Zhang       nz = 0;
133b4fc646aSLois Curfman McInnes       for (i=0; i<c->ncolors; i++) {
13477431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"  Information for color %D\n",i);CHKERRQ(ierr);
13577431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    Number of columns %D\n",c->ncolumns[i]);CHKERRQ(ierr);
136b4fc646aSLois Curfman McInnes         for (j=0; j<c->ncolumns[i]; j++) {
13777431f27SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"      %D\n",c->columns[i][j]);CHKERRQ(ierr);
138639f9d9dSBarry Smith         }
13977431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    Number of rows %D\n",c->nrows[i]);CHKERRQ(ierr);
140b4fc646aSLois Curfman McInnes         for (j=0; j<c->nrows[i]; j++) {
141b312708bSHong Zhang           row  = c->matentry[nz].row;
142b312708bSHong Zhang           col  = c->matentry[nz++].col;
143b312708bSHong Zhang           ierr = PetscViewerASCIIPrintf(viewer,"      %D %D \n",row,col);CHKERRQ(ierr);
144b4fc646aSLois Curfman McInnes         }
145bbf0e169SBarry Smith       }
146bbf0e169SBarry Smith     }
147b0a32e0cSBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
148bbf0e169SBarry Smith   }
1493a40ed3dSBarry Smith   PetscFunctionReturn(0);
150639f9d9dSBarry Smith }
151639f9d9dSBarry Smith 
1524a2ae208SSatish Balay #undef __FUNCT__
1534a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetParameters"
154639f9d9dSBarry Smith /*@
155b4fc646aSLois Curfman McInnes    MatFDColoringSetParameters - Sets the parameters for the sparse approximation of
156b4fc646aSLois Curfman McInnes    a Jacobian matrix using finite differences.
157639f9d9dSBarry Smith 
1583f9fe445SBarry Smith    Logically Collective on MatFDColoring
159ef5ee4d1SLois Curfman McInnes 
160ef5ee4d1SLois Curfman McInnes    The Jacobian is estimated with the differencing approximation
161ef5ee4d1SLois Curfman McInnes .vb
16265f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
163f23b5b22SLois Curfman McInnes        h = error_rel*u[i]                 if  abs(u[i]) > umin
164f23b5b22SLois Curfman McInnes          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
165ef5ee4d1SLois Curfman McInnes        dx_{i} = (0, ... 1, .... 0)
166ef5ee4d1SLois Curfman McInnes .ve
167639f9d9dSBarry Smith 
168639f9d9dSBarry Smith    Input Parameters:
169ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
170639f9d9dSBarry Smith .  error_rel - relative error
171f23b5b22SLois Curfman McInnes -  umin - minimum allowable u-value magnitude
172fee21e36SBarry Smith 
17315091d37SBarry Smith    Level: advanced
17415091d37SBarry Smith 
175b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, parameters
176b4fc646aSLois Curfman McInnes 
17795b89fc3SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringSetFromOptions()
17895b89fc3SBarry Smith 
179639f9d9dSBarry Smith @*/
1807087cfbeSBarry Smith PetscErrorCode  MatFDColoringSetParameters(MatFDColoring matfd,PetscReal error,PetscReal umin)
181639f9d9dSBarry Smith {
1823a40ed3dSBarry Smith   PetscFunctionBegin;
1830700a824SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
184c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(matfd,error,2);
185c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(matfd,umin,3);
186639f9d9dSBarry Smith   if (error != PETSC_DEFAULT) matfd->error_rel = error;
187639f9d9dSBarry Smith   if (umin != PETSC_DEFAULT)  matfd->umin      = umin;
1883a40ed3dSBarry Smith   PetscFunctionReturn(0);
189639f9d9dSBarry Smith }
190639f9d9dSBarry Smith 
191005c665bSBarry Smith 
192ff0cfa39SBarry Smith 
1934a2ae208SSatish Balay #undef __FUNCT__
1944a9d489dSBarry Smith #define __FUNCT__ "MatFDColoringGetFunction"
1954a9d489dSBarry Smith /*@C
1964a9d489dSBarry Smith    MatFDColoringGetFunction - Gets the function to use for computing the Jacobian.
1974a9d489dSBarry Smith 
1983f9fe445SBarry Smith    Not Collective
1994a9d489dSBarry Smith 
2004a9d489dSBarry Smith    Input Parameters:
2014a9d489dSBarry Smith .  coloring - the coloring context
2024a9d489dSBarry Smith 
2034a9d489dSBarry Smith    Output Parameters:
2044a9d489dSBarry Smith +  f - the function
2054a9d489dSBarry Smith -  fctx - the optional user-defined function context
2064a9d489dSBarry Smith 
2074a9d489dSBarry Smith    Level: intermediate
2084a9d489dSBarry Smith 
2094a9d489dSBarry Smith .keywords: Mat, Jacobian, finite differences, set, function
21095b89fc3SBarry Smith 
21195b89fc3SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringSetFunction(), MatFDColoringSetFromOptions()
21295b89fc3SBarry Smith 
2134a9d489dSBarry Smith @*/
2147087cfbeSBarry Smith PetscErrorCode  MatFDColoringGetFunction(MatFDColoring matfd,PetscErrorCode (**f)(void),void **fctx)
2154a9d489dSBarry Smith {
2164a9d489dSBarry Smith   PetscFunctionBegin;
2170700a824SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
2184a9d489dSBarry Smith   if (f) *f = matfd->f;
2194a9d489dSBarry Smith   if (fctx) *fctx = matfd->fctx;
2204a9d489dSBarry Smith   PetscFunctionReturn(0);
2214a9d489dSBarry Smith }
2224a9d489dSBarry Smith 
2234a9d489dSBarry Smith #undef __FUNCT__
2244a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFunction"
225d64ed03dSBarry Smith /*@C
226005c665bSBarry Smith    MatFDColoringSetFunction - Sets the function to use for computing the Jacobian.
227005c665bSBarry Smith 
2283f9fe445SBarry Smith    Logically Collective on MatFDColoring
229fee21e36SBarry Smith 
230ef5ee4d1SLois Curfman McInnes    Input Parameters:
231ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
232ef5ee4d1SLois Curfman McInnes .  f - the function
233ef5ee4d1SLois Curfman McInnes -  fctx - the optional user-defined function context
234ef5ee4d1SLois Curfman McInnes 
2357850c7c0SBarry Smith    Calling sequence of (*f) function:
2367850c7c0SBarry Smith     For SNES:    PetscErrorCode (*f)(SNES,Vec,Vec,void*)
237ab637aeaSJed Brown     If not using SNES: PetscErrorCode (*f)(void *dummy,Vec,Vec,void*) and dummy is ignored
23815091d37SBarry Smith 
2397850c7c0SBarry Smith    Level: advanced
2407850c7c0SBarry Smith 
241ab637aeaSJed Brown    Notes: This function is usually used automatically by SNES (when one uses SNESSetJacobian() with the argument
2428d359177SBarry Smith      SNESComputeJacobianDefaultColor()) and only needs to be used by someone computing a matrix via coloring directly by
243ab637aeaSJed Brown      calling MatFDColoringApply()
2447850c7c0SBarry Smith 
2457850c7c0SBarry Smith    Fortran Notes:
2467850c7c0SBarry Smith     In Fortran you must call MatFDColoringSetFunction() for a coloring object to
247ab637aeaSJed Brown   be used without SNES or within the SNES solvers.
248f881d145SBarry Smith 
249b4fc646aSLois Curfman McInnes .keywords: Mat, Jacobian, finite differences, set, function
25095b89fc3SBarry Smith 
25195b89fc3SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringGetFunction(), MatFDColoringSetFromOptions()
25295b89fc3SBarry Smith 
253005c665bSBarry Smith @*/
2547087cfbeSBarry Smith PetscErrorCode  MatFDColoringSetFunction(MatFDColoring matfd,PetscErrorCode (*f)(void),void *fctx)
255005c665bSBarry Smith {
2563a40ed3dSBarry Smith   PetscFunctionBegin;
2570700a824SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
258005c665bSBarry Smith   matfd->f    = f;
259005c665bSBarry Smith   matfd->fctx = fctx;
2603a40ed3dSBarry Smith   PetscFunctionReturn(0);
261005c665bSBarry Smith }
262005c665bSBarry Smith 
2634a2ae208SSatish Balay #undef __FUNCT__
2644a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFromOptions"
265639f9d9dSBarry Smith /*@
266b4fc646aSLois Curfman McInnes    MatFDColoringSetFromOptions - Sets coloring finite difference parameters from
267639f9d9dSBarry Smith    the options database.
268639f9d9dSBarry Smith 
269fee21e36SBarry Smith    Collective on MatFDColoring
270fee21e36SBarry Smith 
27165f2ba5bSLois Curfman McInnes    The Jacobian, F'(u), is estimated with the differencing approximation
272ef5ee4d1SLois Curfman McInnes .vb
27365f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
274f23b5b22SLois Curfman McInnes        h = error_rel*u[i]                 if  abs(u[i]) > umin
275f23b5b22SLois Curfman McInnes          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
276ef5ee4d1SLois Curfman McInnes        dx_{i} = (0, ... 1, .... 0)
277ef5ee4d1SLois Curfman McInnes .ve
278ef5ee4d1SLois Curfman McInnes 
279ef5ee4d1SLois Curfman McInnes    Input Parameter:
280ef5ee4d1SLois Curfman McInnes .  coloring - the coloring context
281ef5ee4d1SLois Curfman McInnes 
282b4fc646aSLois Curfman McInnes    Options Database Keys:
283d15ffeeaSSatish Balay +  -mat_fd_coloring_err <err> - Sets <err> (square root
284ef5ee4d1SLois Curfman McInnes            of relative error in the function)
285f23b5b22SLois Curfman McInnes .  -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude
2863ec795f1SBarry Smith .  -mat_fd_type - "wp" or "ds" (see MATMFFD_WP or MATMFFD_DS)
287ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_view - Activates basic viewing
288e350db43SBarry Smith .  -mat_fd_coloring_view ::ascii_info - Activates viewing info
289e350db43SBarry Smith -  -mat_fd_coloring_view draw - Activates drawing
290639f9d9dSBarry Smith 
29115091d37SBarry Smith     Level: intermediate
29215091d37SBarry Smith 
293b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, parameters
294d1c39f55SBarry Smith 
295d1c39f55SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringView(), MatFDColoringSetParameters()
296d1c39f55SBarry Smith 
297639f9d9dSBarry Smith @*/
2987087cfbeSBarry Smith PetscErrorCode  MatFDColoringSetFromOptions(MatFDColoring matfd)
299639f9d9dSBarry Smith {
300dfbe8321SBarry Smith   PetscErrorCode ierr;
301ace3abfcSBarry Smith   PetscBool      flg;
302efb30889SBarry Smith   char           value[3];
3033a40ed3dSBarry Smith 
3043a40ed3dSBarry Smith   PetscFunctionBegin;
3050700a824SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
306639f9d9dSBarry Smith 
3073194b578SJed Brown   ierr = PetscObjectOptionsBegin((PetscObject)matfd);CHKERRQ(ierr);
30887828ca2SBarry Smith   ierr = PetscOptionsReal("-mat_fd_coloring_err","Square root of relative error in function","MatFDColoringSetParameters",matfd->error_rel,&matfd->error_rel,0);CHKERRQ(ierr);
30987828ca2SBarry Smith   ierr = PetscOptionsReal("-mat_fd_coloring_umin","Minimum allowable u magnitude","MatFDColoringSetParameters",matfd->umin,&matfd->umin,0);CHKERRQ(ierr);
3101bc75a8dSBarry Smith   ierr = PetscOptionsString("-mat_fd_type","Algorithm to compute h, wp or ds","MatFDColoringCreate",matfd->htype,value,3,&flg);CHKERRQ(ierr);
311efb30889SBarry Smith   if (flg) {
312efb30889SBarry Smith     if (value[0] == 'w' && value[1] == 'p') matfd->htype = "wp";
313efb30889SBarry Smith     else if (value[0] == 'd' && value[1] == 's') matfd->htype = "ds";
314e32f2f54SBarry Smith     else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Unknown finite differencing type %s",value);
315efb30889SBarry Smith   }
3165d973c19SBarry Smith   /* process any options handlers added with PetscObjectAddOptionsHandler() */
3175d973c19SBarry Smith   ierr = PetscObjectProcessOptionsHandlers((PetscObject)matfd);CHKERRQ(ierr);
318b0a32e0cSBarry Smith   PetscOptionsEnd();CHKERRQ(ierr);
3193a40ed3dSBarry Smith   PetscFunctionReturn(0);
320005c665bSBarry Smith }
321005c665bSBarry Smith 
3224a2ae208SSatish Balay #undef __FUNCT__
323e350db43SBarry Smith #define __FUNCT__ "MatFDColoringViewFromOptions"
324146574abSBarry Smith PetscErrorCode MatFDColoringViewFromOptions(MatFDColoring fd,const char prefix[],const char optionname[])
325005c665bSBarry Smith {
326dfbe8321SBarry Smith   PetscErrorCode    ierr;
327e350db43SBarry Smith   PetscBool         flg;
3283050cee2SBarry Smith   PetscViewer       viewer;
329cffb1e40SBarry Smith   PetscViewerFormat format;
330005c665bSBarry Smith 
3313a40ed3dSBarry Smith   PetscFunctionBegin;
332146574abSBarry Smith   if (prefix) {
333146574abSBarry Smith     ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)fd),prefix,optionname,&viewer,&format,&flg);CHKERRQ(ierr);
334146574abSBarry Smith   } else {
335ce94432eSBarry Smith     ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)fd),((PetscObject)fd)->prefix,optionname,&viewer,&format,&flg);CHKERRQ(ierr);
336146574abSBarry Smith   }
337005c665bSBarry Smith   if (flg) {
338cffb1e40SBarry Smith     ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr);
3393050cee2SBarry Smith     ierr = MatFDColoringView(fd,viewer);CHKERRQ(ierr);
340cffb1e40SBarry Smith     ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
341cffb1e40SBarry Smith     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
342005c665bSBarry Smith   }
3433a40ed3dSBarry Smith   PetscFunctionReturn(0);
344bbf0e169SBarry Smith }
345bbf0e169SBarry Smith 
3464a2ae208SSatish Balay #undef __FUNCT__
3474a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringCreate"
34805869f15SSatish Balay /*@
349639f9d9dSBarry Smith    MatFDColoringCreate - Creates a matrix coloring context for finite difference
350639f9d9dSBarry Smith    computation of Jacobians.
351bbf0e169SBarry Smith 
352ef5ee4d1SLois Curfman McInnes    Collective on Mat
353ef5ee4d1SLois Curfman McInnes 
354639f9d9dSBarry Smith    Input Parameters:
355ef5ee4d1SLois Curfman McInnes +  mat - the matrix containing the nonzero structure of the Jacobian
356e727c939SJed Brown -  iscoloring - the coloring of the matrix; usually obtained with MatGetColoring() or DMCreateColoring()
357bbf0e169SBarry Smith 
358bbf0e169SBarry Smith     Output Parameter:
359639f9d9dSBarry Smith .   color - the new coloring context
360bbf0e169SBarry Smith 
36115091d37SBarry Smith     Level: intermediate
36215091d37SBarry Smith 
3638d359177SBarry Smith .seealso: MatFDColoringDestroy(),SNESComputeJacobianDefaultColor(), ISColoringCreate(),
364d1c39f55SBarry Smith           MatFDColoringSetFunction(), MatFDColoringSetFromOptions(), MatFDColoringApply(),
365e727c939SJed Brown           MatFDColoringView(), MatFDColoringSetParameters(), MatGetColoring(), DMCreateColoring()
366bbf0e169SBarry Smith @*/
3677087cfbeSBarry Smith PetscErrorCode  MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
368bbf0e169SBarry Smith {
369639f9d9dSBarry Smith   MatFDColoring  c;
370639f9d9dSBarry Smith   MPI_Comm       comm;
371dfbe8321SBarry Smith   PetscErrorCode ierr;
37238baddfdSBarry Smith   PetscInt       M,N;
373639f9d9dSBarry Smith 
3743a40ed3dSBarry Smith   PetscFunctionBegin;
375f141eedfSHong Zhang   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Matrix must be assembled by calls to MatAssemblyBegin/End();");
376d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
377639f9d9dSBarry Smith   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
378ce94432eSBarry Smith   if (M != N) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Only for square matrices");
379f881d145SBarry Smith   ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
38067c2884eSBarry Smith   ierr = PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_CLASSID,"MatFDColoring","Jacobian computation via finite differences with coloring","Mat",comm,MatFDColoringDestroy,MatFDColoringView);CHKERRQ(ierr);
381639f9d9dSBarry Smith 
382b8f8c88eSHong Zhang   c->ctype = iscoloring->ctype;
383b8f8c88eSHong Zhang 
384f830108cSBarry Smith   if (mat->ops->fdcoloringcreate) {
385f830108cSBarry Smith     ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);CHKERRQ(ierr);
386ce94432eSBarry Smith   } else SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Code not yet written for matrix type %s",((PetscObject)mat)->type_name);
387639f9d9dSBarry Smith 
388f77a5242SKarl Rupp   ierr = MatGetVecs(mat,NULL,&c->w1);CHKERRQ(ierr);
3893bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)c,(PetscObject)c->w1);CHKERRQ(ierr);
390b8f8c88eSHong Zhang   ierr = VecDuplicate(c->w1,&c->w2);CHKERRQ(ierr);
3913bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)c,(PetscObject)c->w2);CHKERRQ(ierr);
392b8f8c88eSHong Zhang 
39377d8c4bbSBarry Smith   c->error_rel    = PETSC_SQRT_MACHINE_EPSILON;
39477d8c4bbSBarry Smith   c->umin         = 100.0*PETSC_SQRT_MACHINE_EPSILON;
39549b058dcSBarry Smith   c->currentcolor = -1;
396efb30889SBarry Smith   c->htype        = "wp";
3974e269d77SPeter Brune   c->fset         = PETSC_FALSE;
398639f9d9dSBarry Smith 
399639f9d9dSBarry Smith   *color = c;
4004e269d77SPeter Brune   ierr   = PetscObjectCompose((PetscObject)mat,"SNESMatFDColoring",(PetscObject)c);CHKERRQ(ierr);
401d5ba7fb7SMatthew Knepley   ierr   = PetscLogEventEnd(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
4023a40ed3dSBarry Smith   PetscFunctionReturn(0);
403bbf0e169SBarry Smith }
404bbf0e169SBarry Smith 
4054a2ae208SSatish Balay #undef __FUNCT__
4064a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringDestroy"
40705869f15SSatish Balay /*@
408639f9d9dSBarry Smith     MatFDColoringDestroy - Destroys a matrix coloring context that was created
409639f9d9dSBarry Smith     via MatFDColoringCreate().
410bbf0e169SBarry Smith 
411ef5ee4d1SLois Curfman McInnes     Collective on MatFDColoring
412ef5ee4d1SLois Curfman McInnes 
413b4fc646aSLois Curfman McInnes     Input Parameter:
414639f9d9dSBarry Smith .   c - coloring context
415bbf0e169SBarry Smith 
41615091d37SBarry Smith     Level: intermediate
41715091d37SBarry Smith 
418639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
419bbf0e169SBarry Smith @*/
4206bf464f9SBarry Smith PetscErrorCode  MatFDColoringDestroy(MatFDColoring *c)
421bbf0e169SBarry Smith {
4226849ba73SBarry Smith   PetscErrorCode ierr;
42338baddfdSBarry Smith   PetscInt       i;
424bbf0e169SBarry Smith 
4253a40ed3dSBarry Smith   PetscFunctionBegin;
4266bf464f9SBarry Smith   if (!*c) PetscFunctionReturn(0);
4276bf464f9SBarry Smith   if (--((PetscObject)(*c))->refct > 0) {*c = 0; PetscFunctionReturn(0);}
428d4bb536fSBarry Smith 
4296bf464f9SBarry Smith   for (i=0; i<(*c)->ncolors; i++) {
4306bf464f9SBarry Smith     ierr = PetscFree((*c)->columns[i]);CHKERRQ(ierr);
431bbf0e169SBarry Smith   }
4326bf464f9SBarry Smith   ierr = PetscFree((*c)->ncolumns);CHKERRQ(ierr);
4336bf464f9SBarry Smith   ierr = PetscFree((*c)->columns);CHKERRQ(ierr);
4346bf464f9SBarry Smith   ierr = PetscFree((*c)->nrows);CHKERRQ(ierr);
435573f477fSHong Zhang   ierr = PetscFree((*c)->matentry);CHKERRQ(ierr);
436b7799381SHong Zhang   ierr = PetscFree((*c)->dy);CHKERRQ(ierr);
4379e917edbSHong Zhang   if ((*c)->vscale) {ierr = VecDestroy(&(*c)->vscale);CHKERRQ(ierr);}
4386bf464f9SBarry Smith   ierr = VecDestroy(&(*c)->w1);CHKERRQ(ierr);
4396bf464f9SBarry Smith   ierr = VecDestroy(&(*c)->w2);CHKERRQ(ierr);
4406bf464f9SBarry Smith   ierr = VecDestroy(&(*c)->w3);CHKERRQ(ierr);
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 @*/
4667087cfbeSBarry Smith PetscErrorCode  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 
47849b058dcSBarry Smith #undef __FUNCT__
4794a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringApply"
48043a90d84SBarry Smith /*@
481e0907662SLois Curfman McInnes     MatFDColoringApply - Given a matrix for which a MatFDColoring context
482e0907662SLois Curfman McInnes     has been created, computes the Jacobian for a function via finite differences.
48343a90d84SBarry Smith 
484fee21e36SBarry Smith     Collective on MatFDColoring
485fee21e36SBarry Smith 
486ef5ee4d1SLois Curfman McInnes     Input Parameters:
487ef5ee4d1SLois Curfman McInnes +   mat - location to store Jacobian
488ef5ee4d1SLois Curfman McInnes .   coloring - coloring context created with MatFDColoringCreate()
489ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
4907850c7c0SBarry Smith -   sctx - context required by function, if this is being used with the SNES solver then it is SNES object, otherwise it is null
491ef5ee4d1SLois Curfman McInnes 
4928bba8e72SBarry Smith     Options Database Keys:
493ebd3b9afSBarry Smith +    -mat_fd_type - "wp" or "ds"  (see MATMFFD_WP or MATMFFD_DS)
494b9722508SLois Curfman McInnes .    -mat_fd_coloring_view - Activates basic viewing or coloring
495e350db43SBarry Smith .    -mat_fd_coloring_view draw - Activates drawing of coloring
496e350db43SBarry Smith -    -mat_fd_coloring_view ::ascii_info - Activates viewing of coloring info
4978bba8e72SBarry Smith 
49815091d37SBarry Smith     Level: intermediate
49998d79826SSatish Balay 
5007850c7c0SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringSetFunction()
50143a90d84SBarry Smith 
50243a90d84SBarry Smith .keywords: coloring, Jacobian, finite differences
50343a90d84SBarry Smith @*/
5047087cfbeSBarry Smith PetscErrorCode  MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
50543a90d84SBarry Smith {
5063acb8795SBarry Smith   PetscErrorCode ierr;
507*684f2004SHong Zhang   PetscBool      flg = PETSC_FALSE;
5083acb8795SBarry Smith 
5093acb8795SBarry Smith   PetscFunctionBegin;
5100700a824SBarry Smith   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
5110700a824SBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_CLASSID,2);
5120700a824SBarry Smith   PetscValidHeaderSpecific(x1,VEC_CLASSID,3);
513ce94432eSBarry Smith   if (!coloring->f) SETERRQ(PetscObjectComm((PetscObject)J),PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetFunction()");
514e32f2f54SBarry Smith   if (!J->ops->fdcoloringapply) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not supported for this matrix type %s",((PetscObject)J)->type_name);
515*684f2004SHong Zhang 
516*684f2004SHong Zhang   ierr = MatSetUnfactored(J);CHKERRQ(ierr);
517*684f2004SHong Zhang   ierr = PetscOptionsGetBool(NULL,"-mat_fd_coloring_dont_rezero",&flg,NULL);CHKERRQ(ierr);
518*684f2004SHong Zhang   if (flg) {
519*684f2004SHong Zhang     ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr);
520*684f2004SHong Zhang   } else {
521*684f2004SHong Zhang     PetscBool assembled;
522*684f2004SHong Zhang     ierr = MatAssembled(J,&assembled);CHKERRQ(ierr);
523*684f2004SHong Zhang     if (assembled) {
524*684f2004SHong Zhang       ierr = MatZeroEntries(J);CHKERRQ(ierr);
525*684f2004SHong Zhang     }
526*684f2004SHong Zhang   }
527*684f2004SHong Zhang 
5285922145eSHong Zhang   ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
5293acb8795SBarry Smith   ierr = (*J->ops->fdcoloringapply)(J,coloring,x1,flag,sctx);CHKERRQ(ierr);
5305922145eSHong Zhang   ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
5313acb8795SBarry Smith   PetscFunctionReturn(0);
5323acb8795SBarry Smith }
533