xref: /petsc/src/mat/impls/mffd/mffd.c (revision 0b4b7b1c20c2ed4ade67e3d50a7710fe0ffbfca5)
1a1f56445SPierre Jolivet #include <../src/mat/impls/shell/shell.h> /*I  "petscmat.h"   I*/
2a1f56445SPierre Jolivet #include <../src/mat/impls/mffd/mffdimpl.h>
3e884886fSBarry Smith 
4f4259b30SLisandro Dalcin PetscFunctionList MatMFFDList              = NULL;
5ace3abfcSBarry Smith PetscBool         MatMFFDRegisterAllCalled = PETSC_FALSE;
6e884886fSBarry Smith 
77087cfbeSBarry Smith PetscClassId  MATMFFD_CLASSID;
8166c7f25SBarry Smith PetscLogEvent MATMFFD_Mult;
9e884886fSBarry Smith 
10ace3abfcSBarry Smith static PetscBool MatMFFDPackageInitialized = PETSC_FALSE;
1166976f2fSJacob Faibussowitsch 
12b022a5c1SBarry Smith /*@C
1311a5261eSBarry Smith   MatMFFDFinalizePackage - This function destroys everything in the MATMFFD` package. It is
1411a5261eSBarry Smith   called from `PetscFinalize()`.
15b022a5c1SBarry Smith 
16b022a5c1SBarry Smith   Level: developer
17b022a5c1SBarry Smith 
181cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `PetscFinalize()`, `MatCreateMFFD()`, `MatCreateSNESMF()`
19b022a5c1SBarry Smith @*/
20d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDFinalizePackage(void)
21d71ae5a4SJacob Faibussowitsch {
22b022a5c1SBarry Smith   PetscFunctionBegin;
239566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListDestroy(&MatMFFDList));
24b022a5c1SBarry Smith   MatMFFDPackageInitialized = PETSC_FALSE;
25b022a5c1SBarry Smith   MatMFFDRegisterAllCalled  = PETSC_FALSE;
263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
27b022a5c1SBarry Smith }
28b022a5c1SBarry Smith 
293ec795f1SBarry Smith /*@C
3011a5261eSBarry Smith   MatMFFDInitializePackage - This function initializes everything in the MATMFFD` package. It is called
3111a5261eSBarry Smith   from `MatInitializePackage()`.
323ec795f1SBarry Smith 
333ec795f1SBarry Smith   Level: developer
343ec795f1SBarry Smith 
351cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `PetscInitialize()`
363ec795f1SBarry Smith @*/
37d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDInitializePackage(void)
38d71ae5a4SJacob Faibussowitsch {
393ec795f1SBarry Smith   char      logList[256];
408e81d068SLisandro Dalcin   PetscBool opt, pkg;
413ec795f1SBarry Smith 
423ec795f1SBarry Smith   PetscFunctionBegin;
433ba16761SJacob Faibussowitsch   if (MatMFFDPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS);
44b022a5c1SBarry Smith   MatMFFDPackageInitialized = PETSC_TRUE;
453ec795f1SBarry Smith   /* Register Classes */
469566063dSJacob Faibussowitsch   PetscCall(PetscClassIdRegister("MatMFFD", &MATMFFD_CLASSID));
473ec795f1SBarry Smith   /* Register Constructors */
489566063dSJacob Faibussowitsch   PetscCall(MatMFFDRegisterAll());
493ec795f1SBarry Smith   /* Register Events */
509566063dSJacob Faibussowitsch   PetscCall(PetscLogEventRegister("MatMult MF", MATMFFD_CLASSID, &MATMFFD_Mult));
51e94e781bSJacob Faibussowitsch   /* Process Info */
52e94e781bSJacob Faibussowitsch   {
53e94e781bSJacob Faibussowitsch     PetscClassId classids[1];
54e94e781bSJacob Faibussowitsch 
55e94e781bSJacob Faibussowitsch     classids[0] = MATMFFD_CLASSID;
569566063dSJacob Faibussowitsch     PetscCall(PetscInfoProcessClass("matmffd", 1, classids));
573ec795f1SBarry Smith   }
583ec795f1SBarry Smith   /* Process summary exclusions */
599566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL, NULL, "-log_exclude", logList, sizeof(logList), &opt));
603ec795f1SBarry Smith   if (opt) {
619566063dSJacob Faibussowitsch     PetscCall(PetscStrInList("matmffd", logList, ',', &pkg));
629566063dSJacob Faibussowitsch     if (pkg) PetscCall(PetscLogEventExcludeClass(MATMFFD_CLASSID));
633ec795f1SBarry Smith   }
648e81d068SLisandro Dalcin   /* Register package finalizer */
659566063dSJacob Faibussowitsch   PetscCall(PetscRegisterFinalize(MatMFFDFinalizePackage));
663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
673ec795f1SBarry Smith }
683ec795f1SBarry Smith 
69d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDSetType_MFFD(Mat mat, MatMFFDType ftype)
70d71ae5a4SJacob Faibussowitsch {
71789d8953SBarry Smith   MatMFFD   ctx;
72789d8953SBarry Smith   PetscBool match;
735f80ce2aSJacob Faibussowitsch   PetscErrorCode (*r)(MatMFFD);
74789d8953SBarry Smith 
75789d8953SBarry Smith   PetscFunctionBegin;
76789d8953SBarry Smith   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
774f572ea9SToby Isaac   PetscAssertPointer(ftype, 2);
789566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
79789d8953SBarry Smith 
80789d8953SBarry Smith   /* already set, so just return */
819566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)ctx, ftype, &match));
823ba16761SJacob Faibussowitsch   if (match) PetscFunctionReturn(PETSC_SUCCESS);
83789d8953SBarry Smith 
84789d8953SBarry Smith   /* destroy the old one if it exists */
85dbbe0bcdSBarry Smith   PetscTryTypeMethod(ctx, destroy);
86789d8953SBarry Smith 
879566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(MatMFFDList, ftype, &r));
886adde796SStefano Zampini   PetscCheck(r, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown MatMFFD type %s given", ftype);
899566063dSJacob Faibussowitsch   PetscCall((*r)(ctx));
909566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)ctx, ftype));
913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
92789d8953SBarry Smith }
93789d8953SBarry Smith 
94cc4c1da9SBarry Smith /*@
95e884886fSBarry Smith   MatMFFDSetType - Sets the method that is used to compute the
96da81f932SPierre Jolivet   differencing parameter for finite difference matrix-free formulations.
97e884886fSBarry Smith 
98e884886fSBarry Smith   Input Parameters:
9911a5261eSBarry Smith + mat   - the "matrix-free" matrix created via `MatCreateSNESMF()`, or `MatCreateMFFD()`
10011a5261eSBarry Smith           or `MatSetType`(mat,`MATMFFD`);
10111a5261eSBarry Smith - ftype - the type requested, either `MATMFFD_WP` or `MATMFFD_DS`
102e884886fSBarry Smith 
103e884886fSBarry Smith   Level: advanced
104e884886fSBarry Smith 
10511a5261eSBarry Smith   Note:
1062ef1f0ffSBarry Smith   For example, such routines can compute `h` for use in
107e884886fSBarry Smith   Jacobian-vector products of the form
1082ef1f0ffSBarry Smith .vb
109e884886fSBarry Smith 
110e884886fSBarry Smith                         F(x+ha) - F(x)
111e884886fSBarry Smith           F'(u)a  ~=  ----------------
112e884886fSBarry Smith                               h
1132ef1f0ffSBarry Smith .ve
114e884886fSBarry Smith 
1151cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MATMFFD_WP`, `MATMFFD_DS`, `MatCreateSNESMF()`, `MatMFFDRegister()`, `MatMFFDSetFunction()`, `MatCreateMFFD()`
116e884886fSBarry Smith @*/
117d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetType(Mat mat, MatMFFDType ftype)
118d71ae5a4SJacob Faibussowitsch {
119e884886fSBarry Smith   PetscFunctionBegin;
1200700a824SBarry Smith   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
1214f572ea9SToby Isaac   PetscAssertPointer(ftype, 2);
122cac4c232SBarry Smith   PetscTryMethod(mat, "MatMFFDSetType_C", (Mat, MatMFFDType), (mat, ftype));
1233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
124e884886fSBarry Smith }
125e884886fSBarry Smith 
1265d21e1f8SStefano Zampini static PetscErrorCode MatGetDiagonal_MFFD(Mat, Vec);
1275d21e1f8SStefano Zampini 
128e884886fSBarry Smith typedef PetscErrorCode (*FCN1)(void *, Vec); /* force argument to next function to not be extern C*/
129d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDSetFunctioniBase_MFFD(Mat mat, FCN1 func)
130d71ae5a4SJacob Faibussowitsch {
131789d8953SBarry Smith   MatMFFD ctx;
132e884886fSBarry Smith 
133e884886fSBarry Smith   PetscFunctionBegin;
1349566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
135e884886fSBarry Smith   ctx->funcisetbase = func;
1363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
137e884886fSBarry Smith }
138e884886fSBarry Smith 
139e884886fSBarry Smith typedef PetscErrorCode (*FCN2)(void *, PetscInt, Vec, PetscScalar *); /* force argument to next function to not be extern C*/
140d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDSetFunctioni_MFFD(Mat mat, FCN2 funci)
141d71ae5a4SJacob Faibussowitsch {
142789d8953SBarry Smith   MatMFFD ctx;
143e884886fSBarry Smith 
144e884886fSBarry Smith   PetscFunctionBegin;
1459566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
146e884886fSBarry Smith   ctx->funci = funci;
1479566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(mat, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiagonal_MFFD));
1483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1495d21e1f8SStefano Zampini }
150789d8953SBarry Smith 
151d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDGetH_MFFD(Mat mat, PetscScalar *h)
152d71ae5a4SJacob Faibussowitsch {
153789d8953SBarry Smith   MatMFFD ctx;
154789d8953SBarry Smith 
155789d8953SBarry Smith   PetscFunctionBegin;
1569566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
157789d8953SBarry Smith   *h = ctx->currenth;
1583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
159e884886fSBarry Smith }
160e884886fSBarry Smith 
161d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDResetHHistory_MFFD(Mat J)
162d71ae5a4SJacob Faibussowitsch {
163789d8953SBarry Smith   MatMFFD ctx;
164bc0f08ceSBarry Smith 
165bc0f08ceSBarry Smith   PetscFunctionBegin;
1669566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J, &ctx));
167bc0f08ceSBarry Smith   ctx->ncurrenth = 0;
1683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
169bc0f08ceSBarry Smith }
170e884886fSBarry Smith 
1711c84c290SBarry Smith /*@C
17220f4b53cSBarry Smith   MatMFFDRegister - Adds a method to the `MATMFFD` registry.
1731c84c290SBarry Smith 
174cc4c1da9SBarry Smith   Not Collective, No Fortran Support
1751c84c290SBarry Smith 
1761c84c290SBarry Smith   Input Parameters:
17720f4b53cSBarry Smith + sname    - name of a new user-defined compute-h module
17820f4b53cSBarry Smith - function - routine to create method context
1791c84c290SBarry Smith 
1801c84c290SBarry Smith   Level: developer
1811c84c290SBarry Smith 
18211a5261eSBarry Smith   Note:
18311a5261eSBarry Smith   `MatMFFDRegister()` may be called multiple times to add several user-defined solvers.
1841c84c290SBarry Smith 
185fe59aa6dSJacob Faibussowitsch   Example Usage:
1861c84c290SBarry Smith .vb
187bdf89e91SBarry Smith    MatMFFDRegister("my_h", MyHCreate);
1881c84c290SBarry Smith .ve
1891c84c290SBarry Smith 
190a3b724e8SBarry Smith   Then, your solver can be chosen with the procedural interface via `MatMFFDSetType`(mfctx, "my_h")` or at runtime via the option
191a3b724e8SBarry Smith   `-mat_mffd_type my_h`
1921c84c290SBarry Smith 
1931cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatMFFDRegisterAll()`, `MatMFFDRegisterDestroy()`
1941c84c290SBarry Smith  @*/
195d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDRegister(const char sname[], PetscErrorCode (*function)(MatMFFD))
196d71ae5a4SJacob Faibussowitsch {
197e884886fSBarry Smith   PetscFunctionBegin;
1989566063dSJacob Faibussowitsch   PetscCall(MatInitializePackage());
1999566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&MatMFFDList, sname, function));
2003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
201e884886fSBarry Smith }
202e884886fSBarry Smith 
203d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDestroy_MFFD(Mat mat)
204d71ae5a4SJacob Faibussowitsch {
205789d8953SBarry Smith   MatMFFD ctx;
206e884886fSBarry Smith 
207e884886fSBarry Smith   PetscFunctionBegin;
2089566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
2099566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ctx->w));
2109566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ctx->current_u));
21148a46eb9SPierre Jolivet   if (ctx->current_f_allocated) PetscCall(VecDestroy(&ctx->current_f));
212dbbe0bcdSBarry Smith   PetscTryTypeMethod(ctx, destroy);
2139566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(&ctx));
214e884886fSBarry Smith 
2159566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetBase_C", NULL));
2169566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetFunctioniBase_C", NULL));
2179566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetFunctioni_C", NULL));
2189566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetFunction_C", NULL));
2199566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetFunctionError_C", NULL));
2209566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetCheckh_C", NULL));
2219566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetPeriod_C", NULL));
2229566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDResetHHistory_C", NULL));
2239566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetHHistory_C", NULL));
2249566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetType_C", NULL));
2259566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDGetH_C", NULL));
2262e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatSNESMFSetReuseBase_C", NULL));
2272e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatSNESMFGetReuseBase_C", NULL));
228ee912ba9SBarry Smith   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatShellSetContext_C", NULL));
2293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
230e884886fSBarry Smith }
231e884886fSBarry Smith 
232e884886fSBarry Smith /*
233e884886fSBarry Smith    MatMFFDView_MFFD - Views matrix-free parameters.
234e884886fSBarry Smith 
235e884886fSBarry Smith */
236d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_MFFD(Mat J, PetscViewer viewer)
237d71ae5a4SJacob Faibussowitsch {
238789d8953SBarry Smith   MatMFFD     ctx;
239a283635eSDmitry Karpeev   PetscBool   iascii, viewbase, viewfunction;
240a283635eSDmitry Karpeev   const char *prefix;
241e884886fSBarry Smith 
242e884886fSBarry Smith   PetscFunctionBegin;
2439566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J, &ctx));
2449566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
245e884886fSBarry Smith   if (iascii) {
2469566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Matrix-free approximation:\n"));
2479566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
2489566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "err=%g (relative error in function evaluation)\n", (double)ctx->error_rel));
2497adad957SLisandro Dalcin     if (!((PetscObject)ctx)->type_name) {
2509566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "The compute h routine has not yet been set\n"));
251e884886fSBarry Smith     } else {
2529566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Using %s compute h routine\n", ((PetscObject)ctx)->type_name));
253e884886fSBarry Smith     }
254c92e3469SBarry Smith #if defined(PETSC_USE_COMPLEX)
25548a46eb9SPierre Jolivet     if (ctx->usecomplex) PetscCall(PetscViewerASCIIPrintf(viewer, "Using Lyness complex number trick to compute the matrix-vector product\n"));
256c92e3469SBarry Smith #endif
257dbbe0bcdSBarry Smith     PetscTryTypeMethod(ctx, view, viewer);
2589566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetOptionsPrefix((PetscObject)J, &prefix));
259a283635eSDmitry Karpeev 
2609566063dSJacob Faibussowitsch     PetscCall(PetscOptionsHasName(((PetscObject)J)->options, prefix, "-mat_mffd_view_base", &viewbase));
261a283635eSDmitry Karpeev     if (viewbase) {
2629566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Base:\n"));
2639566063dSJacob Faibussowitsch       PetscCall(VecView(ctx->current_u, viewer));
264a283635eSDmitry Karpeev     }
2659566063dSJacob Faibussowitsch     PetscCall(PetscOptionsHasName(((PetscObject)J)->options, prefix, "-mat_mffd_view_function", &viewfunction));
266a283635eSDmitry Karpeev     if (viewfunction) {
2679566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Function:\n"));
2689566063dSJacob Faibussowitsch       PetscCall(VecView(ctx->current_f, viewer));
269a283635eSDmitry Karpeev     }
2709566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
27111aeaf0aSBarry Smith   }
2723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
273e884886fSBarry Smith }
274e884886fSBarry Smith 
275e884886fSBarry Smith /*
276e884886fSBarry Smith    MatAssemblyEnd_MFFD - Resets the ctx->ncurrenth to zero. This
277e884886fSBarry Smith    allows the user to indicate the beginning of a new linear solve by calling
27801c1178eSBarry Smith    MatAssemblyXXX() on the matrix-free matrix. This then allows the
2791d0fab5eSBarry Smith    MatCreateMFFD_WP() to properly compute ||U|| only the first time
280e884886fSBarry Smith    in the linear solver rather than every time.
2815a576424SJed Brown 
282cc2e6a90SBarry Smith    This function is referenced directly from MatAssemblyEnd_SNESMF(), which may be in a different shared library hence
283cc2e6a90SBarry Smith    it must be labeled as PETSC_EXTERN
284e884886fSBarry Smith */
285d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatAssemblyEnd_MFFD(Mat J, MatAssemblyType mt)
286d71ae5a4SJacob Faibussowitsch {
287789d8953SBarry Smith   MatMFFD j;
288e884886fSBarry Smith 
289e884886fSBarry Smith   PetscFunctionBegin;
2909566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J, &j));
2919566063dSJacob Faibussowitsch   PetscCall(MatMFFDResetHHistory(J));
2923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
293e884886fSBarry Smith }
294e884886fSBarry Smith 
295e884886fSBarry Smith /*
296e884886fSBarry Smith   MatMult_MFFD - Default matrix-free form for Jacobian-vector product, y = F'(u)*a:
297e884886fSBarry Smith 
298e884886fSBarry Smith         y ~= (F(u + ha) - F(u))/h,
299e884886fSBarry Smith   where F = nonlinear function, as set by SNESSetFunction()
300e884886fSBarry Smith         u = current iterate
301e884886fSBarry Smith         h = difference interval
302e884886fSBarry Smith */
303d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMult_MFFD(Mat mat, Vec a, Vec y)
304d71ae5a4SJacob Faibussowitsch {
305789d8953SBarry Smith   MatMFFD     ctx;
306e884886fSBarry Smith   PetscScalar h;
307e884886fSBarry Smith   Vec         w, U, F;
308ace3abfcSBarry Smith   PetscBool   zeroa;
309e884886fSBarry Smith 
310e884886fSBarry Smith   PetscFunctionBegin;
3119566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
31200045ab3SPierre Jolivet   PetscCheck(ctx->current_u, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "MatMFFDSetBase() has not been called, this is often caused by forgetting to call MatAssemblyBegin/End on the first Mat in the SNES compute function");
313e884886fSBarry Smith   /* We log matrix-free matrix-vector products separately, so that we can
314e884886fSBarry Smith      separate the performance monitoring from the cases that use conventional
315e884886fSBarry Smith      storage.  We may eventually modify event logging to associate events
316e884886fSBarry Smith      with particular objects, hence alleviating the more general problem. */
3179566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(MATMFFD_Mult, a, y, 0, 0));
318e884886fSBarry Smith 
319e884886fSBarry Smith   w = ctx->w;
320e884886fSBarry Smith   U = ctx->current_u;
3213ec795f1SBarry Smith   F = ctx->current_f;
322e884886fSBarry Smith   /*
323e884886fSBarry Smith       Compute differencing parameter
324e884886fSBarry Smith   */
3254a2f8832SBarry Smith   if (!((PetscObject)ctx)->type_name) {
3269566063dSJacob Faibussowitsch     PetscCall(MatMFFDSetType(mat, MATMFFD_WP));
3279566063dSJacob Faibussowitsch     PetscCall(MatSetFromOptions(mat));
328e884886fSBarry Smith   }
329dbbe0bcdSBarry Smith   PetscUseTypeMethod(ctx, compute, U, a, &h, &zeroa);
330e884886fSBarry Smith   if (zeroa) {
3319566063dSJacob Faibussowitsch     PetscCall(VecSet(y, 0.0));
33285d14658SLisandro Dalcin     PetscCall(PetscLogEventEnd(MATMFFD_Mult, a, y, 0, 0));
3333ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
334e884886fSBarry Smith   }
335e884886fSBarry Smith 
33608401ef6SPierre Jolivet   PetscCheck(!mat->erroriffailure || !PetscIsInfOrNanScalar(h), PETSC_COMM_SELF, PETSC_ERR_PLIB, "Computed Nan differencing parameter h");
33748a46eb9SPierre Jolivet   if (ctx->checkh) PetscCall((*ctx->checkh)(ctx->checkhctx, U, a, &h));
338e884886fSBarry Smith 
339e884886fSBarry Smith   /* keep a record of the current differencing parameter h */
340e884886fSBarry Smith   ctx->currenth = h;
341e884886fSBarry Smith #if defined(PETSC_USE_COMPLEX)
3429566063dSJacob Faibussowitsch   PetscCall(PetscInfo(mat, "Current differencing parameter: %g + %g i\n", (double)PetscRealPart(h), (double)PetscImaginaryPart(h)));
343e884886fSBarry Smith #else
3449566063dSJacob Faibussowitsch   PetscCall(PetscInfo(mat, "Current differencing parameter: %15.12e\n", (double)PetscRealPart(h)));
345e884886fSBarry Smith #endif
346ad540459SPierre Jolivet   if (ctx->historyh && ctx->ncurrenth < ctx->maxcurrenth) ctx->historyh[ctx->ncurrenth] = h;
347e884886fSBarry Smith   ctx->ncurrenth++;
348e884886fSBarry Smith 
349c92e3469SBarry Smith #if defined(PETSC_USE_COMPLEX)
350c92e3469SBarry Smith   if (ctx->usecomplex) h = PETSC_i * h;
351c92e3469SBarry Smith #endif
352c92e3469SBarry Smith 
353e884886fSBarry Smith   /* w = u + ha */
3549566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(w, h, a, U));
355e884886fSBarry Smith 
3561b797266SDmitry Karpeev   /* compute func(U) as base for differencing; only needed first time in and not when provided by user */
35748a46eb9SPierre Jolivet   if (ctx->ncurrenth == 1 && ctx->current_f_allocated) PetscCall((*ctx->func)(ctx->funcctx, U, F));
3589566063dSJacob Faibussowitsch   PetscCall((*ctx->func)(ctx->funcctx, w, y));
359e884886fSBarry Smith 
360c92e3469SBarry Smith #if defined(PETSC_USE_COMPLEX)
361c92e3469SBarry Smith   if (ctx->usecomplex) {
3629566063dSJacob Faibussowitsch     PetscCall(VecImaginaryPart(y));
363c92e3469SBarry Smith     h = PetscImaginaryPart(h);
364c92e3469SBarry Smith   } else {
3659566063dSJacob Faibussowitsch     PetscCall(VecAXPY(y, -1.0, F));
366c92e3469SBarry Smith   }
367c92e3469SBarry Smith #else
3689566063dSJacob Faibussowitsch   PetscCall(VecAXPY(y, -1.0, F));
369c92e3469SBarry Smith #endif
3709566063dSJacob Faibussowitsch   PetscCall(VecScale(y, 1.0 / h));
3719566063dSJacob Faibussowitsch   if (mat->nullsp) PetscCall(MatNullSpaceRemove(mat->nullsp, y));
372e884886fSBarry Smith 
3739566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(MATMFFD_Mult, a, y, 0, 0));
3743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
375e884886fSBarry Smith }
376e884886fSBarry Smith 
377e884886fSBarry Smith /*
37801c1178eSBarry Smith   MatGetDiagonal_MFFD - Gets the diagonal for a matrix-free matrix
379e884886fSBarry Smith 
380e884886fSBarry Smith         y ~= (F(u + ha) - F(u))/h,
381e884886fSBarry Smith   where F = nonlinear function, as set by SNESSetFunction()
382e884886fSBarry Smith         u = current iterate
383e884886fSBarry Smith         h = difference interval
384e884886fSBarry Smith */
38566976f2fSJacob Faibussowitsch static PetscErrorCode MatGetDiagonal_MFFD(Mat mat, Vec a)
386d71ae5a4SJacob Faibussowitsch {
387789d8953SBarry Smith   MatMFFD     ctx;
388e884886fSBarry Smith   PetscScalar h, *aa, *ww, v;
389e884886fSBarry Smith   PetscReal   epsilon = PETSC_SQRT_MACHINE_EPSILON, umin = 100.0 * PETSC_SQRT_MACHINE_EPSILON;
390e884886fSBarry Smith   Vec         w, U;
391e884886fSBarry Smith   PetscInt    i, rstart, rend;
392e884886fSBarry Smith 
393e884886fSBarry Smith   PetscFunctionBegin;
3949566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
39528b400f6SJacob Faibussowitsch   PetscCheck(ctx->func, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Requires calling MatMFFDSetFunction() first");
39628b400f6SJacob Faibussowitsch   PetscCheck(ctx->funci, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Requires calling MatMFFDSetFunctioni() first");
397e884886fSBarry Smith   w = ctx->w;
398e884886fSBarry Smith   U = ctx->current_u;
3999566063dSJacob Faibussowitsch   PetscCall((*ctx->func)(ctx->funcctx, U, a));
4001baa6e33SBarry Smith   if (ctx->funcisetbase) PetscCall((*ctx->funcisetbase)(ctx->funcctx, U));
4019566063dSJacob Faibussowitsch   PetscCall(VecCopy(U, w));
402e884886fSBarry Smith 
4039566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(a, &rstart, &rend));
4049566063dSJacob Faibussowitsch   PetscCall(VecGetArray(a, &aa));
405e884886fSBarry Smith   for (i = rstart; i < rend; i++) {
4069566063dSJacob Faibussowitsch     PetscCall(VecGetArray(w, &ww));
407e884886fSBarry Smith     h = ww[i - rstart];
408e884886fSBarry Smith     if (h == 0.0) h = 1.0;
409e884886fSBarry Smith     if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0) h = umin;
410e884886fSBarry Smith     else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin;
411e884886fSBarry Smith     h *= epsilon;
412e884886fSBarry Smith 
413e884886fSBarry Smith     ww[i - rstart] += h;
4149566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(w, &ww));
4159566063dSJacob Faibussowitsch     PetscCall((*ctx->funci)(ctx->funcctx, i, w, &v));
416e884886fSBarry Smith     aa[i - rstart] = (v - aa[i - rstart]) / h;
417e884886fSBarry Smith 
4189566063dSJacob Faibussowitsch     PetscCall(VecGetArray(w, &ww));
419e884886fSBarry Smith     ww[i - rstart] -= h;
4209566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(w, &ww));
421e884886fSBarry Smith   }
4229566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(a, &aa));
4233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
424e884886fSBarry Smith }
425e884886fSBarry Smith 
426d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatMFFDSetBase_MFFD(Mat J, Vec U, Vec F)
427d71ae5a4SJacob Faibussowitsch {
428789d8953SBarry Smith   MatMFFD ctx;
429e884886fSBarry Smith 
430e884886fSBarry Smith   PetscFunctionBegin;
4319566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J, &ctx));
4329566063dSJacob Faibussowitsch   PetscCall(MatMFFDResetHHistory(J));
433c51ad4d4SStefano Zampini   if (!ctx->current_u) {
4349566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(U, &ctx->current_u));
4359566063dSJacob Faibussowitsch     PetscCall(VecLockReadPush(ctx->current_u));
436c51ad4d4SStefano Zampini   }
4379566063dSJacob Faibussowitsch   PetscCall(VecLockReadPop(ctx->current_u));
4389566063dSJacob Faibussowitsch   PetscCall(VecCopy(U, ctx->current_u));
4399566063dSJacob Faibussowitsch   PetscCall(VecLockReadPush(ctx->current_u));
44052121784SBarry Smith   if (F) {
4419566063dSJacob Faibussowitsch     if (ctx->current_f_allocated) PetscCall(VecDestroy(&ctx->current_f));
4423ec795f1SBarry Smith     ctx->current_f           = F;
443cfe22f5eSBarry Smith     ctx->current_f_allocated = PETSC_FALSE;
444cfe22f5eSBarry Smith   } else if (!ctx->current_f_allocated) {
4459566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(J, NULL, &ctx->current_f));
446cfe22f5eSBarry Smith     ctx->current_f_allocated = PETSC_TRUE;
44752121784SBarry Smith   }
44848a46eb9SPierre Jolivet   if (!ctx->w) PetscCall(VecDuplicate(ctx->current_u, &ctx->w));
449e884886fSBarry Smith   J->assembled = PETSC_TRUE;
4503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
451e884886fSBarry Smith }
4525a576424SJed Brown 
453e884886fSBarry Smith typedef PetscErrorCode (*FCN3)(void *, Vec, Vec, PetscScalar *); /* force argument to next function to not be extern C*/
454b2573a8aSBarry Smith 
455d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDSetCheckh_MFFD(Mat J, FCN3 fun, void *ectx)
456d71ae5a4SJacob Faibussowitsch {
457789d8953SBarry Smith   MatMFFD ctx;
458e884886fSBarry Smith 
459e884886fSBarry Smith   PetscFunctionBegin;
4609566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J, &ctx));
461e884886fSBarry Smith   ctx->checkh    = fun;
462e884886fSBarry Smith   ctx->checkhctx = ectx;
4633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
464e884886fSBarry Smith }
465e884886fSBarry Smith 
466cc4c1da9SBarry Smith /*@
4676aa9148fSLisandro Dalcin   MatMFFDSetOptionsPrefix - Sets the prefix used for searching for all
46811a5261eSBarry Smith   MATMFFD` options in the database.
4696aa9148fSLisandro Dalcin 
470c3339decSBarry Smith   Collective
4716aa9148fSLisandro Dalcin 
472d8d19677SJose E. Roman   Input Parameters:
473fe59aa6dSJacob Faibussowitsch + mat    - the `MATMFFD` context
4746aa9148fSLisandro Dalcin - prefix - the prefix to prepend to all option names
4756aa9148fSLisandro Dalcin 
47611a5261eSBarry Smith   Note:
4776aa9148fSLisandro Dalcin   A hyphen (-) must NOT be given at the beginning of the prefix name.
4786aa9148fSLisandro Dalcin   The first character of all runtime options is AUTOMATICALLY the hyphen.
4796aa9148fSLisandro Dalcin 
4806aa9148fSLisandro Dalcin   Level: advanced
4816aa9148fSLisandro Dalcin 
4821cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatSetFromOptions()`, `MatCreateSNESMF()`, `MatCreateMFFD()`
4836aa9148fSLisandro Dalcin @*/
484d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetOptionsPrefix(Mat mat, const char prefix[])
485d71ae5a4SJacob Faibussowitsch {
486789d8953SBarry Smith   MatMFFD mfctx;
4875fd66863SKarl Rupp 
4886aa9148fSLisandro Dalcin   PetscFunctionBegin;
4890700a824SBarry Smith   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
4909566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &mfctx));
4910700a824SBarry Smith   PetscValidHeaderSpecific(mfctx, MATMFFD_CLASSID, 1);
4929566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)mfctx, prefix));
4933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4946aa9148fSLisandro Dalcin }
4956aa9148fSLisandro Dalcin 
496d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSetFromOptions_MFFD(Mat mat, PetscOptionItems *PetscOptionsObject)
497d71ae5a4SJacob Faibussowitsch {
498789d8953SBarry Smith   MatMFFD   mfctx;
499ace3abfcSBarry Smith   PetscBool flg;
500e884886fSBarry Smith   char      ftype[256];
501e884886fSBarry Smith 
502e884886fSBarry Smith   PetscFunctionBegin;
5039566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &mfctx));
504dbbe0bcdSBarry Smith   PetscValidHeaderSpecific(mfctx, MATMFFD_CLASSID, 1);
505d0609cedSBarry Smith   PetscObjectOptionsBegin((PetscObject)mfctx);
5069566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFList("-mat_mffd_type", "Matrix free type", "MatMFFDSetType", MatMFFDList, ((PetscObject)mfctx)->type_name, ftype, 256, &flg));
5071baa6e33SBarry Smith   if (flg) PetscCall(MatMFFDSetType(mat, ftype));
508e884886fSBarry Smith 
5099566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_mffd_err", "set sqrt relative error in function", "MatMFFDSetFunctionError", mfctx->error_rel, &mfctx->error_rel, NULL));
5109566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-mat_mffd_period", "how often h is recomputed", "MatMFFDSetPeriod", mfctx->recomputeperiod, &mfctx->recomputeperiod, NULL));
511e884886fSBarry Smith 
51290d69ab7SBarry Smith   flg = PETSC_FALSE;
5139566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-mat_mffd_check_positivity", "Insure that U + h*a is nonnegative", "MatMFFDSetCheckh", flg, &flg, NULL));
5141baa6e33SBarry Smith   if (flg) PetscCall(MatMFFDSetCheckh(mat, MatMFFDCheckPositivity, NULL));
515c92e3469SBarry Smith #if defined(PETSC_USE_COMPLEX)
5169566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-mat_mffd_complex", "Use Lyness complex number trick to compute the matrix-vector product", "None", mfctx->usecomplex, &mfctx->usecomplex, NULL));
517c92e3469SBarry Smith #endif
518dbbe0bcdSBarry Smith   PetscTryTypeMethod(mfctx, setfromoptions, PetscOptionsObject);
519d0609cedSBarry Smith   PetscOptionsEnd();
5203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
521e884886fSBarry Smith }
522e884886fSBarry Smith 
523d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDSetPeriod_MFFD(Mat mat, PetscInt period)
524d71ae5a4SJacob Faibussowitsch {
525789d8953SBarry Smith   MatMFFD ctx;
526bc0f08ceSBarry Smith 
527bc0f08ceSBarry Smith   PetscFunctionBegin;
5289566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
529bc0f08ceSBarry Smith   ctx->recomputeperiod = period;
5303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
531bc0f08ceSBarry Smith }
532bc0f08ceSBarry Smith 
533d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDSetFunction_MFFD(Mat mat, PetscErrorCode (*func)(void *, Vec, Vec), void *funcctx)
534d71ae5a4SJacob Faibussowitsch {
535789d8953SBarry Smith   MatMFFD ctx;
536bc0f08ceSBarry Smith 
537bc0f08ceSBarry Smith   PetscFunctionBegin;
5389566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
539bc0f08ceSBarry Smith   ctx->func    = func;
540bc0f08ceSBarry Smith   ctx->funcctx = funcctx;
5413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
542bc0f08ceSBarry Smith }
543bc0f08ceSBarry Smith 
544d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDSetFunctionError_MFFD(Mat mat, PetscReal error)
545d71ae5a4SJacob Faibussowitsch {
54613bcc0bdSJacob Faibussowitsch   PetscFunctionBegin;
54713bcc0bdSJacob Faibussowitsch   if (error != (PetscReal)PETSC_DEFAULT) {
548789d8953SBarry Smith     MatMFFD ctx;
549bc0f08ceSBarry Smith 
5509566063dSJacob Faibussowitsch     PetscCall(MatShellGetContext(mat, &ctx));
55113bcc0bdSJacob Faibussowitsch     ctx->error_rel = error;
55213bcc0bdSJacob Faibussowitsch   }
5533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
554bc0f08ceSBarry Smith }
555bc0f08ceSBarry Smith 
55666976f2fSJacob Faibussowitsch static PetscErrorCode MatMFFDSetHHistory_MFFD(Mat J, PetscScalar history[], PetscInt nhistory)
557d71ae5a4SJacob Faibussowitsch {
558789d8953SBarry Smith   MatMFFD ctx;
559789d8953SBarry Smith 
5603b49f96aSBarry Smith   PetscFunctionBegin;
5619566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J, &ctx));
562789d8953SBarry Smith   ctx->historyh    = history;
563789d8953SBarry Smith   ctx->maxcurrenth = nhistory;
564789d8953SBarry Smith   ctx->currenth    = 0.;
5653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5663b49f96aSBarry Smith }
5673b49f96aSBarry Smith 
568e884886fSBarry Smith /*MC
56901c1178eSBarry Smith   MATMFFD - "mffd" - A matrix-free matrix type.
570e884886fSBarry Smith 
571e884886fSBarry Smith   Level: advanced
572e884886fSBarry Smith 
573a1f56445SPierre Jolivet   Developer Notes:
5742ef1f0ffSBarry Smith   This is implemented on top of `MATSHELL` to get support for scaling and shifting without requiring duplicate code
575789d8953SBarry Smith 
576ee912ba9SBarry Smith   Users should not MatShell... operations on this class, there is some error checking for that incorrect usage
577ee912ba9SBarry Smith 
5781cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreateMFFD()`, `MatCreateSNESMF()`, `MatMFFDSetFunction()`, `MatMFFDSetType()`,
579db781477SPatrick Sanan           `MatMFFDSetFunctionError()`, `MatMFFDDSSetUmin()`, `MatMFFDSetFunction()`
580db781477SPatrick Sanan           `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`, `MatCreateSNESMF()`,
581db781477SPatrick Sanan           `MatMFFDGetH()`,
582e884886fSBarry Smith M*/
583d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatCreate_MFFD(Mat A)
584d71ae5a4SJacob Faibussowitsch {
585e884886fSBarry Smith   MatMFFD mfctx;
586e884886fSBarry Smith 
587e884886fSBarry Smith   PetscFunctionBegin;
5889566063dSJacob Faibussowitsch   PetscCall(MatMFFDInitializePackage());
589e884886fSBarry Smith 
5909566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(mfctx, MATMFFD_CLASSID, "MatMFFD", "Matrix-free Finite Differencing", "Mat", PetscObjectComm((PetscObject)A), NULL, NULL));
5912205254eSKarl Rupp 
592e884886fSBarry Smith   mfctx->error_rel                = PETSC_SQRT_MACHINE_EPSILON;
593e884886fSBarry Smith   mfctx->recomputeperiod          = 1;
594e884886fSBarry Smith   mfctx->count                    = 0;
595e884886fSBarry Smith   mfctx->currenth                 = 0.0;
5960298fd71SBarry Smith   mfctx->historyh                 = NULL;
597e884886fSBarry Smith   mfctx->ncurrenth                = 0;
598e884886fSBarry Smith   mfctx->maxcurrenth              = 0;
599f4259b30SLisandro Dalcin   ((PetscObject)mfctx)->type_name = NULL;
600e884886fSBarry Smith 
601e884886fSBarry Smith   /*
602e884886fSBarry Smith      Create the empty data structure to contain compute-h routines.
603e884886fSBarry Smith      These will be filled in below from the command line options or
604e884886fSBarry Smith      a later call with MatMFFDSetType() or if that is not called
605e884886fSBarry Smith      then it will default in the first use of MatMult_MFFD()
606e884886fSBarry Smith   */
607f4259b30SLisandro Dalcin   mfctx->ops->compute        = NULL;
608f4259b30SLisandro Dalcin   mfctx->ops->destroy        = NULL;
609f4259b30SLisandro Dalcin   mfctx->ops->view           = NULL;
610f4259b30SLisandro Dalcin   mfctx->ops->setfromoptions = NULL;
611f4259b30SLisandro Dalcin   mfctx->hctx                = NULL;
612e884886fSBarry Smith 
613f4259b30SLisandro Dalcin   mfctx->func    = NULL;
614f4259b30SLisandro Dalcin   mfctx->funcctx = NULL;
6150298fd71SBarry Smith   mfctx->w       = NULL;
616789d8953SBarry Smith   mfctx->mat     = A;
617e884886fSBarry Smith 
6189566063dSJacob Faibussowitsch   PetscCall(MatSetType(A, MATSHELL));
6199566063dSJacob Faibussowitsch   PetscCall(MatShellSetContext(A, mfctx));
6209566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(A, MATOP_MULT, (void (*)(void))MatMult_MFFD));
6219566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(A, MATOP_DESTROY, (void (*)(void))MatDestroy_MFFD));
6229566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(A, MATOP_VIEW, (void (*)(void))MatView_MFFD));
6239566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(A, MATOP_ASSEMBLY_END, (void (*)(void))MatAssemblyEnd_MFFD));
6249566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(A, MATOP_SET_FROM_OPTIONS, (void (*)(void))MatSetFromOptions_MFFD));
625a1f56445SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatShellSetContext_C", MatShellSetContext_Immutable));
626a1f56445SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatShellSetContextDestroy_C", MatShellSetContextDestroy_Immutable));
627a1f56445SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatShellSetManageScalingShifts_C", MatShellSetManageScalingShifts_Immutable));
628e884886fSBarry Smith 
6299566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetBase_C", MatMFFDSetBase_MFFD));
6309566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetFunctioniBase_C", MatMFFDSetFunctioniBase_MFFD));
6319566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetFunctioni_C", MatMFFDSetFunctioni_MFFD));
6329566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetFunction_C", MatMFFDSetFunction_MFFD));
6339566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetCheckh_C", MatMFFDSetCheckh_MFFD));
6349566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetPeriod_C", MatMFFDSetPeriod_MFFD));
6359566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetFunctionError_C", MatMFFDSetFunctionError_MFFD));
6369566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDResetHHistory_C", MatMFFDResetHHistory_MFFD));
6379566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetHHistory_C", MatMFFDSetHHistory_MFFD));
6389566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetType_C", MatMFFDSetType_MFFD));
6399566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDGetH_C", MatMFFDGetH_MFFD));
6409566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATMFFD));
6413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
642e884886fSBarry Smith }
643e884886fSBarry Smith 
644e884886fSBarry Smith /*@
645*0b4b7b1cSBarry Smith   MatCreateMFFD - Creates a matrix-free matrix of type `MATMFFD` that uses finite differences on a provided function to
646*0b4b7b1cSBarry Smith   approximately multiply a vector by the matrix (Jacobian) . See also `MatCreateSNESMF()`
647e884886fSBarry Smith 
64811a5261eSBarry Smith   Collective
649e884886fSBarry Smith 
650e884886fSBarry Smith   Input Parameters:
651fef1beadSBarry Smith + comm - MPI communicator
6522ef1f0ffSBarry Smith . m    - number of local rows (or `PETSC_DECIDE` to have calculated if `M` is given)
653fef1beadSBarry Smith          This value should be the same as the local size used in creating the
654fef1beadSBarry Smith          y vector for the matrix-vector product y = Ax.
655fef1beadSBarry Smith . n    - This value should be the same as the local size used in creating the
65611a5261eSBarry Smith          x vector for the matrix-vector product y = Ax. (or `PETSC_DECIDE` to have
6572ef1f0ffSBarry Smith          calculated if `N` is given) For square matrices `n` is almost always `m`.
6582ef1f0ffSBarry Smith . M    - number of global rows (or `PETSC_DETERMINE` to have calculated if `m` is given)
6592ef1f0ffSBarry Smith - N    - number of global columns (or `PETSC_DETERMINE` to have calculated if `n` is given)
660fef1beadSBarry Smith 
661e884886fSBarry Smith   Output Parameter:
662e884886fSBarry Smith . J - the matrix-free matrix
663e884886fSBarry Smith 
66411a5261eSBarry Smith   Options Database Keys:
66511a5261eSBarry Smith + -mat_mffd_type             - wp or ds (see `MATMFFD_WP` or `MATMFFD_DS`)
666c92e3469SBarry Smith . -mat_mffd_err              - square root of estimated relative error in function evaluation
667c92e3469SBarry Smith . -mat_mffd_period           - how often h is recomputed, defaults to 1, every time
6682ef1f0ffSBarry Smith . -mat_mffd_check_positivity - possibly decrease `h` until U + h*a has only positive values
66911a5261eSBarry Smith . -mat_mffd_umin <umin>      - Sets umin (for default PETSc routine that computes h only)
670*0b4b7b1cSBarry Smith . -mat_mffd_complex          - use the Lyness trick with complex numbers to compute the matrix-vector product instead of differencing
671c92e3469SBarry Smith                                (requires real valued functions but that PETSc be configured for complex numbers)
672*0b4b7b1cSBarry Smith . -snes_mf                   - use the finite difference based matrix-free matrix with `SNESSolve()` and no preconditioner
673*0b4b7b1cSBarry Smith - -snes_mf_operator          - use the finite difference based matrix-free matrix with `SNESSolve()` but construct a preconditioner
674*0b4b7b1cSBarry Smith                                using the matrix passed as `pmat` to `SNESSetJacobian()`.
6759c6ac3b3SBarry Smith 
676e884886fSBarry Smith   Level: advanced
677e884886fSBarry Smith 
678e884886fSBarry Smith   Notes:
679*0b4b7b1cSBarry Smith   Use `MatMFFDSetFunction()` to provide the function that will be differenced to compute the matrix-vector product.
680*0b4b7b1cSBarry Smith 
681*0b4b7b1cSBarry Smith   The matrix-free matrix context contains the function pointers
682e884886fSBarry Smith   and work space for performing finite difference approximations of
683e884886fSBarry Smith   Jacobian-vector products, F'(u)*a,
684e884886fSBarry Smith 
685e884886fSBarry Smith   The default code uses the following approach to compute h
686e884886fSBarry Smith 
687e884886fSBarry Smith .vb
688e884886fSBarry Smith      F'(u)*a = [F(u+h*a) - F(u)]/h where
689e884886fSBarry Smith      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
690e884886fSBarry Smith        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
691e884886fSBarry Smith  where
692e884886fSBarry Smith      error_rel = square root of relative error in function evaluation
693e884886fSBarry Smith      umin = minimum iterate parameter
694e884886fSBarry Smith .ve
695e884886fSBarry Smith 
696*0b4b7b1cSBarry Smith   To have `SNES` use the matrix-free finite difference matrix-vector product and not provide a separate matrix
697*0b4b7b1cSBarry Smith   from which to compute the preconditioner (the `pmat` argument `SNESSetJacobian()`), then simply call `SNESSetJacobian()`
698*0b4b7b1cSBarry Smith   with `NULL` for the matrices and `MatMFFDComputeJacobian()`. Or use the options database option `-snes_mf`
699ca93e954SBarry Smith 
700*0b4b7b1cSBarry Smith   The user can set `error_rel` via `MatMFFDSetFunctionError()` and `umin` via `MatMFFDDSSetUmin()`.
701e884886fSBarry Smith 
702*0b4b7b1cSBarry Smith   Use `MATSHELL` or `MatCreateShell()` to provide your own custom matrix-vector operation.
703e884886fSBarry Smith 
7041cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatDestroy()`, `MatMFFDSetFunctionError()`, `MatMFFDDSSetUmin()`, `MatMFFDSetFunction()`
705*0b4b7b1cSBarry Smith           `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`, `MatCreateSNESMF()`, `MatCreateShell()`, `MATSHELL`,
706db781477SPatrick Sanan           `MatMFFDGetH()`, `MatMFFDRegister()`, `MatMFFDComputeJacobian()`
707e884886fSBarry Smith @*/
708d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateMFFD(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt M, PetscInt N, Mat *J)
709d71ae5a4SJacob Faibussowitsch {
710e884886fSBarry Smith   PetscFunctionBegin;
7119566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, J));
7129566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*J, m, n, M, N));
7139566063dSJacob Faibussowitsch   PetscCall(MatSetType(*J, MATMFFD));
7149566063dSJacob Faibussowitsch   PetscCall(MatSetUp(*J));
7153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
716e884886fSBarry Smith }
717e884886fSBarry Smith 
718e884886fSBarry Smith /*@
71911a5261eSBarry Smith   MatMFFDGetH - Gets the last value that was used as the differencing for a `MATMFFD` matrix
720e884886fSBarry Smith   parameter.
721e884886fSBarry Smith 
722e884886fSBarry Smith   Not Collective
723e884886fSBarry Smith 
724e884886fSBarry Smith   Input Parameters:
72511a5261eSBarry Smith . mat - the `MATMFFD` matrix
726e884886fSBarry Smith 
727fd292e60Sprj-   Output Parameter:
728e884886fSBarry Smith . h - the differencing step size
729e884886fSBarry Smith 
730e884886fSBarry Smith   Level: advanced
731e884886fSBarry Smith 
732fe59aa6dSJacob Faibussowitsch .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatCreateSNESMF()`, `MatMFFDSetHHistory()`, `MatCreateMFFD()`, `MatMFFDResetHHistory()`
733e884886fSBarry Smith @*/
734d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDGetH(Mat mat, PetscScalar *h)
735d71ae5a4SJacob Faibussowitsch {
736e884886fSBarry Smith   PetscFunctionBegin;
73788b4c220SStefano Zampini   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
7384f572ea9SToby Isaac   PetscAssertPointer(h, 2);
739cac4c232SBarry Smith   PetscUseMethod(mat, "MatMFFDGetH_C", (Mat, PetscScalar *), (mat, h));
7403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
741e884886fSBarry Smith }
742e884886fSBarry Smith 
743e884886fSBarry Smith /*@C
74401c1178eSBarry Smith   MatMFFDSetFunction - Sets the function used in applying the matrix-free `MATMFFD` matrix.
745e884886fSBarry Smith 
746c3339decSBarry Smith   Logically Collective
747e884886fSBarry Smith 
748e884886fSBarry Smith   Input Parameters:
74901c1178eSBarry Smith + mat     - the matrix-free matrix `MATMFFD` created via `MatCreateSNESMF()` or `MatCreateMFFD()`
750e884886fSBarry Smith . func    - the function to use
751e884886fSBarry Smith - funcctx - optional function context passed to function
752e884886fSBarry Smith 
7532920cce0SJacob Faibussowitsch   Calling sequence of `func`:
75414f633e2SBarry Smith + funcctx - user provided context
75514f633e2SBarry Smith . x       - input vector
75614f633e2SBarry Smith - f       - computed output function
75714f633e2SBarry Smith 
758e884886fSBarry Smith   Level: advanced
759e884886fSBarry Smith 
760e884886fSBarry Smith   Notes:
76101c1178eSBarry Smith   If you use this you MUST call `MatAssemblyBegin()` and `MatAssemblyEnd()` on the matrix-free
762e884886fSBarry Smith   matrix inside your compute Jacobian routine
763e884886fSBarry Smith 
76411a5261eSBarry Smith   If this is not set then it will use the function set with `SNESSetFunction()` if `MatCreateSNESMF()` was used.
765e884886fSBarry Smith 
766fe59aa6dSJacob Faibussowitsch .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatCreateSNESMF()`, `MatMFFDGetH()`, `MatCreateMFFD()`,
7672920cce0SJacob Faibussowitsch           `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`, `SNESSetFunction()`
768e884886fSBarry Smith @*/
7692920cce0SJacob Faibussowitsch PetscErrorCode MatMFFDSetFunction(Mat mat, PetscErrorCode (*func)(void *funcctx, Vec x, Vec f), void *funcctx)
770d71ae5a4SJacob Faibussowitsch {
771e884886fSBarry Smith   PetscFunctionBegin;
77288b4c220SStefano Zampini   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
773cac4c232SBarry Smith   PetscTryMethod(mat, "MatMFFDSetFunction_C", (Mat, PetscErrorCode (*)(void *, Vec, Vec), void *), (mat, func, funcctx));
7743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
775e884886fSBarry Smith }
776e884886fSBarry Smith 
777e884886fSBarry Smith /*@C
77811a5261eSBarry Smith   MatMFFDSetFunctioni - Sets the function for a single component for a `MATMFFD` matrix
779e884886fSBarry Smith 
780c3339decSBarry Smith   Logically Collective
781e884886fSBarry Smith 
782e884886fSBarry Smith   Input Parameters:
78301c1178eSBarry Smith + mat   - the matrix-free matrix `MATMFFD`
784e884886fSBarry Smith - funci - the function to use
785e884886fSBarry Smith 
786e884886fSBarry Smith   Level: advanced
787e884886fSBarry Smith 
788e884886fSBarry Smith   Notes:
78901c1178eSBarry Smith   If you use this you MUST call `MatAssemblyBegin()` and `MatAssemblyEnd()` on the matrix-free
79094eb4328SStefano Zampini   matrix inside your compute Jacobian routine.
79111a5261eSBarry Smith 
79294eb4328SStefano Zampini   This function is necessary to compute the diagonal of the matrix.
793486afcceSStefano Zampini   funci must not contain any MPI call as it is called inside a loop on the local portion of the vector.
794e884886fSBarry Smith 
7951cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatCreateSNESMF()`, `MatMFFDGetH()`, `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`, `SNESetFunction()`, `MatGetDiagonal()`
796e884886fSBarry Smith @*/
797d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetFunctioni(Mat mat, PetscErrorCode (*funci)(void *, PetscInt, Vec, PetscScalar *))
798d71ae5a4SJacob Faibussowitsch {
799e884886fSBarry Smith   PetscFunctionBegin;
8000700a824SBarry Smith   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
801cac4c232SBarry Smith   PetscTryMethod(mat, "MatMFFDSetFunctioni_C", (Mat, PetscErrorCode (*)(void *, PetscInt, Vec, PetscScalar *)), (mat, funci));
8023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
803e884886fSBarry Smith }
804e884886fSBarry Smith 
805e884886fSBarry Smith /*@C
80611a5261eSBarry Smith   MatMFFDSetFunctioniBase - Sets the base vector for a single component function evaluation for a `MATMFFD` matrix
807e884886fSBarry Smith 
808c3339decSBarry Smith   Logically Collective
809e884886fSBarry Smith 
810e884886fSBarry Smith   Input Parameters:
81101c1178eSBarry Smith + mat  - the `MATMFFD` matrix-free matrix
812e884886fSBarry Smith - func - the function to use
813e884886fSBarry Smith 
814e884886fSBarry Smith   Level: advanced
815e884886fSBarry Smith 
816e884886fSBarry Smith   Notes:
81701c1178eSBarry Smith   If you use this you MUST call `MatAssemblyBegin()` and `MatAssemblyEnd()` on the matrix-free
81894eb4328SStefano Zampini   matrix inside your compute Jacobian routine.
819e884886fSBarry Smith 
82011a5261eSBarry Smith   This function is necessary to compute the diagonal of the matrix, used for example with `PCJACOBI`
82111a5261eSBarry Smith 
822fe59aa6dSJacob Faibussowitsch .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatCreateSNESMF()`, `MatMFFDGetH()`, `MatCreateMFFD()`,
823db781477SPatrick Sanan           `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`, `SNESetFunction()`, `MatGetDiagonal()`
824e884886fSBarry Smith @*/
825d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetFunctioniBase(Mat mat, PetscErrorCode (*func)(void *, Vec))
826d71ae5a4SJacob Faibussowitsch {
827e884886fSBarry Smith   PetscFunctionBegin;
8280700a824SBarry Smith   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
829cac4c232SBarry Smith   PetscTryMethod(mat, "MatMFFDSetFunctioniBase_C", (Mat, PetscErrorCode (*)(void *, Vec)), (mat, func));
8303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
831e884886fSBarry Smith }
832e884886fSBarry Smith 
833e884886fSBarry Smith /*@
83411a5261eSBarry Smith   MatMFFDSetPeriod - Sets how often h is recomputed for a `MATMFFD` matrix, by default it is every time
835e884886fSBarry Smith 
836c3339decSBarry Smith   Logically Collective
837e884886fSBarry Smith 
838e884886fSBarry Smith   Input Parameters:
83901c1178eSBarry Smith + mat    - the `MATMFFD` matrix-free matrix
840e884886fSBarry Smith - period - 1 for every time, 2 for every second etc
841e884886fSBarry Smith 
8422ef1f0ffSBarry Smith   Options Database Key:
8432ef1f0ffSBarry Smith . -mat_mffd_period <period> - Sets how often `h` is recomputed
844e884886fSBarry Smith 
845e884886fSBarry Smith   Level: advanced
846e884886fSBarry Smith 
8471cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatCreateSNESMF()`, `MatMFFDGetH()`,
848db781477SPatrick Sanan           `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`
849e884886fSBarry Smith @*/
850d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetPeriod(Mat mat, PetscInt period)
851d71ae5a4SJacob Faibussowitsch {
852e884886fSBarry Smith   PetscFunctionBegin;
85388b4c220SStefano Zampini   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
85488b4c220SStefano Zampini   PetscValidLogicalCollectiveInt(mat, period, 2);
855cac4c232SBarry Smith   PetscTryMethod(mat, "MatMFFDSetPeriod_C", (Mat, PetscInt), (mat, period));
8563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
857e884886fSBarry Smith }
858e884886fSBarry Smith 
859e884886fSBarry Smith /*@
86011a5261eSBarry Smith   MatMFFDSetFunctionError - Sets the error_rel for the approximation of matrix-vector products using finite differences with the `MATMFFD` matrix
861e884886fSBarry Smith 
862c3339decSBarry Smith   Logically Collective
863e884886fSBarry Smith 
864e884886fSBarry Smith   Input Parameters:
86501c1178eSBarry Smith + mat   - the `MATMFFD` matrix-free matrix
866fe59aa6dSJacob Faibussowitsch - error - relative error (should be set to the square root of the relative error in the function evaluations)
867e884886fSBarry Smith 
8682ef1f0ffSBarry Smith   Options Database Key:
869a2b725a8SWilliam Gropp . -mat_mffd_err <error_rel> - Sets error_rel
870e884886fSBarry Smith 
871e884886fSBarry Smith   Level: advanced
872e884886fSBarry Smith 
87311a5261eSBarry Smith   Note:
874e884886fSBarry Smith   The default matrix-free matrix-vector product routine computes
875e884886fSBarry Smith .vb
876e884886fSBarry Smith      F'(u)*a = [F(u+h*a) - F(u)]/h where
877e884886fSBarry Smith      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
878e884886fSBarry Smith        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   else
879e884886fSBarry Smith .ve
880e884886fSBarry Smith 
881fe59aa6dSJacob Faibussowitsch .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatCreateSNESMF()`, `MatMFFDGetH()`, `MatCreateMFFD()`,
882db781477SPatrick Sanan           `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`
883e884886fSBarry Smith @*/
884d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetFunctionError(Mat mat, PetscReal error)
885d71ae5a4SJacob Faibussowitsch {
886e884886fSBarry Smith   PetscFunctionBegin;
88788b4c220SStefano Zampini   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
88888b4c220SStefano Zampini   PetscValidLogicalCollectiveReal(mat, error, 2);
889cac4c232SBarry Smith   PetscTryMethod(mat, "MatMFFDSetFunctionError_C", (Mat, PetscReal), (mat, error));
8903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
891e884886fSBarry Smith }
892e884886fSBarry Smith 
893e884886fSBarry Smith /*@
894e884886fSBarry Smith   MatMFFDSetHHistory - Sets an array to collect a history of the
89511a5261eSBarry Smith   differencing values (h) computed for the matrix-free product `MATMFFD` matrix
896e884886fSBarry Smith 
897c3339decSBarry Smith   Logically Collective
898e884886fSBarry Smith 
899e884886fSBarry Smith   Input Parameters:
90011a5261eSBarry Smith + J        - the `MATMFFD` matrix-free matrix
901a5b23f4aSJose E. Roman . history  - space to hold the history
902e884886fSBarry Smith - nhistory - number of entries in history, if more entries are generated than
903e884886fSBarry Smith               nhistory, then the later ones are discarded
904e884886fSBarry Smith 
905e884886fSBarry Smith   Level: advanced
906e884886fSBarry Smith 
90711a5261eSBarry Smith   Note:
90811a5261eSBarry Smith   Use `MatMFFDResetHHistory()` to reset the history counter and collect
909e884886fSBarry Smith   a new batch of differencing parameters, h.
910e884886fSBarry Smith 
9111cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatMFFDGetH()`, `MatCreateSNESMF()`,
912db781477SPatrick Sanan           `MatMFFDResetHHistory()`, `MatMFFDSetFunctionError()`
913e884886fSBarry Smith @*/
914d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetHHistory(Mat J, PetscScalar history[], PetscInt nhistory)
915d71ae5a4SJacob Faibussowitsch {
916e884886fSBarry Smith   PetscFunctionBegin;
91788b4c220SStefano Zampini   PetscValidHeaderSpecific(J, MAT_CLASSID, 1);
9184f572ea9SToby Isaac   if (history) PetscAssertPointer(history, 2);
91988b4c220SStefano Zampini   PetscValidLogicalCollectiveInt(J, nhistory, 3);
920cac4c232SBarry Smith   PetscUseMethod(J, "MatMFFDSetHHistory_C", (Mat, PetscScalar[], PetscInt), (J, history, nhistory));
9213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
922e884886fSBarry Smith }
923e884886fSBarry Smith 
924e884886fSBarry Smith /*@
925e884886fSBarry Smith   MatMFFDResetHHistory - Resets the counter to zero to begin
92611a5261eSBarry Smith   collecting a new set of differencing histories for the `MATMFFD` matrix
927e884886fSBarry Smith 
928c3339decSBarry Smith   Logically Collective
929e884886fSBarry Smith 
9302fe279fdSBarry Smith   Input Parameter:
931e884886fSBarry Smith . J - the matrix-free matrix context
932e884886fSBarry Smith 
933e884886fSBarry Smith   Level: advanced
934e884886fSBarry Smith 
93511a5261eSBarry Smith   Note:
93611a5261eSBarry Smith   Use `MatMFFDSetHHistory()` to create the original history counter.
937e884886fSBarry Smith 
9381cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatMFFDGetH()`, `MatCreateSNESMF()`,
939db781477SPatrick Sanan           `MatMFFDSetHHistory()`, `MatMFFDSetFunctionError()`
940e884886fSBarry Smith @*/
941d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDResetHHistory(Mat J)
942d71ae5a4SJacob Faibussowitsch {
943e884886fSBarry Smith   PetscFunctionBegin;
94488b4c220SStefano Zampini   PetscValidHeaderSpecific(J, MAT_CLASSID, 1);
945cac4c232SBarry Smith   PetscTryMethod(J, "MatMFFDResetHHistory_C", (Mat), (J));
9463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
947e884886fSBarry Smith }
948e884886fSBarry Smith 
949487a658cSBarry Smith /*@
9502ef1f0ffSBarry Smith   MatMFFDSetBase - Sets the vector `U` at which matrix vector products of the
95111a5261eSBarry Smith   Jacobian are computed for the `MATMFFD` matrix
952e884886fSBarry Smith 
953c3339decSBarry Smith   Logically Collective
954e884886fSBarry Smith 
955e884886fSBarry Smith   Input Parameters:
95611a5261eSBarry Smith + J - the `MATMFFD` matrix
9573ec795f1SBarry Smith . U - the vector
958bcddec3dSBarry Smith - F - (optional) vector that contains F(u) if it has been already computed
959e884886fSBarry Smith 
9602ef1f0ffSBarry Smith   Level: advanced
9612ef1f0ffSBarry Smith 
96295452b02SPatrick Sanan   Notes:
96395452b02SPatrick Sanan   This is rarely used directly
964e884886fSBarry Smith 
9652ef1f0ffSBarry Smith   If `F` is provided then it is not recomputed. Otherwise the function is evaluated at the base
96611a5261eSBarry Smith   point during the first `MatMult()` after each call to `MatMFFDSetBase()`.
967dff2f722SBarry Smith 
96842747ad1SJacob Faibussowitsch .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatMult()`
969e884886fSBarry Smith @*/
970d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetBase(Mat J, Vec U, Vec F)
971d71ae5a4SJacob Faibussowitsch {
972e884886fSBarry Smith   PetscFunctionBegin;
9730700a824SBarry Smith   PetscValidHeaderSpecific(J, MAT_CLASSID, 1);
9740700a824SBarry Smith   PetscValidHeaderSpecific(U, VEC_CLASSID, 2);
9750700a824SBarry Smith   if (F) PetscValidHeaderSpecific(F, VEC_CLASSID, 3);
976cac4c232SBarry Smith   PetscTryMethod(J, "MatMFFDSetBase_C", (Mat, Vec, Vec), (J, U, F));
9773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
978e884886fSBarry Smith }
979e884886fSBarry Smith 
980e884886fSBarry Smith /*@C
981e884886fSBarry Smith   MatMFFDSetCheckh - Sets a function that checks the computed h and adjusts
98211a5261eSBarry Smith   it to satisfy some criteria for the `MATMFFD` matrix
983e884886fSBarry Smith 
984c3339decSBarry Smith   Logically Collective
985e884886fSBarry Smith 
986e884886fSBarry Smith   Input Parameters:
98711a5261eSBarry Smith + J   - the `MATMFFD` matrix
9882ef1f0ffSBarry Smith . fun - the function that checks `h`
989e884886fSBarry Smith - ctx - any context needed by the function
990e884886fSBarry Smith 
991e884886fSBarry Smith   Options Database Keys:
99267b8a455SSatish Balay . -mat_mffd_check_positivity <bool> - Insure that U + h*a is non-negative
993e884886fSBarry Smith 
994e884886fSBarry Smith   Level: advanced
995e884886fSBarry Smith 
99695452b02SPatrick Sanan   Notes:
9972ef1f0ffSBarry Smith   For example, `MatMFFDCheckPositivity()` insures that all entries of U + h*a are non-negative
998e884886fSBarry Smith 
9992ef1f0ffSBarry Smith   The function you provide is called after the default `h` has been computed and allows you to
1000755b7f64SBarry Smith   modify it.
1001755b7f64SBarry Smith 
10021cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatMFFDCheckPositivity()`
1003e884886fSBarry Smith @*/
1004d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetCheckh(Mat J, PetscErrorCode (*fun)(void *, Vec, Vec, PetscScalar *), void *ctx)
1005d71ae5a4SJacob Faibussowitsch {
1006e884886fSBarry Smith   PetscFunctionBegin;
10070700a824SBarry Smith   PetscValidHeaderSpecific(J, MAT_CLASSID, 1);
1008cac4c232SBarry Smith   PetscTryMethod(J, "MatMFFDSetCheckh_C", (Mat, PetscErrorCode (*)(void *, Vec, Vec, PetscScalar *), void *), (J, fun, ctx));
10093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1010e884886fSBarry Smith }
1011e884886fSBarry Smith 
1012e884886fSBarry Smith /*@
1013e884886fSBarry Smith   MatMFFDCheckPositivity - Checks that all entries in U + h*a are positive or
101411a5261eSBarry Smith   zero, decreases h until this is satisfied for a `MATMFFD` matrix
1015e884886fSBarry Smith 
1016c3339decSBarry Smith   Logically Collective
1017e884886fSBarry Smith 
1018e884886fSBarry Smith   Input Parameters:
1019e884886fSBarry Smith + U     - base vector that is added to
1020e884886fSBarry Smith . a     - vector that is added
1021e884886fSBarry Smith . h     - scaling factor on a
1022e884886fSBarry Smith - dummy - context variable (unused)
1023e884886fSBarry Smith 
1024e884886fSBarry Smith   Options Database Keys:
102567b8a455SSatish Balay . -mat_mffd_check_positivity <bool> - Insure that U + h*a is nonnegative
1026e884886fSBarry Smith 
1027e884886fSBarry Smith   Level: advanced
1028e884886fSBarry Smith 
102911a5261eSBarry Smith   Note:
10302ef1f0ffSBarry Smith   This is rarely used directly, rather it is passed as an argument to `MatMFFDSetCheckh()`
1031e884886fSBarry Smith 
10321cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatMFFDSetCheckh()`
1033e884886fSBarry Smith @*/
1034d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDCheckPositivity(void *dummy, Vec U, Vec a, PetscScalar *h)
1035d71ae5a4SJacob Faibussowitsch {
1036e884886fSBarry Smith   PetscReal    val, minval;
1037e884886fSBarry Smith   PetscScalar *u_vec, *a_vec;
1038e884886fSBarry Smith   PetscInt     i, n;
1039e884886fSBarry Smith   MPI_Comm     comm;
1040e884886fSBarry Smith 
1041e884886fSBarry Smith   PetscFunctionBegin;
104288b4c220SStefano Zampini   PetscValidHeaderSpecific(U, VEC_CLASSID, 2);
104388b4c220SStefano Zampini   PetscValidHeaderSpecific(a, VEC_CLASSID, 3);
10444f572ea9SToby Isaac   PetscAssertPointer(h, 4);
10459566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)U, &comm));
10469566063dSJacob Faibussowitsch   PetscCall(VecGetArray(U, &u_vec));
10479566063dSJacob Faibussowitsch   PetscCall(VecGetArray(a, &a_vec));
10489566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(U, &n));
1049c068d9bbSLisandro Dalcin   minval = PetscAbsScalar(*h) * PetscRealConstant(1.01);
1050e884886fSBarry Smith   for (i = 0; i < n; i++) {
1051e884886fSBarry Smith     if (PetscRealPart(u_vec[i] + *h * a_vec[i]) <= 0.0) {
1052e884886fSBarry Smith       val = PetscAbsScalar(u_vec[i] / a_vec[i]);
1053e884886fSBarry Smith       if (val < minval) minval = val;
1054e884886fSBarry Smith     }
1055e884886fSBarry Smith   }
10569566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(U, &u_vec));
10579566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(a, &a_vec));
1058462c564dSBarry Smith   PetscCallMPI(MPIU_Allreduce(&minval, &val, 1, MPIU_REAL, MPIU_MIN, comm));
1059e884886fSBarry Smith   if (val <= PetscAbsScalar(*h)) {
10609566063dSJacob Faibussowitsch     PetscCall(PetscInfo(U, "Scaling back h from %g to %g\n", (double)PetscRealPart(*h), (double)(.99 * val)));
1061e884886fSBarry Smith     if (PetscRealPart(*h) > 0.0) *h = 0.99 * val;
1062e884886fSBarry Smith     else *h = -0.99 * val;
1063e884886fSBarry Smith   }
10643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1065e884886fSBarry Smith }
1066