xref: /petsc/src/mat/impls/mffd/mffd.c (revision b45d2f2cb7e031d9c0de5873eca80614ca7b863b)
1e884886fSBarry Smith 
2*b45d2f2cSJed Brown #include <petsc-private/matimpl.h>
3c6db04a5SJed Brown #include <../src/mat/impls/mffd/mffdimpl.h>   /*I  "petscmat.h"   I*/
4e884886fSBarry Smith 
5b022a5c1SBarry Smith PetscFList 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 {
25b022a5c1SBarry Smith   PetscFunctionBegin;
26b022a5c1SBarry Smith   MatMFFDPackageInitialized = PETSC_FALSE;
27b022a5c1SBarry Smith   MatMFFDRegisterAllCalled  = PETSC_FALSE;
28b022a5c1SBarry Smith   MatMFFDList               = PETSC_NULL;
29b022a5c1SBarry Smith   PetscFunctionReturn(0);
30b022a5c1SBarry Smith }
31b022a5c1SBarry Smith 
32e884886fSBarry Smith #undef __FUNCT__
333ec795f1SBarry Smith #define __FUNCT__ "MatMFFDInitializePackage"
343ec795f1SBarry Smith /*@C
353ec795f1SBarry Smith   MatMFFDInitializePackage - This function initializes everything in the MatMFFD package. It is called
363ec795f1SBarry Smith   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to MatCreate_MFFD()
373ec795f1SBarry Smith   when using static libraries.
383ec795f1SBarry Smith 
393ec795f1SBarry Smith   Input Parameter:
403ec795f1SBarry Smith . path - The dynamic library path, or PETSC_NULL
413ec795f1SBarry Smith 
423ec795f1SBarry Smith   Level: developer
433ec795f1SBarry Smith 
443ec795f1SBarry Smith .keywords: Vec, initialize, package
453ec795f1SBarry Smith .seealso: PetscInitialize()
463ec795f1SBarry Smith @*/
477087cfbeSBarry Smith PetscErrorCode  MatMFFDInitializePackage(const char path[])
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 */
603ec795f1SBarry Smith   ierr = MatMFFDRegisterAll(path);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 */
653ec795f1SBarry Smith   ierr = PetscOptionsGetString(PETSC_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 */
733ec795f1SBarry Smith   ierr = PetscOptionsGetString(PETSC_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 
105d2d6cebeSBarry Smith .seealso: MatCreateSNESMF(), MatMFFDRegisterDynamic(), MatMFFDSetFunction()
106e884886fSBarry Smith @*/
1077087cfbeSBarry Smith PetscErrorCode  MatMFFDSetType(Mat mat,const 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 
1179a13eae7SJed Brown   ierr = PetscTypeCompare((PetscObject)mat,MATMFFD,&match);CHKERRQ(ierr);
1189a13eae7SJed Brown   if (!match) PetscFunctionReturn(0);
1199a13eae7SJed Brown 
120e884886fSBarry Smith   /* already set, so just return */
121e884886fSBarry Smith   ierr = PetscTypeCompare((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 
1294b91b6eaSBarry Smith   ierr =  PetscFListFind(MatMFFDList,((PetscObject)ctx)->comm,ftype,PETSC_TRUE,(void (**)(void)) &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 EXTERN_C_BEGIN
138e884886fSBarry Smith #undef __FUNCT__
139c879e56bSBarry Smith #define __FUNCT__ "MatMFFDSetFunctioniBase_MFFD"
1407087cfbeSBarry Smith PetscErrorCode  MatMFFDSetFunctioniBase_MFFD(Mat mat,FCN1 func)
141e884886fSBarry Smith {
142e884886fSBarry Smith   MatMFFD ctx = (MatMFFD)mat->data;
143e884886fSBarry Smith 
144e884886fSBarry Smith   PetscFunctionBegin;
145e884886fSBarry Smith   ctx->funcisetbase = func;
146e884886fSBarry Smith   PetscFunctionReturn(0);
147e884886fSBarry Smith }
148e884886fSBarry Smith EXTERN_C_END
149e884886fSBarry Smith 
150e884886fSBarry Smith typedef PetscErrorCode (*FCN2)(void*,PetscInt,Vec,PetscScalar*); /* force argument to next function to not be extern C*/
151e884886fSBarry Smith EXTERN_C_BEGIN
152e884886fSBarry Smith #undef __FUNCT__
153c879e56bSBarry Smith #define __FUNCT__ "MatMFFDSetFunctioni_MFFD"
1547087cfbeSBarry Smith PetscErrorCode  MatMFFDSetFunctioni_MFFD(Mat mat,FCN2 funci)
155e884886fSBarry Smith {
156e884886fSBarry Smith   MatMFFD ctx = (MatMFFD)mat->data;
157e884886fSBarry Smith 
158e884886fSBarry Smith   PetscFunctionBegin;
159e884886fSBarry Smith   ctx->funci = funci;
160e884886fSBarry Smith   PetscFunctionReturn(0);
161e884886fSBarry Smith }
162e884886fSBarry Smith EXTERN_C_END
163e884886fSBarry Smith 
164bc0f08ceSBarry Smith EXTERN_C_BEGIN
165bc0f08ceSBarry Smith #undef __FUNCT__
166bc0f08ceSBarry Smith #define __FUNCT__ "MatMFFDResetHHistory_MFFD"
167bc0f08ceSBarry Smith PetscErrorCode  MatMFFDResetHHistory_MFFD(Mat J)
168bc0f08ceSBarry Smith {
169bc0f08ceSBarry Smith   MatMFFD ctx = (MatMFFD)J->data;
170bc0f08ceSBarry Smith 
171bc0f08ceSBarry Smith   PetscFunctionBegin;
172bc0f08ceSBarry Smith   ctx->ncurrenth    = 0;
173bc0f08ceSBarry Smith   PetscFunctionReturn(0);
174bc0f08ceSBarry Smith }
175bc0f08ceSBarry Smith EXTERN_C_END
176e884886fSBarry Smith 
177e884886fSBarry Smith #undef __FUNCT__
178e884886fSBarry Smith #define __FUNCT__ "MatMFFDRegister"
1797087cfbeSBarry Smith PetscErrorCode  MatMFFDRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(MatMFFD))
180e884886fSBarry Smith {
181e884886fSBarry Smith   PetscErrorCode ierr;
182e884886fSBarry Smith   char           fullname[PETSC_MAX_PATH_LEN];
183e884886fSBarry Smith 
184e884886fSBarry Smith   PetscFunctionBegin;
185e884886fSBarry Smith   ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
186b022a5c1SBarry Smith   ierr = PetscFListAdd(&MatMFFDList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
187e884886fSBarry Smith   PetscFunctionReturn(0);
188e884886fSBarry Smith }
189e884886fSBarry Smith 
190e884886fSBarry Smith 
191e884886fSBarry Smith #undef __FUNCT__
192e884886fSBarry Smith #define __FUNCT__ "MatMFFDRegisterDestroy"
193e884886fSBarry Smith /*@C
194e884886fSBarry Smith    MatMFFDRegisterDestroy - Frees the list of MatMFFD methods that were
195e884886fSBarry Smith    registered by MatMFFDRegisterDynamic).
196e884886fSBarry Smith 
197e884886fSBarry Smith    Not Collective
198e884886fSBarry Smith 
199e884886fSBarry Smith    Level: developer
200e884886fSBarry Smith 
201e884886fSBarry Smith .keywords: MatMFFD, register, destroy
202e884886fSBarry Smith 
203e884886fSBarry Smith .seealso: MatMFFDRegisterDynamic), MatMFFDRegisterAll()
204e884886fSBarry Smith @*/
2057087cfbeSBarry Smith PetscErrorCode  MatMFFDRegisterDestroy(void)
206e884886fSBarry Smith {
207e884886fSBarry Smith   PetscErrorCode ierr;
208e884886fSBarry Smith 
209e884886fSBarry Smith   PetscFunctionBegin;
210b022a5c1SBarry Smith   ierr = PetscFListDestroy(&MatMFFDList);CHKERRQ(ierr);
211e884886fSBarry Smith   MatMFFDRegisterAllCalled = PETSC_FALSE;
212e884886fSBarry Smith   PetscFunctionReturn(0);
213e884886fSBarry Smith }
214e884886fSBarry Smith 
215bc0f08ceSBarry Smith EXTERN_C_BEGIN
216bc0f08ceSBarry Smith #undef __FUNCT__
217bc0f08ceSBarry Smith #define __FUNCT__ "MatMFFDAddNullSpace_MFFD"
218bc0f08ceSBarry Smith PetscErrorCode  MatMFFDAddNullSpace_MFFD(Mat J,MatNullSpace nullsp)
219bc0f08ceSBarry Smith {
220bc0f08ceSBarry Smith   PetscErrorCode ierr;
221bc0f08ceSBarry Smith   MatMFFD      ctx = (MatMFFD)J->data;
222bc0f08ceSBarry Smith 
223bc0f08ceSBarry Smith   PetscFunctionBegin;
224bc0f08ceSBarry Smith   ierr = PetscObjectReference((PetscObject)nullsp);CHKERRQ(ierr);
225bc0f08ceSBarry Smith   if (ctx->sp) { ierr = MatNullSpaceDestroy(&ctx->sp);CHKERRQ(ierr); }
226bc0f08ceSBarry Smith   ctx->sp = nullsp;
227bc0f08ceSBarry Smith   PetscFunctionReturn(0);
228bc0f08ceSBarry Smith }
229bc0f08ceSBarry Smith EXTERN_C_END
230bc0f08ceSBarry Smith 
231e884886fSBarry Smith /* ----------------------------------------------------------------------------------------*/
232e884886fSBarry Smith #undef __FUNCT__
233e884886fSBarry Smith #define __FUNCT__ "MatDestroy_MFFD"
234e884886fSBarry Smith PetscErrorCode MatDestroy_MFFD(Mat mat)
235e884886fSBarry Smith {
236e884886fSBarry Smith   PetscErrorCode ierr;
237e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)mat->data;
238e884886fSBarry Smith 
239e884886fSBarry Smith   PetscFunctionBegin;
2406bf464f9SBarry Smith   ierr = VecDestroy(&ctx->w);CHKERRQ(ierr);
2416bf464f9SBarry Smith   ierr = VecDestroy(&ctx->drscale);CHKERRQ(ierr);
2426bf464f9SBarry Smith   ierr = VecDestroy(&ctx->dlscale);CHKERRQ(ierr);
2436bf464f9SBarry Smith   ierr = VecDestroy(&ctx->dshift);CHKERRQ(ierr);
244cfe22f5eSBarry Smith   if (ctx->current_f_allocated) {
2456bf464f9SBarry Smith     ierr = VecDestroy(&ctx->current_f);CHKERRQ(ierr);
246cfe22f5eSBarry Smith   }
247e884886fSBarry Smith   if (ctx->ops->destroy) {ierr = (*ctx->ops->destroy)(ctx);CHKERRQ(ierr);}
2486bf464f9SBarry Smith   ierr = MatNullSpaceDestroy(&ctx->sp);CHKERRQ(ierr);
2496bf464f9SBarry Smith   ierr = PetscHeaderDestroy(&ctx);CHKERRQ(ierr);
250bf0cc555SLisandro Dalcin   mat->data = 0;
251e884886fSBarry Smith 
252e884886fSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetBase_C","",PETSC_NULL);CHKERRQ(ierr);
253e884886fSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctioniBase_C","",PETSC_NULL);CHKERRQ(ierr);
254e884886fSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctioni_C","",PETSC_NULL);CHKERRQ(ierr);
255bc0f08ceSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunction_C","",PETSC_NULL);CHKERRQ(ierr);
256bc0f08ceSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctionError_C","",PETSC_NULL);CHKERRQ(ierr);
257e884886fSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetCheckh_C","",PETSC_NULL);CHKERRQ(ierr);
258bc0f08ceSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetPeriod_C","",PETSC_NULL);CHKERRQ(ierr);
259bc0f08ceSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDResetHHistory_C","",PETSC_NULL);CHKERRQ(ierr);
260bc0f08ceSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDAddNullSpace_C","",PETSC_NULL);CHKERRQ(ierr);
261e884886fSBarry Smith 
262e884886fSBarry Smith   PetscFunctionReturn(0);
263e884886fSBarry Smith }
264e884886fSBarry Smith 
265e884886fSBarry Smith #undef __FUNCT__
266e884886fSBarry Smith #define __FUNCT__ "MatView_MFFD"
267e884886fSBarry Smith /*
268e884886fSBarry Smith    MatMFFDView_MFFD - Views matrix-free parameters.
269e884886fSBarry Smith 
270e884886fSBarry Smith */
271e884886fSBarry Smith PetscErrorCode MatView_MFFD(Mat J,PetscViewer viewer)
272e884886fSBarry Smith {
273e884886fSBarry Smith   PetscErrorCode ierr;
274e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)J->data;
275a283635eSDmitry Karpeev   PetscBool      iascii, viewbase, viewfunction;
276a283635eSDmitry Karpeev   const char*    prefix;
277e884886fSBarry Smith 
278e884886fSBarry Smith   PetscFunctionBegin;
2792692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
280e884886fSBarry Smith   if (iascii) {
281a283635eSDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"Matrix-free approximation:\n");CHKERRQ(ierr);
282a283635eSDmitry Karpeev     ierr = PetscViewerASCIIPushTab(viewer); CHKERRQ(ierr);
283e884886fSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"err=%G (relative error in function evaluation)\n",ctx->error_rel);CHKERRQ(ierr);
2847adad957SLisandro Dalcin     if (!((PetscObject)ctx)->type_name) {
285e884886fSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"The compute h routine has not yet been set\n");CHKERRQ(ierr);
286e884886fSBarry Smith     } else {
2877adad957SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"Using %s compute h routine\n",((PetscObject)ctx)->type_name);CHKERRQ(ierr);
288e884886fSBarry Smith     }
289e884886fSBarry Smith     if (ctx->ops->view) {
290e884886fSBarry Smith       ierr = (*ctx->ops->view)(ctx,viewer);CHKERRQ(ierr);
291e884886fSBarry Smith     }
292a283635eSDmitry Karpeev     ierr = PetscObjectGetOptionsPrefix((PetscObject)J, &prefix); CHKERRQ(ierr);
293a283635eSDmitry Karpeev 
294a283635eSDmitry Karpeev     ierr = PetscOptionsHasName(prefix, "-mat_mffd_view_base", &viewbase); CHKERRQ(ierr);
295a283635eSDmitry Karpeev     if(viewbase) {
296a283635eSDmitry Karpeev       ierr = PetscViewerASCIIPrintf(viewer, "Base:\n");     CHKERRQ(ierr);
297a283635eSDmitry Karpeev       ierr = VecView(ctx->current_u, viewer);                 CHKERRQ(ierr);
298a283635eSDmitry Karpeev     }
299a283635eSDmitry Karpeev     ierr = PetscOptionsHasName(prefix, "-mat_mffd_view_function", &viewfunction); CHKERRQ(ierr);
300a283635eSDmitry Karpeev     if(viewfunction) {
301a283635eSDmitry Karpeev       ierr = PetscViewerASCIIPrintf(viewer, "Function:\n"); CHKERRQ(ierr);
302a283635eSDmitry Karpeev       ierr = VecView(ctx->current_f, viewer);                 CHKERRQ(ierr);
303a283635eSDmitry Karpeev     }
304a283635eSDmitry Karpeev     ierr = PetscViewerASCIIPopTab(viewer); CHKERRQ(ierr);
305a8248277SBarry Smith   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for matrix-free matrix",((PetscObject)viewer)->type_name);
306e884886fSBarry Smith   PetscFunctionReturn(0);
307e884886fSBarry Smith }
308e884886fSBarry Smith 
309e884886fSBarry Smith #undef __FUNCT__
310e884886fSBarry Smith #define __FUNCT__ "MatAssemblyEnd_MFFD"
311e884886fSBarry Smith /*
312e884886fSBarry Smith    MatAssemblyEnd_MFFD - Resets the ctx->ncurrenth to zero. This
313e884886fSBarry Smith    allows the user to indicate the beginning of a new linear solve by calling
314e884886fSBarry Smith    MatAssemblyXXX() on the matrix free matrix. This then allows the
3151d0fab5eSBarry Smith    MatCreateMFFD_WP() to properly compute ||U|| only the first time
316e884886fSBarry Smith    in the linear solver rather than every time.
317e884886fSBarry Smith */
318e884886fSBarry Smith PetscErrorCode MatAssemblyEnd_MFFD(Mat J,MatAssemblyType mt)
319e884886fSBarry Smith {
320e884886fSBarry Smith   PetscErrorCode ierr;
321e884886fSBarry Smith   MatMFFD        j = (MatMFFD)J->data;
322e884886fSBarry Smith 
323e884886fSBarry Smith   PetscFunctionBegin;
324e884886fSBarry Smith   ierr      = MatMFFDResetHHistory(J);CHKERRQ(ierr);
325e884886fSBarry Smith   j->vshift = 0.0;
326e884886fSBarry Smith   j->vscale = 1.0;
327e884886fSBarry Smith   PetscFunctionReturn(0);
328e884886fSBarry Smith }
329e884886fSBarry Smith 
330e884886fSBarry Smith #undef __FUNCT__
331e884886fSBarry Smith #define __FUNCT__ "MatMult_MFFD"
332e884886fSBarry Smith /*
333e884886fSBarry Smith   MatMult_MFFD - Default matrix-free form for Jacobian-vector product, y = F'(u)*a:
334e884886fSBarry Smith 
335e884886fSBarry Smith         y ~= (F(u + ha) - F(u))/h,
336e884886fSBarry Smith   where F = nonlinear function, as set by SNESSetFunction()
337e884886fSBarry Smith         u = current iterate
338e884886fSBarry Smith         h = difference interval
339e884886fSBarry Smith */
340e884886fSBarry Smith PetscErrorCode MatMult_MFFD(Mat mat,Vec a,Vec y)
341e884886fSBarry Smith {
342e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)mat->data;
343e884886fSBarry Smith   PetscScalar    h;
344e884886fSBarry Smith   Vec            w,U,F;
345e884886fSBarry Smith   PetscErrorCode ierr;
346ace3abfcSBarry Smith   PetscBool      zeroa;
347e884886fSBarry Smith 
348e884886fSBarry Smith   PetscFunctionBegin;
34984da1bf7SJed Brown   if (!ctx->current_u) SETERRQ(((PetscObject)mat)->comm,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");
350e884886fSBarry Smith   /* We log matrix-free matrix-vector products separately, so that we can
351e884886fSBarry Smith      separate the performance monitoring from the cases that use conventional
352e884886fSBarry Smith      storage.  We may eventually modify event logging to associate events
353e884886fSBarry Smith      with particular objects, hence alleviating the more general problem. */
354e884886fSBarry Smith   ierr = PetscLogEventBegin(MATMFFD_Mult,a,y,0,0);CHKERRQ(ierr);
355e884886fSBarry Smith 
356e884886fSBarry Smith   w    = ctx->w;
357e884886fSBarry Smith   U    = ctx->current_u;
3583ec795f1SBarry Smith   F    = ctx->current_f;
359e884886fSBarry Smith   /*
360e884886fSBarry Smith       Compute differencing parameter
361e884886fSBarry Smith   */
362e884886fSBarry Smith   if (!ctx->ops->compute) {
363e884886fSBarry Smith     ierr = MatMFFDSetType(mat,MATMFFD_WP);CHKERRQ(ierr);
3649c6ac3b3SBarry Smith     ierr = MatSetFromOptions(mat);CHKERRQ(ierr);
365e884886fSBarry Smith   }
366e884886fSBarry Smith   ierr = (*ctx->ops->compute)(ctx,U,a,&h,&zeroa);CHKERRQ(ierr);
367e884886fSBarry Smith   if (zeroa) {
368e884886fSBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
369e884886fSBarry Smith     PetscFunctionReturn(0);
370e884886fSBarry Smith   }
371e884886fSBarry Smith 
372e32f2f54SBarry Smith   if (PetscIsInfOrNanScalar(h)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Computed Nan differencing parameter h");
373e884886fSBarry Smith   if (ctx->checkh) {
374e884886fSBarry Smith     ierr = (*ctx->checkh)(ctx->checkhctx,U,a,&h);CHKERRQ(ierr);
375e884886fSBarry Smith   }
376e884886fSBarry Smith 
377e884886fSBarry Smith   /* keep a record of the current differencing parameter h */
378e884886fSBarry Smith   ctx->currenth = h;
379e884886fSBarry Smith #if defined(PETSC_USE_COMPLEX)
380e884886fSBarry Smith   ierr = PetscInfo2(mat,"Current differencing parameter: %G + %G i\n",PetscRealPart(h),PetscImaginaryPart(h));CHKERRQ(ierr);
381e884886fSBarry Smith #else
382e884886fSBarry Smith   ierr = PetscInfo1(mat,"Current differencing parameter: %15.12e\n",h);CHKERRQ(ierr);
383e884886fSBarry Smith #endif
384e884886fSBarry Smith   if (ctx->historyh && ctx->ncurrenth < ctx->maxcurrenth) {
385e884886fSBarry Smith     ctx->historyh[ctx->ncurrenth] = h;
386e884886fSBarry Smith   }
387e884886fSBarry Smith   ctx->ncurrenth++;
388e884886fSBarry Smith 
389e884886fSBarry Smith   /* w = u + ha */
390c3bb7e23SBarry Smith   if (ctx->drscale) {
391c3bb7e23SBarry Smith     ierr = VecPointwiseMult(ctx->drscale,a,U);CHKERRQ(ierr);
392c3bb7e23SBarry Smith     ierr = VecAYPX(U,h,w);CHKERRQ(ierr);
393c3bb7e23SBarry Smith   } else {
394e884886fSBarry Smith     ierr = VecWAXPY(w,h,a,U);CHKERRQ(ierr);
395c3bb7e23SBarry Smith   }
396e884886fSBarry Smith 
3971b797266SDmitry Karpeev   /* compute func(U) as base for differencing; only needed first time in and not when provided by user */
3981b797266SDmitry Karpeev   if (ctx->ncurrenth == 1 && ctx->current_f_allocated) {
399e884886fSBarry Smith     ierr = (*ctx->func)(ctx->funcctx,U,F);CHKERRQ(ierr);
40052121784SBarry Smith   }
401e884886fSBarry Smith   ierr = (*ctx->func)(ctx->funcctx,w,y);CHKERRQ(ierr);
402e884886fSBarry Smith 
403e884886fSBarry Smith   ierr = VecAXPY(y,-1.0,F);CHKERRQ(ierr);
404e884886fSBarry Smith   ierr = VecScale(y,1.0/h);CHKERRQ(ierr);
405e884886fSBarry Smith 
406c35ec66cSMatthew G Knepley   if ((ctx->vshift != 0.0) || (ctx->vscale != 1.0)) {
407e884886fSBarry Smith     ierr = VecAXPBY(y,ctx->vshift,ctx->vscale,a);CHKERRQ(ierr);
408c3bb7e23SBarry Smith   }
409c3bb7e23SBarry Smith   if (ctx->dlscale) {
410c3bb7e23SBarry Smith     ierr = VecPointwiseMult(y,ctx->dlscale,y);CHKERRQ(ierr);
411c3bb7e23SBarry Smith   }
412c3bb7e23SBarry Smith   if (ctx->dshift) {
413c3bb7e23SBarry Smith     ierr = VecPointwiseMult(ctx->dshift,a,U);CHKERRQ(ierr);
414c3bb7e23SBarry Smith     ierr = VecAXPY(y,1.0,U);CHKERRQ(ierr);
415c3bb7e23SBarry Smith   }
416e884886fSBarry Smith 
417e884886fSBarry Smith   if (ctx->sp) {ierr = MatNullSpaceRemove(ctx->sp,y,PETSC_NULL);CHKERRQ(ierr);}
418e884886fSBarry Smith 
419e884886fSBarry Smith   ierr = PetscLogEventEnd(MATMFFD_Mult,a,y,0,0);CHKERRQ(ierr);
420e884886fSBarry Smith   PetscFunctionReturn(0);
421e884886fSBarry Smith }
422e884886fSBarry Smith 
423e884886fSBarry Smith #undef __FUNCT__
424e884886fSBarry Smith #define __FUNCT__ "MatGetDiagonal_MFFD"
425e884886fSBarry Smith /*
426e884886fSBarry Smith   MatGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix
427e884886fSBarry Smith 
428e884886fSBarry Smith         y ~= (F(u + ha) - F(u))/h,
429e884886fSBarry Smith   where F = nonlinear function, as set by SNESSetFunction()
430e884886fSBarry Smith         u = current iterate
431e884886fSBarry Smith         h = difference interval
432e884886fSBarry Smith */
433e884886fSBarry Smith PetscErrorCode MatGetDiagonal_MFFD(Mat mat,Vec a)
434e884886fSBarry Smith {
435e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)mat->data;
436e884886fSBarry Smith   PetscScalar    h,*aa,*ww,v;
437e884886fSBarry Smith   PetscReal      epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON;
438e884886fSBarry Smith   Vec            w,U;
439e884886fSBarry Smith   PetscErrorCode ierr;
440e884886fSBarry Smith   PetscInt       i,rstart,rend;
441e884886fSBarry Smith 
442e884886fSBarry Smith   PetscFunctionBegin;
443e7e72b3dSBarry Smith   if (!ctx->funci) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Requires calling MatMFFDSetFunctioni() first");
444e884886fSBarry Smith 
445e884886fSBarry Smith   w    = ctx->w;
446e884886fSBarry Smith   U    = ctx->current_u;
447e884886fSBarry Smith   ierr = (*ctx->func)(ctx->funcctx,U,a);CHKERRQ(ierr);
448e884886fSBarry Smith   ierr = (*ctx->funcisetbase)(ctx->funcctx,U);CHKERRQ(ierr);
449e884886fSBarry Smith   ierr = VecCopy(U,w);CHKERRQ(ierr);
450e884886fSBarry Smith 
451e884886fSBarry Smith   ierr = VecGetOwnershipRange(a,&rstart,&rend);CHKERRQ(ierr);
452e884886fSBarry Smith   ierr = VecGetArray(a,&aa);CHKERRQ(ierr);
453e884886fSBarry Smith   for (i=rstart; i<rend; i++) {
454e884886fSBarry Smith     ierr = VecGetArray(w,&ww);CHKERRQ(ierr);
455e884886fSBarry Smith     h  = ww[i-rstart];
456e884886fSBarry Smith     if (h == 0.0) h = 1.0;
457e884886fSBarry Smith #if !defined(PETSC_USE_COMPLEX)
458e884886fSBarry Smith     if (h < umin && h >= 0.0)      h = umin;
459e884886fSBarry Smith     else if (h < 0.0 && h > -umin) h = -umin;
460e884886fSBarry Smith #else
461e884886fSBarry Smith     if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0)     h = umin;
462e884886fSBarry Smith     else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin;
463e884886fSBarry Smith #endif
464e884886fSBarry Smith     h     *= epsilon;
465e884886fSBarry Smith 
466e884886fSBarry Smith     ww[i-rstart] += h;
467e884886fSBarry Smith     ierr = VecRestoreArray(w,&ww);CHKERRQ(ierr);
468e884886fSBarry Smith     ierr          = (*ctx->funci)(ctx->funcctx,i,w,&v);CHKERRQ(ierr);
469e884886fSBarry Smith     aa[i-rstart]  = (v - aa[i-rstart])/h;
470e884886fSBarry Smith 
471e884886fSBarry Smith     /* possibly shift and scale result */
472c35ec66cSMatthew G Knepley     if ((ctx->vshift != 0.0) || (ctx->vscale != 1.0)) {
473e884886fSBarry Smith       aa[i - rstart] = ctx->vshift + ctx->vscale*aa[i-rstart];
474c3bb7e23SBarry Smith     }
475e884886fSBarry Smith 
476e884886fSBarry Smith     ierr = VecGetArray(w,&ww);CHKERRQ(ierr);
477e884886fSBarry Smith     ww[i-rstart] -= h;
478e884886fSBarry Smith     ierr = VecRestoreArray(w,&ww);CHKERRQ(ierr);
479e884886fSBarry Smith   }
480e884886fSBarry Smith   ierr = VecRestoreArray(a,&aa);CHKERRQ(ierr);
481e884886fSBarry Smith   PetscFunctionReturn(0);
482e884886fSBarry Smith }
483e884886fSBarry Smith 
484e884886fSBarry Smith #undef __FUNCT__
485c3bb7e23SBarry Smith #define __FUNCT__ "MatDiagonalScale_MFFD"
486c3bb7e23SBarry Smith PetscErrorCode MatDiagonalScale_MFFD(Mat mat,Vec ll,Vec rr)
487c3bb7e23SBarry Smith {
488c3bb7e23SBarry Smith   MatMFFD        aij = (MatMFFD)mat->data;
489c3bb7e23SBarry Smith   PetscErrorCode ierr;
490c3bb7e23SBarry Smith 
491c3bb7e23SBarry Smith   PetscFunctionBegin;
492c3bb7e23SBarry Smith   if (ll && !aij->dlscale) {
493c3bb7e23SBarry Smith     ierr = VecDuplicate(ll,&aij->dlscale);CHKERRQ(ierr);
494c3bb7e23SBarry Smith   }
495c3bb7e23SBarry Smith   if (rr && !aij->drscale) {
496c3bb7e23SBarry Smith     ierr = VecDuplicate(rr,&aij->drscale);CHKERRQ(ierr);
497c3bb7e23SBarry Smith   }
498c3bb7e23SBarry Smith   if (ll) {
499c3bb7e23SBarry Smith     ierr = VecCopy(ll,aij->dlscale);CHKERRQ(ierr);
500c3bb7e23SBarry Smith   }
501c3bb7e23SBarry Smith   if (rr) {
502c3bb7e23SBarry Smith     ierr = VecCopy(rr,aij->drscale);CHKERRQ(ierr);
503c3bb7e23SBarry Smith   }
504c3bb7e23SBarry Smith   PetscFunctionReturn(0);
505c3bb7e23SBarry Smith }
506c3bb7e23SBarry Smith 
507c3bb7e23SBarry Smith #undef __FUNCT__
508c3bb7e23SBarry Smith #define __FUNCT__ "MatDiagonalSet_MFFD"
509c3bb7e23SBarry Smith PetscErrorCode MatDiagonalSet_MFFD(Mat mat,Vec ll,InsertMode mode)
510c3bb7e23SBarry Smith {
511c3bb7e23SBarry Smith   MatMFFD        aij = (MatMFFD)mat->data;
512c3bb7e23SBarry Smith   PetscErrorCode ierr;
513c3bb7e23SBarry Smith 
514c3bb7e23SBarry Smith   PetscFunctionBegin;
515c3bb7e23SBarry Smith   if (mode == INSERT_VALUES) SETERRQ(((PetscObject)mat)->comm,PETSC_ERR_SUP,"No diagonal set with INSERT_VALUES");
516c3bb7e23SBarry Smith   if (!aij->dshift) {
517c3bb7e23SBarry Smith     ierr = VecDuplicate(ll,&aij->dshift);CHKERRQ(ierr);
518c3bb7e23SBarry Smith   }
519c3bb7e23SBarry Smith   ierr = VecAXPY(aij->dshift,1.0,ll);CHKERRQ(ierr);
520c3bb7e23SBarry Smith   PetscFunctionReturn(0);
521c3bb7e23SBarry Smith }
522c3bb7e23SBarry Smith 
523c3bb7e23SBarry Smith #undef __FUNCT__
524e884886fSBarry Smith #define __FUNCT__ "MatShift_MFFD"
525e884886fSBarry Smith PetscErrorCode MatShift_MFFD(Mat Y,PetscScalar a)
526e884886fSBarry Smith {
527e884886fSBarry Smith   MatMFFD shell = (MatMFFD)Y->data;
528e884886fSBarry Smith   PetscFunctionBegin;
529e884886fSBarry Smith   shell->vshift += a;
530e884886fSBarry Smith   PetscFunctionReturn(0);
531e884886fSBarry Smith }
532e884886fSBarry Smith 
533e884886fSBarry Smith #undef __FUNCT__
534e884886fSBarry Smith #define __FUNCT__ "MatScale_MFFD"
535e884886fSBarry Smith PetscErrorCode MatScale_MFFD(Mat Y,PetscScalar a)
536e884886fSBarry Smith {
537e884886fSBarry Smith   MatMFFD shell = (MatMFFD)Y->data;
538e884886fSBarry Smith   PetscFunctionBegin;
539e884886fSBarry Smith   shell->vscale *= a;
540e884886fSBarry Smith   PetscFunctionReturn(0);
541e884886fSBarry Smith }
542e884886fSBarry Smith 
543e884886fSBarry Smith EXTERN_C_BEGIN
544e884886fSBarry Smith #undef __FUNCT__
545c879e56bSBarry Smith #define __FUNCT__ "MatMFFDSetBase_MFFD"
5467087cfbeSBarry Smith PetscErrorCode  MatMFFDSetBase_MFFD(Mat J,Vec U,Vec F)
547e884886fSBarry Smith {
548e884886fSBarry Smith   PetscErrorCode ierr;
549e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)J->data;
550e884886fSBarry Smith 
551e884886fSBarry Smith   PetscFunctionBegin;
552e884886fSBarry Smith   ierr = MatMFFDResetHHistory(J);CHKERRQ(ierr);
553e884886fSBarry Smith   ctx->current_u = U;
55452121784SBarry Smith   if (F) {
5556bf464f9SBarry Smith     if (ctx->current_f_allocated) {ierr = VecDestroy(&ctx->current_f);CHKERRQ(ierr);}
5563ec795f1SBarry Smith     ctx->current_f           = F;
557cfe22f5eSBarry Smith     ctx->current_f_allocated = PETSC_FALSE;
558cfe22f5eSBarry Smith   } else if (!ctx->current_f_allocated) {
55952121784SBarry Smith     ierr = VecDuplicate(ctx->current_u, &ctx->current_f);CHKERRQ(ierr);
560cfe22f5eSBarry Smith     ctx->current_f_allocated = PETSC_TRUE;
56152121784SBarry Smith   }
562e884886fSBarry Smith   if (!ctx->w) {
563e884886fSBarry Smith     ierr = VecDuplicate(ctx->current_u, &ctx->w);CHKERRQ(ierr);
564e884886fSBarry Smith   }
565e884886fSBarry Smith   J->assembled = PETSC_TRUE;
566e884886fSBarry Smith   PetscFunctionReturn(0);
567e884886fSBarry Smith }
568e884886fSBarry Smith EXTERN_C_END
569e884886fSBarry Smith typedef PetscErrorCode (*FCN3)(void*,Vec,Vec,PetscScalar*); /* force argument to next function to not be extern C*/
570e884886fSBarry Smith EXTERN_C_BEGIN
571e884886fSBarry Smith #undef __FUNCT__
572c879e56bSBarry Smith #define __FUNCT__ "MatMFFDSetCheckh_MFFD"
5737087cfbeSBarry Smith PetscErrorCode  MatMFFDSetCheckh_MFFD(Mat J,FCN3 fun,void*ectx)
574e884886fSBarry Smith {
575e884886fSBarry Smith   MatMFFD ctx = (MatMFFD)J->data;
576e884886fSBarry Smith 
577e884886fSBarry Smith   PetscFunctionBegin;
578e884886fSBarry Smith   ctx->checkh    = fun;
579e884886fSBarry Smith   ctx->checkhctx = ectx;
580e884886fSBarry Smith   PetscFunctionReturn(0);
581e884886fSBarry Smith }
582e884886fSBarry Smith EXTERN_C_END
583e884886fSBarry Smith 
584e884886fSBarry Smith #undef __FUNCT__
5856aa9148fSLisandro Dalcin #define __FUNCT__ "MatMFFDSetOptionsPrefix"
5866aa9148fSLisandro Dalcin /*@C
5876aa9148fSLisandro Dalcin    MatMFFDSetOptionsPrefix - Sets the prefix used for searching for all
5886aa9148fSLisandro Dalcin    MatMFFD options in the database.
5896aa9148fSLisandro Dalcin 
5906aa9148fSLisandro Dalcin    Collective on Mat
5916aa9148fSLisandro Dalcin 
5926aa9148fSLisandro Dalcin    Input Parameter:
5936aa9148fSLisandro Dalcin +  A - the Mat context
5946aa9148fSLisandro Dalcin -  prefix - the prefix to prepend to all option names
5956aa9148fSLisandro Dalcin 
5966aa9148fSLisandro Dalcin    Notes:
5976aa9148fSLisandro Dalcin    A hyphen (-) must NOT be given at the beginning of the prefix name.
5986aa9148fSLisandro Dalcin    The first character of all runtime options is AUTOMATICALLY the hyphen.
5996aa9148fSLisandro Dalcin 
6006aa9148fSLisandro Dalcin    Level: advanced
6016aa9148fSLisandro Dalcin 
6026aa9148fSLisandro Dalcin .keywords: SNES, matrix-free, parameters
6036aa9148fSLisandro Dalcin 
6049c6ac3b3SBarry Smith .seealso: MatSetFromOptions(), MatCreateSNESMF()
6056aa9148fSLisandro Dalcin @*/
6067087cfbeSBarry Smith PetscErrorCode  MatMFFDSetOptionsPrefix(Mat mat,const char prefix[])
6076aa9148fSLisandro Dalcin 
6086aa9148fSLisandro Dalcin {
6091db2f0e5SJed Brown   MatMFFD        mfctx = mat ? (MatMFFD)mat->data : (MatMFFD)PETSC_NULL;
6106aa9148fSLisandro Dalcin   PetscErrorCode ierr;
6116aa9148fSLisandro Dalcin   PetscFunctionBegin;
6120700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
6130700a824SBarry Smith   PetscValidHeaderSpecific(mfctx,MATMFFD_CLASSID,1);
6146aa9148fSLisandro Dalcin   ierr = PetscObjectSetOptionsPrefix((PetscObject)mfctx,prefix);CHKERRQ(ierr);
6156aa9148fSLisandro Dalcin   PetscFunctionReturn(0);
6166aa9148fSLisandro Dalcin }
6176aa9148fSLisandro Dalcin 
6186aa9148fSLisandro Dalcin #undef __FUNCT__
6199c6ac3b3SBarry Smith #define __FUNCT__ "MatSetFromOptions_MFFD"
6209c6ac3b3SBarry Smith PetscErrorCode  MatSetFromOptions_MFFD(Mat mat)
621e884886fSBarry Smith {
62271f08665SBarry Smith   MatMFFD        mfctx = (MatMFFD)mat->data;
623e884886fSBarry Smith   PetscErrorCode ierr;
624ace3abfcSBarry Smith   PetscBool      flg;
625e884886fSBarry Smith   char           ftype[256];
626e884886fSBarry Smith 
627e884886fSBarry Smith   PetscFunctionBegin;
6280700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
6290700a824SBarry Smith   PetscValidHeaderSpecific(mfctx,MATMFFD_CLASSID,1);
6303194b578SJed Brown   ierr = PetscObjectOptionsBegin((PetscObject)mfctx);CHKERRQ(ierr);
631b022a5c1SBarry Smith   ierr = PetscOptionsList("-mat_mffd_type","Matrix free type","MatMFFDSetType",MatMFFDList,((PetscObject)mfctx)->type_name,ftype,256,&flg);CHKERRQ(ierr);
632e884886fSBarry Smith   if (flg) {
633e884886fSBarry Smith     ierr = MatMFFDSetType(mat,ftype);CHKERRQ(ierr);
634e884886fSBarry Smith   }
635e884886fSBarry Smith 
636e884886fSBarry Smith   ierr = PetscOptionsReal("-mat_mffd_err","set sqrt relative error in function","MatMFFDSetFunctionError",mfctx->error_rel,&mfctx->error_rel,0);CHKERRQ(ierr);
637e884886fSBarry Smith   ierr = PetscOptionsInt("-mat_mffd_period","how often h is recomputed","MatMFFDSetPeriod",mfctx->recomputeperiod,&mfctx->recomputeperiod,0);CHKERRQ(ierr);
638e884886fSBarry Smith 
63990d69ab7SBarry Smith   flg  = PETSC_FALSE;
640acfcf0e5SJed Brown   ierr = PetscOptionsBool("-mat_mffd_check_positivity","Insure that U + h*a is nonnegative","MatMFFDSetCheckh",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
641e884886fSBarry Smith   if (flg) {
642e884886fSBarry Smith     ierr = MatMFFDSetCheckh(mat,MatMFFDCheckPositivity,0);CHKERRQ(ierr);
643e884886fSBarry Smith   }
644e884886fSBarry Smith   if (mfctx->ops->setfromoptions) {
645e884886fSBarry Smith     ierr = (*mfctx->ops->setfromoptions)(mfctx);CHKERRQ(ierr);
646e884886fSBarry Smith   }
647e884886fSBarry Smith   ierr = PetscOptionsEnd();CHKERRQ(ierr);
648e884886fSBarry Smith   PetscFunctionReturn(0);
649e884886fSBarry Smith }
650e884886fSBarry Smith 
651bc0f08ceSBarry Smith EXTERN_C_BEGIN
652bc0f08ceSBarry Smith #undef __FUNCT__
653bc0f08ceSBarry Smith #define __FUNCT__ "MatMFFDSetPeriod_MFFD"
654bc0f08ceSBarry Smith PetscErrorCode  MatMFFDSetPeriod_MFFD(Mat mat,PetscInt period)
655bc0f08ceSBarry Smith {
656bc0f08ceSBarry Smith   MatMFFD ctx = (MatMFFD)mat->data;
657bc0f08ceSBarry Smith 
658bc0f08ceSBarry Smith   PetscFunctionBegin;
659bc0f08ceSBarry Smith   PetscValidLogicalCollectiveInt(mat,period,2);
660bc0f08ceSBarry Smith   ctx->recomputeperiod = period;
661bc0f08ceSBarry Smith   PetscFunctionReturn(0);
662bc0f08ceSBarry Smith }
663bc0f08ceSBarry Smith EXTERN_C_END
664bc0f08ceSBarry Smith 
665bc0f08ceSBarry Smith EXTERN_C_BEGIN
666bc0f08ceSBarry Smith #undef __FUNCT__
667bc0f08ceSBarry Smith #define __FUNCT__ "MatMFFDSetFunction_MFFD"
668bc0f08ceSBarry Smith PetscErrorCode  MatMFFDSetFunction_MFFD(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx)
669bc0f08ceSBarry Smith {
670bc0f08ceSBarry Smith   MatMFFD ctx = (MatMFFD)mat->data;
671bc0f08ceSBarry Smith 
672bc0f08ceSBarry Smith   PetscFunctionBegin;
673bc0f08ceSBarry Smith   ctx->func    = func;
674bc0f08ceSBarry Smith   ctx->funcctx = funcctx;
675bc0f08ceSBarry Smith   PetscFunctionReturn(0);
676bc0f08ceSBarry Smith }
677bc0f08ceSBarry Smith EXTERN_C_END
678bc0f08ceSBarry Smith 
679bc0f08ceSBarry Smith EXTERN_C_BEGIN
680bc0f08ceSBarry Smith #undef __FUNCT__
681bc0f08ceSBarry Smith #define __FUNCT__ "MatMFFDSetFunctionError_MFFD"
682bc0f08ceSBarry Smith PetscErrorCode  MatMFFDSetFunctionError_MFFD(Mat mat,PetscReal error)
683bc0f08ceSBarry Smith {
684bc0f08ceSBarry Smith   MatMFFD ctx = (MatMFFD)mat->data;
685bc0f08ceSBarry Smith 
686bc0f08ceSBarry Smith   PetscFunctionBegin;
687bc0f08ceSBarry Smith   PetscValidLogicalCollectiveReal(mat,error,2);
688bc0f08ceSBarry Smith   if (error != PETSC_DEFAULT) ctx->error_rel = error;
689bc0f08ceSBarry Smith   PetscFunctionReturn(0);
690bc0f08ceSBarry Smith }
691bc0f08ceSBarry Smith EXTERN_C_END
692bc0f08ceSBarry Smith 
693e884886fSBarry Smith /*MC
694e884886fSBarry Smith   MATMFFD - MATMFFD = "mffd" - A matrix free matrix type.
695e884886fSBarry Smith 
696e884886fSBarry Smith   Level: advanced
697e884886fSBarry Smith 
698d2d6cebeSBarry Smith .seealso: MatCreateMFFD(), MatCreateSNESMF(), MatMFFDSetFunction()
699e884886fSBarry Smith M*/
700e884886fSBarry Smith EXTERN_C_BEGIN
701e884886fSBarry Smith #undef __FUNCT__
702e884886fSBarry Smith #define __FUNCT__ "MatCreate_MFFD"
7037087cfbeSBarry Smith PetscErrorCode  MatCreate_MFFD(Mat A)
704e884886fSBarry Smith {
705e884886fSBarry Smith   MatMFFD         mfctx;
706e884886fSBarry Smith   PetscErrorCode  ierr;
707e884886fSBarry Smith 
708e884886fSBarry Smith   PetscFunctionBegin;
7093ec795f1SBarry Smith #ifndef PETSC_USE_DYNAMIC_LIBRARIES
7103ec795f1SBarry Smith   ierr = MatMFFDInitializePackage(PETSC_NULL);CHKERRQ(ierr);
7113ec795f1SBarry Smith #endif
712e884886fSBarry Smith 
7133194b578SJed Brown   ierr = PetscHeaderCreate(mfctx,_p_MatMFFD,struct _MFOps,MATMFFD_CLASSID,0,"MatMFFD","Matrix-free Finite Differencing","Mat",((PetscObject)A)->comm,MatDestroy_MFFD,MatView_MFFD);CHKERRQ(ierr);
714e884886fSBarry Smith   mfctx->sp              = 0;
715e884886fSBarry Smith   mfctx->error_rel       = PETSC_SQRT_MACHINE_EPSILON;
716e884886fSBarry Smith   mfctx->recomputeperiod = 1;
717e884886fSBarry Smith   mfctx->count           = 0;
718e884886fSBarry Smith   mfctx->currenth        = 0.0;
719e884886fSBarry Smith   mfctx->historyh        = PETSC_NULL;
720e884886fSBarry Smith   mfctx->ncurrenth       = 0;
721e884886fSBarry Smith   mfctx->maxcurrenth     = 0;
7227adad957SLisandro Dalcin   ((PetscObject)mfctx)->type_name       = 0;
723e884886fSBarry Smith 
724e884886fSBarry Smith   mfctx->vshift          = 0.0;
725e884886fSBarry Smith   mfctx->vscale          = 1.0;
726e884886fSBarry Smith 
727e884886fSBarry Smith   /*
728e884886fSBarry Smith      Create the empty data structure to contain compute-h routines.
729e884886fSBarry Smith      These will be filled in below from the command line options or
730e884886fSBarry Smith      a later call with MatMFFDSetType() or if that is not called
731e884886fSBarry Smith      then it will default in the first use of MatMult_MFFD()
732e884886fSBarry Smith   */
733e884886fSBarry Smith   mfctx->ops->compute        = 0;
734e884886fSBarry Smith   mfctx->ops->destroy        = 0;
735e884886fSBarry Smith   mfctx->ops->view           = 0;
736e884886fSBarry Smith   mfctx->ops->setfromoptions = 0;
737e884886fSBarry Smith   mfctx->hctx                = 0;
738e884886fSBarry Smith 
739e884886fSBarry Smith   mfctx->func                = 0;
740e884886fSBarry Smith   mfctx->funcctx             = 0;
741e884886fSBarry Smith   mfctx->w                   = PETSC_NULL;
742e884886fSBarry Smith 
743e884886fSBarry Smith   A->data                = mfctx;
744e884886fSBarry Smith 
745e884886fSBarry Smith   A->ops->mult           = MatMult_MFFD;
746e884886fSBarry Smith   A->ops->destroy        = MatDestroy_MFFD;
747e884886fSBarry Smith   A->ops->view           = MatView_MFFD;
748e884886fSBarry Smith   A->ops->assemblyend    = MatAssemblyEnd_MFFD;
749e884886fSBarry Smith   A->ops->getdiagonal    = MatGetDiagonal_MFFD;
750e884886fSBarry Smith   A->ops->scale          = MatScale_MFFD;
751e884886fSBarry Smith   A->ops->shift          = MatShift_MFFD;
752c3bb7e23SBarry Smith   A->ops->diagonalscale  = MatDiagonalScale_MFFD;
753c3bb7e23SBarry Smith   A->ops->diagonalset    = MatDiagonalSet_MFFD;
7549c6ac3b3SBarry Smith   A->ops->setfromoptions = MatSetFromOptions_MFFD;
755e884886fSBarry Smith   A->assembled = PETSC_TRUE;
756e884886fSBarry Smith 
75726283091SBarry Smith   ierr = PetscLayoutSetBlockSize(A->rmap,1);CHKERRQ(ierr);
75826283091SBarry Smith   ierr = PetscLayoutSetBlockSize(A->cmap,1);CHKERRQ(ierr);
75926283091SBarry Smith   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
76026283091SBarry Smith   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
761ee1cef2cSJed Brown 
762c879e56bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDSetBase_C","MatMFFDSetBase_MFFD",MatMFFDSetBase_MFFD);CHKERRQ(ierr);
763c879e56bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDSetFunctioniBase_C","MatMFFDSetFunctioniBase_MFFD",MatMFFDSetFunctioniBase_MFFD);CHKERRQ(ierr);
764c879e56bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDSetFunctioni_C","MatMFFDSetFunctioni_MFFD",MatMFFDSetFunctioni_MFFD);CHKERRQ(ierr);
765bc0f08ceSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDSetFunction_C","MatMFFDSetFunction_MFFD",MatMFFDSetFunction_MFFD);CHKERRQ(ierr);
766c879e56bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDSetCheckh_C","MatMFFDSetCheckh_MFFD",MatMFFDSetCheckh_MFFD);CHKERRQ(ierr);
767bc0f08ceSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDSetPeriod_C","MatMFFDSetPeriod_MFFD",MatMFFDSetPeriod_MFFD);CHKERRQ(ierr);
768bc0f08ceSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDSetFunctionError_C","MatMFFDSetFunctionError_MFFD",MatMFFDSetFunctionError_MFFD);CHKERRQ(ierr);
769bc0f08ceSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDResetHHistory_C","MatMFFDResetHHistory_MFFD",MatMFFDResetHHistory_MFFD);CHKERRQ(ierr);
770bc0f08ceSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDAddNullSpace_C","MatMFFDAddNullSpace_MFFD",MatMFFDAddNullSpace_MFFD);CHKERRQ(ierr);
771e884886fSBarry Smith   mfctx->mat = A;
772e884886fSBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,MATMFFD);CHKERRQ(ierr);
773e884886fSBarry Smith   PetscFunctionReturn(0);
774e884886fSBarry Smith }
775e884886fSBarry Smith EXTERN_C_END
776e884886fSBarry Smith 
777e884886fSBarry Smith #undef __FUNCT__
778e884886fSBarry Smith #define __FUNCT__ "MatCreateMFFD"
779e884886fSBarry Smith /*@
780e884886fSBarry Smith    MatCreateMFFD - Creates a matrix-free matrix. See also MatCreateSNESMF()
781e884886fSBarry Smith 
782e884886fSBarry Smith    Collective on Vec
783e884886fSBarry Smith 
784e884886fSBarry Smith    Input Parameters:
785fef1beadSBarry Smith +  comm - MPI communicator
786fef1beadSBarry Smith .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
787fef1beadSBarry Smith            This value should be the same as the local size used in creating the
788fef1beadSBarry Smith            y vector for the matrix-vector product y = Ax.
789fef1beadSBarry Smith .  n - This value should be the same as the local size used in creating the
790fef1beadSBarry Smith        x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have
791fef1beadSBarry Smith        calculated if N is given) For square matrices n is almost always m.
792fef1beadSBarry Smith .  M - number of global rows (or PETSC_DETERMINE to have calculated if m is given)
793fef1beadSBarry Smith -  N - number of global columns (or PETSC_DETERMINE to have calculated if n is given)
794fef1beadSBarry Smith 
795e884886fSBarry Smith 
796e884886fSBarry Smith    Output Parameter:
797e884886fSBarry Smith .  J - the matrix-free matrix
798e884886fSBarry Smith 
7999c6ac3b3SBarry Smith    Options Database Keys: call MatSetFromOptions() to trigger these
8009c6ac3b3SBarry Smith +  -mat_mffd_type - wp or ds (see MATMFFD_WP or MATMFFD_DS)
8019c6ac3b3SBarry Smith -  -mat_mffd_err - square root of estimated relative error in function evaluation
8029c6ac3b3SBarry Smith -  -mat_mffd_period - how often h is recomputed, defaults to 1, everytime
8039c6ac3b3SBarry Smith 
8049c6ac3b3SBarry Smith 
805e884886fSBarry Smith    Level: advanced
806e884886fSBarry Smith 
807e884886fSBarry Smith    Notes:
808e884886fSBarry Smith    The matrix-free matrix context merely contains the function pointers
809e884886fSBarry Smith    and work space for performing finite difference approximations of
810e884886fSBarry Smith    Jacobian-vector products, F'(u)*a,
811e884886fSBarry Smith 
812e884886fSBarry Smith    The default code uses the following approach to compute h
813e884886fSBarry Smith 
814e884886fSBarry Smith .vb
815e884886fSBarry Smith      F'(u)*a = [F(u+h*a) - F(u)]/h where
816e884886fSBarry Smith      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
817e884886fSBarry Smith        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
818e884886fSBarry Smith  where
819e884886fSBarry Smith      error_rel = square root of relative error in function evaluation
820e884886fSBarry Smith      umin = minimum iterate parameter
821e884886fSBarry Smith .ve
822e884886fSBarry Smith 
823e884886fSBarry Smith    The user can set the error_rel via MatMFFDSetFunctionError() and
8240598bfebSBarry Smith    umin via MatMFFDDSSetUmin(); see the <A href="../../docs/manual.pdf#nameddest=ch_snes">SNES chapter of the users manual</A> for details.
825e884886fSBarry Smith 
826e884886fSBarry Smith    The user should call MatDestroy() when finished with the matrix-free
827e884886fSBarry Smith    matrix context.
828e884886fSBarry Smith 
829e884886fSBarry Smith    Options Database Keys:
830e884886fSBarry Smith +  -mat_mffd_err <error_rel> - Sets error_rel
831e884886fSBarry Smith .  -mat_mffd_unim <umin> - Sets umin (for default PETSc routine that computes h only)
832e884886fSBarry Smith -  -mat_mffd_check_positivity
833e884886fSBarry Smith 
834e884886fSBarry Smith .keywords: default, matrix-free, create, matrix
835e884886fSBarry Smith 
8361d0fab5eSBarry Smith .seealso: MatDestroy(), MatMFFDSetFunctionError(), MatMFFDDSSetUmin(), MatMFFDSetFunction()
837e884886fSBarry Smith           MatMFFDSetHHistory(), MatMFFDResetHHistory(), MatCreateSNESMF(),
8381d0fab5eSBarry Smith           MatMFFDGetH(), MatMFFDRegisterDynamic), MatMFFDComputeJacobian()
839e884886fSBarry Smith 
840e884886fSBarry Smith @*/
8417087cfbeSBarry Smith PetscErrorCode  MatCreateMFFD(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *J)
842e884886fSBarry Smith {
843e884886fSBarry Smith   PetscErrorCode ierr;
844e884886fSBarry Smith 
845e884886fSBarry Smith   PetscFunctionBegin;
846e884886fSBarry Smith   ierr = MatCreate(comm,J);CHKERRQ(ierr);
847fef1beadSBarry Smith   ierr = MatSetSizes(*J,m,n,M,N);CHKERRQ(ierr);
848e884886fSBarry Smith   ierr = MatSetType(*J,MATMFFD);CHKERRQ(ierr);
849be7c243fSJed Brown   ierr = MatSetUp(*J);CHKERRQ(ierr);
850e884886fSBarry Smith   PetscFunctionReturn(0);
851e884886fSBarry Smith }
852e884886fSBarry Smith 
853e884886fSBarry Smith 
854e884886fSBarry Smith #undef __FUNCT__
855e884886fSBarry Smith #define __FUNCT__ "MatMFFDGetH"
856e884886fSBarry Smith /*@
857e884886fSBarry Smith    MatMFFDGetH - Gets the last value that was used as the differencing
858e884886fSBarry Smith    parameter.
859e884886fSBarry Smith 
860e884886fSBarry Smith    Not Collective
861e884886fSBarry Smith 
862e884886fSBarry Smith    Input Parameters:
863e884886fSBarry Smith .  mat - the matrix obtained with MatCreateSNESMF()
864e884886fSBarry Smith 
865e884886fSBarry Smith    Output Paramter:
866e884886fSBarry Smith .  h - the differencing step size
867e884886fSBarry Smith 
868e884886fSBarry Smith    Level: advanced
869e884886fSBarry Smith 
870e884886fSBarry Smith .keywords: SNES, matrix-free, parameters
871e884886fSBarry Smith 
8721d0fab5eSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDSetHHistory(), MatCreateMFFD(), MATMFFD, MatMFFDResetHHistory()
873e884886fSBarry Smith @*/
8747087cfbeSBarry Smith PetscErrorCode  MatMFFDGetH(Mat mat,PetscScalar *h)
875e884886fSBarry Smith {
876e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)mat->data;
877bc0f08ceSBarry Smith   PetscErrorCode ierr;
878bc0f08ceSBarry Smith   PetscBool      match;
879e884886fSBarry Smith 
880e884886fSBarry Smith   PetscFunctionBegin;
881bc0f08ceSBarry Smith   ierr = PetscTypeCompare((PetscObject)mat,MATMFFD,&match);CHKERRQ(ierr);
882bc0f08ceSBarry Smith   if (!match) SETERRQ(((PetscObject)mat)->comm,PETSC_ERR_ARG_WRONG,"Not a MFFD matrix");
883bc0f08ceSBarry Smith 
884e884886fSBarry Smith   *h = ctx->currenth;
885e884886fSBarry Smith   PetscFunctionReturn(0);
886e884886fSBarry Smith }
887e884886fSBarry Smith 
888e884886fSBarry Smith #undef __FUNCT__
889e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetFunction"
890e884886fSBarry Smith /*@C
891e884886fSBarry Smith    MatMFFDSetFunction - Sets the function used in applying the matrix free.
892e884886fSBarry Smith 
8933f9fe445SBarry Smith    Logically Collective on Mat
894e884886fSBarry Smith 
895e884886fSBarry Smith    Input Parameters:
896e884886fSBarry Smith +  mat - the matrix free matrix created via MatCreateSNESMF()
897e884886fSBarry Smith .  func - the function to use
898e884886fSBarry Smith -  funcctx - optional function context passed to function
899e884886fSBarry Smith 
900e884886fSBarry Smith    Level: advanced
901e884886fSBarry Smith 
902e884886fSBarry Smith    Notes:
903e884886fSBarry Smith     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
904e884886fSBarry Smith     matrix inside your compute Jacobian routine
905e884886fSBarry Smith 
9063ec795f1SBarry Smith     If this is not set then it will use the function set with SNESSetFunction() if MatCreateSNESMF() was used.
907e884886fSBarry Smith 
908e884886fSBarry Smith .keywords: SNES, matrix-free, function
909e884886fSBarry Smith 
9101d0fab5eSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD,
9111d0fab5eSBarry Smith           MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction()
912e884886fSBarry Smith @*/
9137087cfbeSBarry Smith PetscErrorCode  MatMFFDSetFunction(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx)
914e884886fSBarry Smith {
915bc0f08ceSBarry Smith   PetscErrorCode ierr;
916e884886fSBarry Smith 
917e884886fSBarry Smith   PetscFunctionBegin;
918bc0f08ceSBarry Smith   ierr = PetscTryMethod(mat,"MatMFFDSetFunction_C",(Mat,PetscErrorCode (*)(void*,Vec,Vec),void*),(mat,func,funcctx));CHKERRQ(ierr);
919e884886fSBarry Smith   PetscFunctionReturn(0);
920e884886fSBarry Smith }
921e884886fSBarry Smith 
922e884886fSBarry Smith #undef __FUNCT__
923e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetFunctioni"
924e884886fSBarry Smith /*@C
925e884886fSBarry Smith    MatMFFDSetFunctioni - Sets the function for a single component
926e884886fSBarry Smith 
9273f9fe445SBarry Smith    Logically Collective on Mat
928e884886fSBarry Smith 
929e884886fSBarry Smith    Input Parameters:
930e884886fSBarry Smith +  mat - the matrix free matrix created via MatCreateSNESMF()
931e884886fSBarry Smith -  funci - the function to use
932e884886fSBarry Smith 
933e884886fSBarry Smith    Level: advanced
934e884886fSBarry Smith 
935e884886fSBarry Smith    Notes:
936e884886fSBarry Smith     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
937e884886fSBarry Smith     matrix inside your compute Jacobian routine
938e884886fSBarry Smith 
939e884886fSBarry Smith 
940e884886fSBarry Smith .keywords: SNES, matrix-free, function
941e884886fSBarry Smith 
9421d0fab5eSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction()
9431d0fab5eSBarry Smith 
944e884886fSBarry Smith @*/
9457087cfbeSBarry Smith PetscErrorCode  MatMFFDSetFunctioni(Mat mat,PetscErrorCode (*funci)(void*,PetscInt,Vec,PetscScalar*))
946e884886fSBarry Smith {
9474ac538c5SBarry Smith   PetscErrorCode ierr;
948e884886fSBarry Smith 
949e884886fSBarry Smith   PetscFunctionBegin;
9500700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
9514ac538c5SBarry Smith   ierr = PetscTryMethod(mat,"MatMFFDSetFunctioni_C",(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)),(mat,funci));CHKERRQ(ierr);
952e884886fSBarry Smith   PetscFunctionReturn(0);
953e884886fSBarry Smith }
954e884886fSBarry Smith 
955e884886fSBarry Smith 
956e884886fSBarry Smith #undef __FUNCT__
957e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetFunctioniBase"
958e884886fSBarry Smith /*@C
959e884886fSBarry Smith    MatMFFDSetFunctioniBase - Sets the base vector for a single component function evaluation
960e884886fSBarry Smith 
9613f9fe445SBarry Smith    Logically Collective on Mat
962e884886fSBarry Smith 
963e884886fSBarry Smith    Input Parameters:
964e884886fSBarry Smith +  mat - the matrix free matrix created via MatCreateSNESMF()
965e884886fSBarry Smith -  func - the function to use
966e884886fSBarry Smith 
967e884886fSBarry Smith    Level: advanced
968e884886fSBarry Smith 
969e884886fSBarry Smith    Notes:
970e884886fSBarry Smith     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
971e884886fSBarry Smith     matrix inside your compute Jacobian routine
972e884886fSBarry Smith 
973e884886fSBarry Smith 
974e884886fSBarry Smith .keywords: SNES, matrix-free, function
975e884886fSBarry Smith 
976e884886fSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD
9771d0fab5eSBarry Smith           MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction()
978e884886fSBarry Smith @*/
9797087cfbeSBarry Smith PetscErrorCode  MatMFFDSetFunctioniBase(Mat mat,PetscErrorCode (*func)(void*,Vec))
980e884886fSBarry Smith {
9814ac538c5SBarry Smith   PetscErrorCode ierr;
982e884886fSBarry Smith 
983e884886fSBarry Smith   PetscFunctionBegin;
9840700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
9854ac538c5SBarry Smith   ierr = PetscTryMethod(mat,"MatMFFDSetFunctioniBase_C",(Mat,PetscErrorCode (*)(void*,Vec)),(mat,func));CHKERRQ(ierr);
986e884886fSBarry Smith   PetscFunctionReturn(0);
987e884886fSBarry Smith }
988e884886fSBarry Smith 
989e884886fSBarry Smith #undef __FUNCT__
990e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetPeriod"
991e884886fSBarry Smith /*@
992e884886fSBarry Smith    MatMFFDSetPeriod - Sets how often h is recomputed, by default it is everytime
993e884886fSBarry Smith 
9943f9fe445SBarry Smith    Logically Collective on Mat
995e884886fSBarry Smith 
996e884886fSBarry Smith    Input Parameters:
997e884886fSBarry Smith +  mat - the matrix free matrix created via MatCreateSNESMF()
998e884886fSBarry Smith -  period - 1 for everytime, 2 for every second etc
999e884886fSBarry Smith 
1000e884886fSBarry Smith    Options Database Keys:
1001e884886fSBarry Smith +  -mat_mffd_period <period>
1002e884886fSBarry Smith 
1003e884886fSBarry Smith    Level: advanced
1004e884886fSBarry Smith 
1005e884886fSBarry Smith 
1006e884886fSBarry Smith .keywords: SNES, matrix-free, parameters
1007e884886fSBarry Smith 
1008e884886fSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(),
10091d0fab5eSBarry Smith           MatMFFDSetHHistory(), MatMFFDResetHHistory()
1010e884886fSBarry Smith @*/
10117087cfbeSBarry Smith PetscErrorCode  MatMFFDSetPeriod(Mat mat,PetscInt period)
1012e884886fSBarry Smith {
1013bc0f08ceSBarry Smith   PetscErrorCode ierr;
1014e884886fSBarry Smith 
1015e884886fSBarry Smith   PetscFunctionBegin;
1016bc0f08ceSBarry Smith   ierr = PetscTryMethod(mat,"MatMFFDSetPeriod_C",(Mat,PetscInt),(mat,period));CHKERRQ(ierr);
1017e884886fSBarry Smith   PetscFunctionReturn(0);
1018e884886fSBarry Smith }
1019e884886fSBarry Smith 
1020e884886fSBarry Smith #undef __FUNCT__
1021e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetFunctionError"
1022e884886fSBarry Smith /*@
1023e884886fSBarry Smith    MatMFFDSetFunctionError - Sets the error_rel for the approximation of
1024e884886fSBarry Smith    matrix-vector products using finite differences.
1025e884886fSBarry Smith 
10263f9fe445SBarry Smith    Logically Collective on Mat
1027e884886fSBarry Smith 
1028e884886fSBarry Smith    Input Parameters:
1029e884886fSBarry Smith +  mat - the matrix free matrix created via MatCreateMFFD() or MatCreateSNESMF()
1030e884886fSBarry Smith -  error_rel - relative error (should be set to the square root of
1031e884886fSBarry Smith                the relative error in the function evaluations)
1032e884886fSBarry Smith 
1033e884886fSBarry Smith    Options Database Keys:
1034e884886fSBarry Smith +  -mat_mffd_err <error_rel> - Sets error_rel
1035e884886fSBarry Smith 
1036e884886fSBarry Smith    Level: advanced
1037e884886fSBarry Smith 
1038e884886fSBarry Smith    Notes:
1039e884886fSBarry Smith    The default matrix-free matrix-vector product routine computes
1040e884886fSBarry Smith .vb
1041e884886fSBarry Smith      F'(u)*a = [F(u+h*a) - F(u)]/h where
1042e884886fSBarry Smith      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
1043e884886fSBarry Smith        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   else
1044e884886fSBarry Smith .ve
1045e884886fSBarry Smith 
1046e884886fSBarry Smith .keywords: SNES, matrix-free, parameters
1047e884886fSBarry Smith 
1048e884886fSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD
10491d0fab5eSBarry Smith           MatMFFDSetHHistory(), MatMFFDResetHHistory()
1050e884886fSBarry Smith @*/
10517087cfbeSBarry Smith PetscErrorCode  MatMFFDSetFunctionError(Mat mat,PetscReal error)
1052e884886fSBarry Smith {
1053bc0f08ceSBarry Smith   PetscErrorCode ierr;
1054e884886fSBarry Smith 
1055e884886fSBarry Smith   PetscFunctionBegin;
1056bc0f08ceSBarry Smith   ierr = PetscTryMethod(mat,"MatMFFDSetFunctionError_C",(Mat,PetscReal),(mat,error));CHKERRQ(ierr);
1057e884886fSBarry Smith   PetscFunctionReturn(0);
1058e884886fSBarry Smith }
1059e884886fSBarry Smith 
1060e884886fSBarry Smith #undef __FUNCT__
1061e884886fSBarry Smith #define __FUNCT__ "MatMFFDAddNullSpace"
1062e884886fSBarry Smith /*@
1063e884886fSBarry Smith    MatMFFDAddNullSpace - Provides a null space that an operator is
1064e884886fSBarry Smith    supposed to have.  Since roundoff will create a small component in
1065e884886fSBarry Smith    the null space, if you know the null space you may have it
1066e884886fSBarry Smith    automatically removed.
1067e884886fSBarry Smith 
10683f9fe445SBarry Smith    Logically Collective on Mat
1069e884886fSBarry Smith 
1070e884886fSBarry Smith    Input Parameters:
1071e884886fSBarry Smith +  J - the matrix-free matrix context
1072e884886fSBarry Smith -  nullsp - object created with MatNullSpaceCreate()
1073e884886fSBarry Smith 
1074e884886fSBarry Smith    Level: advanced
1075e884886fSBarry Smith 
1076e884886fSBarry Smith .keywords: SNES, matrix-free, null space
1077e884886fSBarry Smith 
1078e884886fSBarry Smith .seealso: MatNullSpaceCreate(), MatMFFDGetH(), MatCreateSNESMF(), MatCreateMFFD(), MATMFFD
10791d0fab5eSBarry Smith           MatMFFDSetHHistory(), MatMFFDResetHHistory()
1080e884886fSBarry Smith @*/
10817087cfbeSBarry Smith PetscErrorCode  MatMFFDAddNullSpace(Mat J,MatNullSpace nullsp)
1082e884886fSBarry Smith {
1083e884886fSBarry Smith   PetscErrorCode ierr;
1084e884886fSBarry Smith 
1085e884886fSBarry Smith   PetscFunctionBegin;
1086bc0f08ceSBarry Smith   ierr = PetscTryMethod(J,"MatMFFDAddNullSpace_C",(Mat,MatNullSpace),(J,nullsp));CHKERRQ(ierr);
1087e884886fSBarry Smith   PetscFunctionReturn(0);
1088e884886fSBarry Smith }
1089e884886fSBarry Smith 
1090e884886fSBarry Smith #undef __FUNCT__
1091e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetHHistory"
1092e884886fSBarry Smith /*@
1093e884886fSBarry Smith    MatMFFDSetHHistory - Sets an array to collect a history of the
1094e884886fSBarry Smith    differencing values (h) computed for the matrix-free product.
1095e884886fSBarry Smith 
10963f9fe445SBarry Smith    Logically Collective on Mat
1097e884886fSBarry Smith 
1098e884886fSBarry Smith    Input Parameters:
1099e884886fSBarry Smith +  J - the matrix-free matrix context
1100e884886fSBarry Smith .  histroy - space to hold the history
1101e884886fSBarry Smith -  nhistory - number of entries in history, if more entries are generated than
1102e884886fSBarry Smith               nhistory, then the later ones are discarded
1103e884886fSBarry Smith 
1104e884886fSBarry Smith    Level: advanced
1105e884886fSBarry Smith 
1106e884886fSBarry Smith    Notes:
1107e884886fSBarry Smith    Use MatMFFDResetHHistory() to reset the history counter and collect
1108e884886fSBarry Smith    a new batch of differencing parameters, h.
1109e884886fSBarry Smith 
1110e884886fSBarry Smith .keywords: SNES, matrix-free, h history, differencing history
1111e884886fSBarry Smith 
1112e884886fSBarry Smith .seealso: MatMFFDGetH(), MatCreateSNESMF(),
11131d0fab5eSBarry Smith           MatMFFDResetHHistory(), MatMFFDSetFunctionError()
1114e884886fSBarry Smith 
1115e884886fSBarry Smith @*/
11167087cfbeSBarry Smith PetscErrorCode  MatMFFDSetHHistory(Mat J,PetscScalar history[],PetscInt nhistory)
1117e884886fSBarry Smith {
1118e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)J->data;
1119bc0f08ceSBarry Smith   PetscErrorCode ierr;
1120bc0f08ceSBarry Smith   PetscBool      match;
1121e884886fSBarry Smith 
1122e884886fSBarry Smith   PetscFunctionBegin;
1123bc0f08ceSBarry Smith   ierr = PetscTypeCompare((PetscObject)J,MATMFFD,&match);CHKERRQ(ierr);
1124bc0f08ceSBarry Smith   if (!match) SETERRQ(((PetscObject)J)->comm,PETSC_ERR_ARG_WRONG,"Not a MFFD matrix");
1125e884886fSBarry Smith   ctx->historyh    = history;
1126e884886fSBarry Smith   ctx->maxcurrenth = nhistory;
112775567043SBarry Smith   ctx->currenth    = 0.;
1128e884886fSBarry Smith   PetscFunctionReturn(0);
1129e884886fSBarry Smith }
1130e884886fSBarry Smith 
1131bc0f08ceSBarry Smith 
1132e884886fSBarry Smith #undef __FUNCT__
1133e884886fSBarry Smith #define __FUNCT__ "MatMFFDResetHHistory"
1134e884886fSBarry Smith /*@
1135e884886fSBarry Smith    MatMFFDResetHHistory - Resets the counter to zero to begin
1136e884886fSBarry Smith    collecting a new set of differencing histories.
1137e884886fSBarry Smith 
11383f9fe445SBarry Smith    Logically Collective on Mat
1139e884886fSBarry Smith 
1140e884886fSBarry Smith    Input Parameters:
1141e884886fSBarry Smith .  J - the matrix-free matrix context
1142e884886fSBarry Smith 
1143e884886fSBarry Smith    Level: advanced
1144e884886fSBarry Smith 
1145e884886fSBarry Smith    Notes:
1146e884886fSBarry Smith    Use MatMFFDSetHHistory() to create the original history counter.
1147e884886fSBarry Smith 
1148e884886fSBarry Smith .keywords: SNES, matrix-free, h history, differencing history
1149e884886fSBarry Smith 
1150e884886fSBarry Smith .seealso: MatMFFDGetH(), MatCreateSNESMF(),
11511d0fab5eSBarry Smith           MatMFFDSetHHistory(), MatMFFDSetFunctionError()
1152e884886fSBarry Smith 
1153e884886fSBarry Smith @*/
11547087cfbeSBarry Smith PetscErrorCode  MatMFFDResetHHistory(Mat J)
1155e884886fSBarry Smith {
1156bc0f08ceSBarry Smith   PetscErrorCode ierr;
1157e884886fSBarry Smith 
1158e884886fSBarry Smith   PetscFunctionBegin;
1159bc0f08ceSBarry Smith   ierr = PetscTryMethod(J,"MatMFFDResetHHistory_C",(Mat),(J));CHKERRQ(ierr);
1160e884886fSBarry Smith   PetscFunctionReturn(0);
1161e884886fSBarry Smith }
1162e884886fSBarry Smith 
1163e884886fSBarry Smith 
1164e884886fSBarry Smith #undef __FUNCT__
1165e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetBase"
1166e884886fSBarry Smith /*@
1167e884886fSBarry Smith     MatMFFDSetBase - Sets the vector U at which matrix vector products of the
1168e884886fSBarry Smith         Jacobian are computed
1169e884886fSBarry Smith 
11703f9fe445SBarry Smith     Logically Collective on Mat
1171e884886fSBarry Smith 
1172e884886fSBarry Smith     Input Parameters:
1173e884886fSBarry Smith +   J - the MatMFFD matrix
11743ec795f1SBarry Smith .   U - the vector
1175bcddec3dSBarry Smith -   F - (optional) vector that contains F(u) if it has been already computed
1176e884886fSBarry Smith 
1177e884886fSBarry Smith     Notes: This is rarely used directly
1178e884886fSBarry Smith 
11798af5ae88SBarry Smith     If F is provided then it is not recomputed. Otherwise the function is evaluated at the base
11808af5ae88SBarry Smith     point during the first MatMult() after each call to MatMFFDSetBase().
1181dff2f722SBarry Smith 
1182e884886fSBarry Smith     Level: advanced
1183e884886fSBarry Smith 
1184e884886fSBarry Smith @*/
11857087cfbeSBarry Smith PetscErrorCode  MatMFFDSetBase(Mat J,Vec U,Vec F)
1186e884886fSBarry Smith {
11874ac538c5SBarry Smith   PetscErrorCode ierr;
1188e884886fSBarry Smith 
1189e884886fSBarry Smith   PetscFunctionBegin;
11900700a824SBarry Smith   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
11910700a824SBarry Smith   PetscValidHeaderSpecific(U,VEC_CLASSID,2);
11920700a824SBarry Smith   if (F) PetscValidHeaderSpecific(F,VEC_CLASSID,3);
11934ac538c5SBarry Smith   ierr = PetscTryMethod(J,"MatMFFDSetBase_C",(Mat,Vec,Vec),(J,U,F));CHKERRQ(ierr);
1194e884886fSBarry Smith   PetscFunctionReturn(0);
1195e884886fSBarry Smith }
1196e884886fSBarry Smith 
1197e884886fSBarry Smith #undef __FUNCT__
1198e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetCheckh"
1199e884886fSBarry Smith /*@C
1200e884886fSBarry Smith     MatMFFDSetCheckh - Sets a function that checks the computed h and adjusts
1201e884886fSBarry Smith         it to satisfy some criteria
1202e884886fSBarry Smith 
12033f9fe445SBarry Smith     Logically Collective on Mat
1204e884886fSBarry Smith 
1205e884886fSBarry Smith     Input Parameters:
1206e884886fSBarry Smith +   J - the MatMFFD matrix
1207e884886fSBarry Smith .   fun - the function that checks h
1208e884886fSBarry Smith -   ctx - any context needed by the function
1209e884886fSBarry Smith 
1210e884886fSBarry Smith     Options Database Keys:
1211e884886fSBarry Smith .   -mat_mffd_check_positivity
1212e884886fSBarry Smith 
1213e884886fSBarry Smith     Level: advanced
1214e884886fSBarry Smith 
1215e884886fSBarry Smith     Notes: For example, MatMFFDSetCheckPositivity() insures that all entries
1216e884886fSBarry Smith        of U + h*a are non-negative
1217e884886fSBarry Smith 
1218e884886fSBarry Smith .seealso:  MatMFFDSetCheckPositivity()
1219e884886fSBarry Smith @*/
12207087cfbeSBarry Smith PetscErrorCode  MatMFFDSetCheckh(Mat J,PetscErrorCode (*fun)(void*,Vec,Vec,PetscScalar*),void* ctx)
1221e884886fSBarry Smith {
12224ac538c5SBarry Smith   PetscErrorCode ierr;
1223e884886fSBarry Smith 
1224e884886fSBarry Smith   PetscFunctionBegin;
12250700a824SBarry Smith   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
12264ac538c5SBarry Smith   ierr = PetscTryMethod(J,"MatMFFDSetCheckh_C",(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*),(J,fun,ctx));CHKERRQ(ierr);
1227e884886fSBarry Smith   PetscFunctionReturn(0);
1228e884886fSBarry Smith }
1229e884886fSBarry Smith 
1230e884886fSBarry Smith #undef __FUNCT__
1231e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetCheckPositivity"
1232e884886fSBarry Smith /*@
1233e884886fSBarry Smith     MatMFFDCheckPositivity - Checks that all entries in U + h*a are positive or
1234e884886fSBarry Smith         zero, decreases h until this is satisfied.
1235e884886fSBarry Smith 
12363f9fe445SBarry Smith     Logically Collective on Vec
1237e884886fSBarry Smith 
1238e884886fSBarry Smith     Input Parameters:
1239e884886fSBarry Smith +   U - base vector that is added to
1240e884886fSBarry Smith .   a - vector that is added
1241e884886fSBarry Smith .   h - scaling factor on a
1242e884886fSBarry Smith -   dummy - context variable (unused)
1243e884886fSBarry Smith 
1244e884886fSBarry Smith     Options Database Keys:
1245e884886fSBarry Smith .   -mat_mffd_check_positivity
1246e884886fSBarry Smith 
1247e884886fSBarry Smith     Level: advanced
1248e884886fSBarry Smith 
1249e884886fSBarry Smith     Notes: This is rarely used directly, rather it is passed as an argument to
1250e884886fSBarry Smith            MatMFFDSetCheckh()
1251e884886fSBarry Smith 
1252e884886fSBarry Smith .seealso:  MatMFFDSetCheckh()
1253e884886fSBarry Smith @*/
12547087cfbeSBarry Smith PetscErrorCode  MatMFFDCheckPositivity(void* dummy,Vec U,Vec a,PetscScalar *h)
1255e884886fSBarry Smith {
1256e884886fSBarry Smith   PetscReal      val, minval;
1257e884886fSBarry Smith   PetscScalar    *u_vec, *a_vec;
1258e884886fSBarry Smith   PetscErrorCode ierr;
1259e884886fSBarry Smith   PetscInt       i,n;
1260e884886fSBarry Smith   MPI_Comm       comm;
1261e884886fSBarry Smith 
1262e884886fSBarry Smith   PetscFunctionBegin;
1263e884886fSBarry Smith   ierr = PetscObjectGetComm((PetscObject)U,&comm);CHKERRQ(ierr);
1264e884886fSBarry Smith   ierr = VecGetArray(U,&u_vec);CHKERRQ(ierr);
1265e884886fSBarry Smith   ierr = VecGetArray(a,&a_vec);CHKERRQ(ierr);
1266e884886fSBarry Smith   ierr = VecGetLocalSize(U,&n);CHKERRQ(ierr);
1267e884886fSBarry Smith   minval = PetscAbsScalar(*h*1.01);
1268e884886fSBarry Smith   for(i=0;i<n;i++) {
1269e884886fSBarry Smith     if (PetscRealPart(u_vec[i] + *h*a_vec[i]) <= 0.0) {
1270e884886fSBarry Smith       val = PetscAbsScalar(u_vec[i]/a_vec[i]);
1271e884886fSBarry Smith       if (val < minval) minval = val;
1272e884886fSBarry Smith     }
1273e884886fSBarry Smith   }
1274e884886fSBarry Smith   ierr = VecRestoreArray(U,&u_vec);CHKERRQ(ierr);
1275e884886fSBarry Smith   ierr = VecRestoreArray(a,&a_vec);CHKERRQ(ierr);
1276d9822059SBarry Smith   ierr = MPI_Allreduce(&minval,&val,1,MPIU_REAL,MPIU_MIN,comm);CHKERRQ(ierr);
1277e884886fSBarry Smith   if (val <= PetscAbsScalar(*h)) {
1278e884886fSBarry Smith     ierr = PetscInfo2(U,"Scaling back h from %G to %G\n",PetscRealPart(*h),.99*val);CHKERRQ(ierr);
1279e884886fSBarry Smith     if (PetscRealPart(*h) > 0.0) *h =  0.99*val;
1280e884886fSBarry Smith     else                         *h = -0.99*val;
1281e884886fSBarry Smith   }
1282e884886fSBarry Smith   PetscFunctionReturn(0);
1283e884886fSBarry Smith }
1284e884886fSBarry Smith 
1285e884886fSBarry Smith 
1286e884886fSBarry Smith 
1287e884886fSBarry Smith 
1288e884886fSBarry Smith 
1289e884886fSBarry Smith 
1290e884886fSBarry Smith 
1291