xref: /petsc/src/mat/impls/mffd/mffd.c (revision e94e781be4d0de67afa8d29cbcd676556dbc0369)
1e884886fSBarry Smith 
2af0996ceSBarry Smith #include <petsc/private/matimpl.h>
3c6db04a5SJed Brown #include <../src/mat/impls/mffd/mffdimpl.h>   /*I  "petscmat.h"   I*/
4e884886fSBarry Smith 
5140e18c1SBarry Smith PetscFunctionList MatMFFDList              = 0;
6ace3abfcSBarry Smith PetscBool         MatMFFDRegisterAllCalled = PETSC_FALSE;
7e884886fSBarry Smith 
87087cfbeSBarry Smith PetscClassId  MATMFFD_CLASSID;
9166c7f25SBarry Smith PetscLogEvent MATMFFD_Mult;
10e884886fSBarry Smith 
11ace3abfcSBarry Smith static PetscBool MatMFFDPackageInitialized = PETSC_FALSE;
12b022a5c1SBarry Smith /*@C
132390153bSJed Brown   MatMFFDFinalizePackage - This function destroys everything in the MatMFFD package. It is
14b022a5c1SBarry Smith   called from PetscFinalize().
15b022a5c1SBarry Smith 
16b022a5c1SBarry Smith   Level: developer
17b022a5c1SBarry Smith 
18755b7f64SBarry Smith .seealso: PetscFinalize(), MatCreateMFFD(), MatCreateSNESMF()
19b022a5c1SBarry Smith @*/
207087cfbeSBarry Smith PetscErrorCode  MatMFFDFinalizePackage(void)
21b022a5c1SBarry Smith {
22a703d84cSPatrick Lacasse   PetscErrorCode ierr;
23a703d84cSPatrick Lacasse 
24b022a5c1SBarry Smith   PetscFunctionBegin;
2537e93019SBarry Smith   ierr = PetscFunctionListDestroy(&MatMFFDList);CHKERRQ(ierr);
26b022a5c1SBarry Smith   MatMFFDPackageInitialized = PETSC_FALSE;
27b022a5c1SBarry Smith   MatMFFDRegisterAllCalled  = PETSC_FALSE;
28b022a5c1SBarry Smith   PetscFunctionReturn(0);
29b022a5c1SBarry Smith }
30b022a5c1SBarry Smith 
313ec795f1SBarry Smith /*@C
323ec795f1SBarry Smith   MatMFFDInitializePackage - This function initializes everything in the MatMFFD package. It is called
338a690491SBarry Smith   from MatInitializePackage().
343ec795f1SBarry Smith 
353ec795f1SBarry Smith   Level: developer
363ec795f1SBarry Smith 
373ec795f1SBarry Smith .seealso: PetscInitialize()
383ec795f1SBarry Smith @*/
39607a6623SBarry Smith PetscErrorCode  MatMFFDInitializePackage(void)
403ec795f1SBarry Smith {
413ec795f1SBarry Smith   char           logList[256];
428e81d068SLisandro Dalcin   PetscBool      opt,pkg;
433ec795f1SBarry Smith   PetscErrorCode ierr;
443ec795f1SBarry Smith 
453ec795f1SBarry Smith   PetscFunctionBegin;
46b022a5c1SBarry Smith   if (MatMFFDPackageInitialized) PetscFunctionReturn(0);
47b022a5c1SBarry Smith   MatMFFDPackageInitialized = PETSC_TRUE;
483ec795f1SBarry Smith   /* Register Classes */
490700a824SBarry Smith   ierr = PetscClassIdRegister("MatMFFD",&MATMFFD_CLASSID);CHKERRQ(ierr);
503ec795f1SBarry Smith   /* Register Constructors */
51607a6623SBarry Smith   ierr = MatMFFDRegisterAll();CHKERRQ(ierr);
523ec795f1SBarry Smith   /* Register Events */
530700a824SBarry Smith   ierr = PetscLogEventRegister("MatMult MF",MATMFFD_CLASSID,&MATMFFD_Mult);CHKERRQ(ierr);
54*e94e781bSJacob Faibussowitsch  /* Process Info */
55*e94e781bSJacob Faibussowitsch   {
56*e94e781bSJacob Faibussowitsch     PetscClassId  classids[1];
57*e94e781bSJacob Faibussowitsch 
58*e94e781bSJacob Faibussowitsch     classids[0] = MATMFFD_CLASSID;
59*e94e781bSJacob Faibussowitsch     ierr = PetscInfoProcessClass("matmffd", 1, classids);CHKERRQ(ierr);
603ec795f1SBarry Smith   }
613ec795f1SBarry Smith   /* Process summary exclusions */
628e81d068SLisandro Dalcin   ierr = PetscOptionsGetString(NULL,NULL,"-log_exclude",logList,sizeof(logList),&opt);CHKERRQ(ierr);
633ec795f1SBarry Smith   if (opt) {
648e81d068SLisandro Dalcin     ierr = PetscStrInList("matmffd",logList,',',&pkg);CHKERRQ(ierr);
65fa2bb9feSLisandro Dalcin     if (pkg) {ierr = PetscLogEventExcludeClass(MATMFFD_CLASSID);CHKERRQ(ierr);}
663ec795f1SBarry Smith   }
678e81d068SLisandro Dalcin   /* Register package finalizer */
68b022a5c1SBarry Smith   ierr = PetscRegisterFinalize(MatMFFDFinalizePackage);CHKERRQ(ierr);
693ec795f1SBarry Smith   PetscFunctionReturn(0);
703ec795f1SBarry Smith }
713ec795f1SBarry Smith 
72789d8953SBarry Smith static PetscErrorCode  MatMFFDSetType_MFFD(Mat mat,MatMFFDType ftype)
73789d8953SBarry Smith {
74789d8953SBarry Smith   PetscErrorCode ierr,(*r)(MatMFFD);
75789d8953SBarry Smith   MatMFFD        ctx;
76789d8953SBarry Smith   PetscBool      match;
77789d8953SBarry Smith 
78789d8953SBarry Smith   PetscFunctionBegin;
79789d8953SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
80789d8953SBarry Smith   PetscValidCharPointer(ftype,2);
81789d8953SBarry Smith   ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr);
82789d8953SBarry Smith 
83789d8953SBarry Smith   /* already set, so just return */
84789d8953SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)ctx,ftype,&match);CHKERRQ(ierr);
85789d8953SBarry Smith   if (match) PetscFunctionReturn(0);
86789d8953SBarry Smith 
87789d8953SBarry Smith   /* destroy the old one if it exists */
88789d8953SBarry Smith   if (ctx->ops->destroy) {
89789d8953SBarry Smith     ierr = (*ctx->ops->destroy)(ctx);CHKERRQ(ierr);
90789d8953SBarry Smith   }
91789d8953SBarry Smith 
92789d8953SBarry Smith   ierr =  PetscFunctionListFind(MatMFFDList,ftype,&r);CHKERRQ(ierr);
93789d8953SBarry Smith   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown MatMFFD type %s given",ftype);
94789d8953SBarry Smith   ierr = (*r)(ctx);CHKERRQ(ierr);
95789d8953SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)ctx,ftype);CHKERRQ(ierr);
96789d8953SBarry Smith   PetscFunctionReturn(0);
97789d8953SBarry Smith }
98789d8953SBarry Smith 
99e884886fSBarry Smith /*@C
100e884886fSBarry Smith     MatMFFDSetType - Sets the method that is used to compute the
101e884886fSBarry Smith     differencing parameter for finite differene matrix-free formulations.
102e884886fSBarry Smith 
103e884886fSBarry Smith     Input Parameters:
104e884886fSBarry Smith +   mat - the "matrix-free" matrix created via MatCreateSNESMF(), or MatCreateMFFD()
105e884886fSBarry Smith           or MatSetType(mat,MATMFFD);
106e884886fSBarry Smith -   ftype - the type requested, either MATMFFD_WP or MATMFFD_DS
107e884886fSBarry Smith 
108e884886fSBarry Smith     Level: advanced
109e884886fSBarry Smith 
110e884886fSBarry Smith     Notes:
111e884886fSBarry Smith     For example, such routines can compute h for use in
112e884886fSBarry Smith     Jacobian-vector products of the form
113e884886fSBarry Smith 
114e884886fSBarry Smith                         F(x+ha) - F(x)
115e884886fSBarry Smith           F'(u)a  ~=  ----------------
116e884886fSBarry Smith                               h
117e884886fSBarry Smith 
118755b7f64SBarry Smith .seealso: MatCreateSNESMF(), MatMFFDRegister(), MatMFFDSetFunction(), MatCreateMFFD()
119e884886fSBarry Smith @*/
12019fd82e9SBarry Smith PetscErrorCode  MatMFFDSetType(Mat mat,MatMFFDType ftype)
121e884886fSBarry Smith {
122789d8953SBarry Smith   PetscErrorCode ierr;
123e884886fSBarry Smith 
124e884886fSBarry Smith   PetscFunctionBegin;
1250700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
126e884886fSBarry Smith   PetscValidCharPointer(ftype,2);
127789d8953SBarry Smith   ierr = PetscTryMethod(mat,"MatMFFDSetType_C",(Mat,MatMFFDType),(mat,ftype));CHKERRQ(ierr);
128e884886fSBarry Smith   PetscFunctionReturn(0);
129e884886fSBarry Smith }
130e884886fSBarry Smith 
1315d21e1f8SStefano Zampini static PetscErrorCode MatGetDiagonal_MFFD(Mat,Vec);
1325d21e1f8SStefano Zampini 
133e884886fSBarry Smith typedef PetscErrorCode (*FCN1)(void*,Vec); /* force argument to next function to not be extern C*/
13440244768SBarry Smith static PetscErrorCode  MatMFFDSetFunctioniBase_MFFD(Mat mat,FCN1 func)
135e884886fSBarry Smith {
136789d8953SBarry Smith   MatMFFD        ctx;
137789d8953SBarry Smith   PetscErrorCode ierr;
138e884886fSBarry Smith 
139e884886fSBarry Smith   PetscFunctionBegin;
140789d8953SBarry Smith   ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr);
141e884886fSBarry Smith   ctx->funcisetbase = func;
142e884886fSBarry Smith   PetscFunctionReturn(0);
143e884886fSBarry Smith }
144e884886fSBarry Smith 
145e884886fSBarry Smith typedef PetscErrorCode (*FCN2)(void*,PetscInt,Vec,PetscScalar*); /* force argument to next function to not be extern C*/
14640244768SBarry Smith static PetscErrorCode  MatMFFDSetFunctioni_MFFD(Mat mat,FCN2 funci)
147e884886fSBarry Smith {
148789d8953SBarry Smith   MatMFFD        ctx;
149789d8953SBarry Smith   PetscErrorCode ierr;
150e884886fSBarry Smith 
151e884886fSBarry Smith   PetscFunctionBegin;
152789d8953SBarry Smith   ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr);
153e884886fSBarry Smith   ctx->funci = funci;
154789d8953SBarry Smith   ierr = MatShellSetOperation(mat,MATOP_GET_DIAGONAL,(void (*)(void))MatGetDiagonal_MFFD);CHKERRQ(ierr);
155789d8953SBarry Smith   PetscFunctionReturn(0);
1565d21e1f8SStefano Zampini }
157789d8953SBarry Smith 
158789d8953SBarry Smith static PetscErrorCode MatMFFDGetH_MFFD(Mat mat,PetscScalar *h)
159789d8953SBarry Smith {
160789d8953SBarry Smith   MatMFFD        ctx;
161789d8953SBarry Smith   PetscErrorCode ierr;
162789d8953SBarry Smith 
163789d8953SBarry Smith   PetscFunctionBegin;
164789d8953SBarry Smith   ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr);
165789d8953SBarry Smith   *h = ctx->currenth;
166e884886fSBarry Smith   PetscFunctionReturn(0);
167e884886fSBarry Smith }
168e884886fSBarry Smith 
16940244768SBarry Smith static PetscErrorCode  MatMFFDResetHHistory_MFFD(Mat J)
170bc0f08ceSBarry Smith {
171789d8953SBarry Smith   MatMFFD        ctx;
172789d8953SBarry Smith   PetscErrorCode ierr;
173bc0f08ceSBarry Smith 
174bc0f08ceSBarry Smith   PetscFunctionBegin;
175789d8953SBarry Smith   ierr = MatShellGetContext(J,&ctx);CHKERRQ(ierr);
176bc0f08ceSBarry Smith   ctx->ncurrenth = 0;
177bc0f08ceSBarry Smith   PetscFunctionReturn(0);
178bc0f08ceSBarry Smith }
179e884886fSBarry Smith 
1801c84c290SBarry Smith /*@C
1811c84c290SBarry Smith    MatMFFDRegister - Adds a method to the MatMFFD registry.
1821c84c290SBarry Smith 
1831c84c290SBarry Smith    Not Collective
1841c84c290SBarry Smith 
1851c84c290SBarry Smith    Input Parameters:
1861c84c290SBarry Smith +  name_solver - name of a new user-defined compute-h module
1871c84c290SBarry Smith -  routine_create - routine to create method context
1881c84c290SBarry Smith 
1891c84c290SBarry Smith    Level: developer
1901c84c290SBarry Smith 
1911c84c290SBarry Smith    Notes:
1921c84c290SBarry Smith    MatMFFDRegister() may be called multiple times to add several user-defined solvers.
1931c84c290SBarry Smith 
1941c84c290SBarry Smith    Sample usage:
1951c84c290SBarry Smith .vb
196bdf89e91SBarry Smith    MatMFFDRegister("my_h",MyHCreate);
1971c84c290SBarry Smith .ve
1981c84c290SBarry Smith 
1991c84c290SBarry Smith    Then, your solver can be chosen with the procedural interface via
2001c84c290SBarry Smith $     MatMFFDSetType(mfctx,"my_h")
2011c84c290SBarry Smith    or at runtime via the option
202be7a6d03SBarry Smith $     -mat_mffd_type my_h
2031c84c290SBarry Smith 
2041c84c290SBarry Smith .seealso: MatMFFDRegisterAll(), MatMFFDRegisterDestroy()
2051c84c290SBarry Smith  @*/
206bdf89e91SBarry Smith PetscErrorCode  MatMFFDRegister(const char sname[],PetscErrorCode (*function)(MatMFFD))
207e884886fSBarry Smith {
208e884886fSBarry Smith   PetscErrorCode ierr;
209e884886fSBarry Smith 
210e884886fSBarry Smith   PetscFunctionBegin;
2119cc31a68SJed Brown   ierr = MatInitializePackage();CHKERRQ(ierr);
212a240a19fSJed Brown   ierr = PetscFunctionListAdd(&MatMFFDList,sname,function);CHKERRQ(ierr);
213e884886fSBarry Smith   PetscFunctionReturn(0);
214e884886fSBarry Smith }
215e884886fSBarry Smith 
216e884886fSBarry Smith /* ----------------------------------------------------------------------------------------*/
21740244768SBarry Smith static PetscErrorCode MatDestroy_MFFD(Mat mat)
218e884886fSBarry Smith {
219e884886fSBarry Smith   PetscErrorCode ierr;
220789d8953SBarry Smith   MatMFFD        ctx;
221e884886fSBarry Smith 
222e884886fSBarry Smith   PetscFunctionBegin;
223789d8953SBarry Smith   ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr);
2246bf464f9SBarry Smith   ierr = VecDestroy(&ctx->w);CHKERRQ(ierr);
225c51ad4d4SStefano Zampini   ierr = VecDestroy(&ctx->current_u);CHKERRQ(ierr);
226cfe22f5eSBarry Smith   if (ctx->current_f_allocated) {
2276bf464f9SBarry Smith     ierr = VecDestroy(&ctx->current_f);CHKERRQ(ierr);
228cfe22f5eSBarry Smith   }
229e884886fSBarry Smith   if (ctx->ops->destroy) {ierr = (*ctx->ops->destroy)(ctx);CHKERRQ(ierr);}
2306bf464f9SBarry Smith   ierr = PetscHeaderDestroy(&ctx);CHKERRQ(ierr);
231e884886fSBarry Smith 
232bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetBase_C",NULL);CHKERRQ(ierr);
233bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctioniBase_C",NULL);CHKERRQ(ierr);
234bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctioni_C",NULL);CHKERRQ(ierr);
235bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunction_C",NULL);CHKERRQ(ierr);
236bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctionError_C",NULL);CHKERRQ(ierr);
237bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetCheckh_C",NULL);CHKERRQ(ierr);
238bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetPeriod_C",NULL);CHKERRQ(ierr);
239bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDResetHHistory_C",NULL);CHKERRQ(ierr);
240789d8953SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetHHistory_C",NULL);CHKERRQ(ierr);
241789d8953SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetType_C",NULL);CHKERRQ(ierr);
242789d8953SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDGetH_C",NULL);CHKERRQ(ierr);
243e884886fSBarry Smith   PetscFunctionReturn(0);
244e884886fSBarry Smith }
245e884886fSBarry Smith 
246e884886fSBarry Smith /*
247e884886fSBarry Smith    MatMFFDView_MFFD - Views matrix-free parameters.
248e884886fSBarry Smith 
249e884886fSBarry Smith */
25040244768SBarry Smith static PetscErrorCode MatView_MFFD(Mat J,PetscViewer viewer)
251e884886fSBarry Smith {
252e884886fSBarry Smith   PetscErrorCode ierr;
253789d8953SBarry Smith   MatMFFD        ctx;
254a283635eSDmitry Karpeev   PetscBool      iascii, viewbase, viewfunction;
255a283635eSDmitry Karpeev   const char     *prefix;
256e884886fSBarry Smith 
257e884886fSBarry Smith   PetscFunctionBegin;
258789d8953SBarry Smith   ierr = MatShellGetContext(J,&ctx);CHKERRQ(ierr);
259251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
260e884886fSBarry Smith   if (iascii) {
261a283635eSDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"Matrix-free approximation:\n");CHKERRQ(ierr);
262a283635eSDmitry Karpeev     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
26357622a8eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"err=%g (relative error in function evaluation)\n",(double)ctx->error_rel);CHKERRQ(ierr);
2647adad957SLisandro Dalcin     if (!((PetscObject)ctx)->type_name) {
265e884886fSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"The compute h routine has not yet been set\n");CHKERRQ(ierr);
266e884886fSBarry Smith     } else {
2677adad957SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"Using %s compute h routine\n",((PetscObject)ctx)->type_name);CHKERRQ(ierr);
268e884886fSBarry Smith     }
269c92e3469SBarry Smith #if defined(PETSC_USE_COMPLEX)
270c92e3469SBarry Smith     if (ctx->usecomplex) {
271c92e3469SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"Using Lyness complex number trick to compute the matrix-vector product\n");CHKERRQ(ierr);
272c92e3469SBarry Smith     }
273c92e3469SBarry Smith #endif
274e884886fSBarry Smith     if (ctx->ops->view) {
275e884886fSBarry Smith       ierr = (*ctx->ops->view)(ctx,viewer);CHKERRQ(ierr);
276e884886fSBarry Smith     }
277a283635eSDmitry Karpeev     ierr = PetscObjectGetOptionsPrefix((PetscObject)J, &prefix);CHKERRQ(ierr);
278a283635eSDmitry Karpeev 
279c5929fdfSBarry Smith     ierr = PetscOptionsHasName(((PetscObject)J)->options,prefix, "-mat_mffd_view_base", &viewbase);CHKERRQ(ierr);
280a283635eSDmitry Karpeev     if (viewbase) {
281a283635eSDmitry Karpeev       ierr = PetscViewerASCIIPrintf(viewer, "Base:\n");CHKERRQ(ierr);
282a283635eSDmitry Karpeev       ierr = VecView(ctx->current_u, viewer);CHKERRQ(ierr);
283a283635eSDmitry Karpeev     }
284c5929fdfSBarry Smith     ierr = PetscOptionsHasName(((PetscObject)J)->options,prefix, "-mat_mffd_view_function", &viewfunction);CHKERRQ(ierr);
285a283635eSDmitry Karpeev     if (viewfunction) {
286a283635eSDmitry Karpeev       ierr = PetscViewerASCIIPrintf(viewer, "Function:\n");CHKERRQ(ierr);
287a283635eSDmitry Karpeev       ierr = VecView(ctx->current_f, viewer);CHKERRQ(ierr);
288a283635eSDmitry Karpeev     }
289a283635eSDmitry Karpeev     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
29011aeaf0aSBarry Smith   }
291e884886fSBarry Smith   PetscFunctionReturn(0);
292e884886fSBarry Smith }
293e884886fSBarry Smith 
294e884886fSBarry Smith /*
295e884886fSBarry Smith    MatAssemblyEnd_MFFD - Resets the ctx->ncurrenth to zero. This
296e884886fSBarry Smith    allows the user to indicate the beginning of a new linear solve by calling
297e884886fSBarry Smith    MatAssemblyXXX() on the matrix free matrix. This then allows the
2981d0fab5eSBarry Smith    MatCreateMFFD_WP() to properly compute ||U|| only the first time
299e884886fSBarry Smith    in the linear solver rather than every time.
3005a576424SJed Brown 
301cc2e6a90SBarry Smith    This function is referenced directly from MatAssemblyEnd_SNESMF(), which may be in a different shared library hence
302cc2e6a90SBarry Smith    it must be labeled as PETSC_EXTERN
303e884886fSBarry Smith */
3045a576424SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyEnd_MFFD(Mat J,MatAssemblyType mt)
305e884886fSBarry Smith {
306e884886fSBarry Smith   PetscErrorCode ierr;
307789d8953SBarry Smith   MatMFFD        j;
308e884886fSBarry Smith 
309e884886fSBarry Smith   PetscFunctionBegin;
310789d8953SBarry Smith   ierr = MatShellGetContext(J,&j);CHKERRQ(ierr);
311e884886fSBarry Smith   ierr = MatMFFDResetHHistory(J);CHKERRQ(ierr);
312e884886fSBarry Smith   PetscFunctionReturn(0);
313e884886fSBarry Smith }
314e884886fSBarry Smith 
315e884886fSBarry Smith /*
316e884886fSBarry Smith   MatMult_MFFD - Default matrix-free form for Jacobian-vector product, y = F'(u)*a:
317e884886fSBarry Smith 
318e884886fSBarry Smith         y ~= (F(u + ha) - F(u))/h,
319e884886fSBarry Smith   where F = nonlinear function, as set by SNESSetFunction()
320e884886fSBarry Smith         u = current iterate
321e884886fSBarry Smith         h = difference interval
322e884886fSBarry Smith */
32340244768SBarry Smith static PetscErrorCode MatMult_MFFD(Mat mat,Vec a,Vec y)
324e884886fSBarry Smith {
325789d8953SBarry Smith   MatMFFD        ctx;
326e884886fSBarry Smith   PetscScalar    h;
327e884886fSBarry Smith   Vec            w,U,F;
328e884886fSBarry Smith   PetscErrorCode ierr;
329ace3abfcSBarry Smith   PetscBool      zeroa;
330e884886fSBarry Smith 
331e884886fSBarry Smith   PetscFunctionBegin;
332789d8953SBarry Smith   ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr);
333ce94432eSBarry 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");
334e884886fSBarry Smith   /* We log matrix-free matrix-vector products separately, so that we can
335e884886fSBarry Smith      separate the performance monitoring from the cases that use conventional
336e884886fSBarry Smith      storage.  We may eventually modify event logging to associate events
337e884886fSBarry Smith      with particular objects, hence alleviating the more general problem. */
338e884886fSBarry Smith   ierr = PetscLogEventBegin(MATMFFD_Mult,a,y,0,0);CHKERRQ(ierr);
339e884886fSBarry Smith 
340e884886fSBarry Smith   w = ctx->w;
341e884886fSBarry Smith   U = ctx->current_u;
3423ec795f1SBarry Smith   F = ctx->current_f;
343e884886fSBarry Smith   /*
344e884886fSBarry Smith       Compute differencing parameter
345e884886fSBarry Smith   */
3464a2f8832SBarry Smith   if (!((PetscObject)ctx)->type_name) {
347e884886fSBarry Smith     ierr = MatMFFDSetType(mat,MATMFFD_WP);CHKERRQ(ierr);
3489c6ac3b3SBarry Smith     ierr = MatSetFromOptions(mat);CHKERRQ(ierr);
349e884886fSBarry Smith   }
350e884886fSBarry Smith   ierr = (*ctx->ops->compute)(ctx,U,a,&h,&zeroa);CHKERRQ(ierr);
351e884886fSBarry Smith   if (zeroa) {
352e884886fSBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
353e884886fSBarry Smith     PetscFunctionReturn(0);
354e884886fSBarry Smith   }
355e884886fSBarry Smith 
35684d44b13SHong Zhang   if (mat->erroriffailure && PetscIsInfOrNanScalar(h)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Computed Nan differencing parameter h");
357e884886fSBarry Smith   if (ctx->checkh) {
358e884886fSBarry Smith     ierr = (*ctx->checkh)(ctx->checkhctx,U,a,&h);CHKERRQ(ierr);
359e884886fSBarry Smith   }
360e884886fSBarry Smith 
361e884886fSBarry Smith   /* keep a record of the current differencing parameter h */
362e884886fSBarry Smith   ctx->currenth = h;
363e884886fSBarry Smith #if defined(PETSC_USE_COMPLEX)
36457622a8eSBarry Smith   ierr = PetscInfo2(mat,"Current differencing parameter: %g + %g i\n",(double)PetscRealPart(h),(double)PetscImaginaryPart(h));CHKERRQ(ierr);
365e884886fSBarry Smith #else
366e884886fSBarry Smith   ierr = PetscInfo1(mat,"Current differencing parameter: %15.12e\n",h);CHKERRQ(ierr);
367e884886fSBarry Smith #endif
368e884886fSBarry Smith   if (ctx->historyh && ctx->ncurrenth < ctx->maxcurrenth) {
369e884886fSBarry Smith     ctx->historyh[ctx->ncurrenth] = h;
370e884886fSBarry Smith   }
371e884886fSBarry Smith   ctx->ncurrenth++;
372e884886fSBarry Smith 
373c92e3469SBarry Smith #if defined(PETSC_USE_COMPLEX)
374c92e3469SBarry Smith   if (ctx->usecomplex) h = PETSC_i*h;
375c92e3469SBarry Smith #endif
376c92e3469SBarry Smith 
377e884886fSBarry Smith   /* w = u + ha */
378e884886fSBarry Smith   ierr = VecWAXPY(w,h,a,U);CHKERRQ(ierr);
379e884886fSBarry Smith 
3801b797266SDmitry Karpeev   /* compute func(U) as base for differencing; only needed first time in and not when provided by user */
3811b797266SDmitry Karpeev   if (ctx->ncurrenth == 1 && ctx->current_f_allocated) {
382e884886fSBarry Smith     ierr = (*ctx->func)(ctx->funcctx,U,F);CHKERRQ(ierr);
38352121784SBarry Smith   }
384e884886fSBarry Smith   ierr = (*ctx->func)(ctx->funcctx,w,y);CHKERRQ(ierr);
385e884886fSBarry Smith 
386c92e3469SBarry Smith #if defined(PETSC_USE_COMPLEX)
387c92e3469SBarry Smith   if (ctx->usecomplex) {
388c92e3469SBarry Smith     ierr = VecImaginaryPart(y);CHKERRQ(ierr);
389c92e3469SBarry Smith     h    = PetscImaginaryPart(h);
390c92e3469SBarry Smith   } else {
391e884886fSBarry Smith     ierr = VecAXPY(y,-1.0,F);CHKERRQ(ierr);
392c92e3469SBarry Smith   }
393c92e3469SBarry Smith #else
394c92e3469SBarry Smith   ierr = VecAXPY(y,-1.0,F);CHKERRQ(ierr);
395c92e3469SBarry Smith #endif
396e884886fSBarry Smith   ierr = VecScale(y,1.0/h);CHKERRQ(ierr);
39739601f4eSBarry Smith   if (mat->nullsp) {ierr = MatNullSpaceRemove(mat->nullsp,y);CHKERRQ(ierr);}
398e884886fSBarry Smith 
399e884886fSBarry Smith   ierr = PetscLogEventEnd(MATMFFD_Mult,a,y,0,0);CHKERRQ(ierr);
400e884886fSBarry Smith   PetscFunctionReturn(0);
401e884886fSBarry Smith }
402e884886fSBarry Smith 
403e884886fSBarry Smith /*
404e884886fSBarry Smith   MatGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix
405e884886fSBarry Smith 
406e884886fSBarry Smith         y ~= (F(u + ha) - F(u))/h,
407e884886fSBarry Smith   where F = nonlinear function, as set by SNESSetFunction()
408e884886fSBarry Smith         u = current iterate
409e884886fSBarry Smith         h = difference interval
410e884886fSBarry Smith */
4115d21e1f8SStefano Zampini PetscErrorCode MatGetDiagonal_MFFD(Mat mat,Vec a)
412e884886fSBarry Smith {
413789d8953SBarry Smith   MatMFFD        ctx;
414e884886fSBarry Smith   PetscScalar    h,*aa,*ww,v;
415e884886fSBarry Smith   PetscReal      epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON;
416e884886fSBarry Smith   Vec            w,U;
417e884886fSBarry Smith   PetscErrorCode ierr;
418e884886fSBarry Smith   PetscInt       i,rstart,rend;
419e884886fSBarry Smith 
420e884886fSBarry Smith   PetscFunctionBegin;
421789d8953SBarry Smith   ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr);
422486afcceSStefano Zampini   if (!ctx->func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Requires calling MatMFFDSetFunction() first");
423486afcceSStefano Zampini   if (!ctx->funci) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Requires calling MatMFFDSetFunctioni() first");
424e884886fSBarry Smith   w    = ctx->w;
425e884886fSBarry Smith   U    = ctx->current_u;
426e884886fSBarry Smith   ierr = (*ctx->func)(ctx->funcctx,U,a);CHKERRQ(ierr);
427789d8953SBarry Smith   if (ctx->funcisetbase) {
428e884886fSBarry Smith     ierr = (*ctx->funcisetbase)(ctx->funcctx,U);CHKERRQ(ierr);
429789d8953SBarry Smith   }
430e884886fSBarry Smith   ierr = VecCopy(U,w);CHKERRQ(ierr);
431e884886fSBarry Smith 
432e884886fSBarry Smith   ierr = VecGetOwnershipRange(a,&rstart,&rend);CHKERRQ(ierr);
433e884886fSBarry Smith   ierr = VecGetArray(a,&aa);CHKERRQ(ierr);
434e884886fSBarry Smith   for (i=rstart; i<rend; i++) {
435e884886fSBarry Smith     ierr = VecGetArray(w,&ww);CHKERRQ(ierr);
436e884886fSBarry Smith     h    = ww[i-rstart];
437e884886fSBarry Smith     if (h == 0.0) h = 1.0;
438e884886fSBarry Smith     if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0)     h = umin;
439e884886fSBarry Smith     else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin;
440e884886fSBarry Smith     h *= epsilon;
441e884886fSBarry Smith 
442e884886fSBarry Smith     ww[i-rstart] += h;
443e884886fSBarry Smith     ierr          = VecRestoreArray(w,&ww);CHKERRQ(ierr);
444e884886fSBarry Smith     ierr          = (*ctx->funci)(ctx->funcctx,i,w,&v);CHKERRQ(ierr);
445e884886fSBarry Smith     aa[i-rstart]  = (v - aa[i-rstart])/h;
446e884886fSBarry Smith 
447e884886fSBarry Smith     ierr          = VecGetArray(w,&ww);CHKERRQ(ierr);
448e884886fSBarry Smith     ww[i-rstart] -= h;
449e884886fSBarry Smith     ierr          = VecRestoreArray(w,&ww);CHKERRQ(ierr);
450e884886fSBarry Smith   }
451e884886fSBarry Smith   ierr = VecRestoreArray(a,&aa);CHKERRQ(ierr);
452e884886fSBarry Smith   PetscFunctionReturn(0);
453e884886fSBarry Smith }
454e884886fSBarry Smith 
455d8fa6a54SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDSetBase_MFFD(Mat J,Vec U,Vec F)
456e884886fSBarry Smith {
457e884886fSBarry Smith   PetscErrorCode ierr;
458789d8953SBarry Smith   MatMFFD        ctx;
459e884886fSBarry Smith 
460e884886fSBarry Smith   PetscFunctionBegin;
461789d8953SBarry Smith   ierr = MatShellGetContext(J,&ctx);CHKERRQ(ierr);
462e884886fSBarry Smith   ierr = MatMFFDResetHHistory(J);CHKERRQ(ierr);
463c51ad4d4SStefano Zampini   if (!ctx->current_u) {
464c51ad4d4SStefano Zampini     ierr = VecDuplicate(U,&ctx->current_u);CHKERRQ(ierr);
4658860a134SJunchao Zhang     ierr = VecLockReadPush(ctx->current_u);CHKERRQ(ierr);
466c51ad4d4SStefano Zampini   }
4678860a134SJunchao Zhang   ierr = VecLockReadPop(ctx->current_u);CHKERRQ(ierr);
468c51ad4d4SStefano Zampini   ierr = VecCopy(U,ctx->current_u);CHKERRQ(ierr);
4698860a134SJunchao Zhang   ierr = VecLockReadPush(ctx->current_u);CHKERRQ(ierr);
47052121784SBarry Smith   if (F) {
4716bf464f9SBarry Smith     if (ctx->current_f_allocated) {ierr = VecDestroy(&ctx->current_f);CHKERRQ(ierr);}
4723ec795f1SBarry Smith     ctx->current_f           = F;
473cfe22f5eSBarry Smith     ctx->current_f_allocated = PETSC_FALSE;
474cfe22f5eSBarry Smith   } else if (!ctx->current_f_allocated) {
47506c4b6baSBarry Smith     ierr = MatCreateVecs(J,NULL,&ctx->current_f);CHKERRQ(ierr);
4762205254eSKarl Rupp 
477cfe22f5eSBarry Smith     ctx->current_f_allocated = PETSC_TRUE;
47852121784SBarry Smith   }
479e884886fSBarry Smith   if (!ctx->w) {
480e884886fSBarry Smith     ierr = VecDuplicate(ctx->current_u,&ctx->w);CHKERRQ(ierr);
481e884886fSBarry Smith   }
482e884886fSBarry Smith   J->assembled = PETSC_TRUE;
483e884886fSBarry Smith   PetscFunctionReturn(0);
484e884886fSBarry Smith }
4855a576424SJed Brown 
486e884886fSBarry Smith typedef PetscErrorCode (*FCN3)(void*,Vec,Vec,PetscScalar*); /* force argument to next function to not be extern C*/
487b2573a8aSBarry Smith 
48840244768SBarry Smith static PetscErrorCode  MatMFFDSetCheckh_MFFD(Mat J,FCN3 fun,void *ectx)
489e884886fSBarry Smith {
490789d8953SBarry Smith   MatMFFD        ctx;
491789d8953SBarry Smith   PetscErrorCode ierr;
492e884886fSBarry Smith 
493e884886fSBarry Smith   PetscFunctionBegin;
494789d8953SBarry Smith   ierr = MatShellGetContext(J,&ctx);CHKERRQ(ierr);
495e884886fSBarry Smith   ctx->checkh    = fun;
496e884886fSBarry Smith   ctx->checkhctx = ectx;
497e884886fSBarry Smith   PetscFunctionReturn(0);
498e884886fSBarry Smith }
499e884886fSBarry Smith 
5006aa9148fSLisandro Dalcin /*@C
5016aa9148fSLisandro Dalcin    MatMFFDSetOptionsPrefix - Sets the prefix used for searching for all
5026aa9148fSLisandro Dalcin    MatMFFD options in the database.
5036aa9148fSLisandro Dalcin 
5046aa9148fSLisandro Dalcin    Collective on Mat
5056aa9148fSLisandro Dalcin 
5066aa9148fSLisandro Dalcin    Input Parameter:
5076aa9148fSLisandro Dalcin +  A - the Mat context
5086aa9148fSLisandro Dalcin -  prefix - the prefix to prepend to all option names
5096aa9148fSLisandro Dalcin 
5106aa9148fSLisandro Dalcin    Notes:
5116aa9148fSLisandro Dalcin    A hyphen (-) must NOT be given at the beginning of the prefix name.
5126aa9148fSLisandro Dalcin    The first character of all runtime options is AUTOMATICALLY the hyphen.
5136aa9148fSLisandro Dalcin 
5146aa9148fSLisandro Dalcin    Level: advanced
5156aa9148fSLisandro Dalcin 
516755b7f64SBarry Smith .seealso: MatSetFromOptions(), MatCreateSNESMF(), MatCreateMFFD()
5176aa9148fSLisandro Dalcin @*/
5187087cfbeSBarry Smith PetscErrorCode  MatMFFDSetOptionsPrefix(Mat mat,const char prefix[])
5196aa9148fSLisandro Dalcin 
5206aa9148fSLisandro Dalcin {
521789d8953SBarry Smith   MatMFFD        mfctx;
5226aa9148fSLisandro Dalcin   PetscErrorCode ierr;
5235fd66863SKarl Rupp 
5246aa9148fSLisandro Dalcin   PetscFunctionBegin;
5250700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
526789d8953SBarry Smith   ierr = MatShellGetContext(mat,&mfctx);CHKERRQ(ierr);
5270700a824SBarry Smith   PetscValidHeaderSpecific(mfctx,MATMFFD_CLASSID,1);
5286aa9148fSLisandro Dalcin   ierr = PetscObjectSetOptionsPrefix((PetscObject)mfctx,prefix);CHKERRQ(ierr);
5296aa9148fSLisandro Dalcin   PetscFunctionReturn(0);
5306aa9148fSLisandro Dalcin }
5316aa9148fSLisandro Dalcin 
53240244768SBarry Smith static PetscErrorCode  MatSetFromOptions_MFFD(PetscOptionItems *PetscOptionsObject,Mat mat)
533e884886fSBarry Smith {
534789d8953SBarry Smith   MatMFFD        mfctx;
535e884886fSBarry Smith   PetscErrorCode ierr;
536ace3abfcSBarry Smith   PetscBool      flg;
537e884886fSBarry Smith   char           ftype[256];
538e884886fSBarry Smith 
539e884886fSBarry Smith   PetscFunctionBegin;
5400700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
541789d8953SBarry Smith   ierr = MatShellGetContext(mat,&mfctx);CHKERRQ(ierr);
5420700a824SBarry Smith   PetscValidHeaderSpecific(mfctx,MATMFFD_CLASSID,1);
5433194b578SJed Brown   ierr = PetscObjectOptionsBegin((PetscObject)mfctx);CHKERRQ(ierr);
544a264d7a6SBarry Smith   ierr = PetscOptionsFList("-mat_mffd_type","Matrix free type","MatMFFDSetType",MatMFFDList,((PetscObject)mfctx)->type_name,ftype,256,&flg);CHKERRQ(ierr);
545e884886fSBarry Smith   if (flg) {
546e884886fSBarry Smith     ierr = MatMFFDSetType(mat,ftype);CHKERRQ(ierr);
547e884886fSBarry Smith   }
548e884886fSBarry Smith 
549e884886fSBarry Smith   ierr = PetscOptionsReal("-mat_mffd_err","set sqrt relative error in function","MatMFFDSetFunctionError",mfctx->error_rel,&mfctx->error_rel,0);CHKERRQ(ierr);
550e884886fSBarry Smith   ierr = PetscOptionsInt("-mat_mffd_period","how often h is recomputed","MatMFFDSetPeriod",mfctx->recomputeperiod,&mfctx->recomputeperiod,0);CHKERRQ(ierr);
551e884886fSBarry Smith 
55290d69ab7SBarry Smith   flg  = PETSC_FALSE;
5530298fd71SBarry Smith   ierr = PetscOptionsBool("-mat_mffd_check_positivity","Insure that U + h*a is nonnegative","MatMFFDSetCheckh",flg,&flg,NULL);CHKERRQ(ierr);
554e884886fSBarry Smith   if (flg) {
555e884886fSBarry Smith     ierr = MatMFFDSetCheckh(mat,MatMFFDCheckPositivity,0);CHKERRQ(ierr);
556e884886fSBarry Smith   }
557c92e3469SBarry Smith #if defined(PETSC_USE_COMPLEX)
558c92e3469SBarry Smith   ierr = PetscOptionsBool("-mat_mffd_complex","Use Lyness complex number trick to compute the matrix-vector product","None",mfctx->usecomplex,&mfctx->usecomplex,NULL);CHKERRQ(ierr);
559c92e3469SBarry Smith #endif
560e884886fSBarry Smith   if (mfctx->ops->setfromoptions) {
561e55864a3SBarry Smith     ierr = (*mfctx->ops->setfromoptions)(PetscOptionsObject,mfctx);CHKERRQ(ierr);
562e884886fSBarry Smith   }
563e884886fSBarry Smith   ierr = PetscOptionsEnd();CHKERRQ(ierr);
564e884886fSBarry Smith   PetscFunctionReturn(0);
565e884886fSBarry Smith }
566e884886fSBarry Smith 
56740244768SBarry Smith static PetscErrorCode  MatMFFDSetPeriod_MFFD(Mat mat,PetscInt period)
568bc0f08ceSBarry Smith {
569789d8953SBarry Smith   MatMFFD        ctx;
570789d8953SBarry Smith   PetscErrorCode ierr;
571bc0f08ceSBarry Smith 
572bc0f08ceSBarry Smith   PetscFunctionBegin;
573789d8953SBarry Smith   ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr);
574bc0f08ceSBarry Smith   ctx->recomputeperiod = period;
575bc0f08ceSBarry Smith   PetscFunctionReturn(0);
576bc0f08ceSBarry Smith }
577bc0f08ceSBarry Smith 
57840244768SBarry Smith static PetscErrorCode  MatMFFDSetFunction_MFFD(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx)
579bc0f08ceSBarry Smith {
580789d8953SBarry Smith   MatMFFD        ctx;
581789d8953SBarry Smith   PetscErrorCode ierr;
582bc0f08ceSBarry Smith 
583bc0f08ceSBarry Smith   PetscFunctionBegin;
584789d8953SBarry Smith   ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr);
585bc0f08ceSBarry Smith   ctx->func    = func;
586bc0f08ceSBarry Smith   ctx->funcctx = funcctx;
587bc0f08ceSBarry Smith   PetscFunctionReturn(0);
588bc0f08ceSBarry Smith }
589bc0f08ceSBarry Smith 
59040244768SBarry Smith static PetscErrorCode  MatMFFDSetFunctionError_MFFD(Mat mat,PetscReal error)
591bc0f08ceSBarry Smith {
592789d8953SBarry Smith   MatMFFD        ctx;
593789d8953SBarry Smith   PetscErrorCode ierr;
594bc0f08ceSBarry Smith 
595bc0f08ceSBarry Smith   PetscFunctionBegin;
596789d8953SBarry Smith   ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr);
597bc0f08ceSBarry Smith   if (error != PETSC_DEFAULT) ctx->error_rel = error;
598bc0f08ceSBarry Smith   PetscFunctionReturn(0);
599bc0f08ceSBarry Smith }
600bc0f08ceSBarry Smith 
601789d8953SBarry Smith PetscErrorCode  MatMFFDSetHHistory_MFFD(Mat J,PetscScalar history[],PetscInt nhistory)
6023b49f96aSBarry Smith {
603789d8953SBarry Smith   MatMFFD        ctx;
604789d8953SBarry Smith   PetscErrorCode ierr;
605789d8953SBarry Smith 
6063b49f96aSBarry Smith   PetscFunctionBegin;
607789d8953SBarry Smith   ierr = MatShellGetContext(J,&ctx);CHKERRQ(ierr);
608789d8953SBarry Smith   ctx->historyh    = history;
609789d8953SBarry Smith   ctx->maxcurrenth = nhistory;
610789d8953SBarry Smith   ctx->currenth    = 0.;
6113b49f96aSBarry Smith   PetscFunctionReturn(0);
6123b49f96aSBarry Smith }
6133b49f96aSBarry Smith 
614e884886fSBarry Smith /*MC
615e884886fSBarry Smith   MATMFFD - MATMFFD = "mffd" - A matrix free matrix type.
616e884886fSBarry Smith 
617e884886fSBarry Smith   Level: advanced
618e884886fSBarry Smith 
619789d8953SBarry Smith   Developers Note: This is implemented on top of MATSHELL to get support for scaling and shifting without requiring duplicate code
620789d8953SBarry Smith 
621755b7f64SBarry Smith .seealso: MatCreateMFFD(), MatCreateSNESMF(), MatMFFDSetFunction(), MatMFFDSetType(),
622755b7f64SBarry Smith           MatMFFDSetFunctionError(), MatMFFDDSSetUmin(), MatMFFDSetFunction()
623755b7f64SBarry Smith           MatMFFDSetHHistory(), MatMFFDResetHHistory(), MatCreateSNESMF(),
624755b7f64SBarry Smith           MatMFFDGetH(),
625e884886fSBarry Smith M*/
6268cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_MFFD(Mat A)
627e884886fSBarry Smith {
628e884886fSBarry Smith   MatMFFD        mfctx;
629e884886fSBarry Smith   PetscErrorCode ierr;
630e884886fSBarry Smith 
631e884886fSBarry Smith   PetscFunctionBegin;
632607a6623SBarry Smith   ierr = MatMFFDInitializePackage();CHKERRQ(ierr);
633e884886fSBarry Smith 
634789d8953SBarry Smith   ierr = PetscHeaderCreate(mfctx,MATMFFD_CLASSID,"MatMFFD","Matrix-free Finite Differencing","Mat",PetscObjectComm((PetscObject)A),NULL,NULL);CHKERRQ(ierr);
6352205254eSKarl Rupp 
636e884886fSBarry Smith   mfctx->error_rel                = PETSC_SQRT_MACHINE_EPSILON;
637e884886fSBarry Smith   mfctx->recomputeperiod          = 1;
638e884886fSBarry Smith   mfctx->count                    = 0;
639e884886fSBarry Smith   mfctx->currenth                 = 0.0;
6400298fd71SBarry Smith   mfctx->historyh                 = NULL;
641e884886fSBarry Smith   mfctx->ncurrenth                = 0;
642e884886fSBarry Smith   mfctx->maxcurrenth              = 0;
6437adad957SLisandro Dalcin   ((PetscObject)mfctx)->type_name = 0;
644e884886fSBarry Smith 
645e884886fSBarry Smith   /*
646e884886fSBarry Smith      Create the empty data structure to contain compute-h routines.
647e884886fSBarry Smith      These will be filled in below from the command line options or
648e884886fSBarry Smith      a later call with MatMFFDSetType() or if that is not called
649e884886fSBarry Smith      then it will default in the first use of MatMult_MFFD()
650e884886fSBarry Smith   */
651e884886fSBarry Smith   mfctx->ops->compute        = 0;
652e884886fSBarry Smith   mfctx->ops->destroy        = 0;
653e884886fSBarry Smith   mfctx->ops->view           = 0;
654e884886fSBarry Smith   mfctx->ops->setfromoptions = 0;
655e884886fSBarry Smith   mfctx->hctx                = 0;
656e884886fSBarry Smith 
657e884886fSBarry Smith   mfctx->func    = 0;
658e884886fSBarry Smith   mfctx->funcctx = 0;
6590298fd71SBarry Smith   mfctx->w       = NULL;
660789d8953SBarry Smith   mfctx->mat     = A;
661e884886fSBarry Smith 
662789d8953SBarry Smith   ierr = MatSetType(A,MATSHELL);CHKERRQ(ierr);
663789d8953SBarry Smith   ierr = MatShellSetContext(A,mfctx);CHKERRQ(ierr);
664789d8953SBarry Smith   ierr = MatShellSetOperation(A,MATOP_MULT,(void (*)(void))MatMult_MFFD);CHKERRQ(ierr);
665789d8953SBarry Smith   ierr = MatShellSetOperation(A,MATOP_DESTROY,(void (*)(void))MatDestroy_MFFD);CHKERRQ(ierr);
666789d8953SBarry Smith   ierr = MatShellSetOperation(A,MATOP_VIEW,(void (*)(void))MatView_MFFD);CHKERRQ(ierr);
667789d8953SBarry Smith   ierr = MatShellSetOperation(A,MATOP_ASSEMBLY_END,(void (*)(void))MatAssemblyEnd_MFFD);CHKERRQ(ierr);
668789d8953SBarry Smith   ierr = MatShellSetOperation(A,MATOP_SET_FROM_OPTIONS,(void (*)(void))MatSetFromOptions_MFFD);CHKERRQ(ierr);
669e884886fSBarry Smith 
670bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetBase_C",MatMFFDSetBase_MFFD);CHKERRQ(ierr);
671bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctioniBase_C",MatMFFDSetFunctioniBase_MFFD);CHKERRQ(ierr);
672bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctioni_C",MatMFFDSetFunctioni_MFFD);CHKERRQ(ierr);
673bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunction_C",MatMFFDSetFunction_MFFD);CHKERRQ(ierr);
674bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetCheckh_C",MatMFFDSetCheckh_MFFD);CHKERRQ(ierr);
675bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetPeriod_C",MatMFFDSetPeriod_MFFD);CHKERRQ(ierr);
676bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctionError_C",MatMFFDSetFunctionError_MFFD);CHKERRQ(ierr);
677bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDResetHHistory_C",MatMFFDResetHHistory_MFFD);CHKERRQ(ierr);
678789d8953SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetHHistory_C",MatMFFDSetHHistory_MFFD);CHKERRQ(ierr);
679789d8953SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetType_C",MatMFFDSetType_MFFD);CHKERRQ(ierr);
680789d8953SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDGetH_C",MatMFFDGetH_MFFD);CHKERRQ(ierr);
681e884886fSBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,MATMFFD);CHKERRQ(ierr);
682e884886fSBarry Smith   PetscFunctionReturn(0);
683e884886fSBarry Smith }
684e884886fSBarry Smith 
685e884886fSBarry Smith /*@
686e884886fSBarry Smith    MatCreateMFFD - Creates a matrix-free matrix. See also MatCreateSNESMF()
687e884886fSBarry Smith 
688e884886fSBarry Smith    Collective on Vec
689e884886fSBarry Smith 
690e884886fSBarry Smith    Input Parameters:
691fef1beadSBarry Smith +  comm - MPI communicator
692fef1beadSBarry Smith .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
693fef1beadSBarry Smith            This value should be the same as the local size used in creating the
694fef1beadSBarry Smith            y vector for the matrix-vector product y = Ax.
695fef1beadSBarry Smith .  n - This value should be the same as the local size used in creating the
696fef1beadSBarry Smith        x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have
697fef1beadSBarry Smith        calculated if N is given) For square matrices n is almost always m.
698fef1beadSBarry Smith .  M - number of global rows (or PETSC_DETERMINE to have calculated if m is given)
699fef1beadSBarry Smith -  N - number of global columns (or PETSC_DETERMINE to have calculated if n is given)
700fef1beadSBarry Smith 
701e884886fSBarry Smith 
702e884886fSBarry Smith    Output Parameter:
703e884886fSBarry Smith .  J - the matrix-free matrix
704e884886fSBarry Smith 
7059c6ac3b3SBarry Smith    Options Database Keys: call MatSetFromOptions() to trigger these
7069c6ac3b3SBarry Smith +  -mat_mffd_type - wp or ds (see MATMFFD_WP or MATMFFD_DS)
707c92e3469SBarry Smith .  -mat_mffd_err - square root of estimated relative error in function evaluation
708c92e3469SBarry Smith .  -mat_mffd_period - how often h is recomputed, defaults to 1, everytime
709c92e3469SBarry Smith .  -mat_mffd_check_positivity - possibly decrease h until U + h*a has only positive values
710c92e3469SBarry Smith -  -mat_mffd_complex - use the Lyness trick with complex numbers to compute the matrix-vector product instead of differencing
711c92e3469SBarry Smith                        (requires real valued functions but that PETSc be configured for complex numbers)
7129c6ac3b3SBarry Smith 
7139c6ac3b3SBarry Smith 
714e884886fSBarry Smith    Level: advanced
715e884886fSBarry Smith 
716e884886fSBarry Smith    Notes:
717e884886fSBarry Smith    The matrix-free matrix context merely contains the function pointers
718e884886fSBarry Smith    and work space for performing finite difference approximations of
719e884886fSBarry Smith    Jacobian-vector products, F'(u)*a,
720e884886fSBarry Smith 
721e884886fSBarry Smith    The default code uses the following approach to compute h
722e884886fSBarry Smith 
723e884886fSBarry Smith .vb
724e884886fSBarry Smith      F'(u)*a = [F(u+h*a) - F(u)]/h where
725e884886fSBarry Smith      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
726e884886fSBarry Smith        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
727e884886fSBarry Smith  where
728e884886fSBarry Smith      error_rel = square root of relative error in function evaluation
729e884886fSBarry Smith      umin = minimum iterate parameter
730e884886fSBarry Smith .ve
731e884886fSBarry Smith 
732ca93e954SBarry Smith    You can call SNESSetJacobian() with MatMFFDComputeJacobian() if you are using matrix and not a different
733ca93e954SBarry Smith    preconditioner matrix
734ca93e954SBarry Smith 
735e884886fSBarry Smith    The user can set the error_rel via MatMFFDSetFunctionError() and
736a7f22e61SSatish Balay    umin via MatMFFDDSSetUmin(); see Users-Manual: ch_snes for details.
737e884886fSBarry Smith 
738e884886fSBarry Smith    The user should call MatDestroy() when finished with the matrix-free
739e884886fSBarry Smith    matrix context.
740e884886fSBarry Smith 
741e884886fSBarry Smith    Options Database Keys:
742e884886fSBarry Smith +  -mat_mffd_err <error_rel> - Sets error_rel
743e884886fSBarry Smith .  -mat_mffd_unim <umin> - Sets umin (for default PETSc routine that computes h only)
744e884886fSBarry Smith -  -mat_mffd_check_positivity
745e884886fSBarry Smith 
7461d0fab5eSBarry Smith .seealso: MatDestroy(), MatMFFDSetFunctionError(), MatMFFDDSSetUmin(), MatMFFDSetFunction()
747e884886fSBarry Smith           MatMFFDSetHHistory(), MatMFFDResetHHistory(), MatCreateSNESMF(),
74881242352SJed Brown           MatMFFDGetH(), MatMFFDRegister(), MatMFFDComputeJacobian()
749e884886fSBarry Smith 
750e884886fSBarry Smith @*/
7517087cfbeSBarry Smith PetscErrorCode  MatCreateMFFD(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *J)
752e884886fSBarry Smith {
753e884886fSBarry Smith   PetscErrorCode ierr;
754e884886fSBarry Smith 
755e884886fSBarry Smith   PetscFunctionBegin;
756e884886fSBarry Smith   ierr = MatCreate(comm,J);CHKERRQ(ierr);
757fef1beadSBarry Smith   ierr = MatSetSizes(*J,m,n,M,N);CHKERRQ(ierr);
758e884886fSBarry Smith   ierr = MatSetType(*J,MATMFFD);CHKERRQ(ierr);
759be7c243fSJed Brown   ierr = MatSetUp(*J);CHKERRQ(ierr);
760e884886fSBarry Smith   PetscFunctionReturn(0);
761e884886fSBarry Smith }
762e884886fSBarry Smith 
763e884886fSBarry Smith /*@
764e884886fSBarry Smith    MatMFFDGetH - Gets the last value that was used as the differencing
765e884886fSBarry Smith    parameter.
766e884886fSBarry Smith 
767e884886fSBarry Smith    Not Collective
768e884886fSBarry Smith 
769e884886fSBarry Smith    Input Parameters:
770e884886fSBarry Smith .  mat - the matrix obtained with MatCreateSNESMF()
771e884886fSBarry Smith 
772fd292e60Sprj-    Output Parameter:
773e884886fSBarry Smith .  h - the differencing step size
774e884886fSBarry Smith 
775e884886fSBarry Smith    Level: advanced
776e884886fSBarry Smith 
7771d0fab5eSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDSetHHistory(), MatCreateMFFD(), MATMFFD, MatMFFDResetHHistory()
778e884886fSBarry Smith @*/
7797087cfbeSBarry Smith PetscErrorCode  MatMFFDGetH(Mat mat,PetscScalar *h)
780e884886fSBarry Smith {
781bc0f08ceSBarry Smith   PetscErrorCode ierr;
782e884886fSBarry Smith 
783e884886fSBarry Smith   PetscFunctionBegin;
78488b4c220SStefano Zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
78588b4c220SStefano Zampini   PetscValidPointer(h,2);
786789d8953SBarry Smith   ierr = PetscUseMethod(mat,"MatMFFDGetH_C",(Mat,PetscScalar*),(mat,h));CHKERRQ(ierr);
787e884886fSBarry Smith   PetscFunctionReturn(0);
788e884886fSBarry Smith }
789e884886fSBarry Smith 
790e884886fSBarry Smith /*@C
791e884886fSBarry Smith    MatMFFDSetFunction - Sets the function used in applying the matrix free.
792e884886fSBarry Smith 
7933f9fe445SBarry Smith    Logically Collective on Mat
794e884886fSBarry Smith 
795e884886fSBarry Smith    Input Parameters:
79614f633e2SBarry Smith +  mat - the matrix free matrix created via MatCreateSNESMF() or MatCreateMFFD()
797e884886fSBarry Smith .  func - the function to use
798e884886fSBarry Smith -  funcctx - optional function context passed to function
799e884886fSBarry Smith 
80014f633e2SBarry Smith    Calling Sequence of func:
80114f633e2SBarry Smith $     func (void *funcctx, Vec x, Vec f)
80214f633e2SBarry Smith 
80314f633e2SBarry Smith +  funcctx - user provided context
80414f633e2SBarry Smith .  x - input vector
80514f633e2SBarry Smith -  f - computed output function
80614f633e2SBarry Smith 
807e884886fSBarry Smith    Level: advanced
808e884886fSBarry Smith 
809e884886fSBarry Smith    Notes:
810e884886fSBarry Smith     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
811e884886fSBarry Smith     matrix inside your compute Jacobian routine
812e884886fSBarry Smith 
8133ec795f1SBarry Smith     If this is not set then it will use the function set with SNESSetFunction() if MatCreateSNESMF() was used.
814e884886fSBarry Smith 
8151d0fab5eSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD,
8161d0fab5eSBarry Smith           MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction()
817e884886fSBarry Smith @*/
8187087cfbeSBarry Smith PetscErrorCode  MatMFFDSetFunction(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx)
819e884886fSBarry Smith {
820bc0f08ceSBarry Smith   PetscErrorCode ierr;
821e884886fSBarry Smith 
822e884886fSBarry Smith   PetscFunctionBegin;
82388b4c220SStefano Zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
824bc0f08ceSBarry Smith   ierr = PetscTryMethod(mat,"MatMFFDSetFunction_C",(Mat,PetscErrorCode (*)(void*,Vec,Vec),void*),(mat,func,funcctx));CHKERRQ(ierr);
825e884886fSBarry Smith   PetscFunctionReturn(0);
826e884886fSBarry Smith }
827e884886fSBarry Smith 
828e884886fSBarry Smith /*@C
829e884886fSBarry Smith    MatMFFDSetFunctioni - Sets the function for a single component
830e884886fSBarry Smith 
8313f9fe445SBarry Smith    Logically Collective on Mat
832e884886fSBarry Smith 
833e884886fSBarry Smith    Input Parameters:
834e884886fSBarry Smith +  mat - the matrix free matrix created via MatCreateSNESMF()
835e884886fSBarry Smith -  funci - the function to use
836e884886fSBarry Smith 
837e884886fSBarry Smith    Level: advanced
838e884886fSBarry Smith 
839e884886fSBarry Smith    Notes:
840e884886fSBarry Smith     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
84194eb4328SStefano Zampini     matrix inside your compute Jacobian routine.
84294eb4328SStefano Zampini     This function is necessary to compute the diagonal of the matrix.
843486afcceSStefano Zampini     funci must not contain any MPI call as it is called inside a loop on the local portion of the vector.
844e884886fSBarry Smith 
84594eb4328SStefano Zampini .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction(), MatGetDiagonal()
8461d0fab5eSBarry Smith 
847e884886fSBarry Smith @*/
8487087cfbeSBarry Smith PetscErrorCode  MatMFFDSetFunctioni(Mat mat,PetscErrorCode (*funci)(void*,PetscInt,Vec,PetscScalar*))
849e884886fSBarry Smith {
8504ac538c5SBarry Smith   PetscErrorCode ierr;
851e884886fSBarry Smith 
852e884886fSBarry Smith   PetscFunctionBegin;
8530700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
8544ac538c5SBarry Smith   ierr = PetscTryMethod(mat,"MatMFFDSetFunctioni_C",(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)),(mat,funci));CHKERRQ(ierr);
855e884886fSBarry Smith   PetscFunctionReturn(0);
856e884886fSBarry Smith }
857e884886fSBarry Smith 
858e884886fSBarry Smith /*@C
859e884886fSBarry Smith    MatMFFDSetFunctioniBase - Sets the base vector for a single component function evaluation
860e884886fSBarry Smith 
8613f9fe445SBarry Smith    Logically Collective on Mat
862e884886fSBarry Smith 
863e884886fSBarry Smith    Input Parameters:
864e884886fSBarry Smith +  mat - the matrix free matrix created via MatCreateSNESMF()
865e884886fSBarry Smith -  func - the function to use
866e884886fSBarry Smith 
867e884886fSBarry Smith    Level: advanced
868e884886fSBarry Smith 
869e884886fSBarry Smith    Notes:
870e884886fSBarry Smith     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
87194eb4328SStefano Zampini     matrix inside your compute Jacobian routine.
87294eb4328SStefano Zampini     This function is necessary to compute the diagonal of the matrix.
873e884886fSBarry Smith 
874e884886fSBarry Smith 
875e884886fSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD
87694eb4328SStefano Zampini           MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction(), MatGetDiagonal()
877e884886fSBarry Smith @*/
8787087cfbeSBarry Smith PetscErrorCode  MatMFFDSetFunctioniBase(Mat mat,PetscErrorCode (*func)(void*,Vec))
879e884886fSBarry Smith {
8804ac538c5SBarry Smith   PetscErrorCode ierr;
881e884886fSBarry Smith 
882e884886fSBarry Smith   PetscFunctionBegin;
8830700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
8844ac538c5SBarry Smith   ierr = PetscTryMethod(mat,"MatMFFDSetFunctioniBase_C",(Mat,PetscErrorCode (*)(void*,Vec)),(mat,func));CHKERRQ(ierr);
885e884886fSBarry Smith   PetscFunctionReturn(0);
886e884886fSBarry Smith }
887e884886fSBarry Smith 
888e884886fSBarry Smith /*@
889e884886fSBarry Smith    MatMFFDSetPeriod - Sets how often h is recomputed, by default it is everytime
890e884886fSBarry Smith 
8913f9fe445SBarry Smith    Logically Collective on Mat
892e884886fSBarry Smith 
893e884886fSBarry Smith    Input Parameters:
894e884886fSBarry Smith +  mat - the matrix free matrix created via MatCreateSNESMF()
895e884886fSBarry Smith -  period - 1 for everytime, 2 for every second etc
896e884886fSBarry Smith 
897e884886fSBarry Smith    Options Database Keys:
898a2b725a8SWilliam Gropp .  -mat_mffd_period <period>
899e884886fSBarry Smith 
900e884886fSBarry Smith    Level: advanced
901e884886fSBarry Smith 
902e884886fSBarry Smith 
903e884886fSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(),
9041d0fab5eSBarry Smith           MatMFFDSetHHistory(), MatMFFDResetHHistory()
905e884886fSBarry Smith @*/
9067087cfbeSBarry Smith PetscErrorCode  MatMFFDSetPeriod(Mat mat,PetscInt period)
907e884886fSBarry Smith {
908bc0f08ceSBarry Smith   PetscErrorCode ierr;
909e884886fSBarry Smith 
910e884886fSBarry Smith   PetscFunctionBegin;
91188b4c220SStefano Zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
91288b4c220SStefano Zampini   PetscValidLogicalCollectiveInt(mat,period,2);
913bc0f08ceSBarry Smith   ierr = PetscTryMethod(mat,"MatMFFDSetPeriod_C",(Mat,PetscInt),(mat,period));CHKERRQ(ierr);
914e884886fSBarry Smith   PetscFunctionReturn(0);
915e884886fSBarry Smith }
916e884886fSBarry Smith 
917e884886fSBarry Smith /*@
918e884886fSBarry Smith    MatMFFDSetFunctionError - Sets the error_rel for the approximation of
919e884886fSBarry Smith    matrix-vector products using finite differences.
920e884886fSBarry Smith 
9213f9fe445SBarry Smith    Logically Collective on Mat
922e884886fSBarry Smith 
923e884886fSBarry Smith    Input Parameters:
924e884886fSBarry Smith +  mat - the matrix free matrix created via MatCreateMFFD() or MatCreateSNESMF()
925e884886fSBarry Smith -  error_rel - relative error (should be set to the square root of
926e884886fSBarry Smith                the relative error in the function evaluations)
927e884886fSBarry Smith 
928e884886fSBarry Smith    Options Database Keys:
929a2b725a8SWilliam Gropp .  -mat_mffd_err <error_rel> - Sets error_rel
930e884886fSBarry Smith 
931e884886fSBarry Smith    Level: advanced
932e884886fSBarry Smith 
933e884886fSBarry Smith    Notes:
934e884886fSBarry Smith    The default matrix-free matrix-vector product routine computes
935e884886fSBarry Smith .vb
936e884886fSBarry Smith      F'(u)*a = [F(u+h*a) - F(u)]/h where
937e884886fSBarry Smith      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
938e884886fSBarry Smith        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   else
939e884886fSBarry Smith .ve
940e884886fSBarry Smith 
941e884886fSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD
9421d0fab5eSBarry Smith           MatMFFDSetHHistory(), MatMFFDResetHHistory()
943e884886fSBarry Smith @*/
9447087cfbeSBarry Smith PetscErrorCode  MatMFFDSetFunctionError(Mat mat,PetscReal error)
945e884886fSBarry Smith {
946bc0f08ceSBarry Smith   PetscErrorCode ierr;
947e884886fSBarry Smith 
948e884886fSBarry Smith   PetscFunctionBegin;
94988b4c220SStefano Zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
95088b4c220SStefano Zampini   PetscValidLogicalCollectiveReal(mat,error,2);
951bc0f08ceSBarry Smith   ierr = PetscTryMethod(mat,"MatMFFDSetFunctionError_C",(Mat,PetscReal),(mat,error));CHKERRQ(ierr);
952e884886fSBarry Smith   PetscFunctionReturn(0);
953e884886fSBarry Smith }
954e884886fSBarry Smith 
955e884886fSBarry Smith /*@
956e884886fSBarry Smith    MatMFFDSetHHistory - Sets an array to collect a history of the
957e884886fSBarry Smith    differencing values (h) computed for the matrix-free product.
958e884886fSBarry Smith 
9593f9fe445SBarry Smith    Logically Collective on Mat
960e884886fSBarry Smith 
961e884886fSBarry Smith    Input Parameters:
962e884886fSBarry Smith +  J - the matrix-free matrix context
963e884886fSBarry Smith .  histroy - space to hold the history
964e884886fSBarry Smith -  nhistory - number of entries in history, if more entries are generated than
965e884886fSBarry Smith               nhistory, then the later ones are discarded
966e884886fSBarry Smith 
967e884886fSBarry Smith    Level: advanced
968e884886fSBarry Smith 
969e884886fSBarry Smith    Notes:
970e884886fSBarry Smith    Use MatMFFDResetHHistory() to reset the history counter and collect
971e884886fSBarry Smith    a new batch of differencing parameters, h.
972e884886fSBarry Smith 
973e884886fSBarry Smith .seealso: MatMFFDGetH(), MatCreateSNESMF(),
9741d0fab5eSBarry Smith           MatMFFDResetHHistory(), MatMFFDSetFunctionError()
975e884886fSBarry Smith 
976e884886fSBarry Smith @*/
9777087cfbeSBarry Smith PetscErrorCode  MatMFFDSetHHistory(Mat J,PetscScalar history[],PetscInt nhistory)
978e884886fSBarry Smith {
979bc0f08ceSBarry Smith   PetscErrorCode ierr;
980e884886fSBarry Smith 
981e884886fSBarry Smith   PetscFunctionBegin;
98288b4c220SStefano Zampini   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
98388b4c220SStefano Zampini   if (history) PetscValidPointer(history,2);
98488b4c220SStefano Zampini   PetscValidLogicalCollectiveInt(J,nhistory,3);
985789d8953SBarry Smith   ierr = PetscUseMethod(J,"MatMFFDSetHHistory_C",(Mat,PetscScalar[],PetscInt),(J,history,nhistory));CHKERRQ(ierr);
986e884886fSBarry Smith   PetscFunctionReturn(0);
987e884886fSBarry Smith }
988e884886fSBarry Smith 
989e884886fSBarry Smith /*@
990e884886fSBarry Smith    MatMFFDResetHHistory - Resets the counter to zero to begin
991e884886fSBarry Smith    collecting a new set of differencing histories.
992e884886fSBarry Smith 
9933f9fe445SBarry Smith    Logically Collective on Mat
994e884886fSBarry Smith 
995e884886fSBarry Smith    Input Parameters:
996e884886fSBarry Smith .  J - the matrix-free matrix context
997e884886fSBarry Smith 
998e884886fSBarry Smith    Level: advanced
999e884886fSBarry Smith 
1000e884886fSBarry Smith    Notes:
1001e884886fSBarry Smith    Use MatMFFDSetHHistory() to create the original history counter.
1002e884886fSBarry Smith 
1003e884886fSBarry Smith .seealso: MatMFFDGetH(), MatCreateSNESMF(),
10041d0fab5eSBarry Smith           MatMFFDSetHHistory(), MatMFFDSetFunctionError()
1005e884886fSBarry Smith 
1006e884886fSBarry Smith @*/
10077087cfbeSBarry Smith PetscErrorCode  MatMFFDResetHHistory(Mat J)
1008e884886fSBarry Smith {
1009bc0f08ceSBarry Smith   PetscErrorCode ierr;
1010e884886fSBarry Smith 
1011e884886fSBarry Smith   PetscFunctionBegin;
101288b4c220SStefano Zampini   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
1013bc0f08ceSBarry Smith   ierr = PetscTryMethod(J,"MatMFFDResetHHistory_C",(Mat),(J));CHKERRQ(ierr);
1014e884886fSBarry Smith   PetscFunctionReturn(0);
1015e884886fSBarry Smith }
1016e884886fSBarry Smith 
1017487a658cSBarry Smith /*@
1018e884886fSBarry Smith     MatMFFDSetBase - Sets the vector U at which matrix vector products of the
1019e884886fSBarry Smith         Jacobian are computed
1020e884886fSBarry Smith 
10213f9fe445SBarry Smith     Logically Collective on Mat
1022e884886fSBarry Smith 
1023e884886fSBarry Smith     Input Parameters:
1024e884886fSBarry Smith +   J - the MatMFFD matrix
10253ec795f1SBarry Smith .   U - the vector
1026bcddec3dSBarry Smith -   F - (optional) vector that contains F(u) if it has been already computed
1027e884886fSBarry Smith 
102895452b02SPatrick Sanan     Notes:
102995452b02SPatrick Sanan     This is rarely used directly
1030e884886fSBarry Smith 
10318af5ae88SBarry Smith     If F is provided then it is not recomputed. Otherwise the function is evaluated at the base
10328af5ae88SBarry Smith     point during the first MatMult() after each call to MatMFFDSetBase().
1033dff2f722SBarry Smith 
1034e884886fSBarry Smith     Level: advanced
1035e884886fSBarry Smith 
1036e884886fSBarry Smith @*/
10377087cfbeSBarry Smith PetscErrorCode  MatMFFDSetBase(Mat J,Vec U,Vec F)
1038e884886fSBarry Smith {
10394ac538c5SBarry Smith   PetscErrorCode ierr;
1040e884886fSBarry Smith 
1041e884886fSBarry Smith   PetscFunctionBegin;
10420700a824SBarry Smith   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
10430700a824SBarry Smith   PetscValidHeaderSpecific(U,VEC_CLASSID,2);
10440700a824SBarry Smith   if (F) PetscValidHeaderSpecific(F,VEC_CLASSID,3);
10454ac538c5SBarry Smith   ierr = PetscTryMethod(J,"MatMFFDSetBase_C",(Mat,Vec,Vec),(J,U,F));CHKERRQ(ierr);
1046e884886fSBarry Smith   PetscFunctionReturn(0);
1047e884886fSBarry Smith }
1048e884886fSBarry Smith 
1049e884886fSBarry Smith /*@C
1050e884886fSBarry Smith     MatMFFDSetCheckh - Sets a function that checks the computed h and adjusts
1051e884886fSBarry Smith         it to satisfy some criteria
1052e884886fSBarry Smith 
10533f9fe445SBarry Smith     Logically Collective on Mat
1054e884886fSBarry Smith 
1055e884886fSBarry Smith     Input Parameters:
1056e884886fSBarry Smith +   J - the MatMFFD matrix
1057e884886fSBarry Smith .   fun - the function that checks h
1058e884886fSBarry Smith -   ctx - any context needed by the function
1059e884886fSBarry Smith 
1060e884886fSBarry Smith     Options Database Keys:
1061e884886fSBarry Smith .   -mat_mffd_check_positivity
1062e884886fSBarry Smith 
1063e884886fSBarry Smith     Level: advanced
1064e884886fSBarry Smith 
106595452b02SPatrick Sanan     Notes:
106695452b02SPatrick Sanan     For example, MatMFFDCheckPositivity() insures that all entries
1067e884886fSBarry Smith        of U + h*a are non-negative
1068e884886fSBarry Smith 
1069755b7f64SBarry Smith      The function you provide is called after the default h has been computed and allows you to
1070755b7f64SBarry Smith      modify it.
1071755b7f64SBarry Smith 
1072755b7f64SBarry Smith .seealso:  MatMFFDCheckPositivity()
1073e884886fSBarry Smith @*/
10747087cfbeSBarry Smith PetscErrorCode  MatMFFDSetCheckh(Mat J,PetscErrorCode (*fun)(void*,Vec,Vec,PetscScalar*),void *ctx)
1075e884886fSBarry Smith {
10764ac538c5SBarry Smith   PetscErrorCode ierr;
1077e884886fSBarry Smith 
1078e884886fSBarry Smith   PetscFunctionBegin;
10790700a824SBarry Smith   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
10804ac538c5SBarry Smith   ierr = PetscTryMethod(J,"MatMFFDSetCheckh_C",(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*),(J,fun,ctx));CHKERRQ(ierr);
1081e884886fSBarry Smith   PetscFunctionReturn(0);
1082e884886fSBarry Smith }
1083e884886fSBarry Smith 
1084e884886fSBarry Smith /*@
1085e884886fSBarry Smith     MatMFFDCheckPositivity - Checks that all entries in U + h*a are positive or
1086e884886fSBarry Smith         zero, decreases h until this is satisfied.
1087e884886fSBarry Smith 
10883f9fe445SBarry Smith     Logically Collective on Vec
1089e884886fSBarry Smith 
1090e884886fSBarry Smith     Input Parameters:
1091e884886fSBarry Smith +   U - base vector that is added to
1092e884886fSBarry Smith .   a - vector that is added
1093e884886fSBarry Smith .   h - scaling factor on a
1094e884886fSBarry Smith -   dummy - context variable (unused)
1095e884886fSBarry Smith 
1096e884886fSBarry Smith     Options Database Keys:
1097e884886fSBarry Smith .   -mat_mffd_check_positivity
1098e884886fSBarry Smith 
1099e884886fSBarry Smith     Level: advanced
1100e884886fSBarry Smith 
110195452b02SPatrick Sanan     Notes:
110295452b02SPatrick Sanan     This is rarely used directly, rather it is passed as an argument to
1103e884886fSBarry Smith            MatMFFDSetCheckh()
1104e884886fSBarry Smith 
1105e884886fSBarry Smith .seealso:  MatMFFDSetCheckh()
1106e884886fSBarry Smith @*/
11077087cfbeSBarry Smith PetscErrorCode  MatMFFDCheckPositivity(void *dummy,Vec U,Vec a,PetscScalar *h)
1108e884886fSBarry Smith {
1109e884886fSBarry Smith   PetscReal      val, minval;
1110e884886fSBarry Smith   PetscScalar    *u_vec, *a_vec;
1111e884886fSBarry Smith   PetscErrorCode ierr;
1112e884886fSBarry Smith   PetscInt       i,n;
1113e884886fSBarry Smith   MPI_Comm       comm;
1114e884886fSBarry Smith 
1115e884886fSBarry Smith   PetscFunctionBegin;
111688b4c220SStefano Zampini   PetscValidHeaderSpecific(U,VEC_CLASSID,2);
111788b4c220SStefano Zampini   PetscValidHeaderSpecific(a,VEC_CLASSID,3);
111888b4c220SStefano Zampini   PetscValidPointer(h,4);
1119e884886fSBarry Smith   ierr   = PetscObjectGetComm((PetscObject)U,&comm);CHKERRQ(ierr);
1120e884886fSBarry Smith   ierr   = VecGetArray(U,&u_vec);CHKERRQ(ierr);
1121e884886fSBarry Smith   ierr   = VecGetArray(a,&a_vec);CHKERRQ(ierr);
1122e884886fSBarry Smith   ierr   = VecGetLocalSize(U,&n);CHKERRQ(ierr);
1123c068d9bbSLisandro Dalcin   minval = PetscAbsScalar(*h)*PetscRealConstant(1.01);
1124e884886fSBarry Smith   for (i=0; i<n; i++) {
1125e884886fSBarry Smith     if (PetscRealPart(u_vec[i] + *h*a_vec[i]) <= 0.0) {
1126e884886fSBarry Smith       val = PetscAbsScalar(u_vec[i]/a_vec[i]);
1127e884886fSBarry Smith       if (val < minval) minval = val;
1128e884886fSBarry Smith     }
1129e884886fSBarry Smith   }
1130e884886fSBarry Smith   ierr = VecRestoreArray(U,&u_vec);CHKERRQ(ierr);
1131e884886fSBarry Smith   ierr = VecRestoreArray(a,&a_vec);CHKERRQ(ierr);
1132b2566f29SBarry Smith   ierr = MPIU_Allreduce(&minval,&val,1,MPIU_REAL,MPIU_MIN,comm);CHKERRQ(ierr);
1133e884886fSBarry Smith   if (val <= PetscAbsScalar(*h)) {
11346712e2f1SBarry Smith     ierr = PetscInfo2(U,"Scaling back h from %g to %g\n",(double)PetscRealPart(*h),(double)(.99*val));CHKERRQ(ierr);
1135e884886fSBarry Smith     if (PetscRealPart(*h) > 0.0) *h =  0.99*val;
1136e884886fSBarry Smith     else                         *h = -0.99*val;
1137e884886fSBarry Smith   }
1138e884886fSBarry Smith   PetscFunctionReturn(0);
1139e884886fSBarry Smith }
1140