xref: /petsc/src/mat/impls/mffd/mffd.c (revision 37e93019d958c82a17ba35bccd46909a573bf1c3)
1e884886fSBarry Smith 
2b45d2f2cSJed Brown #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 #undef __FUNCT__
13b022a5c1SBarry Smith #define __FUNCT__ "MatMFFDFinalizePackage"
14b022a5c1SBarry Smith /*@C
152390153bSJed Brown   MatMFFDFinalizePackage - This function destroys everything in the MatMFFD package. It is
16b022a5c1SBarry Smith   called from PetscFinalize().
17b022a5c1SBarry Smith 
18b022a5c1SBarry Smith   Level: developer
19b022a5c1SBarry Smith 
202390153bSJed Brown .keywords: Petsc, destroy, package
21b022a5c1SBarry Smith .seealso: PetscFinalize()
22b022a5c1SBarry Smith @*/
237087cfbeSBarry Smith PetscErrorCode  MatMFFDFinalizePackage(void)
24b022a5c1SBarry Smith {
25a703d84cSPatrick Lacasse   PetscErrorCode ierr;
26a703d84cSPatrick Lacasse 
27b022a5c1SBarry Smith   PetscFunctionBegin;
28*37e93019SBarry Smith   ierr = PetscFunctionListDestroy(&MatMFFDList);CHKERRQ(ierr);
29b022a5c1SBarry Smith   MatMFFDPackageInitialized = PETSC_FALSE;
30b022a5c1SBarry Smith   MatMFFDRegisterAllCalled  = PETSC_FALSE;
31b022a5c1SBarry Smith   PetscFunctionReturn(0);
32b022a5c1SBarry Smith }
33b022a5c1SBarry Smith 
34e884886fSBarry Smith #undef __FUNCT__
353ec795f1SBarry Smith #define __FUNCT__ "MatMFFDInitializePackage"
363ec795f1SBarry Smith /*@C
373ec795f1SBarry Smith   MatMFFDInitializePackage - This function initializes everything in the MatMFFD package. It is called
383ec795f1SBarry Smith   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to MatCreate_MFFD()
393ec795f1SBarry Smith   when using static libraries.
403ec795f1SBarry Smith 
413ec795f1SBarry Smith   Level: developer
423ec795f1SBarry Smith 
433ec795f1SBarry Smith .keywords: Vec, initialize, package
443ec795f1SBarry Smith .seealso: PetscInitialize()
453ec795f1SBarry Smith @*/
46607a6623SBarry Smith PetscErrorCode  MatMFFDInitializePackage(void)
473ec795f1SBarry Smith {
483ec795f1SBarry Smith   char           logList[256];
493ec795f1SBarry Smith   char           *className;
50ace3abfcSBarry Smith   PetscBool      opt;
513ec795f1SBarry Smith   PetscErrorCode ierr;
523ec795f1SBarry Smith 
533ec795f1SBarry Smith   PetscFunctionBegin;
54b022a5c1SBarry Smith   if (MatMFFDPackageInitialized) PetscFunctionReturn(0);
55b022a5c1SBarry Smith   MatMFFDPackageInitialized = PETSC_TRUE;
563ec795f1SBarry Smith   /* Register Classes */
570700a824SBarry Smith   ierr = PetscClassIdRegister("MatMFFD",&MATMFFD_CLASSID);CHKERRQ(ierr);
583ec795f1SBarry Smith   /* Register Constructors */
59607a6623SBarry Smith   ierr = MatMFFDRegisterAll();CHKERRQ(ierr);
603ec795f1SBarry Smith   /* Register Events */
610700a824SBarry Smith   ierr = PetscLogEventRegister("MatMult MF",          MATMFFD_CLASSID,&MATMFFD_Mult);CHKERRQ(ierr);
623ec795f1SBarry Smith 
633ec795f1SBarry Smith   /* Process info exclusions */
640298fd71SBarry Smith   ierr = PetscOptionsGetString(NULL, "-info_exclude", logList, 256, &opt);CHKERRQ(ierr);
653ec795f1SBarry Smith   if (opt) {
663ec795f1SBarry Smith     ierr = PetscStrstr(logList, "matmffd", &className);CHKERRQ(ierr);
673ec795f1SBarry Smith     if (className) {
680700a824SBarry Smith       ierr = PetscInfoDeactivateClass(MATMFFD_CLASSID);CHKERRQ(ierr);
693ec795f1SBarry Smith     }
703ec795f1SBarry Smith   }
713ec795f1SBarry Smith   /* Process summary exclusions */
720298fd71SBarry Smith   ierr = PetscOptionsGetString(NULL, "-log_summary_exclude", logList, 256, &opt);CHKERRQ(ierr);
733ec795f1SBarry Smith   if (opt) {
743ec795f1SBarry Smith     ierr = PetscStrstr(logList, "matmffd", &className);CHKERRQ(ierr);
753ec795f1SBarry Smith     if (className) {
760700a824SBarry Smith       ierr = PetscLogEventDeactivateClass(MATMFFD_CLASSID);CHKERRQ(ierr);
773ec795f1SBarry Smith     }
783ec795f1SBarry Smith   }
79b022a5c1SBarry Smith   ierr = PetscRegisterFinalize(MatMFFDFinalizePackage);CHKERRQ(ierr);
803ec795f1SBarry Smith   PetscFunctionReturn(0);
813ec795f1SBarry Smith }
823ec795f1SBarry Smith 
833ec795f1SBarry Smith #undef __FUNCT__
84e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetType"
85e884886fSBarry Smith /*@C
86e884886fSBarry Smith     MatMFFDSetType - Sets the method that is used to compute the
87e884886fSBarry Smith     differencing parameter for finite differene matrix-free formulations.
88e884886fSBarry Smith 
89e884886fSBarry Smith     Input Parameters:
90e884886fSBarry Smith +   mat - the "matrix-free" matrix created via MatCreateSNESMF(), or MatCreateMFFD()
91e884886fSBarry Smith           or MatSetType(mat,MATMFFD);
92e884886fSBarry Smith -   ftype - the type requested, either MATMFFD_WP or MATMFFD_DS
93e884886fSBarry Smith 
94e884886fSBarry Smith     Level: advanced
95e884886fSBarry Smith 
96e884886fSBarry Smith     Notes:
97e884886fSBarry Smith     For example, such routines can compute h for use in
98e884886fSBarry Smith     Jacobian-vector products of the form
99e884886fSBarry Smith 
100e884886fSBarry Smith                         F(x+ha) - F(x)
101e884886fSBarry Smith           F'(u)a  ~=  ----------------
102e884886fSBarry Smith                               h
103e884886fSBarry Smith 
1041c84c290SBarry Smith .seealso: MatCreateSNESMF(), MatMFFDRegister(), MatMFFDSetFunction()
105e884886fSBarry Smith @*/
10619fd82e9SBarry Smith PetscErrorCode  MatMFFDSetType(Mat mat,MatMFFDType ftype)
107e884886fSBarry Smith {
108e884886fSBarry Smith   PetscErrorCode ierr,(*r)(MatMFFD);
109e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)mat->data;
110ace3abfcSBarry Smith   PetscBool      match;
111e884886fSBarry Smith 
112e884886fSBarry Smith   PetscFunctionBegin;
1130700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
114e884886fSBarry Smith   PetscValidCharPointer(ftype,2);
115e884886fSBarry Smith 
116251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)mat,MATMFFD,&match);CHKERRQ(ierr);
1179a13eae7SJed Brown   if (!match) PetscFunctionReturn(0);
1189a13eae7SJed Brown 
119e884886fSBarry Smith   /* already set, so just return */
120251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)ctx,ftype,&match);CHKERRQ(ierr);
121e884886fSBarry Smith   if (match) PetscFunctionReturn(0);
122e884886fSBarry Smith 
123e884886fSBarry Smith   /* destroy the old one if it exists */
124e884886fSBarry Smith   if (ctx->ops->destroy) {
125e884886fSBarry Smith     ierr = (*ctx->ops->destroy)(ctx);CHKERRQ(ierr);
126e884886fSBarry Smith   }
127e884886fSBarry Smith 
1281c9cd337SJed Brown   ierr =  PetscFunctionListFind(MatMFFDList,ftype,&r);CHKERRQ(ierr);
129e32f2f54SBarry Smith   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown MatMFFD type %s given",ftype);
130e884886fSBarry Smith   ierr = (*r)(ctx);CHKERRQ(ierr);
131e884886fSBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)ctx,ftype);CHKERRQ(ierr);
132e884886fSBarry Smith   PetscFunctionReturn(0);
133e884886fSBarry Smith }
134e884886fSBarry Smith 
135e884886fSBarry Smith typedef PetscErrorCode (*FCN1)(void*,Vec); /* force argument to next function to not be extern C*/
136e884886fSBarry Smith #undef __FUNCT__
137c879e56bSBarry Smith #define __FUNCT__ "MatMFFDSetFunctioniBase_MFFD"
1387087cfbeSBarry Smith PetscErrorCode  MatMFFDSetFunctioniBase_MFFD(Mat mat,FCN1 func)
139e884886fSBarry Smith {
140e884886fSBarry Smith   MatMFFD ctx = (MatMFFD)mat->data;
141e884886fSBarry Smith 
142e884886fSBarry Smith   PetscFunctionBegin;
143e884886fSBarry Smith   ctx->funcisetbase = func;
144e884886fSBarry Smith   PetscFunctionReturn(0);
145e884886fSBarry Smith }
146e884886fSBarry Smith 
147e884886fSBarry Smith typedef PetscErrorCode (*FCN2)(void*,PetscInt,Vec,PetscScalar*); /* force argument to next function to not be extern C*/
148e884886fSBarry Smith #undef __FUNCT__
149c879e56bSBarry Smith #define __FUNCT__ "MatMFFDSetFunctioni_MFFD"
1507087cfbeSBarry Smith PetscErrorCode  MatMFFDSetFunctioni_MFFD(Mat mat,FCN2 funci)
151e884886fSBarry Smith {
152e884886fSBarry Smith   MatMFFD ctx = (MatMFFD)mat->data;
153e884886fSBarry Smith 
154e884886fSBarry Smith   PetscFunctionBegin;
155e884886fSBarry Smith   ctx->funci = funci;
156e884886fSBarry Smith   PetscFunctionReturn(0);
157e884886fSBarry Smith }
158e884886fSBarry Smith 
159bc0f08ceSBarry Smith #undef __FUNCT__
160bc0f08ceSBarry Smith #define __FUNCT__ "MatMFFDResetHHistory_MFFD"
161bc0f08ceSBarry Smith PetscErrorCode  MatMFFDResetHHistory_MFFD(Mat J)
162bc0f08ceSBarry Smith {
163bc0f08ceSBarry Smith   MatMFFD ctx = (MatMFFD)J->data;
164bc0f08ceSBarry Smith 
165bc0f08ceSBarry Smith   PetscFunctionBegin;
166bc0f08ceSBarry Smith   ctx->ncurrenth = 0;
167bc0f08ceSBarry Smith   PetscFunctionReturn(0);
168bc0f08ceSBarry Smith }
169e884886fSBarry Smith 
170e884886fSBarry Smith #undef __FUNCT__
171e884886fSBarry Smith #define __FUNCT__ "MatMFFDRegister"
1721c84c290SBarry Smith /*@C
1731c84c290SBarry Smith    MatMFFDRegister - Adds a method to the MatMFFD registry.
1741c84c290SBarry Smith 
1751c84c290SBarry Smith    Not Collective
1761c84c290SBarry Smith 
1771c84c290SBarry Smith    Input Parameters:
1781c84c290SBarry Smith +  name_solver - name of a new user-defined compute-h module
1791c84c290SBarry Smith -  routine_create - routine to create method context
1801c84c290SBarry Smith 
1811c84c290SBarry Smith    Level: developer
1821c84c290SBarry Smith 
1831c84c290SBarry Smith    Notes:
1841c84c290SBarry Smith    MatMFFDRegister() may be called multiple times to add several user-defined solvers.
1851c84c290SBarry Smith 
1861c84c290SBarry Smith    Sample usage:
1871c84c290SBarry Smith .vb
188bdf89e91SBarry Smith    MatMFFDRegister("my_h",MyHCreate);
1891c84c290SBarry Smith .ve
1901c84c290SBarry Smith 
1911c84c290SBarry Smith    Then, your solver can be chosen with the procedural interface via
1921c84c290SBarry Smith $     MatMFFDSetType(mfctx,"my_h")
1931c84c290SBarry Smith    or at runtime via the option
1941c84c290SBarry Smith $     -snes_mf_type my_h
1951c84c290SBarry Smith 
1961c84c290SBarry Smith .keywords: MatMFFD, register
1971c84c290SBarry Smith 
1981c84c290SBarry Smith .seealso: MatMFFDRegisterAll(), MatMFFDRegisterDestroy()
1991c84c290SBarry Smith  @*/
200bdf89e91SBarry Smith PetscErrorCode  MatMFFDRegister(const char sname[],PetscErrorCode (*function)(MatMFFD))
201e884886fSBarry Smith {
202e884886fSBarry Smith   PetscErrorCode ierr;
203e884886fSBarry Smith 
204e884886fSBarry Smith   PetscFunctionBegin;
205a240a19fSJed Brown   ierr = PetscFunctionListAdd(&MatMFFDList,sname,function);CHKERRQ(ierr);
206e884886fSBarry Smith   PetscFunctionReturn(0);
207e884886fSBarry Smith }
208e884886fSBarry Smith 
209e884886fSBarry Smith 
210e884886fSBarry Smith #undef __FUNCT__
211e884886fSBarry Smith #define __FUNCT__ "MatMFFDRegisterDestroy"
212e884886fSBarry Smith /*@C
213e884886fSBarry Smith    MatMFFDRegisterDestroy - Frees the list of MatMFFD methods that were
21481242352SJed Brown    registered by MatMFFDRegister().
215e884886fSBarry Smith 
216e884886fSBarry Smith    Not Collective
217e884886fSBarry Smith 
218e884886fSBarry Smith    Level: developer
219e884886fSBarry Smith 
220e884886fSBarry Smith .keywords: MatMFFD, register, destroy
221e884886fSBarry Smith 
22281242352SJed Brown .seealso: MatMFFDRegister(), MatMFFDRegisterAll()
223e884886fSBarry Smith @*/
2247087cfbeSBarry Smith PetscErrorCode  MatMFFDRegisterDestroy(void)
225e884886fSBarry Smith {
226e884886fSBarry Smith   PetscErrorCode ierr;
227e884886fSBarry Smith 
228e884886fSBarry Smith   PetscFunctionBegin;
229140e18c1SBarry Smith   ierr = PetscFunctionListDestroy(&MatMFFDList);CHKERRQ(ierr);
2302205254eSKarl Rupp 
231e884886fSBarry Smith   MatMFFDRegisterAllCalled = PETSC_FALSE;
232e884886fSBarry Smith   PetscFunctionReturn(0);
233e884886fSBarry Smith }
234e884886fSBarry Smith 
235bc0f08ceSBarry Smith #undef __FUNCT__
236bc0f08ceSBarry Smith #define __FUNCT__ "MatMFFDAddNullSpace_MFFD"
237bc0f08ceSBarry Smith PetscErrorCode  MatMFFDAddNullSpace_MFFD(Mat J,MatNullSpace nullsp)
238bc0f08ceSBarry Smith {
239bc0f08ceSBarry Smith   PetscErrorCode ierr;
240bc0f08ceSBarry Smith   MatMFFD        ctx = (MatMFFD)J->data;
241bc0f08ceSBarry Smith 
242bc0f08ceSBarry Smith   PetscFunctionBegin;
243bc0f08ceSBarry Smith   ierr = PetscObjectReference((PetscObject)nullsp);CHKERRQ(ierr);
244bc0f08ceSBarry Smith   if (ctx->sp) { ierr = MatNullSpaceDestroy(&ctx->sp);CHKERRQ(ierr); }
245bc0f08ceSBarry Smith   ctx->sp = nullsp;
246bc0f08ceSBarry Smith   PetscFunctionReturn(0);
247bc0f08ceSBarry Smith }
248bc0f08ceSBarry Smith 
249e884886fSBarry Smith /* ----------------------------------------------------------------------------------------*/
250e884886fSBarry Smith #undef __FUNCT__
251e884886fSBarry Smith #define __FUNCT__ "MatDestroy_MFFD"
252e884886fSBarry Smith PetscErrorCode MatDestroy_MFFD(Mat mat)
253e884886fSBarry Smith {
254e884886fSBarry Smith   PetscErrorCode ierr;
255e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)mat->data;
256e884886fSBarry Smith 
257e884886fSBarry Smith   PetscFunctionBegin;
2586bf464f9SBarry Smith   ierr = VecDestroy(&ctx->w);CHKERRQ(ierr);
2596bf464f9SBarry Smith   ierr = VecDestroy(&ctx->drscale);CHKERRQ(ierr);
2606bf464f9SBarry Smith   ierr = VecDestroy(&ctx->dlscale);CHKERRQ(ierr);
2616bf464f9SBarry Smith   ierr = VecDestroy(&ctx->dshift);CHKERRQ(ierr);
262cfe22f5eSBarry Smith   if (ctx->current_f_allocated) {
2636bf464f9SBarry Smith     ierr = VecDestroy(&ctx->current_f);CHKERRQ(ierr);
264cfe22f5eSBarry Smith   }
265e884886fSBarry Smith   if (ctx->ops->destroy) {ierr = (*ctx->ops->destroy)(ctx);CHKERRQ(ierr);}
2666bf464f9SBarry Smith   ierr      = MatNullSpaceDestroy(&ctx->sp);CHKERRQ(ierr);
2676bf464f9SBarry Smith   ierr      = PetscHeaderDestroy(&ctx);CHKERRQ(ierr);
268bf0cc555SLisandro Dalcin   mat->data = 0;
269e884886fSBarry Smith 
270bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetBase_C",NULL);CHKERRQ(ierr);
271bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctioniBase_C",NULL);CHKERRQ(ierr);
272bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctioni_C",NULL);CHKERRQ(ierr);
273bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunction_C",NULL);CHKERRQ(ierr);
274bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctionError_C",NULL);CHKERRQ(ierr);
275bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetCheckh_C",NULL);CHKERRQ(ierr);
276bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetPeriod_C",NULL);CHKERRQ(ierr);
277bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDResetHHistory_C",NULL);CHKERRQ(ierr);
278bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDAddNullSpace_C",NULL);CHKERRQ(ierr);
279e884886fSBarry Smith   PetscFunctionReturn(0);
280e884886fSBarry Smith }
281e884886fSBarry Smith 
282e884886fSBarry Smith #undef __FUNCT__
283e884886fSBarry Smith #define __FUNCT__ "MatView_MFFD"
284e884886fSBarry Smith /*
285e884886fSBarry Smith    MatMFFDView_MFFD - Views matrix-free parameters.
286e884886fSBarry Smith 
287e884886fSBarry Smith */
288e884886fSBarry Smith PetscErrorCode MatView_MFFD(Mat J,PetscViewer viewer)
289e884886fSBarry Smith {
290e884886fSBarry Smith   PetscErrorCode ierr;
291e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)J->data;
292a283635eSDmitry Karpeev   PetscBool      iascii, viewbase, viewfunction;
293a283635eSDmitry Karpeev   const char     *prefix;
294e884886fSBarry Smith 
295e884886fSBarry Smith   PetscFunctionBegin;
296251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
297e884886fSBarry Smith   if (iascii) {
298a283635eSDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"Matrix-free approximation:\n");CHKERRQ(ierr);
299a283635eSDmitry Karpeev     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
300e884886fSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"err=%G (relative error in function evaluation)\n",ctx->error_rel);CHKERRQ(ierr);
3017adad957SLisandro Dalcin     if (!((PetscObject)ctx)->type_name) {
302e884886fSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"The compute h routine has not yet been set\n");CHKERRQ(ierr);
303e884886fSBarry Smith     } else {
3047adad957SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"Using %s compute h routine\n",((PetscObject)ctx)->type_name);CHKERRQ(ierr);
305e884886fSBarry Smith     }
306e884886fSBarry Smith     if (ctx->ops->view) {
307e884886fSBarry Smith       ierr = (*ctx->ops->view)(ctx,viewer);CHKERRQ(ierr);
308e884886fSBarry Smith     }
309a283635eSDmitry Karpeev     ierr = PetscObjectGetOptionsPrefix((PetscObject)J, &prefix);CHKERRQ(ierr);
310a283635eSDmitry Karpeev 
311a283635eSDmitry Karpeev     ierr = PetscOptionsHasName(prefix, "-mat_mffd_view_base", &viewbase);CHKERRQ(ierr);
312a283635eSDmitry Karpeev     if (viewbase) {
313a283635eSDmitry Karpeev       ierr = PetscViewerASCIIPrintf(viewer, "Base:\n");CHKERRQ(ierr);
314a283635eSDmitry Karpeev       ierr = VecView(ctx->current_u, viewer);CHKERRQ(ierr);
315a283635eSDmitry Karpeev     }
316a283635eSDmitry Karpeev     ierr = PetscOptionsHasName(prefix, "-mat_mffd_view_function", &viewfunction);CHKERRQ(ierr);
317a283635eSDmitry Karpeev     if (viewfunction) {
318a283635eSDmitry Karpeev       ierr = PetscViewerASCIIPrintf(viewer, "Function:\n");CHKERRQ(ierr);
319a283635eSDmitry Karpeev       ierr = VecView(ctx->current_f, viewer);CHKERRQ(ierr);
320a283635eSDmitry Karpeev     }
321a283635eSDmitry Karpeev     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
32211aeaf0aSBarry Smith   }
323e884886fSBarry Smith   PetscFunctionReturn(0);
324e884886fSBarry Smith }
325e884886fSBarry Smith 
326e884886fSBarry Smith #undef __FUNCT__
327e884886fSBarry Smith #define __FUNCT__ "MatAssemblyEnd_MFFD"
328e884886fSBarry Smith /*
329e884886fSBarry Smith    MatAssemblyEnd_MFFD - Resets the ctx->ncurrenth to zero. This
330e884886fSBarry Smith    allows the user to indicate the beginning of a new linear solve by calling
331e884886fSBarry Smith    MatAssemblyXXX() on the matrix free matrix. This then allows the
3321d0fab5eSBarry Smith    MatCreateMFFD_WP() to properly compute ||U|| only the first time
333e884886fSBarry Smith    in the linear solver rather than every time.
3345a576424SJed Brown 
3355a576424SJed Brown    This function is referenced directly from MatAssemblyEnd_SNESMF(), which may be in a different shared library.
336e884886fSBarry Smith */
3375a576424SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyEnd_MFFD(Mat J,MatAssemblyType mt)
338e884886fSBarry Smith {
339e884886fSBarry Smith   PetscErrorCode ierr;
340e884886fSBarry Smith   MatMFFD        j = (MatMFFD)J->data;
341e884886fSBarry Smith 
342e884886fSBarry Smith   PetscFunctionBegin;
343e884886fSBarry Smith   ierr      = MatMFFDResetHHistory(J);CHKERRQ(ierr);
344e884886fSBarry Smith   j->vshift = 0.0;
345e884886fSBarry Smith   j->vscale = 1.0;
346e884886fSBarry Smith   PetscFunctionReturn(0);
347e884886fSBarry Smith }
348e884886fSBarry Smith 
349e884886fSBarry Smith #undef __FUNCT__
350e884886fSBarry Smith #define __FUNCT__ "MatMult_MFFD"
351e884886fSBarry Smith /*
352e884886fSBarry Smith   MatMult_MFFD - Default matrix-free form for Jacobian-vector product, y = F'(u)*a:
353e884886fSBarry Smith 
354e884886fSBarry Smith         y ~= (F(u + ha) - F(u))/h,
355e884886fSBarry Smith   where F = nonlinear function, as set by SNESSetFunction()
356e884886fSBarry Smith         u = current iterate
357e884886fSBarry Smith         h = difference interval
358e884886fSBarry Smith */
359e884886fSBarry Smith PetscErrorCode MatMult_MFFD(Mat mat,Vec a,Vec y)
360e884886fSBarry Smith {
361e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)mat->data;
362e884886fSBarry Smith   PetscScalar    h;
363e884886fSBarry Smith   Vec            w,U,F;
364e884886fSBarry Smith   PetscErrorCode ierr;
365ace3abfcSBarry Smith   PetscBool      zeroa;
366e884886fSBarry Smith 
367e884886fSBarry Smith   PetscFunctionBegin;
368ce94432eSBarry 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");
369e884886fSBarry Smith   /* We log matrix-free matrix-vector products separately, so that we can
370e884886fSBarry Smith      separate the performance monitoring from the cases that use conventional
371e884886fSBarry Smith      storage.  We may eventually modify event logging to associate events
372e884886fSBarry Smith      with particular objects, hence alleviating the more general problem. */
373e884886fSBarry Smith   ierr = PetscLogEventBegin(MATMFFD_Mult,a,y,0,0);CHKERRQ(ierr);
374e884886fSBarry Smith 
375e884886fSBarry Smith   w = ctx->w;
376e884886fSBarry Smith   U = ctx->current_u;
3773ec795f1SBarry Smith   F = ctx->current_f;
378e884886fSBarry Smith   /*
379e884886fSBarry Smith       Compute differencing parameter
380e884886fSBarry Smith   */
381e884886fSBarry Smith   if (!ctx->ops->compute) {
382e884886fSBarry Smith     ierr = MatMFFDSetType(mat,MATMFFD_WP);CHKERRQ(ierr);
3839c6ac3b3SBarry Smith     ierr = MatSetFromOptions(mat);CHKERRQ(ierr);
384e884886fSBarry Smith   }
385e884886fSBarry Smith   ierr = (*ctx->ops->compute)(ctx,U,a,&h,&zeroa);CHKERRQ(ierr);
386e884886fSBarry Smith   if (zeroa) {
387e884886fSBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
388e884886fSBarry Smith     PetscFunctionReturn(0);
389e884886fSBarry Smith   }
390e884886fSBarry Smith 
391e32f2f54SBarry Smith   if (PetscIsInfOrNanScalar(h)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Computed Nan differencing parameter h");
392e884886fSBarry Smith   if (ctx->checkh) {
393e884886fSBarry Smith     ierr = (*ctx->checkh)(ctx->checkhctx,U,a,&h);CHKERRQ(ierr);
394e884886fSBarry Smith   }
395e884886fSBarry Smith 
396e884886fSBarry Smith   /* keep a record of the current differencing parameter h */
397e884886fSBarry Smith   ctx->currenth = h;
398e884886fSBarry Smith #if defined(PETSC_USE_COMPLEX)
399e884886fSBarry Smith   ierr = PetscInfo2(mat,"Current differencing parameter: %G + %G i\n",PetscRealPart(h),PetscImaginaryPart(h));CHKERRQ(ierr);
400e884886fSBarry Smith #else
401e884886fSBarry Smith   ierr = PetscInfo1(mat,"Current differencing parameter: %15.12e\n",h);CHKERRQ(ierr);
402e884886fSBarry Smith #endif
403e884886fSBarry Smith   if (ctx->historyh && ctx->ncurrenth < ctx->maxcurrenth) {
404e884886fSBarry Smith     ctx->historyh[ctx->ncurrenth] = h;
405e884886fSBarry Smith   }
406e884886fSBarry Smith   ctx->ncurrenth++;
407e884886fSBarry Smith 
408e884886fSBarry Smith   /* w = u + ha */
409c3bb7e23SBarry Smith   if (ctx->drscale) {
410c3bb7e23SBarry Smith     ierr = VecPointwiseMult(ctx->drscale,a,U);CHKERRQ(ierr);
411c3bb7e23SBarry Smith     ierr = VecAYPX(U,h,w);CHKERRQ(ierr);
412c3bb7e23SBarry Smith   } else {
413e884886fSBarry Smith     ierr = VecWAXPY(w,h,a,U);CHKERRQ(ierr);
414c3bb7e23SBarry Smith   }
415e884886fSBarry Smith 
4161b797266SDmitry Karpeev   /* compute func(U) as base for differencing; only needed first time in and not when provided by user */
4171b797266SDmitry Karpeev   if (ctx->ncurrenth == 1 && ctx->current_f_allocated) {
418e884886fSBarry Smith     ierr = (*ctx->func)(ctx->funcctx,U,F);CHKERRQ(ierr);
41952121784SBarry Smith   }
420e884886fSBarry Smith   ierr = (*ctx->func)(ctx->funcctx,w,y);CHKERRQ(ierr);
421e884886fSBarry Smith 
422e884886fSBarry Smith   ierr = VecAXPY(y,-1.0,F);CHKERRQ(ierr);
423e884886fSBarry Smith   ierr = VecScale(y,1.0/h);CHKERRQ(ierr);
424e884886fSBarry Smith 
425c35ec66cSMatthew G Knepley   if ((ctx->vshift != 0.0) || (ctx->vscale != 1.0)) {
426e884886fSBarry Smith     ierr = VecAXPBY(y,ctx->vshift,ctx->vscale,a);CHKERRQ(ierr);
427c3bb7e23SBarry Smith   }
428c3bb7e23SBarry Smith   if (ctx->dlscale) {
429c3bb7e23SBarry Smith     ierr = VecPointwiseMult(y,ctx->dlscale,y);CHKERRQ(ierr);
430c3bb7e23SBarry Smith   }
431c3bb7e23SBarry Smith   if (ctx->dshift) {
432c3bb7e23SBarry Smith     ierr = VecPointwiseMult(ctx->dshift,a,U);CHKERRQ(ierr);
433c3bb7e23SBarry Smith     ierr = VecAXPY(y,1.0,U);CHKERRQ(ierr);
434c3bb7e23SBarry Smith   }
435e884886fSBarry Smith 
4360298fd71SBarry Smith   if (ctx->sp) {ierr = MatNullSpaceRemove(ctx->sp,y,NULL);CHKERRQ(ierr);}
437e884886fSBarry Smith 
438e884886fSBarry Smith   ierr = PetscLogEventEnd(MATMFFD_Mult,a,y,0,0);CHKERRQ(ierr);
439e884886fSBarry Smith   PetscFunctionReturn(0);
440e884886fSBarry Smith }
441e884886fSBarry Smith 
442e884886fSBarry Smith #undef __FUNCT__
443e884886fSBarry Smith #define __FUNCT__ "MatGetDiagonal_MFFD"
444e884886fSBarry Smith /*
445e884886fSBarry Smith   MatGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix
446e884886fSBarry Smith 
447e884886fSBarry Smith         y ~= (F(u + ha) - F(u))/h,
448e884886fSBarry Smith   where F = nonlinear function, as set by SNESSetFunction()
449e884886fSBarry Smith         u = current iterate
450e884886fSBarry Smith         h = difference interval
451e884886fSBarry Smith */
452e884886fSBarry Smith PetscErrorCode MatGetDiagonal_MFFD(Mat mat,Vec a)
453e884886fSBarry Smith {
454e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)mat->data;
455e884886fSBarry Smith   PetscScalar    h,*aa,*ww,v;
456e884886fSBarry Smith   PetscReal      epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON;
457e884886fSBarry Smith   Vec            w,U;
458e884886fSBarry Smith   PetscErrorCode ierr;
459e884886fSBarry Smith   PetscInt       i,rstart,rend;
460e884886fSBarry Smith 
461e884886fSBarry Smith   PetscFunctionBegin;
462e7e72b3dSBarry Smith   if (!ctx->funci) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Requires calling MatMFFDSetFunctioni() first");
463e884886fSBarry Smith 
464e884886fSBarry Smith   w    = ctx->w;
465e884886fSBarry Smith   U    = ctx->current_u;
466e884886fSBarry Smith   ierr = (*ctx->func)(ctx->funcctx,U,a);CHKERRQ(ierr);
467e884886fSBarry Smith   ierr = (*ctx->funcisetbase)(ctx->funcctx,U);CHKERRQ(ierr);
468e884886fSBarry Smith   ierr = VecCopy(U,w);CHKERRQ(ierr);
469e884886fSBarry Smith 
470e884886fSBarry Smith   ierr = VecGetOwnershipRange(a,&rstart,&rend);CHKERRQ(ierr);
471e884886fSBarry Smith   ierr = VecGetArray(a,&aa);CHKERRQ(ierr);
472e884886fSBarry Smith   for (i=rstart; i<rend; i++) {
473e884886fSBarry Smith     ierr = VecGetArray(w,&ww);CHKERRQ(ierr);
474e884886fSBarry Smith     h    = ww[i-rstart];
475e884886fSBarry Smith     if (h == 0.0) h = 1.0;
476e884886fSBarry Smith     if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0)     h = umin;
477e884886fSBarry Smith     else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin;
478e884886fSBarry Smith     h *= epsilon;
479e884886fSBarry Smith 
480e884886fSBarry Smith     ww[i-rstart] += h;
481e884886fSBarry Smith     ierr          = VecRestoreArray(w,&ww);CHKERRQ(ierr);
482e884886fSBarry Smith     ierr          = (*ctx->funci)(ctx->funcctx,i,w,&v);CHKERRQ(ierr);
483e884886fSBarry Smith     aa[i-rstart]  = (v - aa[i-rstart])/h;
484e884886fSBarry Smith 
485e884886fSBarry Smith     /* possibly shift and scale result */
486c35ec66cSMatthew G Knepley     if ((ctx->vshift != 0.0) || (ctx->vscale != 1.0)) {
487e884886fSBarry Smith       aa[i - rstart] = ctx->vshift + ctx->vscale*aa[i-rstart];
488c3bb7e23SBarry Smith     }
489e884886fSBarry Smith 
490e884886fSBarry Smith     ierr          = VecGetArray(w,&ww);CHKERRQ(ierr);
491e884886fSBarry Smith     ww[i-rstart] -= h;
492e884886fSBarry Smith     ierr          = VecRestoreArray(w,&ww);CHKERRQ(ierr);
493e884886fSBarry Smith   }
494e884886fSBarry Smith   ierr = VecRestoreArray(a,&aa);CHKERRQ(ierr);
495e884886fSBarry Smith   PetscFunctionReturn(0);
496e884886fSBarry Smith }
497e884886fSBarry Smith 
498e884886fSBarry Smith #undef __FUNCT__
499c3bb7e23SBarry Smith #define __FUNCT__ "MatDiagonalScale_MFFD"
500c3bb7e23SBarry Smith PetscErrorCode MatDiagonalScale_MFFD(Mat mat,Vec ll,Vec rr)
501c3bb7e23SBarry Smith {
502c3bb7e23SBarry Smith   MatMFFD        aij = (MatMFFD)mat->data;
503c3bb7e23SBarry Smith   PetscErrorCode ierr;
504c3bb7e23SBarry Smith 
505c3bb7e23SBarry Smith   PetscFunctionBegin;
506c3bb7e23SBarry Smith   if (ll && !aij->dlscale) {
507c3bb7e23SBarry Smith     ierr = VecDuplicate(ll,&aij->dlscale);CHKERRQ(ierr);
508c3bb7e23SBarry Smith   }
509c3bb7e23SBarry Smith   if (rr && !aij->drscale) {
510c3bb7e23SBarry Smith     ierr = VecDuplicate(rr,&aij->drscale);CHKERRQ(ierr);
511c3bb7e23SBarry Smith   }
512c3bb7e23SBarry Smith   if (ll) {
513c3bb7e23SBarry Smith     ierr = VecCopy(ll,aij->dlscale);CHKERRQ(ierr);
514c3bb7e23SBarry Smith   }
515c3bb7e23SBarry Smith   if (rr) {
516c3bb7e23SBarry Smith     ierr = VecCopy(rr,aij->drscale);CHKERRQ(ierr);
517c3bb7e23SBarry Smith   }
518c3bb7e23SBarry Smith   PetscFunctionReturn(0);
519c3bb7e23SBarry Smith }
520c3bb7e23SBarry Smith 
521c3bb7e23SBarry Smith #undef __FUNCT__
522c3bb7e23SBarry Smith #define __FUNCT__ "MatDiagonalSet_MFFD"
523c3bb7e23SBarry Smith PetscErrorCode MatDiagonalSet_MFFD(Mat mat,Vec ll,InsertMode mode)
524c3bb7e23SBarry Smith {
525c3bb7e23SBarry Smith   MatMFFD        aij = (MatMFFD)mat->data;
526c3bb7e23SBarry Smith   PetscErrorCode ierr;
527c3bb7e23SBarry Smith 
528c3bb7e23SBarry Smith   PetscFunctionBegin;
529ce94432eSBarry Smith   if (mode == INSERT_VALUES) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"No diagonal set with INSERT_VALUES");
530c3bb7e23SBarry Smith   if (!aij->dshift) {
531c3bb7e23SBarry Smith     ierr = VecDuplicate(ll,&aij->dshift);CHKERRQ(ierr);
532c3bb7e23SBarry Smith   }
533c3bb7e23SBarry Smith   ierr = VecAXPY(aij->dshift,1.0,ll);CHKERRQ(ierr);
534c3bb7e23SBarry Smith   PetscFunctionReturn(0);
535c3bb7e23SBarry Smith }
536c3bb7e23SBarry Smith 
537c3bb7e23SBarry Smith #undef __FUNCT__
538e884886fSBarry Smith #define __FUNCT__ "MatShift_MFFD"
539e884886fSBarry Smith PetscErrorCode MatShift_MFFD(Mat Y,PetscScalar a)
540e884886fSBarry Smith {
541e884886fSBarry Smith   MatMFFD shell = (MatMFFD)Y->data;
5425fd66863SKarl Rupp 
543e884886fSBarry Smith   PetscFunctionBegin;
544e884886fSBarry Smith   shell->vshift += a;
545e884886fSBarry Smith   PetscFunctionReturn(0);
546e884886fSBarry Smith }
547e884886fSBarry Smith 
548e884886fSBarry Smith #undef __FUNCT__
549e884886fSBarry Smith #define __FUNCT__ "MatScale_MFFD"
550e884886fSBarry Smith PetscErrorCode MatScale_MFFD(Mat Y,PetscScalar a)
551e884886fSBarry Smith {
552e884886fSBarry Smith   MatMFFD shell = (MatMFFD)Y->data;
5535fd66863SKarl Rupp 
554e884886fSBarry Smith   PetscFunctionBegin;
555e884886fSBarry Smith   shell->vscale *= a;
556e884886fSBarry Smith   PetscFunctionReturn(0);
557e884886fSBarry Smith }
558e884886fSBarry Smith 
559e884886fSBarry Smith #undef __FUNCT__
560c879e56bSBarry Smith #define __FUNCT__ "MatMFFDSetBase_MFFD"
561d8fa6a54SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDSetBase_MFFD(Mat J,Vec U,Vec F)
562e884886fSBarry Smith {
563e884886fSBarry Smith   PetscErrorCode ierr;
564e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)J->data;
565e884886fSBarry Smith 
566e884886fSBarry Smith   PetscFunctionBegin;
567e884886fSBarry Smith   ierr = MatMFFDResetHHistory(J);CHKERRQ(ierr);
5682205254eSKarl Rupp 
569e884886fSBarry Smith   ctx->current_u = U;
57052121784SBarry Smith   if (F) {
5716bf464f9SBarry Smith     if (ctx->current_f_allocated) {ierr = VecDestroy(&ctx->current_f);CHKERRQ(ierr);}
5723ec795f1SBarry Smith     ctx->current_f           = F;
573cfe22f5eSBarry Smith     ctx->current_f_allocated = PETSC_FALSE;
574cfe22f5eSBarry Smith   } else if (!ctx->current_f_allocated) {
57552121784SBarry Smith     ierr = VecDuplicate(ctx->current_u, &ctx->current_f);CHKERRQ(ierr);
5762205254eSKarl Rupp 
577cfe22f5eSBarry Smith     ctx->current_f_allocated = PETSC_TRUE;
57852121784SBarry Smith   }
579e884886fSBarry Smith   if (!ctx->w) {
580e884886fSBarry Smith     ierr = VecDuplicate(ctx->current_u, &ctx->w);CHKERRQ(ierr);
581e884886fSBarry Smith   }
582e884886fSBarry Smith   J->assembled = PETSC_TRUE;
583e884886fSBarry Smith   PetscFunctionReturn(0);
584e884886fSBarry Smith }
5855a576424SJed Brown 
586e884886fSBarry Smith typedef PetscErrorCode (*FCN3)(void*,Vec,Vec,PetscScalar*); /* force argument to next function to not be extern C*/
587b2573a8aSBarry Smith 
588e884886fSBarry Smith #undef __FUNCT__
589c879e56bSBarry Smith #define __FUNCT__ "MatMFFDSetCheckh_MFFD"
5907087cfbeSBarry Smith PetscErrorCode  MatMFFDSetCheckh_MFFD(Mat J,FCN3 fun,void *ectx)
591e884886fSBarry Smith {
592e884886fSBarry Smith   MatMFFD ctx = (MatMFFD)J->data;
593e884886fSBarry Smith 
594e884886fSBarry Smith   PetscFunctionBegin;
595e884886fSBarry Smith   ctx->checkh    = fun;
596e884886fSBarry Smith   ctx->checkhctx = ectx;
597e884886fSBarry Smith   PetscFunctionReturn(0);
598e884886fSBarry Smith }
599e884886fSBarry Smith 
600e884886fSBarry Smith #undef __FUNCT__
6016aa9148fSLisandro Dalcin #define __FUNCT__ "MatMFFDSetOptionsPrefix"
6026aa9148fSLisandro Dalcin /*@C
6036aa9148fSLisandro Dalcin    MatMFFDSetOptionsPrefix - Sets the prefix used for searching for all
6046aa9148fSLisandro Dalcin    MatMFFD options in the database.
6056aa9148fSLisandro Dalcin 
6066aa9148fSLisandro Dalcin    Collective on Mat
6076aa9148fSLisandro Dalcin 
6086aa9148fSLisandro Dalcin    Input Parameter:
6096aa9148fSLisandro Dalcin +  A - the Mat context
6106aa9148fSLisandro Dalcin -  prefix - the prefix to prepend to all option names
6116aa9148fSLisandro Dalcin 
6126aa9148fSLisandro Dalcin    Notes:
6136aa9148fSLisandro Dalcin    A hyphen (-) must NOT be given at the beginning of the prefix name.
6146aa9148fSLisandro Dalcin    The first character of all runtime options is AUTOMATICALLY the hyphen.
6156aa9148fSLisandro Dalcin 
6166aa9148fSLisandro Dalcin    Level: advanced
6176aa9148fSLisandro Dalcin 
6186aa9148fSLisandro Dalcin .keywords: SNES, matrix-free, parameters
6196aa9148fSLisandro Dalcin 
6209c6ac3b3SBarry Smith .seealso: MatSetFromOptions(), MatCreateSNESMF()
6216aa9148fSLisandro Dalcin @*/
6227087cfbeSBarry Smith PetscErrorCode  MatMFFDSetOptionsPrefix(Mat mat,const char prefix[])
6236aa9148fSLisandro Dalcin 
6246aa9148fSLisandro Dalcin {
6250298fd71SBarry Smith   MatMFFD        mfctx = mat ? (MatMFFD)mat->data : (MatMFFD)NULL;
6266aa9148fSLisandro Dalcin   PetscErrorCode ierr;
6275fd66863SKarl Rupp 
6286aa9148fSLisandro Dalcin   PetscFunctionBegin;
6290700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
6300700a824SBarry Smith   PetscValidHeaderSpecific(mfctx,MATMFFD_CLASSID,1);
6316aa9148fSLisandro Dalcin   ierr = PetscObjectSetOptionsPrefix((PetscObject)mfctx,prefix);CHKERRQ(ierr);
6326aa9148fSLisandro Dalcin   PetscFunctionReturn(0);
6336aa9148fSLisandro Dalcin }
6346aa9148fSLisandro Dalcin 
6356aa9148fSLisandro Dalcin #undef __FUNCT__
6369c6ac3b3SBarry Smith #define __FUNCT__ "MatSetFromOptions_MFFD"
6379c6ac3b3SBarry Smith PetscErrorCode  MatSetFromOptions_MFFD(Mat mat)
638e884886fSBarry Smith {
63971f08665SBarry Smith   MatMFFD        mfctx = (MatMFFD)mat->data;
640e884886fSBarry Smith   PetscErrorCode ierr;
641ace3abfcSBarry Smith   PetscBool      flg;
642e884886fSBarry Smith   char           ftype[256];
643e884886fSBarry Smith 
644e884886fSBarry Smith   PetscFunctionBegin;
6450700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
6460700a824SBarry Smith   PetscValidHeaderSpecific(mfctx,MATMFFD_CLASSID,1);
6473194b578SJed Brown   ierr = PetscObjectOptionsBegin((PetscObject)mfctx);CHKERRQ(ierr);
648b022a5c1SBarry Smith   ierr = PetscOptionsList("-mat_mffd_type","Matrix free type","MatMFFDSetType",MatMFFDList,((PetscObject)mfctx)->type_name,ftype,256,&flg);CHKERRQ(ierr);
649e884886fSBarry Smith   if (flg) {
650e884886fSBarry Smith     ierr = MatMFFDSetType(mat,ftype);CHKERRQ(ierr);
651e884886fSBarry Smith   }
652e884886fSBarry Smith 
653e884886fSBarry Smith   ierr = PetscOptionsReal("-mat_mffd_err","set sqrt relative error in function","MatMFFDSetFunctionError",mfctx->error_rel,&mfctx->error_rel,0);CHKERRQ(ierr);
654e884886fSBarry Smith   ierr = PetscOptionsInt("-mat_mffd_period","how often h is recomputed","MatMFFDSetPeriod",mfctx->recomputeperiod,&mfctx->recomputeperiod,0);CHKERRQ(ierr);
655e884886fSBarry Smith 
65690d69ab7SBarry Smith   flg  = PETSC_FALSE;
6570298fd71SBarry Smith   ierr = PetscOptionsBool("-mat_mffd_check_positivity","Insure that U + h*a is nonnegative","MatMFFDSetCheckh",flg,&flg,NULL);CHKERRQ(ierr);
658e884886fSBarry Smith   if (flg) {
659e884886fSBarry Smith     ierr = MatMFFDSetCheckh(mat,MatMFFDCheckPositivity,0);CHKERRQ(ierr);
660e884886fSBarry Smith   }
661e884886fSBarry Smith   if (mfctx->ops->setfromoptions) {
662e884886fSBarry Smith     ierr = (*mfctx->ops->setfromoptions)(mfctx);CHKERRQ(ierr);
663e884886fSBarry Smith   }
664e884886fSBarry Smith   ierr = PetscOptionsEnd();CHKERRQ(ierr);
665e884886fSBarry Smith   PetscFunctionReturn(0);
666e884886fSBarry Smith }
667e884886fSBarry Smith 
668bc0f08ceSBarry Smith #undef __FUNCT__
669bc0f08ceSBarry Smith #define __FUNCT__ "MatMFFDSetPeriod_MFFD"
670bc0f08ceSBarry Smith PetscErrorCode  MatMFFDSetPeriod_MFFD(Mat mat,PetscInt period)
671bc0f08ceSBarry Smith {
672bc0f08ceSBarry Smith   MatMFFD ctx = (MatMFFD)mat->data;
673bc0f08ceSBarry Smith 
674bc0f08ceSBarry Smith   PetscFunctionBegin;
675bc0f08ceSBarry Smith   PetscValidLogicalCollectiveInt(mat,period,2);
676bc0f08ceSBarry Smith   ctx->recomputeperiod = period;
677bc0f08ceSBarry Smith   PetscFunctionReturn(0);
678bc0f08ceSBarry Smith }
679bc0f08ceSBarry Smith 
680bc0f08ceSBarry Smith #undef __FUNCT__
681bc0f08ceSBarry Smith #define __FUNCT__ "MatMFFDSetFunction_MFFD"
682bc0f08ceSBarry Smith PetscErrorCode  MatMFFDSetFunction_MFFD(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx)
683bc0f08ceSBarry Smith {
684bc0f08ceSBarry Smith   MatMFFD ctx = (MatMFFD)mat->data;
685bc0f08ceSBarry Smith 
686bc0f08ceSBarry Smith   PetscFunctionBegin;
687bc0f08ceSBarry Smith   ctx->func    = func;
688bc0f08ceSBarry Smith   ctx->funcctx = funcctx;
689bc0f08ceSBarry Smith   PetscFunctionReturn(0);
690bc0f08ceSBarry Smith }
691bc0f08ceSBarry Smith 
692bc0f08ceSBarry Smith #undef __FUNCT__
693bc0f08ceSBarry Smith #define __FUNCT__ "MatMFFDSetFunctionError_MFFD"
694bc0f08ceSBarry Smith PetscErrorCode  MatMFFDSetFunctionError_MFFD(Mat mat,PetscReal error)
695bc0f08ceSBarry Smith {
696bc0f08ceSBarry Smith   MatMFFD ctx = (MatMFFD)mat->data;
697bc0f08ceSBarry Smith 
698bc0f08ceSBarry Smith   PetscFunctionBegin;
699bc0f08ceSBarry Smith   PetscValidLogicalCollectiveReal(mat,error,2);
700bc0f08ceSBarry Smith   if (error != PETSC_DEFAULT) ctx->error_rel = error;
701bc0f08ceSBarry Smith   PetscFunctionReturn(0);
702bc0f08ceSBarry Smith }
703bc0f08ceSBarry Smith 
704e884886fSBarry Smith /*MC
705e884886fSBarry Smith   MATMFFD - MATMFFD = "mffd" - A matrix free matrix type.
706e884886fSBarry Smith 
707e884886fSBarry Smith   Level: advanced
708e884886fSBarry Smith 
709d2d6cebeSBarry Smith .seealso: MatCreateMFFD(), MatCreateSNESMF(), MatMFFDSetFunction()
710e884886fSBarry Smith M*/
711e884886fSBarry Smith #undef __FUNCT__
712e884886fSBarry Smith #define __FUNCT__ "MatCreate_MFFD"
7138cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_MFFD(Mat A)
714e884886fSBarry Smith {
715e884886fSBarry Smith   MatMFFD        mfctx;
716e884886fSBarry Smith   PetscErrorCode ierr;
717e884886fSBarry Smith 
718e884886fSBarry Smith   PetscFunctionBegin;
719519f805aSKarl Rupp #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
720607a6623SBarry Smith   ierr = MatMFFDInitializePackage();CHKERRQ(ierr);
7213ec795f1SBarry Smith #endif
722e884886fSBarry Smith 
723ce94432eSBarry Smith   ierr = PetscHeaderCreate(mfctx,_p_MatMFFD,struct _MFOps,MATMFFD_CLASSID,"MatMFFD","Matrix-free Finite Differencing","Mat",PetscObjectComm((PetscObject)A),MatDestroy_MFFD,MatView_MFFD);CHKERRQ(ierr);
7242205254eSKarl Rupp 
725e884886fSBarry Smith   mfctx->sp                       = 0;
726e884886fSBarry Smith   mfctx->error_rel                = PETSC_SQRT_MACHINE_EPSILON;
727e884886fSBarry Smith   mfctx->recomputeperiod          = 1;
728e884886fSBarry Smith   mfctx->count                    = 0;
729e884886fSBarry Smith   mfctx->currenth                 = 0.0;
7300298fd71SBarry Smith   mfctx->historyh                 = NULL;
731e884886fSBarry Smith   mfctx->ncurrenth                = 0;
732e884886fSBarry Smith   mfctx->maxcurrenth              = 0;
7337adad957SLisandro Dalcin   ((PetscObject)mfctx)->type_name = 0;
734e884886fSBarry Smith 
735e884886fSBarry Smith   mfctx->vshift = 0.0;
736e884886fSBarry Smith   mfctx->vscale = 1.0;
737e884886fSBarry Smith 
738e884886fSBarry Smith   /*
739e884886fSBarry Smith      Create the empty data structure to contain compute-h routines.
740e884886fSBarry Smith      These will be filled in below from the command line options or
741e884886fSBarry Smith      a later call with MatMFFDSetType() or if that is not called
742e884886fSBarry Smith      then it will default in the first use of MatMult_MFFD()
743e884886fSBarry Smith   */
744e884886fSBarry Smith   mfctx->ops->compute        = 0;
745e884886fSBarry Smith   mfctx->ops->destroy        = 0;
746e884886fSBarry Smith   mfctx->ops->view           = 0;
747e884886fSBarry Smith   mfctx->ops->setfromoptions = 0;
748e884886fSBarry Smith   mfctx->hctx                = 0;
749e884886fSBarry Smith 
750e884886fSBarry Smith   mfctx->func    = 0;
751e884886fSBarry Smith   mfctx->funcctx = 0;
7520298fd71SBarry Smith   mfctx->w       = NULL;
753e884886fSBarry Smith 
754e884886fSBarry Smith   A->data = mfctx;
755e884886fSBarry Smith 
756e884886fSBarry Smith   A->ops->mult           = MatMult_MFFD;
757e884886fSBarry Smith   A->ops->destroy        = MatDestroy_MFFD;
758e884886fSBarry Smith   A->ops->view           = MatView_MFFD;
759e884886fSBarry Smith   A->ops->assemblyend    = MatAssemblyEnd_MFFD;
760e884886fSBarry Smith   A->ops->getdiagonal    = MatGetDiagonal_MFFD;
761e884886fSBarry Smith   A->ops->scale          = MatScale_MFFD;
762e884886fSBarry Smith   A->ops->shift          = MatShift_MFFD;
763c3bb7e23SBarry Smith   A->ops->diagonalscale  = MatDiagonalScale_MFFD;
764c3bb7e23SBarry Smith   A->ops->diagonalset    = MatDiagonalSet_MFFD;
7659c6ac3b3SBarry Smith   A->ops->setfromoptions = MatSetFromOptions_MFFD;
766e884886fSBarry Smith   A->assembled           = PETSC_TRUE;
767e884886fSBarry Smith 
76826283091SBarry Smith   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
76926283091SBarry Smith   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
770ee1cef2cSJed Brown 
771bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetBase_C",MatMFFDSetBase_MFFD);CHKERRQ(ierr);
772bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctioniBase_C",MatMFFDSetFunctioniBase_MFFD);CHKERRQ(ierr);
773bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctioni_C",MatMFFDSetFunctioni_MFFD);CHKERRQ(ierr);
774bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunction_C",MatMFFDSetFunction_MFFD);CHKERRQ(ierr);
775bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetCheckh_C",MatMFFDSetCheckh_MFFD);CHKERRQ(ierr);
776bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetPeriod_C",MatMFFDSetPeriod_MFFD);CHKERRQ(ierr);
777bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctionError_C",MatMFFDSetFunctionError_MFFD);CHKERRQ(ierr);
778bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDResetHHistory_C",MatMFFDResetHHistory_MFFD);CHKERRQ(ierr);
779bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDAddNullSpace_C",MatMFFDAddNullSpace_MFFD);CHKERRQ(ierr);
7802205254eSKarl Rupp 
781e884886fSBarry Smith   mfctx->mat = A;
7822205254eSKarl Rupp 
783e884886fSBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,MATMFFD);CHKERRQ(ierr);
784e884886fSBarry Smith   PetscFunctionReturn(0);
785e884886fSBarry Smith }
786e884886fSBarry Smith 
787e884886fSBarry Smith #undef __FUNCT__
788e884886fSBarry Smith #define __FUNCT__ "MatCreateMFFD"
789e884886fSBarry Smith /*@
790e884886fSBarry Smith    MatCreateMFFD - Creates a matrix-free matrix. See also MatCreateSNESMF()
791e884886fSBarry Smith 
792e884886fSBarry Smith    Collective on Vec
793e884886fSBarry Smith 
794e884886fSBarry Smith    Input Parameters:
795fef1beadSBarry Smith +  comm - MPI communicator
796fef1beadSBarry Smith .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
797fef1beadSBarry Smith            This value should be the same as the local size used in creating the
798fef1beadSBarry Smith            y vector for the matrix-vector product y = Ax.
799fef1beadSBarry Smith .  n - This value should be the same as the local size used in creating the
800fef1beadSBarry Smith        x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have
801fef1beadSBarry Smith        calculated if N is given) For square matrices n is almost always m.
802fef1beadSBarry Smith .  M - number of global rows (or PETSC_DETERMINE to have calculated if m is given)
803fef1beadSBarry Smith -  N - number of global columns (or PETSC_DETERMINE to have calculated if n is given)
804fef1beadSBarry Smith 
805e884886fSBarry Smith 
806e884886fSBarry Smith    Output Parameter:
807e884886fSBarry Smith .  J - the matrix-free matrix
808e884886fSBarry Smith 
8099c6ac3b3SBarry Smith    Options Database Keys: call MatSetFromOptions() to trigger these
8109c6ac3b3SBarry Smith +  -mat_mffd_type - wp or ds (see MATMFFD_WP or MATMFFD_DS)
8119c6ac3b3SBarry Smith -  -mat_mffd_err - square root of estimated relative error in function evaluation
8129c6ac3b3SBarry Smith -  -mat_mffd_period - how often h is recomputed, defaults to 1, everytime
8139c6ac3b3SBarry Smith 
8149c6ac3b3SBarry Smith 
815e884886fSBarry Smith    Level: advanced
816e884886fSBarry Smith 
817e884886fSBarry Smith    Notes:
818e884886fSBarry Smith    The matrix-free matrix context merely contains the function pointers
819e884886fSBarry Smith    and work space for performing finite difference approximations of
820e884886fSBarry Smith    Jacobian-vector products, F'(u)*a,
821e884886fSBarry Smith 
822e884886fSBarry Smith    The default code uses the following approach to compute h
823e884886fSBarry Smith 
824e884886fSBarry Smith .vb
825e884886fSBarry Smith      F'(u)*a = [F(u+h*a) - F(u)]/h where
826e884886fSBarry Smith      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
827e884886fSBarry Smith        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
828e884886fSBarry Smith  where
829e884886fSBarry Smith      error_rel = square root of relative error in function evaluation
830e884886fSBarry Smith      umin = minimum iterate parameter
831e884886fSBarry Smith .ve
832e884886fSBarry Smith 
833e884886fSBarry Smith    The user can set the error_rel via MatMFFDSetFunctionError() and
8340598bfebSBarry Smith    umin via MatMFFDDSSetUmin(); see the <A href="../../docs/manual.pdf#nameddest=ch_snes">SNES chapter of the users manual</A> for details.
835e884886fSBarry Smith 
836e884886fSBarry Smith    The user should call MatDestroy() when finished with the matrix-free
837e884886fSBarry Smith    matrix context.
838e884886fSBarry Smith 
839e884886fSBarry Smith    Options Database Keys:
840e884886fSBarry Smith +  -mat_mffd_err <error_rel> - Sets error_rel
841e884886fSBarry Smith .  -mat_mffd_unim <umin> - Sets umin (for default PETSc routine that computes h only)
842e884886fSBarry Smith -  -mat_mffd_check_positivity
843e884886fSBarry Smith 
844e884886fSBarry Smith .keywords: default, matrix-free, create, matrix
845e884886fSBarry Smith 
8461d0fab5eSBarry Smith .seealso: MatDestroy(), MatMFFDSetFunctionError(), MatMFFDDSSetUmin(), MatMFFDSetFunction()
847e884886fSBarry Smith           MatMFFDSetHHistory(), MatMFFDResetHHistory(), MatCreateSNESMF(),
84881242352SJed Brown           MatMFFDGetH(), MatMFFDRegister(), MatMFFDComputeJacobian()
849e884886fSBarry Smith 
850e884886fSBarry Smith @*/
8517087cfbeSBarry Smith PetscErrorCode  MatCreateMFFD(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *J)
852e884886fSBarry Smith {
853e884886fSBarry Smith   PetscErrorCode ierr;
854e884886fSBarry Smith 
855e884886fSBarry Smith   PetscFunctionBegin;
856e884886fSBarry Smith   ierr = MatCreate(comm,J);CHKERRQ(ierr);
857fef1beadSBarry Smith   ierr = MatSetSizes(*J,m,n,M,N);CHKERRQ(ierr);
858e884886fSBarry Smith   ierr = MatSetType(*J,MATMFFD);CHKERRQ(ierr);
859be7c243fSJed Brown   ierr = MatSetUp(*J);CHKERRQ(ierr);
860e884886fSBarry Smith   PetscFunctionReturn(0);
861e884886fSBarry Smith }
862e884886fSBarry Smith 
863e884886fSBarry Smith 
864e884886fSBarry Smith #undef __FUNCT__
865e884886fSBarry Smith #define __FUNCT__ "MatMFFDGetH"
866e884886fSBarry Smith /*@
867e884886fSBarry Smith    MatMFFDGetH - Gets the last value that was used as the differencing
868e884886fSBarry Smith    parameter.
869e884886fSBarry Smith 
870e884886fSBarry Smith    Not Collective
871e884886fSBarry Smith 
872e884886fSBarry Smith    Input Parameters:
873e884886fSBarry Smith .  mat - the matrix obtained with MatCreateSNESMF()
874e884886fSBarry Smith 
875e884886fSBarry Smith    Output Paramter:
876e884886fSBarry Smith .  h - the differencing step size
877e884886fSBarry Smith 
878e884886fSBarry Smith    Level: advanced
879e884886fSBarry Smith 
880e884886fSBarry Smith .keywords: SNES, matrix-free, parameters
881e884886fSBarry Smith 
8821d0fab5eSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDSetHHistory(), MatCreateMFFD(), MATMFFD, MatMFFDResetHHistory()
883e884886fSBarry Smith @*/
8847087cfbeSBarry Smith PetscErrorCode  MatMFFDGetH(Mat mat,PetscScalar *h)
885e884886fSBarry Smith {
886e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)mat->data;
887bc0f08ceSBarry Smith   PetscErrorCode ierr;
888bc0f08ceSBarry Smith   PetscBool      match;
889e884886fSBarry Smith 
890e884886fSBarry Smith   PetscFunctionBegin;
891251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)mat,MATMFFD,&match);CHKERRQ(ierr);
892ce94432eSBarry Smith   if (!match) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONG,"Not a MFFD matrix");
893bc0f08ceSBarry Smith 
894e884886fSBarry Smith   *h = ctx->currenth;
895e884886fSBarry Smith   PetscFunctionReturn(0);
896e884886fSBarry Smith }
897e884886fSBarry Smith 
898e884886fSBarry Smith #undef __FUNCT__
899e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetFunction"
900e884886fSBarry Smith /*@C
901e884886fSBarry Smith    MatMFFDSetFunction - Sets the function used in applying the matrix free.
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 -  funcctx - optional function context passed to function
909e884886fSBarry Smith 
910e884886fSBarry Smith    Level: advanced
911e884886fSBarry Smith 
912e884886fSBarry Smith    Notes:
913e884886fSBarry Smith     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
914e884886fSBarry Smith     matrix inside your compute Jacobian routine
915e884886fSBarry Smith 
9163ec795f1SBarry Smith     If this is not set then it will use the function set with SNESSetFunction() if MatCreateSNESMF() was used.
917e884886fSBarry Smith 
918e884886fSBarry Smith .keywords: SNES, matrix-free, function
919e884886fSBarry Smith 
9201d0fab5eSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD,
9211d0fab5eSBarry Smith           MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction()
922e884886fSBarry Smith @*/
9237087cfbeSBarry Smith PetscErrorCode  MatMFFDSetFunction(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx)
924e884886fSBarry Smith {
925bc0f08ceSBarry Smith   PetscErrorCode ierr;
926e884886fSBarry Smith 
927e884886fSBarry Smith   PetscFunctionBegin;
928bc0f08ceSBarry Smith   ierr = PetscTryMethod(mat,"MatMFFDSetFunction_C",(Mat,PetscErrorCode (*)(void*,Vec,Vec),void*),(mat,func,funcctx));CHKERRQ(ierr);
929e884886fSBarry Smith   PetscFunctionReturn(0);
930e884886fSBarry Smith }
931e884886fSBarry Smith 
932e884886fSBarry Smith #undef __FUNCT__
933e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetFunctioni"
934e884886fSBarry Smith /*@C
935e884886fSBarry Smith    MatMFFDSetFunctioni - Sets the function for a single component
936e884886fSBarry Smith 
9373f9fe445SBarry Smith    Logically Collective on Mat
938e884886fSBarry Smith 
939e884886fSBarry Smith    Input Parameters:
940e884886fSBarry Smith +  mat - the matrix free matrix created via MatCreateSNESMF()
941e884886fSBarry Smith -  funci - the function to use
942e884886fSBarry Smith 
943e884886fSBarry Smith    Level: advanced
944e884886fSBarry Smith 
945e884886fSBarry Smith    Notes:
946e884886fSBarry Smith     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
947e884886fSBarry Smith     matrix inside your compute Jacobian routine
948e884886fSBarry Smith 
949e884886fSBarry Smith 
950e884886fSBarry Smith .keywords: SNES, matrix-free, function
951e884886fSBarry Smith 
9521d0fab5eSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction()
9531d0fab5eSBarry Smith 
954e884886fSBarry Smith @*/
9557087cfbeSBarry Smith PetscErrorCode  MatMFFDSetFunctioni(Mat mat,PetscErrorCode (*funci)(void*,PetscInt,Vec,PetscScalar*))
956e884886fSBarry Smith {
9574ac538c5SBarry Smith   PetscErrorCode ierr;
958e884886fSBarry Smith 
959e884886fSBarry Smith   PetscFunctionBegin;
9600700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
9614ac538c5SBarry Smith   ierr = PetscTryMethod(mat,"MatMFFDSetFunctioni_C",(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)),(mat,funci));CHKERRQ(ierr);
962e884886fSBarry Smith   PetscFunctionReturn(0);
963e884886fSBarry Smith }
964e884886fSBarry Smith 
965e884886fSBarry Smith 
966e884886fSBarry Smith #undef __FUNCT__
967e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetFunctioniBase"
968e884886fSBarry Smith /*@C
969e884886fSBarry Smith    MatMFFDSetFunctioniBase - Sets the base vector for a single component function evaluation
970e884886fSBarry Smith 
9713f9fe445SBarry Smith    Logically Collective on Mat
972e884886fSBarry Smith 
973e884886fSBarry Smith    Input Parameters:
974e884886fSBarry Smith +  mat - the matrix free matrix created via MatCreateSNESMF()
975e884886fSBarry Smith -  func - the function to use
976e884886fSBarry Smith 
977e884886fSBarry Smith    Level: advanced
978e884886fSBarry Smith 
979e884886fSBarry Smith    Notes:
980e884886fSBarry Smith     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
981e884886fSBarry Smith     matrix inside your compute Jacobian routine
982e884886fSBarry Smith 
983e884886fSBarry Smith 
984e884886fSBarry Smith .keywords: SNES, matrix-free, function
985e884886fSBarry Smith 
986e884886fSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD
9871d0fab5eSBarry Smith           MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction()
988e884886fSBarry Smith @*/
9897087cfbeSBarry Smith PetscErrorCode  MatMFFDSetFunctioniBase(Mat mat,PetscErrorCode (*func)(void*,Vec))
990e884886fSBarry Smith {
9914ac538c5SBarry Smith   PetscErrorCode ierr;
992e884886fSBarry Smith 
993e884886fSBarry Smith   PetscFunctionBegin;
9940700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
9954ac538c5SBarry Smith   ierr = PetscTryMethod(mat,"MatMFFDSetFunctioniBase_C",(Mat,PetscErrorCode (*)(void*,Vec)),(mat,func));CHKERRQ(ierr);
996e884886fSBarry Smith   PetscFunctionReturn(0);
997e884886fSBarry Smith }
998e884886fSBarry Smith 
999e884886fSBarry Smith #undef __FUNCT__
1000e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetPeriod"
1001e884886fSBarry Smith /*@
1002e884886fSBarry Smith    MatMFFDSetPeriod - Sets how often h is recomputed, by default it is everytime
1003e884886fSBarry Smith 
10043f9fe445SBarry Smith    Logically Collective on Mat
1005e884886fSBarry Smith 
1006e884886fSBarry Smith    Input Parameters:
1007e884886fSBarry Smith +  mat - the matrix free matrix created via MatCreateSNESMF()
1008e884886fSBarry Smith -  period - 1 for everytime, 2 for every second etc
1009e884886fSBarry Smith 
1010e884886fSBarry Smith    Options Database Keys:
1011e884886fSBarry Smith +  -mat_mffd_period <period>
1012e884886fSBarry Smith 
1013e884886fSBarry Smith    Level: advanced
1014e884886fSBarry Smith 
1015e884886fSBarry Smith 
1016e884886fSBarry Smith .keywords: SNES, matrix-free, parameters
1017e884886fSBarry Smith 
1018e884886fSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(),
10191d0fab5eSBarry Smith           MatMFFDSetHHistory(), MatMFFDResetHHistory()
1020e884886fSBarry Smith @*/
10217087cfbeSBarry Smith PetscErrorCode  MatMFFDSetPeriod(Mat mat,PetscInt period)
1022e884886fSBarry Smith {
1023bc0f08ceSBarry Smith   PetscErrorCode ierr;
1024e884886fSBarry Smith 
1025e884886fSBarry Smith   PetscFunctionBegin;
1026bc0f08ceSBarry Smith   ierr = PetscTryMethod(mat,"MatMFFDSetPeriod_C",(Mat,PetscInt),(mat,period));CHKERRQ(ierr);
1027e884886fSBarry Smith   PetscFunctionReturn(0);
1028e884886fSBarry Smith }
1029e884886fSBarry Smith 
1030e884886fSBarry Smith #undef __FUNCT__
1031e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetFunctionError"
1032e884886fSBarry Smith /*@
1033e884886fSBarry Smith    MatMFFDSetFunctionError - Sets the error_rel for the approximation of
1034e884886fSBarry Smith    matrix-vector products using finite differences.
1035e884886fSBarry Smith 
10363f9fe445SBarry Smith    Logically Collective on Mat
1037e884886fSBarry Smith 
1038e884886fSBarry Smith    Input Parameters:
1039e884886fSBarry Smith +  mat - the matrix free matrix created via MatCreateMFFD() or MatCreateSNESMF()
1040e884886fSBarry Smith -  error_rel - relative error (should be set to the square root of
1041e884886fSBarry Smith                the relative error in the function evaluations)
1042e884886fSBarry Smith 
1043e884886fSBarry Smith    Options Database Keys:
1044e884886fSBarry Smith +  -mat_mffd_err <error_rel> - Sets error_rel
1045e884886fSBarry Smith 
1046e884886fSBarry Smith    Level: advanced
1047e884886fSBarry Smith 
1048e884886fSBarry Smith    Notes:
1049e884886fSBarry Smith    The default matrix-free matrix-vector product routine computes
1050e884886fSBarry Smith .vb
1051e884886fSBarry Smith      F'(u)*a = [F(u+h*a) - F(u)]/h where
1052e884886fSBarry Smith      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
1053e884886fSBarry Smith        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   else
1054e884886fSBarry Smith .ve
1055e884886fSBarry Smith 
1056e884886fSBarry Smith .keywords: SNES, matrix-free, parameters
1057e884886fSBarry Smith 
1058e884886fSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD
10591d0fab5eSBarry Smith           MatMFFDSetHHistory(), MatMFFDResetHHistory()
1060e884886fSBarry Smith @*/
10617087cfbeSBarry Smith PetscErrorCode  MatMFFDSetFunctionError(Mat mat,PetscReal error)
1062e884886fSBarry Smith {
1063bc0f08ceSBarry Smith   PetscErrorCode ierr;
1064e884886fSBarry Smith 
1065e884886fSBarry Smith   PetscFunctionBegin;
1066bc0f08ceSBarry Smith   ierr = PetscTryMethod(mat,"MatMFFDSetFunctionError_C",(Mat,PetscReal),(mat,error));CHKERRQ(ierr);
1067e884886fSBarry Smith   PetscFunctionReturn(0);
1068e884886fSBarry Smith }
1069e884886fSBarry Smith 
1070e884886fSBarry Smith #undef __FUNCT__
1071e884886fSBarry Smith #define __FUNCT__ "MatMFFDAddNullSpace"
1072e884886fSBarry Smith /*@
1073e884886fSBarry Smith    MatMFFDAddNullSpace - Provides a null space that an operator is
1074e884886fSBarry Smith    supposed to have.  Since roundoff will create a small component in
1075e884886fSBarry Smith    the null space, if you know the null space you may have it
1076e884886fSBarry Smith    automatically removed.
1077e884886fSBarry Smith 
10783f9fe445SBarry Smith    Logically Collective on Mat
1079e884886fSBarry Smith 
1080e884886fSBarry Smith    Input Parameters:
1081e884886fSBarry Smith +  J - the matrix-free matrix context
1082e884886fSBarry Smith -  nullsp - object created with MatNullSpaceCreate()
1083e884886fSBarry Smith 
1084e884886fSBarry Smith    Level: advanced
1085e884886fSBarry Smith 
1086e884886fSBarry Smith .keywords: SNES, matrix-free, null space
1087e884886fSBarry Smith 
1088e884886fSBarry Smith .seealso: MatNullSpaceCreate(), MatMFFDGetH(), MatCreateSNESMF(), MatCreateMFFD(), MATMFFD
10891d0fab5eSBarry Smith           MatMFFDSetHHistory(), MatMFFDResetHHistory()
1090e884886fSBarry Smith @*/
10917087cfbeSBarry Smith PetscErrorCode  MatMFFDAddNullSpace(Mat J,MatNullSpace nullsp)
1092e884886fSBarry Smith {
1093e884886fSBarry Smith   PetscErrorCode ierr;
1094e884886fSBarry Smith 
1095e884886fSBarry Smith   PetscFunctionBegin;
1096bc0f08ceSBarry Smith   ierr = PetscTryMethod(J,"MatMFFDAddNullSpace_C",(Mat,MatNullSpace),(J,nullsp));CHKERRQ(ierr);
1097e884886fSBarry Smith   PetscFunctionReturn(0);
1098e884886fSBarry Smith }
1099e884886fSBarry Smith 
1100e884886fSBarry Smith #undef __FUNCT__
1101e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetHHistory"
1102e884886fSBarry Smith /*@
1103e884886fSBarry Smith    MatMFFDSetHHistory - Sets an array to collect a history of the
1104e884886fSBarry Smith    differencing values (h) computed for the matrix-free product.
1105e884886fSBarry Smith 
11063f9fe445SBarry Smith    Logically Collective on Mat
1107e884886fSBarry Smith 
1108e884886fSBarry Smith    Input Parameters:
1109e884886fSBarry Smith +  J - the matrix-free matrix context
1110e884886fSBarry Smith .  histroy - space to hold the history
1111e884886fSBarry Smith -  nhistory - number of entries in history, if more entries are generated than
1112e884886fSBarry Smith               nhistory, then the later ones are discarded
1113e884886fSBarry Smith 
1114e884886fSBarry Smith    Level: advanced
1115e884886fSBarry Smith 
1116e884886fSBarry Smith    Notes:
1117e884886fSBarry Smith    Use MatMFFDResetHHistory() to reset the history counter and collect
1118e884886fSBarry Smith    a new batch of differencing parameters, h.
1119e884886fSBarry Smith 
1120e884886fSBarry Smith .keywords: SNES, matrix-free, h history, differencing history
1121e884886fSBarry Smith 
1122e884886fSBarry Smith .seealso: MatMFFDGetH(), MatCreateSNESMF(),
11231d0fab5eSBarry Smith           MatMFFDResetHHistory(), MatMFFDSetFunctionError()
1124e884886fSBarry Smith 
1125e884886fSBarry Smith @*/
11267087cfbeSBarry Smith PetscErrorCode  MatMFFDSetHHistory(Mat J,PetscScalar history[],PetscInt nhistory)
1127e884886fSBarry Smith {
1128e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)J->data;
1129bc0f08ceSBarry Smith   PetscErrorCode ierr;
1130bc0f08ceSBarry Smith   PetscBool      match;
1131e884886fSBarry Smith 
1132e884886fSBarry Smith   PetscFunctionBegin;
1133251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)J,MATMFFD,&match);CHKERRQ(ierr);
1134ce94432eSBarry Smith   if (!match) SETERRQ(PetscObjectComm((PetscObject)J),PETSC_ERR_ARG_WRONG,"Not a MFFD matrix");
1135e884886fSBarry Smith   ctx->historyh    = history;
1136e884886fSBarry Smith   ctx->maxcurrenth = nhistory;
113775567043SBarry Smith   ctx->currenth    = 0.;
1138e884886fSBarry Smith   PetscFunctionReturn(0);
1139e884886fSBarry Smith }
1140e884886fSBarry Smith 
1141bc0f08ceSBarry Smith 
1142e884886fSBarry Smith #undef __FUNCT__
1143e884886fSBarry Smith #define __FUNCT__ "MatMFFDResetHHistory"
1144e884886fSBarry Smith /*@
1145e884886fSBarry Smith    MatMFFDResetHHistory - Resets the counter to zero to begin
1146e884886fSBarry Smith    collecting a new set of differencing histories.
1147e884886fSBarry Smith 
11483f9fe445SBarry Smith    Logically Collective on Mat
1149e884886fSBarry Smith 
1150e884886fSBarry Smith    Input Parameters:
1151e884886fSBarry Smith .  J - the matrix-free matrix context
1152e884886fSBarry Smith 
1153e884886fSBarry Smith    Level: advanced
1154e884886fSBarry Smith 
1155e884886fSBarry Smith    Notes:
1156e884886fSBarry Smith    Use MatMFFDSetHHistory() to create the original history counter.
1157e884886fSBarry Smith 
1158e884886fSBarry Smith .keywords: SNES, matrix-free, h history, differencing history
1159e884886fSBarry Smith 
1160e884886fSBarry Smith .seealso: MatMFFDGetH(), MatCreateSNESMF(),
11611d0fab5eSBarry Smith           MatMFFDSetHHistory(), MatMFFDSetFunctionError()
1162e884886fSBarry Smith 
1163e884886fSBarry Smith @*/
11647087cfbeSBarry Smith PetscErrorCode  MatMFFDResetHHistory(Mat J)
1165e884886fSBarry Smith {
1166bc0f08ceSBarry Smith   PetscErrorCode ierr;
1167e884886fSBarry Smith 
1168e884886fSBarry Smith   PetscFunctionBegin;
1169bc0f08ceSBarry Smith   ierr = PetscTryMethod(J,"MatMFFDResetHHistory_C",(Mat),(J));CHKERRQ(ierr);
1170e884886fSBarry Smith   PetscFunctionReturn(0);
1171e884886fSBarry Smith }
1172e884886fSBarry Smith 
1173e884886fSBarry Smith 
1174e884886fSBarry Smith #undef __FUNCT__
1175e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetBase"
1176e884886fSBarry Smith /*@
1177e884886fSBarry Smith     MatMFFDSetBase - Sets the vector U at which matrix vector products of the
1178e884886fSBarry Smith         Jacobian are computed
1179e884886fSBarry Smith 
11803f9fe445SBarry Smith     Logically Collective on Mat
1181e884886fSBarry Smith 
1182e884886fSBarry Smith     Input Parameters:
1183e884886fSBarry Smith +   J - the MatMFFD matrix
11843ec795f1SBarry Smith .   U - the vector
1185bcddec3dSBarry Smith -   F - (optional) vector that contains F(u) if it has been already computed
1186e884886fSBarry Smith 
1187e884886fSBarry Smith     Notes: This is rarely used directly
1188e884886fSBarry Smith 
11898af5ae88SBarry Smith     If F is provided then it is not recomputed. Otherwise the function is evaluated at the base
11908af5ae88SBarry Smith     point during the first MatMult() after each call to MatMFFDSetBase().
1191dff2f722SBarry Smith 
1192e884886fSBarry Smith     Level: advanced
1193e884886fSBarry Smith 
1194e884886fSBarry Smith @*/
11957087cfbeSBarry Smith PetscErrorCode  MatMFFDSetBase(Mat J,Vec U,Vec F)
1196e884886fSBarry Smith {
11974ac538c5SBarry Smith   PetscErrorCode ierr;
1198e884886fSBarry Smith 
1199e884886fSBarry Smith   PetscFunctionBegin;
12000700a824SBarry Smith   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
12010700a824SBarry Smith   PetscValidHeaderSpecific(U,VEC_CLASSID,2);
12020700a824SBarry Smith   if (F) PetscValidHeaderSpecific(F,VEC_CLASSID,3);
12034ac538c5SBarry Smith   ierr = PetscTryMethod(J,"MatMFFDSetBase_C",(Mat,Vec,Vec),(J,U,F));CHKERRQ(ierr);
1204e884886fSBarry Smith   PetscFunctionReturn(0);
1205e884886fSBarry Smith }
1206e884886fSBarry Smith 
1207e884886fSBarry Smith #undef __FUNCT__
1208e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetCheckh"
1209e884886fSBarry Smith /*@C
1210e884886fSBarry Smith     MatMFFDSetCheckh - Sets a function that checks the computed h and adjusts
1211e884886fSBarry Smith         it to satisfy some criteria
1212e884886fSBarry Smith 
12133f9fe445SBarry Smith     Logically Collective on Mat
1214e884886fSBarry Smith 
1215e884886fSBarry Smith     Input Parameters:
1216e884886fSBarry Smith +   J - the MatMFFD matrix
1217e884886fSBarry Smith .   fun - the function that checks h
1218e884886fSBarry Smith -   ctx - any context needed by the function
1219e884886fSBarry Smith 
1220e884886fSBarry Smith     Options Database Keys:
1221e884886fSBarry Smith .   -mat_mffd_check_positivity
1222e884886fSBarry Smith 
1223e884886fSBarry Smith     Level: advanced
1224e884886fSBarry Smith 
1225e884886fSBarry Smith     Notes: For example, MatMFFDSetCheckPositivity() insures that all entries
1226e884886fSBarry Smith        of U + h*a are non-negative
1227e884886fSBarry Smith 
1228e884886fSBarry Smith .seealso:  MatMFFDSetCheckPositivity()
1229e884886fSBarry Smith @*/
12307087cfbeSBarry Smith PetscErrorCode  MatMFFDSetCheckh(Mat J,PetscErrorCode (*fun)(void*,Vec,Vec,PetscScalar*),void *ctx)
1231e884886fSBarry Smith {
12324ac538c5SBarry Smith   PetscErrorCode ierr;
1233e884886fSBarry Smith 
1234e884886fSBarry Smith   PetscFunctionBegin;
12350700a824SBarry Smith   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
12364ac538c5SBarry Smith   ierr = PetscTryMethod(J,"MatMFFDSetCheckh_C",(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*),(J,fun,ctx));CHKERRQ(ierr);
1237e884886fSBarry Smith   PetscFunctionReturn(0);
1238e884886fSBarry Smith }
1239e884886fSBarry Smith 
1240e884886fSBarry Smith #undef __FUNCT__
1241e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetCheckPositivity"
1242e884886fSBarry Smith /*@
1243e884886fSBarry Smith     MatMFFDCheckPositivity - Checks that all entries in U + h*a are positive or
1244e884886fSBarry Smith         zero, decreases h until this is satisfied.
1245e884886fSBarry Smith 
12463f9fe445SBarry Smith     Logically Collective on Vec
1247e884886fSBarry Smith 
1248e884886fSBarry Smith     Input Parameters:
1249e884886fSBarry Smith +   U - base vector that is added to
1250e884886fSBarry Smith .   a - vector that is added
1251e884886fSBarry Smith .   h - scaling factor on a
1252e884886fSBarry Smith -   dummy - context variable (unused)
1253e884886fSBarry Smith 
1254e884886fSBarry Smith     Options Database Keys:
1255e884886fSBarry Smith .   -mat_mffd_check_positivity
1256e884886fSBarry Smith 
1257e884886fSBarry Smith     Level: advanced
1258e884886fSBarry Smith 
1259e884886fSBarry Smith     Notes: This is rarely used directly, rather it is passed as an argument to
1260e884886fSBarry Smith            MatMFFDSetCheckh()
1261e884886fSBarry Smith 
1262e884886fSBarry Smith .seealso:  MatMFFDSetCheckh()
1263e884886fSBarry Smith @*/
12647087cfbeSBarry Smith PetscErrorCode  MatMFFDCheckPositivity(void *dummy,Vec U,Vec a,PetscScalar *h)
1265e884886fSBarry Smith {
1266e884886fSBarry Smith   PetscReal      val, minval;
1267e884886fSBarry Smith   PetscScalar    *u_vec, *a_vec;
1268e884886fSBarry Smith   PetscErrorCode ierr;
1269e884886fSBarry Smith   PetscInt       i,n;
1270e884886fSBarry Smith   MPI_Comm       comm;
1271e884886fSBarry Smith 
1272e884886fSBarry Smith   PetscFunctionBegin;
1273e884886fSBarry Smith   ierr   = PetscObjectGetComm((PetscObject)U,&comm);CHKERRQ(ierr);
1274e884886fSBarry Smith   ierr   = VecGetArray(U,&u_vec);CHKERRQ(ierr);
1275e884886fSBarry Smith   ierr   = VecGetArray(a,&a_vec);CHKERRQ(ierr);
1276e884886fSBarry Smith   ierr   = VecGetLocalSize(U,&n);CHKERRQ(ierr);
1277e884886fSBarry Smith   minval = PetscAbsScalar(*h*1.01);
1278e884886fSBarry Smith   for (i=0; i<n; i++) {
1279e884886fSBarry Smith     if (PetscRealPart(u_vec[i] + *h*a_vec[i]) <= 0.0) {
1280e884886fSBarry Smith       val = PetscAbsScalar(u_vec[i]/a_vec[i]);
1281e884886fSBarry Smith       if (val < minval) minval = val;
1282e884886fSBarry Smith     }
1283e884886fSBarry Smith   }
1284e884886fSBarry Smith   ierr = VecRestoreArray(U,&u_vec);CHKERRQ(ierr);
1285e884886fSBarry Smith   ierr = VecRestoreArray(a,&a_vec);CHKERRQ(ierr);
1286d9822059SBarry Smith   ierr = MPI_Allreduce(&minval,&val,1,MPIU_REAL,MPIU_MIN,comm);CHKERRQ(ierr);
1287e884886fSBarry Smith   if (val <= PetscAbsScalar(*h)) {
1288e884886fSBarry Smith     ierr = PetscInfo2(U,"Scaling back h from %G to %G\n",PetscRealPart(*h),.99*val);CHKERRQ(ierr);
1289e884886fSBarry Smith     if (PetscRealPart(*h) > 0.0) *h =  0.99*val;
1290e884886fSBarry Smith     else                         *h = -0.99*val;
1291e884886fSBarry Smith   }
1292e884886fSBarry Smith   PetscFunctionReturn(0);
1293e884886fSBarry Smith }
1294e884886fSBarry Smith 
1295e884886fSBarry Smith 
1296e884886fSBarry Smith 
1297e884886fSBarry Smith 
1298e884886fSBarry Smith 
1299e884886fSBarry Smith 
1300e884886fSBarry Smith 
1301