xref: /petsc/src/mat/matfd/fdmatrix.c (revision c3339decea92175325d9368fa13196bcd0e0e58b)
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 
10d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringSetF(MatFDColoring fd, Vec F)
11d71ae5a4SJacob Faibussowitsch {
123a7fca6bSBarry Smith   PetscFunctionBegin;
134e269d77SPeter Brune   if (F) {
149566063dSJacob Faibussowitsch     PetscCall(VecCopy(F, fd->w1));
154e269d77SPeter Brune     fd->fset = PETSC_TRUE;
164e269d77SPeter Brune   } else {
174e269d77SPeter Brune     fd->fset = PETSC_FALSE;
184e269d77SPeter Brune   }
193a7fca6bSBarry Smith   PetscFunctionReturn(0);
203a7fca6bSBarry Smith }
213a7fca6bSBarry Smith 
229804daf3SBarry Smith #include <petscdraw.h>
23d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatFDColoringView_Draw_Zoom(PetscDraw draw, void *Aa)
24d71ae5a4SJacob Faibussowitsch {
259194eea9SBarry Smith   MatFDColoring fd = (MatFDColoring)Aa;
26b312708bSHong Zhang   PetscInt      i, j, nz, row;
279194eea9SBarry Smith   PetscReal     x, y;
28b312708bSHong Zhang   MatEntry     *Jentry = fd->matentry;
299194eea9SBarry Smith 
309194eea9SBarry Smith   PetscFunctionBegin;
319194eea9SBarry Smith   /* loop over colors  */
32b312708bSHong Zhang   nz = 0;
339194eea9SBarry Smith   for (i = 0; i < fd->ncolors; i++) {
349194eea9SBarry Smith     for (j = 0; j < fd->nrows[i]; j++) {
35b312708bSHong Zhang       row = Jentry[nz].row;
36b312708bSHong Zhang       y   = fd->M - row - fd->rstart;
37b312708bSHong Zhang       x   = (PetscReal)Jentry[nz++].col;
389566063dSJacob Faibussowitsch       PetscCall(PetscDrawRectangle(draw, x, y, x + 1, y + 1, i + 1, i + 1, i + 1, i + 1));
399194eea9SBarry Smith     }
409194eea9SBarry Smith   }
419194eea9SBarry Smith   PetscFunctionReturn(0);
429194eea9SBarry Smith }
439194eea9SBarry Smith 
44d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatFDColoringView_Draw(MatFDColoring fd, PetscViewer viewer)
45d71ae5a4SJacob Faibussowitsch {
46ace3abfcSBarry Smith   PetscBool isnull;
47b0a32e0cSBarry Smith   PetscDraw draw;
4854d96fa3SBarry Smith   PetscReal xr, yr, xl, yl, h, w;
49005c665bSBarry Smith 
503a40ed3dSBarry Smith   PetscFunctionBegin;
519566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
529566063dSJacob Faibussowitsch   PetscCall(PetscDrawIsNull(draw, &isnull));
53832b7cebSLisandro Dalcin   if (isnull) PetscFunctionReturn(0);
54005c665bSBarry Smith 
559371c9d4SSatish Balay   xr = fd->N;
569371c9d4SSatish Balay   yr = fd->M;
579371c9d4SSatish Balay   h  = yr / 10.0;
589371c9d4SSatish Balay   w  = xr / 10.0;
599371c9d4SSatish Balay   xr += w;
609371c9d4SSatish Balay   yr += h;
619371c9d4SSatish Balay   xl = -w;
629371c9d4SSatish Balay   yl = -h;
639566063dSJacob Faibussowitsch   PetscCall(PetscDrawSetCoordinates(draw, xl, yl, xr, yr));
649566063dSJacob Faibussowitsch   PetscCall(PetscObjectCompose((PetscObject)fd, "Zoomviewer", (PetscObject)viewer));
659566063dSJacob Faibussowitsch   PetscCall(PetscDrawZoom(draw, MatFDColoringView_Draw_Zoom, fd));
669566063dSJacob Faibussowitsch   PetscCall(PetscObjectCompose((PetscObject)fd, "Zoomviewer", NULL));
679566063dSJacob Faibussowitsch   PetscCall(PetscDrawSave(draw));
683a40ed3dSBarry Smith   PetscFunctionReturn(0);
69005c665bSBarry Smith }
70005c665bSBarry Smith 
71bbf0e169SBarry Smith /*@C
72639f9d9dSBarry Smith    MatFDColoringView - Views a finite difference coloring context.
73bbf0e169SBarry Smith 
7411a5261eSBarry Smith    Collective on c
75fee21e36SBarry Smith 
76ef5ee4d1SLois Curfman McInnes    Input  Parameters:
77ef5ee4d1SLois Curfman McInnes +  c - the coloring context
78ef5ee4d1SLois Curfman McInnes -  viewer - visualization context
79ef5ee4d1SLois Curfman McInnes 
8015091d37SBarry Smith    Level: intermediate
8115091d37SBarry Smith 
82b4fc646aSLois Curfman McInnes    Notes:
83b4fc646aSLois Curfman McInnes    The available visualization contexts include
8411a5261eSBarry Smith +     `PETSC_VIEWER_STDOUT_SELF` - standard output (default)
8511a5261eSBarry Smith .     `PETSC_VIEWER_STDOUT_WORLD` - synchronized standard
86ef5ee4d1SLois Curfman McInnes         output where only the first processor opens
87ef5ee4d1SLois Curfman McInnes         the file.  All other processors send their
88ef5ee4d1SLois Curfman McInnes         data to the first processor to print.
8911a5261eSBarry Smith -     `PETSC_VIEWER_DRAW_WORLD` - graphical display of nonzero structure
90bbf0e169SBarry Smith 
919194eea9SBarry Smith      Since PETSc uses only a small number of basic colors (currently 33), if the coloring
92b0a32e0cSBarry Smith    involves more than 33 then some seemingly identical colors are displayed making it look
939194eea9SBarry Smith    like an illegal coloring. This is just a graphical artifact.
949194eea9SBarry Smith 
9511a5261eSBarry Smith .seealso: `MatFDColoring`, `MatFDColoringCreate()`
96bbf0e169SBarry Smith @*/
97d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringView(MatFDColoring c, PetscViewer viewer)
98d71ae5a4SJacob Faibussowitsch {
9938baddfdSBarry Smith   PetscInt          i, j;
100ace3abfcSBarry Smith   PetscBool         isdraw, iascii;
101f3ef73ceSBarry Smith   PetscViewerFormat format;
102bbf0e169SBarry Smith 
1033a40ed3dSBarry Smith   PetscFunctionBegin;
1040700a824SBarry Smith   PetscValidHeaderSpecific(c, MAT_FDCOLORING_CLASSID, 1);
10548a46eb9SPierre Jolivet   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)c), &viewer));
1060700a824SBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
107c9780b6fSBarry Smith   PetscCheckSameComm(c, 1, viewer, 2);
108bbf0e169SBarry Smith 
1099566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
1109566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
1110f5bd95cSBarry Smith   if (isdraw) {
1129566063dSJacob Faibussowitsch     PetscCall(MatFDColoringView_Draw(c, viewer));
11332077d6dSBarry Smith   } else if (iascii) {
1149566063dSJacob Faibussowitsch     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)c, viewer));
1159566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Error tolerance=%g\n", (double)c->error_rel));
1169566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Umin=%g\n", (double)c->umin));
1179566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Number of colors=%" PetscInt_FMT "\n", c->ncolors));
118ae09f205SBarry Smith 
1199566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetFormat(viewer, &format));
120fb9695e5SSatish Balay     if (format != PETSC_VIEWER_ASCII_INFO) {
121b312708bSHong Zhang       PetscInt row, col, nz;
122b312708bSHong Zhang       nz = 0;
123b4fc646aSLois Curfman McInnes       for (i = 0; i < c->ncolors; i++) {
1249566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  Information for color %" PetscInt_FMT "\n", i));
1259566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "    Number of columns %" PetscInt_FMT "\n", c->ncolumns[i]));
12648a46eb9SPierre Jolivet         for (j = 0; j < c->ncolumns[i]; j++) PetscCall(PetscViewerASCIIPrintf(viewer, "      %" PetscInt_FMT "\n", c->columns[i][j]));
1279566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "    Number of rows %" PetscInt_FMT "\n", c->nrows[i]));
1285bdb020cSBarry Smith         if (c->matentry) {
129b4fc646aSLois Curfman McInnes           for (j = 0; j < c->nrows[i]; j++) {
130b312708bSHong Zhang             row = c->matentry[nz].row;
131b312708bSHong Zhang             col = c->matentry[nz++].col;
1329566063dSJacob Faibussowitsch             PetscCall(PetscViewerASCIIPrintf(viewer, "      %" PetscInt_FMT " %" PetscInt_FMT " \n", row, col));
133b4fc646aSLois Curfman McInnes           }
134bbf0e169SBarry Smith         }
135bbf0e169SBarry Smith       }
1365bdb020cSBarry Smith     }
1379566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(viewer));
138bbf0e169SBarry Smith   }
1393a40ed3dSBarry Smith   PetscFunctionReturn(0);
140639f9d9dSBarry Smith }
141639f9d9dSBarry Smith 
142639f9d9dSBarry Smith /*@
143b4fc646aSLois Curfman McInnes    MatFDColoringSetParameters - Sets the parameters for the sparse approximation of
144b4fc646aSLois Curfman McInnes    a Jacobian matrix using finite differences.
145639f9d9dSBarry Smith 
146*c3339decSBarry Smith    Logically Collective
147ef5ee4d1SLois Curfman McInnes 
148ef5ee4d1SLois Curfman McInnes    The Jacobian is estimated with the differencing approximation
149ef5ee4d1SLois Curfman McInnes .vb
15065f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
151d461c19dSHong Zhang        htype = 'ds':
152f23b5b22SLois Curfman McInnes          h = error_rel*u[i]                 if  abs(u[i]) > umin
153f23b5b22SLois Curfman McInnes            = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
154ef5ee4d1SLois Curfman McInnes          dx_{i} = (0, ... 1, .... 0)
155d461c19dSHong Zhang 
156d461c19dSHong Zhang        htype = 'wp':
157d461c19dSHong Zhang          h = error_rel * sqrt(1 + ||u||)
158ef5ee4d1SLois Curfman McInnes .ve
159639f9d9dSBarry Smith 
160639f9d9dSBarry Smith    Input Parameters:
161ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
162639f9d9dSBarry Smith .  error_rel - relative error
163f23b5b22SLois Curfman McInnes -  umin - minimum allowable u-value magnitude
164fee21e36SBarry Smith 
16515091d37SBarry Smith    Level: advanced
16615091d37SBarry Smith 
16711a5261eSBarry Smith .seealso: `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringSetFromOptions()`
168639f9d9dSBarry Smith @*/
169d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringSetParameters(MatFDColoring matfd, PetscReal error, PetscReal umin)
170d71ae5a4SJacob Faibussowitsch {
1713a40ed3dSBarry Smith   PetscFunctionBegin;
1720700a824SBarry Smith   PetscValidHeaderSpecific(matfd, MAT_FDCOLORING_CLASSID, 1);
173c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(matfd, error, 2);
174c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(matfd, umin, 3);
175639f9d9dSBarry Smith   if (error != PETSC_DEFAULT) matfd->error_rel = error;
176639f9d9dSBarry Smith   if (umin != PETSC_DEFAULT) matfd->umin = umin;
1773a40ed3dSBarry Smith   PetscFunctionReturn(0);
178639f9d9dSBarry Smith }
179639f9d9dSBarry Smith 
180f86b9fbaSHong Zhang /*@
181c8a9c622SHong Zhang    MatFDColoringSetBlockSize - Sets block size for efficient inserting entries of Jacobian matrix.
182005c665bSBarry Smith 
183*c3339decSBarry Smith    Logically Collective
184f86b9fbaSHong Zhang 
185f86b9fbaSHong Zhang    Input Parameters:
186f86b9fbaSHong Zhang +  coloring - the coloring context
187f86b9fbaSHong Zhang .  brows - number of rows in the block
188f86b9fbaSHong Zhang -  bcols - number of columns in the block
189f86b9fbaSHong Zhang 
190f86b9fbaSHong Zhang    Level: intermediate
191f86b9fbaSHong Zhang 
19211a5261eSBarry Smith .seealso: `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringSetFromOptions()`
193f86b9fbaSHong Zhang @*/
194d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring matfd, PetscInt brows, PetscInt bcols)
195d71ae5a4SJacob Faibussowitsch {
196f86b9fbaSHong Zhang   PetscFunctionBegin;
197f86b9fbaSHong Zhang   PetscValidHeaderSpecific(matfd, MAT_FDCOLORING_CLASSID, 1);
198f86b9fbaSHong Zhang   PetscValidLogicalCollectiveInt(matfd, brows, 2);
199f86b9fbaSHong Zhang   PetscValidLogicalCollectiveInt(matfd, bcols, 3);
200f86b9fbaSHong Zhang   if (brows != PETSC_DEFAULT) matfd->brows = brows;
201f86b9fbaSHong Zhang   if (bcols != PETSC_DEFAULT) matfd->bcols = bcols;
202f86b9fbaSHong Zhang   PetscFunctionReturn(0);
203f86b9fbaSHong Zhang }
204f86b9fbaSHong Zhang 
205f86b9fbaSHong Zhang /*@
206f86b9fbaSHong Zhang    MatFDColoringSetUp - Sets up the internal data structures of matrix coloring context for the later use.
207f86b9fbaSHong Zhang 
208*c3339decSBarry Smith    Collective
209f86b9fbaSHong Zhang 
210f86b9fbaSHong Zhang    Input Parameters:
211f86b9fbaSHong Zhang +  mat - the matrix containing the nonzero structure of the Jacobian
21211a5261eSBarry Smith .  iscoloring - the coloring of the matrix; usually obtained with `MatGetColoring()` or `DMCreateColoring()`
213f86b9fbaSHong Zhang -  color - the matrix coloring context
214f86b9fbaSHong Zhang 
215f86b9fbaSHong Zhang    Level: beginner
216f86b9fbaSHong Zhang 
21711a5261eSBarry Smith    Notes:
21811a5261eSBarry Smith    When the coloring type is `IS_COLORING_LOCAL` the coloring is in the local ordering of the unknowns.
219bdaf1daeSBarry Smith 
22011a5261eSBarry Smith .seealso: `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringDestroy()`
221f86b9fbaSHong Zhang @*/
222d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringSetUp(Mat mat, ISColoring iscoloring, MatFDColoring color)
223d71ae5a4SJacob Faibussowitsch {
224bdaf1daeSBarry Smith   PetscBool eq;
225f86b9fbaSHong Zhang 
226f86b9fbaSHong Zhang   PetscFunctionBegin;
227c8a9c622SHong Zhang   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
228c8a9c622SHong Zhang   PetscValidHeaderSpecific(color, MAT_FDCOLORING_CLASSID, 3);
229c8a9c622SHong Zhang   if (color->setupcalled) PetscFunctionReturn(0);
2309566063dSJacob Faibussowitsch   PetscCall(PetscObjectCompareId((PetscObject)mat, color->matid, &eq));
23128b400f6SJacob Faibussowitsch   PetscCheck(eq, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONG, "Matrix used with MatFDColoringSetUp() must be that used with MatFDColoringCreate()");
232c8a9c622SHong Zhang 
2339566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(MAT_FDColoringSetUp, mat, 0, 0, 0));
234dbbe0bcdSBarry Smith   PetscUseTypeMethod(mat, fdcoloringsetup, iscoloring, color);
235c8a9c622SHong Zhang 
236c8a9c622SHong Zhang   color->setupcalled = PETSC_TRUE;
2379566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(MAT_FDColoringSetUp, mat, 0, 0, 0));
238f86b9fbaSHong Zhang   PetscFunctionReturn(0);
239f86b9fbaSHong Zhang }
240ff0cfa39SBarry Smith 
2414a9d489dSBarry Smith /*@C
2424a9d489dSBarry Smith    MatFDColoringGetFunction - Gets the function to use for computing the Jacobian.
2434a9d489dSBarry Smith 
2443f9fe445SBarry Smith    Not Collective
2454a9d489dSBarry Smith 
246f899ff85SJose E. Roman    Input Parameter:
2474a9d489dSBarry Smith .  coloring - the coloring context
2484a9d489dSBarry Smith 
2494a9d489dSBarry Smith    Output Parameters:
2504a9d489dSBarry Smith +  f - the function
2514a9d489dSBarry Smith -  fctx - the optional user-defined function context
2524a9d489dSBarry Smith 
2534a9d489dSBarry Smith    Level: intermediate
2544a9d489dSBarry Smith 
25511a5261eSBarry Smith .seealso: `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringSetFunction()`, `MatFDColoringSetFromOptions()`
2564a9d489dSBarry Smith @*/
257d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringGetFunction(MatFDColoring matfd, PetscErrorCode (**f)(void), void **fctx)
258d71ae5a4SJacob Faibussowitsch {
2594a9d489dSBarry Smith   PetscFunctionBegin;
2600700a824SBarry Smith   PetscValidHeaderSpecific(matfd, MAT_FDCOLORING_CLASSID, 1);
2614a9d489dSBarry Smith   if (f) *f = matfd->f;
2624a9d489dSBarry Smith   if (fctx) *fctx = matfd->fctx;
2634a9d489dSBarry Smith   PetscFunctionReturn(0);
2644a9d489dSBarry Smith }
2654a9d489dSBarry Smith 
266d64ed03dSBarry Smith /*@C
267005c665bSBarry Smith    MatFDColoringSetFunction - Sets the function to use for computing the Jacobian.
268005c665bSBarry Smith 
269*c3339decSBarry Smith    Logically Collective
270fee21e36SBarry Smith 
271ef5ee4d1SLois Curfman McInnes    Input Parameters:
272ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
273ef5ee4d1SLois Curfman McInnes .  f - the function
274ef5ee4d1SLois Curfman McInnes -  fctx - the optional user-defined function context
275ef5ee4d1SLois Curfman McInnes 
2767850c7c0SBarry Smith    Calling sequence of (*f) function:
27711a5261eSBarry Smith     For `SNES`:    PetscErrorCode (*f)(SNES,Vec,Vec,void*)
27811a5261eSBarry Smith     If not using `SNES`: PetscErrorCode (*f)(void *dummy,Vec,Vec,void*) and dummy is ignored
27915091d37SBarry Smith 
2807850c7c0SBarry Smith    Level: advanced
2817850c7c0SBarry Smith 
28211a5261eSBarry Smith    Note:
28311a5261eSBarry Smith     This function is usually used automatically by `SNES` (when one uses `SNESSetJacobian()` with the argument
28411a5261eSBarry Smith      `SNESComputeJacobianDefaultColor()`) and only needs to be used by someone computing a matrix via coloring directly by
28511a5261eSBarry Smith      calling `MatFDColoringApply()`
2867850c7c0SBarry Smith 
28711a5261eSBarry Smith    Fortran Note:
28811a5261eSBarry Smith     In Fortran you must call `MatFDColoringSetFunction()` for a coloring object to
28911a5261eSBarry Smith   be used without `SNES` or within the `SNES` solvers.
290f881d145SBarry Smith 
29111a5261eSBarry Smith .seealso: `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringGetFunction()`, `MatFDColoringSetFromOptions()`
292005c665bSBarry Smith @*/
293d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringSetFunction(MatFDColoring matfd, PetscErrorCode (*f)(void), void *fctx)
294d71ae5a4SJacob Faibussowitsch {
2953a40ed3dSBarry Smith   PetscFunctionBegin;
2960700a824SBarry Smith   PetscValidHeaderSpecific(matfd, MAT_FDCOLORING_CLASSID, 1);
297005c665bSBarry Smith   matfd->f    = f;
298005c665bSBarry Smith   matfd->fctx = fctx;
2993a40ed3dSBarry Smith   PetscFunctionReturn(0);
300005c665bSBarry Smith }
301005c665bSBarry Smith 
302639f9d9dSBarry Smith /*@
303b4fc646aSLois Curfman McInnes    MatFDColoringSetFromOptions - Sets coloring finite difference parameters from
304639f9d9dSBarry Smith    the options database.
305639f9d9dSBarry Smith 
306*c3339decSBarry Smith    Collective
307fee21e36SBarry Smith 
30865f2ba5bSLois Curfman McInnes    The Jacobian, F'(u), is estimated with the differencing approximation
309ef5ee4d1SLois Curfman McInnes .vb
31065f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
311f23b5b22SLois Curfman McInnes        h = error_rel*u[i]                 if  abs(u[i]) > umin
312f23b5b22SLois Curfman McInnes          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
313ef5ee4d1SLois Curfman McInnes        dx_{i} = (0, ... 1, .... 0)
314ef5ee4d1SLois Curfman McInnes .ve
315ef5ee4d1SLois Curfman McInnes 
316ef5ee4d1SLois Curfman McInnes    Input Parameter:
317ef5ee4d1SLois Curfman McInnes .  coloring - the coloring context
318ef5ee4d1SLois Curfman McInnes 
319b4fc646aSLois Curfman McInnes    Options Database Keys:
3205620d6dcSBarry Smith +  -mat_fd_coloring_err <err> - Sets <err> (square root of relative error in the function)
321f23b5b22SLois Curfman McInnes .  -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude
3223ec795f1SBarry Smith .  -mat_fd_type - "wp" or "ds" (see MATMFFD_WP or MATMFFD_DS)
323ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_view - Activates basic viewing
324e350db43SBarry Smith .  -mat_fd_coloring_view ::ascii_info - Activates viewing info
325e350db43SBarry Smith -  -mat_fd_coloring_view draw - Activates drawing
326639f9d9dSBarry Smith 
32715091d37SBarry Smith     Level: intermediate
32815091d37SBarry Smith 
32911a5261eSBarry Smith .seealso: `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringView()`, `MatFDColoringSetParameters()`
330639f9d9dSBarry Smith @*/
331d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring matfd)
332d71ae5a4SJacob Faibussowitsch {
333ace3abfcSBarry Smith   PetscBool flg;
334efb30889SBarry Smith   char      value[3];
3353a40ed3dSBarry Smith 
3363a40ed3dSBarry Smith   PetscFunctionBegin;
3370700a824SBarry Smith   PetscValidHeaderSpecific(matfd, MAT_FDCOLORING_CLASSID, 1);
338639f9d9dSBarry Smith 
339d0609cedSBarry Smith   PetscObjectOptionsBegin((PetscObject)matfd);
3409566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_fd_coloring_err", "Square root of relative error in function", "MatFDColoringSetParameters", matfd->error_rel, &matfd->error_rel, NULL));
3419566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_fd_coloring_umin", "Minimum allowable u magnitude", "MatFDColoringSetParameters", matfd->umin, &matfd->umin, NULL));
3429566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-mat_fd_type", "Algorithm to compute h, wp or ds", "MatFDColoringCreate", matfd->htype, value, sizeof(value), &flg));
343efb30889SBarry Smith   if (flg) {
344efb30889SBarry Smith     if (value[0] == 'w' && value[1] == 'p') matfd->htype = "wp";
345efb30889SBarry Smith     else if (value[0] == 'd' && value[1] == 's') matfd->htype = "ds";
34698921bdaSJacob Faibussowitsch     else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unknown finite differencing type %s", value);
347efb30889SBarry Smith   }
3489566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-mat_fd_coloring_brows", "Number of block rows", "MatFDColoringSetBlockSize", matfd->brows, &matfd->brows, NULL));
3499566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-mat_fd_coloring_bcols", "Number of block columns", "MatFDColoringSetBlockSize", matfd->bcols, &matfd->bcols, &flg));
35093dfae19SHong Zhang   if (flg && matfd->bcols > matfd->ncolors) {
35193dfae19SHong Zhang     /* input bcols cannot be > matfd->ncolors, thus set it as ncolors */
35293dfae19SHong Zhang     matfd->bcols = matfd->ncolors;
35393dfae19SHong Zhang   }
354f86b9fbaSHong Zhang 
3555d973c19SBarry Smith   /* process any options handlers added with PetscObjectAddOptionsHandler() */
356dbbe0bcdSBarry Smith   PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)matfd, PetscOptionsObject));
357d0609cedSBarry Smith   PetscOptionsEnd();
3583a40ed3dSBarry Smith   PetscFunctionReturn(0);
359005c665bSBarry Smith }
360005c665bSBarry Smith 
36161ab5df0SBarry Smith /*@C
36261ab5df0SBarry Smith    MatFDColoringSetType - Sets the approach for computing the finite difference parameter
36361ab5df0SBarry Smith 
364*c3339decSBarry Smith    Collective
36561ab5df0SBarry Smith 
36661ab5df0SBarry Smith    Input Parameters:
36761ab5df0SBarry Smith +  coloring - the coloring context
36811a5261eSBarry Smith -  type - either `MATMFFD_WP` or `MATMFFD_DS`
36961ab5df0SBarry Smith 
37061ab5df0SBarry Smith    Options Database Keys:
37161ab5df0SBarry Smith .  -mat_fd_type - "wp" or "ds"
37261ab5df0SBarry Smith 
37311a5261eSBarry Smith    Note:
37411a5261eSBarry Smith    It is goofy that the argument type is `MatMFFDType` since the `MatFDColoring` actually computes the matrix entries
37511a5261eSBarry Smith          but the process of computing the entries is the same as as with the `MATMFFD` operation so we should reuse the names instead of
37661ab5df0SBarry Smith          introducing another one.
37761ab5df0SBarry Smith 
37861ab5df0SBarry Smith    Level: intermediate
37961ab5df0SBarry Smith 
38011a5261eSBarry Smith .seealso: `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringView()`, `MatFDColoringSetParameters()`
38161ab5df0SBarry Smith @*/
382d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringSetType(MatFDColoring matfd, MatMFFDType type)
383d71ae5a4SJacob Faibussowitsch {
38461ab5df0SBarry Smith   PetscFunctionBegin;
38561ab5df0SBarry Smith   PetscValidHeaderSpecific(matfd, MAT_FDCOLORING_CLASSID, 1);
38661ab5df0SBarry Smith   /*
38761ab5df0SBarry Smith      It is goofy to handle the strings this way but currently there is no code to free a dynamically created matfd->htype
38861ab5df0SBarry Smith      and this function is being provided as patch for a release so it shouldn't change the implementaton
38961ab5df0SBarry Smith   */
39061ab5df0SBarry Smith   if (type[0] == 'w' && type[1] == 'p') matfd->htype = "wp";
39161ab5df0SBarry Smith   else if (type[0] == 'd' && type[1] == 's') matfd->htype = "ds";
39298921bdaSJacob Faibussowitsch   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unknown finite differencing type %s", type);
39361ab5df0SBarry Smith   PetscFunctionReturn(0);
39461ab5df0SBarry Smith }
39561ab5df0SBarry Smith 
396d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringViewFromOptions(MatFDColoring fd, const char prefix[], const char optionname[])
397d71ae5a4SJacob Faibussowitsch {
398e350db43SBarry Smith   PetscBool         flg;
3993050cee2SBarry Smith   PetscViewer       viewer;
400cffb1e40SBarry Smith   PetscViewerFormat format;
401005c665bSBarry Smith 
4023a40ed3dSBarry Smith   PetscFunctionBegin;
403146574abSBarry Smith   if (prefix) {
4049566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)fd), ((PetscObject)fd)->options, prefix, optionname, &viewer, &format, &flg));
405146574abSBarry Smith   } else {
4069566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)fd), ((PetscObject)fd)->options, ((PetscObject)fd)->prefix, optionname, &viewer, &format, &flg));
407146574abSBarry Smith   }
408005c665bSBarry Smith   if (flg) {
4099566063dSJacob Faibussowitsch     PetscCall(PetscViewerPushFormat(viewer, format));
4109566063dSJacob Faibussowitsch     PetscCall(MatFDColoringView(fd, viewer));
4119566063dSJacob Faibussowitsch     PetscCall(PetscViewerPopFormat(viewer));
4129566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
413005c665bSBarry Smith   }
4143a40ed3dSBarry Smith   PetscFunctionReturn(0);
415bbf0e169SBarry Smith }
416bbf0e169SBarry Smith 
41705869f15SSatish Balay /*@
418639f9d9dSBarry Smith    MatFDColoringCreate - Creates a matrix coloring context for finite difference
419639f9d9dSBarry Smith    computation of Jacobians.
420bbf0e169SBarry Smith 
421*c3339decSBarry Smith    Collective
422ef5ee4d1SLois Curfman McInnes 
423639f9d9dSBarry Smith    Input Parameters:
424ef5ee4d1SLois Curfman McInnes +  mat - the matrix containing the nonzero structure of the Jacobian
42511a5261eSBarry Smith -  iscoloring - the coloring of the matrix; usually obtained with `MatColoringCreate()` or `DMCreateColoring()`
426bbf0e169SBarry Smith 
427bbf0e169SBarry Smith     Output Parameter:
428639f9d9dSBarry Smith .   color - the new coloring context
429bbf0e169SBarry Smith 
43015091d37SBarry Smith     Level: intermediate
43115091d37SBarry Smith 
43211a5261eSBarry Smith .seealso: `MatFDColoring`, `MatFDColoringDestroy()`, `SNESComputeJacobianDefaultColor()`, `ISColoringCreate()`,
433db781477SPatrick Sanan           `MatFDColoringSetFunction()`, `MatFDColoringSetFromOptions()`, `MatFDColoringApply()`,
434db781477SPatrick Sanan           `MatFDColoringView()`, `MatFDColoringSetParameters()`, `MatColoringCreate()`, `DMCreateColoring()`, `MatFDColoringSetValues()`
435bbf0e169SBarry Smith @*/
436d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringCreate(Mat mat, ISColoring iscoloring, MatFDColoring *color)
437d71ae5a4SJacob Faibussowitsch {
438639f9d9dSBarry Smith   MatFDColoring c;
439639f9d9dSBarry Smith   MPI_Comm      comm;
44038baddfdSBarry Smith   PetscInt      M, N;
441639f9d9dSBarry Smith 
4423a40ed3dSBarry Smith   PetscFunctionBegin;
443c8a9c622SHong Zhang   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
44428b400f6SJacob Faibussowitsch   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Matrix must be assembled by calls to MatAssemblyBegin/End();");
4459566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(MAT_FDColoringCreate, mat, 0, 0, 0));
4469566063dSJacob Faibussowitsch   PetscCall(MatGetSize(mat, &M, &N));
44708401ef6SPierre Jolivet   PetscCheck(M == N, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Only for square matrices");
4489566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)mat, &comm));
4499566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(c, MAT_FDCOLORING_CLASSID, "MatFDColoring", "Jacobian computation via finite differences with coloring", "Mat", comm, MatFDColoringDestroy, MatFDColoringView));
450639f9d9dSBarry Smith 
451b8f8c88eSHong Zhang   c->ctype = iscoloring->ctype;
4529566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetId((PetscObject)mat, &c->matid));
453b8f8c88eSHong Zhang 
454dbbe0bcdSBarry Smith   PetscUseTypeMethod(mat, fdcoloringcreate, iscoloring, c);
45593dfae19SHong Zhang 
4569566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(mat, NULL, &c->w1));
457ce911d59SRichard Tran Mills   /* Vec is used intensively in particular piece of scalar CPU code; won't benefit from bouncing back and forth to the GPU */
4589566063dSJacob Faibussowitsch   PetscCall(VecBindToCPU(c->w1, PETSC_TRUE));
4599566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(c->w1, &c->w2));
460ce911d59SRichard Tran Mills   /* Vec is used intensively in particular piece of scalar CPU code; won't benefit from bouncing back and forth to the GPU */
4619566063dSJacob Faibussowitsch   PetscCall(VecBindToCPU(c->w2, PETSC_TRUE));
462b8f8c88eSHong Zhang 
46377d8c4bbSBarry Smith   c->error_rel    = PETSC_SQRT_MACHINE_EPSILON;
46477d8c4bbSBarry Smith   c->umin         = 100.0 * PETSC_SQRT_MACHINE_EPSILON;
46549b058dcSBarry Smith   c->currentcolor = -1;
466efb30889SBarry Smith   c->htype        = "wp";
4674e269d77SPeter Brune   c->fset         = PETSC_FALSE;
468c8a9c622SHong Zhang   c->setupcalled  = PETSC_FALSE;
469639f9d9dSBarry Smith 
470639f9d9dSBarry Smith   *color = c;
4719566063dSJacob Faibussowitsch   PetscCall(PetscObjectCompose((PetscObject)mat, "SNESMatFDColoring", (PetscObject)c));
4729566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(MAT_FDColoringCreate, mat, 0, 0, 0));
4733a40ed3dSBarry Smith   PetscFunctionReturn(0);
474bbf0e169SBarry Smith }
475bbf0e169SBarry Smith 
47605869f15SSatish Balay /*@
477639f9d9dSBarry Smith     MatFDColoringDestroy - Destroys a matrix coloring context that was created
47811a5261eSBarry Smith     via `MatFDColoringCreate()`.
479bbf0e169SBarry Smith 
480*c3339decSBarry Smith     Collective
481ef5ee4d1SLois Curfman McInnes 
482b4fc646aSLois Curfman McInnes     Input Parameter:
483639f9d9dSBarry Smith .   c - coloring context
484bbf0e169SBarry Smith 
48515091d37SBarry Smith     Level: intermediate
48615091d37SBarry Smith 
48711a5261eSBarry Smith .seealso: `MatFDColoring`, `MatFDColoringCreate()`
488bbf0e169SBarry Smith @*/
489d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringDestroy(MatFDColoring *c)
490d71ae5a4SJacob Faibussowitsch {
49138baddfdSBarry Smith   PetscInt      i;
4920df34763SHong Zhang   MatFDColoring color = *c;
493bbf0e169SBarry Smith 
4943a40ed3dSBarry Smith   PetscFunctionBegin;
4956bf464f9SBarry Smith   if (!*c) PetscFunctionReturn(0);
4969371c9d4SSatish Balay   if (--((PetscObject)color)->refct > 0) {
4979371c9d4SSatish Balay     *c = NULL;
4989371c9d4SSatish Balay     PetscFunctionReturn(0);
4999371c9d4SSatish Balay   }
500d4bb536fSBarry Smith 
501071fcb05SBarry Smith   /* we do not free the column arrays since their entries are owned by the ISs in color->isa */
50248a46eb9SPierre Jolivet   for (i = 0; i < color->ncolors; i++) PetscCall(ISDestroy(&color->isa[i]));
5039566063dSJacob Faibussowitsch   PetscCall(PetscFree(color->isa));
5049566063dSJacob Faibussowitsch   PetscCall(PetscFree2(color->ncolumns, color->columns));
5059566063dSJacob Faibussowitsch   PetscCall(PetscFree(color->nrows));
5060df34763SHong Zhang   if (color->htype[0] == 'w') {
5079566063dSJacob Faibussowitsch     PetscCall(PetscFree(color->matentry2));
5080df34763SHong Zhang   } else {
5099566063dSJacob Faibussowitsch     PetscCall(PetscFree(color->matentry));
5100df34763SHong Zhang   }
5119566063dSJacob Faibussowitsch   PetscCall(PetscFree(color->dy));
5129566063dSJacob Faibussowitsch   if (color->vscale) PetscCall(VecDestroy(&color->vscale));
5139566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&color->w1));
5149566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&color->w2));
5159566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&color->w3));
5169566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(c));
5173a40ed3dSBarry Smith   PetscFunctionReturn(0);
518bbf0e169SBarry Smith }
51943a90d84SBarry Smith 
52049b058dcSBarry Smith /*@C
52149b058dcSBarry Smith     MatFDColoringGetPerturbedColumns - Returns the indices of the columns that
52249b058dcSBarry Smith       that are currently being perturbed.
52349b058dcSBarry Smith 
52449b058dcSBarry Smith     Not Collective
52549b058dcSBarry Smith 
526f899ff85SJose E. Roman     Input Parameter:
52711a5261eSBarry Smith .   coloring - coloring context created with `MatFDColoringCreate()`
52849b058dcSBarry Smith 
52949b058dcSBarry Smith     Output Parameters:
53049b058dcSBarry Smith +   n - the number of local columns being perturbed
53149b058dcSBarry Smith -   cols - the column indices, in global numbering
53249b058dcSBarry Smith 
53321fcc2ddSBarry Smith    Level: advanced
53449b058dcSBarry Smith 
53511a5261eSBarry Smith    Note:
53611a5261eSBarry Smith    IF the matrix type is `MATBAIJ`, then the block column indices are returned
537e66d0d93SMatthew Knepley 
538edaa7c33SBarry Smith    Fortran Note:
539edaa7c33SBarry Smith    This routine has a different interface for Fortran
54021fcc2ddSBarry Smith $     #include <petsc/finclude/petscmat.h>
54121fcc2ddSBarry Smith $          use petscmat
542edaa7c33SBarry Smith $          PetscInt, pointer :: array(:)
543edaa7c33SBarry Smith $          PetscErrorCode  ierr
544edaa7c33SBarry Smith $          MatFDColoring   i
545edaa7c33SBarry Smith $          call MatFDColoringGetPerturbedColumnsF90(i,array,ierr)
546edaa7c33SBarry Smith $      use the entries of array ...
547edaa7c33SBarry Smith $          call MatFDColoringRestorePerturbedColumnsF90(i,array,ierr)
548edaa7c33SBarry Smith 
54911a5261eSBarry Smith .seealso: `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringDestroy()`, `MatFDColoringView()`, `MatFDColoringApply()`
55049b058dcSBarry Smith @*/
551d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring coloring, PetscInt *n, const PetscInt *cols[])
552d71ae5a4SJacob Faibussowitsch {
55349b058dcSBarry Smith   PetscFunctionBegin;
55449b058dcSBarry Smith   if (coloring->currentcolor >= 0) {
55549b058dcSBarry Smith     *n    = coloring->ncolumns[coloring->currentcolor];
55649b058dcSBarry Smith     *cols = coloring->columns[coloring->currentcolor];
55749b058dcSBarry Smith   } else {
55849b058dcSBarry Smith     *n = 0;
55949b058dcSBarry Smith   }
56049b058dcSBarry Smith   PetscFunctionReturn(0);
56149b058dcSBarry Smith }
56249b058dcSBarry Smith 
56343a90d84SBarry Smith /*@
56411a5261eSBarry Smith     MatFDColoringApply - Given a matrix for which a `MatFDColoring` context
565e0907662SLois Curfman McInnes     has been created, computes the Jacobian for a function via finite differences.
56643a90d84SBarry Smith 
567*c3339decSBarry Smith     Collective
568fee21e36SBarry Smith 
569ef5ee4d1SLois Curfman McInnes     Input Parameters:
570ef5ee4d1SLois Curfman McInnes +   mat - location to store Jacobian
57111a5261eSBarry Smith .   coloring - coloring context created with `MatFDColoringCreate()`
572ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
57311a5261eSBarry Smith -   sctx - context required by function, if this is being used with the SNES solver then it is `SNES` object, otherwise it is null
574ef5ee4d1SLois Curfman McInnes 
5758bba8e72SBarry Smith     Options Database Keys:
57611a5261eSBarry Smith +    -mat_fd_type - "wp" or "ds"  (see `MATMFFD_WP` or `MATMFFD_DS`)
577b9722508SLois Curfman McInnes .    -mat_fd_coloring_view - Activates basic viewing or coloring
578e350db43SBarry Smith .    -mat_fd_coloring_view draw - Activates drawing of coloring
579e350db43SBarry Smith -    -mat_fd_coloring_view ::ascii_info - Activates viewing of coloring info
5808bba8e72SBarry Smith 
58115091d37SBarry Smith     Level: intermediate
58298d79826SSatish Balay 
58311a5261eSBarry Smith .seealso: `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringDestroy()`, `MatFDColoringView()`, `MatFDColoringSetFunction()`, `MatFDColoringSetValues()`
58443a90d84SBarry Smith @*/
585d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringApply(Mat J, MatFDColoring coloring, Vec x1, void *sctx)
586d71ae5a4SJacob Faibussowitsch {
587bdaf1daeSBarry Smith   PetscBool eq;
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);
5939566063dSJacob Faibussowitsch   PetscCall(PetscObjectCompareId((PetscObject)J, coloring->matid, &eq));
59428b400f6SJacob Faibussowitsch   PetscCheck(eq, PetscObjectComm((PetscObject)J), PETSC_ERR_ARG_WRONG, "Matrix used with MatFDColoringApply() must be that used with MatFDColoringCreate()");
59528b400f6SJacob Faibussowitsch   PetscCheck(coloring->f, PetscObjectComm((PetscObject)J), PETSC_ERR_ARG_WRONGSTATE, "Must call MatFDColoringSetFunction()");
59628b400f6SJacob Faibussowitsch   PetscCheck(coloring->setupcalled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call MatFDColoringSetUp()");
597684f2004SHong Zhang 
5989566063dSJacob Faibussowitsch   PetscCall(MatSetUnfactored(J));
5999566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(MAT_FDColoringApply, coloring, J, x1, 0));
600dbbe0bcdSBarry Smith   PetscUseTypeMethod(J, fdcoloringapply, coloring, x1, sctx);
6019566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(MAT_FDColoringApply, coloring, J, x1, 0));
6025bdb020cSBarry Smith   if (!coloring->viewed) {
6039566063dSJacob Faibussowitsch     PetscCall(MatFDColoringViewFromOptions(coloring, NULL, "-mat_fd_coloring_view"));
6045bdb020cSBarry Smith     coloring->viewed = PETSC_TRUE;
6055bdb020cSBarry Smith   }
6063acb8795SBarry Smith   PetscFunctionReturn(0);
6073acb8795SBarry Smith }
608