xref: /petsc/src/mat/matfd/fdmatrix.c (revision 832b7ceb95bd65cfd99922d89fc8cdbc69c04de1)
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 
7af0996ceSBarry Smith #include <petsc/private/matimpl.h>        /*I "petscmat.h" I*/
8af0996ceSBarry Smith #include <petsc/private/isimpl.h>
9bbf0e169SBarry Smith 
104a2ae208SSatish Balay #undef __FUNCT__
113a7fca6bSBarry Smith #define __FUNCT__ "MatFDColoringSetF"
127087cfbeSBarry Smith PetscErrorCode  MatFDColoringSetF(MatFDColoring fd,Vec F)
133a7fca6bSBarry Smith {
144e269d77SPeter Brune   PetscErrorCode ierr;
156e111a19SKarl Rupp 
163a7fca6bSBarry Smith   PetscFunctionBegin;
174e269d77SPeter Brune   if (F) {
184e269d77SPeter Brune     ierr     = VecCopy(F,fd->w1);CHKERRQ(ierr);
194e269d77SPeter Brune     fd->fset = PETSC_TRUE;
204e269d77SPeter Brune   } else {
214e269d77SPeter Brune     fd->fset = PETSC_FALSE;
224e269d77SPeter Brune   }
233a7fca6bSBarry Smith   PetscFunctionReturn(0);
243a7fca6bSBarry Smith }
253a7fca6bSBarry Smith 
269804daf3SBarry Smith #include <petscdraw.h>
273a7fca6bSBarry Smith #undef __FUNCT__
284a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Draw_Zoom"
296849ba73SBarry Smith static PetscErrorCode MatFDColoringView_Draw_Zoom(PetscDraw draw,void *Aa)
309194eea9SBarry Smith {
319194eea9SBarry Smith   MatFDColoring  fd = (MatFDColoring)Aa;
32dfbe8321SBarry Smith   PetscErrorCode ierr;
33b312708bSHong Zhang   PetscInt       i,j,nz,row;
349194eea9SBarry Smith   PetscReal      x,y;
35b312708bSHong Zhang   MatEntry       *Jentry=fd->matentry;
369194eea9SBarry Smith 
379194eea9SBarry Smith   PetscFunctionBegin;
389194eea9SBarry Smith   /* loop over colors  */
39b312708bSHong Zhang   nz = 0;
409194eea9SBarry Smith   for (i=0; i<fd->ncolors; i++) {
419194eea9SBarry Smith     for (j=0; j<fd->nrows[i]; j++) {
42b312708bSHong Zhang       row = Jentry[nz].row;
43b312708bSHong Zhang       y   = fd->M - row - fd->rstart;
44b312708bSHong Zhang       x   = (PetscReal)Jentry[nz++].col;
45b0a32e0cSBarry Smith       ierr = PetscDrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);CHKERRQ(ierr);
469194eea9SBarry Smith     }
479194eea9SBarry Smith   }
489194eea9SBarry Smith   PetscFunctionReturn(0);
499194eea9SBarry Smith }
509194eea9SBarry Smith 
514a2ae208SSatish Balay #undef __FUNCT__
524a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Draw"
536849ba73SBarry Smith static PetscErrorCode MatFDColoringView_Draw(MatFDColoring fd,PetscViewer viewer)
54005c665bSBarry Smith {
55dfbe8321SBarry Smith   PetscErrorCode ierr;
56ace3abfcSBarry Smith   PetscBool      isnull;
57b0a32e0cSBarry Smith   PetscDraw      draw;
5854d96fa3SBarry Smith   PetscReal      xr,yr,xl,yl,h,w;
59005c665bSBarry Smith 
603a40ed3dSBarry Smith   PetscFunctionBegin;
61b0a32e0cSBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
62*832b7cebSLisandro Dalcin   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
63*832b7cebSLisandro Dalcin   if (isnull) PetscFunctionReturn(0);
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);
68*832b7cebSLisandro Dalcin   ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
69b0a32e0cSBarry Smith   ierr = PetscDrawZoom(draw,MatFDColoringView_Draw_Zoom,fd);CHKERRQ(ierr);
70f77a5242SKarl Rupp   ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",NULL);CHKERRQ(ierr);
71*832b7cebSLisandro Dalcin   ierr = PetscDrawSave(draw);CHKERRQ(ierr);
723a40ed3dSBarry Smith   PetscFunctionReturn(0);
73005c665bSBarry Smith }
74005c665bSBarry Smith 
754a2ae208SSatish Balay #undef __FUNCT__
764a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView"
77bbf0e169SBarry Smith /*@C
78639f9d9dSBarry Smith    MatFDColoringView - Views a finite difference coloring context.
79bbf0e169SBarry Smith 
804c49b128SBarry Smith    Collective on MatFDColoring
81fee21e36SBarry Smith 
82ef5ee4d1SLois Curfman McInnes    Input  Parameters:
83ef5ee4d1SLois Curfman McInnes +  c - the coloring context
84ef5ee4d1SLois Curfman McInnes -  viewer - visualization context
85ef5ee4d1SLois Curfman McInnes 
8615091d37SBarry Smith    Level: intermediate
8715091d37SBarry Smith 
88b4fc646aSLois Curfman McInnes    Notes:
89b4fc646aSLois Curfman McInnes    The available visualization contexts include
90b0a32e0cSBarry Smith +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
91b0a32e0cSBarry Smith .     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
92ef5ee4d1SLois Curfman McInnes         output where only the first processor opens
93ef5ee4d1SLois Curfman McInnes         the file.  All other processors send their
94ef5ee4d1SLois Curfman McInnes         data to the first processor to print.
95b0a32e0cSBarry Smith -     PETSC_VIEWER_DRAW_WORLD - graphical display of nonzero structure
96bbf0e169SBarry Smith 
979194eea9SBarry Smith    Notes:
989194eea9SBarry Smith      Since PETSc uses only a small number of basic colors (currently 33), if the coloring
99b0a32e0cSBarry Smith    involves more than 33 then some seemingly identical colors are displayed making it look
1009194eea9SBarry Smith    like an illegal coloring. This is just a graphical artifact.
1019194eea9SBarry Smith 
102639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
103005c665bSBarry Smith 
104b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, view
105bbf0e169SBarry Smith @*/
1067087cfbeSBarry Smith PetscErrorCode  MatFDColoringView(MatFDColoring c,PetscViewer viewer)
107bbf0e169SBarry Smith {
1086849ba73SBarry Smith   PetscErrorCode    ierr;
10938baddfdSBarry Smith   PetscInt          i,j;
110ace3abfcSBarry Smith   PetscBool         isdraw,iascii;
111f3ef73ceSBarry Smith   PetscViewerFormat format;
112bbf0e169SBarry Smith 
1133a40ed3dSBarry Smith   PetscFunctionBegin;
1140700a824SBarry Smith   PetscValidHeaderSpecific(c,MAT_FDCOLORING_CLASSID,1);
1153050cee2SBarry Smith   if (!viewer) {
116ce94432eSBarry Smith     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)c),&viewer);CHKERRQ(ierr);
1173050cee2SBarry Smith   }
1180700a824SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
119c9780b6fSBarry Smith   PetscCheckSameComm(c,1,viewer,2);
120bbf0e169SBarry Smith 
121251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
122251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1230f5bd95cSBarry Smith   if (isdraw) {
124b4fc646aSLois Curfman McInnes     ierr = MatFDColoringView_Draw(c,viewer);CHKERRQ(ierr);
12532077d6dSBarry Smith   } else if (iascii) {
126dae58748SBarry Smith     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)c,viewer);CHKERRQ(ierr);
12757622a8eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Error tolerance=%g\n",(double)c->error_rel);CHKERRQ(ierr);
12857622a8eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Umin=%g\n",(double)c->umin);CHKERRQ(ierr);
12977431f27SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Number of colors=%D\n",c->ncolors);CHKERRQ(ierr);
130ae09f205SBarry Smith 
131b0a32e0cSBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
132fb9695e5SSatish Balay     if (format != PETSC_VIEWER_ASCII_INFO) {
133b312708bSHong Zhang       PetscInt row,col,nz;
134b312708bSHong Zhang       nz = 0;
135b4fc646aSLois Curfman McInnes       for (i=0; i<c->ncolors; i++) {
13677431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"  Information for color %D\n",i);CHKERRQ(ierr);
13777431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    Number of columns %D\n",c->ncolumns[i]);CHKERRQ(ierr);
138b4fc646aSLois Curfman McInnes         for (j=0; j<c->ncolumns[i]; j++) {
13977431f27SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"      %D\n",c->columns[i][j]);CHKERRQ(ierr);
140639f9d9dSBarry Smith         }
14177431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    Number of rows %D\n",c->nrows[i]);CHKERRQ(ierr);
142b4fc646aSLois Curfman McInnes         for (j=0; j<c->nrows[i]; j++) {
143b312708bSHong Zhang           row  = c->matentry[nz].row;
144b312708bSHong Zhang           col  = c->matentry[nz++].col;
145b312708bSHong Zhang           ierr = PetscViewerASCIIPrintf(viewer,"      %D %D \n",row,col);CHKERRQ(ierr);
146b4fc646aSLois Curfman McInnes         }
147bbf0e169SBarry Smith       }
148bbf0e169SBarry Smith     }
149b0a32e0cSBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
150bbf0e169SBarry Smith   }
1513a40ed3dSBarry Smith   PetscFunctionReturn(0);
152639f9d9dSBarry Smith }
153639f9d9dSBarry Smith 
1544a2ae208SSatish Balay #undef __FUNCT__
1554a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetParameters"
156639f9d9dSBarry Smith /*@
157b4fc646aSLois Curfman McInnes    MatFDColoringSetParameters - Sets the parameters for the sparse approximation of
158b4fc646aSLois Curfman McInnes    a Jacobian matrix using finite differences.
159639f9d9dSBarry Smith 
1603f9fe445SBarry Smith    Logically Collective on MatFDColoring
161ef5ee4d1SLois Curfman McInnes 
162ef5ee4d1SLois Curfman McInnes    The Jacobian is estimated with the differencing approximation
163ef5ee4d1SLois Curfman McInnes .vb
16465f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
165f23b5b22SLois Curfman McInnes        h = error_rel*u[i]                 if  abs(u[i]) > umin
166f23b5b22SLois Curfman McInnes          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
167ef5ee4d1SLois Curfman McInnes        dx_{i} = (0, ... 1, .... 0)
168ef5ee4d1SLois Curfman McInnes .ve
169639f9d9dSBarry Smith 
170639f9d9dSBarry Smith    Input Parameters:
171ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
172639f9d9dSBarry Smith .  error_rel - relative error
173f23b5b22SLois Curfman McInnes -  umin - minimum allowable u-value magnitude
174fee21e36SBarry Smith 
17515091d37SBarry Smith    Level: advanced
17615091d37SBarry Smith 
177b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, parameters
178b4fc646aSLois Curfman McInnes 
17995b89fc3SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringSetFromOptions()
18095b89fc3SBarry Smith 
181639f9d9dSBarry Smith @*/
1827087cfbeSBarry Smith PetscErrorCode MatFDColoringSetParameters(MatFDColoring matfd,PetscReal error,PetscReal umin)
183639f9d9dSBarry Smith {
1843a40ed3dSBarry Smith   PetscFunctionBegin;
1850700a824SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
186c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(matfd,error,2);
187c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(matfd,umin,3);
188639f9d9dSBarry Smith   if (error != PETSC_DEFAULT) matfd->error_rel = error;
189639f9d9dSBarry Smith   if (umin != PETSC_DEFAULT)  matfd->umin      = umin;
1903a40ed3dSBarry Smith   PetscFunctionReturn(0);
191639f9d9dSBarry Smith }
192639f9d9dSBarry Smith 
193f86b9fbaSHong Zhang #undef __FUNCT__
194f86b9fbaSHong Zhang #define __FUNCT__ "MatFDColoringSetBlockSize"
195f86b9fbaSHong Zhang /*@
196c8a9c622SHong Zhang    MatFDColoringSetBlockSize - Sets block size for efficient inserting entries of Jacobian matrix.
197005c665bSBarry Smith 
198f86b9fbaSHong Zhang    Logically Collective on MatFDColoring
199f86b9fbaSHong Zhang 
200f86b9fbaSHong Zhang    Input Parameters:
201f86b9fbaSHong Zhang +  coloring - the coloring context
202f86b9fbaSHong Zhang .  brows - number of rows in the block
203f86b9fbaSHong Zhang -  bcols - number of columns in the block
204f86b9fbaSHong Zhang 
205f86b9fbaSHong Zhang    Level: intermediate
206f86b9fbaSHong Zhang 
207f86b9fbaSHong Zhang .keywords: Mat, coloring
208f86b9fbaSHong Zhang 
209f86b9fbaSHong Zhang .seealso: MatFDColoringCreate(), MatFDColoringSetFromOptions()
210f86b9fbaSHong Zhang 
211f86b9fbaSHong Zhang @*/
212f86b9fbaSHong Zhang PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring matfd,PetscInt brows,PetscInt bcols)
213f86b9fbaSHong Zhang {
214f86b9fbaSHong Zhang   PetscFunctionBegin;
215f86b9fbaSHong Zhang   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
216f86b9fbaSHong Zhang   PetscValidLogicalCollectiveInt(matfd,brows,2);
217f86b9fbaSHong Zhang   PetscValidLogicalCollectiveInt(matfd,bcols,3);
218f86b9fbaSHong Zhang   if (brows != PETSC_DEFAULT) matfd->brows = brows;
219f86b9fbaSHong Zhang   if (bcols != PETSC_DEFAULT) matfd->bcols = bcols;
220f86b9fbaSHong Zhang   PetscFunctionReturn(0);
221f86b9fbaSHong Zhang }
222f86b9fbaSHong Zhang 
223f86b9fbaSHong Zhang #undef __FUNCT__
224f86b9fbaSHong Zhang #define __FUNCT__ "MatFDColoringSetUp"
225f86b9fbaSHong Zhang /*@
226f86b9fbaSHong Zhang    MatFDColoringSetUp - Sets up the internal data structures of matrix coloring context for the later use.
227f86b9fbaSHong Zhang 
228f86b9fbaSHong Zhang    Collective on Mat
229f86b9fbaSHong Zhang 
230f86b9fbaSHong Zhang    Input Parameters:
231f86b9fbaSHong Zhang +  mat - the matrix containing the nonzero structure of the Jacobian
232f86b9fbaSHong Zhang .  iscoloring - the coloring of the matrix; usually obtained with MatGetColoring() or DMCreateColoring()
233f86b9fbaSHong Zhang -  color - the matrix coloring context
234f86b9fbaSHong Zhang 
235f86b9fbaSHong Zhang    Level: beginner
236f86b9fbaSHong Zhang 
237f86b9fbaSHong Zhang .keywords: MatFDColoring, setup
238f86b9fbaSHong Zhang 
239f86b9fbaSHong Zhang .seealso: MatFDColoringCreate(), MatFDColoringDestroy()
240f86b9fbaSHong Zhang @*/
241f86b9fbaSHong Zhang PetscErrorCode MatFDColoringSetUp(Mat mat,ISColoring iscoloring,MatFDColoring color)
242f86b9fbaSHong Zhang {
243f86b9fbaSHong Zhang   PetscErrorCode ierr;
244f86b9fbaSHong Zhang 
245f86b9fbaSHong Zhang   PetscFunctionBegin;
246c8a9c622SHong Zhang   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
247c8a9c622SHong Zhang   PetscValidHeaderSpecific(color,MAT_FDCOLORING_CLASSID,3);
248c8a9c622SHong Zhang   if (color->setupcalled) PetscFunctionReturn(0);
249c8a9c622SHong Zhang 
2500df34763SHong Zhang   ierr = PetscLogEventBegin(MAT_FDColoringSetUp,mat,0,0,0);CHKERRQ(ierr);
251f86b9fbaSHong Zhang   if (mat->ops->fdcoloringsetup) {
252f86b9fbaSHong Zhang     ierr = (*mat->ops->fdcoloringsetup)(mat,iscoloring,color);CHKERRQ(ierr);
253f86b9fbaSHong Zhang   } else SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Code not yet written for matrix type %s",((PetscObject)mat)->type_name);
254c8a9c622SHong Zhang 
255c8a9c622SHong Zhang   color->setupcalled = PETSC_TRUE;
2560df34763SHong Zhang    ierr   = PetscLogEventEnd(MAT_FDColoringSetUp,mat,0,0,0);CHKERRQ(ierr);
257f86b9fbaSHong Zhang   PetscFunctionReturn(0);
258f86b9fbaSHong Zhang }
259ff0cfa39SBarry Smith 
2604a2ae208SSatish Balay #undef __FUNCT__
2614a9d489dSBarry Smith #define __FUNCT__ "MatFDColoringGetFunction"
2624a9d489dSBarry Smith /*@C
2634a9d489dSBarry Smith    MatFDColoringGetFunction - Gets the function to use for computing the Jacobian.
2644a9d489dSBarry Smith 
2653f9fe445SBarry Smith    Not Collective
2664a9d489dSBarry Smith 
2674a9d489dSBarry Smith    Input Parameters:
2684a9d489dSBarry Smith .  coloring - the coloring context
2694a9d489dSBarry Smith 
2704a9d489dSBarry Smith    Output Parameters:
2714a9d489dSBarry Smith +  f - the function
2724a9d489dSBarry Smith -  fctx - the optional user-defined function context
2734a9d489dSBarry Smith 
2744a9d489dSBarry Smith    Level: intermediate
2754a9d489dSBarry Smith 
2764a9d489dSBarry Smith .keywords: Mat, Jacobian, finite differences, set, function
27795b89fc3SBarry Smith 
27895b89fc3SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringSetFunction(), MatFDColoringSetFromOptions()
27995b89fc3SBarry Smith 
2804a9d489dSBarry Smith @*/
2817087cfbeSBarry Smith PetscErrorCode  MatFDColoringGetFunction(MatFDColoring matfd,PetscErrorCode (**f)(void),void **fctx)
2824a9d489dSBarry Smith {
2834a9d489dSBarry Smith   PetscFunctionBegin;
2840700a824SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
2854a9d489dSBarry Smith   if (f) *f = matfd->f;
2864a9d489dSBarry Smith   if (fctx) *fctx = matfd->fctx;
2874a9d489dSBarry Smith   PetscFunctionReturn(0);
2884a9d489dSBarry Smith }
2894a9d489dSBarry Smith 
2904a9d489dSBarry Smith #undef __FUNCT__
2914a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFunction"
292d64ed03dSBarry Smith /*@C
293005c665bSBarry Smith    MatFDColoringSetFunction - Sets the function to use for computing the Jacobian.
294005c665bSBarry Smith 
2953f9fe445SBarry Smith    Logically Collective on MatFDColoring
296fee21e36SBarry Smith 
297ef5ee4d1SLois Curfman McInnes    Input Parameters:
298ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
299ef5ee4d1SLois Curfman McInnes .  f - the function
300ef5ee4d1SLois Curfman McInnes -  fctx - the optional user-defined function context
301ef5ee4d1SLois Curfman McInnes 
3027850c7c0SBarry Smith    Calling sequence of (*f) function:
3037850c7c0SBarry Smith     For SNES:    PetscErrorCode (*f)(SNES,Vec,Vec,void*)
304ab637aeaSJed Brown     If not using SNES: PetscErrorCode (*f)(void *dummy,Vec,Vec,void*) and dummy is ignored
30515091d37SBarry Smith 
3067850c7c0SBarry Smith    Level: advanced
3077850c7c0SBarry Smith 
308ab637aeaSJed Brown    Notes: This function is usually used automatically by SNES (when one uses SNESSetJacobian() with the argument
3098d359177SBarry Smith      SNESComputeJacobianDefaultColor()) and only needs to be used by someone computing a matrix via coloring directly by
310ab637aeaSJed Brown      calling MatFDColoringApply()
3117850c7c0SBarry Smith 
3127850c7c0SBarry Smith    Fortran Notes:
3137850c7c0SBarry Smith     In Fortran you must call MatFDColoringSetFunction() for a coloring object to
314ab637aeaSJed Brown   be used without SNES or within the SNES solvers.
315f881d145SBarry Smith 
316b4fc646aSLois Curfman McInnes .keywords: Mat, Jacobian, finite differences, set, function
31795b89fc3SBarry Smith 
31895b89fc3SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringGetFunction(), MatFDColoringSetFromOptions()
31995b89fc3SBarry Smith 
320005c665bSBarry Smith @*/
3217087cfbeSBarry Smith PetscErrorCode  MatFDColoringSetFunction(MatFDColoring matfd,PetscErrorCode (*f)(void),void *fctx)
322005c665bSBarry Smith {
3233a40ed3dSBarry Smith   PetscFunctionBegin;
3240700a824SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
325005c665bSBarry Smith   matfd->f    = f;
326005c665bSBarry Smith   matfd->fctx = fctx;
3273a40ed3dSBarry Smith   PetscFunctionReturn(0);
328005c665bSBarry Smith }
329005c665bSBarry Smith 
3304a2ae208SSatish Balay #undef __FUNCT__
3314a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFromOptions"
332639f9d9dSBarry Smith /*@
333b4fc646aSLois Curfman McInnes    MatFDColoringSetFromOptions - Sets coloring finite difference parameters from
334639f9d9dSBarry Smith    the options database.
335639f9d9dSBarry Smith 
336fee21e36SBarry Smith    Collective on MatFDColoring
337fee21e36SBarry Smith 
33865f2ba5bSLois Curfman McInnes    The Jacobian, F'(u), is estimated with the differencing approximation
339ef5ee4d1SLois Curfman McInnes .vb
34065f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
341f23b5b22SLois Curfman McInnes        h = error_rel*u[i]                 if  abs(u[i]) > umin
342f23b5b22SLois Curfman McInnes          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
343ef5ee4d1SLois Curfman McInnes        dx_{i} = (0, ... 1, .... 0)
344ef5ee4d1SLois Curfman McInnes .ve
345ef5ee4d1SLois Curfman McInnes 
346ef5ee4d1SLois Curfman McInnes    Input Parameter:
347ef5ee4d1SLois Curfman McInnes .  coloring - the coloring context
348ef5ee4d1SLois Curfman McInnes 
349b4fc646aSLois Curfman McInnes    Options Database Keys:
3505620d6dcSBarry Smith +  -mat_fd_coloring_err <err> - Sets <err> (square root of relative error in the function)
351f23b5b22SLois Curfman McInnes .  -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude
3523ec795f1SBarry Smith .  -mat_fd_type - "wp" or "ds" (see MATMFFD_WP or MATMFFD_DS)
353ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_view - Activates basic viewing
354e350db43SBarry Smith .  -mat_fd_coloring_view ::ascii_info - Activates viewing info
355e350db43SBarry Smith -  -mat_fd_coloring_view draw - Activates drawing
356639f9d9dSBarry Smith 
35715091d37SBarry Smith     Level: intermediate
35815091d37SBarry Smith 
359b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, parameters
360d1c39f55SBarry Smith 
361d1c39f55SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringView(), MatFDColoringSetParameters()
362d1c39f55SBarry Smith 
363639f9d9dSBarry Smith @*/
3647087cfbeSBarry Smith PetscErrorCode  MatFDColoringSetFromOptions(MatFDColoring matfd)
365639f9d9dSBarry Smith {
366dfbe8321SBarry Smith   PetscErrorCode ierr;
367ace3abfcSBarry Smith   PetscBool      flg;
368efb30889SBarry Smith   char           value[3];
3693a40ed3dSBarry Smith 
3703a40ed3dSBarry Smith   PetscFunctionBegin;
3710700a824SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
372639f9d9dSBarry Smith 
3733194b578SJed Brown   ierr = PetscObjectOptionsBegin((PetscObject)matfd);CHKERRQ(ierr);
37487828ca2SBarry Smith   ierr = PetscOptionsReal("-mat_fd_coloring_err","Square root of relative error in function","MatFDColoringSetParameters",matfd->error_rel,&matfd->error_rel,0);CHKERRQ(ierr);
37587828ca2SBarry Smith   ierr = PetscOptionsReal("-mat_fd_coloring_umin","Minimum allowable u magnitude","MatFDColoringSetParameters",matfd->umin,&matfd->umin,0);CHKERRQ(ierr);
3761bc75a8dSBarry Smith   ierr = PetscOptionsString("-mat_fd_type","Algorithm to compute h, wp or ds","MatFDColoringCreate",matfd->htype,value,3,&flg);CHKERRQ(ierr);
377efb30889SBarry Smith   if (flg) {
378efb30889SBarry Smith     if (value[0] == 'w' && value[1] == 'p') matfd->htype = "wp";
379efb30889SBarry Smith     else if (value[0] == 'd' && value[1] == 's') matfd->htype = "ds";
380e32f2f54SBarry Smith     else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Unknown finite differencing type %s",value);
381efb30889SBarry Smith   }
382f86b9fbaSHong Zhang   ierr = PetscOptionsInt("-mat_fd_coloring_brows","Number of block rows","MatFDColoringSetBlockSize",matfd->brows,&matfd->brows,NULL);CHKERRQ(ierr);
38393dfae19SHong Zhang   ierr = PetscOptionsInt("-mat_fd_coloring_bcols","Number of block columns","MatFDColoringSetBlockSize",matfd->bcols,&matfd->bcols,&flg);CHKERRQ(ierr);
38493dfae19SHong Zhang   if (flg && matfd->bcols > matfd->ncolors) {
38593dfae19SHong Zhang     /* input bcols cannot be > matfd->ncolors, thus set it as ncolors */
38693dfae19SHong Zhang     matfd->bcols = matfd->ncolors;
38793dfae19SHong Zhang   }
388f86b9fbaSHong Zhang 
3895d973c19SBarry Smith   /* process any options handlers added with PetscObjectAddOptionsHandler() */
3900633abcbSJed Brown   ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)matfd);CHKERRQ(ierr);
391b0a32e0cSBarry Smith   PetscOptionsEnd();CHKERRQ(ierr);
3923a40ed3dSBarry Smith   PetscFunctionReturn(0);
393005c665bSBarry Smith }
394005c665bSBarry Smith 
3954a2ae208SSatish Balay #undef __FUNCT__
396e350db43SBarry Smith #define __FUNCT__ "MatFDColoringViewFromOptions"
397146574abSBarry Smith PetscErrorCode MatFDColoringViewFromOptions(MatFDColoring fd,const char prefix[],const char optionname[])
398005c665bSBarry Smith {
399dfbe8321SBarry Smith   PetscErrorCode    ierr;
400e350db43SBarry Smith   PetscBool         flg;
4013050cee2SBarry Smith   PetscViewer       viewer;
402cffb1e40SBarry Smith   PetscViewerFormat format;
403005c665bSBarry Smith 
4043a40ed3dSBarry Smith   PetscFunctionBegin;
405146574abSBarry Smith   if (prefix) {
406146574abSBarry Smith     ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)fd),prefix,optionname,&viewer,&format,&flg);CHKERRQ(ierr);
407146574abSBarry Smith   } else {
408ce94432eSBarry Smith     ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)fd),((PetscObject)fd)->prefix,optionname,&viewer,&format,&flg);CHKERRQ(ierr);
409146574abSBarry Smith   }
410005c665bSBarry Smith   if (flg) {
411cffb1e40SBarry Smith     ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr);
4123050cee2SBarry Smith     ierr = MatFDColoringView(fd,viewer);CHKERRQ(ierr);
413cffb1e40SBarry Smith     ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
414cffb1e40SBarry Smith     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
415005c665bSBarry Smith   }
4163a40ed3dSBarry Smith   PetscFunctionReturn(0);
417bbf0e169SBarry Smith }
418bbf0e169SBarry Smith 
4194a2ae208SSatish Balay #undef __FUNCT__
4204a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringCreate"
42105869f15SSatish Balay /*@
422639f9d9dSBarry Smith    MatFDColoringCreate - Creates a matrix coloring context for finite difference
423639f9d9dSBarry Smith    computation of Jacobians.
424bbf0e169SBarry Smith 
425ef5ee4d1SLois Curfman McInnes    Collective on Mat
426ef5ee4d1SLois Curfman McInnes 
427639f9d9dSBarry Smith    Input Parameters:
428ef5ee4d1SLois Curfman McInnes +  mat - the matrix containing the nonzero structure of the Jacobian
4297086a01eSPeter Brune -  iscoloring - the coloring of the matrix; usually obtained with MatColoringCreate() or DMCreateColoring()
430bbf0e169SBarry Smith 
431bbf0e169SBarry Smith     Output Parameter:
432639f9d9dSBarry Smith .   color - the new coloring context
433bbf0e169SBarry Smith 
43415091d37SBarry Smith     Level: intermediate
43515091d37SBarry Smith 
4368d359177SBarry Smith .seealso: MatFDColoringDestroy(),SNESComputeJacobianDefaultColor(), ISColoringCreate(),
437d1c39f55SBarry Smith           MatFDColoringSetFunction(), MatFDColoringSetFromOptions(), MatFDColoringApply(),
4387086a01eSPeter Brune           MatFDColoringView(), MatFDColoringSetParameters(), MatColoringCreate(), DMCreateColoring()
439bbf0e169SBarry Smith @*/
4407087cfbeSBarry Smith PetscErrorCode  MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
441bbf0e169SBarry Smith {
442639f9d9dSBarry Smith   MatFDColoring  c;
443639f9d9dSBarry Smith   MPI_Comm       comm;
444dfbe8321SBarry Smith   PetscErrorCode ierr;
44538baddfdSBarry Smith   PetscInt       M,N;
446639f9d9dSBarry Smith 
4473a40ed3dSBarry Smith   PetscFunctionBegin;
448c8a9c622SHong Zhang   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
449f141eedfSHong Zhang   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Matrix must be assembled by calls to MatAssemblyBegin/End();");
450d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
451639f9d9dSBarry Smith   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
452ce94432eSBarry Smith   if (M != N) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Only for square matrices");
453f881d145SBarry Smith   ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
45473107ff1SLisandro Dalcin   ierr = PetscHeaderCreate(c,MAT_FDCOLORING_CLASSID,"MatFDColoring","Jacobian computation via finite differences with coloring","Mat",comm,MatFDColoringDestroy,MatFDColoringView);CHKERRQ(ierr);
455639f9d9dSBarry Smith 
456b8f8c88eSHong Zhang   c->ctype = iscoloring->ctype;
457b8f8c88eSHong Zhang 
45893dfae19SHong Zhang   if (mat->ops->fdcoloringcreate) {
45993dfae19SHong Zhang     ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);CHKERRQ(ierr);
46093dfae19SHong Zhang   } else SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Code not yet written for matrix type %s",((PetscObject)mat)->type_name);
46193dfae19SHong Zhang 
4622a7a6963SBarry Smith   ierr = MatCreateVecs(mat,NULL,&c->w1);CHKERRQ(ierr);
4633bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)c,(PetscObject)c->w1);CHKERRQ(ierr);
464b8f8c88eSHong Zhang   ierr = VecDuplicate(c->w1,&c->w2);CHKERRQ(ierr);
4653bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)c,(PetscObject)c->w2);CHKERRQ(ierr);
466b8f8c88eSHong Zhang 
46777d8c4bbSBarry Smith   c->error_rel    = PETSC_SQRT_MACHINE_EPSILON;
46877d8c4bbSBarry Smith   c->umin         = 100.0*PETSC_SQRT_MACHINE_EPSILON;
46949b058dcSBarry Smith   c->currentcolor = -1;
470efb30889SBarry Smith   c->htype        = "wp";
4714e269d77SPeter Brune   c->fset         = PETSC_FALSE;
472c8a9c622SHong Zhang   c->setupcalled  = PETSC_FALSE;
473639f9d9dSBarry Smith 
474639f9d9dSBarry Smith   *color = c;
4754e269d77SPeter Brune   ierr   = PetscObjectCompose((PetscObject)mat,"SNESMatFDColoring",(PetscObject)c);CHKERRQ(ierr);
476d5ba7fb7SMatthew Knepley   ierr   = PetscLogEventEnd(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
4773a40ed3dSBarry Smith   PetscFunctionReturn(0);
478bbf0e169SBarry Smith }
479bbf0e169SBarry Smith 
4804a2ae208SSatish Balay #undef __FUNCT__
4814a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringDestroy"
48205869f15SSatish Balay /*@
483639f9d9dSBarry Smith     MatFDColoringDestroy - Destroys a matrix coloring context that was created
484639f9d9dSBarry Smith     via MatFDColoringCreate().
485bbf0e169SBarry Smith 
486ef5ee4d1SLois Curfman McInnes     Collective on MatFDColoring
487ef5ee4d1SLois Curfman McInnes 
488b4fc646aSLois Curfman McInnes     Input Parameter:
489639f9d9dSBarry Smith .   c - coloring context
490bbf0e169SBarry Smith 
49115091d37SBarry Smith     Level: intermediate
49215091d37SBarry Smith 
493639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
494bbf0e169SBarry Smith @*/
4956bf464f9SBarry Smith PetscErrorCode  MatFDColoringDestroy(MatFDColoring *c)
496bbf0e169SBarry Smith {
4976849ba73SBarry Smith   PetscErrorCode ierr;
49838baddfdSBarry Smith   PetscInt       i;
4990df34763SHong Zhang   MatFDColoring  color = *c;
500bbf0e169SBarry Smith 
5013a40ed3dSBarry Smith   PetscFunctionBegin;
5026bf464f9SBarry Smith   if (!*c) PetscFunctionReturn(0);
5030df34763SHong Zhang   if (--((PetscObject)color)->refct > 0) {*c = 0; PetscFunctionReturn(0);}
504d4bb536fSBarry Smith 
5050df34763SHong Zhang   for (i=0; i<color->ncolors; i++) {
5060df34763SHong Zhang     ierr = PetscFree(color->columns[i]);CHKERRQ(ierr);
507bbf0e169SBarry Smith   }
5080df34763SHong Zhang   ierr = PetscFree(color->ncolumns);CHKERRQ(ierr);
5090df34763SHong Zhang   ierr = PetscFree(color->columns);CHKERRQ(ierr);
5100df34763SHong Zhang   ierr = PetscFree(color->nrows);CHKERRQ(ierr);
5110df34763SHong Zhang   if (color->htype[0] == 'w') {
5120df34763SHong Zhang     ierr = PetscFree(color->matentry2);CHKERRQ(ierr);
5130df34763SHong Zhang   } else {
5140df34763SHong Zhang     ierr = PetscFree(color->matentry);CHKERRQ(ierr);
5150df34763SHong Zhang   }
5160df34763SHong Zhang   ierr = PetscFree(color->dy);CHKERRQ(ierr);
5170df34763SHong Zhang   if (color->vscale) {ierr = VecDestroy(&color->vscale);CHKERRQ(ierr);}
5180df34763SHong Zhang   ierr = VecDestroy(&color->w1);CHKERRQ(ierr);
5190df34763SHong Zhang   ierr = VecDestroy(&color->w2);CHKERRQ(ierr);
5200df34763SHong Zhang   ierr = VecDestroy(&color->w3);CHKERRQ(ierr);
521d38fa0fbSBarry Smith   ierr = PetscHeaderDestroy(c);CHKERRQ(ierr);
5223a40ed3dSBarry Smith   PetscFunctionReturn(0);
523bbf0e169SBarry Smith }
52443a90d84SBarry Smith 
5254a2ae208SSatish Balay #undef __FUNCT__
52649b058dcSBarry Smith #define __FUNCT__ "MatFDColoringGetPerturbedColumns"
52749b058dcSBarry Smith /*@C
52849b058dcSBarry Smith     MatFDColoringGetPerturbedColumns - Returns the indices of the columns that
52949b058dcSBarry Smith       that are currently being perturbed.
53049b058dcSBarry Smith 
53149b058dcSBarry Smith     Not Collective
53249b058dcSBarry Smith 
53349b058dcSBarry Smith     Input Parameters:
53449b058dcSBarry Smith .   coloring - coloring context created with MatFDColoringCreate()
53549b058dcSBarry Smith 
53649b058dcSBarry Smith     Output Parameters:
53749b058dcSBarry Smith +   n - the number of local columns being perturbed
53849b058dcSBarry Smith -   cols - the column indices, in global numbering
53949b058dcSBarry Smith 
54049b058dcSBarry Smith    Level: intermediate
54149b058dcSBarry Smith 
54249b058dcSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringApply()
54349b058dcSBarry Smith 
54449b058dcSBarry Smith .keywords: coloring, Jacobian, finite differences
54549b058dcSBarry Smith @*/
5467087cfbeSBarry Smith PetscErrorCode  MatFDColoringGetPerturbedColumns(MatFDColoring coloring,PetscInt *n,PetscInt *cols[])
54749b058dcSBarry Smith {
54849b058dcSBarry Smith   PetscFunctionBegin;
54949b058dcSBarry Smith   if (coloring->currentcolor >= 0) {
55049b058dcSBarry Smith     *n    = coloring->ncolumns[coloring->currentcolor];
55149b058dcSBarry Smith     *cols = coloring->columns[coloring->currentcolor];
55249b058dcSBarry Smith   } else {
55349b058dcSBarry Smith     *n = 0;
55449b058dcSBarry Smith   }
55549b058dcSBarry Smith   PetscFunctionReturn(0);
55649b058dcSBarry Smith }
55749b058dcSBarry Smith 
55849b058dcSBarry Smith #undef __FUNCT__
5594a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringApply"
56043a90d84SBarry Smith /*@
561e0907662SLois Curfman McInnes     MatFDColoringApply - Given a matrix for which a MatFDColoring context
562e0907662SLois Curfman McInnes     has been created, computes the Jacobian for a function via finite differences.
56343a90d84SBarry Smith 
564fee21e36SBarry Smith     Collective on MatFDColoring
565fee21e36SBarry Smith 
566ef5ee4d1SLois Curfman McInnes     Input Parameters:
567ef5ee4d1SLois Curfman McInnes +   mat - location to store Jacobian
568ef5ee4d1SLois Curfman McInnes .   coloring - coloring context created with MatFDColoringCreate()
569ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
5707850c7c0SBarry Smith -   sctx - context required by function, if this is being used with the SNES solver then it is SNES object, otherwise it is null
571ef5ee4d1SLois Curfman McInnes 
5728bba8e72SBarry Smith     Options Database Keys:
573ebd3b9afSBarry Smith +    -mat_fd_type - "wp" or "ds"  (see MATMFFD_WP or MATMFFD_DS)
574b9722508SLois Curfman McInnes .    -mat_fd_coloring_view - Activates basic viewing or coloring
575e350db43SBarry Smith .    -mat_fd_coloring_view draw - Activates drawing of coloring
576e350db43SBarry Smith -    -mat_fd_coloring_view ::ascii_info - Activates viewing of coloring info
5778bba8e72SBarry Smith 
57815091d37SBarry Smith     Level: intermediate
57998d79826SSatish Balay 
5807850c7c0SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringSetFunction()
58143a90d84SBarry Smith 
58243a90d84SBarry Smith .keywords: coloring, Jacobian, finite differences
58343a90d84SBarry Smith @*/
584d1e9a80fSBarry Smith PetscErrorCode  MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,void *sctx)
58543a90d84SBarry Smith {
5863acb8795SBarry Smith   PetscErrorCode ierr;
587684f2004SHong Zhang   PetscBool      flg = PETSC_FALSE;
5883acb8795SBarry Smith 
5893acb8795SBarry Smith   PetscFunctionBegin;
5900700a824SBarry Smith   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
5910700a824SBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_CLASSID,2);
5920700a824SBarry Smith   PetscValidHeaderSpecific(x1,VEC_CLASSID,3);
593ce94432eSBarry Smith   if (!coloring->f) SETERRQ(PetscObjectComm((PetscObject)J),PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetFunction()");
594e32f2f54SBarry Smith   if (!J->ops->fdcoloringapply) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not supported for this matrix type %s",((PetscObject)J)->type_name);
595c8a9c622SHong Zhang   if (!coloring->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetUp()");
596684f2004SHong Zhang 
597684f2004SHong Zhang   ierr = MatSetUnfactored(J);CHKERRQ(ierr);
598c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)coloring)->options,NULL,"-mat_fd_coloring_dont_rezero",&flg,NULL);CHKERRQ(ierr);
599684f2004SHong Zhang   if (flg) {
600684f2004SHong Zhang     ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr);
601684f2004SHong Zhang   } else {
602684f2004SHong Zhang     PetscBool assembled;
603684f2004SHong Zhang     ierr = MatAssembled(J,&assembled);CHKERRQ(ierr);
604684f2004SHong Zhang     if (assembled) {
605684f2004SHong Zhang       ierr = MatZeroEntries(J);CHKERRQ(ierr);
606684f2004SHong Zhang     }
607684f2004SHong Zhang   }
608684f2004SHong Zhang 
6095922145eSHong Zhang   ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
610d1e9a80fSBarry Smith   ierr = (*J->ops->fdcoloringapply)(J,coloring,x1,sctx);CHKERRQ(ierr);
6115922145eSHong Zhang   ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
6123acb8795SBarry Smith   PetscFunctionReturn(0);
6133acb8795SBarry Smith }
614