xref: /petsc/src/mat/matfd/fdmatrix.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
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 
109371c9d4SSatish Balay PetscErrorCode MatFDColoringSetF(MatFDColoring fd, Vec F) {
113a7fca6bSBarry Smith   PetscFunctionBegin;
124e269d77SPeter Brune   if (F) {
139566063dSJacob Faibussowitsch     PetscCall(VecCopy(F, fd->w1));
144e269d77SPeter Brune     fd->fset = PETSC_TRUE;
154e269d77SPeter Brune   } else {
164e269d77SPeter Brune     fd->fset = PETSC_FALSE;
174e269d77SPeter Brune   }
183a7fca6bSBarry Smith   PetscFunctionReturn(0);
193a7fca6bSBarry Smith }
203a7fca6bSBarry Smith 
219804daf3SBarry Smith #include <petscdraw.h>
229371c9d4SSatish Balay static PetscErrorCode MatFDColoringView_Draw_Zoom(PetscDraw draw, void *Aa) {
239194eea9SBarry Smith   MatFDColoring fd = (MatFDColoring)Aa;
24b312708bSHong Zhang   PetscInt      i, j, nz, row;
259194eea9SBarry Smith   PetscReal     x, y;
26b312708bSHong Zhang   MatEntry     *Jentry = fd->matentry;
279194eea9SBarry Smith 
289194eea9SBarry Smith   PetscFunctionBegin;
299194eea9SBarry Smith   /* loop over colors  */
30b312708bSHong Zhang   nz = 0;
319194eea9SBarry Smith   for (i = 0; i < fd->ncolors; i++) {
329194eea9SBarry Smith     for (j = 0; j < fd->nrows[i]; j++) {
33b312708bSHong Zhang       row = Jentry[nz].row;
34b312708bSHong Zhang       y   = fd->M - row - fd->rstart;
35b312708bSHong Zhang       x   = (PetscReal)Jentry[nz++].col;
369566063dSJacob Faibussowitsch       PetscCall(PetscDrawRectangle(draw, x, y, x + 1, y + 1, i + 1, i + 1, i + 1, i + 1));
379194eea9SBarry Smith     }
389194eea9SBarry Smith   }
399194eea9SBarry Smith   PetscFunctionReturn(0);
409194eea9SBarry Smith }
419194eea9SBarry Smith 
429371c9d4SSatish Balay static PetscErrorCode MatFDColoringView_Draw(MatFDColoring fd, PetscViewer viewer) {
43ace3abfcSBarry Smith   PetscBool isnull;
44b0a32e0cSBarry Smith   PetscDraw draw;
4554d96fa3SBarry Smith   PetscReal xr, yr, xl, yl, h, w;
46005c665bSBarry Smith 
473a40ed3dSBarry Smith   PetscFunctionBegin;
489566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
499566063dSJacob Faibussowitsch   PetscCall(PetscDrawIsNull(draw, &isnull));
50832b7cebSLisandro Dalcin   if (isnull) PetscFunctionReturn(0);
51005c665bSBarry Smith 
529371c9d4SSatish Balay   xr = fd->N;
539371c9d4SSatish Balay   yr = fd->M;
549371c9d4SSatish Balay   h  = yr / 10.0;
559371c9d4SSatish Balay   w  = xr / 10.0;
569371c9d4SSatish Balay   xr += w;
579371c9d4SSatish Balay   yr += h;
589371c9d4SSatish Balay   xl = -w;
599371c9d4SSatish Balay   yl = -h;
609566063dSJacob Faibussowitsch   PetscCall(PetscDrawSetCoordinates(draw, xl, yl, xr, yr));
619566063dSJacob Faibussowitsch   PetscCall(PetscObjectCompose((PetscObject)fd, "Zoomviewer", (PetscObject)viewer));
629566063dSJacob Faibussowitsch   PetscCall(PetscDrawZoom(draw, MatFDColoringView_Draw_Zoom, fd));
639566063dSJacob Faibussowitsch   PetscCall(PetscObjectCompose((PetscObject)fd, "Zoomviewer", NULL));
649566063dSJacob Faibussowitsch   PetscCall(PetscDrawSave(draw));
653a40ed3dSBarry Smith   PetscFunctionReturn(0);
66005c665bSBarry Smith }
67005c665bSBarry Smith 
68bbf0e169SBarry Smith /*@C
69639f9d9dSBarry Smith    MatFDColoringView - Views a finite difference coloring context.
70bbf0e169SBarry Smith 
714c49b128SBarry Smith    Collective on MatFDColoring
72fee21e36SBarry Smith 
73ef5ee4d1SLois Curfman McInnes    Input  Parameters:
74ef5ee4d1SLois Curfman McInnes +  c - the coloring context
75ef5ee4d1SLois Curfman McInnes -  viewer - visualization context
76ef5ee4d1SLois Curfman McInnes 
7715091d37SBarry Smith    Level: intermediate
7815091d37SBarry Smith 
79b4fc646aSLois Curfman McInnes    Notes:
80b4fc646aSLois Curfman McInnes    The available visualization contexts include
81b0a32e0cSBarry Smith +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
82b0a32e0cSBarry Smith .     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
83ef5ee4d1SLois Curfman McInnes         output where only the first processor opens
84ef5ee4d1SLois Curfman McInnes         the file.  All other processors send their
85ef5ee4d1SLois Curfman McInnes         data to the first processor to print.
86b0a32e0cSBarry Smith -     PETSC_VIEWER_DRAW_WORLD - graphical display of nonzero structure
87bbf0e169SBarry Smith 
889194eea9SBarry Smith    Notes:
899194eea9SBarry Smith      Since PETSc uses only a small number of basic colors (currently 33), if the coloring
90b0a32e0cSBarry Smith    involves more than 33 then some seemingly identical colors are displayed making it look
919194eea9SBarry Smith    like an illegal coloring. This is just a graphical artifact.
929194eea9SBarry Smith 
93db781477SPatrick Sanan .seealso: `MatFDColoringCreate()`
94005c665bSBarry Smith 
95bbf0e169SBarry Smith @*/
969371c9d4SSatish Balay PetscErrorCode MatFDColoringView(MatFDColoring c, PetscViewer viewer) {
9738baddfdSBarry Smith   PetscInt          i, j;
98ace3abfcSBarry Smith   PetscBool         isdraw, iascii;
99f3ef73ceSBarry Smith   PetscViewerFormat format;
100bbf0e169SBarry Smith 
1013a40ed3dSBarry Smith   PetscFunctionBegin;
1020700a824SBarry Smith   PetscValidHeaderSpecific(c, MAT_FDCOLORING_CLASSID, 1);
103*48a46eb9SPierre Jolivet   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)c), &viewer));
1040700a824SBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
105c9780b6fSBarry Smith   PetscCheckSameComm(c, 1, viewer, 2);
106bbf0e169SBarry Smith 
1079566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
1089566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
1090f5bd95cSBarry Smith   if (isdraw) {
1109566063dSJacob Faibussowitsch     PetscCall(MatFDColoringView_Draw(c, viewer));
11132077d6dSBarry Smith   } else if (iascii) {
1129566063dSJacob Faibussowitsch     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)c, viewer));
1139566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Error tolerance=%g\n", (double)c->error_rel));
1149566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Umin=%g\n", (double)c->umin));
1159566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Number of colors=%" PetscInt_FMT "\n", c->ncolors));
116ae09f205SBarry Smith 
1179566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetFormat(viewer, &format));
118fb9695e5SSatish Balay     if (format != PETSC_VIEWER_ASCII_INFO) {
119b312708bSHong Zhang       PetscInt row, col, nz;
120b312708bSHong Zhang       nz = 0;
121b4fc646aSLois Curfman McInnes       for (i = 0; i < c->ncolors; i++) {
1229566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  Information for color %" PetscInt_FMT "\n", i));
1239566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "    Number of columns %" PetscInt_FMT "\n", c->ncolumns[i]));
124*48a46eb9SPierre Jolivet         for (j = 0; j < c->ncolumns[i]; j++) PetscCall(PetscViewerASCIIPrintf(viewer, "      %" PetscInt_FMT "\n", c->columns[i][j]));
1259566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "    Number of rows %" PetscInt_FMT "\n", c->nrows[i]));
1265bdb020cSBarry Smith         if (c->matentry) {
127b4fc646aSLois Curfman McInnes           for (j = 0; j < c->nrows[i]; j++) {
128b312708bSHong Zhang             row = c->matentry[nz].row;
129b312708bSHong Zhang             col = c->matentry[nz++].col;
1309566063dSJacob Faibussowitsch             PetscCall(PetscViewerASCIIPrintf(viewer, "      %" PetscInt_FMT " %" PetscInt_FMT " \n", row, col));
131b4fc646aSLois Curfman McInnes           }
132bbf0e169SBarry Smith         }
133bbf0e169SBarry Smith       }
1345bdb020cSBarry Smith     }
1359566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(viewer));
136bbf0e169SBarry Smith   }
1373a40ed3dSBarry Smith   PetscFunctionReturn(0);
138639f9d9dSBarry Smith }
139639f9d9dSBarry Smith 
140639f9d9dSBarry Smith /*@
141b4fc646aSLois Curfman McInnes    MatFDColoringSetParameters - Sets the parameters for the sparse approximation of
142b4fc646aSLois Curfman McInnes    a Jacobian matrix using finite differences.
143639f9d9dSBarry Smith 
1443f9fe445SBarry Smith    Logically Collective on MatFDColoring
145ef5ee4d1SLois Curfman McInnes 
146ef5ee4d1SLois Curfman McInnes    The Jacobian is estimated with the differencing approximation
147ef5ee4d1SLois Curfman McInnes .vb
14865f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
149d461c19dSHong Zhang        htype = 'ds':
150f23b5b22SLois Curfman McInnes          h = error_rel*u[i]                 if  abs(u[i]) > umin
151f23b5b22SLois Curfman McInnes            = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
152ef5ee4d1SLois Curfman McInnes          dx_{i} = (0, ... 1, .... 0)
153d461c19dSHong Zhang 
154d461c19dSHong Zhang        htype = 'wp':
155d461c19dSHong Zhang          h = error_rel * sqrt(1 + ||u||)
156ef5ee4d1SLois Curfman McInnes .ve
157639f9d9dSBarry Smith 
158639f9d9dSBarry Smith    Input Parameters:
159ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
160639f9d9dSBarry Smith .  error_rel - relative error
161f23b5b22SLois Curfman McInnes -  umin - minimum allowable u-value magnitude
162fee21e36SBarry Smith 
16315091d37SBarry Smith    Level: advanced
16415091d37SBarry Smith 
165db781477SPatrick Sanan .seealso: `MatFDColoringCreate()`, `MatFDColoringSetFromOptions()`
16695b89fc3SBarry Smith 
167639f9d9dSBarry Smith @*/
1689371c9d4SSatish Balay PetscErrorCode MatFDColoringSetParameters(MatFDColoring matfd, PetscReal error, PetscReal umin) {
1693a40ed3dSBarry Smith   PetscFunctionBegin;
1700700a824SBarry Smith   PetscValidHeaderSpecific(matfd, MAT_FDCOLORING_CLASSID, 1);
171c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(matfd, error, 2);
172c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(matfd, umin, 3);
173639f9d9dSBarry Smith   if (error != PETSC_DEFAULT) matfd->error_rel = error;
174639f9d9dSBarry Smith   if (umin != PETSC_DEFAULT) matfd->umin = umin;
1753a40ed3dSBarry Smith   PetscFunctionReturn(0);
176639f9d9dSBarry Smith }
177639f9d9dSBarry Smith 
178f86b9fbaSHong Zhang /*@
179c8a9c622SHong Zhang    MatFDColoringSetBlockSize - Sets block size for efficient inserting entries of Jacobian matrix.
180005c665bSBarry Smith 
181f86b9fbaSHong Zhang    Logically Collective on MatFDColoring
182f86b9fbaSHong Zhang 
183f86b9fbaSHong Zhang    Input Parameters:
184f86b9fbaSHong Zhang +  coloring - the coloring context
185f86b9fbaSHong Zhang .  brows - number of rows in the block
186f86b9fbaSHong Zhang -  bcols - number of columns in the block
187f86b9fbaSHong Zhang 
188f86b9fbaSHong Zhang    Level: intermediate
189f86b9fbaSHong Zhang 
190db781477SPatrick Sanan .seealso: `MatFDColoringCreate()`, `MatFDColoringSetFromOptions()`
191f86b9fbaSHong Zhang 
192f86b9fbaSHong Zhang @*/
1939371c9d4SSatish Balay PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring matfd, PetscInt brows, PetscInt bcols) {
194f86b9fbaSHong Zhang   PetscFunctionBegin;
195f86b9fbaSHong Zhang   PetscValidHeaderSpecific(matfd, MAT_FDCOLORING_CLASSID, 1);
196f86b9fbaSHong Zhang   PetscValidLogicalCollectiveInt(matfd, brows, 2);
197f86b9fbaSHong Zhang   PetscValidLogicalCollectiveInt(matfd, bcols, 3);
198f86b9fbaSHong Zhang   if (brows != PETSC_DEFAULT) matfd->brows = brows;
199f86b9fbaSHong Zhang   if (bcols != PETSC_DEFAULT) matfd->bcols = bcols;
200f86b9fbaSHong Zhang   PetscFunctionReturn(0);
201f86b9fbaSHong Zhang }
202f86b9fbaSHong Zhang 
203f86b9fbaSHong Zhang /*@
204f86b9fbaSHong Zhang    MatFDColoringSetUp - Sets up the internal data structures of matrix coloring context for the later use.
205f86b9fbaSHong Zhang 
206f86b9fbaSHong Zhang    Collective on Mat
207f86b9fbaSHong Zhang 
208f86b9fbaSHong Zhang    Input Parameters:
209f86b9fbaSHong Zhang +  mat - the matrix containing the nonzero structure of the Jacobian
210f86b9fbaSHong Zhang .  iscoloring - the coloring of the matrix; usually obtained with MatGetColoring() or DMCreateColoring()
211f86b9fbaSHong Zhang -  color - the matrix coloring context
212f86b9fbaSHong Zhang 
213f86b9fbaSHong Zhang    Level: beginner
214f86b9fbaSHong Zhang 
215bdaf1daeSBarry Smith    Notes: When the coloring type is IS_COLORING_LOCAL the coloring is in the local ordering of the unknowns.
216bdaf1daeSBarry Smith 
217db781477SPatrick Sanan .seealso: `MatFDColoringCreate()`, `MatFDColoringDestroy()`
218f86b9fbaSHong Zhang @*/
2199371c9d4SSatish Balay PetscErrorCode MatFDColoringSetUp(Mat mat, ISColoring iscoloring, MatFDColoring color) {
220bdaf1daeSBarry Smith   PetscBool eq;
221f86b9fbaSHong Zhang 
222f86b9fbaSHong Zhang   PetscFunctionBegin;
223c8a9c622SHong Zhang   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
224c8a9c622SHong Zhang   PetscValidHeaderSpecific(color, MAT_FDCOLORING_CLASSID, 3);
225c8a9c622SHong Zhang   if (color->setupcalled) PetscFunctionReturn(0);
2269566063dSJacob Faibussowitsch   PetscCall(PetscObjectCompareId((PetscObject)mat, color->matid, &eq));
22728b400f6SJacob Faibussowitsch   PetscCheck(eq, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONG, "Matrix used with MatFDColoringSetUp() must be that used with MatFDColoringCreate()");
228c8a9c622SHong Zhang 
2299566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(MAT_FDColoringSetUp, mat, 0, 0, 0));
230dbbe0bcdSBarry Smith   PetscUseTypeMethod(mat, fdcoloringsetup, iscoloring, color);
231c8a9c622SHong Zhang 
232c8a9c622SHong Zhang   color->setupcalled = PETSC_TRUE;
2339566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(MAT_FDColoringSetUp, mat, 0, 0, 0));
234f86b9fbaSHong Zhang   PetscFunctionReturn(0);
235f86b9fbaSHong Zhang }
236ff0cfa39SBarry Smith 
2374a9d489dSBarry Smith /*@C
2384a9d489dSBarry Smith    MatFDColoringGetFunction - Gets the function to use for computing the Jacobian.
2394a9d489dSBarry Smith 
2403f9fe445SBarry Smith    Not Collective
2414a9d489dSBarry Smith 
242f899ff85SJose E. Roman    Input Parameter:
2434a9d489dSBarry Smith .  coloring - the coloring context
2444a9d489dSBarry Smith 
2454a9d489dSBarry Smith    Output Parameters:
2464a9d489dSBarry Smith +  f - the function
2474a9d489dSBarry Smith -  fctx - the optional user-defined function context
2484a9d489dSBarry Smith 
2494a9d489dSBarry Smith    Level: intermediate
2504a9d489dSBarry Smith 
251db781477SPatrick Sanan .seealso: `MatFDColoringCreate()`, `MatFDColoringSetFunction()`, `MatFDColoringSetFromOptions()`
25295b89fc3SBarry Smith 
2534a9d489dSBarry Smith @*/
2549371c9d4SSatish Balay PetscErrorCode MatFDColoringGetFunction(MatFDColoring matfd, PetscErrorCode (**f)(void), void **fctx) {
2554a9d489dSBarry Smith   PetscFunctionBegin;
2560700a824SBarry Smith   PetscValidHeaderSpecific(matfd, MAT_FDCOLORING_CLASSID, 1);
2574a9d489dSBarry Smith   if (f) *f = matfd->f;
2584a9d489dSBarry Smith   if (fctx) *fctx = matfd->fctx;
2594a9d489dSBarry Smith   PetscFunctionReturn(0);
2604a9d489dSBarry Smith }
2614a9d489dSBarry Smith 
262d64ed03dSBarry Smith /*@C
263005c665bSBarry Smith    MatFDColoringSetFunction - Sets the function to use for computing the Jacobian.
264005c665bSBarry Smith 
2653f9fe445SBarry Smith    Logically Collective on MatFDColoring
266fee21e36SBarry Smith 
267ef5ee4d1SLois Curfman McInnes    Input Parameters:
268ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
269ef5ee4d1SLois Curfman McInnes .  f - the function
270ef5ee4d1SLois Curfman McInnes -  fctx - the optional user-defined function context
271ef5ee4d1SLois Curfman McInnes 
2727850c7c0SBarry Smith    Calling sequence of (*f) function:
2737850c7c0SBarry Smith     For SNES:    PetscErrorCode (*f)(SNES,Vec,Vec,void*)
274ab637aeaSJed Brown     If not using SNES: PetscErrorCode (*f)(void *dummy,Vec,Vec,void*) and dummy is ignored
27515091d37SBarry Smith 
2767850c7c0SBarry Smith    Level: advanced
2777850c7c0SBarry Smith 
27895452b02SPatrick Sanan    Notes:
27995452b02SPatrick Sanan     This function is usually used automatically by SNES (when one uses SNESSetJacobian() with the argument
2808d359177SBarry Smith      SNESComputeJacobianDefaultColor()) and only needs to be used by someone computing a matrix via coloring directly by
281ab637aeaSJed Brown      calling MatFDColoringApply()
2827850c7c0SBarry Smith 
2837850c7c0SBarry Smith    Fortran Notes:
2847850c7c0SBarry Smith     In Fortran you must call MatFDColoringSetFunction() for a coloring object to
285ab637aeaSJed Brown   be used without SNES or within the SNES solvers.
286f881d145SBarry Smith 
287db781477SPatrick Sanan .seealso: `MatFDColoringCreate()`, `MatFDColoringGetFunction()`, `MatFDColoringSetFromOptions()`
28895b89fc3SBarry Smith 
289005c665bSBarry Smith @*/
2909371c9d4SSatish Balay PetscErrorCode MatFDColoringSetFunction(MatFDColoring matfd, PetscErrorCode (*f)(void), void *fctx) {
2913a40ed3dSBarry Smith   PetscFunctionBegin;
2920700a824SBarry Smith   PetscValidHeaderSpecific(matfd, MAT_FDCOLORING_CLASSID, 1);
293005c665bSBarry Smith   matfd->f    = f;
294005c665bSBarry Smith   matfd->fctx = fctx;
2953a40ed3dSBarry Smith   PetscFunctionReturn(0);
296005c665bSBarry Smith }
297005c665bSBarry Smith 
298639f9d9dSBarry Smith /*@
299b4fc646aSLois Curfman McInnes    MatFDColoringSetFromOptions - Sets coloring finite difference parameters from
300639f9d9dSBarry Smith    the options database.
301639f9d9dSBarry Smith 
302fee21e36SBarry Smith    Collective on MatFDColoring
303fee21e36SBarry Smith 
30465f2ba5bSLois Curfman McInnes    The Jacobian, F'(u), is estimated with the differencing approximation
305ef5ee4d1SLois Curfman McInnes .vb
30665f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
307f23b5b22SLois Curfman McInnes        h = error_rel*u[i]                 if  abs(u[i]) > umin
308f23b5b22SLois Curfman McInnes          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
309ef5ee4d1SLois Curfman McInnes        dx_{i} = (0, ... 1, .... 0)
310ef5ee4d1SLois Curfman McInnes .ve
311ef5ee4d1SLois Curfman McInnes 
312ef5ee4d1SLois Curfman McInnes    Input Parameter:
313ef5ee4d1SLois Curfman McInnes .  coloring - the coloring context
314ef5ee4d1SLois Curfman McInnes 
315b4fc646aSLois Curfman McInnes    Options Database Keys:
3165620d6dcSBarry Smith +  -mat_fd_coloring_err <err> - Sets <err> (square root of relative error in the function)
317f23b5b22SLois Curfman McInnes .  -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude
3183ec795f1SBarry Smith .  -mat_fd_type - "wp" or "ds" (see MATMFFD_WP or MATMFFD_DS)
319ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_view - Activates basic viewing
320e350db43SBarry Smith .  -mat_fd_coloring_view ::ascii_info - Activates viewing info
321e350db43SBarry Smith -  -mat_fd_coloring_view draw - Activates drawing
322639f9d9dSBarry Smith 
32315091d37SBarry Smith     Level: intermediate
32415091d37SBarry Smith 
325db781477SPatrick Sanan .seealso: `MatFDColoringCreate()`, `MatFDColoringView()`, `MatFDColoringSetParameters()`
326d1c39f55SBarry Smith 
327639f9d9dSBarry Smith @*/
3289371c9d4SSatish Balay PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring matfd) {
329ace3abfcSBarry Smith   PetscBool flg;
330efb30889SBarry Smith   char      value[3];
3313a40ed3dSBarry Smith 
3323a40ed3dSBarry Smith   PetscFunctionBegin;
3330700a824SBarry Smith   PetscValidHeaderSpecific(matfd, MAT_FDCOLORING_CLASSID, 1);
334639f9d9dSBarry Smith 
335d0609cedSBarry Smith   PetscObjectOptionsBegin((PetscObject)matfd);
3369566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_fd_coloring_err", "Square root of relative error in function", "MatFDColoringSetParameters", matfd->error_rel, &matfd->error_rel, NULL));
3379566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_fd_coloring_umin", "Minimum allowable u magnitude", "MatFDColoringSetParameters", matfd->umin, &matfd->umin, NULL));
3389566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-mat_fd_type", "Algorithm to compute h, wp or ds", "MatFDColoringCreate", matfd->htype, value, sizeof(value), &flg));
339efb30889SBarry Smith   if (flg) {
340efb30889SBarry Smith     if (value[0] == 'w' && value[1] == 'p') matfd->htype = "wp";
341efb30889SBarry Smith     else if (value[0] == 'd' && value[1] == 's') matfd->htype = "ds";
34298921bdaSJacob Faibussowitsch     else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unknown finite differencing type %s", value);
343efb30889SBarry Smith   }
3449566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-mat_fd_coloring_brows", "Number of block rows", "MatFDColoringSetBlockSize", matfd->brows, &matfd->brows, NULL));
3459566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-mat_fd_coloring_bcols", "Number of block columns", "MatFDColoringSetBlockSize", matfd->bcols, &matfd->bcols, &flg));
34693dfae19SHong Zhang   if (flg && matfd->bcols > matfd->ncolors) {
34793dfae19SHong Zhang     /* input bcols cannot be > matfd->ncolors, thus set it as ncolors */
34893dfae19SHong Zhang     matfd->bcols = matfd->ncolors;
34993dfae19SHong Zhang   }
350f86b9fbaSHong Zhang 
3515d973c19SBarry Smith   /* process any options handlers added with PetscObjectAddOptionsHandler() */
352dbbe0bcdSBarry Smith   PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)matfd, PetscOptionsObject));
353d0609cedSBarry Smith   PetscOptionsEnd();
3543a40ed3dSBarry Smith   PetscFunctionReturn(0);
355005c665bSBarry Smith }
356005c665bSBarry Smith 
35761ab5df0SBarry Smith /*@C
35861ab5df0SBarry Smith    MatFDColoringSetType - Sets the approach for computing the finite difference parameter
35961ab5df0SBarry Smith 
36061ab5df0SBarry Smith    Collective on MatFDColoring
36161ab5df0SBarry Smith 
36261ab5df0SBarry Smith    Input Parameters:
36361ab5df0SBarry Smith +  coloring - the coloring context
36461ab5df0SBarry Smith -  type - either MATMFFD_WP or MATMFFD_DS
36561ab5df0SBarry Smith 
36661ab5df0SBarry Smith    Options Database Keys:
36761ab5df0SBarry Smith .  -mat_fd_type - "wp" or "ds"
36861ab5df0SBarry Smith 
36961ab5df0SBarry Smith    Note: It is goofy that the argument type is MatMFFDType since the MatFDColoring actually computes the matrix entries
37061ab5df0SBarry 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
37161ab5df0SBarry Smith          introducing another one.
37261ab5df0SBarry Smith 
37361ab5df0SBarry Smith    Level: intermediate
37461ab5df0SBarry Smith 
375db781477SPatrick Sanan .seealso: `MatFDColoringCreate()`, `MatFDColoringView()`, `MatFDColoringSetParameters()`
37661ab5df0SBarry Smith 
37761ab5df0SBarry Smith @*/
3789371c9d4SSatish Balay PetscErrorCode MatFDColoringSetType(MatFDColoring matfd, MatMFFDType type) {
37961ab5df0SBarry Smith   PetscFunctionBegin;
38061ab5df0SBarry Smith   PetscValidHeaderSpecific(matfd, MAT_FDCOLORING_CLASSID, 1);
38161ab5df0SBarry Smith   /*
38261ab5df0SBarry Smith      It is goofy to handle the strings this way but currently there is no code to free a dynamically created matfd->htype
38361ab5df0SBarry Smith      and this function is being provided as patch for a release so it shouldn't change the implementaton
38461ab5df0SBarry Smith   */
38561ab5df0SBarry Smith   if (type[0] == 'w' && type[1] == 'p') matfd->htype = "wp";
38661ab5df0SBarry Smith   else if (type[0] == 'd' && type[1] == 's') matfd->htype = "ds";
38798921bdaSJacob Faibussowitsch   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unknown finite differencing type %s", type);
38861ab5df0SBarry Smith   PetscFunctionReturn(0);
38961ab5df0SBarry Smith }
39061ab5df0SBarry Smith 
3919371c9d4SSatish Balay PetscErrorCode MatFDColoringViewFromOptions(MatFDColoring fd, const char prefix[], const char optionname[]) {
392e350db43SBarry Smith   PetscBool         flg;
3933050cee2SBarry Smith   PetscViewer       viewer;
394cffb1e40SBarry Smith   PetscViewerFormat format;
395005c665bSBarry Smith 
3963a40ed3dSBarry Smith   PetscFunctionBegin;
397146574abSBarry Smith   if (prefix) {
3989566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)fd), ((PetscObject)fd)->options, prefix, optionname, &viewer, &format, &flg));
399146574abSBarry Smith   } else {
4009566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)fd), ((PetscObject)fd)->options, ((PetscObject)fd)->prefix, optionname, &viewer, &format, &flg));
401146574abSBarry Smith   }
402005c665bSBarry Smith   if (flg) {
4039566063dSJacob Faibussowitsch     PetscCall(PetscViewerPushFormat(viewer, format));
4049566063dSJacob Faibussowitsch     PetscCall(MatFDColoringView(fd, viewer));
4059566063dSJacob Faibussowitsch     PetscCall(PetscViewerPopFormat(viewer));
4069566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
407005c665bSBarry Smith   }
4083a40ed3dSBarry Smith   PetscFunctionReturn(0);
409bbf0e169SBarry Smith }
410bbf0e169SBarry Smith 
41105869f15SSatish Balay /*@
412639f9d9dSBarry Smith    MatFDColoringCreate - Creates a matrix coloring context for finite difference
413639f9d9dSBarry Smith    computation of Jacobians.
414bbf0e169SBarry Smith 
415ef5ee4d1SLois Curfman McInnes    Collective on Mat
416ef5ee4d1SLois Curfman McInnes 
417639f9d9dSBarry Smith    Input Parameters:
418ef5ee4d1SLois Curfman McInnes +  mat - the matrix containing the nonzero structure of the Jacobian
4197086a01eSPeter Brune -  iscoloring - the coloring of the matrix; usually obtained with MatColoringCreate() or DMCreateColoring()
420bbf0e169SBarry Smith 
421bbf0e169SBarry Smith     Output Parameter:
422639f9d9dSBarry Smith .   color - the new coloring context
423bbf0e169SBarry Smith 
42415091d37SBarry Smith     Level: intermediate
42515091d37SBarry Smith 
426c2e3fba1SPatrick Sanan .seealso: `MatFDColoringDestroy()`, `SNESComputeJacobianDefaultColor()`, `ISColoringCreate()`,
427db781477SPatrick Sanan           `MatFDColoringSetFunction()`, `MatFDColoringSetFromOptions()`, `MatFDColoringApply()`,
428db781477SPatrick Sanan           `MatFDColoringView()`, `MatFDColoringSetParameters()`, `MatColoringCreate()`, `DMCreateColoring()`, `MatFDColoringSetValues()`
429bbf0e169SBarry Smith @*/
4309371c9d4SSatish Balay PetscErrorCode MatFDColoringCreate(Mat mat, ISColoring iscoloring, MatFDColoring *color) {
431639f9d9dSBarry Smith   MatFDColoring c;
432639f9d9dSBarry Smith   MPI_Comm      comm;
43338baddfdSBarry Smith   PetscInt      M, N;
434639f9d9dSBarry Smith 
4353a40ed3dSBarry Smith   PetscFunctionBegin;
436c8a9c622SHong Zhang   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
43728b400f6SJacob Faibussowitsch   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Matrix must be assembled by calls to MatAssemblyBegin/End();");
4389566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(MAT_FDColoringCreate, mat, 0, 0, 0));
4399566063dSJacob Faibussowitsch   PetscCall(MatGetSize(mat, &M, &N));
44008401ef6SPierre Jolivet   PetscCheck(M == N, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Only for square matrices");
4419566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)mat, &comm));
4429566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(c, MAT_FDCOLORING_CLASSID, "MatFDColoring", "Jacobian computation via finite differences with coloring", "Mat", comm, MatFDColoringDestroy, MatFDColoringView));
443639f9d9dSBarry Smith 
444b8f8c88eSHong Zhang   c->ctype = iscoloring->ctype;
4459566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetId((PetscObject)mat, &c->matid));
446b8f8c88eSHong Zhang 
447dbbe0bcdSBarry Smith   PetscUseTypeMethod(mat, fdcoloringcreate, iscoloring, c);
44893dfae19SHong Zhang 
4499566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(mat, NULL, &c->w1));
450ce911d59SRichard Tran Mills   /* Vec is used intensively in particular piece of scalar CPU code; won't benefit from bouncing back and forth to the GPU */
4519566063dSJacob Faibussowitsch   PetscCall(VecBindToCPU(c->w1, PETSC_TRUE));
4529566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)c, (PetscObject)c->w1));
4539566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(c->w1, &c->w2));
454ce911d59SRichard Tran Mills   /* Vec is used intensively in particular piece of scalar CPU code; won't benefit from bouncing back and forth to the GPU */
4559566063dSJacob Faibussowitsch   PetscCall(VecBindToCPU(c->w2, PETSC_TRUE));
4569566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)c, (PetscObject)c->w2));
457b8f8c88eSHong Zhang 
45877d8c4bbSBarry Smith   c->error_rel    = PETSC_SQRT_MACHINE_EPSILON;
45977d8c4bbSBarry Smith   c->umin         = 100.0 * PETSC_SQRT_MACHINE_EPSILON;
46049b058dcSBarry Smith   c->currentcolor = -1;
461efb30889SBarry Smith   c->htype        = "wp";
4624e269d77SPeter Brune   c->fset         = PETSC_FALSE;
463c8a9c622SHong Zhang   c->setupcalled  = PETSC_FALSE;
464639f9d9dSBarry Smith 
465639f9d9dSBarry Smith   *color = c;
4669566063dSJacob Faibussowitsch   PetscCall(PetscObjectCompose((PetscObject)mat, "SNESMatFDColoring", (PetscObject)c));
4679566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(MAT_FDColoringCreate, mat, 0, 0, 0));
4683a40ed3dSBarry Smith   PetscFunctionReturn(0);
469bbf0e169SBarry Smith }
470bbf0e169SBarry Smith 
47105869f15SSatish Balay /*@
472639f9d9dSBarry Smith     MatFDColoringDestroy - Destroys a matrix coloring context that was created
473639f9d9dSBarry Smith     via MatFDColoringCreate().
474bbf0e169SBarry Smith 
475ef5ee4d1SLois Curfman McInnes     Collective on MatFDColoring
476ef5ee4d1SLois Curfman McInnes 
477b4fc646aSLois Curfman McInnes     Input Parameter:
478639f9d9dSBarry Smith .   c - coloring context
479bbf0e169SBarry Smith 
48015091d37SBarry Smith     Level: intermediate
48115091d37SBarry Smith 
482db781477SPatrick Sanan .seealso: `MatFDColoringCreate()`
483bbf0e169SBarry Smith @*/
4849371c9d4SSatish Balay PetscErrorCode MatFDColoringDestroy(MatFDColoring *c) {
48538baddfdSBarry Smith   PetscInt      i;
4860df34763SHong Zhang   MatFDColoring color = *c;
487bbf0e169SBarry Smith 
4883a40ed3dSBarry Smith   PetscFunctionBegin;
4896bf464f9SBarry Smith   if (!*c) PetscFunctionReturn(0);
4909371c9d4SSatish Balay   if (--((PetscObject)color)->refct > 0) {
4919371c9d4SSatish Balay     *c = NULL;
4929371c9d4SSatish Balay     PetscFunctionReturn(0);
4939371c9d4SSatish Balay   }
494d4bb536fSBarry Smith 
495071fcb05SBarry Smith   /* we do not free the column arrays since their entries are owned by the ISs in color->isa */
496*48a46eb9SPierre Jolivet   for (i = 0; i < color->ncolors; i++) PetscCall(ISDestroy(&color->isa[i]));
4979566063dSJacob Faibussowitsch   PetscCall(PetscFree(color->isa));
4989566063dSJacob Faibussowitsch   PetscCall(PetscFree2(color->ncolumns, color->columns));
4999566063dSJacob Faibussowitsch   PetscCall(PetscFree(color->nrows));
5000df34763SHong Zhang   if (color->htype[0] == 'w') {
5019566063dSJacob Faibussowitsch     PetscCall(PetscFree(color->matentry2));
5020df34763SHong Zhang   } else {
5039566063dSJacob Faibussowitsch     PetscCall(PetscFree(color->matentry));
5040df34763SHong Zhang   }
5059566063dSJacob Faibussowitsch   PetscCall(PetscFree(color->dy));
5069566063dSJacob Faibussowitsch   if (color->vscale) PetscCall(VecDestroy(&color->vscale));
5079566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&color->w1));
5089566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&color->w2));
5099566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&color->w3));
5109566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(c));
5113a40ed3dSBarry Smith   PetscFunctionReturn(0);
512bbf0e169SBarry Smith }
51343a90d84SBarry Smith 
51449b058dcSBarry Smith /*@C
51549b058dcSBarry Smith     MatFDColoringGetPerturbedColumns - Returns the indices of the columns that
51649b058dcSBarry Smith       that are currently being perturbed.
51749b058dcSBarry Smith 
51849b058dcSBarry Smith     Not Collective
51949b058dcSBarry Smith 
520f899ff85SJose E. Roman     Input Parameter:
52149b058dcSBarry Smith .   coloring - coloring context created with MatFDColoringCreate()
52249b058dcSBarry Smith 
52349b058dcSBarry Smith     Output Parameters:
52449b058dcSBarry Smith +   n - the number of local columns being perturbed
52549b058dcSBarry Smith -   cols - the column indices, in global numbering
52649b058dcSBarry Smith 
52721fcc2ddSBarry Smith    Level: advanced
52849b058dcSBarry Smith 
529e66d0d93SMatthew Knepley    Note: IF the matrix type is BAIJ, then the block column indices are returned
530e66d0d93SMatthew Knepley 
531edaa7c33SBarry Smith    Fortran Note:
532edaa7c33SBarry Smith    This routine has a different interface for Fortran
53321fcc2ddSBarry Smith $     #include <petsc/finclude/petscmat.h>
53421fcc2ddSBarry Smith $          use petscmat
535edaa7c33SBarry Smith $          PetscInt, pointer :: array(:)
536edaa7c33SBarry Smith $          PetscErrorCode  ierr
537edaa7c33SBarry Smith $          MatFDColoring   i
538edaa7c33SBarry Smith $          call MatFDColoringGetPerturbedColumnsF90(i,array,ierr)
539edaa7c33SBarry Smith $      use the entries of array ...
540edaa7c33SBarry Smith $          call MatFDColoringRestorePerturbedColumnsF90(i,array,ierr)
541edaa7c33SBarry Smith 
542db781477SPatrick Sanan .seealso: `MatFDColoringCreate()`, `MatFDColoringDestroy()`, `MatFDColoringView()`, `MatFDColoringApply()`
54349b058dcSBarry Smith 
54449b058dcSBarry Smith @*/
5459371c9d4SSatish Balay PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring coloring, PetscInt *n, const PetscInt *cols[]) {
54649b058dcSBarry Smith   PetscFunctionBegin;
54749b058dcSBarry Smith   if (coloring->currentcolor >= 0) {
54849b058dcSBarry Smith     *n    = coloring->ncolumns[coloring->currentcolor];
54949b058dcSBarry Smith     *cols = coloring->columns[coloring->currentcolor];
55049b058dcSBarry Smith   } else {
55149b058dcSBarry Smith     *n = 0;
55249b058dcSBarry Smith   }
55349b058dcSBarry Smith   PetscFunctionReturn(0);
55449b058dcSBarry Smith }
55549b058dcSBarry Smith 
55643a90d84SBarry Smith /*@
557e0907662SLois Curfman McInnes     MatFDColoringApply - Given a matrix for which a MatFDColoring context
558e0907662SLois Curfman McInnes     has been created, computes the Jacobian for a function via finite differences.
55943a90d84SBarry Smith 
560fee21e36SBarry Smith     Collective on MatFDColoring
561fee21e36SBarry Smith 
562ef5ee4d1SLois Curfman McInnes     Input Parameters:
563ef5ee4d1SLois Curfman McInnes +   mat - location to store Jacobian
564ef5ee4d1SLois Curfman McInnes .   coloring - coloring context created with MatFDColoringCreate()
565ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
5667850c7c0SBarry Smith -   sctx - context required by function, if this is being used with the SNES solver then it is SNES object, otherwise it is null
567ef5ee4d1SLois Curfman McInnes 
5688bba8e72SBarry Smith     Options Database Keys:
569ebd3b9afSBarry Smith +    -mat_fd_type - "wp" or "ds"  (see MATMFFD_WP or MATMFFD_DS)
570b9722508SLois Curfman McInnes .    -mat_fd_coloring_view - Activates basic viewing or coloring
571e350db43SBarry Smith .    -mat_fd_coloring_view draw - Activates drawing of coloring
572e350db43SBarry Smith -    -mat_fd_coloring_view ::ascii_info - Activates viewing of coloring info
5738bba8e72SBarry Smith 
57415091d37SBarry Smith     Level: intermediate
57598d79826SSatish Balay 
576db781477SPatrick Sanan .seealso: `MatFDColoringCreate()`, `MatFDColoringDestroy()`, `MatFDColoringView()`, `MatFDColoringSetFunction()`, `MatFDColoringSetValues()`
57743a90d84SBarry Smith 
57843a90d84SBarry Smith @*/
5799371c9d4SSatish Balay PetscErrorCode MatFDColoringApply(Mat J, MatFDColoring coloring, Vec x1, void *sctx) {
580bdaf1daeSBarry Smith   PetscBool eq;
5813acb8795SBarry Smith 
5823acb8795SBarry Smith   PetscFunctionBegin;
5830700a824SBarry Smith   PetscValidHeaderSpecific(J, MAT_CLASSID, 1);
5840700a824SBarry Smith   PetscValidHeaderSpecific(coloring, MAT_FDCOLORING_CLASSID, 2);
5850700a824SBarry Smith   PetscValidHeaderSpecific(x1, VEC_CLASSID, 3);
5869566063dSJacob Faibussowitsch   PetscCall(PetscObjectCompareId((PetscObject)J, coloring->matid, &eq));
58728b400f6SJacob Faibussowitsch   PetscCheck(eq, PetscObjectComm((PetscObject)J), PETSC_ERR_ARG_WRONG, "Matrix used with MatFDColoringApply() must be that used with MatFDColoringCreate()");
58828b400f6SJacob Faibussowitsch   PetscCheck(coloring->f, PetscObjectComm((PetscObject)J), PETSC_ERR_ARG_WRONGSTATE, "Must call MatFDColoringSetFunction()");
58928b400f6SJacob Faibussowitsch   PetscCheck(coloring->setupcalled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call MatFDColoringSetUp()");
590684f2004SHong Zhang 
5919566063dSJacob Faibussowitsch   PetscCall(MatSetUnfactored(J));
5929566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(MAT_FDColoringApply, coloring, J, x1, 0));
593dbbe0bcdSBarry Smith   PetscUseTypeMethod(J, fdcoloringapply, coloring, x1, sctx);
5949566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(MAT_FDColoringApply, coloring, J, x1, 0));
5955bdb020cSBarry Smith   if (!coloring->viewed) {
5969566063dSJacob Faibussowitsch     PetscCall(MatFDColoringViewFromOptions(coloring, NULL, "-mat_fd_coloring_view"));
5975bdb020cSBarry Smith     coloring->viewed = PETSC_TRUE;
5985bdb020cSBarry Smith   }
5993acb8795SBarry Smith   PetscFunctionReturn(0);
6003acb8795SBarry Smith }
601