xref: /petsc/src/mat/impls/mffd/mffd.c (revision c2e3fba1fe1cda7e6350bbca19c4ed35ce95940a)
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
132390153bSJed Brown   MatMFFDFinalizePackage - This function destroys everything in the MatMFFD package. It is
14b022a5c1SBarry Smith   called from PetscFinalize().
15b022a5c1SBarry Smith 
16b022a5c1SBarry Smith   Level: developer
17b022a5c1SBarry Smith 
18db781477SPatrick Sanan .seealso: `PetscFinalize()`, `MatCreateMFFD()`, `MatCreateSNESMF()`
19b022a5c1SBarry Smith @*/
207087cfbeSBarry Smith PetscErrorCode  MatMFFDFinalizePackage(void)
21b022a5c1SBarry Smith {
22b022a5c1SBarry Smith   PetscFunctionBegin;
239566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListDestroy(&MatMFFDList));
24b022a5c1SBarry Smith   MatMFFDPackageInitialized = PETSC_FALSE;
25b022a5c1SBarry Smith   MatMFFDRegisterAllCalled  = PETSC_FALSE;
26b022a5c1SBarry Smith   PetscFunctionReturn(0);
27b022a5c1SBarry Smith }
28b022a5c1SBarry Smith 
293ec795f1SBarry Smith /*@C
303ec795f1SBarry Smith   MatMFFDInitializePackage - This function initializes everything in the MatMFFD package. It is called
318a690491SBarry Smith   from MatInitializePackage().
323ec795f1SBarry Smith 
333ec795f1SBarry Smith   Level: developer
343ec795f1SBarry Smith 
35db781477SPatrick Sanan .seealso: `PetscInitialize()`
363ec795f1SBarry Smith @*/
37607a6623SBarry Smith PetscErrorCode  MatMFFDInitializePackage(void)
383ec795f1SBarry Smith {
393ec795f1SBarry Smith   char           logList[256];
408e81d068SLisandro Dalcin   PetscBool      opt,pkg;
413ec795f1SBarry Smith 
423ec795f1SBarry Smith   PetscFunctionBegin;
43b022a5c1SBarry Smith   if (MatMFFDPackageInitialized) PetscFunctionReturn(0);
44b022a5c1SBarry Smith   MatMFFDPackageInitialized = PETSC_TRUE;
453ec795f1SBarry Smith   /* Register Classes */
469566063dSJacob Faibussowitsch   PetscCall(PetscClassIdRegister("MatMFFD",&MATMFFD_CLASSID));
473ec795f1SBarry Smith   /* Register Constructors */
489566063dSJacob Faibussowitsch   PetscCall(MatMFFDRegisterAll());
493ec795f1SBarry Smith   /* Register Events */
509566063dSJacob Faibussowitsch   PetscCall(PetscLogEventRegister("MatMult MF",MATMFFD_CLASSID,&MATMFFD_Mult));
51e94e781bSJacob Faibussowitsch  /* Process Info */
52e94e781bSJacob Faibussowitsch   {
53e94e781bSJacob Faibussowitsch     PetscClassId  classids[1];
54e94e781bSJacob Faibussowitsch 
55e94e781bSJacob Faibussowitsch     classids[0] = MATMFFD_CLASSID;
569566063dSJacob Faibussowitsch     PetscCall(PetscInfoProcessClass("matmffd", 1, classids));
573ec795f1SBarry Smith   }
583ec795f1SBarry Smith   /* Process summary exclusions */
599566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL,NULL,"-log_exclude",logList,sizeof(logList),&opt));
603ec795f1SBarry Smith   if (opt) {
619566063dSJacob Faibussowitsch     PetscCall(PetscStrInList("matmffd",logList,',',&pkg));
629566063dSJacob Faibussowitsch     if (pkg) PetscCall(PetscLogEventExcludeClass(MATMFFD_CLASSID));
633ec795f1SBarry Smith   }
648e81d068SLisandro Dalcin   /* Register package finalizer */
659566063dSJacob Faibussowitsch   PetscCall(PetscRegisterFinalize(MatMFFDFinalizePackage));
663ec795f1SBarry Smith   PetscFunctionReturn(0);
673ec795f1SBarry Smith }
683ec795f1SBarry Smith 
69789d8953SBarry Smith static PetscErrorCode  MatMFFDSetType_MFFD(Mat mat,MatMFFDType ftype)
70789d8953SBarry Smith {
71789d8953SBarry Smith   MatMFFD        ctx;
72789d8953SBarry Smith   PetscBool      match;
735f80ce2aSJacob Faibussowitsch   PetscErrorCode (*r)(MatMFFD);
74789d8953SBarry Smith 
75789d8953SBarry Smith   PetscFunctionBegin;
76789d8953SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
77789d8953SBarry Smith   PetscValidCharPointer(ftype,2);
789566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat,&ctx));
79789d8953SBarry Smith 
80789d8953SBarry Smith   /* already set, so just return */
819566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)ctx,ftype,&match));
82789d8953SBarry Smith   if (match) PetscFunctionReturn(0);
83789d8953SBarry Smith 
84789d8953SBarry Smith   /* destroy the old one if it exists */
859566063dSJacob Faibussowitsch   if (ctx->ops->destroy) PetscCall((*ctx->ops->destroy)(ctx));
86789d8953SBarry Smith 
879566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(MatMFFDList,ftype,&r));
885f80ce2aSJacob Faibussowitsch   PetscCheck(r,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown MatMFFD type %s given",ftype);
899566063dSJacob Faibussowitsch   PetscCall((*r)(ctx));
909566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)ctx,ftype));
91789d8953SBarry Smith   PetscFunctionReturn(0);
92789d8953SBarry Smith }
93789d8953SBarry Smith 
94e884886fSBarry Smith /*@C
95e884886fSBarry Smith     MatMFFDSetType - Sets the method that is used to compute the
96e884886fSBarry Smith     differencing parameter for finite differene matrix-free formulations.
97e884886fSBarry Smith 
98e884886fSBarry Smith     Input Parameters:
99e884886fSBarry Smith +   mat - the "matrix-free" matrix created via MatCreateSNESMF(), or MatCreateMFFD()
100e884886fSBarry Smith           or MatSetType(mat,MATMFFD);
101e884886fSBarry Smith -   ftype - the type requested, either MATMFFD_WP or MATMFFD_DS
102e884886fSBarry Smith 
103e884886fSBarry Smith     Level: advanced
104e884886fSBarry Smith 
105e884886fSBarry Smith     Notes:
106e884886fSBarry Smith     For example, such routines can compute h for use in
107e884886fSBarry Smith     Jacobian-vector products of the form
108e884886fSBarry Smith 
109e884886fSBarry Smith                         F(x+ha) - F(x)
110e884886fSBarry Smith           F'(u)a  ~=  ----------------
111e884886fSBarry Smith                               h
112e884886fSBarry Smith 
113db781477SPatrick Sanan .seealso: `MatCreateSNESMF()`, `MatMFFDRegister()`, `MatMFFDSetFunction()`, `MatCreateMFFD()`
114e884886fSBarry Smith @*/
11519fd82e9SBarry Smith PetscErrorCode  MatMFFDSetType(Mat mat,MatMFFDType ftype)
116e884886fSBarry Smith {
117e884886fSBarry Smith   PetscFunctionBegin;
1180700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
119e884886fSBarry Smith   PetscValidCharPointer(ftype,2);
120cac4c232SBarry Smith   PetscTryMethod(mat,"MatMFFDSetType_C",(Mat,MatMFFDType),(mat,ftype));
121e884886fSBarry Smith   PetscFunctionReturn(0);
122e884886fSBarry Smith }
123e884886fSBarry Smith 
1245d21e1f8SStefano Zampini static PetscErrorCode MatGetDiagonal_MFFD(Mat,Vec);
1255d21e1f8SStefano Zampini 
126e884886fSBarry Smith typedef PetscErrorCode (*FCN1)(void*,Vec); /* force argument to next function to not be extern C*/
12740244768SBarry Smith static PetscErrorCode  MatMFFDSetFunctioniBase_MFFD(Mat mat,FCN1 func)
128e884886fSBarry Smith {
129789d8953SBarry Smith   MatMFFD        ctx;
130e884886fSBarry Smith 
131e884886fSBarry Smith   PetscFunctionBegin;
1329566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat,&ctx));
133e884886fSBarry Smith   ctx->funcisetbase = func;
134e884886fSBarry Smith   PetscFunctionReturn(0);
135e884886fSBarry Smith }
136e884886fSBarry Smith 
137e884886fSBarry Smith typedef PetscErrorCode (*FCN2)(void*,PetscInt,Vec,PetscScalar*); /* force argument to next function to not be extern C*/
13840244768SBarry Smith static PetscErrorCode  MatMFFDSetFunctioni_MFFD(Mat mat,FCN2 funci)
139e884886fSBarry Smith {
140789d8953SBarry Smith   MatMFFD        ctx;
141e884886fSBarry Smith 
142e884886fSBarry Smith   PetscFunctionBegin;
1439566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat,&ctx));
144e884886fSBarry Smith   ctx->funci = funci;
1459566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(mat,MATOP_GET_DIAGONAL,(void (*)(void))MatGetDiagonal_MFFD));
146789d8953SBarry Smith   PetscFunctionReturn(0);
1475d21e1f8SStefano Zampini }
148789d8953SBarry Smith 
149789d8953SBarry Smith static PetscErrorCode MatMFFDGetH_MFFD(Mat mat,PetscScalar *h)
150789d8953SBarry Smith {
151789d8953SBarry Smith   MatMFFD        ctx;
152789d8953SBarry Smith 
153789d8953SBarry Smith   PetscFunctionBegin;
1549566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat,&ctx));
155789d8953SBarry Smith   *h = ctx->currenth;
156e884886fSBarry Smith   PetscFunctionReturn(0);
157e884886fSBarry Smith }
158e884886fSBarry Smith 
15940244768SBarry Smith static PetscErrorCode  MatMFFDResetHHistory_MFFD(Mat J)
160bc0f08ceSBarry Smith {
161789d8953SBarry Smith   MatMFFD        ctx;
162bc0f08ceSBarry Smith 
163bc0f08ceSBarry Smith   PetscFunctionBegin;
1649566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J,&ctx));
165bc0f08ceSBarry Smith   ctx->ncurrenth = 0;
166bc0f08ceSBarry Smith   PetscFunctionReturn(0);
167bc0f08ceSBarry Smith }
168e884886fSBarry Smith 
1691c84c290SBarry Smith /*@C
1701c84c290SBarry Smith    MatMFFDRegister - Adds a method to the MatMFFD registry.
1711c84c290SBarry Smith 
1721c84c290SBarry Smith    Not Collective
1731c84c290SBarry Smith 
1741c84c290SBarry Smith    Input Parameters:
1751c84c290SBarry Smith +  name_solver - name of a new user-defined compute-h module
1761c84c290SBarry Smith -  routine_create - routine to create method context
1771c84c290SBarry Smith 
1781c84c290SBarry Smith    Level: developer
1791c84c290SBarry Smith 
1801c84c290SBarry Smith    Notes:
1811c84c290SBarry Smith    MatMFFDRegister() may be called multiple times to add several user-defined solvers.
1821c84c290SBarry Smith 
1831c84c290SBarry Smith    Sample usage:
1841c84c290SBarry Smith .vb
185bdf89e91SBarry Smith    MatMFFDRegister("my_h",MyHCreate);
1861c84c290SBarry Smith .ve
1871c84c290SBarry Smith 
1881c84c290SBarry Smith    Then, your solver can be chosen with the procedural interface via
1891c84c290SBarry Smith $     MatMFFDSetType(mfctx,"my_h")
1901c84c290SBarry Smith    or at runtime via the option
191be7a6d03SBarry Smith $     -mat_mffd_type my_h
1921c84c290SBarry Smith 
193db781477SPatrick Sanan .seealso: `MatMFFDRegisterAll()`, `MatMFFDRegisterDestroy()`
1941c84c290SBarry Smith  @*/
195bdf89e91SBarry Smith PetscErrorCode  MatMFFDRegister(const char sname[],PetscErrorCode (*function)(MatMFFD))
196e884886fSBarry Smith {
197e884886fSBarry Smith   PetscFunctionBegin;
1989566063dSJacob Faibussowitsch   PetscCall(MatInitializePackage());
1999566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&MatMFFDList,sname,function));
200e884886fSBarry Smith   PetscFunctionReturn(0);
201e884886fSBarry Smith }
202e884886fSBarry Smith 
203e884886fSBarry Smith /* ----------------------------------------------------------------------------------------*/
20440244768SBarry Smith static PetscErrorCode MatDestroy_MFFD(Mat mat)
205e884886fSBarry Smith {
206789d8953SBarry Smith   MatMFFD        ctx;
207e884886fSBarry Smith 
208e884886fSBarry Smith   PetscFunctionBegin;
2099566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat,&ctx));
2109566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ctx->w));
2119566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ctx->current_u));
212cfe22f5eSBarry Smith   if (ctx->current_f_allocated) {
2139566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&ctx->current_f));
214cfe22f5eSBarry Smith   }
2159566063dSJacob Faibussowitsch   if (ctx->ops->destroy) PetscCall((*ctx->ops->destroy)(ctx));
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));
229e884886fSBarry Smith   PetscFunctionReturn(0);
230e884886fSBarry Smith }
231e884886fSBarry Smith 
232e884886fSBarry Smith /*
233e884886fSBarry Smith    MatMFFDView_MFFD - Views matrix-free parameters.
234e884886fSBarry Smith 
235e884886fSBarry Smith */
23640244768SBarry Smith static PetscErrorCode MatView_MFFD(Mat J,PetscViewer viewer)
237e884886fSBarry Smith {
238789d8953SBarry Smith   MatMFFD        ctx;
239a283635eSDmitry Karpeev   PetscBool      iascii, viewbase, viewfunction;
240a283635eSDmitry Karpeev   const char     *prefix;
241e884886fSBarry Smith 
242e884886fSBarry Smith   PetscFunctionBegin;
2439566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J,&ctx));
2449566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
245e884886fSBarry Smith   if (iascii) {
2469566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"Matrix-free approximation:\n"));
2479566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
2489566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"err=%g (relative error in function evaluation)\n",(double)ctx->error_rel));
2497adad957SLisandro Dalcin     if (!((PetscObject)ctx)->type_name) {
2509566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"The compute h routine has not yet been set\n"));
251e884886fSBarry Smith     } else {
2529566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"Using %s compute h routine\n",((PetscObject)ctx)->type_name));
253e884886fSBarry Smith     }
254c92e3469SBarry Smith #if defined(PETSC_USE_COMPLEX)
255c92e3469SBarry Smith     if (ctx->usecomplex) {
2569566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"Using Lyness complex number trick to compute the matrix-vector product\n"));
257c92e3469SBarry Smith     }
258c92e3469SBarry Smith #endif
259e884886fSBarry Smith     if (ctx->ops->view) {
2609566063dSJacob Faibussowitsch       PetscCall((*ctx->ops->view)(ctx,viewer));
261e884886fSBarry Smith     }
2629566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetOptionsPrefix((PetscObject)J, &prefix));
263a283635eSDmitry Karpeev 
2649566063dSJacob Faibussowitsch     PetscCall(PetscOptionsHasName(((PetscObject)J)->options,prefix, "-mat_mffd_view_base", &viewbase));
265a283635eSDmitry Karpeev     if (viewbase) {
2669566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Base:\n"));
2679566063dSJacob Faibussowitsch       PetscCall(VecView(ctx->current_u, viewer));
268a283635eSDmitry Karpeev     }
2699566063dSJacob Faibussowitsch     PetscCall(PetscOptionsHasName(((PetscObject)J)->options,prefix, "-mat_mffd_view_function", &viewfunction));
270a283635eSDmitry Karpeev     if (viewfunction) {
2719566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Function:\n"));
2729566063dSJacob Faibussowitsch       PetscCall(VecView(ctx->current_f, viewer));
273a283635eSDmitry Karpeev     }
2749566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
27511aeaf0aSBarry Smith   }
276e884886fSBarry Smith   PetscFunctionReturn(0);
277e884886fSBarry Smith }
278e884886fSBarry Smith 
279e884886fSBarry Smith /*
280e884886fSBarry Smith    MatAssemblyEnd_MFFD - Resets the ctx->ncurrenth to zero. This
281e884886fSBarry Smith    allows the user to indicate the beginning of a new linear solve by calling
282e884886fSBarry Smith    MatAssemblyXXX() on the matrix free matrix. This then allows the
2831d0fab5eSBarry Smith    MatCreateMFFD_WP() to properly compute ||U|| only the first time
284e884886fSBarry Smith    in the linear solver rather than every time.
2855a576424SJed Brown 
286cc2e6a90SBarry Smith    This function is referenced directly from MatAssemblyEnd_SNESMF(), which may be in a different shared library hence
287cc2e6a90SBarry Smith    it must be labeled as PETSC_EXTERN
288e884886fSBarry Smith */
2895a576424SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyEnd_MFFD(Mat J,MatAssemblyType mt)
290e884886fSBarry Smith {
291789d8953SBarry Smith   MatMFFD        j;
292e884886fSBarry Smith 
293e884886fSBarry Smith   PetscFunctionBegin;
2949566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J,&j));
2959566063dSJacob Faibussowitsch   PetscCall(MatMFFDResetHHistory(J));
296e884886fSBarry Smith   PetscFunctionReturn(0);
297e884886fSBarry Smith }
298e884886fSBarry Smith 
299e884886fSBarry Smith /*
300e884886fSBarry Smith   MatMult_MFFD - Default matrix-free form for Jacobian-vector product, y = F'(u)*a:
301e884886fSBarry Smith 
302e884886fSBarry Smith         y ~= (F(u + ha) - F(u))/h,
303e884886fSBarry Smith   where F = nonlinear function, as set by SNESSetFunction()
304e884886fSBarry Smith         u = current iterate
305e884886fSBarry Smith         h = difference interval
306e884886fSBarry Smith */
30740244768SBarry Smith static PetscErrorCode MatMult_MFFD(Mat mat,Vec a,Vec y)
308e884886fSBarry Smith {
309789d8953SBarry Smith   MatMFFD        ctx;
310e884886fSBarry Smith   PetscScalar    h;
311e884886fSBarry Smith   Vec            w,U,F;
312ace3abfcSBarry Smith   PetscBool      zeroa;
313e884886fSBarry Smith 
314e884886fSBarry Smith   PetscFunctionBegin;
3159566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat,&ctx));
31628b400f6SJacob 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");
317e884886fSBarry Smith   /* We log matrix-free matrix-vector products separately, so that we can
318e884886fSBarry Smith      separate the performance monitoring from the cases that use conventional
319e884886fSBarry Smith      storage.  We may eventually modify event logging to associate events
320e884886fSBarry Smith      with particular objects, hence alleviating the more general problem. */
3219566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(MATMFFD_Mult,a,y,0,0));
322e884886fSBarry Smith 
323e884886fSBarry Smith   w = ctx->w;
324e884886fSBarry Smith   U = ctx->current_u;
3253ec795f1SBarry Smith   F = ctx->current_f;
326e884886fSBarry Smith   /*
327e884886fSBarry Smith       Compute differencing parameter
328e884886fSBarry Smith   */
3294a2f8832SBarry Smith   if (!((PetscObject)ctx)->type_name) {
3309566063dSJacob Faibussowitsch     PetscCall(MatMFFDSetType(mat,MATMFFD_WP));
3319566063dSJacob Faibussowitsch     PetscCall(MatSetFromOptions(mat));
332e884886fSBarry Smith   }
3339566063dSJacob Faibussowitsch   PetscCall((*ctx->ops->compute)(ctx,U,a,&h,&zeroa));
334e884886fSBarry Smith   if (zeroa) {
3359566063dSJacob Faibussowitsch     PetscCall(VecSet(y,0.0));
336e884886fSBarry Smith     PetscFunctionReturn(0);
337e884886fSBarry Smith   }
338e884886fSBarry Smith 
33908401ef6SPierre Jolivet   PetscCheck(!mat->erroriffailure || !PetscIsInfOrNanScalar(h),PETSC_COMM_SELF,PETSC_ERR_PLIB,"Computed Nan differencing parameter h");
340e884886fSBarry Smith   if (ctx->checkh) {
3419566063dSJacob Faibussowitsch     PetscCall((*ctx->checkh)(ctx->checkhctx,U,a,&h));
342e884886fSBarry Smith   }
343e884886fSBarry Smith 
344e884886fSBarry Smith   /* keep a record of the current differencing parameter h */
345e884886fSBarry Smith   ctx->currenth = h;
346e884886fSBarry Smith #if defined(PETSC_USE_COMPLEX)
3479566063dSJacob Faibussowitsch   PetscCall(PetscInfo(mat,"Current differencing parameter: %g + %g i\n",(double)PetscRealPart(h),(double)PetscImaginaryPart(h)));
348e884886fSBarry Smith #else
3499566063dSJacob Faibussowitsch   PetscCall(PetscInfo(mat,"Current differencing parameter: %15.12e\n",(double)PetscRealPart(h)));
350e884886fSBarry Smith #endif
351e884886fSBarry Smith   if (ctx->historyh && ctx->ncurrenth < ctx->maxcurrenth) {
352e884886fSBarry Smith     ctx->historyh[ctx->ncurrenth] = h;
353e884886fSBarry Smith   }
354e884886fSBarry Smith   ctx->ncurrenth++;
355e884886fSBarry Smith 
356c92e3469SBarry Smith #if defined(PETSC_USE_COMPLEX)
357c92e3469SBarry Smith   if (ctx->usecomplex) h = PETSC_i*h;
358c92e3469SBarry Smith #endif
359c92e3469SBarry Smith 
360e884886fSBarry Smith   /* w = u + ha */
3619566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(w,h,a,U));
362e884886fSBarry Smith 
3631b797266SDmitry Karpeev   /* compute func(U) as base for differencing; only needed first time in and not when provided by user */
3641b797266SDmitry Karpeev   if (ctx->ncurrenth == 1 && ctx->current_f_allocated) {
3659566063dSJacob Faibussowitsch     PetscCall((*ctx->func)(ctx->funcctx,U,F));
36652121784SBarry Smith   }
3679566063dSJacob Faibussowitsch   PetscCall((*ctx->func)(ctx->funcctx,w,y));
368e884886fSBarry Smith 
369c92e3469SBarry Smith #if defined(PETSC_USE_COMPLEX)
370c92e3469SBarry Smith   if (ctx->usecomplex) {
3719566063dSJacob Faibussowitsch     PetscCall(VecImaginaryPart(y));
372c92e3469SBarry Smith     h    = PetscImaginaryPart(h);
373c92e3469SBarry Smith   } else {
3749566063dSJacob Faibussowitsch     PetscCall(VecAXPY(y,-1.0,F));
375c92e3469SBarry Smith   }
376c92e3469SBarry Smith #else
3779566063dSJacob Faibussowitsch   PetscCall(VecAXPY(y,-1.0,F));
378c92e3469SBarry Smith #endif
3799566063dSJacob Faibussowitsch   PetscCall(VecScale(y,1.0/h));
3809566063dSJacob Faibussowitsch   if (mat->nullsp) PetscCall(MatNullSpaceRemove(mat->nullsp,y));
381e884886fSBarry Smith 
3829566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(MATMFFD_Mult,a,y,0,0));
383e884886fSBarry Smith   PetscFunctionReturn(0);
384e884886fSBarry Smith }
385e884886fSBarry Smith 
386e884886fSBarry Smith /*
387e884886fSBarry Smith   MatGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix
388e884886fSBarry Smith 
389e884886fSBarry Smith         y ~= (F(u + ha) - F(u))/h,
390e884886fSBarry Smith   where F = nonlinear function, as set by SNESSetFunction()
391e884886fSBarry Smith         u = current iterate
392e884886fSBarry Smith         h = difference interval
393e884886fSBarry Smith */
3945d21e1f8SStefano Zampini PetscErrorCode MatGetDiagonal_MFFD(Mat mat,Vec a)
395e884886fSBarry Smith {
396789d8953SBarry Smith   MatMFFD        ctx;
397e884886fSBarry Smith   PetscScalar    h,*aa,*ww,v;
398e884886fSBarry Smith   PetscReal      epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON;
399e884886fSBarry Smith   Vec            w,U;
400e884886fSBarry Smith   PetscInt       i,rstart,rend;
401e884886fSBarry Smith 
402e884886fSBarry Smith   PetscFunctionBegin;
4039566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat,&ctx));
40428b400f6SJacob Faibussowitsch   PetscCheck(ctx->func,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Requires calling MatMFFDSetFunction() first");
40528b400f6SJacob Faibussowitsch   PetscCheck(ctx->funci,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Requires calling MatMFFDSetFunctioni() first");
406e884886fSBarry Smith   w    = ctx->w;
407e884886fSBarry Smith   U    = ctx->current_u;
4089566063dSJacob Faibussowitsch   PetscCall((*ctx->func)(ctx->funcctx,U,a));
409789d8953SBarry Smith   if (ctx->funcisetbase) {
4109566063dSJacob Faibussowitsch     PetscCall((*ctx->funcisetbase)(ctx->funcctx,U));
411789d8953SBarry Smith   }
4129566063dSJacob Faibussowitsch   PetscCall(VecCopy(U,w));
413e884886fSBarry Smith 
4149566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(a,&rstart,&rend));
4159566063dSJacob Faibussowitsch   PetscCall(VecGetArray(a,&aa));
416e884886fSBarry Smith   for (i=rstart; i<rend; i++) {
4179566063dSJacob Faibussowitsch     PetscCall(VecGetArray(w,&ww));
418e884886fSBarry Smith     h    = ww[i-rstart];
419e884886fSBarry Smith     if (h == 0.0) h = 1.0;
420e884886fSBarry Smith     if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0)     h = umin;
421e884886fSBarry Smith     else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin;
422e884886fSBarry Smith     h *= epsilon;
423e884886fSBarry Smith 
424e884886fSBarry Smith     ww[i-rstart] += h;
4259566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(w,&ww));
4269566063dSJacob Faibussowitsch     PetscCall((*ctx->funci)(ctx->funcctx,i,w,&v));
427e884886fSBarry Smith     aa[i-rstart]  = (v - aa[i-rstart])/h;
428e884886fSBarry Smith 
4299566063dSJacob Faibussowitsch     PetscCall(VecGetArray(w,&ww));
430e884886fSBarry Smith     ww[i-rstart] -= h;
4319566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(w,&ww));
432e884886fSBarry Smith   }
4339566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(a,&aa));
434e884886fSBarry Smith   PetscFunctionReturn(0);
435e884886fSBarry Smith }
436e884886fSBarry Smith 
437d8fa6a54SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDSetBase_MFFD(Mat J,Vec U,Vec F)
438e884886fSBarry Smith {
439789d8953SBarry Smith   MatMFFD        ctx;
440e884886fSBarry Smith 
441e884886fSBarry Smith   PetscFunctionBegin;
4429566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J,&ctx));
4439566063dSJacob Faibussowitsch   PetscCall(MatMFFDResetHHistory(J));
444c51ad4d4SStefano Zampini   if (!ctx->current_u) {
4459566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(U,&ctx->current_u));
4469566063dSJacob Faibussowitsch     PetscCall(VecLockReadPush(ctx->current_u));
447c51ad4d4SStefano Zampini   }
4489566063dSJacob Faibussowitsch   PetscCall(VecLockReadPop(ctx->current_u));
4499566063dSJacob Faibussowitsch   PetscCall(VecCopy(U,ctx->current_u));
4509566063dSJacob Faibussowitsch   PetscCall(VecLockReadPush(ctx->current_u));
45152121784SBarry Smith   if (F) {
4529566063dSJacob Faibussowitsch     if (ctx->current_f_allocated) PetscCall(VecDestroy(&ctx->current_f));
4533ec795f1SBarry Smith     ctx->current_f           = F;
454cfe22f5eSBarry Smith     ctx->current_f_allocated = PETSC_FALSE;
455cfe22f5eSBarry Smith   } else if (!ctx->current_f_allocated) {
4569566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(J,NULL,&ctx->current_f));
457cfe22f5eSBarry Smith     ctx->current_f_allocated = PETSC_TRUE;
45852121784SBarry Smith   }
459e884886fSBarry Smith   if (!ctx->w) {
4609566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ctx->current_u,&ctx->w));
461e884886fSBarry Smith   }
462e884886fSBarry Smith   J->assembled = PETSC_TRUE;
463e884886fSBarry Smith   PetscFunctionReturn(0);
464e884886fSBarry Smith }
4655a576424SJed Brown 
466e884886fSBarry Smith typedef PetscErrorCode (*FCN3)(void*,Vec,Vec,PetscScalar*); /* force argument to next function to not be extern C*/
467b2573a8aSBarry Smith 
46840244768SBarry Smith static PetscErrorCode  MatMFFDSetCheckh_MFFD(Mat J,FCN3 fun,void *ectx)
469e884886fSBarry Smith {
470789d8953SBarry Smith   MatMFFD        ctx;
471e884886fSBarry Smith 
472e884886fSBarry Smith   PetscFunctionBegin;
4739566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J,&ctx));
474e884886fSBarry Smith   ctx->checkh    = fun;
475e884886fSBarry Smith   ctx->checkhctx = ectx;
476e884886fSBarry Smith   PetscFunctionReturn(0);
477e884886fSBarry Smith }
478e884886fSBarry Smith 
4796aa9148fSLisandro Dalcin /*@C
4806aa9148fSLisandro Dalcin    MatMFFDSetOptionsPrefix - Sets the prefix used for searching for all
4816aa9148fSLisandro Dalcin    MatMFFD options in the database.
4826aa9148fSLisandro Dalcin 
4836aa9148fSLisandro Dalcin    Collective on Mat
4846aa9148fSLisandro Dalcin 
485d8d19677SJose E. Roman    Input Parameters:
4866aa9148fSLisandro Dalcin +  A - the Mat context
4876aa9148fSLisandro Dalcin -  prefix - the prefix to prepend to all option names
4886aa9148fSLisandro Dalcin 
4896aa9148fSLisandro Dalcin    Notes:
4906aa9148fSLisandro Dalcin    A hyphen (-) must NOT be given at the beginning of the prefix name.
4916aa9148fSLisandro Dalcin    The first character of all runtime options is AUTOMATICALLY the hyphen.
4926aa9148fSLisandro Dalcin 
4936aa9148fSLisandro Dalcin    Level: advanced
4946aa9148fSLisandro Dalcin 
495db781477SPatrick Sanan .seealso: `MatSetFromOptions()`, `MatCreateSNESMF()`, `MatCreateMFFD()`
4966aa9148fSLisandro Dalcin @*/
4977087cfbeSBarry Smith PetscErrorCode  MatMFFDSetOptionsPrefix(Mat mat,const char prefix[])
4986aa9148fSLisandro Dalcin {
499789d8953SBarry Smith   MatMFFD mfctx;
5005fd66863SKarl Rupp 
5016aa9148fSLisandro Dalcin   PetscFunctionBegin;
5020700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
5039566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat,&mfctx));
5040700a824SBarry Smith   PetscValidHeaderSpecific(mfctx,MATMFFD_CLASSID,1);
5059566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)mfctx,prefix));
5066aa9148fSLisandro Dalcin   PetscFunctionReturn(0);
5076aa9148fSLisandro Dalcin }
5086aa9148fSLisandro Dalcin 
50940244768SBarry Smith static PetscErrorCode  MatSetFromOptions_MFFD(PetscOptionItems *PetscOptionsObject,Mat mat)
510e884886fSBarry Smith {
511789d8953SBarry Smith   MatMFFD        mfctx;
512ace3abfcSBarry Smith   PetscBool      flg;
513e884886fSBarry Smith   char           ftype[256];
514e884886fSBarry Smith 
515e884886fSBarry Smith   PetscFunctionBegin;
516064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(mat,MAT_CLASSID,2);
5179566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat,&mfctx));
518064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(mfctx,MATMFFD_CLASSID,2);
519d0609cedSBarry Smith   PetscObjectOptionsBegin((PetscObject)mfctx);
5209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFList("-mat_mffd_type","Matrix free type","MatMFFDSetType",MatMFFDList,((PetscObject)mfctx)->type_name,ftype,256,&flg));
521e884886fSBarry Smith   if (flg) {
5229566063dSJacob Faibussowitsch     PetscCall(MatMFFDSetType(mat,ftype));
523e884886fSBarry Smith   }
524e884886fSBarry Smith 
5259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_mffd_err","set sqrt relative error in function","MatMFFDSetFunctionError",mfctx->error_rel,&mfctx->error_rel,NULL));
5269566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-mat_mffd_period","how often h is recomputed","MatMFFDSetPeriod",mfctx->recomputeperiod,&mfctx->recomputeperiod,NULL));
527e884886fSBarry Smith 
52890d69ab7SBarry Smith   flg  = PETSC_FALSE;
5299566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-mat_mffd_check_positivity","Insure that U + h*a is nonnegative","MatMFFDSetCheckh",flg,&flg,NULL));
530e884886fSBarry Smith   if (flg) {
5319566063dSJacob Faibussowitsch     PetscCall(MatMFFDSetCheckh(mat,MatMFFDCheckPositivity,NULL));
532e884886fSBarry Smith   }
533c92e3469SBarry Smith #if defined(PETSC_USE_COMPLEX)
5349566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-mat_mffd_complex","Use Lyness complex number trick to compute the matrix-vector product","None",mfctx->usecomplex,&mfctx->usecomplex,NULL));
535c92e3469SBarry Smith #endif
536e884886fSBarry Smith   if (mfctx->ops->setfromoptions) {
5379566063dSJacob Faibussowitsch     PetscCall((*mfctx->ops->setfromoptions)(PetscOptionsObject,mfctx));
538e884886fSBarry Smith   }
539d0609cedSBarry Smith   PetscOptionsEnd();
540e884886fSBarry Smith   PetscFunctionReturn(0);
541e884886fSBarry Smith }
542e884886fSBarry Smith 
54340244768SBarry Smith static PetscErrorCode  MatMFFDSetPeriod_MFFD(Mat mat,PetscInt period)
544bc0f08ceSBarry Smith {
545789d8953SBarry Smith   MatMFFD        ctx;
546bc0f08ceSBarry Smith 
547bc0f08ceSBarry Smith   PetscFunctionBegin;
5489566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat,&ctx));
549bc0f08ceSBarry Smith   ctx->recomputeperiod = period;
550bc0f08ceSBarry Smith   PetscFunctionReturn(0);
551bc0f08ceSBarry Smith }
552bc0f08ceSBarry Smith 
55340244768SBarry Smith static PetscErrorCode  MatMFFDSetFunction_MFFD(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx)
554bc0f08ceSBarry Smith {
555789d8953SBarry Smith   MatMFFD        ctx;
556bc0f08ceSBarry Smith 
557bc0f08ceSBarry Smith   PetscFunctionBegin;
5589566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat,&ctx));
559bc0f08ceSBarry Smith   ctx->func    = func;
560bc0f08ceSBarry Smith   ctx->funcctx = funcctx;
561bc0f08ceSBarry Smith   PetscFunctionReturn(0);
562bc0f08ceSBarry Smith }
563bc0f08ceSBarry Smith 
56440244768SBarry Smith static PetscErrorCode  MatMFFDSetFunctionError_MFFD(Mat mat,PetscReal error)
565bc0f08ceSBarry Smith {
566789d8953SBarry Smith   MatMFFD        ctx;
567bc0f08ceSBarry Smith 
568bc0f08ceSBarry Smith   PetscFunctionBegin;
5699566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat,&ctx));
570bc0f08ceSBarry Smith   if (error != PETSC_DEFAULT) ctx->error_rel = error;
571bc0f08ceSBarry Smith   PetscFunctionReturn(0);
572bc0f08ceSBarry Smith }
573bc0f08ceSBarry Smith 
574789d8953SBarry Smith PetscErrorCode  MatMFFDSetHHistory_MFFD(Mat J,PetscScalar history[],PetscInt nhistory)
5753b49f96aSBarry Smith {
576789d8953SBarry Smith   MatMFFD        ctx;
577789d8953SBarry Smith 
5783b49f96aSBarry Smith   PetscFunctionBegin;
5799566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J,&ctx));
580789d8953SBarry Smith   ctx->historyh    = history;
581789d8953SBarry Smith   ctx->maxcurrenth = nhistory;
582789d8953SBarry Smith   ctx->currenth    = 0.;
5833b49f96aSBarry Smith   PetscFunctionReturn(0);
5843b49f96aSBarry Smith }
5853b49f96aSBarry Smith 
586e884886fSBarry Smith /*MC
587e884886fSBarry Smith   MATMFFD - MATMFFD = "mffd" - A matrix free matrix type.
588e884886fSBarry Smith 
589e884886fSBarry Smith   Level: advanced
590e884886fSBarry Smith 
591789d8953SBarry Smith   Developers Note: This is implemented on top of MATSHELL to get support for scaling and shifting without requiring duplicate code
592789d8953SBarry Smith 
593db781477SPatrick Sanan .seealso: `MatCreateMFFD()`, `MatCreateSNESMF()`, `MatMFFDSetFunction()`, `MatMFFDSetType()`,
594db781477SPatrick Sanan           `MatMFFDSetFunctionError()`, `MatMFFDDSSetUmin()`, `MatMFFDSetFunction()`
595db781477SPatrick Sanan           `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`, `MatCreateSNESMF()`,
596db781477SPatrick Sanan           `MatMFFDGetH()`,
597e884886fSBarry Smith M*/
5988cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_MFFD(Mat A)
599e884886fSBarry Smith {
600e884886fSBarry Smith   MatMFFD        mfctx;
601e884886fSBarry Smith 
602e884886fSBarry Smith   PetscFunctionBegin;
6039566063dSJacob Faibussowitsch   PetscCall(MatMFFDInitializePackage());
604e884886fSBarry Smith 
6059566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(mfctx,MATMFFD_CLASSID,"MatMFFD","Matrix-free Finite Differencing","Mat",PetscObjectComm((PetscObject)A),NULL,NULL));
6062205254eSKarl Rupp 
607e884886fSBarry Smith   mfctx->error_rel                = PETSC_SQRT_MACHINE_EPSILON;
608e884886fSBarry Smith   mfctx->recomputeperiod          = 1;
609e884886fSBarry Smith   mfctx->count                    = 0;
610e884886fSBarry Smith   mfctx->currenth                 = 0.0;
6110298fd71SBarry Smith   mfctx->historyh                 = NULL;
612e884886fSBarry Smith   mfctx->ncurrenth                = 0;
613e884886fSBarry Smith   mfctx->maxcurrenth              = 0;
614f4259b30SLisandro Dalcin   ((PetscObject)mfctx)->type_name = NULL;
615e884886fSBarry Smith 
616e884886fSBarry Smith   /*
617e884886fSBarry Smith      Create the empty data structure to contain compute-h routines.
618e884886fSBarry Smith      These will be filled in below from the command line options or
619e884886fSBarry Smith      a later call with MatMFFDSetType() or if that is not called
620e884886fSBarry Smith      then it will default in the first use of MatMult_MFFD()
621e884886fSBarry Smith   */
622f4259b30SLisandro Dalcin   mfctx->ops->compute        = NULL;
623f4259b30SLisandro Dalcin   mfctx->ops->destroy        = NULL;
624f4259b30SLisandro Dalcin   mfctx->ops->view           = NULL;
625f4259b30SLisandro Dalcin   mfctx->ops->setfromoptions = NULL;
626f4259b30SLisandro Dalcin   mfctx->hctx                = NULL;
627e884886fSBarry Smith 
628f4259b30SLisandro Dalcin   mfctx->func    = NULL;
629f4259b30SLisandro Dalcin   mfctx->funcctx = NULL;
6300298fd71SBarry Smith   mfctx->w       = NULL;
631789d8953SBarry Smith   mfctx->mat     = A;
632e884886fSBarry Smith 
6339566063dSJacob Faibussowitsch   PetscCall(MatSetType(A,MATSHELL));
6349566063dSJacob Faibussowitsch   PetscCall(MatShellSetContext(A,mfctx));
6359566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(A,MATOP_MULT,(void (*)(void))MatMult_MFFD));
6369566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(A,MATOP_DESTROY,(void (*)(void))MatDestroy_MFFD));
6379566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(A,MATOP_VIEW,(void (*)(void))MatView_MFFD));
6389566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(A,MATOP_ASSEMBLY_END,(void (*)(void))MatAssemblyEnd_MFFD));
6399566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(A,MATOP_SET_FROM_OPTIONS,(void (*)(void))MatSetFromOptions_MFFD));
640e884886fSBarry Smith 
6419566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetBase_C",MatMFFDSetBase_MFFD));
6429566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctioniBase_C",MatMFFDSetFunctioniBase_MFFD));
6439566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctioni_C",MatMFFDSetFunctioni_MFFD));
6449566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunction_C",MatMFFDSetFunction_MFFD));
6459566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetCheckh_C",MatMFFDSetCheckh_MFFD));
6469566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetPeriod_C",MatMFFDSetPeriod_MFFD));
6479566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctionError_C",MatMFFDSetFunctionError_MFFD));
6489566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatMFFDResetHHistory_C",MatMFFDResetHHistory_MFFD));
6499566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetHHistory_C",MatMFFDSetHHistory_MFFD));
6509566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetType_C",MatMFFDSetType_MFFD));
6519566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatMFFDGetH_C",MatMFFDGetH_MFFD));
6529566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)A,MATMFFD));
653e884886fSBarry Smith   PetscFunctionReturn(0);
654e884886fSBarry Smith }
655e884886fSBarry Smith 
656e884886fSBarry Smith /*@
657e884886fSBarry Smith    MatCreateMFFD - Creates a matrix-free matrix. See also MatCreateSNESMF()
658e884886fSBarry Smith 
659e884886fSBarry Smith    Collective on Vec
660e884886fSBarry Smith 
661e884886fSBarry Smith    Input Parameters:
662fef1beadSBarry Smith +  comm - MPI communicator
663fef1beadSBarry Smith .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
664fef1beadSBarry Smith            This value should be the same as the local size used in creating the
665fef1beadSBarry Smith            y vector for the matrix-vector product y = Ax.
666fef1beadSBarry Smith .  n - This value should be the same as the local size used in creating the
667fef1beadSBarry Smith        x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have
668fef1beadSBarry Smith        calculated if N is given) For square matrices n is almost always m.
669fef1beadSBarry Smith .  M - number of global rows (or PETSC_DETERMINE to have calculated if m is given)
670fef1beadSBarry Smith -  N - number of global columns (or PETSC_DETERMINE to have calculated if n is given)
671fef1beadSBarry Smith 
672e884886fSBarry Smith    Output Parameter:
673e884886fSBarry Smith .  J - the matrix-free matrix
674e884886fSBarry Smith 
6759c6ac3b3SBarry Smith    Options Database Keys: call MatSetFromOptions() to trigger these
6769c6ac3b3SBarry Smith +  -mat_mffd_type - wp or ds (see MATMFFD_WP or MATMFFD_DS)
677c92e3469SBarry Smith .  -mat_mffd_err - square root of estimated relative error in function evaluation
678c92e3469SBarry Smith .  -mat_mffd_period - how often h is recomputed, defaults to 1, everytime
679c92e3469SBarry Smith .  -mat_mffd_check_positivity - possibly decrease h until U + h*a has only positive values
680c92e3469SBarry Smith -  -mat_mffd_complex - use the Lyness trick with complex numbers to compute the matrix-vector product instead of differencing
681c92e3469SBarry Smith                        (requires real valued functions but that PETSc be configured for complex numbers)
6829c6ac3b3SBarry Smith 
683e884886fSBarry Smith    Level: advanced
684e884886fSBarry Smith 
685e884886fSBarry Smith    Notes:
686e884886fSBarry Smith    The matrix-free matrix context merely contains the function pointers
687e884886fSBarry Smith    and work space for performing finite difference approximations of
688e884886fSBarry Smith    Jacobian-vector products, F'(u)*a,
689e884886fSBarry Smith 
690e884886fSBarry Smith    The default code uses the following approach to compute h
691e884886fSBarry Smith 
692e884886fSBarry Smith .vb
693e884886fSBarry Smith      F'(u)*a = [F(u+h*a) - F(u)]/h where
694e884886fSBarry Smith      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
695e884886fSBarry Smith        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
696e884886fSBarry Smith  where
697e884886fSBarry Smith      error_rel = square root of relative error in function evaluation
698e884886fSBarry Smith      umin = minimum iterate parameter
699e884886fSBarry Smith .ve
700e884886fSBarry Smith 
701ca93e954SBarry Smith    You can call SNESSetJacobian() with MatMFFDComputeJacobian() if you are using matrix and not a different
702ca93e954SBarry Smith    preconditioner matrix
703ca93e954SBarry Smith 
704e884886fSBarry Smith    The user can set the error_rel via MatMFFDSetFunctionError() and
705a7f22e61SSatish Balay    umin via MatMFFDDSSetUmin(); see Users-Manual: ch_snes for details.
706e884886fSBarry Smith 
707e884886fSBarry Smith    The user should call MatDestroy() when finished with the matrix-free
708e884886fSBarry Smith    matrix context.
709e884886fSBarry Smith 
710e884886fSBarry Smith    Options Database Keys:
711e884886fSBarry Smith +  -mat_mffd_err <error_rel> - Sets error_rel
7121e2ae407SBarry Smith .  -mat_mffd_umin <umin> - Sets umin (for default PETSc routine that computes h only)
71367b8a455SSatish Balay -  -mat_mffd_check_positivity - check positivity
714e884886fSBarry Smith 
715db781477SPatrick Sanan .seealso: `MatDestroy()`, `MatMFFDSetFunctionError()`, `MatMFFDDSSetUmin()`, `MatMFFDSetFunction()`
716db781477SPatrick Sanan           `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`, `MatCreateSNESMF()`,
717db781477SPatrick Sanan           `MatMFFDGetH()`, `MatMFFDRegister()`, `MatMFFDComputeJacobian()`
718e884886fSBarry Smith 
719e884886fSBarry Smith @*/
7207087cfbeSBarry Smith PetscErrorCode  MatCreateMFFD(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *J)
721e884886fSBarry Smith {
722e884886fSBarry Smith   PetscFunctionBegin;
7239566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm,J));
7249566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*J,m,n,M,N));
7259566063dSJacob Faibussowitsch   PetscCall(MatSetType(*J,MATMFFD));
7269566063dSJacob Faibussowitsch   PetscCall(MatSetUp(*J));
727e884886fSBarry Smith   PetscFunctionReturn(0);
728e884886fSBarry Smith }
729e884886fSBarry Smith 
730e884886fSBarry Smith /*@
731e884886fSBarry Smith    MatMFFDGetH - Gets the last value that was used as the differencing
732e884886fSBarry Smith    parameter.
733e884886fSBarry Smith 
734e884886fSBarry Smith    Not Collective
735e884886fSBarry Smith 
736e884886fSBarry Smith    Input Parameters:
737e884886fSBarry Smith .  mat - the matrix obtained with MatCreateSNESMF()
738e884886fSBarry Smith 
739fd292e60Sprj-    Output Parameter:
740e884886fSBarry Smith .  h - the differencing step size
741e884886fSBarry Smith 
742e884886fSBarry Smith    Level: advanced
743e884886fSBarry Smith 
744*c2e3fba1SPatrick Sanan .seealso: `MatCreateSNESMF()`, `MatMFFDSetHHistory()`, `MatCreateMFFD()`, `MATMFFD`, `MatMFFDResetHHistory()`
745e884886fSBarry Smith @*/
7467087cfbeSBarry Smith PetscErrorCode  MatMFFDGetH(Mat mat,PetscScalar *h)
747e884886fSBarry Smith {
748e884886fSBarry Smith   PetscFunctionBegin;
74988b4c220SStefano Zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
750dadcf809SJacob Faibussowitsch   PetscValidScalarPointer(h,2);
751cac4c232SBarry Smith   PetscUseMethod(mat,"MatMFFDGetH_C",(Mat,PetscScalar*),(mat,h));
752e884886fSBarry Smith   PetscFunctionReturn(0);
753e884886fSBarry Smith }
754e884886fSBarry Smith 
755e884886fSBarry Smith /*@C
756e884886fSBarry Smith    MatMFFDSetFunction - Sets the function used in applying the matrix free.
757e884886fSBarry Smith 
7583f9fe445SBarry Smith    Logically Collective on Mat
759e884886fSBarry Smith 
760e884886fSBarry Smith    Input Parameters:
76114f633e2SBarry Smith +  mat - the matrix free matrix created via MatCreateSNESMF() or MatCreateMFFD()
762e884886fSBarry Smith .  func - the function to use
763e884886fSBarry Smith -  funcctx - optional function context passed to function
764e884886fSBarry Smith 
76514f633e2SBarry Smith    Calling Sequence of func:
76614f633e2SBarry Smith $     func (void *funcctx, Vec x, Vec f)
76714f633e2SBarry Smith 
76814f633e2SBarry Smith +  funcctx - user provided context
76914f633e2SBarry Smith .  x - input vector
77014f633e2SBarry Smith -  f - computed output function
77114f633e2SBarry Smith 
772e884886fSBarry Smith    Level: advanced
773e884886fSBarry Smith 
774e884886fSBarry Smith    Notes:
775e884886fSBarry Smith     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
776e884886fSBarry Smith     matrix inside your compute Jacobian routine
777e884886fSBarry Smith 
7783ec795f1SBarry Smith     If this is not set then it will use the function set with SNESSetFunction() if MatCreateSNESMF() was used.
779e884886fSBarry Smith 
780*c2e3fba1SPatrick Sanan .seealso: `MatCreateSNESMF()`, `MatMFFDGetH()`, `MatCreateMFFD()`, `MATMFFD`,
781db781477SPatrick Sanan           `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`, `SNESetFunction()`
782e884886fSBarry Smith @*/
7837087cfbeSBarry Smith PetscErrorCode  MatMFFDSetFunction(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx)
784e884886fSBarry Smith {
785e884886fSBarry Smith   PetscFunctionBegin;
78688b4c220SStefano Zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
787cac4c232SBarry Smith   PetscTryMethod(mat,"MatMFFDSetFunction_C",(Mat,PetscErrorCode (*)(void*,Vec,Vec),void*),(mat,func,funcctx));
788e884886fSBarry Smith   PetscFunctionReturn(0);
789e884886fSBarry Smith }
790e884886fSBarry Smith 
791e884886fSBarry Smith /*@C
792e884886fSBarry Smith    MatMFFDSetFunctioni - Sets the function for a single component
793e884886fSBarry Smith 
7943f9fe445SBarry Smith    Logically Collective on Mat
795e884886fSBarry Smith 
796e884886fSBarry Smith    Input Parameters:
797e884886fSBarry Smith +  mat - the matrix free matrix created via MatCreateSNESMF()
798e884886fSBarry Smith -  funci - the function to use
799e884886fSBarry Smith 
800e884886fSBarry Smith    Level: advanced
801e884886fSBarry Smith 
802e884886fSBarry Smith    Notes:
803e884886fSBarry Smith     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
80494eb4328SStefano Zampini     matrix inside your compute Jacobian routine.
80594eb4328SStefano Zampini     This function is necessary to compute the diagonal of the matrix.
806486afcceSStefano Zampini     funci must not contain any MPI call as it is called inside a loop on the local portion of the vector.
807e884886fSBarry Smith 
808*c2e3fba1SPatrick Sanan .seealso: `MatCreateSNESMF()`, `MatMFFDGetH()`, `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`, `SNESetFunction()`, `MatGetDiagonal()`
8091d0fab5eSBarry Smith 
810e884886fSBarry Smith @*/
8117087cfbeSBarry Smith PetscErrorCode  MatMFFDSetFunctioni(Mat mat,PetscErrorCode (*funci)(void*,PetscInt,Vec,PetscScalar*))
812e884886fSBarry Smith {
813e884886fSBarry Smith   PetscFunctionBegin;
8140700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
815cac4c232SBarry Smith   PetscTryMethod(mat,"MatMFFDSetFunctioni_C",(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)),(mat,funci));
816e884886fSBarry Smith   PetscFunctionReturn(0);
817e884886fSBarry Smith }
818e884886fSBarry Smith 
819e884886fSBarry Smith /*@C
820e884886fSBarry Smith    MatMFFDSetFunctioniBase - Sets the base vector for a single component function evaluation
821e884886fSBarry Smith 
8223f9fe445SBarry Smith    Logically Collective on Mat
823e884886fSBarry Smith 
824e884886fSBarry Smith    Input Parameters:
825e884886fSBarry Smith +  mat - the matrix free matrix created via MatCreateSNESMF()
826e884886fSBarry Smith -  func - the function to use
827e884886fSBarry Smith 
828e884886fSBarry Smith    Level: advanced
829e884886fSBarry Smith 
830e884886fSBarry Smith    Notes:
831e884886fSBarry Smith     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
83294eb4328SStefano Zampini     matrix inside your compute Jacobian routine.
83394eb4328SStefano Zampini     This function is necessary to compute the diagonal of the matrix.
834e884886fSBarry Smith 
835*c2e3fba1SPatrick Sanan .seealso: `MatCreateSNESMF()`, `MatMFFDGetH()`, `MatCreateMFFD()`, `MATMFFD`
836db781477SPatrick Sanan           `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`, `SNESetFunction()`, `MatGetDiagonal()`
837e884886fSBarry Smith @*/
8387087cfbeSBarry Smith PetscErrorCode  MatMFFDSetFunctioniBase(Mat mat,PetscErrorCode (*func)(void*,Vec))
839e884886fSBarry Smith {
840e884886fSBarry Smith   PetscFunctionBegin;
8410700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
842cac4c232SBarry Smith   PetscTryMethod(mat,"MatMFFDSetFunctioniBase_C",(Mat,PetscErrorCode (*)(void*,Vec)),(mat,func));
843e884886fSBarry Smith   PetscFunctionReturn(0);
844e884886fSBarry Smith }
845e884886fSBarry Smith 
846e884886fSBarry Smith /*@
847e884886fSBarry Smith    MatMFFDSetPeriod - Sets how often h is recomputed, by default it is everytime
848e884886fSBarry Smith 
8493f9fe445SBarry Smith    Logically Collective on Mat
850e884886fSBarry Smith 
851e884886fSBarry Smith    Input Parameters:
852e884886fSBarry Smith +  mat - the matrix free matrix created via MatCreateSNESMF()
853e884886fSBarry Smith -  period - 1 for everytime, 2 for every second etc
854e884886fSBarry Smith 
855e884886fSBarry Smith    Options Database Keys:
85667b8a455SSatish Balay .  -mat_mffd_period <period> - Sets how often h is recomputed
857e884886fSBarry Smith 
858e884886fSBarry Smith    Level: advanced
859e884886fSBarry Smith 
860*c2e3fba1SPatrick Sanan .seealso: `MatCreateSNESMF()`, `MatMFFDGetH()`,
861db781477SPatrick Sanan           `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`
862e884886fSBarry Smith @*/
8637087cfbeSBarry Smith PetscErrorCode  MatMFFDSetPeriod(Mat mat,PetscInt period)
864e884886fSBarry Smith {
865e884886fSBarry Smith   PetscFunctionBegin;
86688b4c220SStefano Zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
86788b4c220SStefano Zampini   PetscValidLogicalCollectiveInt(mat,period,2);
868cac4c232SBarry Smith   PetscTryMethod(mat,"MatMFFDSetPeriod_C",(Mat,PetscInt),(mat,period));
869e884886fSBarry Smith   PetscFunctionReturn(0);
870e884886fSBarry Smith }
871e884886fSBarry Smith 
872e884886fSBarry Smith /*@
873e884886fSBarry Smith    MatMFFDSetFunctionError - Sets the error_rel for the approximation of
874e884886fSBarry Smith    matrix-vector products using finite differences.
875e884886fSBarry Smith 
8763f9fe445SBarry Smith    Logically Collective on Mat
877e884886fSBarry Smith 
878e884886fSBarry Smith    Input Parameters:
879e884886fSBarry Smith +  mat - the matrix free matrix created via MatCreateMFFD() or MatCreateSNESMF()
880e884886fSBarry Smith -  error_rel - relative error (should be set to the square root of
881e884886fSBarry Smith                the relative error in the function evaluations)
882e884886fSBarry Smith 
883e884886fSBarry Smith    Options Database Keys:
884a2b725a8SWilliam Gropp .  -mat_mffd_err <error_rel> - Sets error_rel
885e884886fSBarry Smith 
886e884886fSBarry Smith    Level: advanced
887e884886fSBarry Smith 
888e884886fSBarry Smith    Notes:
889e884886fSBarry Smith    The default matrix-free matrix-vector product routine computes
890e884886fSBarry Smith .vb
891e884886fSBarry Smith      F'(u)*a = [F(u+h*a) - F(u)]/h where
892e884886fSBarry Smith      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
893e884886fSBarry Smith        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   else
894e884886fSBarry Smith .ve
895e884886fSBarry Smith 
896*c2e3fba1SPatrick Sanan .seealso: `MatCreateSNESMF()`, `MatMFFDGetH()`, `MatCreateMFFD()`, `MATMFFD`
897db781477SPatrick Sanan           `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`
898e884886fSBarry Smith @*/
8997087cfbeSBarry Smith PetscErrorCode  MatMFFDSetFunctionError(Mat mat,PetscReal error)
900e884886fSBarry Smith {
901e884886fSBarry Smith   PetscFunctionBegin;
90288b4c220SStefano Zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
90388b4c220SStefano Zampini   PetscValidLogicalCollectiveReal(mat,error,2);
904cac4c232SBarry Smith   PetscTryMethod(mat,"MatMFFDSetFunctionError_C",(Mat,PetscReal),(mat,error));
905e884886fSBarry Smith   PetscFunctionReturn(0);
906e884886fSBarry Smith }
907e884886fSBarry Smith 
908e884886fSBarry Smith /*@
909e884886fSBarry Smith    MatMFFDSetHHistory - Sets an array to collect a history of the
910e884886fSBarry Smith    differencing values (h) computed for the matrix-free product.
911e884886fSBarry Smith 
9123f9fe445SBarry Smith    Logically Collective on Mat
913e884886fSBarry Smith 
914e884886fSBarry Smith    Input Parameters:
915e884886fSBarry Smith +  J - the matrix-free matrix context
916a5b23f4aSJose E. Roman .  history - space to hold the history
917e884886fSBarry Smith -  nhistory - number of entries in history, if more entries are generated than
918e884886fSBarry Smith               nhistory, then the later ones are discarded
919e884886fSBarry Smith 
920e884886fSBarry Smith    Level: advanced
921e884886fSBarry Smith 
922e884886fSBarry Smith    Notes:
923e884886fSBarry Smith    Use MatMFFDResetHHistory() to reset the history counter and collect
924e884886fSBarry Smith    a new batch of differencing parameters, h.
925e884886fSBarry Smith 
926db781477SPatrick Sanan .seealso: `MatMFFDGetH()`, `MatCreateSNESMF()`,
927db781477SPatrick Sanan           `MatMFFDResetHHistory()`, `MatMFFDSetFunctionError()`
928e884886fSBarry Smith 
929e884886fSBarry Smith @*/
9307087cfbeSBarry Smith PetscErrorCode  MatMFFDSetHHistory(Mat J,PetscScalar history[],PetscInt nhistory)
931e884886fSBarry Smith {
932e884886fSBarry Smith   PetscFunctionBegin;
93388b4c220SStefano Zampini   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
934dadcf809SJacob Faibussowitsch   if (history) PetscValidScalarPointer(history,2);
93588b4c220SStefano Zampini   PetscValidLogicalCollectiveInt(J,nhistory,3);
936cac4c232SBarry Smith   PetscUseMethod(J,"MatMFFDSetHHistory_C",(Mat,PetscScalar[],PetscInt),(J,history,nhistory));
937e884886fSBarry Smith   PetscFunctionReturn(0);
938e884886fSBarry Smith }
939e884886fSBarry Smith 
940e884886fSBarry Smith /*@
941e884886fSBarry Smith    MatMFFDResetHHistory - Resets the counter to zero to begin
942e884886fSBarry Smith    collecting a new set of differencing histories.
943e884886fSBarry Smith 
9443f9fe445SBarry Smith    Logically Collective on Mat
945e884886fSBarry Smith 
946e884886fSBarry Smith    Input Parameters:
947e884886fSBarry Smith .  J - the matrix-free matrix context
948e884886fSBarry Smith 
949e884886fSBarry Smith    Level: advanced
950e884886fSBarry Smith 
951e884886fSBarry Smith    Notes:
952e884886fSBarry Smith    Use MatMFFDSetHHistory() to create the original history counter.
953e884886fSBarry Smith 
954db781477SPatrick Sanan .seealso: `MatMFFDGetH()`, `MatCreateSNESMF()`,
955db781477SPatrick Sanan           `MatMFFDSetHHistory()`, `MatMFFDSetFunctionError()`
956e884886fSBarry Smith 
957e884886fSBarry Smith @*/
9587087cfbeSBarry Smith PetscErrorCode  MatMFFDResetHHistory(Mat J)
959e884886fSBarry Smith {
960e884886fSBarry Smith   PetscFunctionBegin;
96188b4c220SStefano Zampini   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
962cac4c232SBarry Smith   PetscTryMethod(J,"MatMFFDResetHHistory_C",(Mat),(J));
963e884886fSBarry Smith   PetscFunctionReturn(0);
964e884886fSBarry Smith }
965e884886fSBarry Smith 
966487a658cSBarry Smith /*@
967e884886fSBarry Smith     MatMFFDSetBase - Sets the vector U at which matrix vector products of the
968e884886fSBarry Smith         Jacobian are computed
969e884886fSBarry Smith 
9703f9fe445SBarry Smith     Logically Collective on Mat
971e884886fSBarry Smith 
972e884886fSBarry Smith     Input Parameters:
973e884886fSBarry Smith +   J - the MatMFFD matrix
9743ec795f1SBarry Smith .   U - the vector
975bcddec3dSBarry Smith -   F - (optional) vector that contains F(u) if it has been already computed
976e884886fSBarry Smith 
97795452b02SPatrick Sanan     Notes:
97895452b02SPatrick Sanan     This is rarely used directly
979e884886fSBarry Smith 
9808af5ae88SBarry Smith     If F is provided then it is not recomputed. Otherwise the function is evaluated at the base
9818af5ae88SBarry Smith     point during the first MatMult() after each call to MatMFFDSetBase().
982dff2f722SBarry Smith 
983e884886fSBarry Smith     Level: advanced
984e884886fSBarry Smith 
985e884886fSBarry Smith @*/
9867087cfbeSBarry Smith PetscErrorCode  MatMFFDSetBase(Mat J,Vec U,Vec F)
987e884886fSBarry Smith {
988e884886fSBarry Smith   PetscFunctionBegin;
9890700a824SBarry Smith   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
9900700a824SBarry Smith   PetscValidHeaderSpecific(U,VEC_CLASSID,2);
9910700a824SBarry Smith   if (F) PetscValidHeaderSpecific(F,VEC_CLASSID,3);
992cac4c232SBarry Smith   PetscTryMethod(J,"MatMFFDSetBase_C",(Mat,Vec,Vec),(J,U,F));
993e884886fSBarry Smith   PetscFunctionReturn(0);
994e884886fSBarry Smith }
995e884886fSBarry Smith 
996e884886fSBarry Smith /*@C
997e884886fSBarry Smith     MatMFFDSetCheckh - Sets a function that checks the computed h and adjusts
998e884886fSBarry Smith         it to satisfy some criteria
999e884886fSBarry Smith 
10003f9fe445SBarry Smith     Logically Collective on Mat
1001e884886fSBarry Smith 
1002e884886fSBarry Smith     Input Parameters:
1003e884886fSBarry Smith +   J - the MatMFFD matrix
1004e884886fSBarry Smith .   fun - the function that checks h
1005e884886fSBarry Smith -   ctx - any context needed by the function
1006e884886fSBarry Smith 
1007e884886fSBarry Smith     Options Database Keys:
100867b8a455SSatish Balay .   -mat_mffd_check_positivity <bool> - Insure that U + h*a is non-negative
1009e884886fSBarry Smith 
1010e884886fSBarry Smith     Level: advanced
1011e884886fSBarry Smith 
101295452b02SPatrick Sanan     Notes:
101395452b02SPatrick Sanan     For example, MatMFFDCheckPositivity() insures that all entries
1014e884886fSBarry Smith        of U + h*a are non-negative
1015e884886fSBarry Smith 
1016755b7f64SBarry Smith      The function you provide is called after the default h has been computed and allows you to
1017755b7f64SBarry Smith      modify it.
1018755b7f64SBarry Smith 
1019db781477SPatrick Sanan .seealso: `MatMFFDCheckPositivity()`
1020e884886fSBarry Smith @*/
10217087cfbeSBarry Smith PetscErrorCode  MatMFFDSetCheckh(Mat J,PetscErrorCode (*fun)(void*,Vec,Vec,PetscScalar*),void *ctx)
1022e884886fSBarry Smith {
1023e884886fSBarry Smith   PetscFunctionBegin;
10240700a824SBarry Smith   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
1025cac4c232SBarry Smith   PetscTryMethod(J,"MatMFFDSetCheckh_C",(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*),(J,fun,ctx));
1026e884886fSBarry Smith   PetscFunctionReturn(0);
1027e884886fSBarry Smith }
1028e884886fSBarry Smith 
1029e884886fSBarry Smith /*@
1030e884886fSBarry Smith     MatMFFDCheckPositivity - Checks that all entries in U + h*a are positive or
1031e884886fSBarry Smith         zero, decreases h until this is satisfied.
1032e884886fSBarry Smith 
10333f9fe445SBarry Smith     Logically Collective on Vec
1034e884886fSBarry Smith 
1035e884886fSBarry Smith     Input Parameters:
1036e884886fSBarry Smith +   U - base vector that is added to
1037e884886fSBarry Smith .   a - vector that is added
1038e884886fSBarry Smith .   h - scaling factor on a
1039e884886fSBarry Smith -   dummy - context variable (unused)
1040e884886fSBarry Smith 
1041e884886fSBarry Smith     Options Database Keys:
104267b8a455SSatish Balay .   -mat_mffd_check_positivity <bool> - Insure that U + h*a is nonnegative
1043e884886fSBarry Smith 
1044e884886fSBarry Smith     Level: advanced
1045e884886fSBarry Smith 
104695452b02SPatrick Sanan     Notes:
104795452b02SPatrick Sanan     This is rarely used directly, rather it is passed as an argument to
1048e884886fSBarry Smith            MatMFFDSetCheckh()
1049e884886fSBarry Smith 
1050db781477SPatrick Sanan .seealso: `MatMFFDSetCheckh()`
1051e884886fSBarry Smith @*/
10527087cfbeSBarry Smith PetscErrorCode  MatMFFDCheckPositivity(void *dummy,Vec U,Vec a,PetscScalar *h)
1053e884886fSBarry Smith {
1054e884886fSBarry Smith   PetscReal      val, minval;
1055e884886fSBarry Smith   PetscScalar    *u_vec, *a_vec;
1056e884886fSBarry Smith   PetscInt       i,n;
1057e884886fSBarry Smith   MPI_Comm       comm;
1058e884886fSBarry Smith 
1059e884886fSBarry Smith   PetscFunctionBegin;
106088b4c220SStefano Zampini   PetscValidHeaderSpecific(U,VEC_CLASSID,2);
106188b4c220SStefano Zampini   PetscValidHeaderSpecific(a,VEC_CLASSID,3);
1062dadcf809SJacob Faibussowitsch   PetscValidScalarPointer(h,4);
10639566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)U,&comm));
10649566063dSJacob Faibussowitsch   PetscCall(VecGetArray(U,&u_vec));
10659566063dSJacob Faibussowitsch   PetscCall(VecGetArray(a,&a_vec));
10669566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(U,&n));
1067c068d9bbSLisandro Dalcin   minval = PetscAbsScalar(*h)*PetscRealConstant(1.01);
1068e884886fSBarry Smith   for (i=0; i<n; i++) {
1069e884886fSBarry Smith     if (PetscRealPart(u_vec[i] + *h*a_vec[i]) <= 0.0) {
1070e884886fSBarry Smith       val = PetscAbsScalar(u_vec[i]/a_vec[i]);
1071e884886fSBarry Smith       if (val < minval) minval = val;
1072e884886fSBarry Smith     }
1073e884886fSBarry Smith   }
10749566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(U,&u_vec));
10759566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(a,&a_vec));
10761c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&minval,&val,1,MPIU_REAL,MPIU_MIN,comm));
1077e884886fSBarry Smith   if (val <= PetscAbsScalar(*h)) {
10789566063dSJacob Faibussowitsch     PetscCall(PetscInfo(U,"Scaling back h from %g to %g\n",(double)PetscRealPart(*h),(double)(.99*val)));
1079e884886fSBarry Smith     if (PetscRealPart(*h) > 0.0) *h =  0.99*val;
1080e884886fSBarry Smith     else                         *h = -0.99*val;
1081e884886fSBarry Smith   }
1082e884886fSBarry Smith   PetscFunctionReturn(0);
1083e884886fSBarry Smith }
1084