xref: /petsc/src/mat/impls/mffd/mffd.c (revision 487a658c8b32ba712a1dc8280daad2fd70c1dcd9)
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];
458e81d068SLisandro Dalcin   PetscBool      opt,pkg;
463ec795f1SBarry Smith   PetscErrorCode ierr;
473ec795f1SBarry Smith 
483ec795f1SBarry Smith   PetscFunctionBegin;
49b022a5c1SBarry Smith   if (MatMFFDPackageInitialized) PetscFunctionReturn(0);
50b022a5c1SBarry Smith   MatMFFDPackageInitialized = PETSC_TRUE;
513ec795f1SBarry Smith   /* Register Classes */
520700a824SBarry Smith   ierr = PetscClassIdRegister("MatMFFD",&MATMFFD_CLASSID);CHKERRQ(ierr);
533ec795f1SBarry Smith   /* Register Constructors */
54607a6623SBarry Smith   ierr = MatMFFDRegisterAll();CHKERRQ(ierr);
553ec795f1SBarry Smith   /* Register Events */
560700a824SBarry Smith   ierr = PetscLogEventRegister("MatMult MF",MATMFFD_CLASSID,&MATMFFD_Mult);CHKERRQ(ierr);
573ec795f1SBarry Smith   /* Process info exclusions */
588e81d068SLisandro Dalcin   ierr = PetscOptionsGetString(NULL,NULL,"-info_exclude",logList,sizeof(logList),&opt);CHKERRQ(ierr);
593ec795f1SBarry Smith   if (opt) {
608e81d068SLisandro Dalcin     ierr = PetscStrInList("matmffd",logList,',',&pkg);CHKERRQ(ierr);
618e81d068SLisandro Dalcin     if (pkg) {ierr = PetscInfoDeactivateClass(MATMFFD_CLASSID);CHKERRQ(ierr);}
623ec795f1SBarry Smith   }
633ec795f1SBarry Smith   /* Process summary exclusions */
648e81d068SLisandro Dalcin   ierr = PetscOptionsGetString(NULL,NULL,"-log_exclude",logList,sizeof(logList),&opt);CHKERRQ(ierr);
653ec795f1SBarry Smith   if (opt) {
668e81d068SLisandro Dalcin     ierr = PetscStrInList("matmffd",logList,',',&pkg);CHKERRQ(ierr);
67fa2bb9feSLisandro Dalcin     if (pkg) {ierr = PetscLogEventExcludeClass(MATMFFD_CLASSID);CHKERRQ(ierr);}
683ec795f1SBarry Smith   }
698e81d068SLisandro Dalcin   /* Register package finalizer */
70b022a5c1SBarry Smith   ierr = PetscRegisterFinalize(MatMFFDFinalizePackage);CHKERRQ(ierr);
713ec795f1SBarry Smith   PetscFunctionReturn(0);
723ec795f1SBarry Smith }
733ec795f1SBarry Smith 
74e884886fSBarry Smith /*@C
75e884886fSBarry Smith     MatMFFDSetType - Sets the method that is used to compute the
76e884886fSBarry Smith     differencing parameter for finite differene matrix-free formulations.
77e884886fSBarry Smith 
78e884886fSBarry Smith     Input Parameters:
79e884886fSBarry Smith +   mat - the "matrix-free" matrix created via MatCreateSNESMF(), or MatCreateMFFD()
80e884886fSBarry Smith           or MatSetType(mat,MATMFFD);
81e884886fSBarry Smith -   ftype - the type requested, either MATMFFD_WP or MATMFFD_DS
82e884886fSBarry Smith 
83e884886fSBarry Smith     Level: advanced
84e884886fSBarry Smith 
85e884886fSBarry Smith     Notes:
86e884886fSBarry Smith     For example, such routines can compute h for use in
87e884886fSBarry Smith     Jacobian-vector products of the form
88e884886fSBarry Smith 
89e884886fSBarry Smith                         F(x+ha) - F(x)
90e884886fSBarry Smith           F'(u)a  ~=  ----------------
91e884886fSBarry Smith                               h
92e884886fSBarry Smith 
93755b7f64SBarry Smith .seealso: MatCreateSNESMF(), MatMFFDRegister(), MatMFFDSetFunction(), MatCreateMFFD()
94e884886fSBarry Smith @*/
9519fd82e9SBarry Smith PetscErrorCode  MatMFFDSetType(Mat mat,MatMFFDType ftype)
96e884886fSBarry Smith {
97e884886fSBarry Smith   PetscErrorCode ierr,(*r)(MatMFFD);
98e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)mat->data;
99ace3abfcSBarry Smith   PetscBool      match;
100e884886fSBarry Smith 
101e884886fSBarry Smith   PetscFunctionBegin;
1020700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
103e884886fSBarry Smith   PetscValidCharPointer(ftype,2);
104e884886fSBarry Smith 
105251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)mat,MATMFFD,&match);CHKERRQ(ierr);
1069a13eae7SJed Brown   if (!match) PetscFunctionReturn(0);
1079a13eae7SJed Brown 
108e884886fSBarry Smith   /* already set, so just return */
109251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)ctx,ftype,&match);CHKERRQ(ierr);
110e884886fSBarry Smith   if (match) PetscFunctionReturn(0);
111e884886fSBarry Smith 
112e884886fSBarry Smith   /* destroy the old one if it exists */
113e884886fSBarry Smith   if (ctx->ops->destroy) {
114e884886fSBarry Smith     ierr = (*ctx->ops->destroy)(ctx);CHKERRQ(ierr);
115e884886fSBarry Smith   }
116e884886fSBarry Smith 
1171c9cd337SJed Brown   ierr =  PetscFunctionListFind(MatMFFDList,ftype,&r);CHKERRQ(ierr);
118e32f2f54SBarry Smith   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown MatMFFD type %s given",ftype);
119e884886fSBarry Smith   ierr = (*r)(ctx);CHKERRQ(ierr);
120e884886fSBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)ctx,ftype);CHKERRQ(ierr);
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 {
129e884886fSBarry Smith   MatMFFD ctx = (MatMFFD)mat->data;
130e884886fSBarry Smith 
131e884886fSBarry Smith   PetscFunctionBegin;
132e884886fSBarry Smith   ctx->funcisetbase = func;
133486afcceSStefano Zampini   /* allow users to compose their own getdiagonal and allow MatHasOperation
134486afcceSStefano Zampini      to return false if the two functions pointers are not set */
135486afcceSStefano Zampini   if (!mat->ops->getdiagonal && func) {
136dbfa1270SStefano Zampini     mat->ops->getdiagonal = MatGetDiagonal_MFFD;
1375d21e1f8SStefano Zampini   }
138e884886fSBarry Smith   PetscFunctionReturn(0);
139e884886fSBarry Smith }
140e884886fSBarry Smith 
141e884886fSBarry Smith typedef PetscErrorCode (*FCN2)(void*,PetscInt,Vec,PetscScalar*); /* force argument to next function to not be extern C*/
14240244768SBarry Smith static PetscErrorCode  MatMFFDSetFunctioni_MFFD(Mat mat,FCN2 funci)
143e884886fSBarry Smith {
144e884886fSBarry Smith   MatMFFD ctx = (MatMFFD)mat->data;
145e884886fSBarry Smith 
146e884886fSBarry Smith   PetscFunctionBegin;
147e884886fSBarry Smith   ctx->funci = funci;
148486afcceSStefano Zampini   /* allow users to compose their own getdiagonal and allow MatHasOperation
149486afcceSStefano Zampini      to return false if the two functions pointers are not set */
150486afcceSStefano Zampini   if (!mat->ops->getdiagonal && funci) {
151dbfa1270SStefano Zampini     mat->ops->getdiagonal = MatGetDiagonal_MFFD;
1525d21e1f8SStefano Zampini   }
153e884886fSBarry Smith   PetscFunctionReturn(0);
154e884886fSBarry Smith }
155e884886fSBarry Smith 
15640244768SBarry Smith static PetscErrorCode  MatMFFDResetHHistory_MFFD(Mat J)
157bc0f08ceSBarry Smith {
158bc0f08ceSBarry Smith   MatMFFD ctx = (MatMFFD)J->data;
159bc0f08ceSBarry Smith 
160bc0f08ceSBarry Smith   PetscFunctionBegin;
161bc0f08ceSBarry Smith   ctx->ncurrenth = 0;
162bc0f08ceSBarry Smith   PetscFunctionReturn(0);
163bc0f08ceSBarry Smith }
164e884886fSBarry Smith 
1651c84c290SBarry Smith /*@C
1661c84c290SBarry Smith    MatMFFDRegister - Adds a method to the MatMFFD registry.
1671c84c290SBarry Smith 
1681c84c290SBarry Smith    Not Collective
1691c84c290SBarry Smith 
1701c84c290SBarry Smith    Input Parameters:
1711c84c290SBarry Smith +  name_solver - name of a new user-defined compute-h module
1721c84c290SBarry Smith -  routine_create - routine to create method context
1731c84c290SBarry Smith 
1741c84c290SBarry Smith    Level: developer
1751c84c290SBarry Smith 
1761c84c290SBarry Smith    Notes:
1771c84c290SBarry Smith    MatMFFDRegister() may be called multiple times to add several user-defined solvers.
1781c84c290SBarry Smith 
1791c84c290SBarry Smith    Sample usage:
1801c84c290SBarry Smith .vb
181bdf89e91SBarry Smith    MatMFFDRegister("my_h",MyHCreate);
1821c84c290SBarry Smith .ve
1831c84c290SBarry Smith 
1841c84c290SBarry Smith    Then, your solver can be chosen with the procedural interface via
1851c84c290SBarry Smith $     MatMFFDSetType(mfctx,"my_h")
1861c84c290SBarry Smith    or at runtime via the option
187be7a6d03SBarry Smith $     -mat_mffd_type my_h
1881c84c290SBarry Smith 
1891c84c290SBarry Smith .keywords: MatMFFD, register
1901c84c290SBarry Smith 
1911c84c290SBarry Smith .seealso: MatMFFDRegisterAll(), MatMFFDRegisterDestroy()
1921c84c290SBarry Smith  @*/
193bdf89e91SBarry Smith PetscErrorCode  MatMFFDRegister(const char sname[],PetscErrorCode (*function)(MatMFFD))
194e884886fSBarry Smith {
195e884886fSBarry Smith   PetscErrorCode ierr;
196e884886fSBarry Smith 
197e884886fSBarry Smith   PetscFunctionBegin;
198a240a19fSJed Brown   ierr = PetscFunctionListAdd(&MatMFFDList,sname,function);CHKERRQ(ierr);
199e884886fSBarry Smith   PetscFunctionReturn(0);
200e884886fSBarry Smith }
201e884886fSBarry Smith 
202e884886fSBarry Smith /* ----------------------------------------------------------------------------------------*/
20340244768SBarry Smith static PetscErrorCode MatDestroy_MFFD(Mat mat)
204e884886fSBarry Smith {
205e884886fSBarry Smith   PetscErrorCode ierr;
206e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)mat->data;
207e884886fSBarry Smith 
208e884886fSBarry Smith   PetscFunctionBegin;
2096bf464f9SBarry Smith   ierr = VecDestroy(&ctx->w);CHKERRQ(ierr);
2106bf464f9SBarry Smith   ierr = VecDestroy(&ctx->drscale);CHKERRQ(ierr);
2116bf464f9SBarry Smith   ierr = VecDestroy(&ctx->dlscale);CHKERRQ(ierr);
2126bf464f9SBarry Smith   ierr = VecDestroy(&ctx->dshift);CHKERRQ(ierr);
213c51ad4d4SStefano Zampini   ierr = VecDestroy(&ctx->dshiftw);CHKERRQ(ierr);
214c51ad4d4SStefano Zampini   ierr = VecDestroy(&ctx->current_u);CHKERRQ(ierr);
215cfe22f5eSBarry Smith   if (ctx->current_f_allocated) {
2166bf464f9SBarry Smith     ierr = VecDestroy(&ctx->current_f);CHKERRQ(ierr);
217cfe22f5eSBarry Smith   }
218e884886fSBarry Smith   if (ctx->ops->destroy) {ierr = (*ctx->ops->destroy)(ctx);CHKERRQ(ierr);}
2196bf464f9SBarry Smith   ierr      = PetscHeaderDestroy(&ctx);CHKERRQ(ierr);
220bf0cc555SLisandro Dalcin   mat->data = 0;
221e884886fSBarry Smith 
222bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetBase_C",NULL);CHKERRQ(ierr);
223bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctioniBase_C",NULL);CHKERRQ(ierr);
224bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctioni_C",NULL);CHKERRQ(ierr);
225bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunction_C",NULL);CHKERRQ(ierr);
226bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctionError_C",NULL);CHKERRQ(ierr);
227bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetCheckh_C",NULL);CHKERRQ(ierr);
228bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetPeriod_C",NULL);CHKERRQ(ierr);
229bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDResetHHistory_C",NULL);CHKERRQ(ierr);
230e884886fSBarry Smith   PetscFunctionReturn(0);
231e884886fSBarry Smith }
232e884886fSBarry Smith 
233e884886fSBarry Smith /*
234e884886fSBarry Smith    MatMFFDView_MFFD - Views matrix-free parameters.
235e884886fSBarry Smith 
236e884886fSBarry Smith */
23740244768SBarry Smith static PetscErrorCode MatView_MFFD(Mat J,PetscViewer viewer)
238e884886fSBarry Smith {
239e884886fSBarry Smith   PetscErrorCode ierr;
240e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)J->data;
241a283635eSDmitry Karpeev   PetscBool      iascii, viewbase, viewfunction;
242a283635eSDmitry Karpeev   const char     *prefix;
243e884886fSBarry Smith 
244e884886fSBarry Smith   PetscFunctionBegin;
245251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
246e884886fSBarry Smith   if (iascii) {
247a283635eSDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"Matrix-free approximation:\n");CHKERRQ(ierr);
248a283635eSDmitry Karpeev     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
24957622a8eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"err=%g (relative error in function evaluation)\n",(double)ctx->error_rel);CHKERRQ(ierr);
2507adad957SLisandro Dalcin     if (!((PetscObject)ctx)->type_name) {
251e884886fSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"The compute h routine has not yet been set\n");CHKERRQ(ierr);
252e884886fSBarry Smith     } else {
2537adad957SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"Using %s compute h routine\n",((PetscObject)ctx)->type_name);CHKERRQ(ierr);
254e884886fSBarry Smith     }
255e884886fSBarry Smith     if (ctx->ops->view) {
256e884886fSBarry Smith       ierr = (*ctx->ops->view)(ctx,viewer);CHKERRQ(ierr);
257e884886fSBarry Smith     }
258a283635eSDmitry Karpeev     ierr = PetscObjectGetOptionsPrefix((PetscObject)J, &prefix);CHKERRQ(ierr);
259a283635eSDmitry Karpeev 
260c5929fdfSBarry Smith     ierr = PetscOptionsHasName(((PetscObject)J)->options,prefix, "-mat_mffd_view_base", &viewbase);CHKERRQ(ierr);
261a283635eSDmitry Karpeev     if (viewbase) {
262a283635eSDmitry Karpeev       ierr = PetscViewerASCIIPrintf(viewer, "Base:\n");CHKERRQ(ierr);
263a283635eSDmitry Karpeev       ierr = VecView(ctx->current_u, viewer);CHKERRQ(ierr);
264a283635eSDmitry Karpeev     }
265c5929fdfSBarry Smith     ierr = PetscOptionsHasName(((PetscObject)J)->options,prefix, "-mat_mffd_view_function", &viewfunction);CHKERRQ(ierr);
266a283635eSDmitry Karpeev     if (viewfunction) {
267a283635eSDmitry Karpeev       ierr = PetscViewerASCIIPrintf(viewer, "Function:\n");CHKERRQ(ierr);
268a283635eSDmitry Karpeev       ierr = VecView(ctx->current_f, viewer);CHKERRQ(ierr);
269a283635eSDmitry Karpeev     }
270a283635eSDmitry Karpeev     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
27111aeaf0aSBarry Smith   }
272e884886fSBarry Smith   PetscFunctionReturn(0);
273e884886fSBarry Smith }
274e884886fSBarry Smith 
275e884886fSBarry Smith /*
276e884886fSBarry Smith    MatAssemblyEnd_MFFD - Resets the ctx->ncurrenth to zero. This
277e884886fSBarry Smith    allows the user to indicate the beginning of a new linear solve by calling
278e884886fSBarry Smith    MatAssemblyXXX() on the matrix free matrix. This then allows the
2791d0fab5eSBarry Smith    MatCreateMFFD_WP() to properly compute ||U|| only the first time
280e884886fSBarry Smith    in the linear solver rather than every time.
2815a576424SJed Brown 
282cc2e6a90SBarry Smith    This function is referenced directly from MatAssemblyEnd_SNESMF(), which may be in a different shared library hence
283cc2e6a90SBarry Smith    it must be labeled as PETSC_EXTERN
284e884886fSBarry Smith */
2855a576424SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyEnd_MFFD(Mat J,MatAssemblyType mt)
286e884886fSBarry Smith {
287e884886fSBarry Smith   PetscErrorCode ierr;
288e884886fSBarry Smith   MatMFFD        j = (MatMFFD)J->data;
289e884886fSBarry Smith 
290e884886fSBarry Smith   PetscFunctionBegin;
291e884886fSBarry Smith   ierr      = MatMFFDResetHHistory(J);CHKERRQ(ierr);
292e884886fSBarry Smith   j->vshift = 0.0;
293e884886fSBarry Smith   j->vscale = 1.0;
294e884886fSBarry Smith   PetscFunctionReturn(0);
295e884886fSBarry Smith }
296e884886fSBarry Smith 
297e884886fSBarry Smith /*
298e884886fSBarry Smith   MatMult_MFFD - Default matrix-free form for Jacobian-vector product, y = F'(u)*a:
299e884886fSBarry Smith 
300e884886fSBarry Smith         y ~= (F(u + ha) - F(u))/h,
301e884886fSBarry Smith   where F = nonlinear function, as set by SNESSetFunction()
302e884886fSBarry Smith         u = current iterate
303e884886fSBarry Smith         h = difference interval
304e884886fSBarry Smith */
30540244768SBarry Smith static PetscErrorCode MatMult_MFFD(Mat mat,Vec a,Vec y)
306e884886fSBarry Smith {
307e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)mat->data;
308e884886fSBarry Smith   PetscScalar    h;
309e884886fSBarry Smith   Vec            w,U,F;
310e884886fSBarry Smith   PetscErrorCode ierr;
311ace3abfcSBarry Smith   PetscBool      zeroa;
312e884886fSBarry Smith 
313e884886fSBarry Smith   PetscFunctionBegin;
314ce94432eSBarry 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");
315e884886fSBarry Smith   /* We log matrix-free matrix-vector products separately, so that we can
316e884886fSBarry Smith      separate the performance monitoring from the cases that use conventional
317e884886fSBarry Smith      storage.  We may eventually modify event logging to associate events
318e884886fSBarry Smith      with particular objects, hence alleviating the more general problem. */
319e884886fSBarry Smith   ierr = PetscLogEventBegin(MATMFFD_Mult,a,y,0,0);CHKERRQ(ierr);
320e884886fSBarry Smith 
321e884886fSBarry Smith   w = ctx->w;
322e884886fSBarry Smith   U = ctx->current_u;
3233ec795f1SBarry Smith   F = ctx->current_f;
324e884886fSBarry Smith   /*
325e884886fSBarry Smith       Compute differencing parameter
326e884886fSBarry Smith   */
3274a2f8832SBarry Smith   if (!((PetscObject)ctx)->type_name) {
328e884886fSBarry Smith     ierr = MatMFFDSetType(mat,MATMFFD_WP);CHKERRQ(ierr);
3299c6ac3b3SBarry Smith     ierr = MatSetFromOptions(mat);CHKERRQ(ierr);
330e884886fSBarry Smith   }
331e884886fSBarry Smith   ierr = (*ctx->ops->compute)(ctx,U,a,&h,&zeroa);CHKERRQ(ierr);
332e884886fSBarry Smith   if (zeroa) {
333e884886fSBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
334e884886fSBarry Smith     PetscFunctionReturn(0);
335e884886fSBarry Smith   }
336e884886fSBarry Smith 
33784d44b13SHong Zhang   if (mat->erroriffailure && PetscIsInfOrNanScalar(h)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Computed Nan differencing parameter h");
338e884886fSBarry Smith   if (ctx->checkh) {
339e884886fSBarry Smith     ierr = (*ctx->checkh)(ctx->checkhctx,U,a,&h);CHKERRQ(ierr);
340e884886fSBarry Smith   }
341e884886fSBarry Smith 
342e884886fSBarry Smith   /* keep a record of the current differencing parameter h */
343e884886fSBarry Smith   ctx->currenth = h;
344e884886fSBarry Smith #if defined(PETSC_USE_COMPLEX)
34557622a8eSBarry Smith   ierr = PetscInfo2(mat,"Current differencing parameter: %g + %g i\n",(double)PetscRealPart(h),(double)PetscImaginaryPart(h));CHKERRQ(ierr);
346e884886fSBarry Smith #else
347e884886fSBarry Smith   ierr = PetscInfo1(mat,"Current differencing parameter: %15.12e\n",h);CHKERRQ(ierr);
348e884886fSBarry Smith #endif
349e884886fSBarry Smith   if (ctx->historyh && ctx->ncurrenth < ctx->maxcurrenth) {
350e884886fSBarry Smith     ctx->historyh[ctx->ncurrenth] = h;
351e884886fSBarry Smith   }
352e884886fSBarry Smith   ctx->ncurrenth++;
353e884886fSBarry Smith 
354e884886fSBarry Smith   /* w = u + ha */
355c3bb7e23SBarry Smith   if (ctx->drscale) {
356c3bb7e23SBarry Smith     ierr = VecPointwiseMult(ctx->drscale,a,U);CHKERRQ(ierr);
357c3bb7e23SBarry Smith     ierr = VecAYPX(U,h,w);CHKERRQ(ierr);
358c3bb7e23SBarry Smith   } else {
359e884886fSBarry Smith     ierr = VecWAXPY(w,h,a,U);CHKERRQ(ierr);
360c3bb7e23SBarry Smith   }
361e884886fSBarry Smith 
3621b797266SDmitry Karpeev   /* compute func(U) as base for differencing; only needed first time in and not when provided by user */
3631b797266SDmitry Karpeev   if (ctx->ncurrenth == 1 && ctx->current_f_allocated) {
364e884886fSBarry Smith     ierr = (*ctx->func)(ctx->funcctx,U,F);CHKERRQ(ierr);
36552121784SBarry Smith   }
366e884886fSBarry Smith   ierr = (*ctx->func)(ctx->funcctx,w,y);CHKERRQ(ierr);
367e884886fSBarry Smith 
368e884886fSBarry Smith   ierr = VecAXPY(y,-1.0,F);CHKERRQ(ierr);
369e884886fSBarry Smith   ierr = VecScale(y,1.0/h);CHKERRQ(ierr);
370e884886fSBarry Smith 
371c35ec66cSMatthew G Knepley   if ((ctx->vshift != 0.0) || (ctx->vscale != 1.0)) {
372e884886fSBarry Smith     ierr = VecAXPBY(y,ctx->vshift,ctx->vscale,a);CHKERRQ(ierr);
373c3bb7e23SBarry Smith   }
374c3bb7e23SBarry Smith   if (ctx->dlscale) {
375c3bb7e23SBarry Smith     ierr = VecPointwiseMult(y,ctx->dlscale,y);CHKERRQ(ierr);
376c3bb7e23SBarry Smith   }
377c3bb7e23SBarry Smith   if (ctx->dshift) {
378c51ad4d4SStefano Zampini     if (!ctx->dshiftw) {
379c51ad4d4SStefano Zampini       ierr = VecDuplicate(y,&ctx->dshiftw);CHKERRQ(ierr);
380c51ad4d4SStefano Zampini     }
381c51ad4d4SStefano Zampini     ierr = VecPointwiseMult(ctx->dshift,a,ctx->dshiftw);CHKERRQ(ierr);
382c51ad4d4SStefano Zampini     ierr = VecAXPY(y,1.0,ctx->dshiftw);CHKERRQ(ierr);
383c3bb7e23SBarry Smith   }
384e884886fSBarry Smith 
38539601f4eSBarry Smith   if (mat->nullsp) {ierr = MatNullSpaceRemove(mat->nullsp,y);CHKERRQ(ierr);}
386e884886fSBarry Smith 
387e884886fSBarry Smith   ierr = PetscLogEventEnd(MATMFFD_Mult,a,y,0,0);CHKERRQ(ierr);
388e884886fSBarry Smith   PetscFunctionReturn(0);
389e884886fSBarry Smith }
390e884886fSBarry Smith 
391e884886fSBarry Smith /*
392e884886fSBarry Smith   MatGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix
393e884886fSBarry Smith 
394e884886fSBarry Smith         y ~= (F(u + ha) - F(u))/h,
395e884886fSBarry Smith   where F = nonlinear function, as set by SNESSetFunction()
396e884886fSBarry Smith         u = current iterate
397e884886fSBarry Smith         h = difference interval
398e884886fSBarry Smith */
3995d21e1f8SStefano Zampini PetscErrorCode MatGetDiagonal_MFFD(Mat mat,Vec a)
400e884886fSBarry Smith {
401e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)mat->data;
402e884886fSBarry Smith   PetscScalar    h,*aa,*ww,v;
403e884886fSBarry Smith   PetscReal      epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON;
404e884886fSBarry Smith   Vec            w,U;
405e884886fSBarry Smith   PetscErrorCode ierr;
406e884886fSBarry Smith   PetscInt       i,rstart,rend;
407e884886fSBarry Smith 
408e884886fSBarry Smith   PetscFunctionBegin;
409486afcceSStefano Zampini   if (!ctx->func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Requires calling MatMFFDSetFunction() first");
410486afcceSStefano Zampini   if (!ctx->funci) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Requires calling MatMFFDSetFunctioni() first");
411486afcceSStefano Zampini   if (!ctx->funcisetbase) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Requires calling MatMFFDSetFunctioniBase() first");
412e884886fSBarry Smith   w    = ctx->w;
413e884886fSBarry Smith   U    = ctx->current_u;
414e884886fSBarry Smith   ierr = (*ctx->func)(ctx->funcctx,U,a);CHKERRQ(ierr);
415e884886fSBarry Smith   ierr = (*ctx->funcisetbase)(ctx->funcctx,U);CHKERRQ(ierr);
416e884886fSBarry Smith   ierr = VecCopy(U,w);CHKERRQ(ierr);
417e884886fSBarry Smith 
418e884886fSBarry Smith   ierr = VecGetOwnershipRange(a,&rstart,&rend);CHKERRQ(ierr);
419e884886fSBarry Smith   ierr = VecGetArray(a,&aa);CHKERRQ(ierr);
420e884886fSBarry Smith   for (i=rstart; i<rend; i++) {
421e884886fSBarry Smith     ierr = VecGetArray(w,&ww);CHKERRQ(ierr);
422e884886fSBarry Smith     h    = ww[i-rstart];
423e884886fSBarry Smith     if (h == 0.0) h = 1.0;
424e884886fSBarry Smith     if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0)     h = umin;
425e884886fSBarry Smith     else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin;
426e884886fSBarry Smith     h *= epsilon;
427e884886fSBarry Smith 
428e884886fSBarry Smith     ww[i-rstart] += h;
429e884886fSBarry Smith     ierr          = VecRestoreArray(w,&ww);CHKERRQ(ierr);
430e884886fSBarry Smith     ierr          = (*ctx->funci)(ctx->funcctx,i,w,&v);CHKERRQ(ierr);
431e884886fSBarry Smith     aa[i-rstart]  = (v - aa[i-rstart])/h;
432e884886fSBarry Smith 
433e884886fSBarry Smith     /* possibly shift and scale result */
434c35ec66cSMatthew G Knepley     if ((ctx->vshift != 0.0) || (ctx->vscale != 1.0)) {
435e884886fSBarry Smith       aa[i - rstart] = ctx->vshift + ctx->vscale*aa[i-rstart];
436c3bb7e23SBarry Smith     }
437e884886fSBarry Smith 
438e884886fSBarry Smith     ierr          = VecGetArray(w,&ww);CHKERRQ(ierr);
439e884886fSBarry Smith     ww[i-rstart] -= h;
440e884886fSBarry Smith     ierr          = VecRestoreArray(w,&ww);CHKERRQ(ierr);
441e884886fSBarry Smith   }
442e884886fSBarry Smith   ierr = VecRestoreArray(a,&aa);CHKERRQ(ierr);
443e884886fSBarry Smith   PetscFunctionReturn(0);
444e884886fSBarry Smith }
445e884886fSBarry Smith 
44640244768SBarry Smith static PetscErrorCode MatDiagonalScale_MFFD(Mat mat,Vec ll,Vec rr)
447c3bb7e23SBarry Smith {
448c3bb7e23SBarry Smith   MatMFFD        aij = (MatMFFD)mat->data;
449c3bb7e23SBarry Smith   PetscErrorCode ierr;
450c3bb7e23SBarry Smith 
451c3bb7e23SBarry Smith   PetscFunctionBegin;
452c3bb7e23SBarry Smith   if (ll && !aij->dlscale) {
453c3bb7e23SBarry Smith     ierr = VecDuplicate(ll,&aij->dlscale);CHKERRQ(ierr);
454c3bb7e23SBarry Smith   }
455c3bb7e23SBarry Smith   if (rr && !aij->drscale) {
456c3bb7e23SBarry Smith     ierr = VecDuplicate(rr,&aij->drscale);CHKERRQ(ierr);
457c3bb7e23SBarry Smith   }
458c3bb7e23SBarry Smith   if (ll) {
459c3bb7e23SBarry Smith     ierr = VecCopy(ll,aij->dlscale);CHKERRQ(ierr);
460c3bb7e23SBarry Smith   }
461c3bb7e23SBarry Smith   if (rr) {
462c3bb7e23SBarry Smith     ierr = VecCopy(rr,aij->drscale);CHKERRQ(ierr);
463c3bb7e23SBarry Smith   }
464c3bb7e23SBarry Smith   PetscFunctionReturn(0);
465c3bb7e23SBarry Smith }
466c3bb7e23SBarry Smith 
46740244768SBarry Smith static PetscErrorCode MatDiagonalSet_MFFD(Mat mat,Vec ll,InsertMode mode)
468c3bb7e23SBarry Smith {
469c3bb7e23SBarry Smith   MatMFFD        aij = (MatMFFD)mat->data;
470c3bb7e23SBarry Smith   PetscErrorCode ierr;
471c3bb7e23SBarry Smith 
472c3bb7e23SBarry Smith   PetscFunctionBegin;
473ce94432eSBarry Smith   if (mode == INSERT_VALUES) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"No diagonal set with INSERT_VALUES");
474c3bb7e23SBarry Smith   if (!aij->dshift) {
475c3bb7e23SBarry Smith     ierr = VecDuplicate(ll,&aij->dshift);CHKERRQ(ierr);
476c3bb7e23SBarry Smith   }
477c3bb7e23SBarry Smith   ierr = VecAXPY(aij->dshift,1.0,ll);CHKERRQ(ierr);
478c3bb7e23SBarry Smith   PetscFunctionReturn(0);
479c3bb7e23SBarry Smith }
480c3bb7e23SBarry Smith 
48140244768SBarry Smith static PetscErrorCode MatShift_MFFD(Mat Y,PetscScalar a)
482e884886fSBarry Smith {
483e884886fSBarry Smith   MatMFFD shell = (MatMFFD)Y->data;
4845fd66863SKarl Rupp 
485e884886fSBarry Smith   PetscFunctionBegin;
486e884886fSBarry Smith   shell->vshift += a;
487e884886fSBarry Smith   PetscFunctionReturn(0);
488e884886fSBarry Smith }
489e884886fSBarry Smith 
49040244768SBarry Smith static PetscErrorCode MatScale_MFFD(Mat Y,PetscScalar a)
491e884886fSBarry Smith {
492e884886fSBarry Smith   MatMFFD shell = (MatMFFD)Y->data;
4935fd66863SKarl Rupp 
494e884886fSBarry Smith   PetscFunctionBegin;
495e884886fSBarry Smith   shell->vscale *= a;
496e884886fSBarry Smith   PetscFunctionReturn(0);
497e884886fSBarry Smith }
498e884886fSBarry Smith 
499d8fa6a54SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDSetBase_MFFD(Mat J,Vec U,Vec F)
500e884886fSBarry Smith {
501e884886fSBarry Smith   PetscErrorCode ierr;
502e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)J->data;
503e884886fSBarry Smith 
504e884886fSBarry Smith   PetscFunctionBegin;
505e884886fSBarry Smith   ierr = MatMFFDResetHHistory(J);CHKERRQ(ierr);
506c51ad4d4SStefano Zampini   if (!ctx->current_u) {
507c51ad4d4SStefano Zampini     ierr = VecDuplicate(U,&ctx->current_u);CHKERRQ(ierr);
508c51ad4d4SStefano Zampini     ierr = VecLockPush(ctx->current_u);CHKERRQ(ierr);
509c51ad4d4SStefano Zampini   }
510c51ad4d4SStefano Zampini   ierr = VecLockPop(ctx->current_u);CHKERRQ(ierr);
511c51ad4d4SStefano Zampini   ierr = VecCopy(U,ctx->current_u);CHKERRQ(ierr);
512c51ad4d4SStefano Zampini   ierr = VecLockPush(ctx->current_u);CHKERRQ(ierr);
51352121784SBarry Smith   if (F) {
5146bf464f9SBarry Smith     if (ctx->current_f_allocated) {ierr = VecDestroy(&ctx->current_f);CHKERRQ(ierr);}
5153ec795f1SBarry Smith     ctx->current_f           = F;
516cfe22f5eSBarry Smith     ctx->current_f_allocated = PETSC_FALSE;
517cfe22f5eSBarry Smith   } else if (!ctx->current_f_allocated) {
51806c4b6baSBarry Smith     ierr = MatCreateVecs(J,NULL,&ctx->current_f);CHKERRQ(ierr);
5192205254eSKarl Rupp 
520cfe22f5eSBarry Smith     ctx->current_f_allocated = PETSC_TRUE;
52152121784SBarry Smith   }
522e884886fSBarry Smith   if (!ctx->w) {
523e884886fSBarry Smith     ierr = VecDuplicate(ctx->current_u,&ctx->w);CHKERRQ(ierr);
524e884886fSBarry Smith   }
525e884886fSBarry Smith   J->assembled = PETSC_TRUE;
526e884886fSBarry Smith   PetscFunctionReturn(0);
527e884886fSBarry Smith }
5285a576424SJed Brown 
529e884886fSBarry Smith typedef PetscErrorCode (*FCN3)(void*,Vec,Vec,PetscScalar*); /* force argument to next function to not be extern C*/
530b2573a8aSBarry Smith 
53140244768SBarry Smith static PetscErrorCode  MatMFFDSetCheckh_MFFD(Mat J,FCN3 fun,void *ectx)
532e884886fSBarry Smith {
533e884886fSBarry Smith   MatMFFD ctx = (MatMFFD)J->data;
534e884886fSBarry Smith 
535e884886fSBarry Smith   PetscFunctionBegin;
536e884886fSBarry Smith   ctx->checkh    = fun;
537e884886fSBarry Smith   ctx->checkhctx = ectx;
538e884886fSBarry Smith   PetscFunctionReturn(0);
539e884886fSBarry Smith }
540e884886fSBarry Smith 
5416aa9148fSLisandro Dalcin /*@C
5426aa9148fSLisandro Dalcin    MatMFFDSetOptionsPrefix - Sets the prefix used for searching for all
5436aa9148fSLisandro Dalcin    MatMFFD options in the database.
5446aa9148fSLisandro Dalcin 
5456aa9148fSLisandro Dalcin    Collective on Mat
5466aa9148fSLisandro Dalcin 
5476aa9148fSLisandro Dalcin    Input Parameter:
5486aa9148fSLisandro Dalcin +  A - the Mat context
5496aa9148fSLisandro Dalcin -  prefix - the prefix to prepend to all option names
5506aa9148fSLisandro Dalcin 
5516aa9148fSLisandro Dalcin    Notes:
5526aa9148fSLisandro Dalcin    A hyphen (-) must NOT be given at the beginning of the prefix name.
5536aa9148fSLisandro Dalcin    The first character of all runtime options is AUTOMATICALLY the hyphen.
5546aa9148fSLisandro Dalcin 
5556aa9148fSLisandro Dalcin    Level: advanced
5566aa9148fSLisandro Dalcin 
5576aa9148fSLisandro Dalcin .keywords: SNES, matrix-free, parameters
5586aa9148fSLisandro Dalcin 
559755b7f64SBarry Smith .seealso: MatSetFromOptions(), MatCreateSNESMF(), MatCreateMFFD()
5606aa9148fSLisandro Dalcin @*/
5617087cfbeSBarry Smith PetscErrorCode  MatMFFDSetOptionsPrefix(Mat mat,const char prefix[])
5626aa9148fSLisandro Dalcin 
5636aa9148fSLisandro Dalcin {
5640298fd71SBarry Smith   MatMFFD        mfctx = mat ? (MatMFFD)mat->data : (MatMFFD)NULL;
5656aa9148fSLisandro Dalcin   PetscErrorCode ierr;
5665fd66863SKarl Rupp 
5676aa9148fSLisandro Dalcin   PetscFunctionBegin;
5680700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
5690700a824SBarry Smith   PetscValidHeaderSpecific(mfctx,MATMFFD_CLASSID,1);
5706aa9148fSLisandro Dalcin   ierr = PetscObjectSetOptionsPrefix((PetscObject)mfctx,prefix);CHKERRQ(ierr);
5716aa9148fSLisandro Dalcin   PetscFunctionReturn(0);
5726aa9148fSLisandro Dalcin }
5736aa9148fSLisandro Dalcin 
57440244768SBarry Smith static PetscErrorCode  MatSetFromOptions_MFFD(PetscOptionItems *PetscOptionsObject,Mat mat)
575e884886fSBarry Smith {
57671f08665SBarry Smith   MatMFFD        mfctx = (MatMFFD)mat->data;
577e884886fSBarry Smith   PetscErrorCode ierr;
578ace3abfcSBarry Smith   PetscBool      flg;
579e884886fSBarry Smith   char           ftype[256];
580e884886fSBarry Smith 
581e884886fSBarry Smith   PetscFunctionBegin;
5820700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
5830700a824SBarry Smith   PetscValidHeaderSpecific(mfctx,MATMFFD_CLASSID,1);
5843194b578SJed Brown   ierr = PetscObjectOptionsBegin((PetscObject)mfctx);CHKERRQ(ierr);
585a264d7a6SBarry Smith   ierr = PetscOptionsFList("-mat_mffd_type","Matrix free type","MatMFFDSetType",MatMFFDList,((PetscObject)mfctx)->type_name,ftype,256,&flg);CHKERRQ(ierr);
586e884886fSBarry Smith   if (flg) {
587e884886fSBarry Smith     ierr = MatMFFDSetType(mat,ftype);CHKERRQ(ierr);
588e884886fSBarry Smith   }
589e884886fSBarry Smith 
590e884886fSBarry Smith   ierr = PetscOptionsReal("-mat_mffd_err","set sqrt relative error in function","MatMFFDSetFunctionError",mfctx->error_rel,&mfctx->error_rel,0);CHKERRQ(ierr);
591e884886fSBarry Smith   ierr = PetscOptionsInt("-mat_mffd_period","how often h is recomputed","MatMFFDSetPeriod",mfctx->recomputeperiod,&mfctx->recomputeperiod,0);CHKERRQ(ierr);
592e884886fSBarry Smith 
59390d69ab7SBarry Smith   flg  = PETSC_FALSE;
5940298fd71SBarry Smith   ierr = PetscOptionsBool("-mat_mffd_check_positivity","Insure that U + h*a is nonnegative","MatMFFDSetCheckh",flg,&flg,NULL);CHKERRQ(ierr);
595e884886fSBarry Smith   if (flg) {
596e884886fSBarry Smith     ierr = MatMFFDSetCheckh(mat,MatMFFDCheckPositivity,0);CHKERRQ(ierr);
597e884886fSBarry Smith   }
598e884886fSBarry Smith   if (mfctx->ops->setfromoptions) {
599e55864a3SBarry Smith     ierr = (*mfctx->ops->setfromoptions)(PetscOptionsObject,mfctx);CHKERRQ(ierr);
600e884886fSBarry Smith   }
601e884886fSBarry Smith   ierr = PetscOptionsEnd();CHKERRQ(ierr);
602e884886fSBarry Smith   PetscFunctionReturn(0);
603e884886fSBarry Smith }
604e884886fSBarry Smith 
60540244768SBarry Smith static PetscErrorCode  MatMFFDSetPeriod_MFFD(Mat mat,PetscInt period)
606bc0f08ceSBarry Smith {
607bc0f08ceSBarry Smith   MatMFFD ctx = (MatMFFD)mat->data;
608bc0f08ceSBarry Smith 
609bc0f08ceSBarry Smith   PetscFunctionBegin;
610bc0f08ceSBarry Smith   ctx->recomputeperiod = period;
611bc0f08ceSBarry Smith   PetscFunctionReturn(0);
612bc0f08ceSBarry Smith }
613bc0f08ceSBarry Smith 
61440244768SBarry Smith static PetscErrorCode  MatMFFDSetFunction_MFFD(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx)
615bc0f08ceSBarry Smith {
616bc0f08ceSBarry Smith   MatMFFD ctx = (MatMFFD)mat->data;
617bc0f08ceSBarry Smith 
618bc0f08ceSBarry Smith   PetscFunctionBegin;
619bc0f08ceSBarry Smith   ctx->func    = func;
620bc0f08ceSBarry Smith   ctx->funcctx = funcctx;
621bc0f08ceSBarry Smith   PetscFunctionReturn(0);
622bc0f08ceSBarry Smith }
623bc0f08ceSBarry Smith 
62440244768SBarry Smith static PetscErrorCode  MatMFFDSetFunctionError_MFFD(Mat mat,PetscReal error)
625bc0f08ceSBarry Smith {
626bc0f08ceSBarry Smith   MatMFFD ctx = (MatMFFD)mat->data;
627bc0f08ceSBarry Smith 
628bc0f08ceSBarry Smith   PetscFunctionBegin;
629bc0f08ceSBarry Smith   if (error != PETSC_DEFAULT) ctx->error_rel = error;
630bc0f08ceSBarry Smith   PetscFunctionReturn(0);
631bc0f08ceSBarry Smith }
632bc0f08ceSBarry Smith 
6333b49f96aSBarry Smith static PetscErrorCode MatMissingDiagonal_MFFD(Mat A,PetscBool  *missing,PetscInt *d)
6343b49f96aSBarry Smith {
6353b49f96aSBarry Smith   PetscFunctionBegin;
6363b49f96aSBarry Smith   *missing = PETSC_FALSE;
6373b49f96aSBarry Smith   PetscFunctionReturn(0);
6383b49f96aSBarry Smith }
6393b49f96aSBarry Smith 
640e884886fSBarry Smith /*MC
641e884886fSBarry Smith   MATMFFD - MATMFFD = "mffd" - A matrix free matrix type.
642e884886fSBarry Smith 
643e884886fSBarry Smith   Level: advanced
644e884886fSBarry Smith 
645755b7f64SBarry Smith .seealso: MatCreateMFFD(), MatCreateSNESMF(), MatMFFDSetFunction(), MatMFFDSetType(),
646755b7f64SBarry Smith           MatMFFDSetFunctionError(), MatMFFDDSSetUmin(), MatMFFDSetFunction()
647755b7f64SBarry Smith           MatMFFDSetHHistory(), MatMFFDResetHHistory(), MatCreateSNESMF(),
648755b7f64SBarry Smith           MatMFFDGetH(),
649e884886fSBarry Smith M*/
6508cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_MFFD(Mat A)
651e884886fSBarry Smith {
652e884886fSBarry Smith   MatMFFD        mfctx;
653e884886fSBarry Smith   PetscErrorCode ierr;
654e884886fSBarry Smith 
655e884886fSBarry Smith   PetscFunctionBegin;
656607a6623SBarry Smith   ierr = MatMFFDInitializePackage();CHKERRQ(ierr);
657e884886fSBarry Smith 
65873107ff1SLisandro Dalcin   ierr = PetscHeaderCreate(mfctx,MATMFFD_CLASSID,"MatMFFD","Matrix-free Finite Differencing","Mat",PetscObjectComm((PetscObject)A),MatDestroy_MFFD,MatView_MFFD);CHKERRQ(ierr);
6592205254eSKarl Rupp 
660e884886fSBarry Smith   mfctx->error_rel                = PETSC_SQRT_MACHINE_EPSILON;
661e884886fSBarry Smith   mfctx->recomputeperiod          = 1;
662e884886fSBarry Smith   mfctx->count                    = 0;
663e884886fSBarry Smith   mfctx->currenth                 = 0.0;
6640298fd71SBarry Smith   mfctx->historyh                 = NULL;
665e884886fSBarry Smith   mfctx->ncurrenth                = 0;
666e884886fSBarry Smith   mfctx->maxcurrenth              = 0;
6677adad957SLisandro Dalcin   ((PetscObject)mfctx)->type_name = 0;
668e884886fSBarry Smith 
669e884886fSBarry Smith   mfctx->vshift = 0.0;
670e884886fSBarry Smith   mfctx->vscale = 1.0;
671e884886fSBarry Smith 
672e884886fSBarry Smith   /*
673e884886fSBarry Smith      Create the empty data structure to contain compute-h routines.
674e884886fSBarry Smith      These will be filled in below from the command line options or
675e884886fSBarry Smith      a later call with MatMFFDSetType() or if that is not called
676e884886fSBarry Smith      then it will default in the first use of MatMult_MFFD()
677e884886fSBarry Smith   */
678e884886fSBarry Smith   mfctx->ops->compute        = 0;
679e884886fSBarry Smith   mfctx->ops->destroy        = 0;
680e884886fSBarry Smith   mfctx->ops->view           = 0;
681e884886fSBarry Smith   mfctx->ops->setfromoptions = 0;
682e884886fSBarry Smith   mfctx->hctx                = 0;
683e884886fSBarry Smith 
684e884886fSBarry Smith   mfctx->func    = 0;
685e884886fSBarry Smith   mfctx->funcctx = 0;
6860298fd71SBarry Smith   mfctx->w       = NULL;
687e884886fSBarry Smith 
688e884886fSBarry Smith   A->data = mfctx;
689e884886fSBarry Smith 
690e884886fSBarry Smith   A->ops->mult            = MatMult_MFFD;
691e884886fSBarry Smith   A->ops->destroy         = MatDestroy_MFFD;
692e884886fSBarry Smith   A->ops->view            = MatView_MFFD;
693e884886fSBarry Smith   A->ops->assemblyend     = MatAssemblyEnd_MFFD;
694e884886fSBarry Smith   A->ops->scale           = MatScale_MFFD;
695e884886fSBarry Smith   A->ops->shift           = MatShift_MFFD;
696c3bb7e23SBarry Smith   A->ops->diagonalscale   = MatDiagonalScale_MFFD;
697c3bb7e23SBarry Smith   A->ops->diagonalset     = MatDiagonalSet_MFFD;
6989c6ac3b3SBarry Smith   A->ops->setfromoptions  = MatSetFromOptions_MFFD;
6993b49f96aSBarry Smith   A->ops->missingdiagonal = MatMissingDiagonal_MFFD;
700e884886fSBarry Smith   A->assembled            = PETSC_TRUE;
701e884886fSBarry Smith 
702bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetBase_C",MatMFFDSetBase_MFFD);CHKERRQ(ierr);
703bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctioniBase_C",MatMFFDSetFunctioniBase_MFFD);CHKERRQ(ierr);
704bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctioni_C",MatMFFDSetFunctioni_MFFD);CHKERRQ(ierr);
705bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunction_C",MatMFFDSetFunction_MFFD);CHKERRQ(ierr);
706bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetCheckh_C",MatMFFDSetCheckh_MFFD);CHKERRQ(ierr);
707bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetPeriod_C",MatMFFDSetPeriod_MFFD);CHKERRQ(ierr);
708bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctionError_C",MatMFFDSetFunctionError_MFFD);CHKERRQ(ierr);
709bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDResetHHistory_C",MatMFFDResetHHistory_MFFD);CHKERRQ(ierr);
7102205254eSKarl Rupp 
711e884886fSBarry Smith   mfctx->mat = A;
7122205254eSKarl Rupp 
713e884886fSBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,MATMFFD);CHKERRQ(ierr);
714e884886fSBarry Smith   PetscFunctionReturn(0);
715e884886fSBarry Smith }
716e884886fSBarry Smith 
717e884886fSBarry Smith /*@
718e884886fSBarry Smith    MatCreateMFFD - Creates a matrix-free matrix. See also MatCreateSNESMF()
719e884886fSBarry Smith 
720e884886fSBarry Smith    Collective on Vec
721e884886fSBarry Smith 
722e884886fSBarry Smith    Input Parameters:
723fef1beadSBarry Smith +  comm - MPI communicator
724fef1beadSBarry Smith .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
725fef1beadSBarry Smith            This value should be the same as the local size used in creating the
726fef1beadSBarry Smith            y vector for the matrix-vector product y = Ax.
727fef1beadSBarry Smith .  n - This value should be the same as the local size used in creating the
728fef1beadSBarry Smith        x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have
729fef1beadSBarry Smith        calculated if N is given) For square matrices n is almost always m.
730fef1beadSBarry Smith .  M - number of global rows (or PETSC_DETERMINE to have calculated if m is given)
731fef1beadSBarry Smith -  N - number of global columns (or PETSC_DETERMINE to have calculated if n is given)
732fef1beadSBarry Smith 
733e884886fSBarry Smith 
734e884886fSBarry Smith    Output Parameter:
735e884886fSBarry Smith .  J - the matrix-free matrix
736e884886fSBarry Smith 
7379c6ac3b3SBarry Smith    Options Database Keys: call MatSetFromOptions() to trigger these
7389c6ac3b3SBarry Smith +  -mat_mffd_type - wp or ds (see MATMFFD_WP or MATMFFD_DS)
7399c6ac3b3SBarry Smith -  -mat_mffd_err - square root of estimated relative error in function evaluation
7409c6ac3b3SBarry Smith -  -mat_mffd_period - how often h is recomputed, defaults to 1, everytime
7419c6ac3b3SBarry Smith 
7429c6ac3b3SBarry Smith 
743e884886fSBarry Smith    Level: advanced
744e884886fSBarry Smith 
745e884886fSBarry Smith    Notes:
746e884886fSBarry Smith    The matrix-free matrix context merely contains the function pointers
747e884886fSBarry Smith    and work space for performing finite difference approximations of
748e884886fSBarry Smith    Jacobian-vector products, F'(u)*a,
749e884886fSBarry Smith 
750e884886fSBarry Smith    The default code uses the following approach to compute h
751e884886fSBarry Smith 
752e884886fSBarry Smith .vb
753e884886fSBarry Smith      F'(u)*a = [F(u+h*a) - F(u)]/h where
754e884886fSBarry Smith      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
755e884886fSBarry Smith        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
756e884886fSBarry Smith  where
757e884886fSBarry Smith      error_rel = square root of relative error in function evaluation
758e884886fSBarry Smith      umin = minimum iterate parameter
759e884886fSBarry Smith .ve
760e884886fSBarry Smith 
761ca93e954SBarry Smith    You can call SNESSetJacobian() with MatMFFDComputeJacobian() if you are using matrix and not a different
762ca93e954SBarry Smith    preconditioner matrix
763ca93e954SBarry Smith 
764e884886fSBarry Smith    The user can set the error_rel via MatMFFDSetFunctionError() and
765a7f22e61SSatish Balay    umin via MatMFFDDSSetUmin(); see Users-Manual: ch_snes for details.
766e884886fSBarry Smith 
767e884886fSBarry Smith    The user should call MatDestroy() when finished with the matrix-free
768e884886fSBarry Smith    matrix context.
769e884886fSBarry Smith 
770e884886fSBarry Smith    Options Database Keys:
771e884886fSBarry Smith +  -mat_mffd_err <error_rel> - Sets error_rel
772e884886fSBarry Smith .  -mat_mffd_unim <umin> - Sets umin (for default PETSc routine that computes h only)
773e884886fSBarry Smith -  -mat_mffd_check_positivity
774e884886fSBarry Smith 
775e884886fSBarry Smith .keywords: default, matrix-free, create, matrix
776e884886fSBarry Smith 
7771d0fab5eSBarry Smith .seealso: MatDestroy(), MatMFFDSetFunctionError(), MatMFFDDSSetUmin(), MatMFFDSetFunction()
778e884886fSBarry Smith           MatMFFDSetHHistory(), MatMFFDResetHHistory(), MatCreateSNESMF(),
77981242352SJed Brown           MatMFFDGetH(), MatMFFDRegister(), MatMFFDComputeJacobian()
780e884886fSBarry Smith 
781e884886fSBarry Smith @*/
7827087cfbeSBarry Smith PetscErrorCode  MatCreateMFFD(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *J)
783e884886fSBarry Smith {
784e884886fSBarry Smith   PetscErrorCode ierr;
785e884886fSBarry Smith 
786e884886fSBarry Smith   PetscFunctionBegin;
787e884886fSBarry Smith   ierr = MatCreate(comm,J);CHKERRQ(ierr);
788fef1beadSBarry Smith   ierr = MatSetSizes(*J,m,n,M,N);CHKERRQ(ierr);
789e884886fSBarry Smith   ierr = MatSetType(*J,MATMFFD);CHKERRQ(ierr);
790be7c243fSJed Brown   ierr = MatSetUp(*J);CHKERRQ(ierr);
791e884886fSBarry Smith   PetscFunctionReturn(0);
792e884886fSBarry Smith }
793e884886fSBarry Smith 
794e884886fSBarry Smith /*@
795e884886fSBarry Smith    MatMFFDGetH - Gets the last value that was used as the differencing
796e884886fSBarry Smith    parameter.
797e884886fSBarry Smith 
798e884886fSBarry Smith    Not Collective
799e884886fSBarry Smith 
800e884886fSBarry Smith    Input Parameters:
801e884886fSBarry Smith .  mat - the matrix obtained with MatCreateSNESMF()
802e884886fSBarry Smith 
803e884886fSBarry Smith    Output Paramter:
804e884886fSBarry Smith .  h - the differencing step size
805e884886fSBarry Smith 
806e884886fSBarry Smith    Level: advanced
807e884886fSBarry Smith 
808e884886fSBarry Smith .keywords: SNES, matrix-free, parameters
809e884886fSBarry Smith 
8101d0fab5eSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDSetHHistory(), MatCreateMFFD(), MATMFFD, MatMFFDResetHHistory()
811e884886fSBarry Smith @*/
8127087cfbeSBarry Smith PetscErrorCode  MatMFFDGetH(Mat mat,PetscScalar *h)
813e884886fSBarry Smith {
814e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)mat->data;
815bc0f08ceSBarry Smith   PetscErrorCode ierr;
816bc0f08ceSBarry Smith   PetscBool      match;
817e884886fSBarry Smith 
818e884886fSBarry Smith   PetscFunctionBegin;
81988b4c220SStefano Zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
82088b4c220SStefano Zampini   PetscValidPointer(h,2);
821251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)mat,MATMFFD,&match);CHKERRQ(ierr);
822ce94432eSBarry Smith   if (!match) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONG,"Not a MFFD matrix");
823bc0f08ceSBarry Smith 
824e884886fSBarry Smith   *h = ctx->currenth;
825e884886fSBarry Smith   PetscFunctionReturn(0);
826e884886fSBarry Smith }
827e884886fSBarry Smith 
828e884886fSBarry Smith /*@C
829e884886fSBarry Smith    MatMFFDSetFunction - Sets the function used in applying the matrix free.
830e884886fSBarry Smith 
8313f9fe445SBarry Smith    Logically Collective on Mat
832e884886fSBarry Smith 
833e884886fSBarry Smith    Input Parameters:
83414f633e2SBarry Smith +  mat - the matrix free matrix created via MatCreateSNESMF() or MatCreateMFFD()
835e884886fSBarry Smith .  func - the function to use
836e884886fSBarry Smith -  funcctx - optional function context passed to function
837e884886fSBarry Smith 
83814f633e2SBarry Smith    Calling Sequence of func:
83914f633e2SBarry Smith $     func (void *funcctx, Vec x, Vec f)
84014f633e2SBarry Smith 
84114f633e2SBarry Smith +  funcctx - user provided context
84214f633e2SBarry Smith .  x - input vector
84314f633e2SBarry Smith -  f - computed output function
84414f633e2SBarry Smith 
845e884886fSBarry Smith    Level: advanced
846e884886fSBarry Smith 
847e884886fSBarry Smith    Notes:
848e884886fSBarry Smith     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
849e884886fSBarry Smith     matrix inside your compute Jacobian routine
850e884886fSBarry Smith 
8513ec795f1SBarry Smith     If this is not set then it will use the function set with SNESSetFunction() if MatCreateSNESMF() was used.
852e884886fSBarry Smith 
853e884886fSBarry Smith .keywords: SNES, matrix-free, function
854e884886fSBarry Smith 
8551d0fab5eSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD,
8561d0fab5eSBarry Smith           MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction()
857e884886fSBarry Smith @*/
8587087cfbeSBarry Smith PetscErrorCode  MatMFFDSetFunction(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx)
859e884886fSBarry Smith {
860bc0f08ceSBarry Smith   PetscErrorCode ierr;
861e884886fSBarry Smith 
862e884886fSBarry Smith   PetscFunctionBegin;
86388b4c220SStefano Zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
864bc0f08ceSBarry Smith   ierr = PetscTryMethod(mat,"MatMFFDSetFunction_C",(Mat,PetscErrorCode (*)(void*,Vec,Vec),void*),(mat,func,funcctx));CHKERRQ(ierr);
865e884886fSBarry Smith   PetscFunctionReturn(0);
866e884886fSBarry Smith }
867e884886fSBarry Smith 
868e884886fSBarry Smith /*@C
869e884886fSBarry Smith    MatMFFDSetFunctioni - Sets the function for a single component
870e884886fSBarry Smith 
8713f9fe445SBarry Smith    Logically Collective on Mat
872e884886fSBarry Smith 
873e884886fSBarry Smith    Input Parameters:
874e884886fSBarry Smith +  mat - the matrix free matrix created via MatCreateSNESMF()
875e884886fSBarry Smith -  funci - the function to use
876e884886fSBarry Smith 
877e884886fSBarry Smith    Level: advanced
878e884886fSBarry Smith 
879e884886fSBarry Smith    Notes:
880e884886fSBarry Smith     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
88194eb4328SStefano Zampini     matrix inside your compute Jacobian routine.
88294eb4328SStefano Zampini     This function is necessary to compute the diagonal of the matrix.
883486afcceSStefano Zampini     funci must not contain any MPI call as it is called inside a loop on the local portion of the vector.
884e884886fSBarry Smith 
885e884886fSBarry Smith .keywords: SNES, matrix-free, function
886e884886fSBarry Smith 
88794eb4328SStefano Zampini .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction(), MatGetDiagonal()
8881d0fab5eSBarry Smith 
889e884886fSBarry Smith @*/
8907087cfbeSBarry Smith PetscErrorCode  MatMFFDSetFunctioni(Mat mat,PetscErrorCode (*funci)(void*,PetscInt,Vec,PetscScalar*))
891e884886fSBarry Smith {
8924ac538c5SBarry Smith   PetscErrorCode ierr;
893e884886fSBarry Smith 
894e884886fSBarry Smith   PetscFunctionBegin;
8950700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
8964ac538c5SBarry Smith   ierr = PetscTryMethod(mat,"MatMFFDSetFunctioni_C",(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)),(mat,funci));CHKERRQ(ierr);
897e884886fSBarry Smith   PetscFunctionReturn(0);
898e884886fSBarry Smith }
899e884886fSBarry Smith 
900e884886fSBarry Smith /*@C
901e884886fSBarry Smith    MatMFFDSetFunctioniBase - Sets the base vector for a single component function evaluation
902e884886fSBarry Smith 
9033f9fe445SBarry Smith    Logically Collective on Mat
904e884886fSBarry Smith 
905e884886fSBarry Smith    Input Parameters:
906e884886fSBarry Smith +  mat - the matrix free matrix created via MatCreateSNESMF()
907e884886fSBarry Smith -  func - the function to use
908e884886fSBarry Smith 
909e884886fSBarry Smith    Level: advanced
910e884886fSBarry Smith 
911e884886fSBarry Smith    Notes:
912e884886fSBarry Smith     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
91394eb4328SStefano Zampini     matrix inside your compute Jacobian routine.
91494eb4328SStefano Zampini     This function is necessary to compute the diagonal of the matrix.
915e884886fSBarry Smith 
916e884886fSBarry Smith 
917e884886fSBarry Smith .keywords: SNES, matrix-free, function
918e884886fSBarry Smith 
919e884886fSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD
92094eb4328SStefano Zampini           MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction(), MatGetDiagonal()
921e884886fSBarry Smith @*/
9227087cfbeSBarry Smith PetscErrorCode  MatMFFDSetFunctioniBase(Mat mat,PetscErrorCode (*func)(void*,Vec))
923e884886fSBarry Smith {
9244ac538c5SBarry Smith   PetscErrorCode ierr;
925e884886fSBarry Smith 
926e884886fSBarry Smith   PetscFunctionBegin;
9270700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
9284ac538c5SBarry Smith   ierr = PetscTryMethod(mat,"MatMFFDSetFunctioniBase_C",(Mat,PetscErrorCode (*)(void*,Vec)),(mat,func));CHKERRQ(ierr);
929e884886fSBarry Smith   PetscFunctionReturn(0);
930e884886fSBarry Smith }
931e884886fSBarry Smith 
932e884886fSBarry Smith /*@
933e884886fSBarry Smith    MatMFFDSetPeriod - Sets how often h is recomputed, by default it is everytime
934e884886fSBarry Smith 
9353f9fe445SBarry Smith    Logically Collective on Mat
936e884886fSBarry Smith 
937e884886fSBarry Smith    Input Parameters:
938e884886fSBarry Smith +  mat - the matrix free matrix created via MatCreateSNESMF()
939e884886fSBarry Smith -  period - 1 for everytime, 2 for every second etc
940e884886fSBarry Smith 
941e884886fSBarry Smith    Options Database Keys:
942e884886fSBarry Smith +  -mat_mffd_period <period>
943e884886fSBarry Smith 
944e884886fSBarry Smith    Level: advanced
945e884886fSBarry Smith 
946e884886fSBarry Smith 
947e884886fSBarry Smith .keywords: SNES, matrix-free, parameters
948e884886fSBarry Smith 
949e884886fSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(),
9501d0fab5eSBarry Smith           MatMFFDSetHHistory(), MatMFFDResetHHistory()
951e884886fSBarry Smith @*/
9527087cfbeSBarry Smith PetscErrorCode  MatMFFDSetPeriod(Mat mat,PetscInt period)
953e884886fSBarry Smith {
954bc0f08ceSBarry Smith   PetscErrorCode ierr;
955e884886fSBarry Smith 
956e884886fSBarry Smith   PetscFunctionBegin;
95788b4c220SStefano Zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
95888b4c220SStefano Zampini   PetscValidLogicalCollectiveInt(mat,period,2);
959bc0f08ceSBarry Smith   ierr = PetscTryMethod(mat,"MatMFFDSetPeriod_C",(Mat,PetscInt),(mat,period));CHKERRQ(ierr);
960e884886fSBarry Smith   PetscFunctionReturn(0);
961e884886fSBarry Smith }
962e884886fSBarry Smith 
963e884886fSBarry Smith /*@
964e884886fSBarry Smith    MatMFFDSetFunctionError - Sets the error_rel for the approximation of
965e884886fSBarry Smith    matrix-vector products using finite differences.
966e884886fSBarry Smith 
9673f9fe445SBarry Smith    Logically Collective on Mat
968e884886fSBarry Smith 
969e884886fSBarry Smith    Input Parameters:
970e884886fSBarry Smith +  mat - the matrix free matrix created via MatCreateMFFD() or MatCreateSNESMF()
971e884886fSBarry Smith -  error_rel - relative error (should be set to the square root of
972e884886fSBarry Smith                the relative error in the function evaluations)
973e884886fSBarry Smith 
974e884886fSBarry Smith    Options Database Keys:
975e884886fSBarry Smith +  -mat_mffd_err <error_rel> - Sets error_rel
976e884886fSBarry Smith 
977e884886fSBarry Smith    Level: advanced
978e884886fSBarry Smith 
979e884886fSBarry Smith    Notes:
980e884886fSBarry Smith    The default matrix-free matrix-vector product routine computes
981e884886fSBarry Smith .vb
982e884886fSBarry Smith      F'(u)*a = [F(u+h*a) - F(u)]/h where
983e884886fSBarry Smith      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
984e884886fSBarry Smith        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   else
985e884886fSBarry Smith .ve
986e884886fSBarry Smith 
987e884886fSBarry Smith .keywords: SNES, matrix-free, parameters
988e884886fSBarry Smith 
989e884886fSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD
9901d0fab5eSBarry Smith           MatMFFDSetHHistory(), MatMFFDResetHHistory()
991e884886fSBarry Smith @*/
9927087cfbeSBarry Smith PetscErrorCode  MatMFFDSetFunctionError(Mat mat,PetscReal error)
993e884886fSBarry Smith {
994bc0f08ceSBarry Smith   PetscErrorCode ierr;
995e884886fSBarry Smith 
996e884886fSBarry Smith   PetscFunctionBegin;
99788b4c220SStefano Zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
99888b4c220SStefano Zampini   PetscValidLogicalCollectiveReal(mat,error,2);
999bc0f08ceSBarry Smith   ierr = PetscTryMethod(mat,"MatMFFDSetFunctionError_C",(Mat,PetscReal),(mat,error));CHKERRQ(ierr);
1000e884886fSBarry Smith   PetscFunctionReturn(0);
1001e884886fSBarry Smith }
1002e884886fSBarry Smith 
1003e884886fSBarry Smith /*@
1004e884886fSBarry Smith    MatMFFDSetHHistory - Sets an array to collect a history of the
1005e884886fSBarry Smith    differencing values (h) computed for the matrix-free product.
1006e884886fSBarry Smith 
10073f9fe445SBarry Smith    Logically Collective on Mat
1008e884886fSBarry Smith 
1009e884886fSBarry Smith    Input Parameters:
1010e884886fSBarry Smith +  J - the matrix-free matrix context
1011e884886fSBarry Smith .  histroy - space to hold the history
1012e884886fSBarry Smith -  nhistory - number of entries in history, if more entries are generated than
1013e884886fSBarry Smith               nhistory, then the later ones are discarded
1014e884886fSBarry Smith 
1015e884886fSBarry Smith    Level: advanced
1016e884886fSBarry Smith 
1017e884886fSBarry Smith    Notes:
1018e884886fSBarry Smith    Use MatMFFDResetHHistory() to reset the history counter and collect
1019e884886fSBarry Smith    a new batch of differencing parameters, h.
1020e884886fSBarry Smith 
1021e884886fSBarry Smith .keywords: SNES, matrix-free, h history, differencing history
1022e884886fSBarry Smith 
1023e884886fSBarry Smith .seealso: MatMFFDGetH(), MatCreateSNESMF(),
10241d0fab5eSBarry Smith           MatMFFDResetHHistory(), MatMFFDSetFunctionError()
1025e884886fSBarry Smith 
1026e884886fSBarry Smith @*/
10277087cfbeSBarry Smith PetscErrorCode  MatMFFDSetHHistory(Mat J,PetscScalar history[],PetscInt nhistory)
1028e884886fSBarry Smith {
1029e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)J->data;
1030bc0f08ceSBarry Smith   PetscErrorCode ierr;
1031bc0f08ceSBarry Smith   PetscBool      match;
1032e884886fSBarry Smith 
1033e884886fSBarry Smith   PetscFunctionBegin;
103488b4c220SStefano Zampini   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
103588b4c220SStefano Zampini   if (history) PetscValidPointer(history,2);
103688b4c220SStefano Zampini   PetscValidLogicalCollectiveInt(J,nhistory,3);
1037251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)J,MATMFFD,&match);CHKERRQ(ierr);
1038ce94432eSBarry Smith   if (!match) SETERRQ(PetscObjectComm((PetscObject)J),PETSC_ERR_ARG_WRONG,"Not a MFFD matrix");
1039e884886fSBarry Smith   ctx->historyh    = history;
1040e884886fSBarry Smith   ctx->maxcurrenth = nhistory;
104175567043SBarry Smith   ctx->currenth    = 0.;
1042e884886fSBarry Smith   PetscFunctionReturn(0);
1043e884886fSBarry Smith }
1044e884886fSBarry Smith 
1045e884886fSBarry Smith /*@
1046e884886fSBarry Smith    MatMFFDResetHHistory - Resets the counter to zero to begin
1047e884886fSBarry Smith    collecting a new set of differencing histories.
1048e884886fSBarry Smith 
10493f9fe445SBarry Smith    Logically Collective on Mat
1050e884886fSBarry Smith 
1051e884886fSBarry Smith    Input Parameters:
1052e884886fSBarry Smith .  J - the matrix-free matrix context
1053e884886fSBarry Smith 
1054e884886fSBarry Smith    Level: advanced
1055e884886fSBarry Smith 
1056e884886fSBarry Smith    Notes:
1057e884886fSBarry Smith    Use MatMFFDSetHHistory() to create the original history counter.
1058e884886fSBarry Smith 
1059e884886fSBarry Smith .keywords: SNES, matrix-free, h history, differencing history
1060e884886fSBarry Smith 
1061e884886fSBarry Smith .seealso: MatMFFDGetH(), MatCreateSNESMF(),
10621d0fab5eSBarry Smith           MatMFFDSetHHistory(), MatMFFDSetFunctionError()
1063e884886fSBarry Smith 
1064e884886fSBarry Smith @*/
10657087cfbeSBarry Smith PetscErrorCode  MatMFFDResetHHistory(Mat J)
1066e884886fSBarry Smith {
1067bc0f08ceSBarry Smith   PetscErrorCode ierr;
1068e884886fSBarry Smith 
1069e884886fSBarry Smith   PetscFunctionBegin;
107088b4c220SStefano Zampini   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
1071bc0f08ceSBarry Smith   ierr = PetscTryMethod(J,"MatMFFDResetHHistory_C",(Mat),(J));CHKERRQ(ierr);
1072e884886fSBarry Smith   PetscFunctionReturn(0);
1073e884886fSBarry Smith }
1074e884886fSBarry Smith 
1075*487a658cSBarry Smith /*@
1076e884886fSBarry Smith     MatMFFDSetBase - Sets the vector U at which matrix vector products of the
1077e884886fSBarry Smith         Jacobian are computed
1078e884886fSBarry Smith 
10793f9fe445SBarry Smith     Logically Collective on Mat
1080e884886fSBarry Smith 
1081e884886fSBarry Smith     Input Parameters:
1082e884886fSBarry Smith +   J - the MatMFFD matrix
10833ec795f1SBarry Smith .   U - the vector
1084bcddec3dSBarry Smith -   F - (optional) vector that contains F(u) if it has been already computed
1085e884886fSBarry Smith 
108695452b02SPatrick Sanan     Notes:
108795452b02SPatrick Sanan     This is rarely used directly
1088e884886fSBarry Smith 
10898af5ae88SBarry Smith     If F is provided then it is not recomputed. Otherwise the function is evaluated at the base
10908af5ae88SBarry Smith     point during the first MatMult() after each call to MatMFFDSetBase().
1091dff2f722SBarry Smith 
1092e884886fSBarry Smith     Level: advanced
1093e884886fSBarry Smith 
1094e884886fSBarry Smith @*/
10957087cfbeSBarry Smith PetscErrorCode  MatMFFDSetBase(Mat J,Vec U,Vec F)
1096e884886fSBarry Smith {
10974ac538c5SBarry Smith   PetscErrorCode ierr;
1098e884886fSBarry Smith 
1099e884886fSBarry Smith   PetscFunctionBegin;
11000700a824SBarry Smith   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
11010700a824SBarry Smith   PetscValidHeaderSpecific(U,VEC_CLASSID,2);
11020700a824SBarry Smith   if (F) PetscValidHeaderSpecific(F,VEC_CLASSID,3);
11034ac538c5SBarry Smith   ierr = PetscTryMethod(J,"MatMFFDSetBase_C",(Mat,Vec,Vec),(J,U,F));CHKERRQ(ierr);
1104e884886fSBarry Smith   PetscFunctionReturn(0);
1105e884886fSBarry Smith }
1106e884886fSBarry Smith 
1107e884886fSBarry Smith /*@C
1108e884886fSBarry Smith     MatMFFDSetCheckh - Sets a function that checks the computed h and adjusts
1109e884886fSBarry Smith         it to satisfy some criteria
1110e884886fSBarry Smith 
11113f9fe445SBarry Smith     Logically Collective on Mat
1112e884886fSBarry Smith 
1113e884886fSBarry Smith     Input Parameters:
1114e884886fSBarry Smith +   J - the MatMFFD matrix
1115e884886fSBarry Smith .   fun - the function that checks h
1116e884886fSBarry Smith -   ctx - any context needed by the function
1117e884886fSBarry Smith 
1118e884886fSBarry Smith     Options Database Keys:
1119e884886fSBarry Smith .   -mat_mffd_check_positivity
1120e884886fSBarry Smith 
1121e884886fSBarry Smith     Level: advanced
1122e884886fSBarry Smith 
112395452b02SPatrick Sanan     Notes:
112495452b02SPatrick Sanan     For example, MatMFFDCheckPositivity() insures that all entries
1125e884886fSBarry Smith        of U + h*a are non-negative
1126e884886fSBarry Smith 
1127755b7f64SBarry Smith      The function you provide is called after the default h has been computed and allows you to
1128755b7f64SBarry Smith      modify it.
1129755b7f64SBarry Smith 
1130755b7f64SBarry Smith .seealso:  MatMFFDCheckPositivity()
1131e884886fSBarry Smith @*/
11327087cfbeSBarry Smith PetscErrorCode  MatMFFDSetCheckh(Mat J,PetscErrorCode (*fun)(void*,Vec,Vec,PetscScalar*),void *ctx)
1133e884886fSBarry Smith {
11344ac538c5SBarry Smith   PetscErrorCode ierr;
1135e884886fSBarry Smith 
1136e884886fSBarry Smith   PetscFunctionBegin;
11370700a824SBarry Smith   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
11384ac538c5SBarry Smith   ierr = PetscTryMethod(J,"MatMFFDSetCheckh_C",(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*),(J,fun,ctx));CHKERRQ(ierr);
1139e884886fSBarry Smith   PetscFunctionReturn(0);
1140e884886fSBarry Smith }
1141e884886fSBarry Smith 
1142e884886fSBarry Smith /*@
1143e884886fSBarry Smith     MatMFFDCheckPositivity - Checks that all entries in U + h*a are positive or
1144e884886fSBarry Smith         zero, decreases h until this is satisfied.
1145e884886fSBarry Smith 
11463f9fe445SBarry Smith     Logically Collective on Vec
1147e884886fSBarry Smith 
1148e884886fSBarry Smith     Input Parameters:
1149e884886fSBarry Smith +   U - base vector that is added to
1150e884886fSBarry Smith .   a - vector that is added
1151e884886fSBarry Smith .   h - scaling factor on a
1152e884886fSBarry Smith -   dummy - context variable (unused)
1153e884886fSBarry Smith 
1154e884886fSBarry Smith     Options Database Keys:
1155e884886fSBarry Smith .   -mat_mffd_check_positivity
1156e884886fSBarry Smith 
1157e884886fSBarry Smith     Level: advanced
1158e884886fSBarry Smith 
115995452b02SPatrick Sanan     Notes:
116095452b02SPatrick Sanan     This is rarely used directly, rather it is passed as an argument to
1161e884886fSBarry Smith            MatMFFDSetCheckh()
1162e884886fSBarry Smith 
1163e884886fSBarry Smith .seealso:  MatMFFDSetCheckh()
1164e884886fSBarry Smith @*/
11657087cfbeSBarry Smith PetscErrorCode  MatMFFDCheckPositivity(void *dummy,Vec U,Vec a,PetscScalar *h)
1166e884886fSBarry Smith {
1167e884886fSBarry Smith   PetscReal      val, minval;
1168e884886fSBarry Smith   PetscScalar    *u_vec, *a_vec;
1169e884886fSBarry Smith   PetscErrorCode ierr;
1170e884886fSBarry Smith   PetscInt       i,n;
1171e884886fSBarry Smith   MPI_Comm       comm;
1172e884886fSBarry Smith 
1173e884886fSBarry Smith   PetscFunctionBegin;
117488b4c220SStefano Zampini   PetscValidHeaderSpecific(U,VEC_CLASSID,2);
117588b4c220SStefano Zampini   PetscValidHeaderSpecific(a,VEC_CLASSID,3);
117688b4c220SStefano Zampini   PetscValidPointer(h,4);
1177e884886fSBarry Smith   ierr   = PetscObjectGetComm((PetscObject)U,&comm);CHKERRQ(ierr);
1178e884886fSBarry Smith   ierr   = VecGetArray(U,&u_vec);CHKERRQ(ierr);
1179e884886fSBarry Smith   ierr   = VecGetArray(a,&a_vec);CHKERRQ(ierr);
1180e884886fSBarry Smith   ierr   = VecGetLocalSize(U,&n);CHKERRQ(ierr);
1181c068d9bbSLisandro Dalcin   minval = PetscAbsScalar(*h)*PetscRealConstant(1.01);
1182e884886fSBarry Smith   for (i=0; i<n; i++) {
1183e884886fSBarry Smith     if (PetscRealPart(u_vec[i] + *h*a_vec[i]) <= 0.0) {
1184e884886fSBarry Smith       val = PetscAbsScalar(u_vec[i]/a_vec[i]);
1185e884886fSBarry Smith       if (val < minval) minval = val;
1186e884886fSBarry Smith     }
1187e884886fSBarry Smith   }
1188e884886fSBarry Smith   ierr = VecRestoreArray(U,&u_vec);CHKERRQ(ierr);
1189e884886fSBarry Smith   ierr = VecRestoreArray(a,&a_vec);CHKERRQ(ierr);
1190b2566f29SBarry Smith   ierr = MPIU_Allreduce(&minval,&val,1,MPIU_REAL,MPIU_MIN,comm);CHKERRQ(ierr);
1191e884886fSBarry Smith   if (val <= PetscAbsScalar(*h)) {
11926712e2f1SBarry Smith     ierr = PetscInfo2(U,"Scaling back h from %g to %g\n",(double)PetscRealPart(*h),(double)(.99*val));CHKERRQ(ierr);
1193e884886fSBarry Smith     if (PetscRealPart(*h) > 0.0) *h =  0.99*val;
1194e884886fSBarry Smith     else                         *h = -0.99*val;
1195e884886fSBarry Smith   }
1196e884886fSBarry Smith   PetscFunctionReturn(0);
1197e884886fSBarry Smith }
1198