xref: /petsc/src/mat/impls/mffd/mffd.c (revision 88b4c22029ec44fe190bdcfb366f38b2676ca873)
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 
129e884886fSBarry Smith typedef PetscErrorCode (*FCN1)(void*,Vec); /* force argument to next function to not be extern C*/
13040244768SBarry Smith static PetscErrorCode  MatMFFDSetFunctioniBase_MFFD(Mat mat,FCN1 func)
131e884886fSBarry Smith {
132e884886fSBarry Smith   MatMFFD ctx = (MatMFFD)mat->data;
133e884886fSBarry Smith 
134e884886fSBarry Smith   PetscFunctionBegin;
135e884886fSBarry Smith   ctx->funcisetbase = func;
136e884886fSBarry Smith   PetscFunctionReturn(0);
137e884886fSBarry Smith }
138e884886fSBarry Smith 
139e884886fSBarry Smith typedef PetscErrorCode (*FCN2)(void*,PetscInt,Vec,PetscScalar*); /* force argument to next function to not be extern C*/
14040244768SBarry Smith static PetscErrorCode  MatMFFDSetFunctioni_MFFD(Mat mat,FCN2 funci)
141e884886fSBarry Smith {
142e884886fSBarry Smith   MatMFFD ctx = (MatMFFD)mat->data;
143e884886fSBarry Smith 
144e884886fSBarry Smith   PetscFunctionBegin;
145e884886fSBarry Smith   ctx->funci = funci;
146e884886fSBarry Smith   PetscFunctionReturn(0);
147e884886fSBarry Smith }
148e884886fSBarry Smith 
14940244768SBarry Smith static PetscErrorCode  MatMFFDResetHHistory_MFFD(Mat J)
150bc0f08ceSBarry Smith {
151bc0f08ceSBarry Smith   MatMFFD ctx = (MatMFFD)J->data;
152bc0f08ceSBarry Smith 
153bc0f08ceSBarry Smith   PetscFunctionBegin;
154bc0f08ceSBarry Smith   ctx->ncurrenth = 0;
155bc0f08ceSBarry Smith   PetscFunctionReturn(0);
156bc0f08ceSBarry Smith }
157e884886fSBarry Smith 
1581c84c290SBarry Smith /*@C
1591c84c290SBarry Smith    MatMFFDRegister - Adds a method to the MatMFFD registry.
1601c84c290SBarry Smith 
1611c84c290SBarry Smith    Not Collective
1621c84c290SBarry Smith 
1631c84c290SBarry Smith    Input Parameters:
1641c84c290SBarry Smith +  name_solver - name of a new user-defined compute-h module
1651c84c290SBarry Smith -  routine_create - routine to create method context
1661c84c290SBarry Smith 
1671c84c290SBarry Smith    Level: developer
1681c84c290SBarry Smith 
1691c84c290SBarry Smith    Notes:
1701c84c290SBarry Smith    MatMFFDRegister() may be called multiple times to add several user-defined solvers.
1711c84c290SBarry Smith 
1721c84c290SBarry Smith    Sample usage:
1731c84c290SBarry Smith .vb
174bdf89e91SBarry Smith    MatMFFDRegister("my_h",MyHCreate);
1751c84c290SBarry Smith .ve
1761c84c290SBarry Smith 
1771c84c290SBarry Smith    Then, your solver can be chosen with the procedural interface via
1781c84c290SBarry Smith $     MatMFFDSetType(mfctx,"my_h")
1791c84c290SBarry Smith    or at runtime via the option
180be7a6d03SBarry Smith $     -mat_mffd_type my_h
1811c84c290SBarry Smith 
1821c84c290SBarry Smith .keywords: MatMFFD, register
1831c84c290SBarry Smith 
1841c84c290SBarry Smith .seealso: MatMFFDRegisterAll(), MatMFFDRegisterDestroy()
1851c84c290SBarry Smith  @*/
186bdf89e91SBarry Smith PetscErrorCode  MatMFFDRegister(const char sname[],PetscErrorCode (*function)(MatMFFD))
187e884886fSBarry Smith {
188e884886fSBarry Smith   PetscErrorCode ierr;
189e884886fSBarry Smith 
190e884886fSBarry Smith   PetscFunctionBegin;
191a240a19fSJed Brown   ierr = PetscFunctionListAdd(&MatMFFDList,sname,function);CHKERRQ(ierr);
192e884886fSBarry Smith   PetscFunctionReturn(0);
193e884886fSBarry Smith }
194e884886fSBarry Smith 
195e884886fSBarry Smith /* ----------------------------------------------------------------------------------------*/
19640244768SBarry Smith static PetscErrorCode MatDestroy_MFFD(Mat mat)
197e884886fSBarry Smith {
198e884886fSBarry Smith   PetscErrorCode ierr;
199e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)mat->data;
200e884886fSBarry Smith 
201e884886fSBarry Smith   PetscFunctionBegin;
2026bf464f9SBarry Smith   ierr = VecDestroy(&ctx->w);CHKERRQ(ierr);
2036bf464f9SBarry Smith   ierr = VecDestroy(&ctx->drscale);CHKERRQ(ierr);
2046bf464f9SBarry Smith   ierr = VecDestroy(&ctx->dlscale);CHKERRQ(ierr);
2056bf464f9SBarry Smith   ierr = VecDestroy(&ctx->dshift);CHKERRQ(ierr);
206cfe22f5eSBarry Smith   if (ctx->current_f_allocated) {
2076bf464f9SBarry Smith     ierr = VecDestroy(&ctx->current_f);CHKERRQ(ierr);
208cfe22f5eSBarry Smith   }
209e884886fSBarry Smith   if (ctx->ops->destroy) {ierr = (*ctx->ops->destroy)(ctx);CHKERRQ(ierr);}
2106bf464f9SBarry Smith   ierr      = PetscHeaderDestroy(&ctx);CHKERRQ(ierr);
211bf0cc555SLisandro Dalcin   mat->data = 0;
212e884886fSBarry Smith 
213bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetBase_C",NULL);CHKERRQ(ierr);
214bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctioniBase_C",NULL);CHKERRQ(ierr);
215bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctioni_C",NULL);CHKERRQ(ierr);
216bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunction_C",NULL);CHKERRQ(ierr);
217bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctionError_C",NULL);CHKERRQ(ierr);
218bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetCheckh_C",NULL);CHKERRQ(ierr);
219bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetPeriod_C",NULL);CHKERRQ(ierr);
220bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDResetHHistory_C",NULL);CHKERRQ(ierr);
221e884886fSBarry Smith   PetscFunctionReturn(0);
222e884886fSBarry Smith }
223e884886fSBarry Smith 
224*88b4c220SStefano Zampini static PetscErrorCode MatSetUp_MFFD(Mat mat)
225*88b4c220SStefano Zampini {
226*88b4c220SStefano Zampini   PetscErrorCode ierr;
227*88b4c220SStefano Zampini 
228*88b4c220SStefano Zampini   PetscFunctionBegin;
229*88b4c220SStefano Zampini   ierr = PetscLayoutSetUp(mat->rmap);CHKERRQ(ierr);
230*88b4c220SStefano Zampini   ierr = PetscLayoutSetUp(mat->cmap);CHKERRQ(ierr);
231*88b4c220SStefano Zampini   PetscFunctionReturn(0);
232*88b4c220SStefano Zampini }
233*88b4c220SStefano Zampini 
234e884886fSBarry Smith /*
235e884886fSBarry Smith    MatMFFDView_MFFD - Views matrix-free parameters.
236e884886fSBarry Smith 
237e884886fSBarry Smith */
23840244768SBarry Smith static PetscErrorCode MatView_MFFD(Mat J,PetscViewer viewer)
239e884886fSBarry Smith {
240e884886fSBarry Smith   PetscErrorCode ierr;
241e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)J->data;
242a283635eSDmitry Karpeev   PetscBool      iascii, viewbase, viewfunction;
243a283635eSDmitry Karpeev   const char     *prefix;
244e884886fSBarry Smith 
245e884886fSBarry Smith   PetscFunctionBegin;
246251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
247e884886fSBarry Smith   if (iascii) {
248a283635eSDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"Matrix-free approximation:\n");CHKERRQ(ierr);
249a283635eSDmitry Karpeev     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
25057622a8eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"err=%g (relative error in function evaluation)\n",(double)ctx->error_rel);CHKERRQ(ierr);
2517adad957SLisandro Dalcin     if (!((PetscObject)ctx)->type_name) {
252e884886fSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"The compute h routine has not yet been set\n");CHKERRQ(ierr);
253e884886fSBarry Smith     } else {
2547adad957SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"Using %s compute h routine\n",((PetscObject)ctx)->type_name);CHKERRQ(ierr);
255e884886fSBarry Smith     }
256e884886fSBarry Smith     if (ctx->ops->view) {
257e884886fSBarry Smith       ierr = (*ctx->ops->view)(ctx,viewer);CHKERRQ(ierr);
258e884886fSBarry Smith     }
259a283635eSDmitry Karpeev     ierr = PetscObjectGetOptionsPrefix((PetscObject)J, &prefix);CHKERRQ(ierr);
260a283635eSDmitry Karpeev 
261c5929fdfSBarry Smith     ierr = PetscOptionsHasName(((PetscObject)J)->options,prefix, "-mat_mffd_view_base", &viewbase);CHKERRQ(ierr);
262a283635eSDmitry Karpeev     if (viewbase) {
263a283635eSDmitry Karpeev       ierr = PetscViewerASCIIPrintf(viewer, "Base:\n");CHKERRQ(ierr);
264a283635eSDmitry Karpeev       ierr = VecView(ctx->current_u, viewer);CHKERRQ(ierr);
265a283635eSDmitry Karpeev     }
266c5929fdfSBarry Smith     ierr = PetscOptionsHasName(((PetscObject)J)->options,prefix, "-mat_mffd_view_function", &viewfunction);CHKERRQ(ierr);
267a283635eSDmitry Karpeev     if (viewfunction) {
268a283635eSDmitry Karpeev       ierr = PetscViewerASCIIPrintf(viewer, "Function:\n");CHKERRQ(ierr);
269a283635eSDmitry Karpeev       ierr = VecView(ctx->current_f, viewer);CHKERRQ(ierr);
270a283635eSDmitry Karpeev     }
271a283635eSDmitry Karpeev     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
27211aeaf0aSBarry Smith   }
273e884886fSBarry Smith   PetscFunctionReturn(0);
274e884886fSBarry Smith }
275e884886fSBarry Smith 
276e884886fSBarry Smith /*
277e884886fSBarry Smith    MatAssemblyEnd_MFFD - Resets the ctx->ncurrenth to zero. This
278e884886fSBarry Smith    allows the user to indicate the beginning of a new linear solve by calling
279e884886fSBarry Smith    MatAssemblyXXX() on the matrix free matrix. This then allows the
2801d0fab5eSBarry Smith    MatCreateMFFD_WP() to properly compute ||U|| only the first time
281e884886fSBarry Smith    in the linear solver rather than every time.
2825a576424SJed Brown 
283cc2e6a90SBarry Smith    This function is referenced directly from MatAssemblyEnd_SNESMF(), which may be in a different shared library hence
284cc2e6a90SBarry Smith    it must be labeled as PETSC_EXTERN
285e884886fSBarry Smith */
2865a576424SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyEnd_MFFD(Mat J,MatAssemblyType mt)
287e884886fSBarry Smith {
288e884886fSBarry Smith   PetscErrorCode ierr;
289e884886fSBarry Smith   MatMFFD        j = (MatMFFD)J->data;
290e884886fSBarry Smith 
291e884886fSBarry Smith   PetscFunctionBegin;
292e884886fSBarry Smith   ierr      = MatMFFDResetHHistory(J);CHKERRQ(ierr);
293e884886fSBarry Smith   j->vshift = 0.0;
294e884886fSBarry Smith   j->vscale = 1.0;
295e884886fSBarry Smith   PetscFunctionReturn(0);
296e884886fSBarry Smith }
297e884886fSBarry Smith 
298e884886fSBarry Smith /*
299e884886fSBarry Smith   MatMult_MFFD - Default matrix-free form for Jacobian-vector product, y = F'(u)*a:
300e884886fSBarry Smith 
301e884886fSBarry Smith         y ~= (F(u + ha) - F(u))/h,
302e884886fSBarry Smith   where F = nonlinear function, as set by SNESSetFunction()
303e884886fSBarry Smith         u = current iterate
304e884886fSBarry Smith         h = difference interval
305e884886fSBarry Smith */
30640244768SBarry Smith static PetscErrorCode MatMult_MFFD(Mat mat,Vec a,Vec y)
307e884886fSBarry Smith {
308e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)mat->data;
309e884886fSBarry Smith   PetscScalar    h;
310e884886fSBarry Smith   Vec            w,U,F;
311e884886fSBarry Smith   PetscErrorCode ierr;
312ace3abfcSBarry Smith   PetscBool      zeroa;
313e884886fSBarry Smith 
314e884886fSBarry Smith   PetscFunctionBegin;
315ce94432eSBarry 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");
316e884886fSBarry Smith   /* We log matrix-free matrix-vector products separately, so that we can
317e884886fSBarry Smith      separate the performance monitoring from the cases that use conventional
318e884886fSBarry Smith      storage.  We may eventually modify event logging to associate events
319e884886fSBarry Smith      with particular objects, hence alleviating the more general problem. */
320e884886fSBarry Smith   ierr = PetscLogEventBegin(MATMFFD_Mult,a,y,0,0);CHKERRQ(ierr);
321e884886fSBarry Smith 
322e884886fSBarry Smith   w = ctx->w;
323e884886fSBarry Smith   U = ctx->current_u;
3243ec795f1SBarry Smith   F = ctx->current_f;
325e884886fSBarry Smith   /*
326e884886fSBarry Smith       Compute differencing parameter
327e884886fSBarry Smith   */
3284a2f8832SBarry Smith   if (!((PetscObject)ctx)->type_name) {
329e884886fSBarry Smith     ierr = MatMFFDSetType(mat,MATMFFD_WP);CHKERRQ(ierr);
3309c6ac3b3SBarry Smith     ierr = MatSetFromOptions(mat);CHKERRQ(ierr);
331e884886fSBarry Smith   }
332e884886fSBarry Smith   ierr = (*ctx->ops->compute)(ctx,U,a,&h,&zeroa);CHKERRQ(ierr);
333e884886fSBarry Smith   if (zeroa) {
334e884886fSBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
335e884886fSBarry Smith     PetscFunctionReturn(0);
336e884886fSBarry Smith   }
337e884886fSBarry Smith 
33884d44b13SHong Zhang   if (mat->erroriffailure && PetscIsInfOrNanScalar(h)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Computed Nan differencing parameter h");
339e884886fSBarry Smith   if (ctx->checkh) {
340e884886fSBarry Smith     ierr = (*ctx->checkh)(ctx->checkhctx,U,a,&h);CHKERRQ(ierr);
341e884886fSBarry Smith   }
342e884886fSBarry Smith 
343e884886fSBarry Smith   /* keep a record of the current differencing parameter h */
344e884886fSBarry Smith   ctx->currenth = h;
345e884886fSBarry Smith #if defined(PETSC_USE_COMPLEX)
34657622a8eSBarry Smith   ierr = PetscInfo2(mat,"Current differencing parameter: %g + %g i\n",(double)PetscRealPart(h),(double)PetscImaginaryPart(h));CHKERRQ(ierr);
347e884886fSBarry Smith #else
348e884886fSBarry Smith   ierr = PetscInfo1(mat,"Current differencing parameter: %15.12e\n",h);CHKERRQ(ierr);
349e884886fSBarry Smith #endif
350e884886fSBarry Smith   if (ctx->historyh && ctx->ncurrenth < ctx->maxcurrenth) {
351e884886fSBarry Smith     ctx->historyh[ctx->ncurrenth] = h;
352e884886fSBarry Smith   }
353e884886fSBarry Smith   ctx->ncurrenth++;
354e884886fSBarry Smith 
355e884886fSBarry Smith   /* w = u + ha */
356c3bb7e23SBarry Smith   if (ctx->drscale) {
357c3bb7e23SBarry Smith     ierr = VecPointwiseMult(ctx->drscale,a,U);CHKERRQ(ierr);
358c3bb7e23SBarry Smith     ierr = VecAYPX(U,h,w);CHKERRQ(ierr);
359c3bb7e23SBarry Smith   } else {
360e884886fSBarry Smith     ierr = VecWAXPY(w,h,a,U);CHKERRQ(ierr);
361c3bb7e23SBarry Smith   }
362e884886fSBarry Smith 
3631b797266SDmitry Karpeev   /* compute func(U) as base for differencing; only needed first time in and not when provided by user */
3641b797266SDmitry Karpeev   if (ctx->ncurrenth == 1 && ctx->current_f_allocated) {
365e884886fSBarry Smith     ierr = (*ctx->func)(ctx->funcctx,U,F);CHKERRQ(ierr);
36652121784SBarry Smith   }
367e884886fSBarry Smith   ierr = (*ctx->func)(ctx->funcctx,w,y);CHKERRQ(ierr);
368e884886fSBarry Smith 
369e884886fSBarry Smith   ierr = VecAXPY(y,-1.0,F);CHKERRQ(ierr);
370e884886fSBarry Smith   ierr = VecScale(y,1.0/h);CHKERRQ(ierr);
371e884886fSBarry Smith 
372c35ec66cSMatthew G Knepley   if ((ctx->vshift != 0.0) || (ctx->vscale != 1.0)) {
373e884886fSBarry Smith     ierr = VecAXPBY(y,ctx->vshift,ctx->vscale,a);CHKERRQ(ierr);
374c3bb7e23SBarry Smith   }
375c3bb7e23SBarry Smith   if (ctx->dlscale) {
376c3bb7e23SBarry Smith     ierr = VecPointwiseMult(y,ctx->dlscale,y);CHKERRQ(ierr);
377c3bb7e23SBarry Smith   }
378c3bb7e23SBarry Smith   if (ctx->dshift) {
379c3bb7e23SBarry Smith     ierr = VecPointwiseMult(ctx->dshift,a,U);CHKERRQ(ierr);
380c3bb7e23SBarry Smith     ierr = VecAXPY(y,1.0,U);CHKERRQ(ierr);
381c3bb7e23SBarry Smith   }
382e884886fSBarry Smith 
38339601f4eSBarry Smith   if (mat->nullsp) {ierr = MatNullSpaceRemove(mat->nullsp,y);CHKERRQ(ierr);}
384e884886fSBarry Smith 
385e884886fSBarry Smith   ierr = PetscLogEventEnd(MATMFFD_Mult,a,y,0,0);CHKERRQ(ierr);
386e884886fSBarry Smith   PetscFunctionReturn(0);
387e884886fSBarry Smith }
388e884886fSBarry Smith 
389e884886fSBarry Smith /*
390e884886fSBarry Smith   MatGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix
391e884886fSBarry Smith 
392e884886fSBarry Smith         y ~= (F(u + ha) - F(u))/h,
393e884886fSBarry Smith   where F = nonlinear function, as set by SNESSetFunction()
394e884886fSBarry Smith         u = current iterate
395e884886fSBarry Smith         h = difference interval
396e884886fSBarry Smith */
39740244768SBarry Smith static PetscErrorCode MatGetDiagonal_MFFD(Mat mat,Vec a)
398e884886fSBarry Smith {
399e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)mat->data;
400e884886fSBarry Smith   PetscScalar    h,*aa,*ww,v;
401e884886fSBarry Smith   PetscReal      epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON;
402e884886fSBarry Smith   Vec            w,U;
403e884886fSBarry Smith   PetscErrorCode ierr;
404e884886fSBarry Smith   PetscInt       i,rstart,rend;
405e884886fSBarry Smith 
406e884886fSBarry Smith   PetscFunctionBegin;
407e7e72b3dSBarry Smith   if (!ctx->funci) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Requires calling MatMFFDSetFunctioni() first");
408e884886fSBarry Smith 
409e884886fSBarry Smith   w    = ctx->w;
410e884886fSBarry Smith   U    = ctx->current_u;
411e884886fSBarry Smith   ierr = (*ctx->func)(ctx->funcctx,U,a);CHKERRQ(ierr);
412e884886fSBarry Smith   ierr = (*ctx->funcisetbase)(ctx->funcctx,U);CHKERRQ(ierr);
413e884886fSBarry Smith   ierr = VecCopy(U,w);CHKERRQ(ierr);
414e884886fSBarry Smith 
415e884886fSBarry Smith   ierr = VecGetOwnershipRange(a,&rstart,&rend);CHKERRQ(ierr);
416e884886fSBarry Smith   ierr = VecGetArray(a,&aa);CHKERRQ(ierr);
417e884886fSBarry Smith   for (i=rstart; i<rend; i++) {
418e884886fSBarry Smith     ierr = VecGetArray(w,&ww);CHKERRQ(ierr);
419e884886fSBarry Smith     h    = ww[i-rstart];
420e884886fSBarry Smith     if (h == 0.0) h = 1.0;
421e884886fSBarry Smith     if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0)     h = umin;
422e884886fSBarry Smith     else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin;
423e884886fSBarry Smith     h *= epsilon;
424e884886fSBarry Smith 
425e884886fSBarry Smith     ww[i-rstart] += h;
426e884886fSBarry Smith     ierr          = VecRestoreArray(w,&ww);CHKERRQ(ierr);
427e884886fSBarry Smith     ierr          = (*ctx->funci)(ctx->funcctx,i,w,&v);CHKERRQ(ierr);
428e884886fSBarry Smith     aa[i-rstart]  = (v - aa[i-rstart])/h;
429e884886fSBarry Smith 
430e884886fSBarry Smith     /* possibly shift and scale result */
431c35ec66cSMatthew G Knepley     if ((ctx->vshift != 0.0) || (ctx->vscale != 1.0)) {
432e884886fSBarry Smith       aa[i - rstart] = ctx->vshift + ctx->vscale*aa[i-rstart];
433c3bb7e23SBarry Smith     }
434e884886fSBarry Smith 
435e884886fSBarry Smith     ierr          = VecGetArray(w,&ww);CHKERRQ(ierr);
436e884886fSBarry Smith     ww[i-rstart] -= h;
437e884886fSBarry Smith     ierr          = VecRestoreArray(w,&ww);CHKERRQ(ierr);
438e884886fSBarry Smith   }
439e884886fSBarry Smith   ierr = VecRestoreArray(a,&aa);CHKERRQ(ierr);
440e884886fSBarry Smith   PetscFunctionReturn(0);
441e884886fSBarry Smith }
442e884886fSBarry Smith 
44340244768SBarry Smith static PetscErrorCode MatDiagonalScale_MFFD(Mat mat,Vec ll,Vec rr)
444c3bb7e23SBarry Smith {
445c3bb7e23SBarry Smith   MatMFFD        aij = (MatMFFD)mat->data;
446c3bb7e23SBarry Smith   PetscErrorCode ierr;
447c3bb7e23SBarry Smith 
448c3bb7e23SBarry Smith   PetscFunctionBegin;
449c3bb7e23SBarry Smith   if (ll && !aij->dlscale) {
450c3bb7e23SBarry Smith     ierr = VecDuplicate(ll,&aij->dlscale);CHKERRQ(ierr);
451c3bb7e23SBarry Smith   }
452c3bb7e23SBarry Smith   if (rr && !aij->drscale) {
453c3bb7e23SBarry Smith     ierr = VecDuplicate(rr,&aij->drscale);CHKERRQ(ierr);
454c3bb7e23SBarry Smith   }
455c3bb7e23SBarry Smith   if (ll) {
456c3bb7e23SBarry Smith     ierr = VecCopy(ll,aij->dlscale);CHKERRQ(ierr);
457c3bb7e23SBarry Smith   }
458c3bb7e23SBarry Smith   if (rr) {
459c3bb7e23SBarry Smith     ierr = VecCopy(rr,aij->drscale);CHKERRQ(ierr);
460c3bb7e23SBarry Smith   }
461c3bb7e23SBarry Smith   PetscFunctionReturn(0);
462c3bb7e23SBarry Smith }
463c3bb7e23SBarry Smith 
46440244768SBarry Smith static PetscErrorCode MatDiagonalSet_MFFD(Mat mat,Vec ll,InsertMode mode)
465c3bb7e23SBarry Smith {
466c3bb7e23SBarry Smith   MatMFFD        aij = (MatMFFD)mat->data;
467c3bb7e23SBarry Smith   PetscErrorCode ierr;
468c3bb7e23SBarry Smith 
469c3bb7e23SBarry Smith   PetscFunctionBegin;
470ce94432eSBarry Smith   if (mode == INSERT_VALUES) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"No diagonal set with INSERT_VALUES");
471c3bb7e23SBarry Smith   if (!aij->dshift) {
472c3bb7e23SBarry Smith     ierr = VecDuplicate(ll,&aij->dshift);CHKERRQ(ierr);
473c3bb7e23SBarry Smith   }
474c3bb7e23SBarry Smith   ierr = VecAXPY(aij->dshift,1.0,ll);CHKERRQ(ierr);
475c3bb7e23SBarry Smith   PetscFunctionReturn(0);
476c3bb7e23SBarry Smith }
477c3bb7e23SBarry Smith 
47840244768SBarry Smith static PetscErrorCode MatShift_MFFD(Mat Y,PetscScalar a)
479e884886fSBarry Smith {
480e884886fSBarry Smith   MatMFFD shell = (MatMFFD)Y->data;
4815fd66863SKarl Rupp 
482e884886fSBarry Smith   PetscFunctionBegin;
483e884886fSBarry Smith   shell->vshift += a;
484e884886fSBarry Smith   PetscFunctionReturn(0);
485e884886fSBarry Smith }
486e884886fSBarry Smith 
48740244768SBarry Smith static PetscErrorCode MatScale_MFFD(Mat Y,PetscScalar a)
488e884886fSBarry Smith {
489e884886fSBarry Smith   MatMFFD shell = (MatMFFD)Y->data;
4905fd66863SKarl Rupp 
491e884886fSBarry Smith   PetscFunctionBegin;
492e884886fSBarry Smith   shell->vscale *= a;
493e884886fSBarry Smith   PetscFunctionReturn(0);
494e884886fSBarry Smith }
495e884886fSBarry Smith 
496d8fa6a54SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDSetBase_MFFD(Mat J,Vec U,Vec F)
497e884886fSBarry Smith {
498e884886fSBarry Smith   PetscErrorCode ierr;
499e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)J->data;
500e884886fSBarry Smith 
501e884886fSBarry Smith   PetscFunctionBegin;
502e884886fSBarry Smith   ierr = MatMFFDResetHHistory(J);CHKERRQ(ierr);
5032205254eSKarl Rupp 
504e884886fSBarry Smith   ctx->current_u = U;
50552121784SBarry Smith   if (F) {
5066bf464f9SBarry Smith     if (ctx->current_f_allocated) {ierr = VecDestroy(&ctx->current_f);CHKERRQ(ierr);}
5073ec795f1SBarry Smith     ctx->current_f           = F;
508cfe22f5eSBarry Smith     ctx->current_f_allocated = PETSC_FALSE;
509cfe22f5eSBarry Smith   } else if (!ctx->current_f_allocated) {
51006c4b6baSBarry Smith     ierr = MatCreateVecs(J,NULL,&ctx->current_f);CHKERRQ(ierr);
5112205254eSKarl Rupp 
512cfe22f5eSBarry Smith     ctx->current_f_allocated = PETSC_TRUE;
51352121784SBarry Smith   }
514e884886fSBarry Smith   if (!ctx->w) {
515e884886fSBarry Smith     ierr = VecDuplicate(ctx->current_u, &ctx->w);CHKERRQ(ierr);
516e884886fSBarry Smith   }
517e884886fSBarry Smith   J->assembled = PETSC_TRUE;
518e884886fSBarry Smith   PetscFunctionReturn(0);
519e884886fSBarry Smith }
5205a576424SJed Brown 
521e884886fSBarry Smith typedef PetscErrorCode (*FCN3)(void*,Vec,Vec,PetscScalar*); /* force argument to next function to not be extern C*/
522b2573a8aSBarry Smith 
52340244768SBarry Smith static PetscErrorCode  MatMFFDSetCheckh_MFFD(Mat J,FCN3 fun,void *ectx)
524e884886fSBarry Smith {
525e884886fSBarry Smith   MatMFFD ctx = (MatMFFD)J->data;
526e884886fSBarry Smith 
527e884886fSBarry Smith   PetscFunctionBegin;
528e884886fSBarry Smith   ctx->checkh    = fun;
529e884886fSBarry Smith   ctx->checkhctx = ectx;
530e884886fSBarry Smith   PetscFunctionReturn(0);
531e884886fSBarry Smith }
532e884886fSBarry Smith 
5336aa9148fSLisandro Dalcin /*@C
5346aa9148fSLisandro Dalcin    MatMFFDSetOptionsPrefix - Sets the prefix used for searching for all
5356aa9148fSLisandro Dalcin    MatMFFD options in the database.
5366aa9148fSLisandro Dalcin 
5376aa9148fSLisandro Dalcin    Collective on Mat
5386aa9148fSLisandro Dalcin 
5396aa9148fSLisandro Dalcin    Input Parameter:
5406aa9148fSLisandro Dalcin +  A - the Mat context
5416aa9148fSLisandro Dalcin -  prefix - the prefix to prepend to all option names
5426aa9148fSLisandro Dalcin 
5436aa9148fSLisandro Dalcin    Notes:
5446aa9148fSLisandro Dalcin    A hyphen (-) must NOT be given at the beginning of the prefix name.
5456aa9148fSLisandro Dalcin    The first character of all runtime options is AUTOMATICALLY the hyphen.
5466aa9148fSLisandro Dalcin 
5476aa9148fSLisandro Dalcin    Level: advanced
5486aa9148fSLisandro Dalcin 
5496aa9148fSLisandro Dalcin .keywords: SNES, matrix-free, parameters
5506aa9148fSLisandro Dalcin 
551755b7f64SBarry Smith .seealso: MatSetFromOptions(), MatCreateSNESMF(), MatCreateMFFD()
5526aa9148fSLisandro Dalcin @*/
5537087cfbeSBarry Smith PetscErrorCode  MatMFFDSetOptionsPrefix(Mat mat,const char prefix[])
5546aa9148fSLisandro Dalcin 
5556aa9148fSLisandro Dalcin {
5560298fd71SBarry Smith   MatMFFD        mfctx = mat ? (MatMFFD)mat->data : (MatMFFD)NULL;
5576aa9148fSLisandro Dalcin   PetscErrorCode ierr;
5585fd66863SKarl Rupp 
5596aa9148fSLisandro Dalcin   PetscFunctionBegin;
5600700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
5610700a824SBarry Smith   PetscValidHeaderSpecific(mfctx,MATMFFD_CLASSID,1);
5626aa9148fSLisandro Dalcin   ierr = PetscObjectSetOptionsPrefix((PetscObject)mfctx,prefix);CHKERRQ(ierr);
5636aa9148fSLisandro Dalcin   PetscFunctionReturn(0);
5646aa9148fSLisandro Dalcin }
5656aa9148fSLisandro Dalcin 
56640244768SBarry Smith static PetscErrorCode  MatSetFromOptions_MFFD(PetscOptionItems *PetscOptionsObject,Mat mat)
567e884886fSBarry Smith {
56871f08665SBarry Smith   MatMFFD        mfctx = (MatMFFD)mat->data;
569e884886fSBarry Smith   PetscErrorCode ierr;
570ace3abfcSBarry Smith   PetscBool      flg;
571e884886fSBarry Smith   char           ftype[256];
572e884886fSBarry Smith 
573e884886fSBarry Smith   PetscFunctionBegin;
5740700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
5750700a824SBarry Smith   PetscValidHeaderSpecific(mfctx,MATMFFD_CLASSID,1);
5763194b578SJed Brown   ierr = PetscObjectOptionsBegin((PetscObject)mfctx);CHKERRQ(ierr);
577a264d7a6SBarry Smith   ierr = PetscOptionsFList("-mat_mffd_type","Matrix free type","MatMFFDSetType",MatMFFDList,((PetscObject)mfctx)->type_name,ftype,256,&flg);CHKERRQ(ierr);
578e884886fSBarry Smith   if (flg) {
579e884886fSBarry Smith     ierr = MatMFFDSetType(mat,ftype);CHKERRQ(ierr);
580e884886fSBarry Smith   }
581e884886fSBarry Smith 
582e884886fSBarry Smith   ierr = PetscOptionsReal("-mat_mffd_err","set sqrt relative error in function","MatMFFDSetFunctionError",mfctx->error_rel,&mfctx->error_rel,0);CHKERRQ(ierr);
583e884886fSBarry Smith   ierr = PetscOptionsInt("-mat_mffd_period","how often h is recomputed","MatMFFDSetPeriod",mfctx->recomputeperiod,&mfctx->recomputeperiod,0);CHKERRQ(ierr);
584e884886fSBarry Smith 
58590d69ab7SBarry Smith   flg  = PETSC_FALSE;
5860298fd71SBarry Smith   ierr = PetscOptionsBool("-mat_mffd_check_positivity","Insure that U + h*a is nonnegative","MatMFFDSetCheckh",flg,&flg,NULL);CHKERRQ(ierr);
587e884886fSBarry Smith   if (flg) {
588e884886fSBarry Smith     ierr = MatMFFDSetCheckh(mat,MatMFFDCheckPositivity,0);CHKERRQ(ierr);
589e884886fSBarry Smith   }
590e884886fSBarry Smith   if (mfctx->ops->setfromoptions) {
591e55864a3SBarry Smith     ierr = (*mfctx->ops->setfromoptions)(PetscOptionsObject,mfctx);CHKERRQ(ierr);
592e884886fSBarry Smith   }
593e884886fSBarry Smith   ierr = PetscOptionsEnd();CHKERRQ(ierr);
594e884886fSBarry Smith   PetscFunctionReturn(0);
595e884886fSBarry Smith }
596e884886fSBarry Smith 
59740244768SBarry Smith static PetscErrorCode  MatMFFDSetPeriod_MFFD(Mat mat,PetscInt period)
598bc0f08ceSBarry Smith {
599bc0f08ceSBarry Smith   MatMFFD ctx = (MatMFFD)mat->data;
600bc0f08ceSBarry Smith 
601bc0f08ceSBarry Smith   PetscFunctionBegin;
602bc0f08ceSBarry Smith   ctx->recomputeperiod = period;
603bc0f08ceSBarry Smith   PetscFunctionReturn(0);
604bc0f08ceSBarry Smith }
605bc0f08ceSBarry Smith 
60640244768SBarry Smith static PetscErrorCode  MatMFFDSetFunction_MFFD(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx)
607bc0f08ceSBarry Smith {
608bc0f08ceSBarry Smith   MatMFFD ctx = (MatMFFD)mat->data;
609bc0f08ceSBarry Smith 
610bc0f08ceSBarry Smith   PetscFunctionBegin;
611bc0f08ceSBarry Smith   ctx->func    = func;
612bc0f08ceSBarry Smith   ctx->funcctx = funcctx;
613bc0f08ceSBarry Smith   PetscFunctionReturn(0);
614bc0f08ceSBarry Smith }
615bc0f08ceSBarry Smith 
61640244768SBarry Smith static PetscErrorCode  MatMFFDSetFunctionError_MFFD(Mat mat,PetscReal error)
617bc0f08ceSBarry Smith {
618bc0f08ceSBarry Smith   MatMFFD ctx = (MatMFFD)mat->data;
619bc0f08ceSBarry Smith 
620bc0f08ceSBarry Smith   PetscFunctionBegin;
621bc0f08ceSBarry Smith   if (error != PETSC_DEFAULT) ctx->error_rel = error;
622bc0f08ceSBarry Smith   PetscFunctionReturn(0);
623bc0f08ceSBarry Smith }
624bc0f08ceSBarry Smith 
6253b49f96aSBarry Smith static PetscErrorCode MatMissingDiagonal_MFFD(Mat A,PetscBool  *missing,PetscInt *d)
6263b49f96aSBarry Smith {
6273b49f96aSBarry Smith   PetscFunctionBegin;
6283b49f96aSBarry Smith   *missing = PETSC_FALSE;
6293b49f96aSBarry Smith   PetscFunctionReturn(0);
6303b49f96aSBarry Smith }
6313b49f96aSBarry Smith 
632e884886fSBarry Smith /*MC
633e884886fSBarry Smith   MATMFFD - MATMFFD = "mffd" - A matrix free matrix type.
634e884886fSBarry Smith 
635e884886fSBarry Smith   Level: advanced
636e884886fSBarry Smith 
637755b7f64SBarry Smith .seealso: MatCreateMFFD(), MatCreateSNESMF(), MatMFFDSetFunction(), MatMFFDSetType(),
638755b7f64SBarry Smith           MatMFFDSetFunctionError(), MatMFFDDSSetUmin(), MatMFFDSetFunction()
639755b7f64SBarry Smith           MatMFFDSetHHistory(), MatMFFDResetHHistory(), MatCreateSNESMF(),
640755b7f64SBarry Smith           MatMFFDGetH(),
641e884886fSBarry Smith M*/
6428cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_MFFD(Mat A)
643e884886fSBarry Smith {
644e884886fSBarry Smith   MatMFFD        mfctx;
645e884886fSBarry Smith   PetscErrorCode ierr;
646e884886fSBarry Smith 
647e884886fSBarry Smith   PetscFunctionBegin;
648607a6623SBarry Smith   ierr = MatMFFDInitializePackage();CHKERRQ(ierr);
649e884886fSBarry Smith 
65073107ff1SLisandro Dalcin   ierr = PetscHeaderCreate(mfctx,MATMFFD_CLASSID,"MatMFFD","Matrix-free Finite Differencing","Mat",PetscObjectComm((PetscObject)A),MatDestroy_MFFD,MatView_MFFD);CHKERRQ(ierr);
6512205254eSKarl Rupp 
652e884886fSBarry Smith   mfctx->error_rel                = PETSC_SQRT_MACHINE_EPSILON;
653e884886fSBarry Smith   mfctx->recomputeperiod          = 1;
654e884886fSBarry Smith   mfctx->count                    = 0;
655e884886fSBarry Smith   mfctx->currenth                 = 0.0;
6560298fd71SBarry Smith   mfctx->historyh                 = NULL;
657e884886fSBarry Smith   mfctx->ncurrenth                = 0;
658e884886fSBarry Smith   mfctx->maxcurrenth              = 0;
6597adad957SLisandro Dalcin   ((PetscObject)mfctx)->type_name = 0;
660e884886fSBarry Smith 
661e884886fSBarry Smith   mfctx->vshift = 0.0;
662e884886fSBarry Smith   mfctx->vscale = 1.0;
663e884886fSBarry Smith 
664e884886fSBarry Smith   /*
665e884886fSBarry Smith      Create the empty data structure to contain compute-h routines.
666e884886fSBarry Smith      These will be filled in below from the command line options or
667e884886fSBarry Smith      a later call with MatMFFDSetType() or if that is not called
668e884886fSBarry Smith      then it will default in the first use of MatMult_MFFD()
669e884886fSBarry Smith   */
670e884886fSBarry Smith   mfctx->ops->compute        = 0;
671e884886fSBarry Smith   mfctx->ops->destroy        = 0;
672e884886fSBarry Smith   mfctx->ops->view           = 0;
673e884886fSBarry Smith   mfctx->ops->setfromoptions = 0;
674e884886fSBarry Smith   mfctx->hctx                = 0;
675e884886fSBarry Smith 
676e884886fSBarry Smith   mfctx->func    = 0;
677e884886fSBarry Smith   mfctx->funcctx = 0;
6780298fd71SBarry Smith   mfctx->w       = NULL;
679e884886fSBarry Smith 
680e884886fSBarry Smith   A->data = mfctx;
681e884886fSBarry Smith 
682e884886fSBarry Smith   A->ops->mult            = MatMult_MFFD;
683e884886fSBarry Smith   A->ops->destroy         = MatDestroy_MFFD;
684*88b4c220SStefano Zampini   A->ops->setup           = MatSetUp_MFFD;
685e884886fSBarry Smith   A->ops->view            = MatView_MFFD;
686e884886fSBarry Smith   A->ops->assemblyend     = MatAssemblyEnd_MFFD;
687e884886fSBarry Smith   A->ops->getdiagonal     = MatGetDiagonal_MFFD;
688e884886fSBarry Smith   A->ops->scale           = MatScale_MFFD;
689e884886fSBarry Smith   A->ops->shift           = MatShift_MFFD;
690c3bb7e23SBarry Smith   A->ops->diagonalscale   = MatDiagonalScale_MFFD;
691c3bb7e23SBarry Smith   A->ops->diagonalset     = MatDiagonalSet_MFFD;
6929c6ac3b3SBarry Smith   A->ops->setfromoptions  = MatSetFromOptions_MFFD;
6933b49f96aSBarry Smith   A->ops->missingdiagonal = MatMissingDiagonal_MFFD;
694e884886fSBarry Smith   A->assembled            = PETSC_TRUE;
695e884886fSBarry Smith 
696bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetBase_C",MatMFFDSetBase_MFFD);CHKERRQ(ierr);
697bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctioniBase_C",MatMFFDSetFunctioniBase_MFFD);CHKERRQ(ierr);
698bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctioni_C",MatMFFDSetFunctioni_MFFD);CHKERRQ(ierr);
699bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunction_C",MatMFFDSetFunction_MFFD);CHKERRQ(ierr);
700bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetCheckh_C",MatMFFDSetCheckh_MFFD);CHKERRQ(ierr);
701bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetPeriod_C",MatMFFDSetPeriod_MFFD);CHKERRQ(ierr);
702bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctionError_C",MatMFFDSetFunctionError_MFFD);CHKERRQ(ierr);
703bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDResetHHistory_C",MatMFFDResetHHistory_MFFD);CHKERRQ(ierr);
7042205254eSKarl Rupp 
705e884886fSBarry Smith   mfctx->mat = A;
7062205254eSKarl Rupp 
707e884886fSBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,MATMFFD);CHKERRQ(ierr);
708e884886fSBarry Smith   PetscFunctionReturn(0);
709e884886fSBarry Smith }
710e884886fSBarry Smith 
711e884886fSBarry Smith /*@
712e884886fSBarry Smith    MatCreateMFFD - Creates a matrix-free matrix. See also MatCreateSNESMF()
713e884886fSBarry Smith 
714e884886fSBarry Smith    Collective on Vec
715e884886fSBarry Smith 
716e884886fSBarry Smith    Input Parameters:
717fef1beadSBarry Smith +  comm - MPI communicator
718fef1beadSBarry Smith .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
719fef1beadSBarry Smith            This value should be the same as the local size used in creating the
720fef1beadSBarry Smith            y vector for the matrix-vector product y = Ax.
721fef1beadSBarry Smith .  n - This value should be the same as the local size used in creating the
722fef1beadSBarry Smith        x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have
723fef1beadSBarry Smith        calculated if N is given) For square matrices n is almost always m.
724fef1beadSBarry Smith .  M - number of global rows (or PETSC_DETERMINE to have calculated if m is given)
725fef1beadSBarry Smith -  N - number of global columns (or PETSC_DETERMINE to have calculated if n is given)
726fef1beadSBarry Smith 
727e884886fSBarry Smith 
728e884886fSBarry Smith    Output Parameter:
729e884886fSBarry Smith .  J - the matrix-free matrix
730e884886fSBarry Smith 
7319c6ac3b3SBarry Smith    Options Database Keys: call MatSetFromOptions() to trigger these
7329c6ac3b3SBarry Smith +  -mat_mffd_type - wp or ds (see MATMFFD_WP or MATMFFD_DS)
7339c6ac3b3SBarry Smith -  -mat_mffd_err - square root of estimated relative error in function evaluation
7349c6ac3b3SBarry Smith -  -mat_mffd_period - how often h is recomputed, defaults to 1, everytime
7359c6ac3b3SBarry Smith 
7369c6ac3b3SBarry Smith 
737e884886fSBarry Smith    Level: advanced
738e884886fSBarry Smith 
739e884886fSBarry Smith    Notes:
740e884886fSBarry Smith    The matrix-free matrix context merely contains the function pointers
741e884886fSBarry Smith    and work space for performing finite difference approximations of
742e884886fSBarry Smith    Jacobian-vector products, F'(u)*a,
743e884886fSBarry Smith 
744e884886fSBarry Smith    The default code uses the following approach to compute h
745e884886fSBarry Smith 
746e884886fSBarry Smith .vb
747e884886fSBarry Smith      F'(u)*a = [F(u+h*a) - F(u)]/h where
748e884886fSBarry Smith      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
749e884886fSBarry Smith        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
750e884886fSBarry Smith  where
751e884886fSBarry Smith      error_rel = square root of relative error in function evaluation
752e884886fSBarry Smith      umin = minimum iterate parameter
753e884886fSBarry Smith .ve
754e884886fSBarry Smith 
755ca93e954SBarry Smith    You can call SNESSetJacobian() with MatMFFDComputeJacobian() if you are using matrix and not a different
756ca93e954SBarry Smith    preconditioner matrix
757ca93e954SBarry Smith 
758e884886fSBarry Smith    The user can set the error_rel via MatMFFDSetFunctionError() and
759a7f22e61SSatish Balay    umin via MatMFFDDSSetUmin(); see Users-Manual: ch_snes for details.
760e884886fSBarry Smith 
761e884886fSBarry Smith    The user should call MatDestroy() when finished with the matrix-free
762e884886fSBarry Smith    matrix context.
763e884886fSBarry Smith 
764e884886fSBarry Smith    Options Database Keys:
765e884886fSBarry Smith +  -mat_mffd_err <error_rel> - Sets error_rel
766e884886fSBarry Smith .  -mat_mffd_unim <umin> - Sets umin (for default PETSc routine that computes h only)
767e884886fSBarry Smith -  -mat_mffd_check_positivity
768e884886fSBarry Smith 
769e884886fSBarry Smith .keywords: default, matrix-free, create, matrix
770e884886fSBarry Smith 
7711d0fab5eSBarry Smith .seealso: MatDestroy(), MatMFFDSetFunctionError(), MatMFFDDSSetUmin(), MatMFFDSetFunction()
772e884886fSBarry Smith           MatMFFDSetHHistory(), MatMFFDResetHHistory(), MatCreateSNESMF(),
77381242352SJed Brown           MatMFFDGetH(), MatMFFDRegister(), MatMFFDComputeJacobian()
774e884886fSBarry Smith 
775e884886fSBarry Smith @*/
7767087cfbeSBarry Smith PetscErrorCode  MatCreateMFFD(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *J)
777e884886fSBarry Smith {
778e884886fSBarry Smith   PetscErrorCode ierr;
779e884886fSBarry Smith 
780e884886fSBarry Smith   PetscFunctionBegin;
781e884886fSBarry Smith   ierr = MatCreate(comm,J);CHKERRQ(ierr);
782fef1beadSBarry Smith   ierr = MatSetSizes(*J,m,n,M,N);CHKERRQ(ierr);
783e884886fSBarry Smith   ierr = MatSetType(*J,MATMFFD);CHKERRQ(ierr);
784be7c243fSJed Brown   ierr = MatSetUp(*J);CHKERRQ(ierr);
785e884886fSBarry Smith   PetscFunctionReturn(0);
786e884886fSBarry Smith }
787e884886fSBarry Smith 
788e884886fSBarry Smith /*@
789e884886fSBarry Smith    MatMFFDGetH - Gets the last value that was used as the differencing
790e884886fSBarry Smith    parameter.
791e884886fSBarry Smith 
792e884886fSBarry Smith    Not Collective
793e884886fSBarry Smith 
794e884886fSBarry Smith    Input Parameters:
795e884886fSBarry Smith .  mat - the matrix obtained with MatCreateSNESMF()
796e884886fSBarry Smith 
797e884886fSBarry Smith    Output Paramter:
798e884886fSBarry Smith .  h - the differencing step size
799e884886fSBarry Smith 
800e884886fSBarry Smith    Level: advanced
801e884886fSBarry Smith 
802e884886fSBarry Smith .keywords: SNES, matrix-free, parameters
803e884886fSBarry Smith 
8041d0fab5eSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDSetHHistory(), MatCreateMFFD(), MATMFFD, MatMFFDResetHHistory()
805e884886fSBarry Smith @*/
8067087cfbeSBarry Smith PetscErrorCode  MatMFFDGetH(Mat mat,PetscScalar *h)
807e884886fSBarry Smith {
808e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)mat->data;
809bc0f08ceSBarry Smith   PetscErrorCode ierr;
810bc0f08ceSBarry Smith   PetscBool      match;
811e884886fSBarry Smith 
812e884886fSBarry Smith   PetscFunctionBegin;
813*88b4c220SStefano Zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
814*88b4c220SStefano Zampini   PetscValidPointer(h,2);
815251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)mat,MATMFFD,&match);CHKERRQ(ierr);
816ce94432eSBarry Smith   if (!match) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONG,"Not a MFFD matrix");
817bc0f08ceSBarry Smith 
818e884886fSBarry Smith   *h = ctx->currenth;
819e884886fSBarry Smith   PetscFunctionReturn(0);
820e884886fSBarry Smith }
821e884886fSBarry Smith 
822e884886fSBarry Smith /*@C
823e884886fSBarry Smith    MatMFFDSetFunction - Sets the function used in applying the matrix free.
824e884886fSBarry Smith 
8253f9fe445SBarry Smith    Logically Collective on Mat
826e884886fSBarry Smith 
827e884886fSBarry Smith    Input Parameters:
82814f633e2SBarry Smith +  mat - the matrix free matrix created via MatCreateSNESMF() or MatCreateMFFD()
829e884886fSBarry Smith .  func - the function to use
830e884886fSBarry Smith -  funcctx - optional function context passed to function
831e884886fSBarry Smith 
83214f633e2SBarry Smith    Calling Sequence of func:
83314f633e2SBarry Smith $     func (void *funcctx, Vec x, Vec f)
83414f633e2SBarry Smith 
83514f633e2SBarry Smith +  funcctx - user provided context
83614f633e2SBarry Smith .  x - input vector
83714f633e2SBarry Smith -  f - computed output function
83814f633e2SBarry Smith 
839e884886fSBarry Smith    Level: advanced
840e884886fSBarry Smith 
841e884886fSBarry Smith    Notes:
842e884886fSBarry Smith     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
843e884886fSBarry Smith     matrix inside your compute Jacobian routine
844e884886fSBarry Smith 
8453ec795f1SBarry Smith     If this is not set then it will use the function set with SNESSetFunction() if MatCreateSNESMF() was used.
846e884886fSBarry Smith 
847e884886fSBarry Smith .keywords: SNES, matrix-free, function
848e884886fSBarry Smith 
8491d0fab5eSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD,
8501d0fab5eSBarry Smith           MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction()
851e884886fSBarry Smith @*/
8527087cfbeSBarry Smith PetscErrorCode  MatMFFDSetFunction(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx)
853e884886fSBarry Smith {
854bc0f08ceSBarry Smith   PetscErrorCode ierr;
855e884886fSBarry Smith 
856e884886fSBarry Smith   PetscFunctionBegin;
857*88b4c220SStefano Zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
858bc0f08ceSBarry Smith   ierr = PetscTryMethod(mat,"MatMFFDSetFunction_C",(Mat,PetscErrorCode (*)(void*,Vec,Vec),void*),(mat,func,funcctx));CHKERRQ(ierr);
859e884886fSBarry Smith   PetscFunctionReturn(0);
860e884886fSBarry Smith }
861e884886fSBarry Smith 
862e884886fSBarry Smith /*@C
863e884886fSBarry Smith    MatMFFDSetFunctioni - Sets the function for a single component
864e884886fSBarry Smith 
8653f9fe445SBarry Smith    Logically Collective on Mat
866e884886fSBarry Smith 
867e884886fSBarry Smith    Input Parameters:
868e884886fSBarry Smith +  mat - the matrix free matrix created via MatCreateSNESMF()
869e884886fSBarry Smith -  funci - the function to use
870e884886fSBarry Smith 
871e884886fSBarry Smith    Level: advanced
872e884886fSBarry Smith 
873e884886fSBarry Smith    Notes:
874e884886fSBarry Smith     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
875e884886fSBarry Smith     matrix inside your compute Jacobian routine
876e884886fSBarry Smith 
877e884886fSBarry Smith 
878e884886fSBarry Smith .keywords: SNES, matrix-free, function
879e884886fSBarry Smith 
8801d0fab5eSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction()
8811d0fab5eSBarry Smith 
882e884886fSBarry Smith @*/
8837087cfbeSBarry Smith PetscErrorCode  MatMFFDSetFunctioni(Mat mat,PetscErrorCode (*funci)(void*,PetscInt,Vec,PetscScalar*))
884e884886fSBarry Smith {
8854ac538c5SBarry Smith   PetscErrorCode ierr;
886e884886fSBarry Smith 
887e884886fSBarry Smith   PetscFunctionBegin;
8880700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
8894ac538c5SBarry Smith   ierr = PetscTryMethod(mat,"MatMFFDSetFunctioni_C",(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)),(mat,funci));CHKERRQ(ierr);
890e884886fSBarry Smith   PetscFunctionReturn(0);
891e884886fSBarry Smith }
892e884886fSBarry Smith 
893e884886fSBarry Smith /*@C
894e884886fSBarry Smith    MatMFFDSetFunctioniBase - Sets the base vector for a single component function evaluation
895e884886fSBarry Smith 
8963f9fe445SBarry Smith    Logically Collective on Mat
897e884886fSBarry Smith 
898e884886fSBarry Smith    Input Parameters:
899e884886fSBarry Smith +  mat - the matrix free matrix created via MatCreateSNESMF()
900e884886fSBarry Smith -  func - the function to use
901e884886fSBarry Smith 
902e884886fSBarry Smith    Level: advanced
903e884886fSBarry Smith 
904e884886fSBarry Smith    Notes:
905e884886fSBarry Smith     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
906e884886fSBarry Smith     matrix inside your compute Jacobian routine
907e884886fSBarry Smith 
908e884886fSBarry Smith 
909e884886fSBarry Smith .keywords: SNES, matrix-free, function
910e884886fSBarry Smith 
911e884886fSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD
9121d0fab5eSBarry Smith           MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction()
913e884886fSBarry Smith @*/
9147087cfbeSBarry Smith PetscErrorCode  MatMFFDSetFunctioniBase(Mat mat,PetscErrorCode (*func)(void*,Vec))
915e884886fSBarry Smith {
9164ac538c5SBarry Smith   PetscErrorCode ierr;
917e884886fSBarry Smith 
918e884886fSBarry Smith   PetscFunctionBegin;
9190700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
9204ac538c5SBarry Smith   ierr = PetscTryMethod(mat,"MatMFFDSetFunctioniBase_C",(Mat,PetscErrorCode (*)(void*,Vec)),(mat,func));CHKERRQ(ierr);
921e884886fSBarry Smith   PetscFunctionReturn(0);
922e884886fSBarry Smith }
923e884886fSBarry Smith 
924e884886fSBarry Smith /*@
925e884886fSBarry Smith    MatMFFDSetPeriod - Sets how often h is recomputed, by default it is everytime
926e884886fSBarry Smith 
9273f9fe445SBarry Smith    Logically Collective on Mat
928e884886fSBarry Smith 
929e884886fSBarry Smith    Input Parameters:
930e884886fSBarry Smith +  mat - the matrix free matrix created via MatCreateSNESMF()
931e884886fSBarry Smith -  period - 1 for everytime, 2 for every second etc
932e884886fSBarry Smith 
933e884886fSBarry Smith    Options Database Keys:
934e884886fSBarry Smith +  -mat_mffd_period <period>
935e884886fSBarry Smith 
936e884886fSBarry Smith    Level: advanced
937e884886fSBarry Smith 
938e884886fSBarry Smith 
939e884886fSBarry Smith .keywords: SNES, matrix-free, parameters
940e884886fSBarry Smith 
941e884886fSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(),
9421d0fab5eSBarry Smith           MatMFFDSetHHistory(), MatMFFDResetHHistory()
943e884886fSBarry Smith @*/
9447087cfbeSBarry Smith PetscErrorCode  MatMFFDSetPeriod(Mat mat,PetscInt period)
945e884886fSBarry Smith {
946bc0f08ceSBarry Smith   PetscErrorCode ierr;
947e884886fSBarry Smith 
948e884886fSBarry Smith   PetscFunctionBegin;
949*88b4c220SStefano Zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
950*88b4c220SStefano Zampini   PetscValidLogicalCollectiveInt(mat,period,2);
951bc0f08ceSBarry Smith   ierr = PetscTryMethod(mat,"MatMFFDSetPeriod_C",(Mat,PetscInt),(mat,period));CHKERRQ(ierr);
952e884886fSBarry Smith   PetscFunctionReturn(0);
953e884886fSBarry Smith }
954e884886fSBarry Smith 
955e884886fSBarry Smith /*@
956e884886fSBarry Smith    MatMFFDSetFunctionError - Sets the error_rel for the approximation of
957e884886fSBarry Smith    matrix-vector products using finite differences.
958e884886fSBarry Smith 
9593f9fe445SBarry Smith    Logically Collective on Mat
960e884886fSBarry Smith 
961e884886fSBarry Smith    Input Parameters:
962e884886fSBarry Smith +  mat - the matrix free matrix created via MatCreateMFFD() or MatCreateSNESMF()
963e884886fSBarry Smith -  error_rel - relative error (should be set to the square root of
964e884886fSBarry Smith                the relative error in the function evaluations)
965e884886fSBarry Smith 
966e884886fSBarry Smith    Options Database Keys:
967e884886fSBarry Smith +  -mat_mffd_err <error_rel> - Sets error_rel
968e884886fSBarry Smith 
969e884886fSBarry Smith    Level: advanced
970e884886fSBarry Smith 
971e884886fSBarry Smith    Notes:
972e884886fSBarry Smith    The default matrix-free matrix-vector product routine computes
973e884886fSBarry Smith .vb
974e884886fSBarry Smith      F'(u)*a = [F(u+h*a) - F(u)]/h where
975e884886fSBarry Smith      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
976e884886fSBarry Smith        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   else
977e884886fSBarry Smith .ve
978e884886fSBarry Smith 
979e884886fSBarry Smith .keywords: SNES, matrix-free, parameters
980e884886fSBarry Smith 
981e884886fSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD
9821d0fab5eSBarry Smith           MatMFFDSetHHistory(), MatMFFDResetHHistory()
983e884886fSBarry Smith @*/
9847087cfbeSBarry Smith PetscErrorCode  MatMFFDSetFunctionError(Mat mat,PetscReal error)
985e884886fSBarry Smith {
986bc0f08ceSBarry Smith   PetscErrorCode ierr;
987e884886fSBarry Smith 
988e884886fSBarry Smith   PetscFunctionBegin;
989*88b4c220SStefano Zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
990*88b4c220SStefano Zampini   PetscValidLogicalCollectiveReal(mat,error,2);
991bc0f08ceSBarry Smith   ierr = PetscTryMethod(mat,"MatMFFDSetFunctionError_C",(Mat,PetscReal),(mat,error));CHKERRQ(ierr);
992e884886fSBarry Smith   PetscFunctionReturn(0);
993e884886fSBarry Smith }
994e884886fSBarry Smith 
995e884886fSBarry Smith /*@
996e884886fSBarry Smith    MatMFFDSetHHistory - Sets an array to collect a history of the
997e884886fSBarry Smith    differencing values (h) computed for the matrix-free product.
998e884886fSBarry Smith 
9993f9fe445SBarry Smith    Logically Collective on Mat
1000e884886fSBarry Smith 
1001e884886fSBarry Smith    Input Parameters:
1002e884886fSBarry Smith +  J - the matrix-free matrix context
1003e884886fSBarry Smith .  histroy - space to hold the history
1004e884886fSBarry Smith -  nhistory - number of entries in history, if more entries are generated than
1005e884886fSBarry Smith               nhistory, then the later ones are discarded
1006e884886fSBarry Smith 
1007e884886fSBarry Smith    Level: advanced
1008e884886fSBarry Smith 
1009e884886fSBarry Smith    Notes:
1010e884886fSBarry Smith    Use MatMFFDResetHHistory() to reset the history counter and collect
1011e884886fSBarry Smith    a new batch of differencing parameters, h.
1012e884886fSBarry Smith 
1013e884886fSBarry Smith .keywords: SNES, matrix-free, h history, differencing history
1014e884886fSBarry Smith 
1015e884886fSBarry Smith .seealso: MatMFFDGetH(), MatCreateSNESMF(),
10161d0fab5eSBarry Smith           MatMFFDResetHHistory(), MatMFFDSetFunctionError()
1017e884886fSBarry Smith 
1018e884886fSBarry Smith @*/
10197087cfbeSBarry Smith PetscErrorCode  MatMFFDSetHHistory(Mat J,PetscScalar history[],PetscInt nhistory)
1020e884886fSBarry Smith {
1021e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)J->data;
1022bc0f08ceSBarry Smith   PetscErrorCode ierr;
1023bc0f08ceSBarry Smith   PetscBool      match;
1024e884886fSBarry Smith 
1025e884886fSBarry Smith   PetscFunctionBegin;
1026*88b4c220SStefano Zampini   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
1027*88b4c220SStefano Zampini   if (history) PetscValidPointer(history,2);
1028*88b4c220SStefano Zampini   PetscValidLogicalCollectiveInt(J,nhistory,3);
1029251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)J,MATMFFD,&match);CHKERRQ(ierr);
1030ce94432eSBarry Smith   if (!match) SETERRQ(PetscObjectComm((PetscObject)J),PETSC_ERR_ARG_WRONG,"Not a MFFD matrix");
1031e884886fSBarry Smith   ctx->historyh    = history;
1032e884886fSBarry Smith   ctx->maxcurrenth = nhistory;
103375567043SBarry Smith   ctx->currenth    = 0.;
1034e884886fSBarry Smith   PetscFunctionReturn(0);
1035e884886fSBarry Smith }
1036e884886fSBarry Smith 
1037e884886fSBarry Smith /*@
1038e884886fSBarry Smith    MatMFFDResetHHistory - Resets the counter to zero to begin
1039e884886fSBarry Smith    collecting a new set of differencing histories.
1040e884886fSBarry Smith 
10413f9fe445SBarry Smith    Logically Collective on Mat
1042e884886fSBarry Smith 
1043e884886fSBarry Smith    Input Parameters:
1044e884886fSBarry Smith .  J - the matrix-free matrix context
1045e884886fSBarry Smith 
1046e884886fSBarry Smith    Level: advanced
1047e884886fSBarry Smith 
1048e884886fSBarry Smith    Notes:
1049e884886fSBarry Smith    Use MatMFFDSetHHistory() to create the original history counter.
1050e884886fSBarry Smith 
1051e884886fSBarry Smith .keywords: SNES, matrix-free, h history, differencing history
1052e884886fSBarry Smith 
1053e884886fSBarry Smith .seealso: MatMFFDGetH(), MatCreateSNESMF(),
10541d0fab5eSBarry Smith           MatMFFDSetHHistory(), MatMFFDSetFunctionError()
1055e884886fSBarry Smith 
1056e884886fSBarry Smith @*/
10577087cfbeSBarry Smith PetscErrorCode  MatMFFDResetHHistory(Mat J)
1058e884886fSBarry Smith {
1059bc0f08ceSBarry Smith   PetscErrorCode ierr;
1060e884886fSBarry Smith 
1061e884886fSBarry Smith   PetscFunctionBegin;
1062*88b4c220SStefano Zampini   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
1063bc0f08ceSBarry Smith   ierr = PetscTryMethod(J,"MatMFFDResetHHistory_C",(Mat),(J));CHKERRQ(ierr);
1064e884886fSBarry Smith   PetscFunctionReturn(0);
1065e884886fSBarry Smith }
1066e884886fSBarry Smith 
1067e884886fSBarry Smith /*@
1068e884886fSBarry Smith     MatMFFDSetBase - Sets the vector U at which matrix vector products of the
1069e884886fSBarry Smith         Jacobian are computed
1070e884886fSBarry Smith 
10713f9fe445SBarry Smith     Logically Collective on Mat
1072e884886fSBarry Smith 
1073e884886fSBarry Smith     Input Parameters:
1074e884886fSBarry Smith +   J - the MatMFFD matrix
10753ec795f1SBarry Smith .   U - the vector
1076bcddec3dSBarry Smith -   F - (optional) vector that contains F(u) if it has been already computed
1077e884886fSBarry Smith 
1078e884886fSBarry Smith     Notes: This is rarely used directly
1079e884886fSBarry Smith 
10808af5ae88SBarry Smith     If F is provided then it is not recomputed. Otherwise the function is evaluated at the base
10818af5ae88SBarry Smith     point during the first MatMult() after each call to MatMFFDSetBase().
1082dff2f722SBarry Smith 
1083e884886fSBarry Smith     Level: advanced
1084e884886fSBarry Smith 
1085e884886fSBarry Smith @*/
10867087cfbeSBarry Smith PetscErrorCode  MatMFFDSetBase(Mat J,Vec U,Vec F)
1087e884886fSBarry Smith {
10884ac538c5SBarry Smith   PetscErrorCode ierr;
1089e884886fSBarry Smith 
1090e884886fSBarry Smith   PetscFunctionBegin;
10910700a824SBarry Smith   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
10920700a824SBarry Smith   PetscValidHeaderSpecific(U,VEC_CLASSID,2);
10930700a824SBarry Smith   if (F) PetscValidHeaderSpecific(F,VEC_CLASSID,3);
10944ac538c5SBarry Smith   ierr = PetscTryMethod(J,"MatMFFDSetBase_C",(Mat,Vec,Vec),(J,U,F));CHKERRQ(ierr);
1095e884886fSBarry Smith   PetscFunctionReturn(0);
1096e884886fSBarry Smith }
1097e884886fSBarry Smith 
1098e884886fSBarry Smith /*@C
1099e884886fSBarry Smith     MatMFFDSetCheckh - Sets a function that checks the computed h and adjusts
1100e884886fSBarry Smith         it to satisfy some criteria
1101e884886fSBarry Smith 
11023f9fe445SBarry Smith     Logically Collective on Mat
1103e884886fSBarry Smith 
1104e884886fSBarry Smith     Input Parameters:
1105e884886fSBarry Smith +   J - the MatMFFD matrix
1106e884886fSBarry Smith .   fun - the function that checks h
1107e884886fSBarry Smith -   ctx - any context needed by the function
1108e884886fSBarry Smith 
1109e884886fSBarry Smith     Options Database Keys:
1110e884886fSBarry Smith .   -mat_mffd_check_positivity
1111e884886fSBarry Smith 
1112e884886fSBarry Smith     Level: advanced
1113e884886fSBarry Smith 
1114755b7f64SBarry Smith     Notes: For example, MatMFFDCheckPositivity() insures that all entries
1115e884886fSBarry Smith        of U + h*a are non-negative
1116e884886fSBarry Smith 
1117755b7f64SBarry Smith      The function you provide is called after the default h has been computed and allows you to
1118755b7f64SBarry Smith      modify it.
1119755b7f64SBarry Smith 
1120755b7f64SBarry Smith .seealso:  MatMFFDCheckPositivity()
1121e884886fSBarry Smith @*/
11227087cfbeSBarry Smith PetscErrorCode  MatMFFDSetCheckh(Mat J,PetscErrorCode (*fun)(void*,Vec,Vec,PetscScalar*),void *ctx)
1123e884886fSBarry Smith {
11244ac538c5SBarry Smith   PetscErrorCode ierr;
1125e884886fSBarry Smith 
1126e884886fSBarry Smith   PetscFunctionBegin;
11270700a824SBarry Smith   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
11284ac538c5SBarry Smith   ierr = PetscTryMethod(J,"MatMFFDSetCheckh_C",(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*),(J,fun,ctx));CHKERRQ(ierr);
1129e884886fSBarry Smith   PetscFunctionReturn(0);
1130e884886fSBarry Smith }
1131e884886fSBarry Smith 
1132e884886fSBarry Smith /*@
1133e884886fSBarry Smith     MatMFFDCheckPositivity - Checks that all entries in U + h*a are positive or
1134e884886fSBarry Smith         zero, decreases h until this is satisfied.
1135e884886fSBarry Smith 
11363f9fe445SBarry Smith     Logically Collective on Vec
1137e884886fSBarry Smith 
1138e884886fSBarry Smith     Input Parameters:
1139e884886fSBarry Smith +   U - base vector that is added to
1140e884886fSBarry Smith .   a - vector that is added
1141e884886fSBarry Smith .   h - scaling factor on a
1142e884886fSBarry Smith -   dummy - context variable (unused)
1143e884886fSBarry Smith 
1144e884886fSBarry Smith     Options Database Keys:
1145e884886fSBarry Smith .   -mat_mffd_check_positivity
1146e884886fSBarry Smith 
1147e884886fSBarry Smith     Level: advanced
1148e884886fSBarry Smith 
1149e884886fSBarry Smith     Notes: This is rarely used directly, rather it is passed as an argument to
1150e884886fSBarry Smith            MatMFFDSetCheckh()
1151e884886fSBarry Smith 
1152e884886fSBarry Smith .seealso:  MatMFFDSetCheckh()
1153e884886fSBarry Smith @*/
11547087cfbeSBarry Smith PetscErrorCode  MatMFFDCheckPositivity(void *dummy,Vec U,Vec a,PetscScalar *h)
1155e884886fSBarry Smith {
1156e884886fSBarry Smith   PetscReal      val, minval;
1157e884886fSBarry Smith   PetscScalar    *u_vec, *a_vec;
1158e884886fSBarry Smith   PetscErrorCode ierr;
1159e884886fSBarry Smith   PetscInt       i,n;
1160e884886fSBarry Smith   MPI_Comm       comm;
1161e884886fSBarry Smith 
1162e884886fSBarry Smith   PetscFunctionBegin;
1163*88b4c220SStefano Zampini   PetscValidHeaderSpecific(U,VEC_CLASSID,2);
1164*88b4c220SStefano Zampini   PetscValidHeaderSpecific(a,VEC_CLASSID,3);
1165*88b4c220SStefano Zampini   PetscValidPointer(h,4);
1166e884886fSBarry Smith   ierr   = PetscObjectGetComm((PetscObject)U,&comm);CHKERRQ(ierr);
1167e884886fSBarry Smith   ierr   = VecGetArray(U,&u_vec);CHKERRQ(ierr);
1168e884886fSBarry Smith   ierr   = VecGetArray(a,&a_vec);CHKERRQ(ierr);
1169e884886fSBarry Smith   ierr   = VecGetLocalSize(U,&n);CHKERRQ(ierr);
1170e884886fSBarry Smith   minval = PetscAbsScalar(*h*1.01);
1171e884886fSBarry Smith   for (i=0; i<n; i++) {
1172e884886fSBarry Smith     if (PetscRealPart(u_vec[i] + *h*a_vec[i]) <= 0.0) {
1173e884886fSBarry Smith       val = PetscAbsScalar(u_vec[i]/a_vec[i]);
1174e884886fSBarry Smith       if (val < minval) minval = val;
1175e884886fSBarry Smith     }
1176e884886fSBarry Smith   }
1177e884886fSBarry Smith   ierr = VecRestoreArray(U,&u_vec);CHKERRQ(ierr);
1178e884886fSBarry Smith   ierr = VecRestoreArray(a,&a_vec);CHKERRQ(ierr);
1179b2566f29SBarry Smith   ierr = MPIU_Allreduce(&minval,&val,1,MPIU_REAL,MPIU_MIN,comm);CHKERRQ(ierr);
1180e884886fSBarry Smith   if (val <= PetscAbsScalar(*h)) {
11816712e2f1SBarry Smith     ierr = PetscInfo2(U,"Scaling back h from %g to %g\n",(double)PetscRealPart(*h),(double)(.99*val));CHKERRQ(ierr);
1182e884886fSBarry Smith     if (PetscRealPart(*h) > 0.0) *h =  0.99*val;
1183e884886fSBarry Smith     else                         *h = -0.99*val;
1184e884886fSBarry Smith   }
1185e884886fSBarry Smith   PetscFunctionReturn(0);
1186e884886fSBarry Smith }
1187