xref: /petsc/src/mat/impls/mffd/mffd.c (revision 11a5261e40035b7c793f2783a2ba6c7cd4f3b077)
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
13*11a5261eSBarry Smith   MatMFFDFinalizePackage - This function destroys everything in the MATMFFD` package. It is
14*11a5261eSBarry Smith   called from `PetscFinalize()`.
15b022a5c1SBarry Smith 
16b022a5c1SBarry Smith   Level: developer
17b022a5c1SBarry Smith 
18*11a5261eSBarry Smith .seealso: `MATMFFD`, `PetscFinalize()`, `MatCreateMFFD()`, `MatCreateSNESMF()`
19b022a5c1SBarry Smith @*/
209371c9d4SSatish Balay PetscErrorCode   MatMFFDFinalizePackage(void) {
21b022a5c1SBarry Smith     PetscFunctionBegin;
229566063dSJacob Faibussowitsch     PetscCall(PetscFunctionListDestroy(&MatMFFDList));
23b022a5c1SBarry Smith     MatMFFDPackageInitialized = PETSC_FALSE;
24b022a5c1SBarry Smith     MatMFFDRegisterAllCalled  = PETSC_FALSE;
25b022a5c1SBarry Smith     PetscFunctionReturn(0);
26b022a5c1SBarry Smith }
27b022a5c1SBarry Smith 
283ec795f1SBarry Smith /*@C
29*11a5261eSBarry Smith   MatMFFDInitializePackage - This function initializes everything in the MATMFFD` package. It is called
30*11a5261eSBarry Smith   from `MatInitializePackage()`.
313ec795f1SBarry Smith 
323ec795f1SBarry Smith   Level: developer
333ec795f1SBarry Smith 
34*11a5261eSBarry Smith .seealso: `MATMFFD`, `PetscInitialize()`
353ec795f1SBarry Smith @*/
369371c9d4SSatish Balay PetscErrorCode MatMFFDInitializePackage(void) {
373ec795f1SBarry Smith   char      logList[256];
388e81d068SLisandro Dalcin   PetscBool opt, pkg;
393ec795f1SBarry Smith 
403ec795f1SBarry Smith   PetscFunctionBegin;
41b022a5c1SBarry Smith   if (MatMFFDPackageInitialized) PetscFunctionReturn(0);
42b022a5c1SBarry Smith   MatMFFDPackageInitialized = PETSC_TRUE;
433ec795f1SBarry Smith   /* Register Classes */
449566063dSJacob Faibussowitsch   PetscCall(PetscClassIdRegister("MatMFFD", &MATMFFD_CLASSID));
453ec795f1SBarry Smith   /* Register Constructors */
469566063dSJacob Faibussowitsch   PetscCall(MatMFFDRegisterAll());
473ec795f1SBarry Smith   /* Register Events */
489566063dSJacob Faibussowitsch   PetscCall(PetscLogEventRegister("MatMult MF", MATMFFD_CLASSID, &MATMFFD_Mult));
49e94e781bSJacob Faibussowitsch   /* Process Info */
50e94e781bSJacob Faibussowitsch   {
51e94e781bSJacob Faibussowitsch     PetscClassId classids[1];
52e94e781bSJacob Faibussowitsch 
53e94e781bSJacob Faibussowitsch     classids[0] = MATMFFD_CLASSID;
549566063dSJacob Faibussowitsch     PetscCall(PetscInfoProcessClass("matmffd", 1, classids));
553ec795f1SBarry Smith   }
563ec795f1SBarry Smith   /* Process summary exclusions */
579566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL, NULL, "-log_exclude", logList, sizeof(logList), &opt));
583ec795f1SBarry Smith   if (opt) {
599566063dSJacob Faibussowitsch     PetscCall(PetscStrInList("matmffd", logList, ',', &pkg));
609566063dSJacob Faibussowitsch     if (pkg) PetscCall(PetscLogEventExcludeClass(MATMFFD_CLASSID));
613ec795f1SBarry Smith   }
628e81d068SLisandro Dalcin   /* Register package finalizer */
639566063dSJacob Faibussowitsch   PetscCall(PetscRegisterFinalize(MatMFFDFinalizePackage));
643ec795f1SBarry Smith   PetscFunctionReturn(0);
653ec795f1SBarry Smith }
663ec795f1SBarry Smith 
679371c9d4SSatish Balay static PetscErrorCode MatMFFDSetType_MFFD(Mat mat, MatMFFDType ftype) {
68789d8953SBarry Smith   MatMFFD   ctx;
69789d8953SBarry Smith   PetscBool match;
705f80ce2aSJacob Faibussowitsch   PetscErrorCode (*r)(MatMFFD);
71789d8953SBarry Smith 
72789d8953SBarry Smith   PetscFunctionBegin;
73789d8953SBarry Smith   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
74789d8953SBarry Smith   PetscValidCharPointer(ftype, 2);
759566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
76789d8953SBarry Smith 
77789d8953SBarry Smith   /* already set, so just return */
789566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)ctx, ftype, &match));
79789d8953SBarry Smith   if (match) PetscFunctionReturn(0);
80789d8953SBarry Smith 
81789d8953SBarry Smith   /* destroy the old one if it exists */
82dbbe0bcdSBarry Smith   PetscTryTypeMethod(ctx, destroy);
83789d8953SBarry Smith 
849566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(MatMFFDList, ftype, &r));
855f80ce2aSJacob Faibussowitsch   PetscCheck(r, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown MatMFFD type %s given", ftype);
869566063dSJacob Faibussowitsch   PetscCall((*r)(ctx));
879566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)ctx, ftype));
88789d8953SBarry Smith   PetscFunctionReturn(0);
89789d8953SBarry Smith }
90789d8953SBarry Smith 
91e884886fSBarry Smith /*@C
92e884886fSBarry Smith     MatMFFDSetType - Sets the method that is used to compute the
93e884886fSBarry Smith     differencing parameter for finite differene matrix-free formulations.
94e884886fSBarry Smith 
95e884886fSBarry Smith     Input Parameters:
96*11a5261eSBarry Smith +   mat - the "matrix-free" matrix created via `MatCreateSNESMF()`, or `MatCreateMFFD()`
97*11a5261eSBarry Smith           or `MatSetType`(mat,`MATMFFD`);
98*11a5261eSBarry Smith -   ftype - the type requested, either `MATMFFD_WP` or `MATMFFD_DS`
99e884886fSBarry Smith 
100e884886fSBarry Smith     Level: advanced
101e884886fSBarry Smith 
102*11a5261eSBarry Smith     Note:
103e884886fSBarry Smith     For example, such routines can compute h for use in
104e884886fSBarry Smith     Jacobian-vector products of the form
105e884886fSBarry Smith 
106e884886fSBarry Smith                         F(x+ha) - F(x)
107e884886fSBarry Smith           F'(u)a  ~=  ----------------
108e884886fSBarry Smith                               h
109e884886fSBarry Smith 
110*11a5261eSBarry Smith .seealso: `MATMFFD`, `MATMFFD_WP`, `MATMFFD_DS`, `MatCreateSNESMF()`, `MatMFFDRegister()`, `MatMFFDSetFunction()`, `MatCreateMFFD()`
111e884886fSBarry Smith @*/
1129371c9d4SSatish Balay PetscErrorCode MatMFFDSetType(Mat mat, MatMFFDType ftype) {
113e884886fSBarry Smith   PetscFunctionBegin;
1140700a824SBarry Smith   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
115e884886fSBarry Smith   PetscValidCharPointer(ftype, 2);
116cac4c232SBarry Smith   PetscTryMethod(mat, "MatMFFDSetType_C", (Mat, MatMFFDType), (mat, ftype));
117e884886fSBarry Smith   PetscFunctionReturn(0);
118e884886fSBarry Smith }
119e884886fSBarry Smith 
1205d21e1f8SStefano Zampini static PetscErrorCode MatGetDiagonal_MFFD(Mat, Vec);
1215d21e1f8SStefano Zampini 
122e884886fSBarry Smith typedef PetscErrorCode (*FCN1)(void *, Vec); /* force argument to next function to not be extern C*/
1239371c9d4SSatish Balay static PetscErrorCode MatMFFDSetFunctioniBase_MFFD(Mat mat, FCN1 func) {
124789d8953SBarry Smith   MatMFFD ctx;
125e884886fSBarry Smith 
126e884886fSBarry Smith   PetscFunctionBegin;
1279566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
128e884886fSBarry Smith   ctx->funcisetbase = func;
129e884886fSBarry Smith   PetscFunctionReturn(0);
130e884886fSBarry Smith }
131e884886fSBarry Smith 
132e884886fSBarry Smith typedef PetscErrorCode (*FCN2)(void *, PetscInt, Vec, PetscScalar *); /* force argument to next function to not be extern C*/
1339371c9d4SSatish Balay static PetscErrorCode MatMFFDSetFunctioni_MFFD(Mat mat, FCN2 funci) {
134789d8953SBarry Smith   MatMFFD ctx;
135e884886fSBarry Smith 
136e884886fSBarry Smith   PetscFunctionBegin;
1379566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
138e884886fSBarry Smith   ctx->funci = funci;
1399566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(mat, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiagonal_MFFD));
140789d8953SBarry Smith   PetscFunctionReturn(0);
1415d21e1f8SStefano Zampini }
142789d8953SBarry Smith 
1439371c9d4SSatish Balay static PetscErrorCode MatMFFDGetH_MFFD(Mat mat, PetscScalar *h) {
144789d8953SBarry Smith   MatMFFD ctx;
145789d8953SBarry Smith 
146789d8953SBarry Smith   PetscFunctionBegin;
1479566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
148789d8953SBarry Smith   *h = ctx->currenth;
149e884886fSBarry Smith   PetscFunctionReturn(0);
150e884886fSBarry Smith }
151e884886fSBarry Smith 
1529371c9d4SSatish Balay static PetscErrorCode MatMFFDResetHHistory_MFFD(Mat J) {
153789d8953SBarry Smith   MatMFFD ctx;
154bc0f08ceSBarry Smith 
155bc0f08ceSBarry Smith   PetscFunctionBegin;
1569566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J, &ctx));
157bc0f08ceSBarry Smith   ctx->ncurrenth = 0;
158bc0f08ceSBarry Smith   PetscFunctionReturn(0);
159bc0f08ceSBarry Smith }
160e884886fSBarry Smith 
1611c84c290SBarry Smith /*@C
162*11a5261eSBarry Smith    MatMFFDRegister - Adds a method to the MATMFFD` registry.
1631c84c290SBarry Smith 
1641c84c290SBarry Smith    Not Collective
1651c84c290SBarry Smith 
1661c84c290SBarry Smith    Input Parameters:
1671c84c290SBarry Smith +  name_solver - name of a new user-defined compute-h module
1681c84c290SBarry Smith -  routine_create - routine to create method context
1691c84c290SBarry Smith 
1701c84c290SBarry Smith    Level: developer
1711c84c290SBarry Smith 
172*11a5261eSBarry Smith    Note:
173*11a5261eSBarry Smith    `MatMFFDRegister()` may be called multiple times to add several user-defined solvers.
1741c84c290SBarry Smith 
1751c84c290SBarry Smith    Sample usage:
1761c84c290SBarry Smith .vb
177bdf89e91SBarry Smith    MatMFFDRegister("my_h",MyHCreate);
1781c84c290SBarry Smith .ve
1791c84c290SBarry Smith 
1801c84c290SBarry Smith    Then, your solver can be chosen with the procedural interface via
181*11a5261eSBarry Smith $     `MatMFFDSetType`(mfctx,"my_h")
1821c84c290SBarry Smith    or at runtime via the option
183be7a6d03SBarry Smith $     -mat_mffd_type my_h
1841c84c290SBarry Smith 
185*11a5261eSBarry Smith .seealso: `MATMFFD`, `MatMFFDRegisterAll()`, `MatMFFDRegisterDestroy()`
1861c84c290SBarry Smith  @*/
1879371c9d4SSatish Balay PetscErrorCode MatMFFDRegister(const char sname[], PetscErrorCode (*function)(MatMFFD)) {
188e884886fSBarry Smith   PetscFunctionBegin;
1899566063dSJacob Faibussowitsch   PetscCall(MatInitializePackage());
1909566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&MatMFFDList, sname, function));
191e884886fSBarry Smith   PetscFunctionReturn(0);
192e884886fSBarry Smith }
193e884886fSBarry Smith 
194e884886fSBarry Smith /* ----------------------------------------------------------------------------------------*/
1959371c9d4SSatish Balay static PetscErrorCode MatDestroy_MFFD(Mat mat) {
196789d8953SBarry Smith   MatMFFD ctx;
197e884886fSBarry Smith 
198e884886fSBarry Smith   PetscFunctionBegin;
1999566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
2009566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ctx->w));
2019566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ctx->current_u));
20248a46eb9SPierre Jolivet   if (ctx->current_f_allocated) PetscCall(VecDestroy(&ctx->current_f));
203dbbe0bcdSBarry Smith   PetscTryTypeMethod(ctx, destroy);
2049566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(&ctx));
205e884886fSBarry Smith 
2069566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetBase_C", NULL));
2079566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetFunctioniBase_C", NULL));
2089566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetFunctioni_C", NULL));
2099566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetFunction_C", NULL));
2109566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetFunctionError_C", NULL));
2119566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetCheckh_C", NULL));
2129566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetPeriod_C", NULL));
2139566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDResetHHistory_C", NULL));
2149566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetHHistory_C", NULL));
2159566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetType_C", NULL));
2169566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDGetH_C", NULL));
2172e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatSNESMFSetReuseBase_C", NULL));
2182e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatSNESMFGetReuseBase_C", NULL));
219e884886fSBarry Smith   PetscFunctionReturn(0);
220e884886fSBarry Smith }
221e884886fSBarry Smith 
222e884886fSBarry Smith /*
223e884886fSBarry Smith    MatMFFDView_MFFD - Views matrix-free parameters.
224e884886fSBarry Smith 
225e884886fSBarry Smith */
2269371c9d4SSatish Balay static PetscErrorCode MatView_MFFD(Mat J, PetscViewer viewer) {
227789d8953SBarry Smith   MatMFFD     ctx;
228a283635eSDmitry Karpeev   PetscBool   iascii, viewbase, viewfunction;
229a283635eSDmitry Karpeev   const char *prefix;
230e884886fSBarry Smith 
231e884886fSBarry Smith   PetscFunctionBegin;
2329566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J, &ctx));
2339566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
234e884886fSBarry Smith   if (iascii) {
2359566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Matrix-free approximation:\n"));
2369566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
2379566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "err=%g (relative error in function evaluation)\n", (double)ctx->error_rel));
2387adad957SLisandro Dalcin     if (!((PetscObject)ctx)->type_name) {
2399566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "The compute h routine has not yet been set\n"));
240e884886fSBarry Smith     } else {
2419566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Using %s compute h routine\n", ((PetscObject)ctx)->type_name));
242e884886fSBarry Smith     }
243c92e3469SBarry Smith #if defined(PETSC_USE_COMPLEX)
24448a46eb9SPierre Jolivet     if (ctx->usecomplex) PetscCall(PetscViewerASCIIPrintf(viewer, "Using Lyness complex number trick to compute the matrix-vector product\n"));
245c92e3469SBarry Smith #endif
246dbbe0bcdSBarry Smith     PetscTryTypeMethod(ctx, view, viewer);
2479566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetOptionsPrefix((PetscObject)J, &prefix));
248a283635eSDmitry Karpeev 
2499566063dSJacob Faibussowitsch     PetscCall(PetscOptionsHasName(((PetscObject)J)->options, prefix, "-mat_mffd_view_base", &viewbase));
250a283635eSDmitry Karpeev     if (viewbase) {
2519566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Base:\n"));
2529566063dSJacob Faibussowitsch       PetscCall(VecView(ctx->current_u, viewer));
253a283635eSDmitry Karpeev     }
2549566063dSJacob Faibussowitsch     PetscCall(PetscOptionsHasName(((PetscObject)J)->options, prefix, "-mat_mffd_view_function", &viewfunction));
255a283635eSDmitry Karpeev     if (viewfunction) {
2569566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Function:\n"));
2579566063dSJacob Faibussowitsch       PetscCall(VecView(ctx->current_f, viewer));
258a283635eSDmitry Karpeev     }
2599566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
26011aeaf0aSBarry Smith   }
261e884886fSBarry Smith   PetscFunctionReturn(0);
262e884886fSBarry Smith }
263e884886fSBarry Smith 
264e884886fSBarry Smith /*
265e884886fSBarry Smith    MatAssemblyEnd_MFFD - Resets the ctx->ncurrenth to zero. This
266e884886fSBarry Smith    allows the user to indicate the beginning of a new linear solve by calling
267e884886fSBarry Smith    MatAssemblyXXX() on the matrix free matrix. This then allows the
2681d0fab5eSBarry Smith    MatCreateMFFD_WP() to properly compute ||U|| only the first time
269e884886fSBarry Smith    in the linear solver rather than every time.
2705a576424SJed Brown 
271cc2e6a90SBarry Smith    This function is referenced directly from MatAssemblyEnd_SNESMF(), which may be in a different shared library hence
272cc2e6a90SBarry Smith    it must be labeled as PETSC_EXTERN
273e884886fSBarry Smith */
2749371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode MatAssemblyEnd_MFFD(Mat J, MatAssemblyType mt) {
275789d8953SBarry Smith   MatMFFD j;
276e884886fSBarry Smith 
277e884886fSBarry Smith   PetscFunctionBegin;
2789566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J, &j));
2799566063dSJacob Faibussowitsch   PetscCall(MatMFFDResetHHistory(J));
280e884886fSBarry Smith   PetscFunctionReturn(0);
281e884886fSBarry Smith }
282e884886fSBarry Smith 
283e884886fSBarry Smith /*
284e884886fSBarry Smith   MatMult_MFFD - Default matrix-free form for Jacobian-vector product, y = F'(u)*a:
285e884886fSBarry Smith 
286e884886fSBarry Smith         y ~= (F(u + ha) - F(u))/h,
287e884886fSBarry Smith   where F = nonlinear function, as set by SNESSetFunction()
288e884886fSBarry Smith         u = current iterate
289e884886fSBarry Smith         h = difference interval
290e884886fSBarry Smith */
2919371c9d4SSatish Balay static PetscErrorCode MatMult_MFFD(Mat mat, Vec a, Vec y) {
292789d8953SBarry Smith   MatMFFD     ctx;
293e884886fSBarry Smith   PetscScalar h;
294e884886fSBarry Smith   Vec         w, U, F;
295ace3abfcSBarry Smith   PetscBool   zeroa;
296e884886fSBarry Smith 
297e884886fSBarry Smith   PetscFunctionBegin;
2989566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
29928b400f6SJacob 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");
300e884886fSBarry Smith   /* We log matrix-free matrix-vector products separately, so that we can
301e884886fSBarry Smith      separate the performance monitoring from the cases that use conventional
302e884886fSBarry Smith      storage.  We may eventually modify event logging to associate events
303e884886fSBarry Smith      with particular objects, hence alleviating the more general problem. */
3049566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(MATMFFD_Mult, a, y, 0, 0));
305e884886fSBarry Smith 
306e884886fSBarry Smith   w = ctx->w;
307e884886fSBarry Smith   U = ctx->current_u;
3083ec795f1SBarry Smith   F = ctx->current_f;
309e884886fSBarry Smith   /*
310e884886fSBarry Smith       Compute differencing parameter
311e884886fSBarry Smith   */
3124a2f8832SBarry Smith   if (!((PetscObject)ctx)->type_name) {
3139566063dSJacob Faibussowitsch     PetscCall(MatMFFDSetType(mat, MATMFFD_WP));
3149566063dSJacob Faibussowitsch     PetscCall(MatSetFromOptions(mat));
315e884886fSBarry Smith   }
316dbbe0bcdSBarry Smith   PetscUseTypeMethod(ctx, compute, U, a, &h, &zeroa);
317e884886fSBarry Smith   if (zeroa) {
3189566063dSJacob Faibussowitsch     PetscCall(VecSet(y, 0.0));
31985d14658SLisandro Dalcin     PetscCall(PetscLogEventEnd(MATMFFD_Mult, a, y, 0, 0));
320e884886fSBarry Smith     PetscFunctionReturn(0);
321e884886fSBarry Smith   }
322e884886fSBarry Smith 
32308401ef6SPierre Jolivet   PetscCheck(!mat->erroriffailure || !PetscIsInfOrNanScalar(h), PETSC_COMM_SELF, PETSC_ERR_PLIB, "Computed Nan differencing parameter h");
32448a46eb9SPierre Jolivet   if (ctx->checkh) PetscCall((*ctx->checkh)(ctx->checkhctx, U, a, &h));
325e884886fSBarry Smith 
326e884886fSBarry Smith   /* keep a record of the current differencing parameter h */
327e884886fSBarry Smith   ctx->currenth = h;
328e884886fSBarry Smith #if defined(PETSC_USE_COMPLEX)
3299566063dSJacob Faibussowitsch   PetscCall(PetscInfo(mat, "Current differencing parameter: %g + %g i\n", (double)PetscRealPart(h), (double)PetscImaginaryPart(h)));
330e884886fSBarry Smith #else
3319566063dSJacob Faibussowitsch   PetscCall(PetscInfo(mat, "Current differencing parameter: %15.12e\n", (double)PetscRealPart(h)));
332e884886fSBarry Smith #endif
333ad540459SPierre Jolivet   if (ctx->historyh && ctx->ncurrenth < ctx->maxcurrenth) ctx->historyh[ctx->ncurrenth] = h;
334e884886fSBarry Smith   ctx->ncurrenth++;
335e884886fSBarry Smith 
336c92e3469SBarry Smith #if defined(PETSC_USE_COMPLEX)
337c92e3469SBarry Smith   if (ctx->usecomplex) h = PETSC_i * h;
338c92e3469SBarry Smith #endif
339c92e3469SBarry Smith 
340e884886fSBarry Smith   /* w = u + ha */
3419566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(w, h, a, U));
342e884886fSBarry Smith 
3431b797266SDmitry Karpeev   /* compute func(U) as base for differencing; only needed first time in and not when provided by user */
34448a46eb9SPierre Jolivet   if (ctx->ncurrenth == 1 && ctx->current_f_allocated) PetscCall((*ctx->func)(ctx->funcctx, U, F));
3459566063dSJacob Faibussowitsch   PetscCall((*ctx->func)(ctx->funcctx, w, y));
346e884886fSBarry Smith 
347c92e3469SBarry Smith #if defined(PETSC_USE_COMPLEX)
348c92e3469SBarry Smith   if (ctx->usecomplex) {
3499566063dSJacob Faibussowitsch     PetscCall(VecImaginaryPart(y));
350c92e3469SBarry Smith     h = PetscImaginaryPart(h);
351c92e3469SBarry Smith   } else {
3529566063dSJacob Faibussowitsch     PetscCall(VecAXPY(y, -1.0, F));
353c92e3469SBarry Smith   }
354c92e3469SBarry Smith #else
3559566063dSJacob Faibussowitsch   PetscCall(VecAXPY(y, -1.0, F));
356c92e3469SBarry Smith #endif
3579566063dSJacob Faibussowitsch   PetscCall(VecScale(y, 1.0 / h));
3589566063dSJacob Faibussowitsch   if (mat->nullsp) PetscCall(MatNullSpaceRemove(mat->nullsp, y));
359e884886fSBarry Smith 
3609566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(MATMFFD_Mult, a, y, 0, 0));
361e884886fSBarry Smith   PetscFunctionReturn(0);
362e884886fSBarry Smith }
363e884886fSBarry Smith 
364e884886fSBarry Smith /*
365e884886fSBarry Smith   MatGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix
366e884886fSBarry Smith 
367e884886fSBarry Smith         y ~= (F(u + ha) - F(u))/h,
368e884886fSBarry Smith   where F = nonlinear function, as set by SNESSetFunction()
369e884886fSBarry Smith         u = current iterate
370e884886fSBarry Smith         h = difference interval
371e884886fSBarry Smith */
3729371c9d4SSatish Balay PetscErrorCode MatGetDiagonal_MFFD(Mat mat, Vec a) {
373789d8953SBarry Smith   MatMFFD     ctx;
374e884886fSBarry Smith   PetscScalar h, *aa, *ww, v;
375e884886fSBarry Smith   PetscReal   epsilon = PETSC_SQRT_MACHINE_EPSILON, umin = 100.0 * PETSC_SQRT_MACHINE_EPSILON;
376e884886fSBarry Smith   Vec         w, U;
377e884886fSBarry Smith   PetscInt    i, rstart, rend;
378e884886fSBarry Smith 
379e884886fSBarry Smith   PetscFunctionBegin;
3809566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
38128b400f6SJacob Faibussowitsch   PetscCheck(ctx->func, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Requires calling MatMFFDSetFunction() first");
38228b400f6SJacob Faibussowitsch   PetscCheck(ctx->funci, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Requires calling MatMFFDSetFunctioni() first");
383e884886fSBarry Smith   w = ctx->w;
384e884886fSBarry Smith   U = ctx->current_u;
3859566063dSJacob Faibussowitsch   PetscCall((*ctx->func)(ctx->funcctx, U, a));
3861baa6e33SBarry Smith   if (ctx->funcisetbase) PetscCall((*ctx->funcisetbase)(ctx->funcctx, U));
3879566063dSJacob Faibussowitsch   PetscCall(VecCopy(U, w));
388e884886fSBarry Smith 
3899566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(a, &rstart, &rend));
3909566063dSJacob Faibussowitsch   PetscCall(VecGetArray(a, &aa));
391e884886fSBarry Smith   for (i = rstart; i < rend; i++) {
3929566063dSJacob Faibussowitsch     PetscCall(VecGetArray(w, &ww));
393e884886fSBarry Smith     h = ww[i - rstart];
394e884886fSBarry Smith     if (h == 0.0) h = 1.0;
395e884886fSBarry Smith     if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0) h = umin;
396e884886fSBarry Smith     else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin;
397e884886fSBarry Smith     h *= epsilon;
398e884886fSBarry Smith 
399e884886fSBarry Smith     ww[i - rstart] += h;
4009566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(w, &ww));
4019566063dSJacob Faibussowitsch     PetscCall((*ctx->funci)(ctx->funcctx, i, w, &v));
402e884886fSBarry Smith     aa[i - rstart] = (v - aa[i - rstart]) / h;
403e884886fSBarry Smith 
4049566063dSJacob Faibussowitsch     PetscCall(VecGetArray(w, &ww));
405e884886fSBarry Smith     ww[i - rstart] -= h;
4069566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(w, &ww));
407e884886fSBarry Smith   }
4089566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(a, &aa));
409e884886fSBarry Smith   PetscFunctionReturn(0);
410e884886fSBarry Smith }
411e884886fSBarry Smith 
4129371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode MatMFFDSetBase_MFFD(Mat J, Vec U, Vec F) {
413789d8953SBarry Smith   MatMFFD ctx;
414e884886fSBarry Smith 
415e884886fSBarry Smith   PetscFunctionBegin;
4169566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J, &ctx));
4179566063dSJacob Faibussowitsch   PetscCall(MatMFFDResetHHistory(J));
418c51ad4d4SStefano Zampini   if (!ctx->current_u) {
4199566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(U, &ctx->current_u));
4209566063dSJacob Faibussowitsch     PetscCall(VecLockReadPush(ctx->current_u));
421c51ad4d4SStefano Zampini   }
4229566063dSJacob Faibussowitsch   PetscCall(VecLockReadPop(ctx->current_u));
4239566063dSJacob Faibussowitsch   PetscCall(VecCopy(U, ctx->current_u));
4249566063dSJacob Faibussowitsch   PetscCall(VecLockReadPush(ctx->current_u));
42552121784SBarry Smith   if (F) {
4269566063dSJacob Faibussowitsch     if (ctx->current_f_allocated) PetscCall(VecDestroy(&ctx->current_f));
4273ec795f1SBarry Smith     ctx->current_f           = F;
428cfe22f5eSBarry Smith     ctx->current_f_allocated = PETSC_FALSE;
429cfe22f5eSBarry Smith   } else if (!ctx->current_f_allocated) {
4309566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(J, NULL, &ctx->current_f));
431cfe22f5eSBarry Smith     ctx->current_f_allocated = PETSC_TRUE;
43252121784SBarry Smith   }
43348a46eb9SPierre Jolivet   if (!ctx->w) PetscCall(VecDuplicate(ctx->current_u, &ctx->w));
434e884886fSBarry Smith   J->assembled = PETSC_TRUE;
435e884886fSBarry Smith   PetscFunctionReturn(0);
436e884886fSBarry Smith }
4375a576424SJed Brown 
438e884886fSBarry Smith typedef PetscErrorCode (*FCN3)(void *, Vec, Vec, PetscScalar *); /* force argument to next function to not be extern C*/
439b2573a8aSBarry Smith 
4409371c9d4SSatish Balay static PetscErrorCode MatMFFDSetCheckh_MFFD(Mat J, FCN3 fun, void *ectx) {
441789d8953SBarry Smith   MatMFFD ctx;
442e884886fSBarry Smith 
443e884886fSBarry Smith   PetscFunctionBegin;
4449566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J, &ctx));
445e884886fSBarry Smith   ctx->checkh    = fun;
446e884886fSBarry Smith   ctx->checkhctx = ectx;
447e884886fSBarry Smith   PetscFunctionReturn(0);
448e884886fSBarry Smith }
449e884886fSBarry Smith 
4506aa9148fSLisandro Dalcin /*@C
4516aa9148fSLisandro Dalcin    MatMFFDSetOptionsPrefix - Sets the prefix used for searching for all
452*11a5261eSBarry Smith    MATMFFD` options in the database.
4536aa9148fSLisandro Dalcin 
454*11a5261eSBarry Smith    Collective on mat
4556aa9148fSLisandro Dalcin 
456d8d19677SJose E. Roman    Input Parameters:
457*11a5261eSBarry Smith +  A - the `MATMFFD` context
4586aa9148fSLisandro Dalcin -  prefix - the prefix to prepend to all option names
4596aa9148fSLisandro Dalcin 
460*11a5261eSBarry Smith    Note:
4616aa9148fSLisandro Dalcin    A hyphen (-) must NOT be given at the beginning of the prefix name.
4626aa9148fSLisandro Dalcin    The first character of all runtime options is AUTOMATICALLY the hyphen.
4636aa9148fSLisandro Dalcin 
4646aa9148fSLisandro Dalcin    Level: advanced
4656aa9148fSLisandro Dalcin 
466*11a5261eSBarry Smith .seealso: `MATMFFD`, `MatSetFromOptions()`, `MatCreateSNESMF()`, `MatCreateMFFD()`
4676aa9148fSLisandro Dalcin @*/
4689371c9d4SSatish Balay PetscErrorCode MatMFFDSetOptionsPrefix(Mat mat, const char prefix[]) {
469789d8953SBarry Smith   MatMFFD mfctx;
4705fd66863SKarl Rupp 
4716aa9148fSLisandro Dalcin   PetscFunctionBegin;
4720700a824SBarry Smith   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
4739566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &mfctx));
4740700a824SBarry Smith   PetscValidHeaderSpecific(mfctx, MATMFFD_CLASSID, 1);
4759566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)mfctx, prefix));
4766aa9148fSLisandro Dalcin   PetscFunctionReturn(0);
4776aa9148fSLisandro Dalcin }
4786aa9148fSLisandro Dalcin 
4799371c9d4SSatish Balay static PetscErrorCode MatSetFromOptions_MFFD(Mat mat, PetscOptionItems *PetscOptionsObject) {
480789d8953SBarry Smith   MatMFFD   mfctx;
481ace3abfcSBarry Smith   PetscBool flg;
482e884886fSBarry Smith   char      ftype[256];
483e884886fSBarry Smith 
484e884886fSBarry Smith   PetscFunctionBegin;
4859566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &mfctx));
486dbbe0bcdSBarry Smith   PetscValidHeaderSpecific(mfctx, MATMFFD_CLASSID, 1);
487d0609cedSBarry Smith   PetscObjectOptionsBegin((PetscObject)mfctx);
4889566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFList("-mat_mffd_type", "Matrix free type", "MatMFFDSetType", MatMFFDList, ((PetscObject)mfctx)->type_name, ftype, 256, &flg));
4891baa6e33SBarry Smith   if (flg) PetscCall(MatMFFDSetType(mat, ftype));
490e884886fSBarry Smith 
4919566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_mffd_err", "set sqrt relative error in function", "MatMFFDSetFunctionError", mfctx->error_rel, &mfctx->error_rel, NULL));
4929566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-mat_mffd_period", "how often h is recomputed", "MatMFFDSetPeriod", mfctx->recomputeperiod, &mfctx->recomputeperiod, NULL));
493e884886fSBarry Smith 
49490d69ab7SBarry Smith   flg = PETSC_FALSE;
4959566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-mat_mffd_check_positivity", "Insure that U + h*a is nonnegative", "MatMFFDSetCheckh", flg, &flg, NULL));
4961baa6e33SBarry Smith   if (flg) PetscCall(MatMFFDSetCheckh(mat, MatMFFDCheckPositivity, NULL));
497c92e3469SBarry Smith #if defined(PETSC_USE_COMPLEX)
4989566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-mat_mffd_complex", "Use Lyness complex number trick to compute the matrix-vector product", "None", mfctx->usecomplex, &mfctx->usecomplex, NULL));
499c92e3469SBarry Smith #endif
500dbbe0bcdSBarry Smith   PetscTryTypeMethod(mfctx, setfromoptions, PetscOptionsObject);
501d0609cedSBarry Smith   PetscOptionsEnd();
502e884886fSBarry Smith   PetscFunctionReturn(0);
503e884886fSBarry Smith }
504e884886fSBarry Smith 
5059371c9d4SSatish Balay static PetscErrorCode MatMFFDSetPeriod_MFFD(Mat mat, PetscInt period) {
506789d8953SBarry Smith   MatMFFD ctx;
507bc0f08ceSBarry Smith 
508bc0f08ceSBarry Smith   PetscFunctionBegin;
5099566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
510bc0f08ceSBarry Smith   ctx->recomputeperiod = period;
511bc0f08ceSBarry Smith   PetscFunctionReturn(0);
512bc0f08ceSBarry Smith }
513bc0f08ceSBarry Smith 
5149371c9d4SSatish Balay static PetscErrorCode MatMFFDSetFunction_MFFD(Mat mat, PetscErrorCode (*func)(void *, Vec, Vec), void *funcctx) {
515789d8953SBarry Smith   MatMFFD ctx;
516bc0f08ceSBarry Smith 
517bc0f08ceSBarry Smith   PetscFunctionBegin;
5189566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
519bc0f08ceSBarry Smith   ctx->func    = func;
520bc0f08ceSBarry Smith   ctx->funcctx = funcctx;
521bc0f08ceSBarry Smith   PetscFunctionReturn(0);
522bc0f08ceSBarry Smith }
523bc0f08ceSBarry Smith 
5249371c9d4SSatish Balay static PetscErrorCode MatMFFDSetFunctionError_MFFD(Mat mat, PetscReal error) {
525789d8953SBarry Smith   MatMFFD ctx;
526bc0f08ceSBarry Smith 
527bc0f08ceSBarry Smith   PetscFunctionBegin;
5289566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
529bc0f08ceSBarry Smith   if (error != PETSC_DEFAULT) ctx->error_rel = error;
530bc0f08ceSBarry Smith   PetscFunctionReturn(0);
531bc0f08ceSBarry Smith }
532bc0f08ceSBarry Smith 
5339371c9d4SSatish Balay PetscErrorCode MatMFFDSetHHistory_MFFD(Mat J, PetscScalar history[], PetscInt nhistory) {
534789d8953SBarry Smith   MatMFFD ctx;
535789d8953SBarry Smith 
5363b49f96aSBarry Smith   PetscFunctionBegin;
5379566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J, &ctx));
538789d8953SBarry Smith   ctx->historyh    = history;
539789d8953SBarry Smith   ctx->maxcurrenth = nhistory;
540789d8953SBarry Smith   ctx->currenth    = 0.;
5413b49f96aSBarry Smith   PetscFunctionReturn(0);
5423b49f96aSBarry Smith }
5433b49f96aSBarry Smith 
544e884886fSBarry Smith /*MC
545e884886fSBarry Smith   MATMFFD - MATMFFD = "mffd" - A matrix free matrix type.
546e884886fSBarry Smith 
547e884886fSBarry Smith   Level: advanced
548e884886fSBarry Smith 
549*11a5261eSBarry Smith   Developers Note: This is implemented on top of `MATSHELL` to get support for scaling and shifting without requiring duplicate code
550789d8953SBarry Smith 
551db781477SPatrick Sanan .seealso: `MatCreateMFFD()`, `MatCreateSNESMF()`, `MatMFFDSetFunction()`, `MatMFFDSetType()`,
552db781477SPatrick Sanan           `MatMFFDSetFunctionError()`, `MatMFFDDSSetUmin()`, `MatMFFDSetFunction()`
553db781477SPatrick Sanan           `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`, `MatCreateSNESMF()`,
554db781477SPatrick Sanan           `MatMFFDGetH()`,
555e884886fSBarry Smith M*/
5569371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode MatCreate_MFFD(Mat A) {
557e884886fSBarry Smith   MatMFFD mfctx;
558e884886fSBarry Smith 
559e884886fSBarry Smith   PetscFunctionBegin;
5609566063dSJacob Faibussowitsch   PetscCall(MatMFFDInitializePackage());
561e884886fSBarry Smith 
5629566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(mfctx, MATMFFD_CLASSID, "MatMFFD", "Matrix-free Finite Differencing", "Mat", PetscObjectComm((PetscObject)A), NULL, NULL));
5632205254eSKarl Rupp 
564e884886fSBarry Smith   mfctx->error_rel                = PETSC_SQRT_MACHINE_EPSILON;
565e884886fSBarry Smith   mfctx->recomputeperiod          = 1;
566e884886fSBarry Smith   mfctx->count                    = 0;
567e884886fSBarry Smith   mfctx->currenth                 = 0.0;
5680298fd71SBarry Smith   mfctx->historyh                 = NULL;
569e884886fSBarry Smith   mfctx->ncurrenth                = 0;
570e884886fSBarry Smith   mfctx->maxcurrenth              = 0;
571f4259b30SLisandro Dalcin   ((PetscObject)mfctx)->type_name = NULL;
572e884886fSBarry Smith 
573e884886fSBarry Smith   /*
574e884886fSBarry Smith      Create the empty data structure to contain compute-h routines.
575e884886fSBarry Smith      These will be filled in below from the command line options or
576e884886fSBarry Smith      a later call with MatMFFDSetType() or if that is not called
577e884886fSBarry Smith      then it will default in the first use of MatMult_MFFD()
578e884886fSBarry Smith   */
579f4259b30SLisandro Dalcin   mfctx->ops->compute        = NULL;
580f4259b30SLisandro Dalcin   mfctx->ops->destroy        = NULL;
581f4259b30SLisandro Dalcin   mfctx->ops->view           = NULL;
582f4259b30SLisandro Dalcin   mfctx->ops->setfromoptions = NULL;
583f4259b30SLisandro Dalcin   mfctx->hctx                = NULL;
584e884886fSBarry Smith 
585f4259b30SLisandro Dalcin   mfctx->func    = NULL;
586f4259b30SLisandro Dalcin   mfctx->funcctx = NULL;
5870298fd71SBarry Smith   mfctx->w       = NULL;
588789d8953SBarry Smith   mfctx->mat     = A;
589e884886fSBarry Smith 
5909566063dSJacob Faibussowitsch   PetscCall(MatSetType(A, MATSHELL));
5919566063dSJacob Faibussowitsch   PetscCall(MatShellSetContext(A, mfctx));
5929566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(A, MATOP_MULT, (void (*)(void))MatMult_MFFD));
5939566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(A, MATOP_DESTROY, (void (*)(void))MatDestroy_MFFD));
5949566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(A, MATOP_VIEW, (void (*)(void))MatView_MFFD));
5959566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(A, MATOP_ASSEMBLY_END, (void (*)(void))MatAssemblyEnd_MFFD));
5969566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(A, MATOP_SET_FROM_OPTIONS, (void (*)(void))MatSetFromOptions_MFFD));
597e884886fSBarry Smith 
5989566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetBase_C", MatMFFDSetBase_MFFD));
5999566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetFunctioniBase_C", MatMFFDSetFunctioniBase_MFFD));
6009566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetFunctioni_C", MatMFFDSetFunctioni_MFFD));
6019566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetFunction_C", MatMFFDSetFunction_MFFD));
6029566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetCheckh_C", MatMFFDSetCheckh_MFFD));
6039566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetPeriod_C", MatMFFDSetPeriod_MFFD));
6049566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetFunctionError_C", MatMFFDSetFunctionError_MFFD));
6059566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDResetHHistory_C", MatMFFDResetHHistory_MFFD));
6069566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetHHistory_C", MatMFFDSetHHistory_MFFD));
6079566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetType_C", MatMFFDSetType_MFFD));
6089566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDGetH_C", MatMFFDGetH_MFFD));
6099566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATMFFD));
610e884886fSBarry Smith   PetscFunctionReturn(0);
611e884886fSBarry Smith }
612e884886fSBarry Smith 
613e884886fSBarry Smith /*@
614*11a5261eSBarry Smith    MatCreateMFFD - Creates a matrix-free matrix of type `MATMFFD`. See also `MatCreateSNESMF()`
615e884886fSBarry Smith 
616*11a5261eSBarry Smith    Collective
617e884886fSBarry Smith 
618e884886fSBarry Smith    Input Parameters:
619fef1beadSBarry Smith +  comm - MPI communicator
620*11a5261eSBarry Smith .  m - number of local rows (or `PETSC_DECIDE` to have calculated if M is given)
621fef1beadSBarry Smith            This value should be the same as the local size used in creating the
622fef1beadSBarry Smith            y vector for the matrix-vector product y = Ax.
623fef1beadSBarry Smith .  n - This value should be the same as the local size used in creating the
624*11a5261eSBarry Smith        x vector for the matrix-vector product y = Ax. (or `PETSC_DECIDE` to have
625fef1beadSBarry Smith        calculated if N is given) For square matrices n is almost always m.
626*11a5261eSBarry Smith .  M - number of global rows (or `PETSC_DETERMINE` to have calculated if m is given)
627*11a5261eSBarry Smith -  N - number of global columns (or `PETSC_DETERMINE` to have calculated if n is given)
628fef1beadSBarry Smith 
629e884886fSBarry Smith    Output Parameter:
630e884886fSBarry Smith .  J - the matrix-free matrix
631e884886fSBarry Smith 
632*11a5261eSBarry Smith    Options Database Keys:
633*11a5261eSBarry Smith +  -mat_mffd_type - wp or ds (see `MATMFFD_WP` or `MATMFFD_DS`)
634c92e3469SBarry Smith .  -mat_mffd_err - square root of estimated relative error in function evaluation
635c92e3469SBarry Smith .  -mat_mffd_period - how often h is recomputed, defaults to 1, every time
636c92e3469SBarry Smith .  -mat_mffd_check_positivity - possibly decrease h until U + h*a has only positive values
637*11a5261eSBarry Smith .  -mat_mffd_umin <umin> - Sets umin (for default PETSc routine that computes h only)
638c92e3469SBarry Smith -  -mat_mffd_complex - use the Lyness trick with complex numbers to compute the matrix-vector product instead of differencing
639c92e3469SBarry Smith                        (requires real valued functions but that PETSc be configured for complex numbers)
6409c6ac3b3SBarry Smith 
641e884886fSBarry Smith    Level: advanced
642e884886fSBarry Smith 
643e884886fSBarry Smith    Notes:
644e884886fSBarry Smith    The matrix-free matrix context merely contains the function pointers
645e884886fSBarry Smith    and work space for performing finite difference approximations of
646e884886fSBarry Smith    Jacobian-vector products, F'(u)*a,
647e884886fSBarry Smith 
648e884886fSBarry Smith    The default code uses the following approach to compute h
649e884886fSBarry Smith 
650e884886fSBarry Smith .vb
651e884886fSBarry Smith      F'(u)*a = [F(u+h*a) - F(u)]/h where
652e884886fSBarry Smith      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
653e884886fSBarry Smith        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
654e884886fSBarry Smith  where
655e884886fSBarry Smith      error_rel = square root of relative error in function evaluation
656e884886fSBarry Smith      umin = minimum iterate parameter
657e884886fSBarry Smith .ve
658e884886fSBarry Smith 
659*11a5261eSBarry Smith    You can call `SNESSetJacobian()` with `MatMFFDComputeJacobian()` if you are using matrix and not a different
660ca93e954SBarry Smith    preconditioner matrix
661ca93e954SBarry Smith 
662*11a5261eSBarry Smith    The user can set the error_rel via `MatMFFDSetFunctionError()` and
663*11a5261eSBarry Smith    umin via `MatMFFDDSSetUmin()`; see Users-Manual: ch_snes for details.
664e884886fSBarry Smith 
665*11a5261eSBarry Smith    The user should call `MatDestroy()` when finished with the matrix-free
666e884886fSBarry Smith    matrix context.
667e884886fSBarry Smith 
668*11a5261eSBarry Smith .seealso: `MATMFFD`, `MatDestroy()`, `MatMFFDSetFunctionError()`, `MatMFFDDSSetUmin()`, `MatMFFDSetFunction()`
669db781477SPatrick Sanan           `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`, `MatCreateSNESMF()`,
670db781477SPatrick Sanan           `MatMFFDGetH()`, `MatMFFDRegister()`, `MatMFFDComputeJacobian()`
671e884886fSBarry Smith @*/
6729371c9d4SSatish Balay PetscErrorCode MatCreateMFFD(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt M, PetscInt N, Mat *J) {
673e884886fSBarry Smith   PetscFunctionBegin;
6749566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, J));
6759566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*J, m, n, M, N));
6769566063dSJacob Faibussowitsch   PetscCall(MatSetType(*J, MATMFFD));
6779566063dSJacob Faibussowitsch   PetscCall(MatSetUp(*J));
678e884886fSBarry Smith   PetscFunctionReturn(0);
679e884886fSBarry Smith }
680e884886fSBarry Smith 
681e884886fSBarry Smith /*@
682*11a5261eSBarry Smith    MatMFFDGetH - Gets the last value that was used as the differencing for a `MATMFFD` matrix
683e884886fSBarry Smith    parameter.
684e884886fSBarry Smith 
685e884886fSBarry Smith    Not Collective
686e884886fSBarry Smith 
687e884886fSBarry Smith    Input Parameters:
688*11a5261eSBarry Smith .  mat - the `MATMFFD` matrix
689e884886fSBarry Smith 
690fd292e60Sprj-    Output Parameter:
691e884886fSBarry Smith .  h - the differencing step size
692e884886fSBarry Smith 
693e884886fSBarry Smith    Level: advanced
694e884886fSBarry Smith 
695*11a5261eSBarry Smith .seealso: `MATMFFD`, `MatCreateSNESMF()`, `MatMFFDSetHHistory()`, `MatCreateMFFD()`, `MATMFFD`, `MatMFFDResetHHistory()`
696e884886fSBarry Smith @*/
6979371c9d4SSatish Balay PetscErrorCode MatMFFDGetH(Mat mat, PetscScalar *h) {
698e884886fSBarry Smith   PetscFunctionBegin;
69988b4c220SStefano Zampini   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
700dadcf809SJacob Faibussowitsch   PetscValidScalarPointer(h, 2);
701cac4c232SBarry Smith   PetscUseMethod(mat, "MatMFFDGetH_C", (Mat, PetscScalar *), (mat, h));
702e884886fSBarry Smith   PetscFunctionReturn(0);
703e884886fSBarry Smith }
704e884886fSBarry Smith 
705e884886fSBarry Smith /*@C
706*11a5261eSBarry Smith    MatMFFDSetFunction - Sets the function used in applying the matrix free `MATMFFD` matrix.
707e884886fSBarry Smith 
708*11a5261eSBarry Smith    Logically Collective on mat
709e884886fSBarry Smith 
710e884886fSBarry Smith    Input Parameters:
711*11a5261eSBarry Smith +  mat - the matrix free matrix `MATMFFD` created via `MatCreateSNESMF()` or `MatCreateMFFD()`
712e884886fSBarry Smith .  func - the function to use
713e884886fSBarry Smith -  funcctx - optional function context passed to function
714e884886fSBarry Smith 
71514f633e2SBarry Smith    Calling Sequence of func:
71614f633e2SBarry Smith $     func (void *funcctx, Vec x, Vec f)
71714f633e2SBarry Smith 
71814f633e2SBarry Smith +  funcctx - user provided context
71914f633e2SBarry Smith .  x - input vector
72014f633e2SBarry Smith -  f - computed output function
72114f633e2SBarry Smith 
722e884886fSBarry Smith    Level: advanced
723e884886fSBarry Smith 
724e884886fSBarry Smith    Notes:
725*11a5261eSBarry Smith     If you use this you MUST call `MatAssemblyBegin()` and `MatAssemblyEnd()` on the matrix free
726e884886fSBarry Smith     matrix inside your compute Jacobian routine
727e884886fSBarry Smith 
728*11a5261eSBarry Smith     If this is not set then it will use the function set with `SNESSetFunction()` if `MatCreateSNESMF()` was used.
729e884886fSBarry Smith 
730*11a5261eSBarry Smith .seealso: `MATMFFD`, `MatCreateSNESMF()`, `MatMFFDGetH()`, `MatCreateMFFD()`, `MATMFFD`,
731db781477SPatrick Sanan           `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`, `SNESetFunction()`
732e884886fSBarry Smith @*/
7339371c9d4SSatish Balay PetscErrorCode MatMFFDSetFunction(Mat mat, PetscErrorCode (*func)(void *, Vec, Vec), void *funcctx) {
734e884886fSBarry Smith   PetscFunctionBegin;
73588b4c220SStefano Zampini   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
736cac4c232SBarry Smith   PetscTryMethod(mat, "MatMFFDSetFunction_C", (Mat, PetscErrorCode(*)(void *, Vec, Vec), void *), (mat, func, funcctx));
737e884886fSBarry Smith   PetscFunctionReturn(0);
738e884886fSBarry Smith }
739e884886fSBarry Smith 
740e884886fSBarry Smith /*@C
741*11a5261eSBarry Smith    MatMFFDSetFunctioni - Sets the function for a single component for a `MATMFFD` matrix
742e884886fSBarry Smith 
743*11a5261eSBarry Smith    Logically Collective on mat
744e884886fSBarry Smith 
745e884886fSBarry Smith    Input Parameters:
746*11a5261eSBarry Smith +  mat - the matrix free matrix `MATMFFD`
747e884886fSBarry Smith -  funci - the function to use
748e884886fSBarry Smith 
749e884886fSBarry Smith    Level: advanced
750e884886fSBarry Smith 
751e884886fSBarry Smith    Notes:
752*11a5261eSBarry Smith     If you use this you MUST call `MatAssemblyBegin()` and `MatAssemblyEnd()` on the matrix free
75394eb4328SStefano Zampini     matrix inside your compute Jacobian routine.
754*11a5261eSBarry Smith 
75594eb4328SStefano Zampini     This function is necessary to compute the diagonal of the matrix.
756486afcceSStefano Zampini     funci must not contain any MPI call as it is called inside a loop on the local portion of the vector.
757e884886fSBarry Smith 
758*11a5261eSBarry Smith .seealso: `MATMFFD`, `MatCreateSNESMF()`, `MatMFFDGetH()`, `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`, `SNESetFunction()`, `MatGetDiagonal()`
759e884886fSBarry Smith @*/
7609371c9d4SSatish Balay PetscErrorCode MatMFFDSetFunctioni(Mat mat, PetscErrorCode (*funci)(void *, PetscInt, Vec, PetscScalar *)) {
761e884886fSBarry Smith   PetscFunctionBegin;
7620700a824SBarry Smith   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
763cac4c232SBarry Smith   PetscTryMethod(mat, "MatMFFDSetFunctioni_C", (Mat, PetscErrorCode(*)(void *, PetscInt, Vec, PetscScalar *)), (mat, funci));
764e884886fSBarry Smith   PetscFunctionReturn(0);
765e884886fSBarry Smith }
766e884886fSBarry Smith 
767e884886fSBarry Smith /*@C
768*11a5261eSBarry Smith    MatMFFDSetFunctioniBase - Sets the base vector for a single component function evaluation for a `MATMFFD` matrix
769e884886fSBarry Smith 
770*11a5261eSBarry Smith    Logically Collective on mat
771e884886fSBarry Smith 
772e884886fSBarry Smith    Input Parameters:
773*11a5261eSBarry Smith +  mat - the `MATMFFD` matrix free matrix
774e884886fSBarry Smith -  func - the function to use
775e884886fSBarry Smith 
776e884886fSBarry Smith    Level: advanced
777e884886fSBarry Smith 
778e884886fSBarry Smith    Notes:
779*11a5261eSBarry Smith     If you use this you MUST call `MatAssemblyBegin()` and `MatAssemblyEnd()` on the matrix free
78094eb4328SStefano Zampini     matrix inside your compute Jacobian routine.
781e884886fSBarry Smith 
782*11a5261eSBarry Smith     This function is necessary to compute the diagonal of the matrix, used for example with `PCJACOBI`
783*11a5261eSBarry Smith 
784*11a5261eSBarry Smith .seealso: `MATMFFD`, `MatCreateSNESMF()`, `MatMFFDGetH()`, `MatCreateMFFD()`, `MATMFFD`
785db781477SPatrick Sanan           `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`, `SNESetFunction()`, `MatGetDiagonal()`
786e884886fSBarry Smith @*/
7879371c9d4SSatish Balay PetscErrorCode MatMFFDSetFunctioniBase(Mat mat, PetscErrorCode (*func)(void *, Vec)) {
788e884886fSBarry Smith   PetscFunctionBegin;
7890700a824SBarry Smith   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
790cac4c232SBarry Smith   PetscTryMethod(mat, "MatMFFDSetFunctioniBase_C", (Mat, PetscErrorCode(*)(void *, Vec)), (mat, func));
791e884886fSBarry Smith   PetscFunctionReturn(0);
792e884886fSBarry Smith }
793e884886fSBarry Smith 
794e884886fSBarry Smith /*@
795*11a5261eSBarry Smith    MatMFFDSetPeriod - Sets how often h is recomputed for a `MATMFFD` matrix, by default it is every time
796e884886fSBarry Smith 
797*11a5261eSBarry Smith    Logically Collective on mat
798e884886fSBarry Smith 
799e884886fSBarry Smith    Input Parameters:
800*11a5261eSBarry Smith +  mat - the `MATMFFD` matrix free matrix
801e884886fSBarry Smith -  period - 1 for every time, 2 for every second etc
802e884886fSBarry Smith 
803e884886fSBarry Smith    Options Database Keys:
80467b8a455SSatish Balay .  -mat_mffd_period <period> - Sets how often h is recomputed
805e884886fSBarry Smith 
806e884886fSBarry Smith    Level: advanced
807e884886fSBarry Smith 
808*11a5261eSBarry Smith .seealso: `MATMFFD`, `MatCreateSNESMF()`, `MatMFFDGetH()`,
809db781477SPatrick Sanan           `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`
810e884886fSBarry Smith @*/
8119371c9d4SSatish Balay PetscErrorCode MatMFFDSetPeriod(Mat mat, PetscInt period) {
812e884886fSBarry Smith   PetscFunctionBegin;
81388b4c220SStefano Zampini   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
81488b4c220SStefano Zampini   PetscValidLogicalCollectiveInt(mat, period, 2);
815cac4c232SBarry Smith   PetscTryMethod(mat, "MatMFFDSetPeriod_C", (Mat, PetscInt), (mat, period));
816e884886fSBarry Smith   PetscFunctionReturn(0);
817e884886fSBarry Smith }
818e884886fSBarry Smith 
819e884886fSBarry Smith /*@
820*11a5261eSBarry Smith    MatMFFDSetFunctionError - Sets the error_rel for the approximation of matrix-vector products using finite differences with the `MATMFFD` matrix
821e884886fSBarry Smith 
822*11a5261eSBarry Smith    Logically Collective on mat
823e884886fSBarry Smith 
824e884886fSBarry Smith    Input Parameters:
825*11a5261eSBarry Smith +  mat - the `MATMFFD` matrix free matrix
826*11a5261eSBarry Smith -  error_rel - relative error (should be set to the square root of the relative error in the function evaluations)
827e884886fSBarry Smith 
828e884886fSBarry Smith    Options Database Keys:
829a2b725a8SWilliam Gropp .  -mat_mffd_err <error_rel> - Sets error_rel
830e884886fSBarry Smith 
831e884886fSBarry Smith    Level: advanced
832e884886fSBarry Smith 
833*11a5261eSBarry Smith    Note:
834e884886fSBarry Smith    The default matrix-free matrix-vector product routine computes
835e884886fSBarry Smith .vb
836e884886fSBarry Smith      F'(u)*a = [F(u+h*a) - F(u)]/h where
837e884886fSBarry Smith      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
838e884886fSBarry Smith        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   else
839e884886fSBarry Smith .ve
840e884886fSBarry Smith 
841*11a5261eSBarry Smith .seealso: `MATMFFD`, `MatCreateSNESMF()`, `MatMFFDGetH()`, `MatCreateMFFD()`, `MATMFFD`
842db781477SPatrick Sanan           `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`
843e884886fSBarry Smith @*/
8449371c9d4SSatish Balay PetscErrorCode MatMFFDSetFunctionError(Mat mat, PetscReal error) {
845e884886fSBarry Smith   PetscFunctionBegin;
84688b4c220SStefano Zampini   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
84788b4c220SStefano Zampini   PetscValidLogicalCollectiveReal(mat, error, 2);
848cac4c232SBarry Smith   PetscTryMethod(mat, "MatMFFDSetFunctionError_C", (Mat, PetscReal), (mat, error));
849e884886fSBarry Smith   PetscFunctionReturn(0);
850e884886fSBarry Smith }
851e884886fSBarry Smith 
852e884886fSBarry Smith /*@
853e884886fSBarry Smith    MatMFFDSetHHistory - Sets an array to collect a history of the
854*11a5261eSBarry Smith    differencing values (h) computed for the matrix-free product `MATMFFD` matrix
855e884886fSBarry Smith 
856*11a5261eSBarry Smith    Logically Collective on mat
857e884886fSBarry Smith 
858e884886fSBarry Smith    Input Parameters:
859*11a5261eSBarry Smith +  J - the `MATMFFD` matrix-free matrix
860a5b23f4aSJose E. Roman .  history - space to hold the history
861e884886fSBarry Smith -  nhistory - number of entries in history, if more entries are generated than
862e884886fSBarry Smith               nhistory, then the later ones are discarded
863e884886fSBarry Smith 
864e884886fSBarry Smith    Level: advanced
865e884886fSBarry Smith 
866*11a5261eSBarry Smith    Note:
867*11a5261eSBarry Smith    Use `MatMFFDResetHHistory()` to reset the history counter and collect
868e884886fSBarry Smith    a new batch of differencing parameters, h.
869e884886fSBarry Smith 
870*11a5261eSBarry Smith .seealso: `MATMFFD`, `MatMFFDGetH()`, `MatCreateSNESMF()`,
871db781477SPatrick Sanan           `MatMFFDResetHHistory()`, `MatMFFDSetFunctionError()`
872e884886fSBarry Smith @*/
8739371c9d4SSatish Balay PetscErrorCode MatMFFDSetHHistory(Mat J, PetscScalar history[], PetscInt nhistory) {
874e884886fSBarry Smith   PetscFunctionBegin;
87588b4c220SStefano Zampini   PetscValidHeaderSpecific(J, MAT_CLASSID, 1);
876dadcf809SJacob Faibussowitsch   if (history) PetscValidScalarPointer(history, 2);
87788b4c220SStefano Zampini   PetscValidLogicalCollectiveInt(J, nhistory, 3);
878cac4c232SBarry Smith   PetscUseMethod(J, "MatMFFDSetHHistory_C", (Mat, PetscScalar[], PetscInt), (J, history, nhistory));
879e884886fSBarry Smith   PetscFunctionReturn(0);
880e884886fSBarry Smith }
881e884886fSBarry Smith 
882e884886fSBarry Smith /*@
883e884886fSBarry Smith    MatMFFDResetHHistory - Resets the counter to zero to begin
884*11a5261eSBarry Smith    collecting a new set of differencing histories for the `MATMFFD` matrix
885e884886fSBarry Smith 
886*11a5261eSBarry Smith    Logically Collective on mat
887e884886fSBarry Smith 
888e884886fSBarry Smith    Input Parameters:
889e884886fSBarry Smith .  J - the matrix-free matrix context
890e884886fSBarry Smith 
891e884886fSBarry Smith    Level: advanced
892e884886fSBarry Smith 
893*11a5261eSBarry Smith    Note:
894*11a5261eSBarry Smith    Use `MatMFFDSetHHistory()` to create the original history counter.
895e884886fSBarry Smith 
896*11a5261eSBarry Smith .seealso: `MATMFFD`, `MatMFFDGetH()`, `MatCreateSNESMF()`,
897db781477SPatrick Sanan           `MatMFFDSetHHistory()`, `MatMFFDSetFunctionError()`
898e884886fSBarry Smith @*/
8999371c9d4SSatish Balay PetscErrorCode MatMFFDResetHHistory(Mat J) {
900e884886fSBarry Smith   PetscFunctionBegin;
90188b4c220SStefano Zampini   PetscValidHeaderSpecific(J, MAT_CLASSID, 1);
902cac4c232SBarry Smith   PetscTryMethod(J, "MatMFFDResetHHistory_C", (Mat), (J));
903e884886fSBarry Smith   PetscFunctionReturn(0);
904e884886fSBarry Smith }
905e884886fSBarry Smith 
906487a658cSBarry Smith /*@
907e884886fSBarry Smith     MatMFFDSetBase - Sets the vector U at which matrix vector products of the
908*11a5261eSBarry Smith         Jacobian are computed for the `MATMFFD` matrix
909e884886fSBarry Smith 
910*11a5261eSBarry Smith     Logically Collective on mat
911e884886fSBarry Smith 
912e884886fSBarry Smith     Input Parameters:
913*11a5261eSBarry Smith +   J - the `MATMFFD` matrix
9143ec795f1SBarry Smith .   U - the vector
915bcddec3dSBarry Smith -   F - (optional) vector that contains F(u) if it has been already computed
916e884886fSBarry Smith 
91795452b02SPatrick Sanan     Notes:
91895452b02SPatrick Sanan     This is rarely used directly
919e884886fSBarry Smith 
9208af5ae88SBarry Smith     If F is provided then it is not recomputed. Otherwise the function is evaluated at the base
921*11a5261eSBarry Smith     point during the first `MatMult()` after each call to `MatMFFDSetBase()`.
922dff2f722SBarry Smith 
923e884886fSBarry Smith     Level: advanced
924e884886fSBarry Smith 
925*11a5261eSBarry Smith .seealso: `MATMFFD`, `MatMult()`, `MatMFFDSetBase()`
926e884886fSBarry Smith @*/
9279371c9d4SSatish Balay PetscErrorCode MatMFFDSetBase(Mat J, Vec U, Vec F) {
928e884886fSBarry Smith   PetscFunctionBegin;
9290700a824SBarry Smith   PetscValidHeaderSpecific(J, MAT_CLASSID, 1);
9300700a824SBarry Smith   PetscValidHeaderSpecific(U, VEC_CLASSID, 2);
9310700a824SBarry Smith   if (F) PetscValidHeaderSpecific(F, VEC_CLASSID, 3);
932cac4c232SBarry Smith   PetscTryMethod(J, "MatMFFDSetBase_C", (Mat, Vec, Vec), (J, U, F));
933e884886fSBarry Smith   PetscFunctionReturn(0);
934e884886fSBarry Smith }
935e884886fSBarry Smith 
936e884886fSBarry Smith /*@C
937e884886fSBarry Smith     MatMFFDSetCheckh - Sets a function that checks the computed h and adjusts
938*11a5261eSBarry Smith         it to satisfy some criteria for the `MATMFFD` matrix
939e884886fSBarry Smith 
940*11a5261eSBarry Smith     Logically Collective on mat
941e884886fSBarry Smith 
942e884886fSBarry Smith     Input Parameters:
943*11a5261eSBarry Smith +   J - the `MATMFFD` matrix
944e884886fSBarry Smith .   fun - the function that checks h
945e884886fSBarry Smith -   ctx - any context needed by the function
946e884886fSBarry Smith 
947e884886fSBarry Smith     Options Database Keys:
94867b8a455SSatish Balay .   -mat_mffd_check_positivity <bool> - Insure that U + h*a is non-negative
949e884886fSBarry Smith 
950e884886fSBarry Smith     Level: advanced
951e884886fSBarry Smith 
95295452b02SPatrick Sanan     Notes:
953*11a5261eSBarry Smith     For example, `MatMFFDCheckPositivity()` insures that all entries
954e884886fSBarry Smith        of U + h*a are non-negative
955e884886fSBarry Smith 
956755b7f64SBarry Smith      The function you provide is called after the default h has been computed and allows you to
957755b7f64SBarry Smith      modify it.
958755b7f64SBarry Smith 
959*11a5261eSBarry Smith .seealso: `MATMFFD`, `MatMFFDCheckPositivity()`
960e884886fSBarry Smith @*/
9619371c9d4SSatish Balay PetscErrorCode MatMFFDSetCheckh(Mat J, PetscErrorCode (*fun)(void *, Vec, Vec, PetscScalar *), void *ctx) {
962e884886fSBarry Smith   PetscFunctionBegin;
9630700a824SBarry Smith   PetscValidHeaderSpecific(J, MAT_CLASSID, 1);
964cac4c232SBarry Smith   PetscTryMethod(J, "MatMFFDSetCheckh_C", (Mat, PetscErrorCode(*)(void *, Vec, Vec, PetscScalar *), void *), (J, fun, ctx));
965e884886fSBarry Smith   PetscFunctionReturn(0);
966e884886fSBarry Smith }
967e884886fSBarry Smith 
968e884886fSBarry Smith /*@
969e884886fSBarry Smith     MatMFFDCheckPositivity - Checks that all entries in U + h*a are positive or
970*11a5261eSBarry Smith         zero, decreases h until this is satisfied for a `MATMFFD` matrix
971e884886fSBarry Smith 
972*11a5261eSBarry Smith     Logically Collective on U
973e884886fSBarry Smith 
974e884886fSBarry Smith     Input Parameters:
975e884886fSBarry Smith +   U - base vector that is added to
976e884886fSBarry Smith .   a - vector that is added
977e884886fSBarry Smith .   h - scaling factor on a
978e884886fSBarry Smith -   dummy - context variable (unused)
979e884886fSBarry Smith 
980e884886fSBarry Smith     Options Database Keys:
98167b8a455SSatish Balay .   -mat_mffd_check_positivity <bool> - Insure that U + h*a is nonnegative
982e884886fSBarry Smith 
983e884886fSBarry Smith     Level: advanced
984e884886fSBarry Smith 
985*11a5261eSBarry Smith     Note:
98695452b02SPatrick Sanan     This is rarely used directly, rather it is passed as an argument to
987*11a5261eSBarry Smith            `MatMFFDSetCheckh()`
988e884886fSBarry Smith 
989*11a5261eSBarry Smith .seealso: `MATMFFD`, `MatMFFDSetCheckh()`
990e884886fSBarry Smith @*/
9919371c9d4SSatish Balay PetscErrorCode MatMFFDCheckPositivity(void *dummy, Vec U, Vec a, PetscScalar *h) {
992e884886fSBarry Smith   PetscReal    val, minval;
993e884886fSBarry Smith   PetscScalar *u_vec, *a_vec;
994e884886fSBarry Smith   PetscInt     i, n;
995e884886fSBarry Smith   MPI_Comm     comm;
996e884886fSBarry Smith 
997e884886fSBarry Smith   PetscFunctionBegin;
99888b4c220SStefano Zampini   PetscValidHeaderSpecific(U, VEC_CLASSID, 2);
99988b4c220SStefano Zampini   PetscValidHeaderSpecific(a, VEC_CLASSID, 3);
1000dadcf809SJacob Faibussowitsch   PetscValidScalarPointer(h, 4);
10019566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)U, &comm));
10029566063dSJacob Faibussowitsch   PetscCall(VecGetArray(U, &u_vec));
10039566063dSJacob Faibussowitsch   PetscCall(VecGetArray(a, &a_vec));
10049566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(U, &n));
1005c068d9bbSLisandro Dalcin   minval = PetscAbsScalar(*h) * PetscRealConstant(1.01);
1006e884886fSBarry Smith   for (i = 0; i < n; i++) {
1007e884886fSBarry Smith     if (PetscRealPart(u_vec[i] + *h * a_vec[i]) <= 0.0) {
1008e884886fSBarry Smith       val = PetscAbsScalar(u_vec[i] / a_vec[i]);
1009e884886fSBarry Smith       if (val < minval) minval = val;
1010e884886fSBarry Smith     }
1011e884886fSBarry Smith   }
10129566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(U, &u_vec));
10139566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(a, &a_vec));
10141c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&minval, &val, 1, MPIU_REAL, MPIU_MIN, comm));
1015e884886fSBarry Smith   if (val <= PetscAbsScalar(*h)) {
10169566063dSJacob Faibussowitsch     PetscCall(PetscInfo(U, "Scaling back h from %g to %g\n", (double)PetscRealPart(*h), (double)(.99 * val)));
1017e884886fSBarry Smith     if (PetscRealPart(*h) > 0.0) *h = 0.99 * val;
1018e884886fSBarry Smith     else *h = -0.99 * val;
1019e884886fSBarry Smith   }
1020e884886fSBarry Smith   PetscFunctionReturn(0);
1021e884886fSBarry Smith }
1022