xref: /petsc/src/mat/impls/mffd/mffd.c (revision 94eb4328b200eda81c2188f851bfa7fa70497842)
1e884886fSBarry Smith 
2af0996ceSBarry Smith #include <petsc/private/matimpl.h>
3c6db04a5SJed Brown #include <../src/mat/impls/mffd/mffdimpl.h>   /*I  "petscmat.h"   I*/
4e884886fSBarry Smith 
5140e18c1SBarry Smith PetscFunctionList MatMFFDList              = 0;
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 
182390153bSJed Brown .keywords: Petsc, destroy, package
19755b7f64SBarry Smith .seealso: PetscFinalize(), MatCreateMFFD(), MatCreateSNESMF()
20b022a5c1SBarry Smith @*/
217087cfbeSBarry Smith PetscErrorCode  MatMFFDFinalizePackage(void)
22b022a5c1SBarry Smith {
23a703d84cSPatrick Lacasse   PetscErrorCode ierr;
24a703d84cSPatrick Lacasse 
25b022a5c1SBarry Smith   PetscFunctionBegin;
2637e93019SBarry Smith   ierr = PetscFunctionListDestroy(&MatMFFDList);CHKERRQ(ierr);
27b022a5c1SBarry Smith   MatMFFDPackageInitialized = PETSC_FALSE;
28b022a5c1SBarry Smith   MatMFFDRegisterAllCalled  = PETSC_FALSE;
29b022a5c1SBarry Smith   PetscFunctionReturn(0);
30b022a5c1SBarry Smith }
31b022a5c1SBarry Smith 
323ec795f1SBarry Smith /*@C
333ec795f1SBarry Smith   MatMFFDInitializePackage - This function initializes everything in the MatMFFD package. It is called
343ec795f1SBarry Smith   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to MatCreate_MFFD()
353ec795f1SBarry Smith   when using static libraries.
363ec795f1SBarry Smith 
373ec795f1SBarry Smith   Level: developer
383ec795f1SBarry Smith 
393ec795f1SBarry Smith .keywords: Vec, initialize, package
403ec795f1SBarry Smith .seealso: PetscInitialize()
413ec795f1SBarry Smith @*/
42607a6623SBarry Smith PetscErrorCode  MatMFFDInitializePackage(void)
433ec795f1SBarry Smith {
443ec795f1SBarry Smith   char           logList[256];
453ec795f1SBarry Smith   char           *className;
46ace3abfcSBarry Smith   PetscBool      opt;
473ec795f1SBarry Smith   PetscErrorCode ierr;
483ec795f1SBarry Smith 
493ec795f1SBarry Smith   PetscFunctionBegin;
50b022a5c1SBarry Smith   if (MatMFFDPackageInitialized) PetscFunctionReturn(0);
51b022a5c1SBarry Smith   MatMFFDPackageInitialized = PETSC_TRUE;
523ec795f1SBarry Smith   /* Register Classes */
530700a824SBarry Smith   ierr = PetscClassIdRegister("MatMFFD",&MATMFFD_CLASSID);CHKERRQ(ierr);
543ec795f1SBarry Smith   /* Register Constructors */
55607a6623SBarry Smith   ierr = MatMFFDRegisterAll();CHKERRQ(ierr);
563ec795f1SBarry Smith   /* Register Events */
570700a824SBarry Smith   ierr = PetscLogEventRegister("MatMult MF",          MATMFFD_CLASSID,&MATMFFD_Mult);CHKERRQ(ierr);
583ec795f1SBarry Smith 
593ec795f1SBarry Smith   /* Process info exclusions */
60c5929fdfSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL, "-info_exclude", logList, 256, &opt);CHKERRQ(ierr);
613ec795f1SBarry Smith   if (opt) {
623ec795f1SBarry Smith     ierr = PetscStrstr(logList, "matmffd", &className);CHKERRQ(ierr);
633ec795f1SBarry Smith     if (className) {
640700a824SBarry Smith       ierr = PetscInfoDeactivateClass(MATMFFD_CLASSID);CHKERRQ(ierr);
653ec795f1SBarry Smith     }
663ec795f1SBarry Smith   }
673ec795f1SBarry Smith   /* Process summary exclusions */
687bf5a629SBarry Smith   ierr = PetscOptionsGetString(NULL,NULL, "-log_exclude", logList, 256, &opt);CHKERRQ(ierr);
693ec795f1SBarry Smith   if (opt) {
703ec795f1SBarry Smith     ierr = PetscStrstr(logList, "matmffd", &className);CHKERRQ(ierr);
713ec795f1SBarry Smith     if (className) {
720700a824SBarry Smith       ierr = PetscLogEventDeactivateClass(MATMFFD_CLASSID);CHKERRQ(ierr);
733ec795f1SBarry Smith     }
743ec795f1SBarry Smith   }
75b022a5c1SBarry Smith   ierr = PetscRegisterFinalize(MatMFFDFinalizePackage);CHKERRQ(ierr);
763ec795f1SBarry Smith   PetscFunctionReturn(0);
773ec795f1SBarry Smith }
783ec795f1SBarry Smith 
79e884886fSBarry Smith /*@C
80e884886fSBarry Smith     MatMFFDSetType - Sets the method that is used to compute the
81e884886fSBarry Smith     differencing parameter for finite differene matrix-free formulations.
82e884886fSBarry Smith 
83e884886fSBarry Smith     Input Parameters:
84e884886fSBarry Smith +   mat - the "matrix-free" matrix created via MatCreateSNESMF(), or MatCreateMFFD()
85e884886fSBarry Smith           or MatSetType(mat,MATMFFD);
86e884886fSBarry Smith -   ftype - the type requested, either MATMFFD_WP or MATMFFD_DS
87e884886fSBarry Smith 
88e884886fSBarry Smith     Level: advanced
89e884886fSBarry Smith 
90e884886fSBarry Smith     Notes:
91e884886fSBarry Smith     For example, such routines can compute h for use in
92e884886fSBarry Smith     Jacobian-vector products of the form
93e884886fSBarry Smith 
94e884886fSBarry Smith                         F(x+ha) - F(x)
95e884886fSBarry Smith           F'(u)a  ~=  ----------------
96e884886fSBarry Smith                               h
97e884886fSBarry Smith 
98755b7f64SBarry Smith .seealso: MatCreateSNESMF(), MatMFFDRegister(), MatMFFDSetFunction(), MatCreateMFFD()
99e884886fSBarry Smith @*/
10019fd82e9SBarry Smith PetscErrorCode  MatMFFDSetType(Mat mat,MatMFFDType ftype)
101e884886fSBarry Smith {
102e884886fSBarry Smith   PetscErrorCode ierr,(*r)(MatMFFD);
103e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)mat->data;
104ace3abfcSBarry Smith   PetscBool      match;
105e884886fSBarry Smith 
106e884886fSBarry Smith   PetscFunctionBegin;
1070700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
108e884886fSBarry Smith   PetscValidCharPointer(ftype,2);
109e884886fSBarry Smith 
110251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)mat,MATMFFD,&match);CHKERRQ(ierr);
1119a13eae7SJed Brown   if (!match) PetscFunctionReturn(0);
1129a13eae7SJed Brown 
113e884886fSBarry Smith   /* already set, so just return */
114251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)ctx,ftype,&match);CHKERRQ(ierr);
115e884886fSBarry Smith   if (match) PetscFunctionReturn(0);
116e884886fSBarry Smith 
117e884886fSBarry Smith   /* destroy the old one if it exists */
118e884886fSBarry Smith   if (ctx->ops->destroy) {
119e884886fSBarry Smith     ierr = (*ctx->ops->destroy)(ctx);CHKERRQ(ierr);
120e884886fSBarry Smith   }
121e884886fSBarry Smith 
1221c9cd337SJed Brown   ierr =  PetscFunctionListFind(MatMFFDList,ftype,&r);CHKERRQ(ierr);
123e32f2f54SBarry Smith   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown MatMFFD type %s given",ftype);
124e884886fSBarry Smith   ierr = (*r)(ctx);CHKERRQ(ierr);
125e884886fSBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)ctx,ftype);CHKERRQ(ierr);
126e884886fSBarry Smith   PetscFunctionReturn(0);
127e884886fSBarry Smith }
128e884886fSBarry Smith 
1295d21e1f8SStefano Zampini static PetscErrorCode MatGetDiagonal_MFFD(Mat,Vec);
1305d21e1f8SStefano Zampini 
131e884886fSBarry Smith typedef PetscErrorCode (*FCN1)(void*,Vec); /* force argument to next function to not be extern C*/
13240244768SBarry Smith static PetscErrorCode  MatMFFDSetFunctioniBase_MFFD(Mat mat,FCN1 func)
133e884886fSBarry Smith {
134e884886fSBarry Smith   MatMFFD ctx = (MatMFFD)mat->data;
135e884886fSBarry Smith 
136e884886fSBarry Smith   PetscFunctionBegin;
137e884886fSBarry Smith   ctx->funcisetbase = func;
1385d21e1f8SStefano Zampini   /* allow users to compose their own getdiagonal and allow MatHasOperation return false if the two functions pointers are not set */
1395d21e1f8SStefano Zampini   if (!mat->ops->getdiagonal || mat->ops->getdiagonal == MatGetDiagonal_MFFD) {
1405d21e1f8SStefano Zampini     mat->ops->getdiagonal = mat->ops->getdiagonal ? ( func ? mat->ops->getdiagonal : NULL) :
1415d21e1f8SStefano Zampini                                                     ( (func && ctx->funci) ? MatGetDiagonal_MFFD : NULL);
1425d21e1f8SStefano Zampini   }
143e884886fSBarry Smith   PetscFunctionReturn(0);
144e884886fSBarry Smith }
145e884886fSBarry Smith 
146e884886fSBarry Smith typedef PetscErrorCode (*FCN2)(void*,PetscInt,Vec,PetscScalar*); /* force argument to next function to not be extern C*/
14740244768SBarry Smith static PetscErrorCode  MatMFFDSetFunctioni_MFFD(Mat mat,FCN2 funci)
148e884886fSBarry Smith {
149e884886fSBarry Smith   MatMFFD ctx = (MatMFFD)mat->data;
150e884886fSBarry Smith 
151e884886fSBarry Smith   PetscFunctionBegin;
152e884886fSBarry Smith   ctx->funci = funci;
1535d21e1f8SStefano Zampini   /* allow users to compose their own getdiagonal and allow MatHasOperation return false if the two functions pointers are not set */
1545d21e1f8SStefano Zampini   if (!mat->ops->getdiagonal || mat->ops->getdiagonal == MatGetDiagonal_MFFD) {
1555d21e1f8SStefano Zampini     mat->ops->getdiagonal = mat->ops->getdiagonal ? ( funci ? mat->ops->getdiagonal : NULL) :
1565d21e1f8SStefano Zampini                                                     ( (funci && ctx->funcisetbase) ? MatGetDiagonal_MFFD : NULL);
1575d21e1f8SStefano Zampini   }
158e884886fSBarry Smith   PetscFunctionReturn(0);
159e884886fSBarry Smith }
160e884886fSBarry Smith 
16140244768SBarry Smith static PetscErrorCode  MatMFFDResetHHistory_MFFD(Mat J)
162bc0f08ceSBarry Smith {
163bc0f08ceSBarry Smith   MatMFFD ctx = (MatMFFD)J->data;
164bc0f08ceSBarry Smith 
165bc0f08ceSBarry Smith   PetscFunctionBegin;
166bc0f08ceSBarry Smith   ctx->ncurrenth = 0;
167bc0f08ceSBarry Smith   PetscFunctionReturn(0);
168bc0f08ceSBarry Smith }
169e884886fSBarry Smith 
1701c84c290SBarry Smith /*@C
1711c84c290SBarry Smith    MatMFFDRegister - Adds a method to the MatMFFD registry.
1721c84c290SBarry Smith 
1731c84c290SBarry Smith    Not Collective
1741c84c290SBarry Smith 
1751c84c290SBarry Smith    Input Parameters:
1761c84c290SBarry Smith +  name_solver - name of a new user-defined compute-h module
1771c84c290SBarry Smith -  routine_create - routine to create method context
1781c84c290SBarry Smith 
1791c84c290SBarry Smith    Level: developer
1801c84c290SBarry Smith 
1811c84c290SBarry Smith    Notes:
1821c84c290SBarry Smith    MatMFFDRegister() may be called multiple times to add several user-defined solvers.
1831c84c290SBarry Smith 
1841c84c290SBarry Smith    Sample usage:
1851c84c290SBarry Smith .vb
186bdf89e91SBarry Smith    MatMFFDRegister("my_h",MyHCreate);
1871c84c290SBarry Smith .ve
1881c84c290SBarry Smith 
1891c84c290SBarry Smith    Then, your solver can be chosen with the procedural interface via
1901c84c290SBarry Smith $     MatMFFDSetType(mfctx,"my_h")
1911c84c290SBarry Smith    or at runtime via the option
192be7a6d03SBarry Smith $     -mat_mffd_type my_h
1931c84c290SBarry Smith 
1941c84c290SBarry Smith .keywords: MatMFFD, register
1951c84c290SBarry Smith 
1961c84c290SBarry Smith .seealso: MatMFFDRegisterAll(), MatMFFDRegisterDestroy()
1971c84c290SBarry Smith  @*/
198bdf89e91SBarry Smith PetscErrorCode  MatMFFDRegister(const char sname[],PetscErrorCode (*function)(MatMFFD))
199e884886fSBarry Smith {
200e884886fSBarry Smith   PetscErrorCode ierr;
201e884886fSBarry Smith 
202e884886fSBarry Smith   PetscFunctionBegin;
203a240a19fSJed Brown   ierr = PetscFunctionListAdd(&MatMFFDList,sname,function);CHKERRQ(ierr);
204e884886fSBarry Smith   PetscFunctionReturn(0);
205e884886fSBarry Smith }
206e884886fSBarry Smith 
207e884886fSBarry Smith /* ----------------------------------------------------------------------------------------*/
20840244768SBarry Smith static PetscErrorCode MatDestroy_MFFD(Mat mat)
209e884886fSBarry Smith {
210e884886fSBarry Smith   PetscErrorCode ierr;
211e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)mat->data;
212e884886fSBarry Smith 
213e884886fSBarry Smith   PetscFunctionBegin;
2146bf464f9SBarry Smith   ierr = VecDestroy(&ctx->w);CHKERRQ(ierr);
2156bf464f9SBarry Smith   ierr = VecDestroy(&ctx->drscale);CHKERRQ(ierr);
2166bf464f9SBarry Smith   ierr = VecDestroy(&ctx->dlscale);CHKERRQ(ierr);
2176bf464f9SBarry Smith   ierr = VecDestroy(&ctx->dshift);CHKERRQ(ierr);
218c51ad4d4SStefano Zampini   ierr = VecDestroy(&ctx->dshiftw);CHKERRQ(ierr);
219c51ad4d4SStefano Zampini   ierr = VecDestroy(&ctx->current_u);CHKERRQ(ierr);
220cfe22f5eSBarry Smith   if (ctx->current_f_allocated) {
2216bf464f9SBarry Smith     ierr = VecDestroy(&ctx->current_f);CHKERRQ(ierr);
222cfe22f5eSBarry Smith   }
223e884886fSBarry Smith   if (ctx->ops->destroy) {ierr = (*ctx->ops->destroy)(ctx);CHKERRQ(ierr);}
2246bf464f9SBarry Smith   ierr      = PetscHeaderDestroy(&ctx);CHKERRQ(ierr);
225bf0cc555SLisandro Dalcin   mat->data = 0;
226e884886fSBarry Smith 
227bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetBase_C",NULL);CHKERRQ(ierr);
228bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctioniBase_C",NULL);CHKERRQ(ierr);
229bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctioni_C",NULL);CHKERRQ(ierr);
230bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunction_C",NULL);CHKERRQ(ierr);
231bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctionError_C",NULL);CHKERRQ(ierr);
232bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetCheckh_C",NULL);CHKERRQ(ierr);
233bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetPeriod_C",NULL);CHKERRQ(ierr);
234bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDResetHHistory_C",NULL);CHKERRQ(ierr);
235e884886fSBarry Smith   PetscFunctionReturn(0);
236e884886fSBarry Smith }
237e884886fSBarry Smith 
23888b4c220SStefano Zampini static PetscErrorCode MatSetUp_MFFD(Mat mat)
23988b4c220SStefano Zampini {
24088b4c220SStefano Zampini   PetscErrorCode ierr;
24188b4c220SStefano Zampini 
24288b4c220SStefano Zampini   PetscFunctionBegin;
24388b4c220SStefano Zampini   ierr = PetscLayoutSetUp(mat->rmap);CHKERRQ(ierr);
24488b4c220SStefano Zampini   ierr = PetscLayoutSetUp(mat->cmap);CHKERRQ(ierr);
24588b4c220SStefano Zampini   PetscFunctionReturn(0);
24688b4c220SStefano Zampini }
24788b4c220SStefano Zampini 
248e884886fSBarry Smith /*
249e884886fSBarry Smith    MatMFFDView_MFFD - Views matrix-free parameters.
250e884886fSBarry Smith 
251e884886fSBarry Smith */
25240244768SBarry Smith static PetscErrorCode MatView_MFFD(Mat J,PetscViewer viewer)
253e884886fSBarry Smith {
254e884886fSBarry Smith   PetscErrorCode ierr;
255e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)J->data;
256a283635eSDmitry Karpeev   PetscBool      iascii, viewbase, viewfunction;
257a283635eSDmitry Karpeev   const char     *prefix;
258e884886fSBarry Smith 
259e884886fSBarry Smith   PetscFunctionBegin;
260251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
261e884886fSBarry Smith   if (iascii) {
262a283635eSDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"Matrix-free approximation:\n");CHKERRQ(ierr);
263a283635eSDmitry Karpeev     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
26457622a8eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"err=%g (relative error in function evaluation)\n",(double)ctx->error_rel);CHKERRQ(ierr);
2657adad957SLisandro Dalcin     if (!((PetscObject)ctx)->type_name) {
266e884886fSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"The compute h routine has not yet been set\n");CHKERRQ(ierr);
267e884886fSBarry Smith     } else {
2687adad957SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"Using %s compute h routine\n",((PetscObject)ctx)->type_name);CHKERRQ(ierr);
269e884886fSBarry Smith     }
270e884886fSBarry Smith     if (ctx->ops->view) {
271e884886fSBarry Smith       ierr = (*ctx->ops->view)(ctx,viewer);CHKERRQ(ierr);
272e884886fSBarry Smith     }
273a283635eSDmitry Karpeev     ierr = PetscObjectGetOptionsPrefix((PetscObject)J, &prefix);CHKERRQ(ierr);
274a283635eSDmitry Karpeev 
275c5929fdfSBarry Smith     ierr = PetscOptionsHasName(((PetscObject)J)->options,prefix, "-mat_mffd_view_base", &viewbase);CHKERRQ(ierr);
276a283635eSDmitry Karpeev     if (viewbase) {
277a283635eSDmitry Karpeev       ierr = PetscViewerASCIIPrintf(viewer, "Base:\n");CHKERRQ(ierr);
278a283635eSDmitry Karpeev       ierr = VecView(ctx->current_u, viewer);CHKERRQ(ierr);
279a283635eSDmitry Karpeev     }
280c5929fdfSBarry Smith     ierr = PetscOptionsHasName(((PetscObject)J)->options,prefix, "-mat_mffd_view_function", &viewfunction);CHKERRQ(ierr);
281a283635eSDmitry Karpeev     if (viewfunction) {
282a283635eSDmitry Karpeev       ierr = PetscViewerASCIIPrintf(viewer, "Function:\n");CHKERRQ(ierr);
283a283635eSDmitry Karpeev       ierr = VecView(ctx->current_f, viewer);CHKERRQ(ierr);
284a283635eSDmitry Karpeev     }
285a283635eSDmitry Karpeev     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
28611aeaf0aSBarry Smith   }
287e884886fSBarry Smith   PetscFunctionReturn(0);
288e884886fSBarry Smith }
289e884886fSBarry Smith 
290e884886fSBarry Smith /*
291e884886fSBarry Smith    MatAssemblyEnd_MFFD - Resets the ctx->ncurrenth to zero. This
292e884886fSBarry Smith    allows the user to indicate the beginning of a new linear solve by calling
293e884886fSBarry Smith    MatAssemblyXXX() on the matrix free matrix. This then allows the
2941d0fab5eSBarry Smith    MatCreateMFFD_WP() to properly compute ||U|| only the first time
295e884886fSBarry Smith    in the linear solver rather than every time.
2965a576424SJed Brown 
297cc2e6a90SBarry Smith    This function is referenced directly from MatAssemblyEnd_SNESMF(), which may be in a different shared library hence
298cc2e6a90SBarry Smith    it must be labeled as PETSC_EXTERN
299e884886fSBarry Smith */
3005a576424SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyEnd_MFFD(Mat J,MatAssemblyType mt)
301e884886fSBarry Smith {
302e884886fSBarry Smith   PetscErrorCode ierr;
303e884886fSBarry Smith   MatMFFD        j = (MatMFFD)J->data;
304e884886fSBarry Smith 
305e884886fSBarry Smith   PetscFunctionBegin;
306e884886fSBarry Smith   ierr      = MatMFFDResetHHistory(J);CHKERRQ(ierr);
307e884886fSBarry Smith   j->vshift = 0.0;
308e884886fSBarry Smith   j->vscale = 1.0;
309e884886fSBarry Smith   PetscFunctionReturn(0);
310e884886fSBarry Smith }
311e884886fSBarry Smith 
312e884886fSBarry Smith /*
313e884886fSBarry Smith   MatMult_MFFD - Default matrix-free form for Jacobian-vector product, y = F'(u)*a:
314e884886fSBarry Smith 
315e884886fSBarry Smith         y ~= (F(u + ha) - F(u))/h,
316e884886fSBarry Smith   where F = nonlinear function, as set by SNESSetFunction()
317e884886fSBarry Smith         u = current iterate
318e884886fSBarry Smith         h = difference interval
319e884886fSBarry Smith */
32040244768SBarry Smith static PetscErrorCode MatMult_MFFD(Mat mat,Vec a,Vec y)
321e884886fSBarry Smith {
322e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)mat->data;
323e884886fSBarry Smith   PetscScalar    h;
324e884886fSBarry Smith   Vec            w,U,F;
325e884886fSBarry Smith   PetscErrorCode ierr;
326ace3abfcSBarry Smith   PetscBool      zeroa;
327e884886fSBarry Smith 
328e884886fSBarry Smith   PetscFunctionBegin;
329ce94432eSBarry Smith   if (!ctx->current_u) SETERRQ(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");
330e884886fSBarry Smith   /* We log matrix-free matrix-vector products separately, so that we can
331e884886fSBarry Smith      separate the performance monitoring from the cases that use conventional
332e884886fSBarry Smith      storage.  We may eventually modify event logging to associate events
333e884886fSBarry Smith      with particular objects, hence alleviating the more general problem. */
334e884886fSBarry Smith   ierr = PetscLogEventBegin(MATMFFD_Mult,a,y,0,0);CHKERRQ(ierr);
335e884886fSBarry Smith 
336e884886fSBarry Smith   w = ctx->w;
337e884886fSBarry Smith   U = ctx->current_u;
3383ec795f1SBarry Smith   F = ctx->current_f;
339e884886fSBarry Smith   /*
340e884886fSBarry Smith       Compute differencing parameter
341e884886fSBarry Smith   */
3424a2f8832SBarry Smith   if (!((PetscObject)ctx)->type_name) {
343e884886fSBarry Smith     ierr = MatMFFDSetType(mat,MATMFFD_WP);CHKERRQ(ierr);
3449c6ac3b3SBarry Smith     ierr = MatSetFromOptions(mat);CHKERRQ(ierr);
345e884886fSBarry Smith   }
346e884886fSBarry Smith   ierr = (*ctx->ops->compute)(ctx,U,a,&h,&zeroa);CHKERRQ(ierr);
347e884886fSBarry Smith   if (zeroa) {
348e884886fSBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
349e884886fSBarry Smith     PetscFunctionReturn(0);
350e884886fSBarry Smith   }
351e884886fSBarry Smith 
35284d44b13SHong Zhang   if (mat->erroriffailure && PetscIsInfOrNanScalar(h)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Computed Nan differencing parameter h");
353e884886fSBarry Smith   if (ctx->checkh) {
354e884886fSBarry Smith     ierr = (*ctx->checkh)(ctx->checkhctx,U,a,&h);CHKERRQ(ierr);
355e884886fSBarry Smith   }
356e884886fSBarry Smith 
357e884886fSBarry Smith   /* keep a record of the current differencing parameter h */
358e884886fSBarry Smith   ctx->currenth = h;
359e884886fSBarry Smith #if defined(PETSC_USE_COMPLEX)
36057622a8eSBarry Smith   ierr = PetscInfo2(mat,"Current differencing parameter: %g + %g i\n",(double)PetscRealPart(h),(double)PetscImaginaryPart(h));CHKERRQ(ierr);
361e884886fSBarry Smith #else
362e884886fSBarry Smith   ierr = PetscInfo1(mat,"Current differencing parameter: %15.12e\n",h);CHKERRQ(ierr);
363e884886fSBarry Smith #endif
364e884886fSBarry Smith   if (ctx->historyh && ctx->ncurrenth < ctx->maxcurrenth) {
365e884886fSBarry Smith     ctx->historyh[ctx->ncurrenth] = h;
366e884886fSBarry Smith   }
367e884886fSBarry Smith   ctx->ncurrenth++;
368e884886fSBarry Smith 
369e884886fSBarry Smith   /* w = u + ha */
370c3bb7e23SBarry Smith   if (ctx->drscale) {
371c3bb7e23SBarry Smith     ierr = VecPointwiseMult(ctx->drscale,a,U);CHKERRQ(ierr);
372c3bb7e23SBarry Smith     ierr = VecAYPX(U,h,w);CHKERRQ(ierr);
373c3bb7e23SBarry Smith   } else {
374e884886fSBarry Smith     ierr = VecWAXPY(w,h,a,U);CHKERRQ(ierr);
375c3bb7e23SBarry Smith   }
376e884886fSBarry Smith 
3771b797266SDmitry Karpeev   /* compute func(U) as base for differencing; only needed first time in and not when provided by user */
3781b797266SDmitry Karpeev   if (ctx->ncurrenth == 1 && ctx->current_f_allocated) {
379e884886fSBarry Smith     ierr = (*ctx->func)(ctx->funcctx,U,F);CHKERRQ(ierr);
38052121784SBarry Smith   }
381e884886fSBarry Smith   ierr = (*ctx->func)(ctx->funcctx,w,y);CHKERRQ(ierr);
382e884886fSBarry Smith 
383e884886fSBarry Smith   ierr = VecAXPY(y,-1.0,F);CHKERRQ(ierr);
384e884886fSBarry Smith   ierr = VecScale(y,1.0/h);CHKERRQ(ierr);
385e884886fSBarry Smith 
386c35ec66cSMatthew G Knepley   if ((ctx->vshift != 0.0) || (ctx->vscale != 1.0)) {
387e884886fSBarry Smith     ierr = VecAXPBY(y,ctx->vshift,ctx->vscale,a);CHKERRQ(ierr);
388c3bb7e23SBarry Smith   }
389c3bb7e23SBarry Smith   if (ctx->dlscale) {
390c3bb7e23SBarry Smith     ierr = VecPointwiseMult(y,ctx->dlscale,y);CHKERRQ(ierr);
391c3bb7e23SBarry Smith   }
392c3bb7e23SBarry Smith   if (ctx->dshift) {
393c51ad4d4SStefano Zampini     if (!ctx->dshiftw) {
394c51ad4d4SStefano Zampini       ierr = VecDuplicate(y,&ctx->dshiftw);CHKERRQ(ierr);
395c51ad4d4SStefano Zampini     }
396c51ad4d4SStefano Zampini     ierr = VecPointwiseMult(ctx->dshift,a,ctx->dshiftw);CHKERRQ(ierr);
397c51ad4d4SStefano Zampini     ierr = VecAXPY(y,1.0,ctx->dshiftw);CHKERRQ(ierr);
398c3bb7e23SBarry Smith   }
399e884886fSBarry Smith 
40039601f4eSBarry Smith   if (mat->nullsp) {ierr = MatNullSpaceRemove(mat->nullsp,y);CHKERRQ(ierr);}
401e884886fSBarry Smith 
402e884886fSBarry Smith   ierr = PetscLogEventEnd(MATMFFD_Mult,a,y,0,0);CHKERRQ(ierr);
403e884886fSBarry Smith   PetscFunctionReturn(0);
404e884886fSBarry Smith }
405e884886fSBarry Smith 
406e884886fSBarry Smith /*
407e884886fSBarry Smith   MatGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix
408e884886fSBarry Smith 
409e884886fSBarry Smith         y ~= (F(u + ha) - F(u))/h,
410e884886fSBarry Smith   where F = nonlinear function, as set by SNESSetFunction()
411e884886fSBarry Smith         u = current iterate
412e884886fSBarry Smith         h = difference interval
413e884886fSBarry Smith */
4145d21e1f8SStefano Zampini PetscErrorCode MatGetDiagonal_MFFD(Mat mat,Vec a)
415e884886fSBarry Smith {
416e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)mat->data;
417e884886fSBarry Smith   PetscScalar    h,*aa,*ww,v;
418e884886fSBarry Smith   PetscReal      epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON;
419e884886fSBarry Smith   Vec            w,U;
420e884886fSBarry Smith   PetscErrorCode ierr;
421e884886fSBarry Smith   PetscInt       i,rstart,rend;
422e884886fSBarry Smith 
423e884886fSBarry Smith   PetscFunctionBegin;
424*94eb4328SStefano Zampini   if (!ctx->funci) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing Functioni function");
425*94eb4328SStefano Zampini   if (!ctx->funcisetbase) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing FunctioniBase function");
426e884886fSBarry Smith   w    = ctx->w;
427e884886fSBarry Smith   U    = ctx->current_u;
428e884886fSBarry Smith   ierr = (*ctx->func)(ctx->funcctx,U,a);CHKERRQ(ierr);
429e884886fSBarry Smith   ierr = (*ctx->funcisetbase)(ctx->funcctx,U);CHKERRQ(ierr);
430e884886fSBarry Smith   ierr = VecCopy(U,w);CHKERRQ(ierr);
431e884886fSBarry Smith 
432e884886fSBarry Smith   ierr = VecGetOwnershipRange(a,&rstart,&rend);CHKERRQ(ierr);
433e884886fSBarry Smith   ierr = VecGetArray(a,&aa);CHKERRQ(ierr);
434e884886fSBarry Smith   for (i=rstart; i<rend; i++) {
435e884886fSBarry Smith     ierr = VecGetArray(w,&ww);CHKERRQ(ierr);
436e884886fSBarry Smith     h    = ww[i-rstart];
437e884886fSBarry Smith     if (h == 0.0) h = 1.0;
438e884886fSBarry Smith     if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0)     h = umin;
439e884886fSBarry Smith     else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin;
440e884886fSBarry Smith     h *= epsilon;
441e884886fSBarry Smith 
442e884886fSBarry Smith     ww[i-rstart] += h;
443e884886fSBarry Smith     ierr          = VecRestoreArray(w,&ww);CHKERRQ(ierr);
444e884886fSBarry Smith     ierr          = (*ctx->funci)(ctx->funcctx,i,w,&v);CHKERRQ(ierr);
445e884886fSBarry Smith     aa[i-rstart]  = (v - aa[i-rstart])/h;
446e884886fSBarry Smith 
447e884886fSBarry Smith     /* possibly shift and scale result */
448c35ec66cSMatthew G Knepley     if ((ctx->vshift != 0.0) || (ctx->vscale != 1.0)) {
449e884886fSBarry Smith       aa[i - rstart] = ctx->vshift + ctx->vscale*aa[i-rstart];
450c3bb7e23SBarry Smith     }
451e884886fSBarry Smith 
452e884886fSBarry Smith     ierr          = VecGetArray(w,&ww);CHKERRQ(ierr);
453e884886fSBarry Smith     ww[i-rstart] -= h;
454e884886fSBarry Smith     ierr          = VecRestoreArray(w,&ww);CHKERRQ(ierr);
455e884886fSBarry Smith   }
456e884886fSBarry Smith   ierr = VecRestoreArray(a,&aa);CHKERRQ(ierr);
457e884886fSBarry Smith   PetscFunctionReturn(0);
458e884886fSBarry Smith }
459e884886fSBarry Smith 
46040244768SBarry Smith static PetscErrorCode MatDiagonalScale_MFFD(Mat mat,Vec ll,Vec rr)
461c3bb7e23SBarry Smith {
462c3bb7e23SBarry Smith   MatMFFD        aij = (MatMFFD)mat->data;
463c3bb7e23SBarry Smith   PetscErrorCode ierr;
464c3bb7e23SBarry Smith 
465c3bb7e23SBarry Smith   PetscFunctionBegin;
466c3bb7e23SBarry Smith   if (ll && !aij->dlscale) {
467c3bb7e23SBarry Smith     ierr = VecDuplicate(ll,&aij->dlscale);CHKERRQ(ierr);
468c3bb7e23SBarry Smith   }
469c3bb7e23SBarry Smith   if (rr && !aij->drscale) {
470c3bb7e23SBarry Smith     ierr = VecDuplicate(rr,&aij->drscale);CHKERRQ(ierr);
471c3bb7e23SBarry Smith   }
472c3bb7e23SBarry Smith   if (ll) {
473c3bb7e23SBarry Smith     ierr = VecCopy(ll,aij->dlscale);CHKERRQ(ierr);
474c3bb7e23SBarry Smith   }
475c3bb7e23SBarry Smith   if (rr) {
476c3bb7e23SBarry Smith     ierr = VecCopy(rr,aij->drscale);CHKERRQ(ierr);
477c3bb7e23SBarry Smith   }
478c3bb7e23SBarry Smith   PetscFunctionReturn(0);
479c3bb7e23SBarry Smith }
480c3bb7e23SBarry Smith 
48140244768SBarry Smith static PetscErrorCode MatDiagonalSet_MFFD(Mat mat,Vec ll,InsertMode mode)
482c3bb7e23SBarry Smith {
483c3bb7e23SBarry Smith   MatMFFD        aij = (MatMFFD)mat->data;
484c3bb7e23SBarry Smith   PetscErrorCode ierr;
485c3bb7e23SBarry Smith 
486c3bb7e23SBarry Smith   PetscFunctionBegin;
487ce94432eSBarry Smith   if (mode == INSERT_VALUES) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"No diagonal set with INSERT_VALUES");
488c3bb7e23SBarry Smith   if (!aij->dshift) {
489c3bb7e23SBarry Smith     ierr = VecDuplicate(ll,&aij->dshift);CHKERRQ(ierr);
490c3bb7e23SBarry Smith   }
491c3bb7e23SBarry Smith   ierr = VecAXPY(aij->dshift,1.0,ll);CHKERRQ(ierr);
492c3bb7e23SBarry Smith   PetscFunctionReturn(0);
493c3bb7e23SBarry Smith }
494c3bb7e23SBarry Smith 
49540244768SBarry Smith static PetscErrorCode MatShift_MFFD(Mat Y,PetscScalar a)
496e884886fSBarry Smith {
497e884886fSBarry Smith   MatMFFD shell = (MatMFFD)Y->data;
4985fd66863SKarl Rupp 
499e884886fSBarry Smith   PetscFunctionBegin;
500e884886fSBarry Smith   shell->vshift += a;
501e884886fSBarry Smith   PetscFunctionReturn(0);
502e884886fSBarry Smith }
503e884886fSBarry Smith 
50440244768SBarry Smith static PetscErrorCode MatScale_MFFD(Mat Y,PetscScalar a)
505e884886fSBarry Smith {
506e884886fSBarry Smith   MatMFFD shell = (MatMFFD)Y->data;
5075fd66863SKarl Rupp 
508e884886fSBarry Smith   PetscFunctionBegin;
509e884886fSBarry Smith   shell->vscale *= a;
510e884886fSBarry Smith   PetscFunctionReturn(0);
511e884886fSBarry Smith }
512e884886fSBarry Smith 
513d8fa6a54SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDSetBase_MFFD(Mat J,Vec U,Vec F)
514e884886fSBarry Smith {
515e884886fSBarry Smith   PetscErrorCode ierr;
516e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)J->data;
517e884886fSBarry Smith 
518e884886fSBarry Smith   PetscFunctionBegin;
519e884886fSBarry Smith   ierr = MatMFFDResetHHistory(J);CHKERRQ(ierr);
520c51ad4d4SStefano Zampini   if (!ctx->current_u) {
521c51ad4d4SStefano Zampini     ierr = VecDuplicate(U,&ctx->current_u);CHKERRQ(ierr);
522c51ad4d4SStefano Zampini     ierr = VecLockPush(ctx->current_u);CHKERRQ(ierr);
523c51ad4d4SStefano Zampini   }
524c51ad4d4SStefano Zampini   ierr = VecLockPop(ctx->current_u);CHKERRQ(ierr);
525c51ad4d4SStefano Zampini   ierr = VecCopy(U,ctx->current_u);CHKERRQ(ierr);
526c51ad4d4SStefano Zampini   ierr = VecLockPush(ctx->current_u);CHKERRQ(ierr);
52752121784SBarry Smith   if (F) {
5286bf464f9SBarry Smith     if (ctx->current_f_allocated) {ierr = VecDestroy(&ctx->current_f);CHKERRQ(ierr);}
5293ec795f1SBarry Smith     ctx->current_f           = F;
530cfe22f5eSBarry Smith     ctx->current_f_allocated = PETSC_FALSE;
531cfe22f5eSBarry Smith   } else if (!ctx->current_f_allocated) {
53206c4b6baSBarry Smith     ierr = MatCreateVecs(J,NULL,&ctx->current_f);CHKERRQ(ierr);
5332205254eSKarl Rupp 
534cfe22f5eSBarry Smith     ctx->current_f_allocated = PETSC_TRUE;
53552121784SBarry Smith   }
536e884886fSBarry Smith   if (!ctx->w) {
537e884886fSBarry Smith     ierr = VecDuplicate(ctx->current_u,&ctx->w);CHKERRQ(ierr);
538e884886fSBarry Smith   }
539e884886fSBarry Smith   J->assembled = PETSC_TRUE;
540e884886fSBarry Smith   PetscFunctionReturn(0);
541e884886fSBarry Smith }
5425a576424SJed Brown 
543e884886fSBarry Smith typedef PetscErrorCode (*FCN3)(void*,Vec,Vec,PetscScalar*); /* force argument to next function to not be extern C*/
544b2573a8aSBarry Smith 
54540244768SBarry Smith static PetscErrorCode  MatMFFDSetCheckh_MFFD(Mat J,FCN3 fun,void *ectx)
546e884886fSBarry Smith {
547e884886fSBarry Smith   MatMFFD ctx = (MatMFFD)J->data;
548e884886fSBarry Smith 
549e884886fSBarry Smith   PetscFunctionBegin;
550e884886fSBarry Smith   ctx->checkh    = fun;
551e884886fSBarry Smith   ctx->checkhctx = ectx;
552e884886fSBarry Smith   PetscFunctionReturn(0);
553e884886fSBarry Smith }
554e884886fSBarry Smith 
5556aa9148fSLisandro Dalcin /*@C
5566aa9148fSLisandro Dalcin    MatMFFDSetOptionsPrefix - Sets the prefix used for searching for all
5576aa9148fSLisandro Dalcin    MatMFFD options in the database.
5586aa9148fSLisandro Dalcin 
5596aa9148fSLisandro Dalcin    Collective on Mat
5606aa9148fSLisandro Dalcin 
5616aa9148fSLisandro Dalcin    Input Parameter:
5626aa9148fSLisandro Dalcin +  A - the Mat context
5636aa9148fSLisandro Dalcin -  prefix - the prefix to prepend to all option names
5646aa9148fSLisandro Dalcin 
5656aa9148fSLisandro Dalcin    Notes:
5666aa9148fSLisandro Dalcin    A hyphen (-) must NOT be given at the beginning of the prefix name.
5676aa9148fSLisandro Dalcin    The first character of all runtime options is AUTOMATICALLY the hyphen.
5686aa9148fSLisandro Dalcin 
5696aa9148fSLisandro Dalcin    Level: advanced
5706aa9148fSLisandro Dalcin 
5716aa9148fSLisandro Dalcin .keywords: SNES, matrix-free, parameters
5726aa9148fSLisandro Dalcin 
573755b7f64SBarry Smith .seealso: MatSetFromOptions(), MatCreateSNESMF(), MatCreateMFFD()
5746aa9148fSLisandro Dalcin @*/
5757087cfbeSBarry Smith PetscErrorCode  MatMFFDSetOptionsPrefix(Mat mat,const char prefix[])
5766aa9148fSLisandro Dalcin 
5776aa9148fSLisandro Dalcin {
5780298fd71SBarry Smith   MatMFFD        mfctx = mat ? (MatMFFD)mat->data : (MatMFFD)NULL;
5796aa9148fSLisandro Dalcin   PetscErrorCode ierr;
5805fd66863SKarl Rupp 
5816aa9148fSLisandro Dalcin   PetscFunctionBegin;
5820700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
5830700a824SBarry Smith   PetscValidHeaderSpecific(mfctx,MATMFFD_CLASSID,1);
5846aa9148fSLisandro Dalcin   ierr = PetscObjectSetOptionsPrefix((PetscObject)mfctx,prefix);CHKERRQ(ierr);
5856aa9148fSLisandro Dalcin   PetscFunctionReturn(0);
5866aa9148fSLisandro Dalcin }
5876aa9148fSLisandro Dalcin 
58840244768SBarry Smith static PetscErrorCode  MatSetFromOptions_MFFD(PetscOptionItems *PetscOptionsObject,Mat mat)
589e884886fSBarry Smith {
59071f08665SBarry Smith   MatMFFD        mfctx = (MatMFFD)mat->data;
591e884886fSBarry Smith   PetscErrorCode ierr;
592ace3abfcSBarry Smith   PetscBool      flg;
593e884886fSBarry Smith   char           ftype[256];
594e884886fSBarry Smith 
595e884886fSBarry Smith   PetscFunctionBegin;
5960700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
5970700a824SBarry Smith   PetscValidHeaderSpecific(mfctx,MATMFFD_CLASSID,1);
5983194b578SJed Brown   ierr = PetscObjectOptionsBegin((PetscObject)mfctx);CHKERRQ(ierr);
599a264d7a6SBarry Smith   ierr = PetscOptionsFList("-mat_mffd_type","Matrix free type","MatMFFDSetType",MatMFFDList,((PetscObject)mfctx)->type_name,ftype,256,&flg);CHKERRQ(ierr);
600e884886fSBarry Smith   if (flg) {
601e884886fSBarry Smith     ierr = MatMFFDSetType(mat,ftype);CHKERRQ(ierr);
602e884886fSBarry Smith   }
603e884886fSBarry Smith 
604e884886fSBarry Smith   ierr = PetscOptionsReal("-mat_mffd_err","set sqrt relative error in function","MatMFFDSetFunctionError",mfctx->error_rel,&mfctx->error_rel,0);CHKERRQ(ierr);
605e884886fSBarry Smith   ierr = PetscOptionsInt("-mat_mffd_period","how often h is recomputed","MatMFFDSetPeriod",mfctx->recomputeperiod,&mfctx->recomputeperiod,0);CHKERRQ(ierr);
606e884886fSBarry Smith 
60790d69ab7SBarry Smith   flg  = PETSC_FALSE;
6080298fd71SBarry Smith   ierr = PetscOptionsBool("-mat_mffd_check_positivity","Insure that U + h*a is nonnegative","MatMFFDSetCheckh",flg,&flg,NULL);CHKERRQ(ierr);
609e884886fSBarry Smith   if (flg) {
610e884886fSBarry Smith     ierr = MatMFFDSetCheckh(mat,MatMFFDCheckPositivity,0);CHKERRQ(ierr);
611e884886fSBarry Smith   }
612e884886fSBarry Smith   if (mfctx->ops->setfromoptions) {
613e55864a3SBarry Smith     ierr = (*mfctx->ops->setfromoptions)(PetscOptionsObject,mfctx);CHKERRQ(ierr);
614e884886fSBarry Smith   }
615e884886fSBarry Smith   ierr = PetscOptionsEnd();CHKERRQ(ierr);
616e884886fSBarry Smith   PetscFunctionReturn(0);
617e884886fSBarry Smith }
618e884886fSBarry Smith 
61940244768SBarry Smith static PetscErrorCode  MatMFFDSetPeriod_MFFD(Mat mat,PetscInt period)
620bc0f08ceSBarry Smith {
621bc0f08ceSBarry Smith   MatMFFD ctx = (MatMFFD)mat->data;
622bc0f08ceSBarry Smith 
623bc0f08ceSBarry Smith   PetscFunctionBegin;
624bc0f08ceSBarry Smith   ctx->recomputeperiod = period;
625bc0f08ceSBarry Smith   PetscFunctionReturn(0);
626bc0f08ceSBarry Smith }
627bc0f08ceSBarry Smith 
62840244768SBarry Smith static PetscErrorCode  MatMFFDSetFunction_MFFD(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx)
629bc0f08ceSBarry Smith {
630bc0f08ceSBarry Smith   MatMFFD ctx = (MatMFFD)mat->data;
631bc0f08ceSBarry Smith 
632bc0f08ceSBarry Smith   PetscFunctionBegin;
633bc0f08ceSBarry Smith   ctx->func    = func;
634bc0f08ceSBarry Smith   ctx->funcctx = funcctx;
635bc0f08ceSBarry Smith   PetscFunctionReturn(0);
636bc0f08ceSBarry Smith }
637bc0f08ceSBarry Smith 
63840244768SBarry Smith static PetscErrorCode  MatMFFDSetFunctionError_MFFD(Mat mat,PetscReal error)
639bc0f08ceSBarry Smith {
640bc0f08ceSBarry Smith   MatMFFD ctx = (MatMFFD)mat->data;
641bc0f08ceSBarry Smith 
642bc0f08ceSBarry Smith   PetscFunctionBegin;
643bc0f08ceSBarry Smith   if (error != PETSC_DEFAULT) ctx->error_rel = error;
644bc0f08ceSBarry Smith   PetscFunctionReturn(0);
645bc0f08ceSBarry Smith }
646bc0f08ceSBarry Smith 
6473b49f96aSBarry Smith static PetscErrorCode MatMissingDiagonal_MFFD(Mat A,PetscBool  *missing,PetscInt *d)
6483b49f96aSBarry Smith {
6493b49f96aSBarry Smith   PetscFunctionBegin;
6503b49f96aSBarry Smith   *missing = PETSC_FALSE;
6513b49f96aSBarry Smith   PetscFunctionReturn(0);
6523b49f96aSBarry Smith }
6533b49f96aSBarry Smith 
654e884886fSBarry Smith /*MC
655e884886fSBarry Smith   MATMFFD - MATMFFD = "mffd" - A matrix free matrix type.
656e884886fSBarry Smith 
657e884886fSBarry Smith   Level: advanced
658e884886fSBarry Smith 
659755b7f64SBarry Smith .seealso: MatCreateMFFD(), MatCreateSNESMF(), MatMFFDSetFunction(), MatMFFDSetType(),
660755b7f64SBarry Smith           MatMFFDSetFunctionError(), MatMFFDDSSetUmin(), MatMFFDSetFunction()
661755b7f64SBarry Smith           MatMFFDSetHHistory(), MatMFFDResetHHistory(), MatCreateSNESMF(),
662755b7f64SBarry Smith           MatMFFDGetH(),
663e884886fSBarry Smith M*/
6648cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_MFFD(Mat A)
665e884886fSBarry Smith {
666e884886fSBarry Smith   MatMFFD        mfctx;
667e884886fSBarry Smith   PetscErrorCode ierr;
668e884886fSBarry Smith 
669e884886fSBarry Smith   PetscFunctionBegin;
670607a6623SBarry Smith   ierr = MatMFFDInitializePackage();CHKERRQ(ierr);
671e884886fSBarry Smith 
67273107ff1SLisandro Dalcin   ierr = PetscHeaderCreate(mfctx,MATMFFD_CLASSID,"MatMFFD","Matrix-free Finite Differencing","Mat",PetscObjectComm((PetscObject)A),MatDestroy_MFFD,MatView_MFFD);CHKERRQ(ierr);
6732205254eSKarl Rupp 
674e884886fSBarry Smith   mfctx->error_rel                = PETSC_SQRT_MACHINE_EPSILON;
675e884886fSBarry Smith   mfctx->recomputeperiod          = 1;
676e884886fSBarry Smith   mfctx->count                    = 0;
677e884886fSBarry Smith   mfctx->currenth                 = 0.0;
6780298fd71SBarry Smith   mfctx->historyh                 = NULL;
679e884886fSBarry Smith   mfctx->ncurrenth                = 0;
680e884886fSBarry Smith   mfctx->maxcurrenth              = 0;
6817adad957SLisandro Dalcin   ((PetscObject)mfctx)->type_name = 0;
682e884886fSBarry Smith 
683e884886fSBarry Smith   mfctx->vshift = 0.0;
684e884886fSBarry Smith   mfctx->vscale = 1.0;
685e884886fSBarry Smith 
686e884886fSBarry Smith   /*
687e884886fSBarry Smith      Create the empty data structure to contain compute-h routines.
688e884886fSBarry Smith      These will be filled in below from the command line options or
689e884886fSBarry Smith      a later call with MatMFFDSetType() or if that is not called
690e884886fSBarry Smith      then it will default in the first use of MatMult_MFFD()
691e884886fSBarry Smith   */
692e884886fSBarry Smith   mfctx->ops->compute        = 0;
693e884886fSBarry Smith   mfctx->ops->destroy        = 0;
694e884886fSBarry Smith   mfctx->ops->view           = 0;
695e884886fSBarry Smith   mfctx->ops->setfromoptions = 0;
696e884886fSBarry Smith   mfctx->hctx                = 0;
697e884886fSBarry Smith 
698e884886fSBarry Smith   mfctx->func    = 0;
699e884886fSBarry Smith   mfctx->funcctx = 0;
7000298fd71SBarry Smith   mfctx->w       = NULL;
701e884886fSBarry Smith 
702e884886fSBarry Smith   A->data = mfctx;
703e884886fSBarry Smith 
704e884886fSBarry Smith   A->ops->mult            = MatMult_MFFD;
705e884886fSBarry Smith   A->ops->destroy         = MatDestroy_MFFD;
70688b4c220SStefano Zampini   A->ops->setup           = MatSetUp_MFFD;
707e884886fSBarry Smith   A->ops->view            = MatView_MFFD;
708e884886fSBarry Smith   A->ops->assemblyend     = MatAssemblyEnd_MFFD;
709e884886fSBarry Smith   A->ops->scale           = MatScale_MFFD;
710e884886fSBarry Smith   A->ops->shift           = MatShift_MFFD;
711c3bb7e23SBarry Smith   A->ops->diagonalscale   = MatDiagonalScale_MFFD;
712c3bb7e23SBarry Smith   A->ops->diagonalset     = MatDiagonalSet_MFFD;
7139c6ac3b3SBarry Smith   A->ops->setfromoptions  = MatSetFromOptions_MFFD;
7143b49f96aSBarry Smith   A->ops->missingdiagonal = MatMissingDiagonal_MFFD;
715e884886fSBarry Smith   A->assembled            = PETSC_TRUE;
716e884886fSBarry Smith 
717bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetBase_C",MatMFFDSetBase_MFFD);CHKERRQ(ierr);
718bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctioniBase_C",MatMFFDSetFunctioniBase_MFFD);CHKERRQ(ierr);
719bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctioni_C",MatMFFDSetFunctioni_MFFD);CHKERRQ(ierr);
720bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunction_C",MatMFFDSetFunction_MFFD);CHKERRQ(ierr);
721bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetCheckh_C",MatMFFDSetCheckh_MFFD);CHKERRQ(ierr);
722bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetPeriod_C",MatMFFDSetPeriod_MFFD);CHKERRQ(ierr);
723bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctionError_C",MatMFFDSetFunctionError_MFFD);CHKERRQ(ierr);
724bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDResetHHistory_C",MatMFFDResetHHistory_MFFD);CHKERRQ(ierr);
7252205254eSKarl Rupp 
726e884886fSBarry Smith   mfctx->mat = A;
7272205254eSKarl Rupp 
728e884886fSBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,MATMFFD);CHKERRQ(ierr);
729e884886fSBarry Smith   PetscFunctionReturn(0);
730e884886fSBarry Smith }
731e884886fSBarry Smith 
732e884886fSBarry Smith /*@
733e884886fSBarry Smith    MatCreateMFFD - Creates a matrix-free matrix. See also MatCreateSNESMF()
734e884886fSBarry Smith 
735e884886fSBarry Smith    Collective on Vec
736e884886fSBarry Smith 
737e884886fSBarry Smith    Input Parameters:
738fef1beadSBarry Smith +  comm - MPI communicator
739fef1beadSBarry Smith .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
740fef1beadSBarry Smith            This value should be the same as the local size used in creating the
741fef1beadSBarry Smith            y vector for the matrix-vector product y = Ax.
742fef1beadSBarry Smith .  n - This value should be the same as the local size used in creating the
743fef1beadSBarry Smith        x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have
744fef1beadSBarry Smith        calculated if N is given) For square matrices n is almost always m.
745fef1beadSBarry Smith .  M - number of global rows (or PETSC_DETERMINE to have calculated if m is given)
746fef1beadSBarry Smith -  N - number of global columns (or PETSC_DETERMINE to have calculated if n is given)
747fef1beadSBarry Smith 
748e884886fSBarry Smith 
749e884886fSBarry Smith    Output Parameter:
750e884886fSBarry Smith .  J - the matrix-free matrix
751e884886fSBarry Smith 
7529c6ac3b3SBarry Smith    Options Database Keys: call MatSetFromOptions() to trigger these
7539c6ac3b3SBarry Smith +  -mat_mffd_type - wp or ds (see MATMFFD_WP or MATMFFD_DS)
7549c6ac3b3SBarry Smith -  -mat_mffd_err - square root of estimated relative error in function evaluation
7559c6ac3b3SBarry Smith -  -mat_mffd_period - how often h is recomputed, defaults to 1, everytime
7569c6ac3b3SBarry Smith 
7579c6ac3b3SBarry Smith 
758e884886fSBarry Smith    Level: advanced
759e884886fSBarry Smith 
760e884886fSBarry Smith    Notes:
761e884886fSBarry Smith    The matrix-free matrix context merely contains the function pointers
762e884886fSBarry Smith    and work space for performing finite difference approximations of
763e884886fSBarry Smith    Jacobian-vector products, F'(u)*a,
764e884886fSBarry Smith 
765e884886fSBarry Smith    The default code uses the following approach to compute h
766e884886fSBarry Smith 
767e884886fSBarry Smith .vb
768e884886fSBarry Smith      F'(u)*a = [F(u+h*a) - F(u)]/h where
769e884886fSBarry Smith      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
770e884886fSBarry Smith        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
771e884886fSBarry Smith  where
772e884886fSBarry Smith      error_rel = square root of relative error in function evaluation
773e884886fSBarry Smith      umin = minimum iterate parameter
774e884886fSBarry Smith .ve
775e884886fSBarry Smith 
776ca93e954SBarry Smith    You can call SNESSetJacobian() with MatMFFDComputeJacobian() if you are using matrix and not a different
777ca93e954SBarry Smith    preconditioner matrix
778ca93e954SBarry Smith 
779e884886fSBarry Smith    The user can set the error_rel via MatMFFDSetFunctionError() and
780a7f22e61SSatish Balay    umin via MatMFFDDSSetUmin(); see Users-Manual: ch_snes for details.
781e884886fSBarry Smith 
782e884886fSBarry Smith    The user should call MatDestroy() when finished with the matrix-free
783e884886fSBarry Smith    matrix context.
784e884886fSBarry Smith 
785e884886fSBarry Smith    Options Database Keys:
786e884886fSBarry Smith +  -mat_mffd_err <error_rel> - Sets error_rel
787e884886fSBarry Smith .  -mat_mffd_unim <umin> - Sets umin (for default PETSc routine that computes h only)
788e884886fSBarry Smith -  -mat_mffd_check_positivity
789e884886fSBarry Smith 
790e884886fSBarry Smith .keywords: default, matrix-free, create, matrix
791e884886fSBarry Smith 
7921d0fab5eSBarry Smith .seealso: MatDestroy(), MatMFFDSetFunctionError(), MatMFFDDSSetUmin(), MatMFFDSetFunction()
793e884886fSBarry Smith           MatMFFDSetHHistory(), MatMFFDResetHHistory(), MatCreateSNESMF(),
79481242352SJed Brown           MatMFFDGetH(), MatMFFDRegister(), MatMFFDComputeJacobian()
795e884886fSBarry Smith 
796e884886fSBarry Smith @*/
7977087cfbeSBarry Smith PetscErrorCode  MatCreateMFFD(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *J)
798e884886fSBarry Smith {
799e884886fSBarry Smith   PetscErrorCode ierr;
800e884886fSBarry Smith 
801e884886fSBarry Smith   PetscFunctionBegin;
802e884886fSBarry Smith   ierr = MatCreate(comm,J);CHKERRQ(ierr);
803fef1beadSBarry Smith   ierr = MatSetSizes(*J,m,n,M,N);CHKERRQ(ierr);
804e884886fSBarry Smith   ierr = MatSetType(*J,MATMFFD);CHKERRQ(ierr);
805be7c243fSJed Brown   ierr = MatSetUp(*J);CHKERRQ(ierr);
806e884886fSBarry Smith   PetscFunctionReturn(0);
807e884886fSBarry Smith }
808e884886fSBarry Smith 
809e884886fSBarry Smith /*@
810e884886fSBarry Smith    MatMFFDGetH - Gets the last value that was used as the differencing
811e884886fSBarry Smith    parameter.
812e884886fSBarry Smith 
813e884886fSBarry Smith    Not Collective
814e884886fSBarry Smith 
815e884886fSBarry Smith    Input Parameters:
816e884886fSBarry Smith .  mat - the matrix obtained with MatCreateSNESMF()
817e884886fSBarry Smith 
818e884886fSBarry Smith    Output Paramter:
819e884886fSBarry Smith .  h - the differencing step size
820e884886fSBarry Smith 
821e884886fSBarry Smith    Level: advanced
822e884886fSBarry Smith 
823e884886fSBarry Smith .keywords: SNES, matrix-free, parameters
824e884886fSBarry Smith 
8251d0fab5eSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDSetHHistory(), MatCreateMFFD(), MATMFFD, MatMFFDResetHHistory()
826e884886fSBarry Smith @*/
8277087cfbeSBarry Smith PetscErrorCode  MatMFFDGetH(Mat mat,PetscScalar *h)
828e884886fSBarry Smith {
829e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)mat->data;
830bc0f08ceSBarry Smith   PetscErrorCode ierr;
831bc0f08ceSBarry Smith   PetscBool      match;
832e884886fSBarry Smith 
833e884886fSBarry Smith   PetscFunctionBegin;
83488b4c220SStefano Zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
83588b4c220SStefano Zampini   PetscValidPointer(h,2);
836251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)mat,MATMFFD,&match);CHKERRQ(ierr);
837ce94432eSBarry Smith   if (!match) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONG,"Not a MFFD matrix");
838bc0f08ceSBarry Smith 
839e884886fSBarry Smith   *h = ctx->currenth;
840e884886fSBarry Smith   PetscFunctionReturn(0);
841e884886fSBarry Smith }
842e884886fSBarry Smith 
843e884886fSBarry Smith /*@C
844e884886fSBarry Smith    MatMFFDSetFunction - Sets the function used in applying the matrix free.
845e884886fSBarry Smith 
8463f9fe445SBarry Smith    Logically Collective on Mat
847e884886fSBarry Smith 
848e884886fSBarry Smith    Input Parameters:
84914f633e2SBarry Smith +  mat - the matrix free matrix created via MatCreateSNESMF() or MatCreateMFFD()
850e884886fSBarry Smith .  func - the function to use
851e884886fSBarry Smith -  funcctx - optional function context passed to function
852e884886fSBarry Smith 
85314f633e2SBarry Smith    Calling Sequence of func:
85414f633e2SBarry Smith $     func (void *funcctx, Vec x, Vec f)
85514f633e2SBarry Smith 
85614f633e2SBarry Smith +  funcctx - user provided context
85714f633e2SBarry Smith .  x - input vector
85814f633e2SBarry Smith -  f - computed output function
85914f633e2SBarry Smith 
860e884886fSBarry Smith    Level: advanced
861e884886fSBarry Smith 
862e884886fSBarry Smith    Notes:
863e884886fSBarry Smith     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
864e884886fSBarry Smith     matrix inside your compute Jacobian routine
865e884886fSBarry Smith 
8663ec795f1SBarry Smith     If this is not set then it will use the function set with SNESSetFunction() if MatCreateSNESMF() was used.
867e884886fSBarry Smith 
868e884886fSBarry Smith .keywords: SNES, matrix-free, function
869e884886fSBarry Smith 
8701d0fab5eSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD,
8711d0fab5eSBarry Smith           MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction()
872e884886fSBarry Smith @*/
8737087cfbeSBarry Smith PetscErrorCode  MatMFFDSetFunction(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx)
874e884886fSBarry Smith {
875bc0f08ceSBarry Smith   PetscErrorCode ierr;
876e884886fSBarry Smith 
877e884886fSBarry Smith   PetscFunctionBegin;
87888b4c220SStefano Zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
879bc0f08ceSBarry Smith   ierr = PetscTryMethod(mat,"MatMFFDSetFunction_C",(Mat,PetscErrorCode (*)(void*,Vec,Vec),void*),(mat,func,funcctx));CHKERRQ(ierr);
880e884886fSBarry Smith   PetscFunctionReturn(0);
881e884886fSBarry Smith }
882e884886fSBarry Smith 
883e884886fSBarry Smith /*@C
884e884886fSBarry Smith    MatMFFDSetFunctioni - Sets the function for a single component
885e884886fSBarry Smith 
8863f9fe445SBarry Smith    Logically Collective on Mat
887e884886fSBarry Smith 
888e884886fSBarry Smith    Input Parameters:
889e884886fSBarry Smith +  mat - the matrix free matrix created via MatCreateSNESMF()
890e884886fSBarry Smith -  funci - the function to use
891e884886fSBarry Smith 
892e884886fSBarry Smith    Level: advanced
893e884886fSBarry Smith 
894e884886fSBarry Smith    Notes:
895e884886fSBarry Smith     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
896*94eb4328SStefano Zampini     matrix inside your compute Jacobian routine.
897*94eb4328SStefano Zampini     This function is necessary to compute the diagonal of the matrix.
898e884886fSBarry Smith 
899e884886fSBarry Smith 
900e884886fSBarry Smith .keywords: SNES, matrix-free, function
901e884886fSBarry Smith 
902*94eb4328SStefano Zampini .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction(), MatGetDiagonal()
9031d0fab5eSBarry Smith 
904e884886fSBarry Smith @*/
9057087cfbeSBarry Smith PetscErrorCode  MatMFFDSetFunctioni(Mat mat,PetscErrorCode (*funci)(void*,PetscInt,Vec,PetscScalar*))
906e884886fSBarry Smith {
9074ac538c5SBarry Smith   PetscErrorCode ierr;
908e884886fSBarry Smith 
909e884886fSBarry Smith   PetscFunctionBegin;
9100700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
9114ac538c5SBarry Smith   ierr = PetscTryMethod(mat,"MatMFFDSetFunctioni_C",(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)),(mat,funci));CHKERRQ(ierr);
912e884886fSBarry Smith   PetscFunctionReturn(0);
913e884886fSBarry Smith }
914e884886fSBarry Smith 
915e884886fSBarry Smith /*@C
916e884886fSBarry Smith    MatMFFDSetFunctioniBase - Sets the base vector for a single component function evaluation
917e884886fSBarry Smith 
9183f9fe445SBarry Smith    Logically Collective on Mat
919e884886fSBarry Smith 
920e884886fSBarry Smith    Input Parameters:
921e884886fSBarry Smith +  mat - the matrix free matrix created via MatCreateSNESMF()
922e884886fSBarry Smith -  func - the function to use
923e884886fSBarry Smith 
924e884886fSBarry Smith    Level: advanced
925e884886fSBarry Smith 
926e884886fSBarry Smith    Notes:
927e884886fSBarry Smith     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
928*94eb4328SStefano Zampini     matrix inside your compute Jacobian routine.
929*94eb4328SStefano Zampini     This function is necessary to compute the diagonal of the matrix.
930e884886fSBarry Smith 
931e884886fSBarry Smith 
932e884886fSBarry Smith .keywords: SNES, matrix-free, function
933e884886fSBarry Smith 
934e884886fSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD
935*94eb4328SStefano Zampini           MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction(), MatGetDiagonal()
936e884886fSBarry Smith @*/
9377087cfbeSBarry Smith PetscErrorCode  MatMFFDSetFunctioniBase(Mat mat,PetscErrorCode (*func)(void*,Vec))
938e884886fSBarry Smith {
9394ac538c5SBarry Smith   PetscErrorCode ierr;
940e884886fSBarry Smith 
941e884886fSBarry Smith   PetscFunctionBegin;
9420700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
9434ac538c5SBarry Smith   ierr = PetscTryMethod(mat,"MatMFFDSetFunctioniBase_C",(Mat,PetscErrorCode (*)(void*,Vec)),(mat,func));CHKERRQ(ierr);
944e884886fSBarry Smith   PetscFunctionReturn(0);
945e884886fSBarry Smith }
946e884886fSBarry Smith 
947e884886fSBarry Smith /*@
948e884886fSBarry Smith    MatMFFDSetPeriod - Sets how often h is recomputed, by default it is everytime
949e884886fSBarry Smith 
9503f9fe445SBarry Smith    Logically Collective on Mat
951e884886fSBarry Smith 
952e884886fSBarry Smith    Input Parameters:
953e884886fSBarry Smith +  mat - the matrix free matrix created via MatCreateSNESMF()
954e884886fSBarry Smith -  period - 1 for everytime, 2 for every second etc
955e884886fSBarry Smith 
956e884886fSBarry Smith    Options Database Keys:
957e884886fSBarry Smith +  -mat_mffd_period <period>
958e884886fSBarry Smith 
959e884886fSBarry Smith    Level: advanced
960e884886fSBarry Smith 
961e884886fSBarry Smith 
962e884886fSBarry Smith .keywords: SNES, matrix-free, parameters
963e884886fSBarry Smith 
964e884886fSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(),
9651d0fab5eSBarry Smith           MatMFFDSetHHistory(), MatMFFDResetHHistory()
966e884886fSBarry Smith @*/
9677087cfbeSBarry Smith PetscErrorCode  MatMFFDSetPeriod(Mat mat,PetscInt period)
968e884886fSBarry Smith {
969bc0f08ceSBarry Smith   PetscErrorCode ierr;
970e884886fSBarry Smith 
971e884886fSBarry Smith   PetscFunctionBegin;
97288b4c220SStefano Zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
97388b4c220SStefano Zampini   PetscValidLogicalCollectiveInt(mat,period,2);
974bc0f08ceSBarry Smith   ierr = PetscTryMethod(mat,"MatMFFDSetPeriod_C",(Mat,PetscInt),(mat,period));CHKERRQ(ierr);
975e884886fSBarry Smith   PetscFunctionReturn(0);
976e884886fSBarry Smith }
977e884886fSBarry Smith 
978e884886fSBarry Smith /*@
979e884886fSBarry Smith    MatMFFDSetFunctionError - Sets the error_rel for the approximation of
980e884886fSBarry Smith    matrix-vector products using finite differences.
981e884886fSBarry Smith 
9823f9fe445SBarry Smith    Logically Collective on Mat
983e884886fSBarry Smith 
984e884886fSBarry Smith    Input Parameters:
985e884886fSBarry Smith +  mat - the matrix free matrix created via MatCreateMFFD() or MatCreateSNESMF()
986e884886fSBarry Smith -  error_rel - relative error (should be set to the square root of
987e884886fSBarry Smith                the relative error in the function evaluations)
988e884886fSBarry Smith 
989e884886fSBarry Smith    Options Database Keys:
990e884886fSBarry Smith +  -mat_mffd_err <error_rel> - Sets error_rel
991e884886fSBarry Smith 
992e884886fSBarry Smith    Level: advanced
993e884886fSBarry Smith 
994e884886fSBarry Smith    Notes:
995e884886fSBarry Smith    The default matrix-free matrix-vector product routine computes
996e884886fSBarry Smith .vb
997e884886fSBarry Smith      F'(u)*a = [F(u+h*a) - F(u)]/h where
998e884886fSBarry Smith      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
999e884886fSBarry Smith        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   else
1000e884886fSBarry Smith .ve
1001e884886fSBarry Smith 
1002e884886fSBarry Smith .keywords: SNES, matrix-free, parameters
1003e884886fSBarry Smith 
1004e884886fSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD
10051d0fab5eSBarry Smith           MatMFFDSetHHistory(), MatMFFDResetHHistory()
1006e884886fSBarry Smith @*/
10077087cfbeSBarry Smith PetscErrorCode  MatMFFDSetFunctionError(Mat mat,PetscReal error)
1008e884886fSBarry Smith {
1009bc0f08ceSBarry Smith   PetscErrorCode ierr;
1010e884886fSBarry Smith 
1011e884886fSBarry Smith   PetscFunctionBegin;
101288b4c220SStefano Zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
101388b4c220SStefano Zampini   PetscValidLogicalCollectiveReal(mat,error,2);
1014bc0f08ceSBarry Smith   ierr = PetscTryMethod(mat,"MatMFFDSetFunctionError_C",(Mat,PetscReal),(mat,error));CHKERRQ(ierr);
1015e884886fSBarry Smith   PetscFunctionReturn(0);
1016e884886fSBarry Smith }
1017e884886fSBarry Smith 
1018e884886fSBarry Smith /*@
1019e884886fSBarry Smith    MatMFFDSetHHistory - Sets an array to collect a history of the
1020e884886fSBarry Smith    differencing values (h) computed for the matrix-free product.
1021e884886fSBarry Smith 
10223f9fe445SBarry Smith    Logically Collective on Mat
1023e884886fSBarry Smith 
1024e884886fSBarry Smith    Input Parameters:
1025e884886fSBarry Smith +  J - the matrix-free matrix context
1026e884886fSBarry Smith .  histroy - space to hold the history
1027e884886fSBarry Smith -  nhistory - number of entries in history, if more entries are generated than
1028e884886fSBarry Smith               nhistory, then the later ones are discarded
1029e884886fSBarry Smith 
1030e884886fSBarry Smith    Level: advanced
1031e884886fSBarry Smith 
1032e884886fSBarry Smith    Notes:
1033e884886fSBarry Smith    Use MatMFFDResetHHistory() to reset the history counter and collect
1034e884886fSBarry Smith    a new batch of differencing parameters, h.
1035e884886fSBarry Smith 
1036e884886fSBarry Smith .keywords: SNES, matrix-free, h history, differencing history
1037e884886fSBarry Smith 
1038e884886fSBarry Smith .seealso: MatMFFDGetH(), MatCreateSNESMF(),
10391d0fab5eSBarry Smith           MatMFFDResetHHistory(), MatMFFDSetFunctionError()
1040e884886fSBarry Smith 
1041e884886fSBarry Smith @*/
10427087cfbeSBarry Smith PetscErrorCode  MatMFFDSetHHistory(Mat J,PetscScalar history[],PetscInt nhistory)
1043e884886fSBarry Smith {
1044e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)J->data;
1045bc0f08ceSBarry Smith   PetscErrorCode ierr;
1046bc0f08ceSBarry Smith   PetscBool      match;
1047e884886fSBarry Smith 
1048e884886fSBarry Smith   PetscFunctionBegin;
104988b4c220SStefano Zampini   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
105088b4c220SStefano Zampini   if (history) PetscValidPointer(history,2);
105188b4c220SStefano Zampini   PetscValidLogicalCollectiveInt(J,nhistory,3);
1052251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)J,MATMFFD,&match);CHKERRQ(ierr);
1053ce94432eSBarry Smith   if (!match) SETERRQ(PetscObjectComm((PetscObject)J),PETSC_ERR_ARG_WRONG,"Not a MFFD matrix");
1054e884886fSBarry Smith   ctx->historyh    = history;
1055e884886fSBarry Smith   ctx->maxcurrenth = nhistory;
105675567043SBarry Smith   ctx->currenth    = 0.;
1057e884886fSBarry Smith   PetscFunctionReturn(0);
1058e884886fSBarry Smith }
1059e884886fSBarry Smith 
1060e884886fSBarry Smith /*@
1061e884886fSBarry Smith    MatMFFDResetHHistory - Resets the counter to zero to begin
1062e884886fSBarry Smith    collecting a new set of differencing histories.
1063e884886fSBarry Smith 
10643f9fe445SBarry Smith    Logically Collective on Mat
1065e884886fSBarry Smith 
1066e884886fSBarry Smith    Input Parameters:
1067e884886fSBarry Smith .  J - the matrix-free matrix context
1068e884886fSBarry Smith 
1069e884886fSBarry Smith    Level: advanced
1070e884886fSBarry Smith 
1071e884886fSBarry Smith    Notes:
1072e884886fSBarry Smith    Use MatMFFDSetHHistory() to create the original history counter.
1073e884886fSBarry Smith 
1074e884886fSBarry Smith .keywords: SNES, matrix-free, h history, differencing history
1075e884886fSBarry Smith 
1076e884886fSBarry Smith .seealso: MatMFFDGetH(), MatCreateSNESMF(),
10771d0fab5eSBarry Smith           MatMFFDSetHHistory(), MatMFFDSetFunctionError()
1078e884886fSBarry Smith 
1079e884886fSBarry Smith @*/
10807087cfbeSBarry Smith PetscErrorCode  MatMFFDResetHHistory(Mat J)
1081e884886fSBarry Smith {
1082bc0f08ceSBarry Smith   PetscErrorCode ierr;
1083e884886fSBarry Smith 
1084e884886fSBarry Smith   PetscFunctionBegin;
108588b4c220SStefano Zampini   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
1086bc0f08ceSBarry Smith   ierr = PetscTryMethod(J,"MatMFFDResetHHistory_C",(Mat),(J));CHKERRQ(ierr);
1087e884886fSBarry Smith   PetscFunctionReturn(0);
1088e884886fSBarry Smith }
1089e884886fSBarry Smith 
1090e884886fSBarry Smith /*@
1091e884886fSBarry Smith     MatMFFDSetBase - Sets the vector U at which matrix vector products of the
1092e884886fSBarry Smith         Jacobian are computed
1093e884886fSBarry Smith 
10943f9fe445SBarry Smith     Logically Collective on Mat
1095e884886fSBarry Smith 
1096e884886fSBarry Smith     Input Parameters:
1097e884886fSBarry Smith +   J - the MatMFFD matrix
10983ec795f1SBarry Smith .   U - the vector
1099bcddec3dSBarry Smith -   F - (optional) vector that contains F(u) if it has been already computed
1100e884886fSBarry Smith 
1101e884886fSBarry Smith     Notes: This is rarely used directly
1102e884886fSBarry Smith 
11038af5ae88SBarry Smith     If F is provided then it is not recomputed. Otherwise the function is evaluated at the base
11048af5ae88SBarry Smith     point during the first MatMult() after each call to MatMFFDSetBase().
1105dff2f722SBarry Smith 
1106e884886fSBarry Smith     Level: advanced
1107e884886fSBarry Smith 
1108e884886fSBarry Smith @*/
11097087cfbeSBarry Smith PetscErrorCode  MatMFFDSetBase(Mat J,Vec U,Vec F)
1110e884886fSBarry Smith {
11114ac538c5SBarry Smith   PetscErrorCode ierr;
1112e884886fSBarry Smith 
1113e884886fSBarry Smith   PetscFunctionBegin;
11140700a824SBarry Smith   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
11150700a824SBarry Smith   PetscValidHeaderSpecific(U,VEC_CLASSID,2);
11160700a824SBarry Smith   if (F) PetscValidHeaderSpecific(F,VEC_CLASSID,3);
11174ac538c5SBarry Smith   ierr = PetscTryMethod(J,"MatMFFDSetBase_C",(Mat,Vec,Vec),(J,U,F));CHKERRQ(ierr);
1118e884886fSBarry Smith   PetscFunctionReturn(0);
1119e884886fSBarry Smith }
1120e884886fSBarry Smith 
1121e884886fSBarry Smith /*@C
1122e884886fSBarry Smith     MatMFFDSetCheckh - Sets a function that checks the computed h and adjusts
1123e884886fSBarry Smith         it to satisfy some criteria
1124e884886fSBarry Smith 
11253f9fe445SBarry Smith     Logically Collective on Mat
1126e884886fSBarry Smith 
1127e884886fSBarry Smith     Input Parameters:
1128e884886fSBarry Smith +   J - the MatMFFD matrix
1129e884886fSBarry Smith .   fun - the function that checks h
1130e884886fSBarry Smith -   ctx - any context needed by the function
1131e884886fSBarry Smith 
1132e884886fSBarry Smith     Options Database Keys:
1133e884886fSBarry Smith .   -mat_mffd_check_positivity
1134e884886fSBarry Smith 
1135e884886fSBarry Smith     Level: advanced
1136e884886fSBarry Smith 
1137755b7f64SBarry Smith     Notes: For example, MatMFFDCheckPositivity() insures that all entries
1138e884886fSBarry Smith        of U + h*a are non-negative
1139e884886fSBarry Smith 
1140755b7f64SBarry Smith      The function you provide is called after the default h has been computed and allows you to
1141755b7f64SBarry Smith      modify it.
1142755b7f64SBarry Smith 
1143755b7f64SBarry Smith .seealso:  MatMFFDCheckPositivity()
1144e884886fSBarry Smith @*/
11457087cfbeSBarry Smith PetscErrorCode  MatMFFDSetCheckh(Mat J,PetscErrorCode (*fun)(void*,Vec,Vec,PetscScalar*),void *ctx)
1146e884886fSBarry Smith {
11474ac538c5SBarry Smith   PetscErrorCode ierr;
1148e884886fSBarry Smith 
1149e884886fSBarry Smith   PetscFunctionBegin;
11500700a824SBarry Smith   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
11514ac538c5SBarry Smith   ierr = PetscTryMethod(J,"MatMFFDSetCheckh_C",(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*),(J,fun,ctx));CHKERRQ(ierr);
1152e884886fSBarry Smith   PetscFunctionReturn(0);
1153e884886fSBarry Smith }
1154e884886fSBarry Smith 
1155e884886fSBarry Smith /*@
1156e884886fSBarry Smith     MatMFFDCheckPositivity - Checks that all entries in U + h*a are positive or
1157e884886fSBarry Smith         zero, decreases h until this is satisfied.
1158e884886fSBarry Smith 
11593f9fe445SBarry Smith     Logically Collective on Vec
1160e884886fSBarry Smith 
1161e884886fSBarry Smith     Input Parameters:
1162e884886fSBarry Smith +   U - base vector that is added to
1163e884886fSBarry Smith .   a - vector that is added
1164e884886fSBarry Smith .   h - scaling factor on a
1165e884886fSBarry Smith -   dummy - context variable (unused)
1166e884886fSBarry Smith 
1167e884886fSBarry Smith     Options Database Keys:
1168e884886fSBarry Smith .   -mat_mffd_check_positivity
1169e884886fSBarry Smith 
1170e884886fSBarry Smith     Level: advanced
1171e884886fSBarry Smith 
1172e884886fSBarry Smith     Notes: This is rarely used directly, rather it is passed as an argument to
1173e884886fSBarry Smith            MatMFFDSetCheckh()
1174e884886fSBarry Smith 
1175e884886fSBarry Smith .seealso:  MatMFFDSetCheckh()
1176e884886fSBarry Smith @*/
11777087cfbeSBarry Smith PetscErrorCode  MatMFFDCheckPositivity(void *dummy,Vec U,Vec a,PetscScalar *h)
1178e884886fSBarry Smith {
1179e884886fSBarry Smith   PetscReal      val, minval;
1180e884886fSBarry Smith   PetscScalar    *u_vec, *a_vec;
1181e884886fSBarry Smith   PetscErrorCode ierr;
1182e884886fSBarry Smith   PetscInt       i,n;
1183e884886fSBarry Smith   MPI_Comm       comm;
1184e884886fSBarry Smith 
1185e884886fSBarry Smith   PetscFunctionBegin;
118688b4c220SStefano Zampini   PetscValidHeaderSpecific(U,VEC_CLASSID,2);
118788b4c220SStefano Zampini   PetscValidHeaderSpecific(a,VEC_CLASSID,3);
118888b4c220SStefano Zampini   PetscValidPointer(h,4);
1189e884886fSBarry Smith   ierr   = PetscObjectGetComm((PetscObject)U,&comm);CHKERRQ(ierr);
1190e884886fSBarry Smith   ierr   = VecGetArray(U,&u_vec);CHKERRQ(ierr);
1191e884886fSBarry Smith   ierr   = VecGetArray(a,&a_vec);CHKERRQ(ierr);
1192e884886fSBarry Smith   ierr   = VecGetLocalSize(U,&n);CHKERRQ(ierr);
1193e884886fSBarry Smith   minval = PetscAbsScalar(*h*1.01);
1194e884886fSBarry Smith   for (i=0; i<n; i++) {
1195e884886fSBarry Smith     if (PetscRealPart(u_vec[i] + *h*a_vec[i]) <= 0.0) {
1196e884886fSBarry Smith       val = PetscAbsScalar(u_vec[i]/a_vec[i]);
1197e884886fSBarry Smith       if (val < minval) minval = val;
1198e884886fSBarry Smith     }
1199e884886fSBarry Smith   }
1200e884886fSBarry Smith   ierr = VecRestoreArray(U,&u_vec);CHKERRQ(ierr);
1201e884886fSBarry Smith   ierr = VecRestoreArray(a,&a_vec);CHKERRQ(ierr);
1202b2566f29SBarry Smith   ierr = MPIU_Allreduce(&minval,&val,1,MPIU_REAL,MPIU_MIN,comm);CHKERRQ(ierr);
1203e884886fSBarry Smith   if (val <= PetscAbsScalar(*h)) {
12046712e2f1SBarry Smith     ierr = PetscInfo2(U,"Scaling back h from %g to %g\n",(double)PetscRealPart(*h),(double)(.99*val));CHKERRQ(ierr);
1205e884886fSBarry Smith     if (PetscRealPart(*h) > 0.0) *h =  0.99*val;
1206e884886fSBarry Smith     else                         *h = -0.99*val;
1207e884886fSBarry Smith   }
1208e884886fSBarry Smith   PetscFunctionReturn(0);
1209e884886fSBarry Smith }
1210