xref: /petsc/src/mat/impls/mffd/mffd.c (revision 66976f2f44dcc61d86a452a70219fb23b45d00f0)
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;
12*66976f2fSJacob Faibussowitsch 
13b022a5c1SBarry Smith /*@C
1411a5261eSBarry Smith   MatMFFDFinalizePackage - This function destroys everything in the MATMFFD` package. It is
1511a5261eSBarry Smith   called from `PetscFinalize()`.
16b022a5c1SBarry Smith 
17b022a5c1SBarry Smith   Level: developer
18b022a5c1SBarry Smith 
191cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `PetscFinalize()`, `MatCreateMFFD()`, `MatCreateSNESMF()`
20b022a5c1SBarry Smith @*/
21d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDFinalizePackage(void)
22d71ae5a4SJacob Faibussowitsch {
23b022a5c1SBarry Smith   PetscFunctionBegin;
249566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListDestroy(&MatMFFDList));
25b022a5c1SBarry Smith   MatMFFDPackageInitialized = PETSC_FALSE;
26b022a5c1SBarry Smith   MatMFFDRegisterAllCalled  = PETSC_FALSE;
273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
28b022a5c1SBarry Smith }
29b022a5c1SBarry Smith 
303ec795f1SBarry Smith /*@C
3111a5261eSBarry Smith   MatMFFDInitializePackage - This function initializes everything in the MATMFFD` package. It is called
3211a5261eSBarry Smith   from `MatInitializePackage()`.
333ec795f1SBarry Smith 
343ec795f1SBarry Smith   Level: developer
353ec795f1SBarry Smith 
361cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `PetscInitialize()`
373ec795f1SBarry Smith @*/
38d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDInitializePackage(void)
39d71ae5a4SJacob Faibussowitsch {
403ec795f1SBarry Smith   char      logList[256];
418e81d068SLisandro Dalcin   PetscBool opt, pkg;
423ec795f1SBarry Smith 
433ec795f1SBarry Smith   PetscFunctionBegin;
443ba16761SJacob Faibussowitsch   if (MatMFFDPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS);
45b022a5c1SBarry Smith   MatMFFDPackageInitialized = PETSC_TRUE;
463ec795f1SBarry Smith   /* Register Classes */
479566063dSJacob Faibussowitsch   PetscCall(PetscClassIdRegister("MatMFFD", &MATMFFD_CLASSID));
483ec795f1SBarry Smith   /* Register Constructors */
499566063dSJacob Faibussowitsch   PetscCall(MatMFFDRegisterAll());
503ec795f1SBarry Smith   /* Register Events */
519566063dSJacob Faibussowitsch   PetscCall(PetscLogEventRegister("MatMult MF", MATMFFD_CLASSID, &MATMFFD_Mult));
52e94e781bSJacob Faibussowitsch   /* Process Info */
53e94e781bSJacob Faibussowitsch   {
54e94e781bSJacob Faibussowitsch     PetscClassId classids[1];
55e94e781bSJacob Faibussowitsch 
56e94e781bSJacob Faibussowitsch     classids[0] = MATMFFD_CLASSID;
579566063dSJacob Faibussowitsch     PetscCall(PetscInfoProcessClass("matmffd", 1, classids));
583ec795f1SBarry Smith   }
593ec795f1SBarry Smith   /* Process summary exclusions */
609566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL, NULL, "-log_exclude", logList, sizeof(logList), &opt));
613ec795f1SBarry Smith   if (opt) {
629566063dSJacob Faibussowitsch     PetscCall(PetscStrInList("matmffd", logList, ',', &pkg));
639566063dSJacob Faibussowitsch     if (pkg) PetscCall(PetscLogEventExcludeClass(MATMFFD_CLASSID));
643ec795f1SBarry Smith   }
658e81d068SLisandro Dalcin   /* Register package finalizer */
669566063dSJacob Faibussowitsch   PetscCall(PetscRegisterFinalize(MatMFFDFinalizePackage));
673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
683ec795f1SBarry Smith }
693ec795f1SBarry Smith 
70d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDSetType_MFFD(Mat mat, MatMFFDType ftype)
71d71ae5a4SJacob Faibussowitsch {
72789d8953SBarry Smith   MatMFFD   ctx;
73789d8953SBarry Smith   PetscBool match;
745f80ce2aSJacob Faibussowitsch   PetscErrorCode (*r)(MatMFFD);
75789d8953SBarry Smith 
76789d8953SBarry Smith   PetscFunctionBegin;
77789d8953SBarry Smith   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
784f572ea9SToby Isaac   PetscAssertPointer(ftype, 2);
799566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
80789d8953SBarry Smith 
81789d8953SBarry Smith   /* already set, so just return */
829566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)ctx, ftype, &match));
833ba16761SJacob Faibussowitsch   if (match) PetscFunctionReturn(PETSC_SUCCESS);
84789d8953SBarry Smith 
85789d8953SBarry Smith   /* destroy the old one if it exists */
86dbbe0bcdSBarry Smith   PetscTryTypeMethod(ctx, destroy);
87789d8953SBarry Smith 
889566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(MatMFFDList, ftype, &r));
896adde796SStefano Zampini   PetscCheck(r, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown MatMFFD type %s given", ftype);
909566063dSJacob Faibussowitsch   PetscCall((*r)(ctx));
919566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)ctx, ftype));
923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
93789d8953SBarry Smith }
94789d8953SBarry Smith 
95e884886fSBarry Smith /*@C
96e884886fSBarry Smith   MatMFFDSetType - Sets the method that is used to compute the
97da81f932SPierre Jolivet   differencing parameter for finite difference matrix-free formulations.
98e884886fSBarry Smith 
99e884886fSBarry Smith   Input Parameters:
10011a5261eSBarry Smith + mat   - the "matrix-free" matrix created via `MatCreateSNESMF()`, or `MatCreateMFFD()`
10111a5261eSBarry Smith           or `MatSetType`(mat,`MATMFFD`);
10211a5261eSBarry Smith - ftype - the type requested, either `MATMFFD_WP` or `MATMFFD_DS`
103e884886fSBarry Smith 
104e884886fSBarry Smith   Level: advanced
105e884886fSBarry Smith 
10611a5261eSBarry Smith   Note:
1072ef1f0ffSBarry Smith   For example, such routines can compute `h` for use in
108e884886fSBarry Smith   Jacobian-vector products of the form
1092ef1f0ffSBarry Smith .vb
110e884886fSBarry Smith 
111e884886fSBarry Smith                         F(x+ha) - F(x)
112e884886fSBarry Smith           F'(u)a  ~=  ----------------
113e884886fSBarry Smith                               h
1142ef1f0ffSBarry Smith .ve
115e884886fSBarry Smith 
1161cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MATMFFD_WP`, `MATMFFD_DS`, `MatCreateSNESMF()`, `MatMFFDRegister()`, `MatMFFDSetFunction()`, `MatCreateMFFD()`
117e884886fSBarry Smith @*/
118d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetType(Mat mat, MatMFFDType ftype)
119d71ae5a4SJacob Faibussowitsch {
120e884886fSBarry Smith   PetscFunctionBegin;
1210700a824SBarry Smith   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
1224f572ea9SToby Isaac   PetscAssertPointer(ftype, 2);
123cac4c232SBarry Smith   PetscTryMethod(mat, "MatMFFDSetType_C", (Mat, MatMFFDType), (mat, ftype));
1243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
125e884886fSBarry Smith }
126e884886fSBarry Smith 
1275d21e1f8SStefano Zampini static PetscErrorCode MatGetDiagonal_MFFD(Mat, Vec);
1285d21e1f8SStefano Zampini 
129e884886fSBarry Smith typedef PetscErrorCode (*FCN1)(void *, Vec); /* force argument to next function to not be extern C*/
130d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDSetFunctioniBase_MFFD(Mat mat, FCN1 func)
131d71ae5a4SJacob Faibussowitsch {
132789d8953SBarry Smith   MatMFFD ctx;
133e884886fSBarry Smith 
134e884886fSBarry Smith   PetscFunctionBegin;
1359566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
136e884886fSBarry Smith   ctx->funcisetbase = func;
1373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
138e884886fSBarry Smith }
139e884886fSBarry Smith 
140e884886fSBarry Smith typedef PetscErrorCode (*FCN2)(void *, PetscInt, Vec, PetscScalar *); /* force argument to next function to not be extern C*/
141d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDSetFunctioni_MFFD(Mat mat, FCN2 funci)
142d71ae5a4SJacob Faibussowitsch {
143789d8953SBarry Smith   MatMFFD ctx;
144e884886fSBarry Smith 
145e884886fSBarry Smith   PetscFunctionBegin;
1469566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
147e884886fSBarry Smith   ctx->funci = funci;
1489566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(mat, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiagonal_MFFD));
1493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1505d21e1f8SStefano Zampini }
151789d8953SBarry Smith 
152d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDGetH_MFFD(Mat mat, PetscScalar *h)
153d71ae5a4SJacob Faibussowitsch {
154789d8953SBarry Smith   MatMFFD ctx;
155789d8953SBarry Smith 
156789d8953SBarry Smith   PetscFunctionBegin;
1579566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
158789d8953SBarry Smith   *h = ctx->currenth;
1593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
160e884886fSBarry Smith }
161e884886fSBarry Smith 
162d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDResetHHistory_MFFD(Mat J)
163d71ae5a4SJacob Faibussowitsch {
164789d8953SBarry Smith   MatMFFD ctx;
165bc0f08ceSBarry Smith 
166bc0f08ceSBarry Smith   PetscFunctionBegin;
1679566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J, &ctx));
168bc0f08ceSBarry Smith   ctx->ncurrenth = 0;
1693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
170bc0f08ceSBarry Smith }
171e884886fSBarry Smith 
1721c84c290SBarry Smith /*@C
17320f4b53cSBarry Smith   MatMFFDRegister - Adds a method to the `MATMFFD` registry.
1741c84c290SBarry Smith 
1751c84c290SBarry Smith   Not Collective
1761c84c290SBarry Smith 
1771c84c290SBarry Smith   Input Parameters:
17820f4b53cSBarry Smith + sname    - name of a new user-defined compute-h module
17920f4b53cSBarry Smith - function - routine to create method context
1801c84c290SBarry Smith 
1811c84c290SBarry Smith   Level: developer
1821c84c290SBarry Smith 
18311a5261eSBarry Smith   Note:
18411a5261eSBarry Smith   `MatMFFDRegister()` may be called multiple times to add several user-defined solvers.
1851c84c290SBarry Smith 
186fe59aa6dSJacob Faibussowitsch   Example Usage:
1871c84c290SBarry Smith .vb
188bdf89e91SBarry Smith    MatMFFDRegister("my_h", MyHCreate);
1891c84c290SBarry Smith .ve
1901c84c290SBarry Smith 
1911c84c290SBarry Smith   Then, your solver can be chosen with the procedural interface via
19211a5261eSBarry Smith $     `MatMFFDSetType`(mfctx, "my_h")
1931c84c290SBarry Smith   or at runtime via the option
194be7a6d03SBarry Smith $     -mat_mffd_type my_h
1951c84c290SBarry Smith 
1961cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatMFFDRegisterAll()`, `MatMFFDRegisterDestroy()`
1971c84c290SBarry Smith  @*/
198d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDRegister(const char sname[], PetscErrorCode (*function)(MatMFFD))
199d71ae5a4SJacob Faibussowitsch {
200e884886fSBarry Smith   PetscFunctionBegin;
2019566063dSJacob Faibussowitsch   PetscCall(MatInitializePackage());
2029566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&MatMFFDList, sname, function));
2033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
204e884886fSBarry Smith }
205e884886fSBarry Smith 
206d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDestroy_MFFD(Mat mat)
207d71ae5a4SJacob Faibussowitsch {
208789d8953SBarry Smith   MatMFFD ctx;
209e884886fSBarry Smith 
210e884886fSBarry Smith   PetscFunctionBegin;
2119566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
2129566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ctx->w));
2139566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ctx->current_u));
21448a46eb9SPierre Jolivet   if (ctx->current_f_allocated) PetscCall(VecDestroy(&ctx->current_f));
215dbbe0bcdSBarry Smith   PetscTryTypeMethod(ctx, destroy);
2169566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(&ctx));
217e884886fSBarry Smith 
2189566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetBase_C", NULL));
2199566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetFunctioniBase_C", NULL));
2209566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetFunctioni_C", NULL));
2219566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetFunction_C", NULL));
2229566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetFunctionError_C", NULL));
2239566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetCheckh_C", NULL));
2249566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetPeriod_C", NULL));
2259566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDResetHHistory_C", NULL));
2269566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetHHistory_C", NULL));
2279566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetType_C", NULL));
2289566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDGetH_C", NULL));
2292e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatSNESMFSetReuseBase_C", NULL));
2302e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatSNESMFGetReuseBase_C", NULL));
2313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
232e884886fSBarry Smith }
233e884886fSBarry Smith 
234e884886fSBarry Smith /*
235e884886fSBarry Smith    MatMFFDView_MFFD - Views matrix-free parameters.
236e884886fSBarry Smith 
237e884886fSBarry Smith */
238d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_MFFD(Mat J, PetscViewer viewer)
239d71ae5a4SJacob Faibussowitsch {
240789d8953SBarry Smith   MatMFFD     ctx;
241a283635eSDmitry Karpeev   PetscBool   iascii, viewbase, viewfunction;
242a283635eSDmitry Karpeev   const char *prefix;
243e884886fSBarry Smith 
244e884886fSBarry Smith   PetscFunctionBegin;
2459566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J, &ctx));
2469566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
247e884886fSBarry Smith   if (iascii) {
2489566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Matrix-free approximation:\n"));
2499566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
2509566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "err=%g (relative error in function evaluation)\n", (double)ctx->error_rel));
2517adad957SLisandro Dalcin     if (!((PetscObject)ctx)->type_name) {
2529566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "The compute h routine has not yet been set\n"));
253e884886fSBarry Smith     } else {
2549566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Using %s compute h routine\n", ((PetscObject)ctx)->type_name));
255e884886fSBarry Smith     }
256c92e3469SBarry Smith #if defined(PETSC_USE_COMPLEX)
25748a46eb9SPierre Jolivet     if (ctx->usecomplex) PetscCall(PetscViewerASCIIPrintf(viewer, "Using Lyness complex number trick to compute the matrix-vector product\n"));
258c92e3469SBarry Smith #endif
259dbbe0bcdSBarry Smith     PetscTryTypeMethod(ctx, view, viewer);
2609566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetOptionsPrefix((PetscObject)J, &prefix));
261a283635eSDmitry Karpeev 
2629566063dSJacob Faibussowitsch     PetscCall(PetscOptionsHasName(((PetscObject)J)->options, prefix, "-mat_mffd_view_base", &viewbase));
263a283635eSDmitry Karpeev     if (viewbase) {
2649566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Base:\n"));
2659566063dSJacob Faibussowitsch       PetscCall(VecView(ctx->current_u, viewer));
266a283635eSDmitry Karpeev     }
2679566063dSJacob Faibussowitsch     PetscCall(PetscOptionsHasName(((PetscObject)J)->options, prefix, "-mat_mffd_view_function", &viewfunction));
268a283635eSDmitry Karpeev     if (viewfunction) {
2699566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Function:\n"));
2709566063dSJacob Faibussowitsch       PetscCall(VecView(ctx->current_f, viewer));
271a283635eSDmitry Karpeev     }
2729566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
27311aeaf0aSBarry Smith   }
2743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
275e884886fSBarry Smith }
276e884886fSBarry Smith 
277e884886fSBarry Smith /*
278e884886fSBarry Smith    MatAssemblyEnd_MFFD - Resets the ctx->ncurrenth to zero. This
279e884886fSBarry Smith    allows the user to indicate the beginning of a new linear solve by calling
28001c1178eSBarry Smith    MatAssemblyXXX() on the matrix-free matrix. This then allows the
2811d0fab5eSBarry Smith    MatCreateMFFD_WP() to properly compute ||U|| only the first time
282e884886fSBarry Smith    in the linear solver rather than every time.
2835a576424SJed Brown 
284cc2e6a90SBarry Smith    This function is referenced directly from MatAssemblyEnd_SNESMF(), which may be in a different shared library hence
285cc2e6a90SBarry Smith    it must be labeled as PETSC_EXTERN
286e884886fSBarry Smith */
287d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatAssemblyEnd_MFFD(Mat J, MatAssemblyType mt)
288d71ae5a4SJacob Faibussowitsch {
289789d8953SBarry Smith   MatMFFD j;
290e884886fSBarry Smith 
291e884886fSBarry Smith   PetscFunctionBegin;
2929566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J, &j));
2939566063dSJacob Faibussowitsch   PetscCall(MatMFFDResetHHistory(J));
2943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
295e884886fSBarry Smith }
296e884886fSBarry Smith 
297e884886fSBarry Smith /*
298e884886fSBarry Smith   MatMult_MFFD - Default matrix-free form for Jacobian-vector product, y = F'(u)*a:
299e884886fSBarry Smith 
300e884886fSBarry Smith         y ~= (F(u + ha) - F(u))/h,
301e884886fSBarry Smith   where F = nonlinear function, as set by SNESSetFunction()
302e884886fSBarry Smith         u = current iterate
303e884886fSBarry Smith         h = difference interval
304e884886fSBarry Smith */
305d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMult_MFFD(Mat mat, Vec a, Vec y)
306d71ae5a4SJacob Faibussowitsch {
307789d8953SBarry Smith   MatMFFD     ctx;
308e884886fSBarry Smith   PetscScalar h;
309e884886fSBarry Smith   Vec         w, U, F;
310ace3abfcSBarry Smith   PetscBool   zeroa;
311e884886fSBarry Smith 
312e884886fSBarry Smith   PetscFunctionBegin;
3139566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
31428b400f6SJacob 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");
315e884886fSBarry Smith   /* We log matrix-free matrix-vector products separately, so that we can
316e884886fSBarry Smith      separate the performance monitoring from the cases that use conventional
317e884886fSBarry Smith      storage.  We may eventually modify event logging to associate events
318e884886fSBarry Smith      with particular objects, hence alleviating the more general problem. */
3199566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(MATMFFD_Mult, a, y, 0, 0));
320e884886fSBarry Smith 
321e884886fSBarry Smith   w = ctx->w;
322e884886fSBarry Smith   U = ctx->current_u;
3233ec795f1SBarry Smith   F = ctx->current_f;
324e884886fSBarry Smith   /*
325e884886fSBarry Smith       Compute differencing parameter
326e884886fSBarry Smith   */
3274a2f8832SBarry Smith   if (!((PetscObject)ctx)->type_name) {
3289566063dSJacob Faibussowitsch     PetscCall(MatMFFDSetType(mat, MATMFFD_WP));
3299566063dSJacob Faibussowitsch     PetscCall(MatSetFromOptions(mat));
330e884886fSBarry Smith   }
331dbbe0bcdSBarry Smith   PetscUseTypeMethod(ctx, compute, U, a, &h, &zeroa);
332e884886fSBarry Smith   if (zeroa) {
3339566063dSJacob Faibussowitsch     PetscCall(VecSet(y, 0.0));
33485d14658SLisandro Dalcin     PetscCall(PetscLogEventEnd(MATMFFD_Mult, a, y, 0, 0));
3353ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
336e884886fSBarry Smith   }
337e884886fSBarry Smith 
33808401ef6SPierre Jolivet   PetscCheck(!mat->erroriffailure || !PetscIsInfOrNanScalar(h), PETSC_COMM_SELF, PETSC_ERR_PLIB, "Computed Nan differencing parameter h");
33948a46eb9SPierre Jolivet   if (ctx->checkh) PetscCall((*ctx->checkh)(ctx->checkhctx, U, a, &h));
340e884886fSBarry Smith 
341e884886fSBarry Smith   /* keep a record of the current differencing parameter h */
342e884886fSBarry Smith   ctx->currenth = h;
343e884886fSBarry Smith #if defined(PETSC_USE_COMPLEX)
3449566063dSJacob Faibussowitsch   PetscCall(PetscInfo(mat, "Current differencing parameter: %g + %g i\n", (double)PetscRealPart(h), (double)PetscImaginaryPart(h)));
345e884886fSBarry Smith #else
3469566063dSJacob Faibussowitsch   PetscCall(PetscInfo(mat, "Current differencing parameter: %15.12e\n", (double)PetscRealPart(h)));
347e884886fSBarry Smith #endif
348ad540459SPierre Jolivet   if (ctx->historyh && ctx->ncurrenth < ctx->maxcurrenth) ctx->historyh[ctx->ncurrenth] = h;
349e884886fSBarry Smith   ctx->ncurrenth++;
350e884886fSBarry Smith 
351c92e3469SBarry Smith #if defined(PETSC_USE_COMPLEX)
352c92e3469SBarry Smith   if (ctx->usecomplex) h = PETSC_i * h;
353c92e3469SBarry Smith #endif
354c92e3469SBarry Smith 
355e884886fSBarry Smith   /* w = u + ha */
3569566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(w, h, a, U));
357e884886fSBarry Smith 
3581b797266SDmitry Karpeev   /* compute func(U) as base for differencing; only needed first time in and not when provided by user */
35948a46eb9SPierre Jolivet   if (ctx->ncurrenth == 1 && ctx->current_f_allocated) PetscCall((*ctx->func)(ctx->funcctx, U, F));
3609566063dSJacob Faibussowitsch   PetscCall((*ctx->func)(ctx->funcctx, w, y));
361e884886fSBarry Smith 
362c92e3469SBarry Smith #if defined(PETSC_USE_COMPLEX)
363c92e3469SBarry Smith   if (ctx->usecomplex) {
3649566063dSJacob Faibussowitsch     PetscCall(VecImaginaryPart(y));
365c92e3469SBarry Smith     h = PetscImaginaryPart(h);
366c92e3469SBarry Smith   } else {
3679566063dSJacob Faibussowitsch     PetscCall(VecAXPY(y, -1.0, F));
368c92e3469SBarry Smith   }
369c92e3469SBarry Smith #else
3709566063dSJacob Faibussowitsch   PetscCall(VecAXPY(y, -1.0, F));
371c92e3469SBarry Smith #endif
3729566063dSJacob Faibussowitsch   PetscCall(VecScale(y, 1.0 / h));
3739566063dSJacob Faibussowitsch   if (mat->nullsp) PetscCall(MatNullSpaceRemove(mat->nullsp, y));
374e884886fSBarry Smith 
3759566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(MATMFFD_Mult, a, y, 0, 0));
3763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
377e884886fSBarry Smith }
378e884886fSBarry Smith 
379e884886fSBarry Smith /*
38001c1178eSBarry Smith   MatGetDiagonal_MFFD - Gets the diagonal for a matrix-free matrix
381e884886fSBarry Smith 
382e884886fSBarry Smith         y ~= (F(u + ha) - F(u))/h,
383e884886fSBarry Smith   where F = nonlinear function, as set by SNESSetFunction()
384e884886fSBarry Smith         u = current iterate
385e884886fSBarry Smith         h = difference interval
386e884886fSBarry Smith */
387*66976f2fSJacob Faibussowitsch static PetscErrorCode MatGetDiagonal_MFFD(Mat mat, Vec a)
388d71ae5a4SJacob Faibussowitsch {
389789d8953SBarry Smith   MatMFFD     ctx;
390e884886fSBarry Smith   PetscScalar h, *aa, *ww, v;
391e884886fSBarry Smith   PetscReal   epsilon = PETSC_SQRT_MACHINE_EPSILON, umin = 100.0 * PETSC_SQRT_MACHINE_EPSILON;
392e884886fSBarry Smith   Vec         w, U;
393e884886fSBarry Smith   PetscInt    i, rstart, rend;
394e884886fSBarry Smith 
395e884886fSBarry Smith   PetscFunctionBegin;
3969566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
39728b400f6SJacob Faibussowitsch   PetscCheck(ctx->func, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Requires calling MatMFFDSetFunction() first");
39828b400f6SJacob Faibussowitsch   PetscCheck(ctx->funci, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Requires calling MatMFFDSetFunctioni() first");
399e884886fSBarry Smith   w = ctx->w;
400e884886fSBarry Smith   U = ctx->current_u;
4019566063dSJacob Faibussowitsch   PetscCall((*ctx->func)(ctx->funcctx, U, a));
4021baa6e33SBarry Smith   if (ctx->funcisetbase) PetscCall((*ctx->funcisetbase)(ctx->funcctx, U));
4039566063dSJacob Faibussowitsch   PetscCall(VecCopy(U, w));
404e884886fSBarry Smith 
4059566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(a, &rstart, &rend));
4069566063dSJacob Faibussowitsch   PetscCall(VecGetArray(a, &aa));
407e884886fSBarry Smith   for (i = rstart; i < rend; i++) {
4089566063dSJacob Faibussowitsch     PetscCall(VecGetArray(w, &ww));
409e884886fSBarry Smith     h = ww[i - rstart];
410e884886fSBarry Smith     if (h == 0.0) h = 1.0;
411e884886fSBarry Smith     if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0) h = umin;
412e884886fSBarry Smith     else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin;
413e884886fSBarry Smith     h *= epsilon;
414e884886fSBarry Smith 
415e884886fSBarry Smith     ww[i - rstart] += h;
4169566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(w, &ww));
4179566063dSJacob Faibussowitsch     PetscCall((*ctx->funci)(ctx->funcctx, i, w, &v));
418e884886fSBarry Smith     aa[i - rstart] = (v - aa[i - rstart]) / h;
419e884886fSBarry Smith 
4209566063dSJacob Faibussowitsch     PetscCall(VecGetArray(w, &ww));
421e884886fSBarry Smith     ww[i - rstart] -= h;
4229566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(w, &ww));
423e884886fSBarry Smith   }
4249566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(a, &aa));
4253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
426e884886fSBarry Smith }
427e884886fSBarry Smith 
428d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatMFFDSetBase_MFFD(Mat J, Vec U, Vec F)
429d71ae5a4SJacob Faibussowitsch {
430789d8953SBarry Smith   MatMFFD ctx;
431e884886fSBarry Smith 
432e884886fSBarry Smith   PetscFunctionBegin;
4339566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J, &ctx));
4349566063dSJacob Faibussowitsch   PetscCall(MatMFFDResetHHistory(J));
435c51ad4d4SStefano Zampini   if (!ctx->current_u) {
4369566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(U, &ctx->current_u));
4379566063dSJacob Faibussowitsch     PetscCall(VecLockReadPush(ctx->current_u));
438c51ad4d4SStefano Zampini   }
4399566063dSJacob Faibussowitsch   PetscCall(VecLockReadPop(ctx->current_u));
4409566063dSJacob Faibussowitsch   PetscCall(VecCopy(U, ctx->current_u));
4419566063dSJacob Faibussowitsch   PetscCall(VecLockReadPush(ctx->current_u));
44252121784SBarry Smith   if (F) {
4439566063dSJacob Faibussowitsch     if (ctx->current_f_allocated) PetscCall(VecDestroy(&ctx->current_f));
4443ec795f1SBarry Smith     ctx->current_f           = F;
445cfe22f5eSBarry Smith     ctx->current_f_allocated = PETSC_FALSE;
446cfe22f5eSBarry Smith   } else if (!ctx->current_f_allocated) {
4479566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(J, NULL, &ctx->current_f));
448cfe22f5eSBarry Smith     ctx->current_f_allocated = PETSC_TRUE;
44952121784SBarry Smith   }
45048a46eb9SPierre Jolivet   if (!ctx->w) PetscCall(VecDuplicate(ctx->current_u, &ctx->w));
451e884886fSBarry Smith   J->assembled = PETSC_TRUE;
4523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
453e884886fSBarry Smith }
4545a576424SJed Brown 
455e884886fSBarry Smith typedef PetscErrorCode (*FCN3)(void *, Vec, Vec, PetscScalar *); /* force argument to next function to not be extern C*/
456b2573a8aSBarry Smith 
457d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDSetCheckh_MFFD(Mat J, FCN3 fun, void *ectx)
458d71ae5a4SJacob Faibussowitsch {
459789d8953SBarry Smith   MatMFFD ctx;
460e884886fSBarry Smith 
461e884886fSBarry Smith   PetscFunctionBegin;
4629566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J, &ctx));
463e884886fSBarry Smith   ctx->checkh    = fun;
464e884886fSBarry Smith   ctx->checkhctx = ectx;
4653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
466e884886fSBarry Smith }
467e884886fSBarry Smith 
4686aa9148fSLisandro Dalcin /*@C
4696aa9148fSLisandro Dalcin   MatMFFDSetOptionsPrefix - Sets the prefix used for searching for all
47011a5261eSBarry Smith   MATMFFD` options in the database.
4716aa9148fSLisandro Dalcin 
472c3339decSBarry Smith   Collective
4736aa9148fSLisandro Dalcin 
474d8d19677SJose E. Roman   Input Parameters:
475fe59aa6dSJacob Faibussowitsch + mat    - the `MATMFFD` context
4766aa9148fSLisandro Dalcin - prefix - the prefix to prepend to all option names
4776aa9148fSLisandro Dalcin 
47811a5261eSBarry Smith   Note:
4796aa9148fSLisandro Dalcin   A hyphen (-) must NOT be given at the beginning of the prefix name.
4806aa9148fSLisandro Dalcin   The first character of all runtime options is AUTOMATICALLY the hyphen.
4816aa9148fSLisandro Dalcin 
4826aa9148fSLisandro Dalcin   Level: advanced
4836aa9148fSLisandro Dalcin 
4841cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatSetFromOptions()`, `MatCreateSNESMF()`, `MatCreateMFFD()`
4856aa9148fSLisandro Dalcin @*/
486d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetOptionsPrefix(Mat mat, const char prefix[])
487d71ae5a4SJacob Faibussowitsch {
488789d8953SBarry Smith   MatMFFD mfctx;
4895fd66863SKarl Rupp 
4906aa9148fSLisandro Dalcin   PetscFunctionBegin;
4910700a824SBarry Smith   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
4929566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &mfctx));
4930700a824SBarry Smith   PetscValidHeaderSpecific(mfctx, MATMFFD_CLASSID, 1);
4949566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)mfctx, prefix));
4953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4966aa9148fSLisandro Dalcin }
4976aa9148fSLisandro Dalcin 
498d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSetFromOptions_MFFD(Mat mat, PetscOptionItems *PetscOptionsObject)
499d71ae5a4SJacob Faibussowitsch {
500789d8953SBarry Smith   MatMFFD   mfctx;
501ace3abfcSBarry Smith   PetscBool flg;
502e884886fSBarry Smith   char      ftype[256];
503e884886fSBarry Smith 
504e884886fSBarry Smith   PetscFunctionBegin;
5059566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &mfctx));
506dbbe0bcdSBarry Smith   PetscValidHeaderSpecific(mfctx, MATMFFD_CLASSID, 1);
507d0609cedSBarry Smith   PetscObjectOptionsBegin((PetscObject)mfctx);
5089566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFList("-mat_mffd_type", "Matrix free type", "MatMFFDSetType", MatMFFDList, ((PetscObject)mfctx)->type_name, ftype, 256, &flg));
5091baa6e33SBarry Smith   if (flg) PetscCall(MatMFFDSetType(mat, ftype));
510e884886fSBarry Smith 
5119566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_mffd_err", "set sqrt relative error in function", "MatMFFDSetFunctionError", mfctx->error_rel, &mfctx->error_rel, NULL));
5129566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-mat_mffd_period", "how often h is recomputed", "MatMFFDSetPeriod", mfctx->recomputeperiod, &mfctx->recomputeperiod, NULL));
513e884886fSBarry Smith 
51490d69ab7SBarry Smith   flg = PETSC_FALSE;
5159566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-mat_mffd_check_positivity", "Insure that U + h*a is nonnegative", "MatMFFDSetCheckh", flg, &flg, NULL));
5161baa6e33SBarry Smith   if (flg) PetscCall(MatMFFDSetCheckh(mat, MatMFFDCheckPositivity, NULL));
517c92e3469SBarry Smith #if defined(PETSC_USE_COMPLEX)
5189566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-mat_mffd_complex", "Use Lyness complex number trick to compute the matrix-vector product", "None", mfctx->usecomplex, &mfctx->usecomplex, NULL));
519c92e3469SBarry Smith #endif
520dbbe0bcdSBarry Smith   PetscTryTypeMethod(mfctx, setfromoptions, PetscOptionsObject);
521d0609cedSBarry Smith   PetscOptionsEnd();
5223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
523e884886fSBarry Smith }
524e884886fSBarry Smith 
525d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDSetPeriod_MFFD(Mat mat, PetscInt period)
526d71ae5a4SJacob Faibussowitsch {
527789d8953SBarry Smith   MatMFFD ctx;
528bc0f08ceSBarry Smith 
529bc0f08ceSBarry Smith   PetscFunctionBegin;
5309566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
531bc0f08ceSBarry Smith   ctx->recomputeperiod = period;
5323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
533bc0f08ceSBarry Smith }
534bc0f08ceSBarry Smith 
535d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDSetFunction_MFFD(Mat mat, PetscErrorCode (*func)(void *, Vec, Vec), void *funcctx)
536d71ae5a4SJacob Faibussowitsch {
537789d8953SBarry Smith   MatMFFD ctx;
538bc0f08ceSBarry Smith 
539bc0f08ceSBarry Smith   PetscFunctionBegin;
5409566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
541bc0f08ceSBarry Smith   ctx->func    = func;
542bc0f08ceSBarry Smith   ctx->funcctx = funcctx;
5433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
544bc0f08ceSBarry Smith }
545bc0f08ceSBarry Smith 
546d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDSetFunctionError_MFFD(Mat mat, PetscReal error)
547d71ae5a4SJacob Faibussowitsch {
54813bcc0bdSJacob Faibussowitsch   PetscFunctionBegin;
54913bcc0bdSJacob Faibussowitsch   if (error != (PetscReal)PETSC_DEFAULT) {
550789d8953SBarry Smith     MatMFFD ctx;
551bc0f08ceSBarry Smith 
5529566063dSJacob Faibussowitsch     PetscCall(MatShellGetContext(mat, &ctx));
55313bcc0bdSJacob Faibussowitsch     ctx->error_rel = error;
55413bcc0bdSJacob Faibussowitsch   }
5553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
556bc0f08ceSBarry Smith }
557bc0f08ceSBarry Smith 
558*66976f2fSJacob Faibussowitsch static PetscErrorCode MatMFFDSetHHistory_MFFD(Mat J, PetscScalar history[], PetscInt nhistory)
559d71ae5a4SJacob Faibussowitsch {
560789d8953SBarry Smith   MatMFFD ctx;
561789d8953SBarry Smith 
5623b49f96aSBarry Smith   PetscFunctionBegin;
5639566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J, &ctx));
564789d8953SBarry Smith   ctx->historyh    = history;
565789d8953SBarry Smith   ctx->maxcurrenth = nhistory;
566789d8953SBarry Smith   ctx->currenth    = 0.;
5673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5683b49f96aSBarry Smith }
5693b49f96aSBarry Smith 
570e884886fSBarry Smith /*MC
57101c1178eSBarry Smith   MATMFFD - "mffd" - A matrix-free matrix type.
572e884886fSBarry Smith 
573e884886fSBarry Smith   Level: advanced
574e884886fSBarry Smith 
5752ef1f0ffSBarry Smith   Developers Note:
5762ef1f0ffSBarry Smith   This is implemented on top of `MATSHELL` to get support for scaling and shifting without requiring duplicate code
577789d8953SBarry 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));
625e884886fSBarry Smith 
6269566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetBase_C", MatMFFDSetBase_MFFD));
6279566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetFunctioniBase_C", MatMFFDSetFunctioniBase_MFFD));
6289566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetFunctioni_C", MatMFFDSetFunctioni_MFFD));
6299566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetFunction_C", MatMFFDSetFunction_MFFD));
6309566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetCheckh_C", MatMFFDSetCheckh_MFFD));
6319566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetPeriod_C", MatMFFDSetPeriod_MFFD));
6329566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetFunctionError_C", MatMFFDSetFunctionError_MFFD));
6339566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDResetHHistory_C", MatMFFDResetHHistory_MFFD));
6349566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetHHistory_C", MatMFFDSetHHistory_MFFD));
6359566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetType_C", MatMFFDSetType_MFFD));
6369566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDGetH_C", MatMFFDGetH_MFFD));
6379566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATMFFD));
6383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
639e884886fSBarry Smith }
640e884886fSBarry Smith 
641e884886fSBarry Smith /*@
64211a5261eSBarry Smith   MatCreateMFFD - Creates a matrix-free matrix of type `MATMFFD`. See also `MatCreateSNESMF()`
643e884886fSBarry Smith 
64411a5261eSBarry Smith   Collective
645e884886fSBarry Smith 
646e884886fSBarry Smith   Input Parameters:
647fef1beadSBarry Smith + comm - MPI communicator
6482ef1f0ffSBarry Smith . m    - number of local rows (or `PETSC_DECIDE` to have calculated if `M` is given)
649fef1beadSBarry Smith            This value should be the same as the local size used in creating the
650fef1beadSBarry Smith            y vector for the matrix-vector product y = Ax.
651fef1beadSBarry Smith . n    - This value should be the same as the local size used in creating the
65211a5261eSBarry Smith        x vector for the matrix-vector product y = Ax. (or `PETSC_DECIDE` to have
6532ef1f0ffSBarry Smith        calculated if `N` is given) For square matrices `n` is almost always `m`.
6542ef1f0ffSBarry Smith . M    - number of global rows (or `PETSC_DETERMINE` to have calculated if `m` is given)
6552ef1f0ffSBarry Smith - N    - number of global columns (or `PETSC_DETERMINE` to have calculated if `n` is given)
656fef1beadSBarry Smith 
657e884886fSBarry Smith   Output Parameter:
658e884886fSBarry Smith . J - the matrix-free matrix
659e884886fSBarry Smith 
66011a5261eSBarry Smith   Options Database Keys:
66111a5261eSBarry Smith + -mat_mffd_type             - wp or ds (see `MATMFFD_WP` or `MATMFFD_DS`)
662c92e3469SBarry Smith . -mat_mffd_err              - square root of estimated relative error in function evaluation
663c92e3469SBarry Smith . -mat_mffd_period           - how often h is recomputed, defaults to 1, every time
6642ef1f0ffSBarry Smith . -mat_mffd_check_positivity - possibly decrease `h` until U + h*a has only positive values
66511a5261eSBarry Smith . -mat_mffd_umin <umin>      - Sets umin (for default PETSc routine that computes h only)
666c92e3469SBarry Smith - -mat_mffd_complex          - use the Lyness trick with complex numbers to compute the matrix-vector product instead of differencing
667c92e3469SBarry Smith                        (requires real valued functions but that PETSc be configured for complex numbers)
6689c6ac3b3SBarry Smith 
669e884886fSBarry Smith   Level: advanced
670e884886fSBarry Smith 
671e884886fSBarry Smith   Notes:
672e884886fSBarry Smith   The matrix-free matrix context merely contains the function pointers
673e884886fSBarry Smith   and work space for performing finite difference approximations of
674e884886fSBarry Smith   Jacobian-vector products, F'(u)*a,
675e884886fSBarry Smith 
676e884886fSBarry Smith   The default code uses the following approach to compute h
677e884886fSBarry Smith 
678e884886fSBarry Smith .vb
679e884886fSBarry Smith      F'(u)*a = [F(u+h*a) - F(u)]/h where
680e884886fSBarry Smith      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
681e884886fSBarry Smith        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
682e884886fSBarry Smith  where
683e884886fSBarry Smith      error_rel = square root of relative error in function evaluation
684e884886fSBarry Smith      umin = minimum iterate parameter
685e884886fSBarry Smith .ve
686e884886fSBarry Smith 
68711a5261eSBarry Smith   You can call `SNESSetJacobian()` with `MatMFFDComputeJacobian()` if you are using matrix and not a different
688ca93e954SBarry Smith   preconditioner matrix
689ca93e954SBarry Smith 
69011a5261eSBarry Smith   The user can set the error_rel via `MatMFFDSetFunctionError()` and
6912ef1f0ffSBarry Smith   umin via `MatMFFDDSSetUmin()`.
692e884886fSBarry Smith 
69311a5261eSBarry Smith   The user should call `MatDestroy()` when finished with the matrix-free
694e884886fSBarry Smith   matrix context.
695e884886fSBarry Smith 
6961cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatDestroy()`, `MatMFFDSetFunctionError()`, `MatMFFDDSSetUmin()`, `MatMFFDSetFunction()`
697db781477SPatrick Sanan           `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`, `MatCreateSNESMF()`,
698db781477SPatrick Sanan           `MatMFFDGetH()`, `MatMFFDRegister()`, `MatMFFDComputeJacobian()`
699e884886fSBarry Smith @*/
700d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateMFFD(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt M, PetscInt N, Mat *J)
701d71ae5a4SJacob Faibussowitsch {
702e884886fSBarry Smith   PetscFunctionBegin;
7039566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, J));
7049566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*J, m, n, M, N));
7059566063dSJacob Faibussowitsch   PetscCall(MatSetType(*J, MATMFFD));
7069566063dSJacob Faibussowitsch   PetscCall(MatSetUp(*J));
7073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
708e884886fSBarry Smith }
709e884886fSBarry Smith 
710e884886fSBarry Smith /*@
71111a5261eSBarry Smith   MatMFFDGetH - Gets the last value that was used as the differencing for a `MATMFFD` matrix
712e884886fSBarry Smith   parameter.
713e884886fSBarry Smith 
714e884886fSBarry Smith   Not Collective
715e884886fSBarry Smith 
716e884886fSBarry Smith   Input Parameters:
71711a5261eSBarry Smith . mat - the `MATMFFD` matrix
718e884886fSBarry Smith 
719fd292e60Sprj-   Output Parameter:
720e884886fSBarry Smith . h - the differencing step size
721e884886fSBarry Smith 
722e884886fSBarry Smith   Level: advanced
723e884886fSBarry Smith 
724fe59aa6dSJacob Faibussowitsch .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatCreateSNESMF()`, `MatMFFDSetHHistory()`, `MatCreateMFFD()`, `MatMFFDResetHHistory()`
725e884886fSBarry Smith @*/
726d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDGetH(Mat mat, PetscScalar *h)
727d71ae5a4SJacob Faibussowitsch {
728e884886fSBarry Smith   PetscFunctionBegin;
72988b4c220SStefano Zampini   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
7304f572ea9SToby Isaac   PetscAssertPointer(h, 2);
731cac4c232SBarry Smith   PetscUseMethod(mat, "MatMFFDGetH_C", (Mat, PetscScalar *), (mat, h));
7323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
733e884886fSBarry Smith }
734e884886fSBarry Smith 
735e884886fSBarry Smith /*@C
73601c1178eSBarry Smith   MatMFFDSetFunction - Sets the function used in applying the matrix-free `MATMFFD` matrix.
737e884886fSBarry Smith 
738c3339decSBarry Smith   Logically Collective
739e884886fSBarry Smith 
740e884886fSBarry Smith   Input Parameters:
74101c1178eSBarry Smith + mat     - the matrix-free matrix `MATMFFD` created via `MatCreateSNESMF()` or `MatCreateMFFD()`
742e884886fSBarry Smith . func    - the function to use
743e884886fSBarry Smith - funcctx - optional function context passed to function
744e884886fSBarry Smith 
7452920cce0SJacob Faibussowitsch   Calling sequence of `func`:
74614f633e2SBarry Smith + funcctx - user provided context
74714f633e2SBarry Smith . x       - input vector
74814f633e2SBarry Smith - f       - computed output function
74914f633e2SBarry Smith 
750e884886fSBarry Smith   Level: advanced
751e884886fSBarry Smith 
752e884886fSBarry Smith   Notes:
75301c1178eSBarry Smith   If you use this you MUST call `MatAssemblyBegin()` and `MatAssemblyEnd()` on the matrix-free
754e884886fSBarry Smith   matrix inside your compute Jacobian routine
755e884886fSBarry Smith 
75611a5261eSBarry Smith   If this is not set then it will use the function set with `SNESSetFunction()` if `MatCreateSNESMF()` was used.
757e884886fSBarry Smith 
758fe59aa6dSJacob Faibussowitsch .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatCreateSNESMF()`, `MatMFFDGetH()`, `MatCreateMFFD()`,
7592920cce0SJacob Faibussowitsch           `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`, `SNESSetFunction()`
760e884886fSBarry Smith @*/
7612920cce0SJacob Faibussowitsch PetscErrorCode MatMFFDSetFunction(Mat mat, PetscErrorCode (*func)(void *funcctx, Vec x, Vec f), void *funcctx)
762d71ae5a4SJacob Faibussowitsch {
763e884886fSBarry Smith   PetscFunctionBegin;
76488b4c220SStefano Zampini   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
765cac4c232SBarry Smith   PetscTryMethod(mat, "MatMFFDSetFunction_C", (Mat, PetscErrorCode(*)(void *, Vec, Vec), void *), (mat, func, funcctx));
7663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
767e884886fSBarry Smith }
768e884886fSBarry Smith 
769e884886fSBarry Smith /*@C
77011a5261eSBarry Smith   MatMFFDSetFunctioni - Sets the function for a single component for a `MATMFFD` matrix
771e884886fSBarry Smith 
772c3339decSBarry Smith   Logically Collective
773e884886fSBarry Smith 
774e884886fSBarry Smith   Input Parameters:
77501c1178eSBarry Smith + mat   - the matrix-free matrix `MATMFFD`
776e884886fSBarry Smith - funci - the function to use
777e884886fSBarry Smith 
778e884886fSBarry Smith   Level: advanced
779e884886fSBarry Smith 
780e884886fSBarry Smith   Notes:
78101c1178eSBarry Smith   If you use this you MUST call `MatAssemblyBegin()` and `MatAssemblyEnd()` on the matrix-free
78294eb4328SStefano Zampini   matrix inside your compute Jacobian routine.
78311a5261eSBarry Smith 
78494eb4328SStefano Zampini   This function is necessary to compute the diagonal of the matrix.
785486afcceSStefano Zampini   funci must not contain any MPI call as it is called inside a loop on the local portion of the vector.
786e884886fSBarry Smith 
7871cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatCreateSNESMF()`, `MatMFFDGetH()`, `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`, `SNESetFunction()`, `MatGetDiagonal()`
788e884886fSBarry Smith @*/
789d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetFunctioni(Mat mat, PetscErrorCode (*funci)(void *, PetscInt, Vec, PetscScalar *))
790d71ae5a4SJacob Faibussowitsch {
791e884886fSBarry Smith   PetscFunctionBegin;
7920700a824SBarry Smith   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
793cac4c232SBarry Smith   PetscTryMethod(mat, "MatMFFDSetFunctioni_C", (Mat, PetscErrorCode(*)(void *, PetscInt, Vec, PetscScalar *)), (mat, funci));
7943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
795e884886fSBarry Smith }
796e884886fSBarry Smith 
797e884886fSBarry Smith /*@C
79811a5261eSBarry Smith   MatMFFDSetFunctioniBase - Sets the base vector for a single component function evaluation for a `MATMFFD` matrix
799e884886fSBarry Smith 
800c3339decSBarry Smith   Logically Collective
801e884886fSBarry Smith 
802e884886fSBarry Smith   Input Parameters:
80301c1178eSBarry Smith + mat  - the `MATMFFD` matrix-free matrix
804e884886fSBarry Smith - func - the function to use
805e884886fSBarry Smith 
806e884886fSBarry Smith   Level: advanced
807e884886fSBarry Smith 
808e884886fSBarry Smith   Notes:
80901c1178eSBarry Smith   If you use this you MUST call `MatAssemblyBegin()` and `MatAssemblyEnd()` on the matrix-free
81094eb4328SStefano Zampini   matrix inside your compute Jacobian routine.
811e884886fSBarry Smith 
81211a5261eSBarry Smith   This function is necessary to compute the diagonal of the matrix, used for example with `PCJACOBI`
81311a5261eSBarry Smith 
814fe59aa6dSJacob Faibussowitsch .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatCreateSNESMF()`, `MatMFFDGetH()`, `MatCreateMFFD()`,
815db781477SPatrick Sanan           `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`, `SNESetFunction()`, `MatGetDiagonal()`
816e884886fSBarry Smith @*/
817d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetFunctioniBase(Mat mat, PetscErrorCode (*func)(void *, Vec))
818d71ae5a4SJacob Faibussowitsch {
819e884886fSBarry Smith   PetscFunctionBegin;
8200700a824SBarry Smith   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
821cac4c232SBarry Smith   PetscTryMethod(mat, "MatMFFDSetFunctioniBase_C", (Mat, PetscErrorCode(*)(void *, Vec)), (mat, func));
8223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
823e884886fSBarry Smith }
824e884886fSBarry Smith 
825e884886fSBarry Smith /*@
82611a5261eSBarry Smith   MatMFFDSetPeriod - Sets how often h is recomputed for a `MATMFFD` matrix, by default it is every time
827e884886fSBarry Smith 
828c3339decSBarry Smith   Logically Collective
829e884886fSBarry Smith 
830e884886fSBarry Smith   Input Parameters:
83101c1178eSBarry Smith + mat    - the `MATMFFD` matrix-free matrix
832e884886fSBarry Smith - period - 1 for every time, 2 for every second etc
833e884886fSBarry Smith 
8342ef1f0ffSBarry Smith   Options Database Key:
8352ef1f0ffSBarry Smith . -mat_mffd_period <period> - Sets how often `h` is recomputed
836e884886fSBarry Smith 
837e884886fSBarry Smith   Level: advanced
838e884886fSBarry Smith 
8391cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatCreateSNESMF()`, `MatMFFDGetH()`,
840db781477SPatrick Sanan           `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`
841e884886fSBarry Smith @*/
842d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetPeriod(Mat mat, PetscInt period)
843d71ae5a4SJacob Faibussowitsch {
844e884886fSBarry Smith   PetscFunctionBegin;
84588b4c220SStefano Zampini   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
84688b4c220SStefano Zampini   PetscValidLogicalCollectiveInt(mat, period, 2);
847cac4c232SBarry Smith   PetscTryMethod(mat, "MatMFFDSetPeriod_C", (Mat, PetscInt), (mat, period));
8483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
849e884886fSBarry Smith }
850e884886fSBarry Smith 
851e884886fSBarry Smith /*@
85211a5261eSBarry Smith   MatMFFDSetFunctionError - Sets the error_rel for the approximation of matrix-vector products using finite differences with the `MATMFFD` matrix
853e884886fSBarry Smith 
854c3339decSBarry Smith   Logically Collective
855e884886fSBarry Smith 
856e884886fSBarry Smith   Input Parameters:
85701c1178eSBarry Smith + mat   - the `MATMFFD` matrix-free matrix
858fe59aa6dSJacob Faibussowitsch - error - relative error (should be set to the square root of the relative error in the function evaluations)
859e884886fSBarry Smith 
8602ef1f0ffSBarry Smith   Options Database Key:
861a2b725a8SWilliam Gropp . -mat_mffd_err <error_rel> - Sets error_rel
862e884886fSBarry Smith 
863e884886fSBarry Smith   Level: advanced
864e884886fSBarry Smith 
86511a5261eSBarry Smith   Note:
866e884886fSBarry Smith   The default matrix-free matrix-vector product routine computes
867e884886fSBarry Smith .vb
868e884886fSBarry Smith      F'(u)*a = [F(u+h*a) - F(u)]/h where
869e884886fSBarry Smith      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
870e884886fSBarry Smith        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   else
871e884886fSBarry Smith .ve
872e884886fSBarry Smith 
873fe59aa6dSJacob Faibussowitsch .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatCreateSNESMF()`, `MatMFFDGetH()`, `MatCreateMFFD()`,
874db781477SPatrick Sanan           `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`
875e884886fSBarry Smith @*/
876d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetFunctionError(Mat mat, PetscReal error)
877d71ae5a4SJacob Faibussowitsch {
878e884886fSBarry Smith   PetscFunctionBegin;
87988b4c220SStefano Zampini   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
88088b4c220SStefano Zampini   PetscValidLogicalCollectiveReal(mat, error, 2);
881cac4c232SBarry Smith   PetscTryMethod(mat, "MatMFFDSetFunctionError_C", (Mat, PetscReal), (mat, error));
8823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
883e884886fSBarry Smith }
884e884886fSBarry Smith 
885e884886fSBarry Smith /*@
886e884886fSBarry Smith   MatMFFDSetHHistory - Sets an array to collect a history of the
88711a5261eSBarry Smith   differencing values (h) computed for the matrix-free product `MATMFFD` matrix
888e884886fSBarry Smith 
889c3339decSBarry Smith   Logically Collective
890e884886fSBarry Smith 
891e884886fSBarry Smith   Input Parameters:
89211a5261eSBarry Smith + J        - the `MATMFFD` matrix-free matrix
893a5b23f4aSJose E. Roman . history  - space to hold the history
894e884886fSBarry Smith - nhistory - number of entries in history, if more entries are generated than
895e884886fSBarry Smith               nhistory, then the later ones are discarded
896e884886fSBarry Smith 
897e884886fSBarry Smith   Level: advanced
898e884886fSBarry Smith 
89911a5261eSBarry Smith   Note:
90011a5261eSBarry Smith   Use `MatMFFDResetHHistory()` to reset the history counter and collect
901e884886fSBarry Smith   a new batch of differencing parameters, h.
902e884886fSBarry Smith 
9031cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatMFFDGetH()`, `MatCreateSNESMF()`,
904db781477SPatrick Sanan           `MatMFFDResetHHistory()`, `MatMFFDSetFunctionError()`
905e884886fSBarry Smith @*/
906d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetHHistory(Mat J, PetscScalar history[], PetscInt nhistory)
907d71ae5a4SJacob Faibussowitsch {
908e884886fSBarry Smith   PetscFunctionBegin;
90988b4c220SStefano Zampini   PetscValidHeaderSpecific(J, MAT_CLASSID, 1);
9104f572ea9SToby Isaac   if (history) PetscAssertPointer(history, 2);
91188b4c220SStefano Zampini   PetscValidLogicalCollectiveInt(J, nhistory, 3);
912cac4c232SBarry Smith   PetscUseMethod(J, "MatMFFDSetHHistory_C", (Mat, PetscScalar[], PetscInt), (J, history, nhistory));
9133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
914e884886fSBarry Smith }
915e884886fSBarry Smith 
916e884886fSBarry Smith /*@
917e884886fSBarry Smith   MatMFFDResetHHistory - Resets the counter to zero to begin
91811a5261eSBarry Smith   collecting a new set of differencing histories for the `MATMFFD` matrix
919e884886fSBarry Smith 
920c3339decSBarry Smith   Logically Collective
921e884886fSBarry Smith 
9222fe279fdSBarry Smith   Input Parameter:
923e884886fSBarry Smith . J - the matrix-free matrix context
924e884886fSBarry Smith 
925e884886fSBarry Smith   Level: advanced
926e884886fSBarry Smith 
92711a5261eSBarry Smith   Note:
92811a5261eSBarry Smith   Use `MatMFFDSetHHistory()` to create the original history counter.
929e884886fSBarry Smith 
9301cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatMFFDGetH()`, `MatCreateSNESMF()`,
931db781477SPatrick Sanan           `MatMFFDSetHHistory()`, `MatMFFDSetFunctionError()`
932e884886fSBarry Smith @*/
933d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDResetHHistory(Mat J)
934d71ae5a4SJacob Faibussowitsch {
935e884886fSBarry Smith   PetscFunctionBegin;
93688b4c220SStefano Zampini   PetscValidHeaderSpecific(J, MAT_CLASSID, 1);
937cac4c232SBarry Smith   PetscTryMethod(J, "MatMFFDResetHHistory_C", (Mat), (J));
9383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
939e884886fSBarry Smith }
940e884886fSBarry Smith 
941487a658cSBarry Smith /*@
9422ef1f0ffSBarry Smith   MatMFFDSetBase - Sets the vector `U` at which matrix vector products of the
94311a5261eSBarry Smith   Jacobian are computed for the `MATMFFD` matrix
944e884886fSBarry Smith 
945c3339decSBarry Smith   Logically Collective
946e884886fSBarry Smith 
947e884886fSBarry Smith   Input Parameters:
94811a5261eSBarry Smith + J - the `MATMFFD` matrix
9493ec795f1SBarry Smith . U - the vector
950bcddec3dSBarry Smith - F - (optional) vector that contains F(u) if it has been already computed
951e884886fSBarry Smith 
9522ef1f0ffSBarry Smith   Level: advanced
9532ef1f0ffSBarry Smith 
95495452b02SPatrick Sanan   Notes:
95595452b02SPatrick Sanan   This is rarely used directly
956e884886fSBarry Smith 
9572ef1f0ffSBarry Smith   If `F` is provided then it is not recomputed. Otherwise the function is evaluated at the base
95811a5261eSBarry Smith   point during the first `MatMult()` after each call to `MatMFFDSetBase()`.
959dff2f722SBarry Smith 
96042747ad1SJacob Faibussowitsch .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatMult()`
961e884886fSBarry Smith @*/
962d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetBase(Mat J, Vec U, Vec F)
963d71ae5a4SJacob Faibussowitsch {
964e884886fSBarry Smith   PetscFunctionBegin;
9650700a824SBarry Smith   PetscValidHeaderSpecific(J, MAT_CLASSID, 1);
9660700a824SBarry Smith   PetscValidHeaderSpecific(U, VEC_CLASSID, 2);
9670700a824SBarry Smith   if (F) PetscValidHeaderSpecific(F, VEC_CLASSID, 3);
968cac4c232SBarry Smith   PetscTryMethod(J, "MatMFFDSetBase_C", (Mat, Vec, Vec), (J, U, F));
9693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
970e884886fSBarry Smith }
971e884886fSBarry Smith 
972e884886fSBarry Smith /*@C
973e884886fSBarry Smith   MatMFFDSetCheckh - Sets a function that checks the computed h and adjusts
97411a5261eSBarry Smith   it to satisfy some criteria for the `MATMFFD` matrix
975e884886fSBarry Smith 
976c3339decSBarry Smith   Logically Collective
977e884886fSBarry Smith 
978e884886fSBarry Smith   Input Parameters:
97911a5261eSBarry Smith + J   - the `MATMFFD` matrix
9802ef1f0ffSBarry Smith . fun - the function that checks `h`
981e884886fSBarry Smith - ctx - any context needed by the function
982e884886fSBarry Smith 
983e884886fSBarry Smith   Options Database Keys:
98467b8a455SSatish Balay . -mat_mffd_check_positivity <bool> - Insure that U + h*a is non-negative
985e884886fSBarry Smith 
986e884886fSBarry Smith   Level: advanced
987e884886fSBarry Smith 
98895452b02SPatrick Sanan   Notes:
9892ef1f0ffSBarry Smith   For example, `MatMFFDCheckPositivity()` insures that all entries of U + h*a are non-negative
990e884886fSBarry Smith 
9912ef1f0ffSBarry Smith   The function you provide is called after the default `h` has been computed and allows you to
992755b7f64SBarry Smith   modify it.
993755b7f64SBarry Smith 
9941cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatMFFDCheckPositivity()`
995e884886fSBarry Smith @*/
996d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetCheckh(Mat J, PetscErrorCode (*fun)(void *, Vec, Vec, PetscScalar *), void *ctx)
997d71ae5a4SJacob Faibussowitsch {
998e884886fSBarry Smith   PetscFunctionBegin;
9990700a824SBarry Smith   PetscValidHeaderSpecific(J, MAT_CLASSID, 1);
1000cac4c232SBarry Smith   PetscTryMethod(J, "MatMFFDSetCheckh_C", (Mat, PetscErrorCode(*)(void *, Vec, Vec, PetscScalar *), void *), (J, fun, ctx));
10013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1002e884886fSBarry Smith }
1003e884886fSBarry Smith 
1004e884886fSBarry Smith /*@
1005e884886fSBarry Smith   MatMFFDCheckPositivity - Checks that all entries in U + h*a are positive or
100611a5261eSBarry Smith   zero, decreases h until this is satisfied for a `MATMFFD` matrix
1007e884886fSBarry Smith 
1008c3339decSBarry Smith   Logically Collective
1009e884886fSBarry Smith 
1010e884886fSBarry Smith   Input Parameters:
1011e884886fSBarry Smith + U     - base vector that is added to
1012e884886fSBarry Smith . a     - vector that is added
1013e884886fSBarry Smith . h     - scaling factor on a
1014e884886fSBarry Smith - dummy - context variable (unused)
1015e884886fSBarry Smith 
1016e884886fSBarry Smith   Options Database Keys:
101767b8a455SSatish Balay . -mat_mffd_check_positivity <bool> - Insure that U + h*a is nonnegative
1018e884886fSBarry Smith 
1019e884886fSBarry Smith   Level: advanced
1020e884886fSBarry Smith 
102111a5261eSBarry Smith   Note:
10222ef1f0ffSBarry Smith   This is rarely used directly, rather it is passed as an argument to `MatMFFDSetCheckh()`
1023e884886fSBarry Smith 
10241cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatMFFDSetCheckh()`
1025e884886fSBarry Smith @*/
1026d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDCheckPositivity(void *dummy, Vec U, Vec a, PetscScalar *h)
1027d71ae5a4SJacob Faibussowitsch {
1028e884886fSBarry Smith   PetscReal    val, minval;
1029e884886fSBarry Smith   PetscScalar *u_vec, *a_vec;
1030e884886fSBarry Smith   PetscInt     i, n;
1031e884886fSBarry Smith   MPI_Comm     comm;
1032e884886fSBarry Smith 
1033e884886fSBarry Smith   PetscFunctionBegin;
103488b4c220SStefano Zampini   PetscValidHeaderSpecific(U, VEC_CLASSID, 2);
103588b4c220SStefano Zampini   PetscValidHeaderSpecific(a, VEC_CLASSID, 3);
10364f572ea9SToby Isaac   PetscAssertPointer(h, 4);
10379566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)U, &comm));
10389566063dSJacob Faibussowitsch   PetscCall(VecGetArray(U, &u_vec));
10399566063dSJacob Faibussowitsch   PetscCall(VecGetArray(a, &a_vec));
10409566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(U, &n));
1041c068d9bbSLisandro Dalcin   minval = PetscAbsScalar(*h) * PetscRealConstant(1.01);
1042e884886fSBarry Smith   for (i = 0; i < n; i++) {
1043e884886fSBarry Smith     if (PetscRealPart(u_vec[i] + *h * a_vec[i]) <= 0.0) {
1044e884886fSBarry Smith       val = PetscAbsScalar(u_vec[i] / a_vec[i]);
1045e884886fSBarry Smith       if (val < minval) minval = val;
1046e884886fSBarry Smith     }
1047e884886fSBarry Smith   }
10489566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(U, &u_vec));
10499566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(a, &a_vec));
10501c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&minval, &val, 1, MPIU_REAL, MPIU_MIN, comm));
1051e884886fSBarry Smith   if (val <= PetscAbsScalar(*h)) {
10529566063dSJacob Faibussowitsch     PetscCall(PetscInfo(U, "Scaling back h from %g to %g\n", (double)PetscRealPart(*h), (double)(.99 * val)));
1053e884886fSBarry Smith     if (PetscRealPart(*h) > 0.0) *h = 0.99 * val;
1054e884886fSBarry Smith     else *h = -0.99 * val;
1055e884886fSBarry Smith   }
10563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1057e884886fSBarry Smith }
1058