xref: /petsc/src/mat/impls/mffd/mffd.c (revision 2ef1f0ff6e3530e8731eb06ad663081f5844f49f)
1e884886fSBarry Smith 
2af0996ceSBarry Smith #include <petsc/private/matimpl.h>
3c6db04a5SJed Brown #include <../src/mat/impls/mffd/mffdimpl.h> /*I  "petscmat.h"   I*/
4e884886fSBarry Smith 
5f4259b30SLisandro Dalcin PetscFunctionList MatMFFDList              = NULL;
6ace3abfcSBarry Smith PetscBool         MatMFFDRegisterAllCalled = PETSC_FALSE;
7e884886fSBarry Smith 
87087cfbeSBarry Smith PetscClassId  MATMFFD_CLASSID;
9166c7f25SBarry Smith PetscLogEvent MATMFFD_Mult;
10e884886fSBarry Smith 
11ace3abfcSBarry Smith static PetscBool MatMFFDPackageInitialized = PETSC_FALSE;
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 
18*2ef1f0ffSBarry Smith .seealso: [](chapter_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 
35*2ef1f0ffSBarry Smith .seealso: [](chapter_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);
77789d8953SBarry Smith   PetscValidCharPointer(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));
885f80ce2aSJacob Faibussowitsch   PetscCheck(r, PETSC_COMM_SELF, 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 
94e884886fSBarry Smith /*@C
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:
106*2ef1f0ffSBarry Smith     For example, such routines can compute `h` for use in
107e884886fSBarry Smith     Jacobian-vector products of the form
108*2ef1f0ffSBarry Smith .vb
109e884886fSBarry Smith 
110e884886fSBarry Smith                         F(x+ha) - F(x)
111e884886fSBarry Smith           F'(u)a  ~=  ----------------
112e884886fSBarry Smith                               h
113*2ef1f0ffSBarry Smith .ve
114e884886fSBarry Smith 
115*2ef1f0ffSBarry Smith .seealso: [](chapter_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);
121e884886fSBarry Smith   PetscValidCharPointer(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
17211a5261eSBarry Smith    MatMFFDRegister - Adds a method to the MATMFFD` registry.
1731c84c290SBarry Smith 
1741c84c290SBarry Smith    Not Collective
1751c84c290SBarry Smith 
1761c84c290SBarry Smith    Input Parameters:
1771c84c290SBarry Smith +  name_solver - name of a new user-defined compute-h module
1781c84c290SBarry Smith -  routine_create - 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 
1851c84c290SBarry Smith    Sample usage:
1861c84c290SBarry Smith .vb
187bdf89e91SBarry Smith    MatMFFDRegister("my_h",MyHCreate);
1881c84c290SBarry Smith .ve
1891c84c290SBarry Smith 
1901c84c290SBarry Smith    Then, your solver can be chosen with the procedural interface via
19111a5261eSBarry Smith $     `MatMFFDSetType`(mfctx,"my_h")
1921c84c290SBarry Smith    or at runtime via the option
193be7a6d03SBarry Smith $     -mat_mffd_type my_h
1941c84c290SBarry Smith 
195*2ef1f0ffSBarry Smith .seealso: [](chapter_matrices), `Mat`, `MATMFFD`, `MatMFFDRegisterAll()`, `MatMFFDRegisterDestroy()`
1961c84c290SBarry Smith  @*/
197d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDRegister(const char sname[], PetscErrorCode (*function)(MatMFFD))
198d71ae5a4SJacob Faibussowitsch {
199e884886fSBarry Smith   PetscFunctionBegin;
2009566063dSJacob Faibussowitsch   PetscCall(MatInitializePackage());
2019566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&MatMFFDList, sname, function));
2023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
203e884886fSBarry Smith }
204e884886fSBarry Smith 
205d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDestroy_MFFD(Mat mat)
206d71ae5a4SJacob Faibussowitsch {
207789d8953SBarry Smith   MatMFFD ctx;
208e884886fSBarry Smith 
209e884886fSBarry Smith   PetscFunctionBegin;
2109566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
2119566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ctx->w));
2129566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ctx->current_u));
21348a46eb9SPierre Jolivet   if (ctx->current_f_allocated) PetscCall(VecDestroy(&ctx->current_f));
214dbbe0bcdSBarry Smith   PetscTryTypeMethod(ctx, destroy);
2159566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(&ctx));
216e884886fSBarry Smith 
2179566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetBase_C", NULL));
2189566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetFunctioniBase_C", NULL));
2199566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetFunctioni_C", NULL));
2209566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetFunction_C", NULL));
2219566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetFunctionError_C", NULL));
2229566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetCheckh_C", NULL));
2239566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetPeriod_C", NULL));
2249566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDResetHHistory_C", NULL));
2259566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetHHistory_C", NULL));
2269566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetType_C", NULL));
2279566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDGetH_C", NULL));
2282e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatSNESMFSetReuseBase_C", NULL));
2292e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatSNESMFGetReuseBase_C", NULL));
2303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
231e884886fSBarry Smith }
232e884886fSBarry Smith 
233e884886fSBarry Smith /*
234e884886fSBarry Smith    MatMFFDView_MFFD - Views matrix-free parameters.
235e884886fSBarry Smith 
236e884886fSBarry Smith */
237d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_MFFD(Mat J, PetscViewer viewer)
238d71ae5a4SJacob Faibussowitsch {
239789d8953SBarry Smith   MatMFFD     ctx;
240a283635eSDmitry Karpeev   PetscBool   iascii, viewbase, viewfunction;
241a283635eSDmitry Karpeev   const char *prefix;
242e884886fSBarry Smith 
243e884886fSBarry Smith   PetscFunctionBegin;
2449566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J, &ctx));
2459566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
246e884886fSBarry Smith   if (iascii) {
2479566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Matrix-free approximation:\n"));
2489566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
2499566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "err=%g (relative error in function evaluation)\n", (double)ctx->error_rel));
2507adad957SLisandro Dalcin     if (!((PetscObject)ctx)->type_name) {
2519566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "The compute h routine has not yet been set\n"));
252e884886fSBarry Smith     } else {
2539566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Using %s compute h routine\n", ((PetscObject)ctx)->type_name));
254e884886fSBarry Smith     }
255c92e3469SBarry Smith #if defined(PETSC_USE_COMPLEX)
25648a46eb9SPierre Jolivet     if (ctx->usecomplex) PetscCall(PetscViewerASCIIPrintf(viewer, "Using Lyness complex number trick to compute the matrix-vector product\n"));
257c92e3469SBarry Smith #endif
258dbbe0bcdSBarry Smith     PetscTryTypeMethod(ctx, view, viewer);
2599566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetOptionsPrefix((PetscObject)J, &prefix));
260a283635eSDmitry Karpeev 
2619566063dSJacob Faibussowitsch     PetscCall(PetscOptionsHasName(((PetscObject)J)->options, prefix, "-mat_mffd_view_base", &viewbase));
262a283635eSDmitry Karpeev     if (viewbase) {
2639566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Base:\n"));
2649566063dSJacob Faibussowitsch       PetscCall(VecView(ctx->current_u, viewer));
265a283635eSDmitry Karpeev     }
2669566063dSJacob Faibussowitsch     PetscCall(PetscOptionsHasName(((PetscObject)J)->options, prefix, "-mat_mffd_view_function", &viewfunction));
267a283635eSDmitry Karpeev     if (viewfunction) {
2689566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Function:\n"));
2699566063dSJacob Faibussowitsch       PetscCall(VecView(ctx->current_f, viewer));
270a283635eSDmitry Karpeev     }
2719566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
27211aeaf0aSBarry Smith   }
2733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
274e884886fSBarry Smith }
275e884886fSBarry Smith 
276e884886fSBarry Smith /*
277e884886fSBarry Smith    MatAssemblyEnd_MFFD - Resets the ctx->ncurrenth to zero. This
278e884886fSBarry Smith    allows the user to indicate the beginning of a new linear solve by calling
279e884886fSBarry Smith    MatAssemblyXXX() on the matrix free matrix. This then allows the
2801d0fab5eSBarry Smith    MatCreateMFFD_WP() to properly compute ||U|| only the first time
281e884886fSBarry Smith    in the linear solver rather than every time.
2825a576424SJed Brown 
283cc2e6a90SBarry Smith    This function is referenced directly from MatAssemblyEnd_SNESMF(), which may be in a different shared library hence
284cc2e6a90SBarry Smith    it must be labeled as PETSC_EXTERN
285e884886fSBarry Smith */
286d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatAssemblyEnd_MFFD(Mat J, MatAssemblyType mt)
287d71ae5a4SJacob Faibussowitsch {
288789d8953SBarry Smith   MatMFFD j;
289e884886fSBarry Smith 
290e884886fSBarry Smith   PetscFunctionBegin;
2919566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J, &j));
2929566063dSJacob Faibussowitsch   PetscCall(MatMFFDResetHHistory(J));
2933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
294e884886fSBarry Smith }
295e884886fSBarry Smith 
296e884886fSBarry Smith /*
297e884886fSBarry Smith   MatMult_MFFD - Default matrix-free form for Jacobian-vector product, y = F'(u)*a:
298e884886fSBarry Smith 
299e884886fSBarry Smith         y ~= (F(u + ha) - F(u))/h,
300e884886fSBarry Smith   where F = nonlinear function, as set by SNESSetFunction()
301e884886fSBarry Smith         u = current iterate
302e884886fSBarry Smith         h = difference interval
303e884886fSBarry Smith */
304d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMult_MFFD(Mat mat, Vec a, Vec y)
305d71ae5a4SJacob Faibussowitsch {
306789d8953SBarry Smith   MatMFFD     ctx;
307e884886fSBarry Smith   PetscScalar h;
308e884886fSBarry Smith   Vec         w, U, F;
309ace3abfcSBarry Smith   PetscBool   zeroa;
310e884886fSBarry Smith 
311e884886fSBarry Smith   PetscFunctionBegin;
3129566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
31328b400f6SJacob Faibussowitsch   PetscCheck(ctx->current_u, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "MatMFFDSetBase() has not been called, this is often caused by forgetting to call \n\t\tMatAssemblyBegin/End on the first Mat in the SNES compute function");
314e884886fSBarry Smith   /* We log matrix-free matrix-vector products separately, so that we can
315e884886fSBarry Smith      separate the performance monitoring from the cases that use conventional
316e884886fSBarry Smith      storage.  We may eventually modify event logging to associate events
317e884886fSBarry Smith      with particular objects, hence alleviating the more general problem. */
3189566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(MATMFFD_Mult, a, y, 0, 0));
319e884886fSBarry Smith 
320e884886fSBarry Smith   w = ctx->w;
321e884886fSBarry Smith   U = ctx->current_u;
3223ec795f1SBarry Smith   F = ctx->current_f;
323e884886fSBarry Smith   /*
324e884886fSBarry Smith       Compute differencing parameter
325e884886fSBarry Smith   */
3264a2f8832SBarry Smith   if (!((PetscObject)ctx)->type_name) {
3279566063dSJacob Faibussowitsch     PetscCall(MatMFFDSetType(mat, MATMFFD_WP));
3289566063dSJacob Faibussowitsch     PetscCall(MatSetFromOptions(mat));
329e884886fSBarry Smith   }
330dbbe0bcdSBarry Smith   PetscUseTypeMethod(ctx, compute, U, a, &h, &zeroa);
331e884886fSBarry Smith   if (zeroa) {
3329566063dSJacob Faibussowitsch     PetscCall(VecSet(y, 0.0));
33385d14658SLisandro Dalcin     PetscCall(PetscLogEventEnd(MATMFFD_Mult, a, y, 0, 0));
3343ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
335e884886fSBarry Smith   }
336e884886fSBarry Smith 
33708401ef6SPierre Jolivet   PetscCheck(!mat->erroriffailure || !PetscIsInfOrNanScalar(h), PETSC_COMM_SELF, PETSC_ERR_PLIB, "Computed Nan differencing parameter h");
33848a46eb9SPierre Jolivet   if (ctx->checkh) PetscCall((*ctx->checkh)(ctx->checkhctx, U, a, &h));
339e884886fSBarry Smith 
340e884886fSBarry Smith   /* keep a record of the current differencing parameter h */
341e884886fSBarry Smith   ctx->currenth = h;
342e884886fSBarry Smith #if defined(PETSC_USE_COMPLEX)
3439566063dSJacob Faibussowitsch   PetscCall(PetscInfo(mat, "Current differencing parameter: %g + %g i\n", (double)PetscRealPart(h), (double)PetscImaginaryPart(h)));
344e884886fSBarry Smith #else
3459566063dSJacob Faibussowitsch   PetscCall(PetscInfo(mat, "Current differencing parameter: %15.12e\n", (double)PetscRealPart(h)));
346e884886fSBarry Smith #endif
347ad540459SPierre Jolivet   if (ctx->historyh && ctx->ncurrenth < ctx->maxcurrenth) ctx->historyh[ctx->ncurrenth] = h;
348e884886fSBarry Smith   ctx->ncurrenth++;
349e884886fSBarry Smith 
350c92e3469SBarry Smith #if defined(PETSC_USE_COMPLEX)
351c92e3469SBarry Smith   if (ctx->usecomplex) h = PETSC_i * h;
352c92e3469SBarry Smith #endif
353c92e3469SBarry Smith 
354e884886fSBarry Smith   /* w = u + ha */
3559566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(w, h, a, U));
356e884886fSBarry Smith 
3571b797266SDmitry Karpeev   /* compute func(U) as base for differencing; only needed first time in and not when provided by user */
35848a46eb9SPierre Jolivet   if (ctx->ncurrenth == 1 && ctx->current_f_allocated) PetscCall((*ctx->func)(ctx->funcctx, U, F));
3599566063dSJacob Faibussowitsch   PetscCall((*ctx->func)(ctx->funcctx, w, y));
360e884886fSBarry Smith 
361c92e3469SBarry Smith #if defined(PETSC_USE_COMPLEX)
362c92e3469SBarry Smith   if (ctx->usecomplex) {
3639566063dSJacob Faibussowitsch     PetscCall(VecImaginaryPart(y));
364c92e3469SBarry Smith     h = PetscImaginaryPart(h);
365c92e3469SBarry Smith   } else {
3669566063dSJacob Faibussowitsch     PetscCall(VecAXPY(y, -1.0, F));
367c92e3469SBarry Smith   }
368c92e3469SBarry Smith #else
3699566063dSJacob Faibussowitsch   PetscCall(VecAXPY(y, -1.0, F));
370c92e3469SBarry Smith #endif
3719566063dSJacob Faibussowitsch   PetscCall(VecScale(y, 1.0 / h));
3729566063dSJacob Faibussowitsch   if (mat->nullsp) PetscCall(MatNullSpaceRemove(mat->nullsp, y));
373e884886fSBarry Smith 
3749566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(MATMFFD_Mult, a, y, 0, 0));
3753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
376e884886fSBarry Smith }
377e884886fSBarry Smith 
378e884886fSBarry Smith /*
379e884886fSBarry Smith   MatGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix
380e884886fSBarry Smith 
381e884886fSBarry Smith         y ~= (F(u + ha) - F(u))/h,
382e884886fSBarry Smith   where F = nonlinear function, as set by SNESSetFunction()
383e884886fSBarry Smith         u = current iterate
384e884886fSBarry Smith         h = difference interval
385e884886fSBarry Smith */
386d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetDiagonal_MFFD(Mat mat, Vec a)
387d71ae5a4SJacob Faibussowitsch {
388789d8953SBarry Smith   MatMFFD     ctx;
389e884886fSBarry Smith   PetscScalar h, *aa, *ww, v;
390e884886fSBarry Smith   PetscReal   epsilon = PETSC_SQRT_MACHINE_EPSILON, umin = 100.0 * PETSC_SQRT_MACHINE_EPSILON;
391e884886fSBarry Smith   Vec         w, U;
392e884886fSBarry Smith   PetscInt    i, rstart, rend;
393e884886fSBarry Smith 
394e884886fSBarry Smith   PetscFunctionBegin;
3959566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
39628b400f6SJacob Faibussowitsch   PetscCheck(ctx->func, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Requires calling MatMFFDSetFunction() first");
39728b400f6SJacob Faibussowitsch   PetscCheck(ctx->funci, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Requires calling MatMFFDSetFunctioni() first");
398e884886fSBarry Smith   w = ctx->w;
399e884886fSBarry Smith   U = ctx->current_u;
4009566063dSJacob Faibussowitsch   PetscCall((*ctx->func)(ctx->funcctx, U, a));
4011baa6e33SBarry Smith   if (ctx->funcisetbase) PetscCall((*ctx->funcisetbase)(ctx->funcctx, U));
4029566063dSJacob Faibussowitsch   PetscCall(VecCopy(U, w));
403e884886fSBarry Smith 
4049566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(a, &rstart, &rend));
4059566063dSJacob Faibussowitsch   PetscCall(VecGetArray(a, &aa));
406e884886fSBarry Smith   for (i = rstart; i < rend; i++) {
4079566063dSJacob Faibussowitsch     PetscCall(VecGetArray(w, &ww));
408e884886fSBarry Smith     h = ww[i - rstart];
409e884886fSBarry Smith     if (h == 0.0) h = 1.0;
410e884886fSBarry Smith     if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0) h = umin;
411e884886fSBarry Smith     else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin;
412e884886fSBarry Smith     h *= epsilon;
413e884886fSBarry Smith 
414e884886fSBarry Smith     ww[i - rstart] += h;
4159566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(w, &ww));
4169566063dSJacob Faibussowitsch     PetscCall((*ctx->funci)(ctx->funcctx, i, w, &v));
417e884886fSBarry Smith     aa[i - rstart] = (v - aa[i - rstart]) / h;
418e884886fSBarry Smith 
4199566063dSJacob Faibussowitsch     PetscCall(VecGetArray(w, &ww));
420e884886fSBarry Smith     ww[i - rstart] -= h;
4219566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(w, &ww));
422e884886fSBarry Smith   }
4239566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(a, &aa));
4243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
425e884886fSBarry Smith }
426e884886fSBarry Smith 
427d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatMFFDSetBase_MFFD(Mat J, Vec U, Vec F)
428d71ae5a4SJacob Faibussowitsch {
429789d8953SBarry Smith   MatMFFD ctx;
430e884886fSBarry Smith 
431e884886fSBarry Smith   PetscFunctionBegin;
4329566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J, &ctx));
4339566063dSJacob Faibussowitsch   PetscCall(MatMFFDResetHHistory(J));
434c51ad4d4SStefano Zampini   if (!ctx->current_u) {
4359566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(U, &ctx->current_u));
4369566063dSJacob Faibussowitsch     PetscCall(VecLockReadPush(ctx->current_u));
437c51ad4d4SStefano Zampini   }
4389566063dSJacob Faibussowitsch   PetscCall(VecLockReadPop(ctx->current_u));
4399566063dSJacob Faibussowitsch   PetscCall(VecCopy(U, ctx->current_u));
4409566063dSJacob Faibussowitsch   PetscCall(VecLockReadPush(ctx->current_u));
44152121784SBarry Smith   if (F) {
4429566063dSJacob Faibussowitsch     if (ctx->current_f_allocated) PetscCall(VecDestroy(&ctx->current_f));
4433ec795f1SBarry Smith     ctx->current_f           = F;
444cfe22f5eSBarry Smith     ctx->current_f_allocated = PETSC_FALSE;
445cfe22f5eSBarry Smith   } else if (!ctx->current_f_allocated) {
4469566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(J, NULL, &ctx->current_f));
447cfe22f5eSBarry Smith     ctx->current_f_allocated = PETSC_TRUE;
44852121784SBarry Smith   }
44948a46eb9SPierre Jolivet   if (!ctx->w) PetscCall(VecDuplicate(ctx->current_u, &ctx->w));
450e884886fSBarry Smith   J->assembled = PETSC_TRUE;
4513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
452e884886fSBarry Smith }
4535a576424SJed Brown 
454e884886fSBarry Smith typedef PetscErrorCode (*FCN3)(void *, Vec, Vec, PetscScalar *); /* force argument to next function to not be extern C*/
455b2573a8aSBarry Smith 
456d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDSetCheckh_MFFD(Mat J, FCN3 fun, void *ectx)
457d71ae5a4SJacob Faibussowitsch {
458789d8953SBarry Smith   MatMFFD ctx;
459e884886fSBarry Smith 
460e884886fSBarry Smith   PetscFunctionBegin;
4619566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J, &ctx));
462e884886fSBarry Smith   ctx->checkh    = fun;
463e884886fSBarry Smith   ctx->checkhctx = ectx;
4643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
465e884886fSBarry Smith }
466e884886fSBarry Smith 
4676aa9148fSLisandro Dalcin /*@C
4686aa9148fSLisandro Dalcin    MatMFFDSetOptionsPrefix - Sets the prefix used for searching for all
46911a5261eSBarry Smith    MATMFFD` options in the database.
4706aa9148fSLisandro Dalcin 
471c3339decSBarry Smith    Collective
4726aa9148fSLisandro Dalcin 
473d8d19677SJose E. Roman    Input Parameters:
47411a5261eSBarry Smith +  A - the `MATMFFD` context
4756aa9148fSLisandro Dalcin -  prefix - the prefix to prepend to all option names
4766aa9148fSLisandro Dalcin 
47711a5261eSBarry Smith    Note:
4786aa9148fSLisandro Dalcin    A hyphen (-) must NOT be given at the beginning of the prefix name.
4796aa9148fSLisandro Dalcin    The first character of all runtime options is AUTOMATICALLY the hyphen.
4806aa9148fSLisandro Dalcin 
4816aa9148fSLisandro Dalcin    Level: advanced
4826aa9148fSLisandro Dalcin 
483*2ef1f0ffSBarry Smith .seealso: [](chapter_matrices), `Mat`, `MATMFFD`, `MatSetFromOptions()`, `MatCreateSNESMF()`, `MatCreateMFFD()`
4846aa9148fSLisandro Dalcin @*/
485d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetOptionsPrefix(Mat mat, const char prefix[])
486d71ae5a4SJacob Faibussowitsch {
487789d8953SBarry Smith   MatMFFD mfctx;
4885fd66863SKarl Rupp 
4896aa9148fSLisandro Dalcin   PetscFunctionBegin;
4900700a824SBarry Smith   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
4919566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &mfctx));
4920700a824SBarry Smith   PetscValidHeaderSpecific(mfctx, MATMFFD_CLASSID, 1);
4939566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)mfctx, prefix));
4943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4956aa9148fSLisandro Dalcin }
4966aa9148fSLisandro Dalcin 
497d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSetFromOptions_MFFD(Mat mat, PetscOptionItems *PetscOptionsObject)
498d71ae5a4SJacob Faibussowitsch {
499789d8953SBarry Smith   MatMFFD   mfctx;
500ace3abfcSBarry Smith   PetscBool flg;
501e884886fSBarry Smith   char      ftype[256];
502e884886fSBarry Smith 
503e884886fSBarry Smith   PetscFunctionBegin;
5049566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &mfctx));
505dbbe0bcdSBarry Smith   PetscValidHeaderSpecific(mfctx, MATMFFD_CLASSID, 1);
506d0609cedSBarry Smith   PetscObjectOptionsBegin((PetscObject)mfctx);
5079566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFList("-mat_mffd_type", "Matrix free type", "MatMFFDSetType", MatMFFDList, ((PetscObject)mfctx)->type_name, ftype, 256, &flg));
5081baa6e33SBarry Smith   if (flg) PetscCall(MatMFFDSetType(mat, ftype));
509e884886fSBarry Smith 
5109566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_mffd_err", "set sqrt relative error in function", "MatMFFDSetFunctionError", mfctx->error_rel, &mfctx->error_rel, NULL));
5119566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-mat_mffd_period", "how often h is recomputed", "MatMFFDSetPeriod", mfctx->recomputeperiod, &mfctx->recomputeperiod, NULL));
512e884886fSBarry Smith 
51390d69ab7SBarry Smith   flg = PETSC_FALSE;
5149566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-mat_mffd_check_positivity", "Insure that U + h*a is nonnegative", "MatMFFDSetCheckh", flg, &flg, NULL));
5151baa6e33SBarry Smith   if (flg) PetscCall(MatMFFDSetCheckh(mat, MatMFFDCheckPositivity, NULL));
516c92e3469SBarry Smith #if defined(PETSC_USE_COMPLEX)
5179566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-mat_mffd_complex", "Use Lyness complex number trick to compute the matrix-vector product", "None", mfctx->usecomplex, &mfctx->usecomplex, NULL));
518c92e3469SBarry Smith #endif
519dbbe0bcdSBarry Smith   PetscTryTypeMethod(mfctx, setfromoptions, PetscOptionsObject);
520d0609cedSBarry Smith   PetscOptionsEnd();
5213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
522e884886fSBarry Smith }
523e884886fSBarry Smith 
524d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDSetPeriod_MFFD(Mat mat, PetscInt period)
525d71ae5a4SJacob Faibussowitsch {
526789d8953SBarry Smith   MatMFFD ctx;
527bc0f08ceSBarry Smith 
528bc0f08ceSBarry Smith   PetscFunctionBegin;
5299566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
530bc0f08ceSBarry Smith   ctx->recomputeperiod = period;
5313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
532bc0f08ceSBarry Smith }
533bc0f08ceSBarry Smith 
534d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDSetFunction_MFFD(Mat mat, PetscErrorCode (*func)(void *, Vec, Vec), void *funcctx)
535d71ae5a4SJacob Faibussowitsch {
536789d8953SBarry Smith   MatMFFD ctx;
537bc0f08ceSBarry Smith 
538bc0f08ceSBarry Smith   PetscFunctionBegin;
5399566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
540bc0f08ceSBarry Smith   ctx->func    = func;
541bc0f08ceSBarry Smith   ctx->funcctx = funcctx;
5423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
543bc0f08ceSBarry Smith }
544bc0f08ceSBarry Smith 
545d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDSetFunctionError_MFFD(Mat mat, PetscReal error)
546d71ae5a4SJacob Faibussowitsch {
547789d8953SBarry Smith   MatMFFD ctx;
548bc0f08ceSBarry Smith 
549bc0f08ceSBarry Smith   PetscFunctionBegin;
5509566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
551bc0f08ceSBarry Smith   if (error != PETSC_DEFAULT) ctx->error_rel = error;
5523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
553bc0f08ceSBarry Smith }
554bc0f08ceSBarry Smith 
555d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetHHistory_MFFD(Mat J, PetscScalar history[], PetscInt nhistory)
556d71ae5a4SJacob Faibussowitsch {
557789d8953SBarry Smith   MatMFFD ctx;
558789d8953SBarry Smith 
5593b49f96aSBarry Smith   PetscFunctionBegin;
5609566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J, &ctx));
561789d8953SBarry Smith   ctx->historyh    = history;
562789d8953SBarry Smith   ctx->maxcurrenth = nhistory;
563789d8953SBarry Smith   ctx->currenth    = 0.;
5643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5653b49f96aSBarry Smith }
5663b49f96aSBarry Smith 
567e884886fSBarry Smith /*MC
568*2ef1f0ffSBarry Smith   MATMFFD - "mffd" - A matrix free matrix type.
569e884886fSBarry Smith 
570e884886fSBarry Smith   Level: advanced
571e884886fSBarry Smith 
572*2ef1f0ffSBarry Smith   Developers Note:
573*2ef1f0ffSBarry Smith   This is implemented on top of `MATSHELL` to get support for scaling and shifting without requiring duplicate code
574789d8953SBarry Smith 
575*2ef1f0ffSBarry Smith .seealso: [](chapter_matrices), `Mat`, `MatCreateMFFD()`, `MatCreateSNESMF()`, `MatMFFDSetFunction()`, `MatMFFDSetType()`,
576db781477SPatrick Sanan           `MatMFFDSetFunctionError()`, `MatMFFDDSSetUmin()`, `MatMFFDSetFunction()`
577db781477SPatrick Sanan           `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`, `MatCreateSNESMF()`,
578db781477SPatrick Sanan           `MatMFFDGetH()`,
579e884886fSBarry Smith M*/
580d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatCreate_MFFD(Mat A)
581d71ae5a4SJacob Faibussowitsch {
582e884886fSBarry Smith   MatMFFD mfctx;
583e884886fSBarry Smith 
584e884886fSBarry Smith   PetscFunctionBegin;
5859566063dSJacob Faibussowitsch   PetscCall(MatMFFDInitializePackage());
586e884886fSBarry Smith 
5879566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(mfctx, MATMFFD_CLASSID, "MatMFFD", "Matrix-free Finite Differencing", "Mat", PetscObjectComm((PetscObject)A), NULL, NULL));
5882205254eSKarl Rupp 
589e884886fSBarry Smith   mfctx->error_rel                = PETSC_SQRT_MACHINE_EPSILON;
590e884886fSBarry Smith   mfctx->recomputeperiod          = 1;
591e884886fSBarry Smith   mfctx->count                    = 0;
592e884886fSBarry Smith   mfctx->currenth                 = 0.0;
5930298fd71SBarry Smith   mfctx->historyh                 = NULL;
594e884886fSBarry Smith   mfctx->ncurrenth                = 0;
595e884886fSBarry Smith   mfctx->maxcurrenth              = 0;
596f4259b30SLisandro Dalcin   ((PetscObject)mfctx)->type_name = NULL;
597e884886fSBarry Smith 
598e884886fSBarry Smith   /*
599e884886fSBarry Smith      Create the empty data structure to contain compute-h routines.
600e884886fSBarry Smith      These will be filled in below from the command line options or
601e884886fSBarry Smith      a later call with MatMFFDSetType() or if that is not called
602e884886fSBarry Smith      then it will default in the first use of MatMult_MFFD()
603e884886fSBarry Smith   */
604f4259b30SLisandro Dalcin   mfctx->ops->compute        = NULL;
605f4259b30SLisandro Dalcin   mfctx->ops->destroy        = NULL;
606f4259b30SLisandro Dalcin   mfctx->ops->view           = NULL;
607f4259b30SLisandro Dalcin   mfctx->ops->setfromoptions = NULL;
608f4259b30SLisandro Dalcin   mfctx->hctx                = NULL;
609e884886fSBarry Smith 
610f4259b30SLisandro Dalcin   mfctx->func    = NULL;
611f4259b30SLisandro Dalcin   mfctx->funcctx = NULL;
6120298fd71SBarry Smith   mfctx->w       = NULL;
613789d8953SBarry Smith   mfctx->mat     = A;
614e884886fSBarry Smith 
6159566063dSJacob Faibussowitsch   PetscCall(MatSetType(A, MATSHELL));
6169566063dSJacob Faibussowitsch   PetscCall(MatShellSetContext(A, mfctx));
6179566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(A, MATOP_MULT, (void (*)(void))MatMult_MFFD));
6189566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(A, MATOP_DESTROY, (void (*)(void))MatDestroy_MFFD));
6199566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(A, MATOP_VIEW, (void (*)(void))MatView_MFFD));
6209566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(A, MATOP_ASSEMBLY_END, (void (*)(void))MatAssemblyEnd_MFFD));
6219566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(A, MATOP_SET_FROM_OPTIONS, (void (*)(void))MatSetFromOptions_MFFD));
622e884886fSBarry Smith 
6239566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetBase_C", MatMFFDSetBase_MFFD));
6249566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetFunctioniBase_C", MatMFFDSetFunctioniBase_MFFD));
6259566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetFunctioni_C", MatMFFDSetFunctioni_MFFD));
6269566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetFunction_C", MatMFFDSetFunction_MFFD));
6279566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetCheckh_C", MatMFFDSetCheckh_MFFD));
6289566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetPeriod_C", MatMFFDSetPeriod_MFFD));
6299566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetFunctionError_C", MatMFFDSetFunctionError_MFFD));
6309566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDResetHHistory_C", MatMFFDResetHHistory_MFFD));
6319566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetHHistory_C", MatMFFDSetHHistory_MFFD));
6329566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetType_C", MatMFFDSetType_MFFD));
6339566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDGetH_C", MatMFFDGetH_MFFD));
6349566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATMFFD));
6353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
636e884886fSBarry Smith }
637e884886fSBarry Smith 
638e884886fSBarry Smith /*@
63911a5261eSBarry Smith    MatCreateMFFD - Creates a matrix-free matrix of type `MATMFFD`. See also `MatCreateSNESMF()`
640e884886fSBarry Smith 
64111a5261eSBarry Smith    Collective
642e884886fSBarry Smith 
643e884886fSBarry Smith    Input Parameters:
644fef1beadSBarry Smith +  comm - MPI communicator
645*2ef1f0ffSBarry Smith .  m - number of local rows (or `PETSC_DECIDE` to have calculated if `M` is given)
646fef1beadSBarry Smith            This value should be the same as the local size used in creating the
647fef1beadSBarry Smith            y vector for the matrix-vector product y = Ax.
648fef1beadSBarry Smith .  n - This value should be the same as the local size used in creating the
64911a5261eSBarry Smith        x vector for the matrix-vector product y = Ax. (or `PETSC_DECIDE` to have
650*2ef1f0ffSBarry Smith        calculated if `N` is given) For square matrices `n` is almost always `m`.
651*2ef1f0ffSBarry Smith .  M - number of global rows (or `PETSC_DETERMINE` to have calculated if `m` is given)
652*2ef1f0ffSBarry Smith -  N - number of global columns (or `PETSC_DETERMINE` to have calculated if `n` is given)
653fef1beadSBarry Smith 
654e884886fSBarry Smith    Output Parameter:
655e884886fSBarry Smith .  J - the matrix-free matrix
656e884886fSBarry Smith 
65711a5261eSBarry Smith    Options Database Keys:
65811a5261eSBarry Smith +  -mat_mffd_type - wp or ds (see `MATMFFD_WP` or `MATMFFD_DS`)
659c92e3469SBarry Smith .  -mat_mffd_err - square root of estimated relative error in function evaluation
660c92e3469SBarry Smith .  -mat_mffd_period - how often h is recomputed, defaults to 1, every time
661*2ef1f0ffSBarry Smith .  -mat_mffd_check_positivity - possibly decrease `h` until U + h*a has only positive values
66211a5261eSBarry Smith .  -mat_mffd_umin <umin> - Sets umin (for default PETSc routine that computes h only)
663c92e3469SBarry Smith -  -mat_mffd_complex - use the Lyness trick with complex numbers to compute the matrix-vector product instead of differencing
664c92e3469SBarry Smith                        (requires real valued functions but that PETSc be configured for complex numbers)
6659c6ac3b3SBarry Smith 
666e884886fSBarry Smith    Level: advanced
667e884886fSBarry Smith 
668e884886fSBarry Smith    Notes:
669e884886fSBarry Smith    The matrix-free matrix context merely contains the function pointers
670e884886fSBarry Smith    and work space for performing finite difference approximations of
671e884886fSBarry Smith    Jacobian-vector products, F'(u)*a,
672e884886fSBarry Smith 
673e884886fSBarry Smith    The default code uses the following approach to compute h
674e884886fSBarry Smith 
675e884886fSBarry Smith .vb
676e884886fSBarry Smith      F'(u)*a = [F(u+h*a) - F(u)]/h where
677e884886fSBarry Smith      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
678e884886fSBarry Smith        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
679e884886fSBarry Smith  where
680e884886fSBarry Smith      error_rel = square root of relative error in function evaluation
681e884886fSBarry Smith      umin = minimum iterate parameter
682e884886fSBarry Smith .ve
683e884886fSBarry Smith 
68411a5261eSBarry Smith    You can call `SNESSetJacobian()` with `MatMFFDComputeJacobian()` if you are using matrix and not a different
685ca93e954SBarry Smith    preconditioner matrix
686ca93e954SBarry Smith 
68711a5261eSBarry Smith    The user can set the error_rel via `MatMFFDSetFunctionError()` and
688*2ef1f0ffSBarry Smith    umin via `MatMFFDDSSetUmin()`.
689e884886fSBarry Smith 
69011a5261eSBarry Smith    The user should call `MatDestroy()` when finished with the matrix-free
691e884886fSBarry Smith    matrix context.
692e884886fSBarry Smith 
693*2ef1f0ffSBarry Smith .seealso: [](chapter_matrices), `Mat`, `MATMFFD`, `MatDestroy()`, `MatMFFDSetFunctionError()`, `MatMFFDDSSetUmin()`, `MatMFFDSetFunction()`
694db781477SPatrick Sanan           `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`, `MatCreateSNESMF()`,
695db781477SPatrick Sanan           `MatMFFDGetH()`, `MatMFFDRegister()`, `MatMFFDComputeJacobian()`
696e884886fSBarry Smith @*/
697d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateMFFD(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt M, PetscInt N, Mat *J)
698d71ae5a4SJacob Faibussowitsch {
699e884886fSBarry Smith   PetscFunctionBegin;
7009566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, J));
7019566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*J, m, n, M, N));
7029566063dSJacob Faibussowitsch   PetscCall(MatSetType(*J, MATMFFD));
7039566063dSJacob Faibussowitsch   PetscCall(MatSetUp(*J));
7043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
705e884886fSBarry Smith }
706e884886fSBarry Smith 
707e884886fSBarry Smith /*@
70811a5261eSBarry Smith    MatMFFDGetH - Gets the last value that was used as the differencing for a `MATMFFD` matrix
709e884886fSBarry Smith    parameter.
710e884886fSBarry Smith 
711e884886fSBarry Smith    Not Collective
712e884886fSBarry Smith 
713e884886fSBarry Smith    Input Parameters:
71411a5261eSBarry Smith .  mat - the `MATMFFD` matrix
715e884886fSBarry Smith 
716fd292e60Sprj-    Output Parameter:
717e884886fSBarry Smith .  h - the differencing step size
718e884886fSBarry Smith 
719e884886fSBarry Smith    Level: advanced
720e884886fSBarry Smith 
721*2ef1f0ffSBarry Smith .seealso: [](chapter_matrices), `Mat`, `MATMFFD`, `MatCreateSNESMF()`, `MatMFFDSetHHistory()`, `MatCreateMFFD()`, `MATMFFD`, `MatMFFDResetHHistory()`
722e884886fSBarry Smith @*/
723d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDGetH(Mat mat, PetscScalar *h)
724d71ae5a4SJacob Faibussowitsch {
725e884886fSBarry Smith   PetscFunctionBegin;
72688b4c220SStefano Zampini   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
727dadcf809SJacob Faibussowitsch   PetscValidScalarPointer(h, 2);
728cac4c232SBarry Smith   PetscUseMethod(mat, "MatMFFDGetH_C", (Mat, PetscScalar *), (mat, h));
7293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
730e884886fSBarry Smith }
731e884886fSBarry Smith 
732e884886fSBarry Smith /*@C
73311a5261eSBarry Smith    MatMFFDSetFunction - Sets the function used in applying the matrix free `MATMFFD` matrix.
734e884886fSBarry Smith 
735c3339decSBarry Smith    Logically Collective
736e884886fSBarry Smith 
737e884886fSBarry Smith    Input Parameters:
73811a5261eSBarry Smith +  mat - the matrix free matrix `MATMFFD` created via `MatCreateSNESMF()` or `MatCreateMFFD()`
739e884886fSBarry Smith .  func - the function to use
740e884886fSBarry Smith -  funcctx - optional function context passed to function
741e884886fSBarry Smith 
74214f633e2SBarry Smith    Calling Sequence of func:
74314f633e2SBarry Smith $     func (void *funcctx, Vec x, Vec f)
74414f633e2SBarry Smith 
74514f633e2SBarry Smith +  funcctx - user provided context
74614f633e2SBarry Smith .  x - input vector
74714f633e2SBarry Smith -  f - computed output function
74814f633e2SBarry Smith 
749e884886fSBarry Smith    Level: advanced
750e884886fSBarry Smith 
751e884886fSBarry Smith    Notes:
75211a5261eSBarry Smith     If you use this you MUST call `MatAssemblyBegin()` and `MatAssemblyEnd()` on the matrix free
753e884886fSBarry Smith     matrix inside your compute Jacobian routine
754e884886fSBarry Smith 
75511a5261eSBarry Smith     If this is not set then it will use the function set with `SNESSetFunction()` if `MatCreateSNESMF()` was used.
756e884886fSBarry Smith 
757*2ef1f0ffSBarry Smith .seealso: [](chapter_matrices), `Mat`, `MATMFFD`, `MatCreateSNESMF()`, `MatMFFDGetH()`, `MatCreateMFFD()`, `MATMFFD`,
758db781477SPatrick Sanan           `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`, `SNESetFunction()`
759e884886fSBarry Smith @*/
760d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetFunction(Mat mat, PetscErrorCode (*func)(void *, Vec, Vec), void *funcctx)
761d71ae5a4SJacob Faibussowitsch {
762e884886fSBarry Smith   PetscFunctionBegin;
76388b4c220SStefano Zampini   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
764cac4c232SBarry Smith   PetscTryMethod(mat, "MatMFFDSetFunction_C", (Mat, PetscErrorCode(*)(void *, Vec, Vec), void *), (mat, func, funcctx));
7653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
766e884886fSBarry Smith }
767e884886fSBarry Smith 
768e884886fSBarry Smith /*@C
76911a5261eSBarry Smith    MatMFFDSetFunctioni - Sets the function for a single component for a `MATMFFD` matrix
770e884886fSBarry Smith 
771c3339decSBarry Smith    Logically Collective
772e884886fSBarry Smith 
773e884886fSBarry Smith    Input Parameters:
77411a5261eSBarry Smith +  mat - the matrix free matrix `MATMFFD`
775e884886fSBarry Smith -  funci - the function to use
776e884886fSBarry Smith 
777e884886fSBarry Smith    Level: advanced
778e884886fSBarry Smith 
779e884886fSBarry Smith    Notes:
78011a5261eSBarry Smith     If you use this you MUST call `MatAssemblyBegin()` and `MatAssemblyEnd()` on the matrix free
78194eb4328SStefano Zampini     matrix inside your compute Jacobian routine.
78211a5261eSBarry Smith 
78394eb4328SStefano Zampini     This function is necessary to compute the diagonal of the matrix.
784486afcceSStefano Zampini     funci must not contain any MPI call as it is called inside a loop on the local portion of the vector.
785e884886fSBarry Smith 
786*2ef1f0ffSBarry Smith .seealso: [](chapter_matrices), `Mat`, `MATMFFD`, `MatCreateSNESMF()`, `MatMFFDGetH()`, `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`, `SNESetFunction()`, `MatGetDiagonal()`
787e884886fSBarry Smith @*/
788d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetFunctioni(Mat mat, PetscErrorCode (*funci)(void *, PetscInt, Vec, PetscScalar *))
789d71ae5a4SJacob Faibussowitsch {
790e884886fSBarry Smith   PetscFunctionBegin;
7910700a824SBarry Smith   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
792cac4c232SBarry Smith   PetscTryMethod(mat, "MatMFFDSetFunctioni_C", (Mat, PetscErrorCode(*)(void *, PetscInt, Vec, PetscScalar *)), (mat, funci));
7933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
794e884886fSBarry Smith }
795e884886fSBarry Smith 
796e884886fSBarry Smith /*@C
79711a5261eSBarry Smith    MatMFFDSetFunctioniBase - Sets the base vector for a single component function evaluation for a `MATMFFD` matrix
798e884886fSBarry Smith 
799c3339decSBarry Smith    Logically Collective
800e884886fSBarry Smith 
801e884886fSBarry Smith    Input Parameters:
80211a5261eSBarry Smith +  mat - the `MATMFFD` matrix free matrix
803e884886fSBarry Smith -  func - the function to use
804e884886fSBarry Smith 
805e884886fSBarry Smith    Level: advanced
806e884886fSBarry Smith 
807e884886fSBarry Smith    Notes:
80811a5261eSBarry Smith     If you use this you MUST call `MatAssemblyBegin()` and `MatAssemblyEnd()` on the matrix free
80994eb4328SStefano Zampini     matrix inside your compute Jacobian routine.
810e884886fSBarry Smith 
81111a5261eSBarry Smith     This function is necessary to compute the diagonal of the matrix, used for example with `PCJACOBI`
81211a5261eSBarry Smith 
813*2ef1f0ffSBarry Smith .seealso: [](chapter_matrices), `Mat`, `MATMFFD`, `MatCreateSNESMF()`, `MatMFFDGetH()`, `MatCreateMFFD()`, `MATMFFD`
814db781477SPatrick Sanan           `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`, `SNESetFunction()`, `MatGetDiagonal()`
815e884886fSBarry Smith @*/
816d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetFunctioniBase(Mat mat, PetscErrorCode (*func)(void *, Vec))
817d71ae5a4SJacob Faibussowitsch {
818e884886fSBarry Smith   PetscFunctionBegin;
8190700a824SBarry Smith   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
820cac4c232SBarry Smith   PetscTryMethod(mat, "MatMFFDSetFunctioniBase_C", (Mat, PetscErrorCode(*)(void *, Vec)), (mat, func));
8213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
822e884886fSBarry Smith }
823e884886fSBarry Smith 
824e884886fSBarry Smith /*@
82511a5261eSBarry Smith    MatMFFDSetPeriod - Sets how often h is recomputed for a `MATMFFD` matrix, by default it is every time
826e884886fSBarry Smith 
827c3339decSBarry Smith    Logically Collective
828e884886fSBarry Smith 
829e884886fSBarry Smith    Input Parameters:
83011a5261eSBarry Smith +  mat - the `MATMFFD` matrix free matrix
831e884886fSBarry Smith -  period - 1 for every time, 2 for every second etc
832e884886fSBarry Smith 
833*2ef1f0ffSBarry Smith    Options Database Key:
834*2ef1f0ffSBarry Smith .  -mat_mffd_period <period> - Sets how often `h` is recomputed
835e884886fSBarry Smith 
836e884886fSBarry Smith    Level: advanced
837e884886fSBarry Smith 
838*2ef1f0ffSBarry Smith .seealso: [](chapter_matrices), `Mat`, `MATMFFD`, `MatCreateSNESMF()`, `MatMFFDGetH()`,
839db781477SPatrick Sanan           `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`
840e884886fSBarry Smith @*/
841d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetPeriod(Mat mat, PetscInt period)
842d71ae5a4SJacob Faibussowitsch {
843e884886fSBarry Smith   PetscFunctionBegin;
84488b4c220SStefano Zampini   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
84588b4c220SStefano Zampini   PetscValidLogicalCollectiveInt(mat, period, 2);
846cac4c232SBarry Smith   PetscTryMethod(mat, "MatMFFDSetPeriod_C", (Mat, PetscInt), (mat, period));
8473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
848e884886fSBarry Smith }
849e884886fSBarry Smith 
850e884886fSBarry Smith /*@
85111a5261eSBarry Smith    MatMFFDSetFunctionError - Sets the error_rel for the approximation of matrix-vector products using finite differences with the `MATMFFD` matrix
852e884886fSBarry Smith 
853c3339decSBarry Smith    Logically Collective
854e884886fSBarry Smith 
855e884886fSBarry Smith    Input Parameters:
85611a5261eSBarry Smith +  mat - the `MATMFFD` matrix free matrix
85711a5261eSBarry Smith -  error_rel - relative error (should be set to the square root of the relative error in the function evaluations)
858e884886fSBarry Smith 
859*2ef1f0ffSBarry Smith    Options Database Key:
860a2b725a8SWilliam Gropp .  -mat_mffd_err <error_rel> - Sets error_rel
861e884886fSBarry Smith 
862e884886fSBarry Smith    Level: advanced
863e884886fSBarry Smith 
86411a5261eSBarry Smith    Note:
865e884886fSBarry Smith    The default matrix-free matrix-vector product routine computes
866e884886fSBarry Smith .vb
867e884886fSBarry Smith      F'(u)*a = [F(u+h*a) - F(u)]/h where
868e884886fSBarry Smith      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
869e884886fSBarry Smith        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   else
870e884886fSBarry Smith .ve
871e884886fSBarry Smith 
872*2ef1f0ffSBarry Smith .seealso: [](chapter_matrices), `Mat`, `MATMFFD`, `MatCreateSNESMF()`, `MatMFFDGetH()`, `MatCreateMFFD()`, `MATMFFD`
873db781477SPatrick Sanan           `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`
874e884886fSBarry Smith @*/
875d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetFunctionError(Mat mat, PetscReal error)
876d71ae5a4SJacob Faibussowitsch {
877e884886fSBarry Smith   PetscFunctionBegin;
87888b4c220SStefano Zampini   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
87988b4c220SStefano Zampini   PetscValidLogicalCollectiveReal(mat, error, 2);
880cac4c232SBarry Smith   PetscTryMethod(mat, "MatMFFDSetFunctionError_C", (Mat, PetscReal), (mat, error));
8813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
882e884886fSBarry Smith }
883e884886fSBarry Smith 
884e884886fSBarry Smith /*@
885e884886fSBarry Smith    MatMFFDSetHHistory - Sets an array to collect a history of the
88611a5261eSBarry Smith    differencing values (h) computed for the matrix-free product `MATMFFD` matrix
887e884886fSBarry Smith 
888c3339decSBarry Smith    Logically Collective
889e884886fSBarry Smith 
890e884886fSBarry Smith    Input Parameters:
89111a5261eSBarry Smith +  J - the `MATMFFD` matrix-free matrix
892a5b23f4aSJose E. Roman .  history - space to hold the history
893e884886fSBarry Smith -  nhistory - number of entries in history, if more entries are generated than
894e884886fSBarry Smith               nhistory, then the later ones are discarded
895e884886fSBarry Smith 
896e884886fSBarry Smith    Level: advanced
897e884886fSBarry Smith 
89811a5261eSBarry Smith    Note:
89911a5261eSBarry Smith    Use `MatMFFDResetHHistory()` to reset the history counter and collect
900e884886fSBarry Smith    a new batch of differencing parameters, h.
901e884886fSBarry Smith 
902*2ef1f0ffSBarry Smith .seealso: [](chapter_matrices), `Mat`, `MATMFFD`, `MatMFFDGetH()`, `MatCreateSNESMF()`,
903db781477SPatrick Sanan           `MatMFFDResetHHistory()`, `MatMFFDSetFunctionError()`
904e884886fSBarry Smith @*/
905d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetHHistory(Mat J, PetscScalar history[], PetscInt nhistory)
906d71ae5a4SJacob Faibussowitsch {
907e884886fSBarry Smith   PetscFunctionBegin;
90888b4c220SStefano Zampini   PetscValidHeaderSpecific(J, MAT_CLASSID, 1);
909dadcf809SJacob Faibussowitsch   if (history) PetscValidScalarPointer(history, 2);
91088b4c220SStefano Zampini   PetscValidLogicalCollectiveInt(J, nhistory, 3);
911cac4c232SBarry Smith   PetscUseMethod(J, "MatMFFDSetHHistory_C", (Mat, PetscScalar[], PetscInt), (J, history, nhistory));
9123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
913e884886fSBarry Smith }
914e884886fSBarry Smith 
915e884886fSBarry Smith /*@
916e884886fSBarry Smith    MatMFFDResetHHistory - Resets the counter to zero to begin
91711a5261eSBarry Smith    collecting a new set of differencing histories for the `MATMFFD` matrix
918e884886fSBarry Smith 
919c3339decSBarry Smith    Logically Collective
920e884886fSBarry Smith 
921e884886fSBarry Smith    Input Parameters:
922e884886fSBarry Smith .  J - the matrix-free matrix context
923e884886fSBarry Smith 
924e884886fSBarry Smith    Level: advanced
925e884886fSBarry Smith 
92611a5261eSBarry Smith    Note:
92711a5261eSBarry Smith    Use `MatMFFDSetHHistory()` to create the original history counter.
928e884886fSBarry Smith 
929*2ef1f0ffSBarry Smith .seealso: [](chapter_matrices), `Mat`, `MATMFFD`, `MatMFFDGetH()`, `MatCreateSNESMF()`,
930db781477SPatrick Sanan           `MatMFFDSetHHistory()`, `MatMFFDSetFunctionError()`
931e884886fSBarry Smith @*/
932d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDResetHHistory(Mat J)
933d71ae5a4SJacob Faibussowitsch {
934e884886fSBarry Smith   PetscFunctionBegin;
93588b4c220SStefano Zampini   PetscValidHeaderSpecific(J, MAT_CLASSID, 1);
936cac4c232SBarry Smith   PetscTryMethod(J, "MatMFFDResetHHistory_C", (Mat), (J));
9373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
938e884886fSBarry Smith }
939e884886fSBarry Smith 
940487a658cSBarry Smith /*@
941*2ef1f0ffSBarry Smith     MatMFFDSetBase - Sets the vector `U` at which matrix vector products of the
94211a5261eSBarry Smith         Jacobian are computed for the `MATMFFD` matrix
943e884886fSBarry Smith 
944c3339decSBarry Smith     Logically Collective
945e884886fSBarry Smith 
946e884886fSBarry Smith     Input Parameters:
94711a5261eSBarry Smith +   J - the `MATMFFD` matrix
9483ec795f1SBarry Smith .   U - the vector
949bcddec3dSBarry Smith -   F - (optional) vector that contains F(u) if it has been already computed
950e884886fSBarry Smith 
951*2ef1f0ffSBarry Smith     Level: advanced
952*2ef1f0ffSBarry Smith 
95395452b02SPatrick Sanan     Notes:
95495452b02SPatrick Sanan     This is rarely used directly
955e884886fSBarry Smith 
956*2ef1f0ffSBarry Smith     If `F` is provided then it is not recomputed. Otherwise the function is evaluated at the base
95711a5261eSBarry Smith     point during the first `MatMult()` after each call to `MatMFFDSetBase()`.
958dff2f722SBarry Smith 
959*2ef1f0ffSBarry Smith .seealso: [](chapter_matrices), `Mat`, `MATMFFD`, `MatMult()`, `MatMFFDSetBase()`
960e884886fSBarry Smith @*/
961d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetBase(Mat J, Vec U, Vec F)
962d71ae5a4SJacob Faibussowitsch {
963e884886fSBarry Smith   PetscFunctionBegin;
9640700a824SBarry Smith   PetscValidHeaderSpecific(J, MAT_CLASSID, 1);
9650700a824SBarry Smith   PetscValidHeaderSpecific(U, VEC_CLASSID, 2);
9660700a824SBarry Smith   if (F) PetscValidHeaderSpecific(F, VEC_CLASSID, 3);
967cac4c232SBarry Smith   PetscTryMethod(J, "MatMFFDSetBase_C", (Mat, Vec, Vec), (J, U, F));
9683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
969e884886fSBarry Smith }
970e884886fSBarry Smith 
971e884886fSBarry Smith /*@C
972e884886fSBarry Smith     MatMFFDSetCheckh - Sets a function that checks the computed h and adjusts
97311a5261eSBarry Smith         it to satisfy some criteria for the `MATMFFD` matrix
974e884886fSBarry Smith 
975c3339decSBarry Smith     Logically Collective
976e884886fSBarry Smith 
977e884886fSBarry Smith     Input Parameters:
97811a5261eSBarry Smith +   J - the `MATMFFD` matrix
979*2ef1f0ffSBarry Smith .   fun - the function that checks `h`
980e884886fSBarry Smith -   ctx - any context needed by the function
981e884886fSBarry Smith 
982e884886fSBarry Smith     Options Database Keys:
98367b8a455SSatish Balay .   -mat_mffd_check_positivity <bool> - Insure that U + h*a is non-negative
984e884886fSBarry Smith 
985e884886fSBarry Smith     Level: advanced
986e884886fSBarry Smith 
98795452b02SPatrick Sanan     Notes:
988*2ef1f0ffSBarry Smith     For example, `MatMFFDCheckPositivity()` insures that all entries of U + h*a are non-negative
989e884886fSBarry Smith 
990*2ef1f0ffSBarry Smith      The function you provide is called after the default `h` has been computed and allows you to
991755b7f64SBarry Smith      modify it.
992755b7f64SBarry Smith 
993*2ef1f0ffSBarry Smith .seealso: [](chapter_matrices), `Mat`, `MATMFFD`, `MatMFFDCheckPositivity()`
994e884886fSBarry Smith @*/
995d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetCheckh(Mat J, PetscErrorCode (*fun)(void *, Vec, Vec, PetscScalar *), void *ctx)
996d71ae5a4SJacob Faibussowitsch {
997e884886fSBarry Smith   PetscFunctionBegin;
9980700a824SBarry Smith   PetscValidHeaderSpecific(J, MAT_CLASSID, 1);
999cac4c232SBarry Smith   PetscTryMethod(J, "MatMFFDSetCheckh_C", (Mat, PetscErrorCode(*)(void *, Vec, Vec, PetscScalar *), void *), (J, fun, ctx));
10003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1001e884886fSBarry Smith }
1002e884886fSBarry Smith 
1003e884886fSBarry Smith /*@
1004e884886fSBarry Smith     MatMFFDCheckPositivity - Checks that all entries in U + h*a are positive or
100511a5261eSBarry Smith         zero, decreases h until this is satisfied for a `MATMFFD` matrix
1006e884886fSBarry Smith 
1007c3339decSBarry Smith     Logically Collective
1008e884886fSBarry Smith 
1009e884886fSBarry Smith     Input Parameters:
1010e884886fSBarry Smith +   U - base vector that is added to
1011e884886fSBarry Smith .   a - vector that is added
1012e884886fSBarry Smith .   h - scaling factor on a
1013e884886fSBarry Smith -   dummy - context variable (unused)
1014e884886fSBarry Smith 
1015e884886fSBarry Smith     Options Database Keys:
101667b8a455SSatish Balay .   -mat_mffd_check_positivity <bool> - Insure that U + h*a is nonnegative
1017e884886fSBarry Smith 
1018e884886fSBarry Smith     Level: advanced
1019e884886fSBarry Smith 
102011a5261eSBarry Smith     Note:
1021*2ef1f0ffSBarry Smith     This is rarely used directly, rather it is passed as an argument to `MatMFFDSetCheckh()`
1022e884886fSBarry Smith 
1023*2ef1f0ffSBarry Smith .seealso: [](chapter_matrices), `Mat`, `MATMFFD`, `MatMFFDSetCheckh()`
1024e884886fSBarry Smith @*/
1025d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDCheckPositivity(void *dummy, Vec U, Vec a, PetscScalar *h)
1026d71ae5a4SJacob Faibussowitsch {
1027e884886fSBarry Smith   PetscReal    val, minval;
1028e884886fSBarry Smith   PetscScalar *u_vec, *a_vec;
1029e884886fSBarry Smith   PetscInt     i, n;
1030e884886fSBarry Smith   MPI_Comm     comm;
1031e884886fSBarry Smith 
1032e884886fSBarry Smith   PetscFunctionBegin;
103388b4c220SStefano Zampini   PetscValidHeaderSpecific(U, VEC_CLASSID, 2);
103488b4c220SStefano Zampini   PetscValidHeaderSpecific(a, VEC_CLASSID, 3);
1035dadcf809SJacob Faibussowitsch   PetscValidScalarPointer(h, 4);
10369566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)U, &comm));
10379566063dSJacob Faibussowitsch   PetscCall(VecGetArray(U, &u_vec));
10389566063dSJacob Faibussowitsch   PetscCall(VecGetArray(a, &a_vec));
10399566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(U, &n));
1040c068d9bbSLisandro Dalcin   minval = PetscAbsScalar(*h) * PetscRealConstant(1.01);
1041e884886fSBarry Smith   for (i = 0; i < n; i++) {
1042e884886fSBarry Smith     if (PetscRealPart(u_vec[i] + *h * a_vec[i]) <= 0.0) {
1043e884886fSBarry Smith       val = PetscAbsScalar(u_vec[i] / a_vec[i]);
1044e884886fSBarry Smith       if (val < minval) minval = val;
1045e884886fSBarry Smith     }
1046e884886fSBarry Smith   }
10479566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(U, &u_vec));
10489566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(a, &a_vec));
10491c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&minval, &val, 1, MPIU_REAL, MPIU_MIN, comm));
1050e884886fSBarry Smith   if (val <= PetscAbsScalar(*h)) {
10519566063dSJacob Faibussowitsch     PetscCall(PetscInfo(U, "Scaling back h from %g to %g\n", (double)PetscRealPart(*h), (double)(.99 * val)));
1052e884886fSBarry Smith     if (PetscRealPart(*h) > 0.0) *h = 0.99 * val;
1053e884886fSBarry Smith     else *h = -0.99 * val;
1054e884886fSBarry Smith   }
10553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1056e884886fSBarry Smith }
1057