xref: /petsc/src/mat/matfd/fdmatrix.c (revision c8a9c62294a3909d4cb0ca2bea407dacfcbcf516)
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 
191f86b9fbaSHong Zhang #undef __FUNCT__
192f86b9fbaSHong Zhang #define __FUNCT__ "MatFDColoringSetBlockSize"
193f86b9fbaSHong Zhang /*@
194*c8a9c622SHong Zhang    MatFDColoringSetBlockSize - Sets block size for efficient inserting entries of Jacobian matrix.
195005c665bSBarry Smith 
196f86b9fbaSHong Zhang    Logically Collective on MatFDColoring
197f86b9fbaSHong Zhang 
198f86b9fbaSHong Zhang    Input Parameters:
199f86b9fbaSHong Zhang +  coloring - the coloring context
200f86b9fbaSHong Zhang .  brows - number of rows in the block
201f86b9fbaSHong Zhang -  bcols - number of columns in the block
202f86b9fbaSHong Zhang 
203f86b9fbaSHong Zhang    Level: intermediate
204f86b9fbaSHong Zhang 
205f86b9fbaSHong Zhang .keywords: Mat, coloring
206f86b9fbaSHong Zhang 
207f86b9fbaSHong Zhang .seealso: MatFDColoringCreate(), MatFDColoringSetFromOptions()
208f86b9fbaSHong Zhang 
209f86b9fbaSHong Zhang @*/
210f86b9fbaSHong Zhang PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring matfd,PetscInt brows,PetscInt bcols)
211f86b9fbaSHong Zhang {
212f86b9fbaSHong Zhang   PetscFunctionBegin;
213f86b9fbaSHong Zhang   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
214f86b9fbaSHong Zhang   PetscValidLogicalCollectiveInt(matfd,brows,2);
215f86b9fbaSHong Zhang   PetscValidLogicalCollectiveInt(matfd,bcols,3);
216f86b9fbaSHong Zhang   if (brows != PETSC_DEFAULT) matfd->brows = brows;
217f86b9fbaSHong Zhang   if (bcols != PETSC_DEFAULT) matfd->bcols = bcols;
218f86b9fbaSHong Zhang   PetscFunctionReturn(0);
219f86b9fbaSHong Zhang }
220f86b9fbaSHong Zhang 
221f86b9fbaSHong Zhang #undef __FUNCT__
222f86b9fbaSHong Zhang #define __FUNCT__ "MatFDColoringSetUp"
223f86b9fbaSHong Zhang /*@
224f86b9fbaSHong Zhang    MatFDColoringSetUp - Sets up the internal data structures of matrix coloring context for the later use.
225f86b9fbaSHong Zhang 
226f86b9fbaSHong Zhang    Collective on Mat
227f86b9fbaSHong Zhang 
228f86b9fbaSHong Zhang    Input Parameters:
229f86b9fbaSHong Zhang +  mat - the matrix containing the nonzero structure of the Jacobian
230f86b9fbaSHong Zhang .  iscoloring - the coloring of the matrix; usually obtained with MatGetColoring() or DMCreateColoring()
231f86b9fbaSHong Zhang -  color - the matrix coloring context
232f86b9fbaSHong Zhang 
233f86b9fbaSHong Zhang    Level: beginner
234f86b9fbaSHong Zhang 
235f86b9fbaSHong Zhang .keywords: MatFDColoring, setup
236f86b9fbaSHong Zhang 
237f86b9fbaSHong Zhang .seealso: MatFDColoringCreate(), MatFDColoringDestroy()
238f86b9fbaSHong Zhang @*/
239f86b9fbaSHong Zhang PetscErrorCode MatFDColoringSetUp(Mat mat,ISColoring iscoloring,MatFDColoring color)
240f86b9fbaSHong Zhang {
241f86b9fbaSHong Zhang   PetscErrorCode ierr;
242f86b9fbaSHong Zhang 
243f86b9fbaSHong Zhang   PetscFunctionBegin;
244*c8a9c622SHong Zhang   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
245*c8a9c622SHong Zhang   PetscValidHeaderSpecific(color,MAT_FDCOLORING_CLASSID,3);
246*c8a9c622SHong Zhang   if (color->setupcalled) PetscFunctionReturn(0);
247*c8a9c622SHong Zhang 
248f86b9fbaSHong Zhang   if (mat->ops->fdcoloringsetup) {
249f86b9fbaSHong Zhang     ierr = (*mat->ops->fdcoloringsetup)(mat,iscoloring,color);CHKERRQ(ierr);
250f86b9fbaSHong Zhang   } else SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Code not yet written for matrix type %s",((PetscObject)mat)->type_name);
251*c8a9c622SHong Zhang 
252*c8a9c622SHong Zhang   color->setupcalled = PETSC_TRUE;
253f86b9fbaSHong Zhang   PetscFunctionReturn(0);
254f86b9fbaSHong Zhang }
255ff0cfa39SBarry Smith 
2564a2ae208SSatish Balay #undef __FUNCT__
2574a9d489dSBarry Smith #define __FUNCT__ "MatFDColoringGetFunction"
2584a9d489dSBarry Smith /*@C
2594a9d489dSBarry Smith    MatFDColoringGetFunction - Gets the function to use for computing the Jacobian.
2604a9d489dSBarry Smith 
2613f9fe445SBarry Smith    Not Collective
2624a9d489dSBarry Smith 
2634a9d489dSBarry Smith    Input Parameters:
2644a9d489dSBarry Smith .  coloring - the coloring context
2654a9d489dSBarry Smith 
2664a9d489dSBarry Smith    Output Parameters:
2674a9d489dSBarry Smith +  f - the function
2684a9d489dSBarry Smith -  fctx - the optional user-defined function context
2694a9d489dSBarry Smith 
2704a9d489dSBarry Smith    Level: intermediate
2714a9d489dSBarry Smith 
2724a9d489dSBarry Smith .keywords: Mat, Jacobian, finite differences, set, function
27395b89fc3SBarry Smith 
27495b89fc3SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringSetFunction(), MatFDColoringSetFromOptions()
27595b89fc3SBarry Smith 
2764a9d489dSBarry Smith @*/
2777087cfbeSBarry Smith PetscErrorCode  MatFDColoringGetFunction(MatFDColoring matfd,PetscErrorCode (**f)(void),void **fctx)
2784a9d489dSBarry Smith {
2794a9d489dSBarry Smith   PetscFunctionBegin;
2800700a824SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
2814a9d489dSBarry Smith   if (f) *f = matfd->f;
2824a9d489dSBarry Smith   if (fctx) *fctx = matfd->fctx;
2834a9d489dSBarry Smith   PetscFunctionReturn(0);
2844a9d489dSBarry Smith }
2854a9d489dSBarry Smith 
2864a9d489dSBarry Smith #undef __FUNCT__
2874a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFunction"
288d64ed03dSBarry Smith /*@C
289005c665bSBarry Smith    MatFDColoringSetFunction - Sets the function to use for computing the Jacobian.
290005c665bSBarry Smith 
2913f9fe445SBarry Smith    Logically Collective on MatFDColoring
292fee21e36SBarry Smith 
293ef5ee4d1SLois Curfman McInnes    Input Parameters:
294ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
295ef5ee4d1SLois Curfman McInnes .  f - the function
296ef5ee4d1SLois Curfman McInnes -  fctx - the optional user-defined function context
297ef5ee4d1SLois Curfman McInnes 
2987850c7c0SBarry Smith    Calling sequence of (*f) function:
2997850c7c0SBarry Smith     For SNES:    PetscErrorCode (*f)(SNES,Vec,Vec,void*)
300ab637aeaSJed Brown     If not using SNES: PetscErrorCode (*f)(void *dummy,Vec,Vec,void*) and dummy is ignored
30115091d37SBarry Smith 
3027850c7c0SBarry Smith    Level: advanced
3037850c7c0SBarry Smith 
304ab637aeaSJed Brown    Notes: This function is usually used automatically by SNES (when one uses SNESSetJacobian() with the argument
3058d359177SBarry Smith      SNESComputeJacobianDefaultColor()) and only needs to be used by someone computing a matrix via coloring directly by
306ab637aeaSJed Brown      calling MatFDColoringApply()
3077850c7c0SBarry Smith 
3087850c7c0SBarry Smith    Fortran Notes:
3097850c7c0SBarry Smith     In Fortran you must call MatFDColoringSetFunction() for a coloring object to
310ab637aeaSJed Brown   be used without SNES or within the SNES solvers.
311f881d145SBarry Smith 
312b4fc646aSLois Curfman McInnes .keywords: Mat, Jacobian, finite differences, set, function
31395b89fc3SBarry Smith 
31495b89fc3SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringGetFunction(), MatFDColoringSetFromOptions()
31595b89fc3SBarry Smith 
316005c665bSBarry Smith @*/
3177087cfbeSBarry Smith PetscErrorCode  MatFDColoringSetFunction(MatFDColoring matfd,PetscErrorCode (*f)(void),void *fctx)
318005c665bSBarry Smith {
3193a40ed3dSBarry Smith   PetscFunctionBegin;
3200700a824SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
321005c665bSBarry Smith   matfd->f    = f;
322005c665bSBarry Smith   matfd->fctx = fctx;
3233a40ed3dSBarry Smith   PetscFunctionReturn(0);
324005c665bSBarry Smith }
325005c665bSBarry Smith 
3264a2ae208SSatish Balay #undef __FUNCT__
3274a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFromOptions"
328639f9d9dSBarry Smith /*@
329b4fc646aSLois Curfman McInnes    MatFDColoringSetFromOptions - Sets coloring finite difference parameters from
330639f9d9dSBarry Smith    the options database.
331639f9d9dSBarry Smith 
332fee21e36SBarry Smith    Collective on MatFDColoring
333fee21e36SBarry Smith 
33465f2ba5bSLois Curfman McInnes    The Jacobian, F'(u), is estimated with the differencing approximation
335ef5ee4d1SLois Curfman McInnes .vb
33665f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
337f23b5b22SLois Curfman McInnes        h = error_rel*u[i]                 if  abs(u[i]) > umin
338f23b5b22SLois Curfman McInnes          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
339ef5ee4d1SLois Curfman McInnes        dx_{i} = (0, ... 1, .... 0)
340ef5ee4d1SLois Curfman McInnes .ve
341ef5ee4d1SLois Curfman McInnes 
342ef5ee4d1SLois Curfman McInnes    Input Parameter:
343ef5ee4d1SLois Curfman McInnes .  coloring - the coloring context
344ef5ee4d1SLois Curfman McInnes 
345b4fc646aSLois Curfman McInnes    Options Database Keys:
346d15ffeeaSSatish Balay +  -mat_fd_coloring_err <err> - Sets <err> (square root
347ef5ee4d1SLois Curfman McInnes            of relative error in the function)
348f23b5b22SLois Curfman McInnes .  -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude
3493ec795f1SBarry Smith .  -mat_fd_type - "wp" or "ds" (see MATMFFD_WP or MATMFFD_DS)
350ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_view - Activates basic viewing
351e350db43SBarry Smith .  -mat_fd_coloring_view ::ascii_info - Activates viewing info
352e350db43SBarry Smith -  -mat_fd_coloring_view draw - Activates drawing
353639f9d9dSBarry Smith 
35415091d37SBarry Smith     Level: intermediate
35515091d37SBarry Smith 
356b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, parameters
357d1c39f55SBarry Smith 
358d1c39f55SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringView(), MatFDColoringSetParameters()
359d1c39f55SBarry Smith 
360639f9d9dSBarry Smith @*/
3617087cfbeSBarry Smith PetscErrorCode  MatFDColoringSetFromOptions(MatFDColoring matfd)
362639f9d9dSBarry Smith {
363dfbe8321SBarry Smith   PetscErrorCode ierr;
364ace3abfcSBarry Smith   PetscBool      flg;
365efb30889SBarry Smith   char           value[3];
3663a40ed3dSBarry Smith 
3673a40ed3dSBarry Smith   PetscFunctionBegin;
3680700a824SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
369639f9d9dSBarry Smith 
3703194b578SJed Brown   ierr = PetscObjectOptionsBegin((PetscObject)matfd);CHKERRQ(ierr);
37187828ca2SBarry Smith   ierr = PetscOptionsReal("-mat_fd_coloring_err","Square root of relative error in function","MatFDColoringSetParameters",matfd->error_rel,&matfd->error_rel,0);CHKERRQ(ierr);
37287828ca2SBarry Smith   ierr = PetscOptionsReal("-mat_fd_coloring_umin","Minimum allowable u magnitude","MatFDColoringSetParameters",matfd->umin,&matfd->umin,0);CHKERRQ(ierr);
3731bc75a8dSBarry Smith   ierr = PetscOptionsString("-mat_fd_type","Algorithm to compute h, wp or ds","MatFDColoringCreate",matfd->htype,value,3,&flg);CHKERRQ(ierr);
374efb30889SBarry Smith   if (flg) {
375efb30889SBarry Smith     if (value[0] == 'w' && value[1] == 'p') matfd->htype = "wp";
376efb30889SBarry Smith     else if (value[0] == 'd' && value[1] == 's') matfd->htype = "ds";
377e32f2f54SBarry Smith     else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Unknown finite differencing type %s",value);
378efb30889SBarry Smith   }
379f86b9fbaSHong Zhang   ierr = PetscOptionsInt("-mat_fd_coloring_brows","Number of block rows","MatFDColoringSetBlockSize",matfd->brows,&matfd->brows,NULL);CHKERRQ(ierr);
38093dfae19SHong Zhang   ierr = PetscOptionsInt("-mat_fd_coloring_bcols","Number of block columns","MatFDColoringSetBlockSize",matfd->bcols,&matfd->bcols,&flg);CHKERRQ(ierr);
38193dfae19SHong Zhang   if (flg && matfd->bcols > matfd->ncolors) {
38293dfae19SHong Zhang     /* input bcols cannot be > matfd->ncolors, thus set it as ncolors */
38393dfae19SHong Zhang     matfd->bcols = matfd->ncolors;
38493dfae19SHong Zhang   }
385f86b9fbaSHong Zhang 
3865d973c19SBarry Smith   /* process any options handlers added with PetscObjectAddOptionsHandler() */
3875d973c19SBarry Smith   ierr = PetscObjectProcessOptionsHandlers((PetscObject)matfd);CHKERRQ(ierr);
388b0a32e0cSBarry Smith   PetscOptionsEnd();CHKERRQ(ierr);
3893a40ed3dSBarry Smith   PetscFunctionReturn(0);
390005c665bSBarry Smith }
391005c665bSBarry Smith 
3924a2ae208SSatish Balay #undef __FUNCT__
393e350db43SBarry Smith #define __FUNCT__ "MatFDColoringViewFromOptions"
394146574abSBarry Smith PetscErrorCode MatFDColoringViewFromOptions(MatFDColoring fd,const char prefix[],const char optionname[])
395005c665bSBarry Smith {
396dfbe8321SBarry Smith   PetscErrorCode    ierr;
397e350db43SBarry Smith   PetscBool         flg;
3983050cee2SBarry Smith   PetscViewer       viewer;
399cffb1e40SBarry Smith   PetscViewerFormat format;
400005c665bSBarry Smith 
4013a40ed3dSBarry Smith   PetscFunctionBegin;
402146574abSBarry Smith   if (prefix) {
403146574abSBarry Smith     ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)fd),prefix,optionname,&viewer,&format,&flg);CHKERRQ(ierr);
404146574abSBarry Smith   } else {
405ce94432eSBarry Smith     ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)fd),((PetscObject)fd)->prefix,optionname,&viewer,&format,&flg);CHKERRQ(ierr);
406146574abSBarry Smith   }
407005c665bSBarry Smith   if (flg) {
408cffb1e40SBarry Smith     ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr);
4093050cee2SBarry Smith     ierr = MatFDColoringView(fd,viewer);CHKERRQ(ierr);
410cffb1e40SBarry Smith     ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
411cffb1e40SBarry Smith     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
412005c665bSBarry Smith   }
4133a40ed3dSBarry Smith   PetscFunctionReturn(0);
414bbf0e169SBarry Smith }
415bbf0e169SBarry Smith 
4164a2ae208SSatish Balay #undef __FUNCT__
4174a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringCreate"
41805869f15SSatish Balay /*@
419639f9d9dSBarry Smith    MatFDColoringCreate - Creates a matrix coloring context for finite difference
420639f9d9dSBarry Smith    computation of Jacobians.
421bbf0e169SBarry Smith 
422ef5ee4d1SLois Curfman McInnes    Collective on Mat
423ef5ee4d1SLois Curfman McInnes 
424639f9d9dSBarry Smith    Input Parameters:
425ef5ee4d1SLois Curfman McInnes +  mat - the matrix containing the nonzero structure of the Jacobian
426e727c939SJed Brown -  iscoloring - the coloring of the matrix; usually obtained with MatGetColoring() or DMCreateColoring()
427bbf0e169SBarry Smith 
428bbf0e169SBarry Smith     Output Parameter:
429639f9d9dSBarry Smith .   color - the new coloring context
430bbf0e169SBarry Smith 
43115091d37SBarry Smith     Level: intermediate
43215091d37SBarry Smith 
4338d359177SBarry Smith .seealso: MatFDColoringDestroy(),SNESComputeJacobianDefaultColor(), ISColoringCreate(),
434d1c39f55SBarry Smith           MatFDColoringSetFunction(), MatFDColoringSetFromOptions(), MatFDColoringApply(),
435e727c939SJed Brown           MatFDColoringView(), MatFDColoringSetParameters(), MatGetColoring(), DMCreateColoring()
436bbf0e169SBarry Smith @*/
4377087cfbeSBarry Smith PetscErrorCode  MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
438bbf0e169SBarry Smith {
439639f9d9dSBarry Smith   MatFDColoring  c;
440639f9d9dSBarry Smith   MPI_Comm       comm;
441dfbe8321SBarry Smith   PetscErrorCode ierr;
44238baddfdSBarry Smith   PetscInt       M,N;
443639f9d9dSBarry Smith 
4443a40ed3dSBarry Smith   PetscFunctionBegin;
445*c8a9c622SHong Zhang   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
446f141eedfSHong Zhang   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Matrix must be assembled by calls to MatAssemblyBegin/End();");
447d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
448639f9d9dSBarry Smith   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
449ce94432eSBarry Smith   if (M != N) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Only for square matrices");
450f881d145SBarry Smith   ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
45167c2884eSBarry Smith   ierr = PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_CLASSID,"MatFDColoring","Jacobian computation via finite differences with coloring","Mat",comm,MatFDColoringDestroy,MatFDColoringView);CHKERRQ(ierr);
452639f9d9dSBarry Smith 
453b8f8c88eSHong Zhang   c->ctype = iscoloring->ctype;
454b8f8c88eSHong Zhang 
45593dfae19SHong Zhang   if (mat->ops->fdcoloringcreate) {
45693dfae19SHong Zhang     ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);CHKERRQ(ierr);
45793dfae19SHong Zhang   } else SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Code not yet written for matrix type %s",((PetscObject)mat)->type_name);
45893dfae19SHong Zhang 
459f77a5242SKarl Rupp   ierr = MatGetVecs(mat,NULL,&c->w1);CHKERRQ(ierr);
4603bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)c,(PetscObject)c->w1);CHKERRQ(ierr);
461b8f8c88eSHong Zhang   ierr = VecDuplicate(c->w1,&c->w2);CHKERRQ(ierr);
4623bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)c,(PetscObject)c->w2);CHKERRQ(ierr);
463b8f8c88eSHong Zhang 
46477d8c4bbSBarry Smith   c->error_rel    = PETSC_SQRT_MACHINE_EPSILON;
46577d8c4bbSBarry Smith   c->umin         = 100.0*PETSC_SQRT_MACHINE_EPSILON;
46649b058dcSBarry Smith   c->currentcolor = -1;
467efb30889SBarry Smith   c->htype        = "wp";
4684e269d77SPeter Brune   c->fset         = PETSC_FALSE;
469*c8a9c622SHong Zhang   c->setupcalled  = PETSC_FALSE;
470639f9d9dSBarry Smith 
471639f9d9dSBarry Smith   *color = c;
4724e269d77SPeter Brune   ierr   = PetscObjectCompose((PetscObject)mat,"SNESMatFDColoring",(PetscObject)c);CHKERRQ(ierr);
473d5ba7fb7SMatthew Knepley   ierr   = PetscLogEventEnd(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
4743a40ed3dSBarry Smith   PetscFunctionReturn(0);
475bbf0e169SBarry Smith }
476bbf0e169SBarry Smith 
4774a2ae208SSatish Balay #undef __FUNCT__
4784a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringDestroy"
47905869f15SSatish Balay /*@
480639f9d9dSBarry Smith     MatFDColoringDestroy - Destroys a matrix coloring context that was created
481639f9d9dSBarry Smith     via MatFDColoringCreate().
482bbf0e169SBarry Smith 
483ef5ee4d1SLois Curfman McInnes     Collective on MatFDColoring
484ef5ee4d1SLois Curfman McInnes 
485b4fc646aSLois Curfman McInnes     Input Parameter:
486639f9d9dSBarry Smith .   c - coloring context
487bbf0e169SBarry Smith 
48815091d37SBarry Smith     Level: intermediate
48915091d37SBarry Smith 
490639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
491bbf0e169SBarry Smith @*/
4926bf464f9SBarry Smith PetscErrorCode  MatFDColoringDestroy(MatFDColoring *c)
493bbf0e169SBarry Smith {
4946849ba73SBarry Smith   PetscErrorCode ierr;
49538baddfdSBarry Smith   PetscInt       i;
496bbf0e169SBarry Smith 
4973a40ed3dSBarry Smith   PetscFunctionBegin;
4986bf464f9SBarry Smith   if (!*c) PetscFunctionReturn(0);
4996bf464f9SBarry Smith   if (--((PetscObject)(*c))->refct > 0) {*c = 0; PetscFunctionReturn(0);}
500d4bb536fSBarry Smith 
5016bf464f9SBarry Smith   for (i=0; i<(*c)->ncolors; i++) {
5026bf464f9SBarry Smith     ierr = PetscFree((*c)->columns[i]);CHKERRQ(ierr);
503bbf0e169SBarry Smith   }
5046bf464f9SBarry Smith   ierr = PetscFree((*c)->ncolumns);CHKERRQ(ierr);
5056bf464f9SBarry Smith   ierr = PetscFree((*c)->columns);CHKERRQ(ierr);
5066bf464f9SBarry Smith   ierr = PetscFree((*c)->nrows);CHKERRQ(ierr);
507573f477fSHong Zhang   ierr = PetscFree((*c)->matentry);CHKERRQ(ierr);
508b7799381SHong Zhang   ierr = PetscFree((*c)->dy);CHKERRQ(ierr);
5099e917edbSHong Zhang   if ((*c)->vscale) {ierr = VecDestroy(&(*c)->vscale);CHKERRQ(ierr);}
5106bf464f9SBarry Smith   ierr = VecDestroy(&(*c)->w1);CHKERRQ(ierr);
5116bf464f9SBarry Smith   ierr = VecDestroy(&(*c)->w2);CHKERRQ(ierr);
5126bf464f9SBarry Smith   ierr = VecDestroy(&(*c)->w3);CHKERRQ(ierr);
513d38fa0fbSBarry Smith   ierr = PetscHeaderDestroy(c);CHKERRQ(ierr);
5143a40ed3dSBarry Smith   PetscFunctionReturn(0);
515bbf0e169SBarry Smith }
51643a90d84SBarry Smith 
5174a2ae208SSatish Balay #undef __FUNCT__
51849b058dcSBarry Smith #define __FUNCT__ "MatFDColoringGetPerturbedColumns"
51949b058dcSBarry Smith /*@C
52049b058dcSBarry Smith     MatFDColoringGetPerturbedColumns - Returns the indices of the columns that
52149b058dcSBarry Smith       that are currently being perturbed.
52249b058dcSBarry Smith 
52349b058dcSBarry Smith     Not Collective
52449b058dcSBarry Smith 
52549b058dcSBarry Smith     Input Parameters:
52649b058dcSBarry Smith .   coloring - coloring context created with MatFDColoringCreate()
52749b058dcSBarry Smith 
52849b058dcSBarry Smith     Output Parameters:
52949b058dcSBarry Smith +   n - the number of local columns being perturbed
53049b058dcSBarry Smith -   cols - the column indices, in global numbering
53149b058dcSBarry Smith 
53249b058dcSBarry Smith    Level: intermediate
53349b058dcSBarry Smith 
53449b058dcSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringApply()
53549b058dcSBarry Smith 
53649b058dcSBarry Smith .keywords: coloring, Jacobian, finite differences
53749b058dcSBarry Smith @*/
5387087cfbeSBarry Smith PetscErrorCode  MatFDColoringGetPerturbedColumns(MatFDColoring coloring,PetscInt *n,PetscInt *cols[])
53949b058dcSBarry Smith {
54049b058dcSBarry Smith   PetscFunctionBegin;
54149b058dcSBarry Smith   if (coloring->currentcolor >= 0) {
54249b058dcSBarry Smith     *n    = coloring->ncolumns[coloring->currentcolor];
54349b058dcSBarry Smith     *cols = coloring->columns[coloring->currentcolor];
54449b058dcSBarry Smith   } else {
54549b058dcSBarry Smith     *n = 0;
54649b058dcSBarry Smith   }
54749b058dcSBarry Smith   PetscFunctionReturn(0);
54849b058dcSBarry Smith }
54949b058dcSBarry Smith 
55049b058dcSBarry Smith #undef __FUNCT__
5514a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringApply"
55243a90d84SBarry Smith /*@
553e0907662SLois Curfman McInnes     MatFDColoringApply - Given a matrix for which a MatFDColoring context
554e0907662SLois Curfman McInnes     has been created, computes the Jacobian for a function via finite differences.
55543a90d84SBarry Smith 
556fee21e36SBarry Smith     Collective on MatFDColoring
557fee21e36SBarry Smith 
558ef5ee4d1SLois Curfman McInnes     Input Parameters:
559ef5ee4d1SLois Curfman McInnes +   mat - location to store Jacobian
560ef5ee4d1SLois Curfman McInnes .   coloring - coloring context created with MatFDColoringCreate()
561ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
5627850c7c0SBarry Smith -   sctx - context required by function, if this is being used with the SNES solver then it is SNES object, otherwise it is null
563ef5ee4d1SLois Curfman McInnes 
5648bba8e72SBarry Smith     Options Database Keys:
565ebd3b9afSBarry Smith +    -mat_fd_type - "wp" or "ds"  (see MATMFFD_WP or MATMFFD_DS)
566b9722508SLois Curfman McInnes .    -mat_fd_coloring_view - Activates basic viewing or coloring
567e350db43SBarry Smith .    -mat_fd_coloring_view draw - Activates drawing of coloring
568e350db43SBarry Smith -    -mat_fd_coloring_view ::ascii_info - Activates viewing of coloring info
5698bba8e72SBarry Smith 
57015091d37SBarry Smith     Level: intermediate
57198d79826SSatish Balay 
5727850c7c0SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringSetFunction()
57343a90d84SBarry Smith 
57443a90d84SBarry Smith .keywords: coloring, Jacobian, finite differences
57543a90d84SBarry Smith @*/
5767087cfbeSBarry Smith PetscErrorCode  MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
57743a90d84SBarry Smith {
5783acb8795SBarry Smith   PetscErrorCode ierr;
579684f2004SHong Zhang   PetscBool      flg = PETSC_FALSE;
5803acb8795SBarry Smith 
5813acb8795SBarry Smith   PetscFunctionBegin;
5820700a824SBarry Smith   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
5830700a824SBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_CLASSID,2);
5840700a824SBarry Smith   PetscValidHeaderSpecific(x1,VEC_CLASSID,3);
585ce94432eSBarry Smith   if (!coloring->f) SETERRQ(PetscObjectComm((PetscObject)J),PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetFunction()");
586e32f2f54SBarry Smith   if (!J->ops->fdcoloringapply) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not supported for this matrix type %s",((PetscObject)J)->type_name);
587*c8a9c622SHong Zhang   if (!coloring->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetUp()");
588684f2004SHong Zhang 
589684f2004SHong Zhang   ierr = MatSetUnfactored(J);CHKERRQ(ierr);
590684f2004SHong Zhang   ierr = PetscOptionsGetBool(NULL,"-mat_fd_coloring_dont_rezero",&flg,NULL);CHKERRQ(ierr);
591684f2004SHong Zhang   if (flg) {
592684f2004SHong Zhang     ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr);
593684f2004SHong Zhang   } else {
594684f2004SHong Zhang     PetscBool assembled;
595684f2004SHong Zhang     ierr = MatAssembled(J,&assembled);CHKERRQ(ierr);
596684f2004SHong Zhang     if (assembled) {
597684f2004SHong Zhang       ierr = MatZeroEntries(J);CHKERRQ(ierr);
598684f2004SHong Zhang     }
599684f2004SHong Zhang   }
600684f2004SHong Zhang 
6015922145eSHong Zhang   ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
6023acb8795SBarry Smith   ierr = (*J->ops->fdcoloringapply)(J,coloring,x1,flag,sctx);CHKERRQ(ierr);
6035922145eSHong Zhang   ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
6043acb8795SBarry Smith   PetscFunctionReturn(0);
6053acb8795SBarry Smith }
606