xref: /petsc/src/mat/impls/mffd/mffd.c (revision 1c9cd33768b1fd01403e37595b64fe66efc857ab)
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;
28b022a5c1SBarry Smith   MatMFFDPackageInitialized = PETSC_FALSE;
29b022a5c1SBarry Smith   MatMFFDRegisterAllCalled  = PETSC_FALSE;
30a703d84cSPatrick Lacasse   ierr = MatMFFDRegisterDestroy();CHKERRQ(ierr);
310298fd71SBarry Smith   MatMFFDList               = NULL;
32b022a5c1SBarry Smith   PetscFunctionReturn(0);
33b022a5c1SBarry Smith }
34b022a5c1SBarry Smith 
35e884886fSBarry Smith #undef __FUNCT__
363ec795f1SBarry Smith #define __FUNCT__ "MatMFFDInitializePackage"
373ec795f1SBarry Smith /*@C
383ec795f1SBarry Smith   MatMFFDInitializePackage - This function initializes everything in the MatMFFD package. It is called
393ec795f1SBarry Smith   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to MatCreate_MFFD()
403ec795f1SBarry Smith   when using static libraries.
413ec795f1SBarry Smith 
423ec795f1SBarry Smith   Level: developer
433ec795f1SBarry Smith 
443ec795f1SBarry Smith .keywords: Vec, initialize, package
453ec795f1SBarry Smith .seealso: PetscInitialize()
463ec795f1SBarry Smith @*/
47607a6623SBarry Smith PetscErrorCode  MatMFFDInitializePackage(void)
483ec795f1SBarry Smith {
493ec795f1SBarry Smith   char           logList[256];
503ec795f1SBarry Smith   char           *className;
51ace3abfcSBarry Smith   PetscBool      opt;
523ec795f1SBarry Smith   PetscErrorCode ierr;
533ec795f1SBarry Smith 
543ec795f1SBarry Smith   PetscFunctionBegin;
55b022a5c1SBarry Smith   if (MatMFFDPackageInitialized) PetscFunctionReturn(0);
56b022a5c1SBarry Smith   MatMFFDPackageInitialized = PETSC_TRUE;
573ec795f1SBarry Smith   /* Register Classes */
580700a824SBarry Smith   ierr = PetscClassIdRegister("MatMFFD",&MATMFFD_CLASSID);CHKERRQ(ierr);
593ec795f1SBarry Smith   /* Register Constructors */
60607a6623SBarry Smith   ierr = MatMFFDRegisterAll();CHKERRQ(ierr);
613ec795f1SBarry Smith   /* Register Events */
620700a824SBarry Smith   ierr = PetscLogEventRegister("MatMult MF",          MATMFFD_CLASSID,&MATMFFD_Mult);CHKERRQ(ierr);
633ec795f1SBarry Smith 
643ec795f1SBarry Smith   /* Process info exclusions */
650298fd71SBarry Smith   ierr = PetscOptionsGetString(NULL, "-info_exclude", logList, 256, &opt);CHKERRQ(ierr);
663ec795f1SBarry Smith   if (opt) {
673ec795f1SBarry Smith     ierr = PetscStrstr(logList, "matmffd", &className);CHKERRQ(ierr);
683ec795f1SBarry Smith     if (className) {
690700a824SBarry Smith       ierr = PetscInfoDeactivateClass(MATMFFD_CLASSID);CHKERRQ(ierr);
703ec795f1SBarry Smith     }
713ec795f1SBarry Smith   }
723ec795f1SBarry Smith   /* Process summary exclusions */
730298fd71SBarry Smith   ierr = PetscOptionsGetString(NULL, "-log_summary_exclude", logList, 256, &opt);CHKERRQ(ierr);
743ec795f1SBarry Smith   if (opt) {
753ec795f1SBarry Smith     ierr = PetscStrstr(logList, "matmffd", &className);CHKERRQ(ierr);
763ec795f1SBarry Smith     if (className) {
770700a824SBarry Smith       ierr = PetscLogEventDeactivateClass(MATMFFD_CLASSID);CHKERRQ(ierr);
783ec795f1SBarry Smith     }
793ec795f1SBarry Smith   }
80b022a5c1SBarry Smith   ierr = PetscRegisterFinalize(MatMFFDFinalizePackage);CHKERRQ(ierr);
813ec795f1SBarry Smith   PetscFunctionReturn(0);
823ec795f1SBarry Smith }
833ec795f1SBarry Smith 
843ec795f1SBarry Smith #undef __FUNCT__
85e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetType"
86e884886fSBarry Smith /*@C
87e884886fSBarry Smith     MatMFFDSetType - Sets the method that is used to compute the
88e884886fSBarry Smith     differencing parameter for finite differene matrix-free formulations.
89e884886fSBarry Smith 
90e884886fSBarry Smith     Input Parameters:
91e884886fSBarry Smith +   mat - the "matrix-free" matrix created via MatCreateSNESMF(), or MatCreateMFFD()
92e884886fSBarry Smith           or MatSetType(mat,MATMFFD);
93e884886fSBarry Smith -   ftype - the type requested, either MATMFFD_WP or MATMFFD_DS
94e884886fSBarry Smith 
95e884886fSBarry Smith     Level: advanced
96e884886fSBarry Smith 
97e884886fSBarry Smith     Notes:
98e884886fSBarry Smith     For example, such routines can compute h for use in
99e884886fSBarry Smith     Jacobian-vector products of the form
100e884886fSBarry Smith 
101e884886fSBarry Smith                         F(x+ha) - F(x)
102e884886fSBarry Smith           F'(u)a  ~=  ----------------
103e884886fSBarry Smith                               h
104e884886fSBarry Smith 
1051c84c290SBarry Smith .seealso: MatCreateSNESMF(), MatMFFDRegister(), MatMFFDSetFunction()
106e884886fSBarry Smith @*/
10719fd82e9SBarry Smith PetscErrorCode  MatMFFDSetType(Mat mat,MatMFFDType ftype)
108e884886fSBarry Smith {
109e884886fSBarry Smith   PetscErrorCode ierr,(*r)(MatMFFD);
110e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)mat->data;
111ace3abfcSBarry Smith   PetscBool      match;
112e884886fSBarry Smith 
113e884886fSBarry Smith   PetscFunctionBegin;
1140700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
115e884886fSBarry Smith   PetscValidCharPointer(ftype,2);
116e884886fSBarry Smith 
117251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)mat,MATMFFD,&match);CHKERRQ(ierr);
1189a13eae7SJed Brown   if (!match) PetscFunctionReturn(0);
1199a13eae7SJed Brown 
120e884886fSBarry Smith   /* already set, so just return */
121251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)ctx,ftype,&match);CHKERRQ(ierr);
122e884886fSBarry Smith   if (match) PetscFunctionReturn(0);
123e884886fSBarry Smith 
124e884886fSBarry Smith   /* destroy the old one if it exists */
125e884886fSBarry Smith   if (ctx->ops->destroy) {
126e884886fSBarry Smith     ierr = (*ctx->ops->destroy)(ctx);CHKERRQ(ierr);
127e884886fSBarry Smith   }
128e884886fSBarry Smith 
129*1c9cd337SJed Brown   ierr =  PetscFunctionListFind(MatMFFDList,ftype,&r);CHKERRQ(ierr);
130e32f2f54SBarry Smith   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown MatMFFD type %s given",ftype);
131e884886fSBarry Smith   ierr = (*r)(ctx);CHKERRQ(ierr);
132e884886fSBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)ctx,ftype);CHKERRQ(ierr);
133e884886fSBarry Smith   PetscFunctionReturn(0);
134e884886fSBarry Smith }
135e884886fSBarry Smith 
136e884886fSBarry Smith typedef PetscErrorCode (*FCN1)(void*,Vec); /* force argument to next function to not be extern C*/
137e884886fSBarry Smith #undef __FUNCT__
138c879e56bSBarry Smith #define __FUNCT__ "MatMFFDSetFunctioniBase_MFFD"
1397087cfbeSBarry Smith PetscErrorCode  MatMFFDSetFunctioniBase_MFFD(Mat mat,FCN1 func)
140e884886fSBarry Smith {
141e884886fSBarry Smith   MatMFFD ctx = (MatMFFD)mat->data;
142e884886fSBarry Smith 
143e884886fSBarry Smith   PetscFunctionBegin;
144e884886fSBarry Smith   ctx->funcisetbase = func;
145e884886fSBarry Smith   PetscFunctionReturn(0);
146e884886fSBarry Smith }
147e884886fSBarry Smith 
148e884886fSBarry Smith typedef PetscErrorCode (*FCN2)(void*,PetscInt,Vec,PetscScalar*); /* force argument to next function to not be extern C*/
149e884886fSBarry Smith #undef __FUNCT__
150c879e56bSBarry Smith #define __FUNCT__ "MatMFFDSetFunctioni_MFFD"
1517087cfbeSBarry Smith PetscErrorCode  MatMFFDSetFunctioni_MFFD(Mat mat,FCN2 funci)
152e884886fSBarry Smith {
153e884886fSBarry Smith   MatMFFD ctx = (MatMFFD)mat->data;
154e884886fSBarry Smith 
155e884886fSBarry Smith   PetscFunctionBegin;
156e884886fSBarry Smith   ctx->funci = funci;
157e884886fSBarry Smith   PetscFunctionReturn(0);
158e884886fSBarry Smith }
159e884886fSBarry Smith 
160bc0f08ceSBarry Smith #undef __FUNCT__
161bc0f08ceSBarry Smith #define __FUNCT__ "MatMFFDResetHHistory_MFFD"
162bc0f08ceSBarry Smith PetscErrorCode  MatMFFDResetHHistory_MFFD(Mat J)
163bc0f08ceSBarry Smith {
164bc0f08ceSBarry Smith   MatMFFD ctx = (MatMFFD)J->data;
165bc0f08ceSBarry Smith 
166bc0f08ceSBarry Smith   PetscFunctionBegin;
167bc0f08ceSBarry Smith   ctx->ncurrenth = 0;
168bc0f08ceSBarry Smith   PetscFunctionReturn(0);
169bc0f08ceSBarry Smith }
170e884886fSBarry Smith 
171e884886fSBarry Smith #undef __FUNCT__
172e884886fSBarry Smith #define __FUNCT__ "MatMFFDRegister"
1731c84c290SBarry Smith /*@C
1741c84c290SBarry Smith    MatMFFDRegister - Adds a method to the MatMFFD registry.
1751c84c290SBarry Smith 
1761c84c290SBarry Smith    Not Collective
1771c84c290SBarry Smith 
1781c84c290SBarry Smith    Input Parameters:
1791c84c290SBarry Smith +  name_solver - name of a new user-defined compute-h module
1801c84c290SBarry Smith -  routine_create - routine to create method context
1811c84c290SBarry Smith 
1821c84c290SBarry Smith    Level: developer
1831c84c290SBarry Smith 
1841c84c290SBarry Smith    Notes:
1851c84c290SBarry Smith    MatMFFDRegister() may be called multiple times to add several user-defined solvers.
1861c84c290SBarry Smith 
1871c84c290SBarry Smith    Sample usage:
1881c84c290SBarry Smith .vb
189bdf89e91SBarry Smith    MatMFFDRegister("my_h",MyHCreate);
1901c84c290SBarry Smith .ve
1911c84c290SBarry Smith 
1921c84c290SBarry Smith    Then, your solver can be chosen with the procedural interface via
1931c84c290SBarry Smith $     MatMFFDSetType(mfctx,"my_h")
1941c84c290SBarry Smith    or at runtime via the option
1951c84c290SBarry Smith $     -snes_mf_type my_h
1961c84c290SBarry Smith 
1971c84c290SBarry Smith .keywords: MatMFFD, register
1981c84c290SBarry Smith 
1991c84c290SBarry Smith .seealso: MatMFFDRegisterAll(), MatMFFDRegisterDestroy()
2001c84c290SBarry Smith  @*/
201bdf89e91SBarry Smith PetscErrorCode  MatMFFDRegister(const char sname[],PetscErrorCode (*function)(MatMFFD))
202e884886fSBarry Smith {
203e884886fSBarry Smith   PetscErrorCode ierr;
204e884886fSBarry Smith 
205e884886fSBarry Smith   PetscFunctionBegin;
206bdf89e91SBarry Smith   ierr = PetscFunctionListAdd(&MatMFFDList,sname,(void (*)(void))function);CHKERRQ(ierr);
207e884886fSBarry Smith   PetscFunctionReturn(0);
208e884886fSBarry Smith }
209e884886fSBarry Smith 
210e884886fSBarry Smith 
211e884886fSBarry Smith #undef __FUNCT__
212e884886fSBarry Smith #define __FUNCT__ "MatMFFDRegisterDestroy"
213e884886fSBarry Smith /*@C
214e884886fSBarry Smith    MatMFFDRegisterDestroy - Frees the list of MatMFFD methods that were
215e884886fSBarry Smith    registered by MatMFFDRegisterDynamic).
216e884886fSBarry Smith 
217e884886fSBarry Smith    Not Collective
218e884886fSBarry Smith 
219e884886fSBarry Smith    Level: developer
220e884886fSBarry Smith 
221e884886fSBarry Smith .keywords: MatMFFD, register, destroy
222e884886fSBarry Smith 
223e884886fSBarry Smith .seealso: MatMFFDRegisterDynamic), MatMFFDRegisterAll()
224e884886fSBarry Smith @*/
2257087cfbeSBarry Smith PetscErrorCode  MatMFFDRegisterDestroy(void)
226e884886fSBarry Smith {
227e884886fSBarry Smith   PetscErrorCode ierr;
228e884886fSBarry Smith 
229e884886fSBarry Smith   PetscFunctionBegin;
230140e18c1SBarry Smith   ierr = PetscFunctionListDestroy(&MatMFFDList);CHKERRQ(ierr);
2312205254eSKarl Rupp 
232e884886fSBarry Smith   MatMFFDRegisterAllCalled = PETSC_FALSE;
233e884886fSBarry Smith   PetscFunctionReturn(0);
234e884886fSBarry Smith }
235e884886fSBarry Smith 
236bc0f08ceSBarry Smith #undef __FUNCT__
237bc0f08ceSBarry Smith #define __FUNCT__ "MatMFFDAddNullSpace_MFFD"
238bc0f08ceSBarry Smith PetscErrorCode  MatMFFDAddNullSpace_MFFD(Mat J,MatNullSpace nullsp)
239bc0f08ceSBarry Smith {
240bc0f08ceSBarry Smith   PetscErrorCode ierr;
241bc0f08ceSBarry Smith   MatMFFD        ctx = (MatMFFD)J->data;
242bc0f08ceSBarry Smith 
243bc0f08ceSBarry Smith   PetscFunctionBegin;
244bc0f08ceSBarry Smith   ierr = PetscObjectReference((PetscObject)nullsp);CHKERRQ(ierr);
245bc0f08ceSBarry Smith   if (ctx->sp) { ierr = MatNullSpaceDestroy(&ctx->sp);CHKERRQ(ierr); }
246bc0f08ceSBarry Smith   ctx->sp = nullsp;
247bc0f08ceSBarry Smith   PetscFunctionReturn(0);
248bc0f08ceSBarry Smith }
249bc0f08ceSBarry Smith 
250e884886fSBarry Smith /* ----------------------------------------------------------------------------------------*/
251e884886fSBarry Smith #undef __FUNCT__
252e884886fSBarry Smith #define __FUNCT__ "MatDestroy_MFFD"
253e884886fSBarry Smith PetscErrorCode MatDestroy_MFFD(Mat mat)
254e884886fSBarry Smith {
255e884886fSBarry Smith   PetscErrorCode ierr;
256e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)mat->data;
257e884886fSBarry Smith 
258e884886fSBarry Smith   PetscFunctionBegin;
2596bf464f9SBarry Smith   ierr = VecDestroy(&ctx->w);CHKERRQ(ierr);
2606bf464f9SBarry Smith   ierr = VecDestroy(&ctx->drscale);CHKERRQ(ierr);
2616bf464f9SBarry Smith   ierr = VecDestroy(&ctx->dlscale);CHKERRQ(ierr);
2626bf464f9SBarry Smith   ierr = VecDestroy(&ctx->dshift);CHKERRQ(ierr);
263cfe22f5eSBarry Smith   if (ctx->current_f_allocated) {
2646bf464f9SBarry Smith     ierr = VecDestroy(&ctx->current_f);CHKERRQ(ierr);
265cfe22f5eSBarry Smith   }
266e884886fSBarry Smith   if (ctx->ops->destroy) {ierr = (*ctx->ops->destroy)(ctx);CHKERRQ(ierr);}
2676bf464f9SBarry Smith   ierr      = MatNullSpaceDestroy(&ctx->sp);CHKERRQ(ierr);
2686bf464f9SBarry Smith   ierr      = PetscHeaderDestroy(&ctx);CHKERRQ(ierr);
269bf0cc555SLisandro Dalcin   mat->data = 0;
270e884886fSBarry Smith 
271bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetBase_C",NULL);CHKERRQ(ierr);
272bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctioniBase_C",NULL);CHKERRQ(ierr);
273bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctioni_C",NULL);CHKERRQ(ierr);
274bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunction_C",NULL);CHKERRQ(ierr);
275bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctionError_C",NULL);CHKERRQ(ierr);
276bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetCheckh_C",NULL);CHKERRQ(ierr);
277bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetPeriod_C",NULL);CHKERRQ(ierr);
278bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDResetHHistory_C",NULL);CHKERRQ(ierr);
279bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDAddNullSpace_C",NULL);CHKERRQ(ierr);
280e884886fSBarry Smith   PetscFunctionReturn(0);
281e884886fSBarry Smith }
282e884886fSBarry Smith 
283e884886fSBarry Smith #undef __FUNCT__
284e884886fSBarry Smith #define __FUNCT__ "MatView_MFFD"
285e884886fSBarry Smith /*
286e884886fSBarry Smith    MatMFFDView_MFFD - Views matrix-free parameters.
287e884886fSBarry Smith 
288e884886fSBarry Smith */
289e884886fSBarry Smith PetscErrorCode MatView_MFFD(Mat J,PetscViewer viewer)
290e884886fSBarry Smith {
291e884886fSBarry Smith   PetscErrorCode ierr;
292e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)J->data;
293a283635eSDmitry Karpeev   PetscBool      iascii, viewbase, viewfunction;
294a283635eSDmitry Karpeev   const char     *prefix;
295e884886fSBarry Smith 
296e884886fSBarry Smith   PetscFunctionBegin;
297251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
298e884886fSBarry Smith   if (iascii) {
299a283635eSDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"Matrix-free approximation:\n");CHKERRQ(ierr);
300a283635eSDmitry Karpeev     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
301e884886fSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"err=%G (relative error in function evaluation)\n",ctx->error_rel);CHKERRQ(ierr);
3027adad957SLisandro Dalcin     if (!((PetscObject)ctx)->type_name) {
303e884886fSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"The compute h routine has not yet been set\n");CHKERRQ(ierr);
304e884886fSBarry Smith     } else {
3057adad957SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"Using %s compute h routine\n",((PetscObject)ctx)->type_name);CHKERRQ(ierr);
306e884886fSBarry Smith     }
307e884886fSBarry Smith     if (ctx->ops->view) {
308e884886fSBarry Smith       ierr = (*ctx->ops->view)(ctx,viewer);CHKERRQ(ierr);
309e884886fSBarry Smith     }
310a283635eSDmitry Karpeev     ierr = PetscObjectGetOptionsPrefix((PetscObject)J, &prefix);CHKERRQ(ierr);
311a283635eSDmitry Karpeev 
312a283635eSDmitry Karpeev     ierr = PetscOptionsHasName(prefix, "-mat_mffd_view_base", &viewbase);CHKERRQ(ierr);
313a283635eSDmitry Karpeev     if (viewbase) {
314a283635eSDmitry Karpeev       ierr = PetscViewerASCIIPrintf(viewer, "Base:\n");CHKERRQ(ierr);
315a283635eSDmitry Karpeev       ierr = VecView(ctx->current_u, viewer);CHKERRQ(ierr);
316a283635eSDmitry Karpeev     }
317a283635eSDmitry Karpeev     ierr = PetscOptionsHasName(prefix, "-mat_mffd_view_function", &viewfunction);CHKERRQ(ierr);
318a283635eSDmitry Karpeev     if (viewfunction) {
319a283635eSDmitry Karpeev       ierr = PetscViewerASCIIPrintf(viewer, "Function:\n");CHKERRQ(ierr);
320a283635eSDmitry Karpeev       ierr = VecView(ctx->current_f, viewer);CHKERRQ(ierr);
321a283635eSDmitry Karpeev     }
322a283635eSDmitry Karpeev     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
32311aeaf0aSBarry Smith   }
324e884886fSBarry Smith   PetscFunctionReturn(0);
325e884886fSBarry Smith }
326e884886fSBarry Smith 
327e884886fSBarry Smith #undef __FUNCT__
328e884886fSBarry Smith #define __FUNCT__ "MatAssemblyEnd_MFFD"
329e884886fSBarry Smith /*
330e884886fSBarry Smith    MatAssemblyEnd_MFFD - Resets the ctx->ncurrenth to zero. This
331e884886fSBarry Smith    allows the user to indicate the beginning of a new linear solve by calling
332e884886fSBarry Smith    MatAssemblyXXX() on the matrix free matrix. This then allows the
3331d0fab5eSBarry Smith    MatCreateMFFD_WP() to properly compute ||U|| only the first time
334e884886fSBarry Smith    in the linear solver rather than every time.
3355a576424SJed Brown 
3365a576424SJed Brown    This function is referenced directly from MatAssemblyEnd_SNESMF(), which may be in a different shared library.
337e884886fSBarry Smith */
3385a576424SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyEnd_MFFD(Mat J,MatAssemblyType mt)
339e884886fSBarry Smith {
340e884886fSBarry Smith   PetscErrorCode ierr;
341e884886fSBarry Smith   MatMFFD        j = (MatMFFD)J->data;
342e884886fSBarry Smith 
343e884886fSBarry Smith   PetscFunctionBegin;
344e884886fSBarry Smith   ierr      = MatMFFDResetHHistory(J);CHKERRQ(ierr);
345e884886fSBarry Smith   j->vshift = 0.0;
346e884886fSBarry Smith   j->vscale = 1.0;
347e884886fSBarry Smith   PetscFunctionReturn(0);
348e884886fSBarry Smith }
349e884886fSBarry Smith 
350e884886fSBarry Smith #undef __FUNCT__
351e884886fSBarry Smith #define __FUNCT__ "MatMult_MFFD"
352e884886fSBarry Smith /*
353e884886fSBarry Smith   MatMult_MFFD - Default matrix-free form for Jacobian-vector product, y = F'(u)*a:
354e884886fSBarry Smith 
355e884886fSBarry Smith         y ~= (F(u + ha) - F(u))/h,
356e884886fSBarry Smith   where F = nonlinear function, as set by SNESSetFunction()
357e884886fSBarry Smith         u = current iterate
358e884886fSBarry Smith         h = difference interval
359e884886fSBarry Smith */
360e884886fSBarry Smith PetscErrorCode MatMult_MFFD(Mat mat,Vec a,Vec y)
361e884886fSBarry Smith {
362e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)mat->data;
363e884886fSBarry Smith   PetscScalar    h;
364e884886fSBarry Smith   Vec            w,U,F;
365e884886fSBarry Smith   PetscErrorCode ierr;
366ace3abfcSBarry Smith   PetscBool      zeroa;
367e884886fSBarry Smith 
368e884886fSBarry Smith   PetscFunctionBegin;
369ce94432eSBarry 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");
370e884886fSBarry Smith   /* We log matrix-free matrix-vector products separately, so that we can
371e884886fSBarry Smith      separate the performance monitoring from the cases that use conventional
372e884886fSBarry Smith      storage.  We may eventually modify event logging to associate events
373e884886fSBarry Smith      with particular objects, hence alleviating the more general problem. */
374e884886fSBarry Smith   ierr = PetscLogEventBegin(MATMFFD_Mult,a,y,0,0);CHKERRQ(ierr);
375e884886fSBarry Smith 
376e884886fSBarry Smith   w = ctx->w;
377e884886fSBarry Smith   U = ctx->current_u;
3783ec795f1SBarry Smith   F = ctx->current_f;
379e884886fSBarry Smith   /*
380e884886fSBarry Smith       Compute differencing parameter
381e884886fSBarry Smith   */
382e884886fSBarry Smith   if (!ctx->ops->compute) {
383e884886fSBarry Smith     ierr = MatMFFDSetType(mat,MATMFFD_WP);CHKERRQ(ierr);
3849c6ac3b3SBarry Smith     ierr = MatSetFromOptions(mat);CHKERRQ(ierr);
385e884886fSBarry Smith   }
386e884886fSBarry Smith   ierr = (*ctx->ops->compute)(ctx,U,a,&h,&zeroa);CHKERRQ(ierr);
387e884886fSBarry Smith   if (zeroa) {
388e884886fSBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
389e884886fSBarry Smith     PetscFunctionReturn(0);
390e884886fSBarry Smith   }
391e884886fSBarry Smith 
392e32f2f54SBarry Smith   if (PetscIsInfOrNanScalar(h)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Computed Nan differencing parameter h");
393e884886fSBarry Smith   if (ctx->checkh) {
394e884886fSBarry Smith     ierr = (*ctx->checkh)(ctx->checkhctx,U,a,&h);CHKERRQ(ierr);
395e884886fSBarry Smith   }
396e884886fSBarry Smith 
397e884886fSBarry Smith   /* keep a record of the current differencing parameter h */
398e884886fSBarry Smith   ctx->currenth = h;
399e884886fSBarry Smith #if defined(PETSC_USE_COMPLEX)
400e884886fSBarry Smith   ierr = PetscInfo2(mat,"Current differencing parameter: %G + %G i\n",PetscRealPart(h),PetscImaginaryPart(h));CHKERRQ(ierr);
401e884886fSBarry Smith #else
402e884886fSBarry Smith   ierr = PetscInfo1(mat,"Current differencing parameter: %15.12e\n",h);CHKERRQ(ierr);
403e884886fSBarry Smith #endif
404e884886fSBarry Smith   if (ctx->historyh && ctx->ncurrenth < ctx->maxcurrenth) {
405e884886fSBarry Smith     ctx->historyh[ctx->ncurrenth] = h;
406e884886fSBarry Smith   }
407e884886fSBarry Smith   ctx->ncurrenth++;
408e884886fSBarry Smith 
409e884886fSBarry Smith   /* w = u + ha */
410c3bb7e23SBarry Smith   if (ctx->drscale) {
411c3bb7e23SBarry Smith     ierr = VecPointwiseMult(ctx->drscale,a,U);CHKERRQ(ierr);
412c3bb7e23SBarry Smith     ierr = VecAYPX(U,h,w);CHKERRQ(ierr);
413c3bb7e23SBarry Smith   } else {
414e884886fSBarry Smith     ierr = VecWAXPY(w,h,a,U);CHKERRQ(ierr);
415c3bb7e23SBarry Smith   }
416e884886fSBarry Smith 
4171b797266SDmitry Karpeev   /* compute func(U) as base for differencing; only needed first time in and not when provided by user */
4181b797266SDmitry Karpeev   if (ctx->ncurrenth == 1 && ctx->current_f_allocated) {
419e884886fSBarry Smith     ierr = (*ctx->func)(ctx->funcctx,U,F);CHKERRQ(ierr);
42052121784SBarry Smith   }
421e884886fSBarry Smith   ierr = (*ctx->func)(ctx->funcctx,w,y);CHKERRQ(ierr);
422e884886fSBarry Smith 
423e884886fSBarry Smith   ierr = VecAXPY(y,-1.0,F);CHKERRQ(ierr);
424e884886fSBarry Smith   ierr = VecScale(y,1.0/h);CHKERRQ(ierr);
425e884886fSBarry Smith 
426c35ec66cSMatthew G Knepley   if ((ctx->vshift != 0.0) || (ctx->vscale != 1.0)) {
427e884886fSBarry Smith     ierr = VecAXPBY(y,ctx->vshift,ctx->vscale,a);CHKERRQ(ierr);
428c3bb7e23SBarry Smith   }
429c3bb7e23SBarry Smith   if (ctx->dlscale) {
430c3bb7e23SBarry Smith     ierr = VecPointwiseMult(y,ctx->dlscale,y);CHKERRQ(ierr);
431c3bb7e23SBarry Smith   }
432c3bb7e23SBarry Smith   if (ctx->dshift) {
433c3bb7e23SBarry Smith     ierr = VecPointwiseMult(ctx->dshift,a,U);CHKERRQ(ierr);
434c3bb7e23SBarry Smith     ierr = VecAXPY(y,1.0,U);CHKERRQ(ierr);
435c3bb7e23SBarry Smith   }
436e884886fSBarry Smith 
4370298fd71SBarry Smith   if (ctx->sp) {ierr = MatNullSpaceRemove(ctx->sp,y,NULL);CHKERRQ(ierr);}
438e884886fSBarry Smith 
439e884886fSBarry Smith   ierr = PetscLogEventEnd(MATMFFD_Mult,a,y,0,0);CHKERRQ(ierr);
440e884886fSBarry Smith   PetscFunctionReturn(0);
441e884886fSBarry Smith }
442e884886fSBarry Smith 
443e884886fSBarry Smith #undef __FUNCT__
444e884886fSBarry Smith #define __FUNCT__ "MatGetDiagonal_MFFD"
445e884886fSBarry Smith /*
446e884886fSBarry Smith   MatGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix
447e884886fSBarry Smith 
448e884886fSBarry Smith         y ~= (F(u + ha) - F(u))/h,
449e884886fSBarry Smith   where F = nonlinear function, as set by SNESSetFunction()
450e884886fSBarry Smith         u = current iterate
451e884886fSBarry Smith         h = difference interval
452e884886fSBarry Smith */
453e884886fSBarry Smith PetscErrorCode MatGetDiagonal_MFFD(Mat mat,Vec a)
454e884886fSBarry Smith {
455e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)mat->data;
456e884886fSBarry Smith   PetscScalar    h,*aa,*ww,v;
457e884886fSBarry Smith   PetscReal      epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON;
458e884886fSBarry Smith   Vec            w,U;
459e884886fSBarry Smith   PetscErrorCode ierr;
460e884886fSBarry Smith   PetscInt       i,rstart,rend;
461e884886fSBarry Smith 
462e884886fSBarry Smith   PetscFunctionBegin;
463e7e72b3dSBarry Smith   if (!ctx->funci) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Requires calling MatMFFDSetFunctioni() first");
464e884886fSBarry Smith 
465e884886fSBarry Smith   w    = ctx->w;
466e884886fSBarry Smith   U    = ctx->current_u;
467e884886fSBarry Smith   ierr = (*ctx->func)(ctx->funcctx,U,a);CHKERRQ(ierr);
468e884886fSBarry Smith   ierr = (*ctx->funcisetbase)(ctx->funcctx,U);CHKERRQ(ierr);
469e884886fSBarry Smith   ierr = VecCopy(U,w);CHKERRQ(ierr);
470e884886fSBarry Smith 
471e884886fSBarry Smith   ierr = VecGetOwnershipRange(a,&rstart,&rend);CHKERRQ(ierr);
472e884886fSBarry Smith   ierr = VecGetArray(a,&aa);CHKERRQ(ierr);
473e884886fSBarry Smith   for (i=rstart; i<rend; i++) {
474e884886fSBarry Smith     ierr = VecGetArray(w,&ww);CHKERRQ(ierr);
475e884886fSBarry Smith     h    = ww[i-rstart];
476e884886fSBarry Smith     if (h == 0.0) h = 1.0;
477e884886fSBarry Smith     if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0)     h = umin;
478e884886fSBarry Smith     else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin;
479e884886fSBarry Smith     h *= epsilon;
480e884886fSBarry Smith 
481e884886fSBarry Smith     ww[i-rstart] += h;
482e884886fSBarry Smith     ierr          = VecRestoreArray(w,&ww);CHKERRQ(ierr);
483e884886fSBarry Smith     ierr          = (*ctx->funci)(ctx->funcctx,i,w,&v);CHKERRQ(ierr);
484e884886fSBarry Smith     aa[i-rstart]  = (v - aa[i-rstart])/h;
485e884886fSBarry Smith 
486e884886fSBarry Smith     /* possibly shift and scale result */
487c35ec66cSMatthew G Knepley     if ((ctx->vshift != 0.0) || (ctx->vscale != 1.0)) {
488e884886fSBarry Smith       aa[i - rstart] = ctx->vshift + ctx->vscale*aa[i-rstart];
489c3bb7e23SBarry Smith     }
490e884886fSBarry Smith 
491e884886fSBarry Smith     ierr          = VecGetArray(w,&ww);CHKERRQ(ierr);
492e884886fSBarry Smith     ww[i-rstart] -= h;
493e884886fSBarry Smith     ierr          = VecRestoreArray(w,&ww);CHKERRQ(ierr);
494e884886fSBarry Smith   }
495e884886fSBarry Smith   ierr = VecRestoreArray(a,&aa);CHKERRQ(ierr);
496e884886fSBarry Smith   PetscFunctionReturn(0);
497e884886fSBarry Smith }
498e884886fSBarry Smith 
499e884886fSBarry Smith #undef __FUNCT__
500c3bb7e23SBarry Smith #define __FUNCT__ "MatDiagonalScale_MFFD"
501c3bb7e23SBarry Smith PetscErrorCode MatDiagonalScale_MFFD(Mat mat,Vec ll,Vec rr)
502c3bb7e23SBarry Smith {
503c3bb7e23SBarry Smith   MatMFFD        aij = (MatMFFD)mat->data;
504c3bb7e23SBarry Smith   PetscErrorCode ierr;
505c3bb7e23SBarry Smith 
506c3bb7e23SBarry Smith   PetscFunctionBegin;
507c3bb7e23SBarry Smith   if (ll && !aij->dlscale) {
508c3bb7e23SBarry Smith     ierr = VecDuplicate(ll,&aij->dlscale);CHKERRQ(ierr);
509c3bb7e23SBarry Smith   }
510c3bb7e23SBarry Smith   if (rr && !aij->drscale) {
511c3bb7e23SBarry Smith     ierr = VecDuplicate(rr,&aij->drscale);CHKERRQ(ierr);
512c3bb7e23SBarry Smith   }
513c3bb7e23SBarry Smith   if (ll) {
514c3bb7e23SBarry Smith     ierr = VecCopy(ll,aij->dlscale);CHKERRQ(ierr);
515c3bb7e23SBarry Smith   }
516c3bb7e23SBarry Smith   if (rr) {
517c3bb7e23SBarry Smith     ierr = VecCopy(rr,aij->drscale);CHKERRQ(ierr);
518c3bb7e23SBarry Smith   }
519c3bb7e23SBarry Smith   PetscFunctionReturn(0);
520c3bb7e23SBarry Smith }
521c3bb7e23SBarry Smith 
522c3bb7e23SBarry Smith #undef __FUNCT__
523c3bb7e23SBarry Smith #define __FUNCT__ "MatDiagonalSet_MFFD"
524c3bb7e23SBarry Smith PetscErrorCode MatDiagonalSet_MFFD(Mat mat,Vec ll,InsertMode mode)
525c3bb7e23SBarry Smith {
526c3bb7e23SBarry Smith   MatMFFD        aij = (MatMFFD)mat->data;
527c3bb7e23SBarry Smith   PetscErrorCode ierr;
528c3bb7e23SBarry Smith 
529c3bb7e23SBarry Smith   PetscFunctionBegin;
530ce94432eSBarry Smith   if (mode == INSERT_VALUES) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"No diagonal set with INSERT_VALUES");
531c3bb7e23SBarry Smith   if (!aij->dshift) {
532c3bb7e23SBarry Smith     ierr = VecDuplicate(ll,&aij->dshift);CHKERRQ(ierr);
533c3bb7e23SBarry Smith   }
534c3bb7e23SBarry Smith   ierr = VecAXPY(aij->dshift,1.0,ll);CHKERRQ(ierr);
535c3bb7e23SBarry Smith   PetscFunctionReturn(0);
536c3bb7e23SBarry Smith }
537c3bb7e23SBarry Smith 
538c3bb7e23SBarry Smith #undef __FUNCT__
539e884886fSBarry Smith #define __FUNCT__ "MatShift_MFFD"
540e884886fSBarry Smith PetscErrorCode MatShift_MFFD(Mat Y,PetscScalar a)
541e884886fSBarry Smith {
542e884886fSBarry Smith   MatMFFD shell = (MatMFFD)Y->data;
5435fd66863SKarl Rupp 
544e884886fSBarry Smith   PetscFunctionBegin;
545e884886fSBarry Smith   shell->vshift += a;
546e884886fSBarry Smith   PetscFunctionReturn(0);
547e884886fSBarry Smith }
548e884886fSBarry Smith 
549e884886fSBarry Smith #undef __FUNCT__
550e884886fSBarry Smith #define __FUNCT__ "MatScale_MFFD"
551e884886fSBarry Smith PetscErrorCode MatScale_MFFD(Mat Y,PetscScalar a)
552e884886fSBarry Smith {
553e884886fSBarry Smith   MatMFFD shell = (MatMFFD)Y->data;
5545fd66863SKarl Rupp 
555e884886fSBarry Smith   PetscFunctionBegin;
556e884886fSBarry Smith   shell->vscale *= a;
557e884886fSBarry Smith   PetscFunctionReturn(0);
558e884886fSBarry Smith }
559e884886fSBarry Smith 
560e884886fSBarry Smith #undef __FUNCT__
561c879e56bSBarry Smith #define __FUNCT__ "MatMFFDSetBase_MFFD"
562d8fa6a54SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDSetBase_MFFD(Mat J,Vec U,Vec F)
563e884886fSBarry Smith {
564e884886fSBarry Smith   PetscErrorCode ierr;
565e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)J->data;
566e884886fSBarry Smith 
567e884886fSBarry Smith   PetscFunctionBegin;
568e884886fSBarry Smith   ierr = MatMFFDResetHHistory(J);CHKERRQ(ierr);
5692205254eSKarl Rupp 
570e884886fSBarry Smith   ctx->current_u = U;
57152121784SBarry Smith   if (F) {
5726bf464f9SBarry Smith     if (ctx->current_f_allocated) {ierr = VecDestroy(&ctx->current_f);CHKERRQ(ierr);}
5733ec795f1SBarry Smith     ctx->current_f           = F;
574cfe22f5eSBarry Smith     ctx->current_f_allocated = PETSC_FALSE;
575cfe22f5eSBarry Smith   } else if (!ctx->current_f_allocated) {
57652121784SBarry Smith     ierr = VecDuplicate(ctx->current_u, &ctx->current_f);CHKERRQ(ierr);
5772205254eSKarl Rupp 
578cfe22f5eSBarry Smith     ctx->current_f_allocated = PETSC_TRUE;
57952121784SBarry Smith   }
580e884886fSBarry Smith   if (!ctx->w) {
581e884886fSBarry Smith     ierr = VecDuplicate(ctx->current_u, &ctx->w);CHKERRQ(ierr);
582e884886fSBarry Smith   }
583e884886fSBarry Smith   J->assembled = PETSC_TRUE;
584e884886fSBarry Smith   PetscFunctionReturn(0);
585e884886fSBarry Smith }
5865a576424SJed Brown 
587e884886fSBarry Smith typedef PetscErrorCode (*FCN3)(void*,Vec,Vec,PetscScalar*); /* force argument to next function to not be extern C*/
588b2573a8aSBarry Smith 
589e884886fSBarry Smith #undef __FUNCT__
590c879e56bSBarry Smith #define __FUNCT__ "MatMFFDSetCheckh_MFFD"
5917087cfbeSBarry Smith PetscErrorCode  MatMFFDSetCheckh_MFFD(Mat J,FCN3 fun,void *ectx)
592e884886fSBarry Smith {
593e884886fSBarry Smith   MatMFFD ctx = (MatMFFD)J->data;
594e884886fSBarry Smith 
595e884886fSBarry Smith   PetscFunctionBegin;
596e884886fSBarry Smith   ctx->checkh    = fun;
597e884886fSBarry Smith   ctx->checkhctx = ectx;
598e884886fSBarry Smith   PetscFunctionReturn(0);
599e884886fSBarry Smith }
600e884886fSBarry Smith 
601e884886fSBarry Smith #undef __FUNCT__
6026aa9148fSLisandro Dalcin #define __FUNCT__ "MatMFFDSetOptionsPrefix"
6036aa9148fSLisandro Dalcin /*@C
6046aa9148fSLisandro Dalcin    MatMFFDSetOptionsPrefix - Sets the prefix used for searching for all
6056aa9148fSLisandro Dalcin    MatMFFD options in the database.
6066aa9148fSLisandro Dalcin 
6076aa9148fSLisandro Dalcin    Collective on Mat
6086aa9148fSLisandro Dalcin 
6096aa9148fSLisandro Dalcin    Input Parameter:
6106aa9148fSLisandro Dalcin +  A - the Mat context
6116aa9148fSLisandro Dalcin -  prefix - the prefix to prepend to all option names
6126aa9148fSLisandro Dalcin 
6136aa9148fSLisandro Dalcin    Notes:
6146aa9148fSLisandro Dalcin    A hyphen (-) must NOT be given at the beginning of the prefix name.
6156aa9148fSLisandro Dalcin    The first character of all runtime options is AUTOMATICALLY the hyphen.
6166aa9148fSLisandro Dalcin 
6176aa9148fSLisandro Dalcin    Level: advanced
6186aa9148fSLisandro Dalcin 
6196aa9148fSLisandro Dalcin .keywords: SNES, matrix-free, parameters
6206aa9148fSLisandro Dalcin 
6219c6ac3b3SBarry Smith .seealso: MatSetFromOptions(), MatCreateSNESMF()
6226aa9148fSLisandro Dalcin @*/
6237087cfbeSBarry Smith PetscErrorCode  MatMFFDSetOptionsPrefix(Mat mat,const char prefix[])
6246aa9148fSLisandro Dalcin 
6256aa9148fSLisandro Dalcin {
6260298fd71SBarry Smith   MatMFFD        mfctx = mat ? (MatMFFD)mat->data : (MatMFFD)NULL;
6276aa9148fSLisandro Dalcin   PetscErrorCode ierr;
6285fd66863SKarl Rupp 
6296aa9148fSLisandro Dalcin   PetscFunctionBegin;
6300700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
6310700a824SBarry Smith   PetscValidHeaderSpecific(mfctx,MATMFFD_CLASSID,1);
6326aa9148fSLisandro Dalcin   ierr = PetscObjectSetOptionsPrefix((PetscObject)mfctx,prefix);CHKERRQ(ierr);
6336aa9148fSLisandro Dalcin   PetscFunctionReturn(0);
6346aa9148fSLisandro Dalcin }
6356aa9148fSLisandro Dalcin 
6366aa9148fSLisandro Dalcin #undef __FUNCT__
6379c6ac3b3SBarry Smith #define __FUNCT__ "MatSetFromOptions_MFFD"
6389c6ac3b3SBarry Smith PetscErrorCode  MatSetFromOptions_MFFD(Mat mat)
639e884886fSBarry Smith {
64071f08665SBarry Smith   MatMFFD        mfctx = (MatMFFD)mat->data;
641e884886fSBarry Smith   PetscErrorCode ierr;
642ace3abfcSBarry Smith   PetscBool      flg;
643e884886fSBarry Smith   char           ftype[256];
644e884886fSBarry Smith 
645e884886fSBarry Smith   PetscFunctionBegin;
6460700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
6470700a824SBarry Smith   PetscValidHeaderSpecific(mfctx,MATMFFD_CLASSID,1);
6483194b578SJed Brown   ierr = PetscObjectOptionsBegin((PetscObject)mfctx);CHKERRQ(ierr);
649b022a5c1SBarry Smith   ierr = PetscOptionsList("-mat_mffd_type","Matrix free type","MatMFFDSetType",MatMFFDList,((PetscObject)mfctx)->type_name,ftype,256,&flg);CHKERRQ(ierr);
650e884886fSBarry Smith   if (flg) {
651e884886fSBarry Smith     ierr = MatMFFDSetType(mat,ftype);CHKERRQ(ierr);
652e884886fSBarry Smith   }
653e884886fSBarry Smith 
654e884886fSBarry Smith   ierr = PetscOptionsReal("-mat_mffd_err","set sqrt relative error in function","MatMFFDSetFunctionError",mfctx->error_rel,&mfctx->error_rel,0);CHKERRQ(ierr);
655e884886fSBarry Smith   ierr = PetscOptionsInt("-mat_mffd_period","how often h is recomputed","MatMFFDSetPeriod",mfctx->recomputeperiod,&mfctx->recomputeperiod,0);CHKERRQ(ierr);
656e884886fSBarry Smith 
65790d69ab7SBarry Smith   flg  = PETSC_FALSE;
6580298fd71SBarry Smith   ierr = PetscOptionsBool("-mat_mffd_check_positivity","Insure that U + h*a is nonnegative","MatMFFDSetCheckh",flg,&flg,NULL);CHKERRQ(ierr);
659e884886fSBarry Smith   if (flg) {
660e884886fSBarry Smith     ierr = MatMFFDSetCheckh(mat,MatMFFDCheckPositivity,0);CHKERRQ(ierr);
661e884886fSBarry Smith   }
662e884886fSBarry Smith   if (mfctx->ops->setfromoptions) {
663e884886fSBarry Smith     ierr = (*mfctx->ops->setfromoptions)(mfctx);CHKERRQ(ierr);
664e884886fSBarry Smith   }
665e884886fSBarry Smith   ierr = PetscOptionsEnd();CHKERRQ(ierr);
666e884886fSBarry Smith   PetscFunctionReturn(0);
667e884886fSBarry Smith }
668e884886fSBarry Smith 
669bc0f08ceSBarry Smith #undef __FUNCT__
670bc0f08ceSBarry Smith #define __FUNCT__ "MatMFFDSetPeriod_MFFD"
671bc0f08ceSBarry Smith PetscErrorCode  MatMFFDSetPeriod_MFFD(Mat mat,PetscInt period)
672bc0f08ceSBarry Smith {
673bc0f08ceSBarry Smith   MatMFFD ctx = (MatMFFD)mat->data;
674bc0f08ceSBarry Smith 
675bc0f08ceSBarry Smith   PetscFunctionBegin;
676bc0f08ceSBarry Smith   PetscValidLogicalCollectiveInt(mat,period,2);
677bc0f08ceSBarry Smith   ctx->recomputeperiod = period;
678bc0f08ceSBarry Smith   PetscFunctionReturn(0);
679bc0f08ceSBarry Smith }
680bc0f08ceSBarry Smith 
681bc0f08ceSBarry Smith #undef __FUNCT__
682bc0f08ceSBarry Smith #define __FUNCT__ "MatMFFDSetFunction_MFFD"
683bc0f08ceSBarry Smith PetscErrorCode  MatMFFDSetFunction_MFFD(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx)
684bc0f08ceSBarry Smith {
685bc0f08ceSBarry Smith   MatMFFD ctx = (MatMFFD)mat->data;
686bc0f08ceSBarry Smith 
687bc0f08ceSBarry Smith   PetscFunctionBegin;
688bc0f08ceSBarry Smith   ctx->func    = func;
689bc0f08ceSBarry Smith   ctx->funcctx = funcctx;
690bc0f08ceSBarry Smith   PetscFunctionReturn(0);
691bc0f08ceSBarry Smith }
692bc0f08ceSBarry Smith 
693bc0f08ceSBarry Smith #undef __FUNCT__
694bc0f08ceSBarry Smith #define __FUNCT__ "MatMFFDSetFunctionError_MFFD"
695bc0f08ceSBarry Smith PetscErrorCode  MatMFFDSetFunctionError_MFFD(Mat mat,PetscReal error)
696bc0f08ceSBarry Smith {
697bc0f08ceSBarry Smith   MatMFFD ctx = (MatMFFD)mat->data;
698bc0f08ceSBarry Smith 
699bc0f08ceSBarry Smith   PetscFunctionBegin;
700bc0f08ceSBarry Smith   PetscValidLogicalCollectiveReal(mat,error,2);
701bc0f08ceSBarry Smith   if (error != PETSC_DEFAULT) ctx->error_rel = error;
702bc0f08ceSBarry Smith   PetscFunctionReturn(0);
703bc0f08ceSBarry Smith }
704bc0f08ceSBarry Smith 
705e884886fSBarry Smith /*MC
706e884886fSBarry Smith   MATMFFD - MATMFFD = "mffd" - A matrix free matrix type.
707e884886fSBarry Smith 
708e884886fSBarry Smith   Level: advanced
709e884886fSBarry Smith 
710d2d6cebeSBarry Smith .seealso: MatCreateMFFD(), MatCreateSNESMF(), MatMFFDSetFunction()
711e884886fSBarry Smith M*/
712e884886fSBarry Smith #undef __FUNCT__
713e884886fSBarry Smith #define __FUNCT__ "MatCreate_MFFD"
7148cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_MFFD(Mat A)
715e884886fSBarry Smith {
716e884886fSBarry Smith   MatMFFD        mfctx;
717e884886fSBarry Smith   PetscErrorCode ierr;
718e884886fSBarry Smith 
719e884886fSBarry Smith   PetscFunctionBegin;
720519f805aSKarl Rupp #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
721607a6623SBarry Smith   ierr = MatMFFDInitializePackage();CHKERRQ(ierr);
7223ec795f1SBarry Smith #endif
723e884886fSBarry Smith 
724ce94432eSBarry 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);
7252205254eSKarl Rupp 
726e884886fSBarry Smith   mfctx->sp                       = 0;
727e884886fSBarry Smith   mfctx->error_rel                = PETSC_SQRT_MACHINE_EPSILON;
728e884886fSBarry Smith   mfctx->recomputeperiod          = 1;
729e884886fSBarry Smith   mfctx->count                    = 0;
730e884886fSBarry Smith   mfctx->currenth                 = 0.0;
7310298fd71SBarry Smith   mfctx->historyh                 = NULL;
732e884886fSBarry Smith   mfctx->ncurrenth                = 0;
733e884886fSBarry Smith   mfctx->maxcurrenth              = 0;
7347adad957SLisandro Dalcin   ((PetscObject)mfctx)->type_name = 0;
735e884886fSBarry Smith 
736e884886fSBarry Smith   mfctx->vshift = 0.0;
737e884886fSBarry Smith   mfctx->vscale = 1.0;
738e884886fSBarry Smith 
739e884886fSBarry Smith   /*
740e884886fSBarry Smith      Create the empty data structure to contain compute-h routines.
741e884886fSBarry Smith      These will be filled in below from the command line options or
742e884886fSBarry Smith      a later call with MatMFFDSetType() or if that is not called
743e884886fSBarry Smith      then it will default in the first use of MatMult_MFFD()
744e884886fSBarry Smith   */
745e884886fSBarry Smith   mfctx->ops->compute        = 0;
746e884886fSBarry Smith   mfctx->ops->destroy        = 0;
747e884886fSBarry Smith   mfctx->ops->view           = 0;
748e884886fSBarry Smith   mfctx->ops->setfromoptions = 0;
749e884886fSBarry Smith   mfctx->hctx                = 0;
750e884886fSBarry Smith 
751e884886fSBarry Smith   mfctx->func    = 0;
752e884886fSBarry Smith   mfctx->funcctx = 0;
7530298fd71SBarry Smith   mfctx->w       = NULL;
754e884886fSBarry Smith 
755e884886fSBarry Smith   A->data = mfctx;
756e884886fSBarry Smith 
757e884886fSBarry Smith   A->ops->mult           = MatMult_MFFD;
758e884886fSBarry Smith   A->ops->destroy        = MatDestroy_MFFD;
759e884886fSBarry Smith   A->ops->view           = MatView_MFFD;
760e884886fSBarry Smith   A->ops->assemblyend    = MatAssemblyEnd_MFFD;
761e884886fSBarry Smith   A->ops->getdiagonal    = MatGetDiagonal_MFFD;
762e884886fSBarry Smith   A->ops->scale          = MatScale_MFFD;
763e884886fSBarry Smith   A->ops->shift          = MatShift_MFFD;
764c3bb7e23SBarry Smith   A->ops->diagonalscale  = MatDiagonalScale_MFFD;
765c3bb7e23SBarry Smith   A->ops->diagonalset    = MatDiagonalSet_MFFD;
7669c6ac3b3SBarry Smith   A->ops->setfromoptions = MatSetFromOptions_MFFD;
767e884886fSBarry Smith   A->assembled           = PETSC_TRUE;
768e884886fSBarry Smith 
76926283091SBarry Smith   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
77026283091SBarry Smith   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
771ee1cef2cSJed Brown 
772bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetBase_C",MatMFFDSetBase_MFFD);CHKERRQ(ierr);
773bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctioniBase_C",MatMFFDSetFunctioniBase_MFFD);CHKERRQ(ierr);
774bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctioni_C",MatMFFDSetFunctioni_MFFD);CHKERRQ(ierr);
775bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunction_C",MatMFFDSetFunction_MFFD);CHKERRQ(ierr);
776bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetCheckh_C",MatMFFDSetCheckh_MFFD);CHKERRQ(ierr);
777bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetPeriod_C",MatMFFDSetPeriod_MFFD);CHKERRQ(ierr);
778bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctionError_C",MatMFFDSetFunctionError_MFFD);CHKERRQ(ierr);
779bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDResetHHistory_C",MatMFFDResetHHistory_MFFD);CHKERRQ(ierr);
780bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDAddNullSpace_C",MatMFFDAddNullSpace_MFFD);CHKERRQ(ierr);
7812205254eSKarl Rupp 
782e884886fSBarry Smith   mfctx->mat = A;
7832205254eSKarl Rupp 
784e884886fSBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,MATMFFD);CHKERRQ(ierr);
785e884886fSBarry Smith   PetscFunctionReturn(0);
786e884886fSBarry Smith }
787e884886fSBarry Smith 
788e884886fSBarry Smith #undef __FUNCT__
789e884886fSBarry Smith #define __FUNCT__ "MatCreateMFFD"
790e884886fSBarry Smith /*@
791e884886fSBarry Smith    MatCreateMFFD - Creates a matrix-free matrix. See also MatCreateSNESMF()
792e884886fSBarry Smith 
793e884886fSBarry Smith    Collective on Vec
794e884886fSBarry Smith 
795e884886fSBarry Smith    Input Parameters:
796fef1beadSBarry Smith +  comm - MPI communicator
797fef1beadSBarry Smith .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
798fef1beadSBarry Smith            This value should be the same as the local size used in creating the
799fef1beadSBarry Smith            y vector for the matrix-vector product y = Ax.
800fef1beadSBarry Smith .  n - This value should be the same as the local size used in creating the
801fef1beadSBarry Smith        x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have
802fef1beadSBarry Smith        calculated if N is given) For square matrices n is almost always m.
803fef1beadSBarry Smith .  M - number of global rows (or PETSC_DETERMINE to have calculated if m is given)
804fef1beadSBarry Smith -  N - number of global columns (or PETSC_DETERMINE to have calculated if n is given)
805fef1beadSBarry Smith 
806e884886fSBarry Smith 
807e884886fSBarry Smith    Output Parameter:
808e884886fSBarry Smith .  J - the matrix-free matrix
809e884886fSBarry Smith 
8109c6ac3b3SBarry Smith    Options Database Keys: call MatSetFromOptions() to trigger these
8119c6ac3b3SBarry Smith +  -mat_mffd_type - wp or ds (see MATMFFD_WP or MATMFFD_DS)
8129c6ac3b3SBarry Smith -  -mat_mffd_err - square root of estimated relative error in function evaluation
8139c6ac3b3SBarry Smith -  -mat_mffd_period - how often h is recomputed, defaults to 1, everytime
8149c6ac3b3SBarry Smith 
8159c6ac3b3SBarry Smith 
816e884886fSBarry Smith    Level: advanced
817e884886fSBarry Smith 
818e884886fSBarry Smith    Notes:
819e884886fSBarry Smith    The matrix-free matrix context merely contains the function pointers
820e884886fSBarry Smith    and work space for performing finite difference approximations of
821e884886fSBarry Smith    Jacobian-vector products, F'(u)*a,
822e884886fSBarry Smith 
823e884886fSBarry Smith    The default code uses the following approach to compute h
824e884886fSBarry Smith 
825e884886fSBarry Smith .vb
826e884886fSBarry Smith      F'(u)*a = [F(u+h*a) - F(u)]/h where
827e884886fSBarry Smith      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
828e884886fSBarry Smith        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
829e884886fSBarry Smith  where
830e884886fSBarry Smith      error_rel = square root of relative error in function evaluation
831e884886fSBarry Smith      umin = minimum iterate parameter
832e884886fSBarry Smith .ve
833e884886fSBarry Smith 
834e884886fSBarry Smith    The user can set the error_rel via MatMFFDSetFunctionError() and
8350598bfebSBarry Smith    umin via MatMFFDDSSetUmin(); see the <A href="../../docs/manual.pdf#nameddest=ch_snes">SNES chapter of the users manual</A> for details.
836e884886fSBarry Smith 
837e884886fSBarry Smith    The user should call MatDestroy() when finished with the matrix-free
838e884886fSBarry Smith    matrix context.
839e884886fSBarry Smith 
840e884886fSBarry Smith    Options Database Keys:
841e884886fSBarry Smith +  -mat_mffd_err <error_rel> - Sets error_rel
842e884886fSBarry Smith .  -mat_mffd_unim <umin> - Sets umin (for default PETSc routine that computes h only)
843e884886fSBarry Smith -  -mat_mffd_check_positivity
844e884886fSBarry Smith 
845e884886fSBarry Smith .keywords: default, matrix-free, create, matrix
846e884886fSBarry Smith 
8471d0fab5eSBarry Smith .seealso: MatDestroy(), MatMFFDSetFunctionError(), MatMFFDDSSetUmin(), MatMFFDSetFunction()
848e884886fSBarry Smith           MatMFFDSetHHistory(), MatMFFDResetHHistory(), MatCreateSNESMF(),
8491d0fab5eSBarry Smith           MatMFFDGetH(), MatMFFDRegisterDynamic), MatMFFDComputeJacobian()
850e884886fSBarry Smith 
851e884886fSBarry Smith @*/
8527087cfbeSBarry Smith PetscErrorCode  MatCreateMFFD(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *J)
853e884886fSBarry Smith {
854e884886fSBarry Smith   PetscErrorCode ierr;
855e884886fSBarry Smith 
856e884886fSBarry Smith   PetscFunctionBegin;
857e884886fSBarry Smith   ierr = MatCreate(comm,J);CHKERRQ(ierr);
858fef1beadSBarry Smith   ierr = MatSetSizes(*J,m,n,M,N);CHKERRQ(ierr);
859e884886fSBarry Smith   ierr = MatSetType(*J,MATMFFD);CHKERRQ(ierr);
860be7c243fSJed Brown   ierr = MatSetUp(*J);CHKERRQ(ierr);
861e884886fSBarry Smith   PetscFunctionReturn(0);
862e884886fSBarry Smith }
863e884886fSBarry Smith 
864e884886fSBarry Smith 
865e884886fSBarry Smith #undef __FUNCT__
866e884886fSBarry Smith #define __FUNCT__ "MatMFFDGetH"
867e884886fSBarry Smith /*@
868e884886fSBarry Smith    MatMFFDGetH - Gets the last value that was used as the differencing
869e884886fSBarry Smith    parameter.
870e884886fSBarry Smith 
871e884886fSBarry Smith    Not Collective
872e884886fSBarry Smith 
873e884886fSBarry Smith    Input Parameters:
874e884886fSBarry Smith .  mat - the matrix obtained with MatCreateSNESMF()
875e884886fSBarry Smith 
876e884886fSBarry Smith    Output Paramter:
877e884886fSBarry Smith .  h - the differencing step size
878e884886fSBarry Smith 
879e884886fSBarry Smith    Level: advanced
880e884886fSBarry Smith 
881e884886fSBarry Smith .keywords: SNES, matrix-free, parameters
882e884886fSBarry Smith 
8831d0fab5eSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDSetHHistory(), MatCreateMFFD(), MATMFFD, MatMFFDResetHHistory()
884e884886fSBarry Smith @*/
8857087cfbeSBarry Smith PetscErrorCode  MatMFFDGetH(Mat mat,PetscScalar *h)
886e884886fSBarry Smith {
887e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)mat->data;
888bc0f08ceSBarry Smith   PetscErrorCode ierr;
889bc0f08ceSBarry Smith   PetscBool      match;
890e884886fSBarry Smith 
891e884886fSBarry Smith   PetscFunctionBegin;
892251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)mat,MATMFFD,&match);CHKERRQ(ierr);
893ce94432eSBarry Smith   if (!match) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONG,"Not a MFFD matrix");
894bc0f08ceSBarry Smith 
895e884886fSBarry Smith   *h = ctx->currenth;
896e884886fSBarry Smith   PetscFunctionReturn(0);
897e884886fSBarry Smith }
898e884886fSBarry Smith 
899e884886fSBarry Smith #undef __FUNCT__
900e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetFunction"
901e884886fSBarry Smith /*@C
902e884886fSBarry Smith    MatMFFDSetFunction - Sets the function used in applying the matrix free.
903e884886fSBarry Smith 
9043f9fe445SBarry Smith    Logically Collective on Mat
905e884886fSBarry Smith 
906e884886fSBarry Smith    Input Parameters:
907e884886fSBarry Smith +  mat - the matrix free matrix created via MatCreateSNESMF()
908e884886fSBarry Smith .  func - the function to use
909e884886fSBarry Smith -  funcctx - optional function context passed to function
910e884886fSBarry Smith 
911e884886fSBarry Smith    Level: advanced
912e884886fSBarry Smith 
913e884886fSBarry Smith    Notes:
914e884886fSBarry Smith     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
915e884886fSBarry Smith     matrix inside your compute Jacobian routine
916e884886fSBarry Smith 
9173ec795f1SBarry Smith     If this is not set then it will use the function set with SNESSetFunction() if MatCreateSNESMF() was used.
918e884886fSBarry Smith 
919e884886fSBarry Smith .keywords: SNES, matrix-free, function
920e884886fSBarry Smith 
9211d0fab5eSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD,
9221d0fab5eSBarry Smith           MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction()
923e884886fSBarry Smith @*/
9247087cfbeSBarry Smith PetscErrorCode  MatMFFDSetFunction(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx)
925e884886fSBarry Smith {
926bc0f08ceSBarry Smith   PetscErrorCode ierr;
927e884886fSBarry Smith 
928e884886fSBarry Smith   PetscFunctionBegin;
929bc0f08ceSBarry Smith   ierr = PetscTryMethod(mat,"MatMFFDSetFunction_C",(Mat,PetscErrorCode (*)(void*,Vec,Vec),void*),(mat,func,funcctx));CHKERRQ(ierr);
930e884886fSBarry Smith   PetscFunctionReturn(0);
931e884886fSBarry Smith }
932e884886fSBarry Smith 
933e884886fSBarry Smith #undef __FUNCT__
934e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetFunctioni"
935e884886fSBarry Smith /*@C
936e884886fSBarry Smith    MatMFFDSetFunctioni - Sets the function for a single component
937e884886fSBarry Smith 
9383f9fe445SBarry Smith    Logically Collective on Mat
939e884886fSBarry Smith 
940e884886fSBarry Smith    Input Parameters:
941e884886fSBarry Smith +  mat - the matrix free matrix created via MatCreateSNESMF()
942e884886fSBarry Smith -  funci - the function to use
943e884886fSBarry Smith 
944e884886fSBarry Smith    Level: advanced
945e884886fSBarry Smith 
946e884886fSBarry Smith    Notes:
947e884886fSBarry Smith     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
948e884886fSBarry Smith     matrix inside your compute Jacobian routine
949e884886fSBarry Smith 
950e884886fSBarry Smith 
951e884886fSBarry Smith .keywords: SNES, matrix-free, function
952e884886fSBarry Smith 
9531d0fab5eSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction()
9541d0fab5eSBarry Smith 
955e884886fSBarry Smith @*/
9567087cfbeSBarry Smith PetscErrorCode  MatMFFDSetFunctioni(Mat mat,PetscErrorCode (*funci)(void*,PetscInt,Vec,PetscScalar*))
957e884886fSBarry Smith {
9584ac538c5SBarry Smith   PetscErrorCode ierr;
959e884886fSBarry Smith 
960e884886fSBarry Smith   PetscFunctionBegin;
9610700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
9624ac538c5SBarry Smith   ierr = PetscTryMethod(mat,"MatMFFDSetFunctioni_C",(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)),(mat,funci));CHKERRQ(ierr);
963e884886fSBarry Smith   PetscFunctionReturn(0);
964e884886fSBarry Smith }
965e884886fSBarry Smith 
966e884886fSBarry Smith 
967e884886fSBarry Smith #undef __FUNCT__
968e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetFunctioniBase"
969e884886fSBarry Smith /*@C
970e884886fSBarry Smith    MatMFFDSetFunctioniBase - Sets the base vector for a single component function evaluation
971e884886fSBarry Smith 
9723f9fe445SBarry Smith    Logically Collective on Mat
973e884886fSBarry Smith 
974e884886fSBarry Smith    Input Parameters:
975e884886fSBarry Smith +  mat - the matrix free matrix created via MatCreateSNESMF()
976e884886fSBarry Smith -  func - the function to use
977e884886fSBarry Smith 
978e884886fSBarry Smith    Level: advanced
979e884886fSBarry Smith 
980e884886fSBarry Smith    Notes:
981e884886fSBarry Smith     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
982e884886fSBarry Smith     matrix inside your compute Jacobian routine
983e884886fSBarry Smith 
984e884886fSBarry Smith 
985e884886fSBarry Smith .keywords: SNES, matrix-free, function
986e884886fSBarry Smith 
987e884886fSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD
9881d0fab5eSBarry Smith           MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction()
989e884886fSBarry Smith @*/
9907087cfbeSBarry Smith PetscErrorCode  MatMFFDSetFunctioniBase(Mat mat,PetscErrorCode (*func)(void*,Vec))
991e884886fSBarry Smith {
9924ac538c5SBarry Smith   PetscErrorCode ierr;
993e884886fSBarry Smith 
994e884886fSBarry Smith   PetscFunctionBegin;
9950700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
9964ac538c5SBarry Smith   ierr = PetscTryMethod(mat,"MatMFFDSetFunctioniBase_C",(Mat,PetscErrorCode (*)(void*,Vec)),(mat,func));CHKERRQ(ierr);
997e884886fSBarry Smith   PetscFunctionReturn(0);
998e884886fSBarry Smith }
999e884886fSBarry Smith 
1000e884886fSBarry Smith #undef __FUNCT__
1001e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetPeriod"
1002e884886fSBarry Smith /*@
1003e884886fSBarry Smith    MatMFFDSetPeriod - Sets how often h is recomputed, by default it is everytime
1004e884886fSBarry Smith 
10053f9fe445SBarry Smith    Logically Collective on Mat
1006e884886fSBarry Smith 
1007e884886fSBarry Smith    Input Parameters:
1008e884886fSBarry Smith +  mat - the matrix free matrix created via MatCreateSNESMF()
1009e884886fSBarry Smith -  period - 1 for everytime, 2 for every second etc
1010e884886fSBarry Smith 
1011e884886fSBarry Smith    Options Database Keys:
1012e884886fSBarry Smith +  -mat_mffd_period <period>
1013e884886fSBarry Smith 
1014e884886fSBarry Smith    Level: advanced
1015e884886fSBarry Smith 
1016e884886fSBarry Smith 
1017e884886fSBarry Smith .keywords: SNES, matrix-free, parameters
1018e884886fSBarry Smith 
1019e884886fSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(),
10201d0fab5eSBarry Smith           MatMFFDSetHHistory(), MatMFFDResetHHistory()
1021e884886fSBarry Smith @*/
10227087cfbeSBarry Smith PetscErrorCode  MatMFFDSetPeriod(Mat mat,PetscInt period)
1023e884886fSBarry Smith {
1024bc0f08ceSBarry Smith   PetscErrorCode ierr;
1025e884886fSBarry Smith 
1026e884886fSBarry Smith   PetscFunctionBegin;
1027bc0f08ceSBarry Smith   ierr = PetscTryMethod(mat,"MatMFFDSetPeriod_C",(Mat,PetscInt),(mat,period));CHKERRQ(ierr);
1028e884886fSBarry Smith   PetscFunctionReturn(0);
1029e884886fSBarry Smith }
1030e884886fSBarry Smith 
1031e884886fSBarry Smith #undef __FUNCT__
1032e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetFunctionError"
1033e884886fSBarry Smith /*@
1034e884886fSBarry Smith    MatMFFDSetFunctionError - Sets the error_rel for the approximation of
1035e884886fSBarry Smith    matrix-vector products using finite differences.
1036e884886fSBarry Smith 
10373f9fe445SBarry Smith    Logically Collective on Mat
1038e884886fSBarry Smith 
1039e884886fSBarry Smith    Input Parameters:
1040e884886fSBarry Smith +  mat - the matrix free matrix created via MatCreateMFFD() or MatCreateSNESMF()
1041e884886fSBarry Smith -  error_rel - relative error (should be set to the square root of
1042e884886fSBarry Smith                the relative error in the function evaluations)
1043e884886fSBarry Smith 
1044e884886fSBarry Smith    Options Database Keys:
1045e884886fSBarry Smith +  -mat_mffd_err <error_rel> - Sets error_rel
1046e884886fSBarry Smith 
1047e884886fSBarry Smith    Level: advanced
1048e884886fSBarry Smith 
1049e884886fSBarry Smith    Notes:
1050e884886fSBarry Smith    The default matrix-free matrix-vector product routine computes
1051e884886fSBarry Smith .vb
1052e884886fSBarry Smith      F'(u)*a = [F(u+h*a) - F(u)]/h where
1053e884886fSBarry Smith      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
1054e884886fSBarry Smith        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   else
1055e884886fSBarry Smith .ve
1056e884886fSBarry Smith 
1057e884886fSBarry Smith .keywords: SNES, matrix-free, parameters
1058e884886fSBarry Smith 
1059e884886fSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD
10601d0fab5eSBarry Smith           MatMFFDSetHHistory(), MatMFFDResetHHistory()
1061e884886fSBarry Smith @*/
10627087cfbeSBarry Smith PetscErrorCode  MatMFFDSetFunctionError(Mat mat,PetscReal error)
1063e884886fSBarry Smith {
1064bc0f08ceSBarry Smith   PetscErrorCode ierr;
1065e884886fSBarry Smith 
1066e884886fSBarry Smith   PetscFunctionBegin;
1067bc0f08ceSBarry Smith   ierr = PetscTryMethod(mat,"MatMFFDSetFunctionError_C",(Mat,PetscReal),(mat,error));CHKERRQ(ierr);
1068e884886fSBarry Smith   PetscFunctionReturn(0);
1069e884886fSBarry Smith }
1070e884886fSBarry Smith 
1071e884886fSBarry Smith #undef __FUNCT__
1072e884886fSBarry Smith #define __FUNCT__ "MatMFFDAddNullSpace"
1073e884886fSBarry Smith /*@
1074e884886fSBarry Smith    MatMFFDAddNullSpace - Provides a null space that an operator is
1075e884886fSBarry Smith    supposed to have.  Since roundoff will create a small component in
1076e884886fSBarry Smith    the null space, if you know the null space you may have it
1077e884886fSBarry Smith    automatically removed.
1078e884886fSBarry Smith 
10793f9fe445SBarry Smith    Logically Collective on Mat
1080e884886fSBarry Smith 
1081e884886fSBarry Smith    Input Parameters:
1082e884886fSBarry Smith +  J - the matrix-free matrix context
1083e884886fSBarry Smith -  nullsp - object created with MatNullSpaceCreate()
1084e884886fSBarry Smith 
1085e884886fSBarry Smith    Level: advanced
1086e884886fSBarry Smith 
1087e884886fSBarry Smith .keywords: SNES, matrix-free, null space
1088e884886fSBarry Smith 
1089e884886fSBarry Smith .seealso: MatNullSpaceCreate(), MatMFFDGetH(), MatCreateSNESMF(), MatCreateMFFD(), MATMFFD
10901d0fab5eSBarry Smith           MatMFFDSetHHistory(), MatMFFDResetHHistory()
1091e884886fSBarry Smith @*/
10927087cfbeSBarry Smith PetscErrorCode  MatMFFDAddNullSpace(Mat J,MatNullSpace nullsp)
1093e884886fSBarry Smith {
1094e884886fSBarry Smith   PetscErrorCode ierr;
1095e884886fSBarry Smith 
1096e884886fSBarry Smith   PetscFunctionBegin;
1097bc0f08ceSBarry Smith   ierr = PetscTryMethod(J,"MatMFFDAddNullSpace_C",(Mat,MatNullSpace),(J,nullsp));CHKERRQ(ierr);
1098e884886fSBarry Smith   PetscFunctionReturn(0);
1099e884886fSBarry Smith }
1100e884886fSBarry Smith 
1101e884886fSBarry Smith #undef __FUNCT__
1102e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetHHistory"
1103e884886fSBarry Smith /*@
1104e884886fSBarry Smith    MatMFFDSetHHistory - Sets an array to collect a history of the
1105e884886fSBarry Smith    differencing values (h) computed for the matrix-free product.
1106e884886fSBarry Smith 
11073f9fe445SBarry Smith    Logically Collective on Mat
1108e884886fSBarry Smith 
1109e884886fSBarry Smith    Input Parameters:
1110e884886fSBarry Smith +  J - the matrix-free matrix context
1111e884886fSBarry Smith .  histroy - space to hold the history
1112e884886fSBarry Smith -  nhistory - number of entries in history, if more entries are generated than
1113e884886fSBarry Smith               nhistory, then the later ones are discarded
1114e884886fSBarry Smith 
1115e884886fSBarry Smith    Level: advanced
1116e884886fSBarry Smith 
1117e884886fSBarry Smith    Notes:
1118e884886fSBarry Smith    Use MatMFFDResetHHistory() to reset the history counter and collect
1119e884886fSBarry Smith    a new batch of differencing parameters, h.
1120e884886fSBarry Smith 
1121e884886fSBarry Smith .keywords: SNES, matrix-free, h history, differencing history
1122e884886fSBarry Smith 
1123e884886fSBarry Smith .seealso: MatMFFDGetH(), MatCreateSNESMF(),
11241d0fab5eSBarry Smith           MatMFFDResetHHistory(), MatMFFDSetFunctionError()
1125e884886fSBarry Smith 
1126e884886fSBarry Smith @*/
11277087cfbeSBarry Smith PetscErrorCode  MatMFFDSetHHistory(Mat J,PetscScalar history[],PetscInt nhistory)
1128e884886fSBarry Smith {
1129e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)J->data;
1130bc0f08ceSBarry Smith   PetscErrorCode ierr;
1131bc0f08ceSBarry Smith   PetscBool      match;
1132e884886fSBarry Smith 
1133e884886fSBarry Smith   PetscFunctionBegin;
1134251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)J,MATMFFD,&match);CHKERRQ(ierr);
1135ce94432eSBarry Smith   if (!match) SETERRQ(PetscObjectComm((PetscObject)J),PETSC_ERR_ARG_WRONG,"Not a MFFD matrix");
1136e884886fSBarry Smith   ctx->historyh    = history;
1137e884886fSBarry Smith   ctx->maxcurrenth = nhistory;
113875567043SBarry Smith   ctx->currenth    = 0.;
1139e884886fSBarry Smith   PetscFunctionReturn(0);
1140e884886fSBarry Smith }
1141e884886fSBarry Smith 
1142bc0f08ceSBarry Smith 
1143e884886fSBarry Smith #undef __FUNCT__
1144e884886fSBarry Smith #define __FUNCT__ "MatMFFDResetHHistory"
1145e884886fSBarry Smith /*@
1146e884886fSBarry Smith    MatMFFDResetHHistory - Resets the counter to zero to begin
1147e884886fSBarry Smith    collecting a new set of differencing histories.
1148e884886fSBarry Smith 
11493f9fe445SBarry Smith    Logically Collective on Mat
1150e884886fSBarry Smith 
1151e884886fSBarry Smith    Input Parameters:
1152e884886fSBarry Smith .  J - the matrix-free matrix context
1153e884886fSBarry Smith 
1154e884886fSBarry Smith    Level: advanced
1155e884886fSBarry Smith 
1156e884886fSBarry Smith    Notes:
1157e884886fSBarry Smith    Use MatMFFDSetHHistory() to create the original history counter.
1158e884886fSBarry Smith 
1159e884886fSBarry Smith .keywords: SNES, matrix-free, h history, differencing history
1160e884886fSBarry Smith 
1161e884886fSBarry Smith .seealso: MatMFFDGetH(), MatCreateSNESMF(),
11621d0fab5eSBarry Smith           MatMFFDSetHHistory(), MatMFFDSetFunctionError()
1163e884886fSBarry Smith 
1164e884886fSBarry Smith @*/
11657087cfbeSBarry Smith PetscErrorCode  MatMFFDResetHHistory(Mat J)
1166e884886fSBarry Smith {
1167bc0f08ceSBarry Smith   PetscErrorCode ierr;
1168e884886fSBarry Smith 
1169e884886fSBarry Smith   PetscFunctionBegin;
1170bc0f08ceSBarry Smith   ierr = PetscTryMethod(J,"MatMFFDResetHHistory_C",(Mat),(J));CHKERRQ(ierr);
1171e884886fSBarry Smith   PetscFunctionReturn(0);
1172e884886fSBarry Smith }
1173e884886fSBarry Smith 
1174e884886fSBarry Smith 
1175e884886fSBarry Smith #undef __FUNCT__
1176e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetBase"
1177e884886fSBarry Smith /*@
1178e884886fSBarry Smith     MatMFFDSetBase - Sets the vector U at which matrix vector products of the
1179e884886fSBarry Smith         Jacobian are computed
1180e884886fSBarry Smith 
11813f9fe445SBarry Smith     Logically Collective on Mat
1182e884886fSBarry Smith 
1183e884886fSBarry Smith     Input Parameters:
1184e884886fSBarry Smith +   J - the MatMFFD matrix
11853ec795f1SBarry Smith .   U - the vector
1186bcddec3dSBarry Smith -   F - (optional) vector that contains F(u) if it has been already computed
1187e884886fSBarry Smith 
1188e884886fSBarry Smith     Notes: This is rarely used directly
1189e884886fSBarry Smith 
11908af5ae88SBarry Smith     If F is provided then it is not recomputed. Otherwise the function is evaluated at the base
11918af5ae88SBarry Smith     point during the first MatMult() after each call to MatMFFDSetBase().
1192dff2f722SBarry Smith 
1193e884886fSBarry Smith     Level: advanced
1194e884886fSBarry Smith 
1195e884886fSBarry Smith @*/
11967087cfbeSBarry Smith PetscErrorCode  MatMFFDSetBase(Mat J,Vec U,Vec F)
1197e884886fSBarry Smith {
11984ac538c5SBarry Smith   PetscErrorCode ierr;
1199e884886fSBarry Smith 
1200e884886fSBarry Smith   PetscFunctionBegin;
12010700a824SBarry Smith   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
12020700a824SBarry Smith   PetscValidHeaderSpecific(U,VEC_CLASSID,2);
12030700a824SBarry Smith   if (F) PetscValidHeaderSpecific(F,VEC_CLASSID,3);
12044ac538c5SBarry Smith   ierr = PetscTryMethod(J,"MatMFFDSetBase_C",(Mat,Vec,Vec),(J,U,F));CHKERRQ(ierr);
1205e884886fSBarry Smith   PetscFunctionReturn(0);
1206e884886fSBarry Smith }
1207e884886fSBarry Smith 
1208e884886fSBarry Smith #undef __FUNCT__
1209e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetCheckh"
1210e884886fSBarry Smith /*@C
1211e884886fSBarry Smith     MatMFFDSetCheckh - Sets a function that checks the computed h and adjusts
1212e884886fSBarry Smith         it to satisfy some criteria
1213e884886fSBarry Smith 
12143f9fe445SBarry Smith     Logically Collective on Mat
1215e884886fSBarry Smith 
1216e884886fSBarry Smith     Input Parameters:
1217e884886fSBarry Smith +   J - the MatMFFD matrix
1218e884886fSBarry Smith .   fun - the function that checks h
1219e884886fSBarry Smith -   ctx - any context needed by the function
1220e884886fSBarry Smith 
1221e884886fSBarry Smith     Options Database Keys:
1222e884886fSBarry Smith .   -mat_mffd_check_positivity
1223e884886fSBarry Smith 
1224e884886fSBarry Smith     Level: advanced
1225e884886fSBarry Smith 
1226e884886fSBarry Smith     Notes: For example, MatMFFDSetCheckPositivity() insures that all entries
1227e884886fSBarry Smith        of U + h*a are non-negative
1228e884886fSBarry Smith 
1229e884886fSBarry Smith .seealso:  MatMFFDSetCheckPositivity()
1230e884886fSBarry Smith @*/
12317087cfbeSBarry Smith PetscErrorCode  MatMFFDSetCheckh(Mat J,PetscErrorCode (*fun)(void*,Vec,Vec,PetscScalar*),void *ctx)
1232e884886fSBarry Smith {
12334ac538c5SBarry Smith   PetscErrorCode ierr;
1234e884886fSBarry Smith 
1235e884886fSBarry Smith   PetscFunctionBegin;
12360700a824SBarry Smith   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
12374ac538c5SBarry Smith   ierr = PetscTryMethod(J,"MatMFFDSetCheckh_C",(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*),(J,fun,ctx));CHKERRQ(ierr);
1238e884886fSBarry Smith   PetscFunctionReturn(0);
1239e884886fSBarry Smith }
1240e884886fSBarry Smith 
1241e884886fSBarry Smith #undef __FUNCT__
1242e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetCheckPositivity"
1243e884886fSBarry Smith /*@
1244e884886fSBarry Smith     MatMFFDCheckPositivity - Checks that all entries in U + h*a are positive or
1245e884886fSBarry Smith         zero, decreases h until this is satisfied.
1246e884886fSBarry Smith 
12473f9fe445SBarry Smith     Logically Collective on Vec
1248e884886fSBarry Smith 
1249e884886fSBarry Smith     Input Parameters:
1250e884886fSBarry Smith +   U - base vector that is added to
1251e884886fSBarry Smith .   a - vector that is added
1252e884886fSBarry Smith .   h - scaling factor on a
1253e884886fSBarry Smith -   dummy - context variable (unused)
1254e884886fSBarry Smith 
1255e884886fSBarry Smith     Options Database Keys:
1256e884886fSBarry Smith .   -mat_mffd_check_positivity
1257e884886fSBarry Smith 
1258e884886fSBarry Smith     Level: advanced
1259e884886fSBarry Smith 
1260e884886fSBarry Smith     Notes: This is rarely used directly, rather it is passed as an argument to
1261e884886fSBarry Smith            MatMFFDSetCheckh()
1262e884886fSBarry Smith 
1263e884886fSBarry Smith .seealso:  MatMFFDSetCheckh()
1264e884886fSBarry Smith @*/
12657087cfbeSBarry Smith PetscErrorCode  MatMFFDCheckPositivity(void *dummy,Vec U,Vec a,PetscScalar *h)
1266e884886fSBarry Smith {
1267e884886fSBarry Smith   PetscReal      val, minval;
1268e884886fSBarry Smith   PetscScalar    *u_vec, *a_vec;
1269e884886fSBarry Smith   PetscErrorCode ierr;
1270e884886fSBarry Smith   PetscInt       i,n;
1271e884886fSBarry Smith   MPI_Comm       comm;
1272e884886fSBarry Smith 
1273e884886fSBarry Smith   PetscFunctionBegin;
1274e884886fSBarry Smith   ierr   = PetscObjectGetComm((PetscObject)U,&comm);CHKERRQ(ierr);
1275e884886fSBarry Smith   ierr   = VecGetArray(U,&u_vec);CHKERRQ(ierr);
1276e884886fSBarry Smith   ierr   = VecGetArray(a,&a_vec);CHKERRQ(ierr);
1277e884886fSBarry Smith   ierr   = VecGetLocalSize(U,&n);CHKERRQ(ierr);
1278e884886fSBarry Smith   minval = PetscAbsScalar(*h*1.01);
1279e884886fSBarry Smith   for (i=0; i<n; i++) {
1280e884886fSBarry Smith     if (PetscRealPart(u_vec[i] + *h*a_vec[i]) <= 0.0) {
1281e884886fSBarry Smith       val = PetscAbsScalar(u_vec[i]/a_vec[i]);
1282e884886fSBarry Smith       if (val < minval) minval = val;
1283e884886fSBarry Smith     }
1284e884886fSBarry Smith   }
1285e884886fSBarry Smith   ierr = VecRestoreArray(U,&u_vec);CHKERRQ(ierr);
1286e884886fSBarry Smith   ierr = VecRestoreArray(a,&a_vec);CHKERRQ(ierr);
1287d9822059SBarry Smith   ierr = MPI_Allreduce(&minval,&val,1,MPIU_REAL,MPIU_MIN,comm);CHKERRQ(ierr);
1288e884886fSBarry Smith   if (val <= PetscAbsScalar(*h)) {
1289e884886fSBarry Smith     ierr = PetscInfo2(U,"Scaling back h from %G to %G\n",PetscRealPart(*h),.99*val);CHKERRQ(ierr);
1290e884886fSBarry Smith     if (PetscRealPart(*h) > 0.0) *h =  0.99*val;
1291e884886fSBarry Smith     else                         *h = -0.99*val;
1292e884886fSBarry Smith   }
1293e884886fSBarry Smith   PetscFunctionReturn(0);
1294e884886fSBarry Smith }
1295e884886fSBarry Smith 
1296e884886fSBarry Smith 
1297e884886fSBarry Smith 
1298e884886fSBarry Smith 
1299e884886fSBarry Smith 
1300e884886fSBarry Smith 
1301e884886fSBarry Smith 
1302