xref: /petsc/src/mat/impls/mffd/mffd.c (revision 88b4c22029ec44fe190bdcfb366f38b2676ca873)
1 
2 #include <petsc/private/matimpl.h>
3 #include <../src/mat/impls/mffd/mffdimpl.h>   /*I  "petscmat.h"   I*/
4 
5 PetscFunctionList MatMFFDList              = 0;
6 PetscBool         MatMFFDRegisterAllCalled = PETSC_FALSE;
7 
8 PetscClassId  MATMFFD_CLASSID;
9 PetscLogEvent MATMFFD_Mult;
10 
11 static PetscBool MatMFFDPackageInitialized = PETSC_FALSE;
12 /*@C
13   MatMFFDFinalizePackage - This function destroys everything in the MatMFFD package. It is
14   called from PetscFinalize().
15 
16   Level: developer
17 
18 .keywords: Petsc, destroy, package
19 .seealso: PetscFinalize(), MatCreateMFFD(), MatCreateSNESMF()
20 @*/
21 PetscErrorCode  MatMFFDFinalizePackage(void)
22 {
23   PetscErrorCode ierr;
24 
25   PetscFunctionBegin;
26   ierr = PetscFunctionListDestroy(&MatMFFDList);CHKERRQ(ierr);
27   MatMFFDPackageInitialized = PETSC_FALSE;
28   MatMFFDRegisterAllCalled  = PETSC_FALSE;
29   PetscFunctionReturn(0);
30 }
31 
32 /*@C
33   MatMFFDInitializePackage - This function initializes everything in the MatMFFD package. It is called
34   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to MatCreate_MFFD()
35   when using static libraries.
36 
37   Level: developer
38 
39 .keywords: Vec, initialize, package
40 .seealso: PetscInitialize()
41 @*/
42 PetscErrorCode  MatMFFDInitializePackage(void)
43 {
44   char           logList[256];
45   char           *className;
46   PetscBool      opt;
47   PetscErrorCode ierr;
48 
49   PetscFunctionBegin;
50   if (MatMFFDPackageInitialized) PetscFunctionReturn(0);
51   MatMFFDPackageInitialized = PETSC_TRUE;
52   /* Register Classes */
53   ierr = PetscClassIdRegister("MatMFFD",&MATMFFD_CLASSID);CHKERRQ(ierr);
54   /* Register Constructors */
55   ierr = MatMFFDRegisterAll();CHKERRQ(ierr);
56   /* Register Events */
57   ierr = PetscLogEventRegister("MatMult MF",          MATMFFD_CLASSID,&MATMFFD_Mult);CHKERRQ(ierr);
58 
59   /* Process info exclusions */
60   ierr = PetscOptionsGetString(NULL,NULL, "-info_exclude", logList, 256, &opt);CHKERRQ(ierr);
61   if (opt) {
62     ierr = PetscStrstr(logList, "matmffd", &className);CHKERRQ(ierr);
63     if (className) {
64       ierr = PetscInfoDeactivateClass(MATMFFD_CLASSID);CHKERRQ(ierr);
65     }
66   }
67   /* Process summary exclusions */
68   ierr = PetscOptionsGetString(NULL,NULL, "-log_exclude", logList, 256, &opt);CHKERRQ(ierr);
69   if (opt) {
70     ierr = PetscStrstr(logList, "matmffd", &className);CHKERRQ(ierr);
71     if (className) {
72       ierr = PetscLogEventDeactivateClass(MATMFFD_CLASSID);CHKERRQ(ierr);
73     }
74   }
75   ierr = PetscRegisterFinalize(MatMFFDFinalizePackage);CHKERRQ(ierr);
76   PetscFunctionReturn(0);
77 }
78 
79 /*@C
80     MatMFFDSetType - Sets the method that is used to compute the
81     differencing parameter for finite differene matrix-free formulations.
82 
83     Input Parameters:
84 +   mat - the "matrix-free" matrix created via MatCreateSNESMF(), or MatCreateMFFD()
85           or MatSetType(mat,MATMFFD);
86 -   ftype - the type requested, either MATMFFD_WP or MATMFFD_DS
87 
88     Level: advanced
89 
90     Notes:
91     For example, such routines can compute h for use in
92     Jacobian-vector products of the form
93 
94                         F(x+ha) - F(x)
95           F'(u)a  ~=  ----------------
96                               h
97 
98 .seealso: MatCreateSNESMF(), MatMFFDRegister(), MatMFFDSetFunction(), MatCreateMFFD()
99 @*/
100 PetscErrorCode  MatMFFDSetType(Mat mat,MatMFFDType ftype)
101 {
102   PetscErrorCode ierr,(*r)(MatMFFD);
103   MatMFFD        ctx = (MatMFFD)mat->data;
104   PetscBool      match;
105 
106   PetscFunctionBegin;
107   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
108   PetscValidCharPointer(ftype,2);
109 
110   ierr = PetscObjectTypeCompare((PetscObject)mat,MATMFFD,&match);CHKERRQ(ierr);
111   if (!match) PetscFunctionReturn(0);
112 
113   /* already set, so just return */
114   ierr = PetscObjectTypeCompare((PetscObject)ctx,ftype,&match);CHKERRQ(ierr);
115   if (match) PetscFunctionReturn(0);
116 
117   /* destroy the old one if it exists */
118   if (ctx->ops->destroy) {
119     ierr = (*ctx->ops->destroy)(ctx);CHKERRQ(ierr);
120   }
121 
122   ierr =  PetscFunctionListFind(MatMFFDList,ftype,&r);CHKERRQ(ierr);
123   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown MatMFFD type %s given",ftype);
124   ierr = (*r)(ctx);CHKERRQ(ierr);
125   ierr = PetscObjectChangeTypeName((PetscObject)ctx,ftype);CHKERRQ(ierr);
126   PetscFunctionReturn(0);
127 }
128 
129 typedef PetscErrorCode (*FCN1)(void*,Vec); /* force argument to next function to not be extern C*/
130 static PetscErrorCode  MatMFFDSetFunctioniBase_MFFD(Mat mat,FCN1 func)
131 {
132   MatMFFD ctx = (MatMFFD)mat->data;
133 
134   PetscFunctionBegin;
135   ctx->funcisetbase = func;
136   PetscFunctionReturn(0);
137 }
138 
139 typedef PetscErrorCode (*FCN2)(void*,PetscInt,Vec,PetscScalar*); /* force argument to next function to not be extern C*/
140 static PetscErrorCode  MatMFFDSetFunctioni_MFFD(Mat mat,FCN2 funci)
141 {
142   MatMFFD ctx = (MatMFFD)mat->data;
143 
144   PetscFunctionBegin;
145   ctx->funci = funci;
146   PetscFunctionReturn(0);
147 }
148 
149 static PetscErrorCode  MatMFFDResetHHistory_MFFD(Mat J)
150 {
151   MatMFFD ctx = (MatMFFD)J->data;
152 
153   PetscFunctionBegin;
154   ctx->ncurrenth = 0;
155   PetscFunctionReturn(0);
156 }
157 
158 /*@C
159    MatMFFDRegister - Adds a method to the MatMFFD registry.
160 
161    Not Collective
162 
163    Input Parameters:
164 +  name_solver - name of a new user-defined compute-h module
165 -  routine_create - routine to create method context
166 
167    Level: developer
168 
169    Notes:
170    MatMFFDRegister() may be called multiple times to add several user-defined solvers.
171 
172    Sample usage:
173 .vb
174    MatMFFDRegister("my_h",MyHCreate);
175 .ve
176 
177    Then, your solver can be chosen with the procedural interface via
178 $     MatMFFDSetType(mfctx,"my_h")
179    or at runtime via the option
180 $     -mat_mffd_type my_h
181 
182 .keywords: MatMFFD, register
183 
184 .seealso: MatMFFDRegisterAll(), MatMFFDRegisterDestroy()
185  @*/
186 PetscErrorCode  MatMFFDRegister(const char sname[],PetscErrorCode (*function)(MatMFFD))
187 {
188   PetscErrorCode ierr;
189 
190   PetscFunctionBegin;
191   ierr = PetscFunctionListAdd(&MatMFFDList,sname,function);CHKERRQ(ierr);
192   PetscFunctionReturn(0);
193 }
194 
195 /* ----------------------------------------------------------------------------------------*/
196 static PetscErrorCode MatDestroy_MFFD(Mat mat)
197 {
198   PetscErrorCode ierr;
199   MatMFFD        ctx = (MatMFFD)mat->data;
200 
201   PetscFunctionBegin;
202   ierr = VecDestroy(&ctx->w);CHKERRQ(ierr);
203   ierr = VecDestroy(&ctx->drscale);CHKERRQ(ierr);
204   ierr = VecDestroy(&ctx->dlscale);CHKERRQ(ierr);
205   ierr = VecDestroy(&ctx->dshift);CHKERRQ(ierr);
206   if (ctx->current_f_allocated) {
207     ierr = VecDestroy(&ctx->current_f);CHKERRQ(ierr);
208   }
209   if (ctx->ops->destroy) {ierr = (*ctx->ops->destroy)(ctx);CHKERRQ(ierr);}
210   ierr      = PetscHeaderDestroy(&ctx);CHKERRQ(ierr);
211   mat->data = 0;
212 
213   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetBase_C",NULL);CHKERRQ(ierr);
214   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctioniBase_C",NULL);CHKERRQ(ierr);
215   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctioni_C",NULL);CHKERRQ(ierr);
216   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunction_C",NULL);CHKERRQ(ierr);
217   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctionError_C",NULL);CHKERRQ(ierr);
218   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetCheckh_C",NULL);CHKERRQ(ierr);
219   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetPeriod_C",NULL);CHKERRQ(ierr);
220   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDResetHHistory_C",NULL);CHKERRQ(ierr);
221   PetscFunctionReturn(0);
222 }
223 
224 static PetscErrorCode MatSetUp_MFFD(Mat mat)
225 {
226   PetscErrorCode ierr;
227 
228   PetscFunctionBegin;
229   ierr = PetscLayoutSetUp(mat->rmap);CHKERRQ(ierr);
230   ierr = PetscLayoutSetUp(mat->cmap);CHKERRQ(ierr);
231   PetscFunctionReturn(0);
232 }
233 
234 /*
235    MatMFFDView_MFFD - Views matrix-free parameters.
236 
237 */
238 static PetscErrorCode MatView_MFFD(Mat J,PetscViewer viewer)
239 {
240   PetscErrorCode ierr;
241   MatMFFD        ctx = (MatMFFD)J->data;
242   PetscBool      iascii, viewbase, viewfunction;
243   const char     *prefix;
244 
245   PetscFunctionBegin;
246   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
247   if (iascii) {
248     ierr = PetscViewerASCIIPrintf(viewer,"Matrix-free approximation:\n");CHKERRQ(ierr);
249     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
250     ierr = PetscViewerASCIIPrintf(viewer,"err=%g (relative error in function evaluation)\n",(double)ctx->error_rel);CHKERRQ(ierr);
251     if (!((PetscObject)ctx)->type_name) {
252       ierr = PetscViewerASCIIPrintf(viewer,"The compute h routine has not yet been set\n");CHKERRQ(ierr);
253     } else {
254       ierr = PetscViewerASCIIPrintf(viewer,"Using %s compute h routine\n",((PetscObject)ctx)->type_name);CHKERRQ(ierr);
255     }
256     if (ctx->ops->view) {
257       ierr = (*ctx->ops->view)(ctx,viewer);CHKERRQ(ierr);
258     }
259     ierr = PetscObjectGetOptionsPrefix((PetscObject)J, &prefix);CHKERRQ(ierr);
260 
261     ierr = PetscOptionsHasName(((PetscObject)J)->options,prefix, "-mat_mffd_view_base", &viewbase);CHKERRQ(ierr);
262     if (viewbase) {
263       ierr = PetscViewerASCIIPrintf(viewer, "Base:\n");CHKERRQ(ierr);
264       ierr = VecView(ctx->current_u, viewer);CHKERRQ(ierr);
265     }
266     ierr = PetscOptionsHasName(((PetscObject)J)->options,prefix, "-mat_mffd_view_function", &viewfunction);CHKERRQ(ierr);
267     if (viewfunction) {
268       ierr = PetscViewerASCIIPrintf(viewer, "Function:\n");CHKERRQ(ierr);
269       ierr = VecView(ctx->current_f, viewer);CHKERRQ(ierr);
270     }
271     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
272   }
273   PetscFunctionReturn(0);
274 }
275 
276 /*
277    MatAssemblyEnd_MFFD - Resets the ctx->ncurrenth to zero. This
278    allows the user to indicate the beginning of a new linear solve by calling
279    MatAssemblyXXX() on the matrix free matrix. This then allows the
280    MatCreateMFFD_WP() to properly compute ||U|| only the first time
281    in the linear solver rather than every time.
282 
283    This function is referenced directly from MatAssemblyEnd_SNESMF(), which may be in a different shared library hence
284    it must be labeled as PETSC_EXTERN
285 */
286 PETSC_EXTERN PetscErrorCode MatAssemblyEnd_MFFD(Mat J,MatAssemblyType mt)
287 {
288   PetscErrorCode ierr;
289   MatMFFD        j = (MatMFFD)J->data;
290 
291   PetscFunctionBegin;
292   ierr      = MatMFFDResetHHistory(J);CHKERRQ(ierr);
293   j->vshift = 0.0;
294   j->vscale = 1.0;
295   PetscFunctionReturn(0);
296 }
297 
298 /*
299   MatMult_MFFD - Default matrix-free form for Jacobian-vector product, y = F'(u)*a:
300 
301         y ~= (F(u + ha) - F(u))/h,
302   where F = nonlinear function, as set by SNESSetFunction()
303         u = current iterate
304         h = difference interval
305 */
306 static PetscErrorCode MatMult_MFFD(Mat mat,Vec a,Vec y)
307 {
308   MatMFFD        ctx = (MatMFFD)mat->data;
309   PetscScalar    h;
310   Vec            w,U,F;
311   PetscErrorCode ierr;
312   PetscBool      zeroa;
313 
314   PetscFunctionBegin;
315   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");
316   /* We log matrix-free matrix-vector products separately, so that we can
317      separate the performance monitoring from the cases that use conventional
318      storage.  We may eventually modify event logging to associate events
319      with particular objects, hence alleviating the more general problem. */
320   ierr = PetscLogEventBegin(MATMFFD_Mult,a,y,0,0);CHKERRQ(ierr);
321 
322   w = ctx->w;
323   U = ctx->current_u;
324   F = ctx->current_f;
325   /*
326       Compute differencing parameter
327   */
328   if (!((PetscObject)ctx)->type_name) {
329     ierr = MatMFFDSetType(mat,MATMFFD_WP);CHKERRQ(ierr);
330     ierr = MatSetFromOptions(mat);CHKERRQ(ierr);
331   }
332   ierr = (*ctx->ops->compute)(ctx,U,a,&h,&zeroa);CHKERRQ(ierr);
333   if (zeroa) {
334     ierr = VecSet(y,0.0);CHKERRQ(ierr);
335     PetscFunctionReturn(0);
336   }
337 
338   if (mat->erroriffailure && PetscIsInfOrNanScalar(h)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Computed Nan differencing parameter h");
339   if (ctx->checkh) {
340     ierr = (*ctx->checkh)(ctx->checkhctx,U,a,&h);CHKERRQ(ierr);
341   }
342 
343   /* keep a record of the current differencing parameter h */
344   ctx->currenth = h;
345 #if defined(PETSC_USE_COMPLEX)
346   ierr = PetscInfo2(mat,"Current differencing parameter: %g + %g i\n",(double)PetscRealPart(h),(double)PetscImaginaryPart(h));CHKERRQ(ierr);
347 #else
348   ierr = PetscInfo1(mat,"Current differencing parameter: %15.12e\n",h);CHKERRQ(ierr);
349 #endif
350   if (ctx->historyh && ctx->ncurrenth < ctx->maxcurrenth) {
351     ctx->historyh[ctx->ncurrenth] = h;
352   }
353   ctx->ncurrenth++;
354 
355   /* w = u + ha */
356   if (ctx->drscale) {
357     ierr = VecPointwiseMult(ctx->drscale,a,U);CHKERRQ(ierr);
358     ierr = VecAYPX(U,h,w);CHKERRQ(ierr);
359   } else {
360     ierr = VecWAXPY(w,h,a,U);CHKERRQ(ierr);
361   }
362 
363   /* compute func(U) as base for differencing; only needed first time in and not when provided by user */
364   if (ctx->ncurrenth == 1 && ctx->current_f_allocated) {
365     ierr = (*ctx->func)(ctx->funcctx,U,F);CHKERRQ(ierr);
366   }
367   ierr = (*ctx->func)(ctx->funcctx,w,y);CHKERRQ(ierr);
368 
369   ierr = VecAXPY(y,-1.0,F);CHKERRQ(ierr);
370   ierr = VecScale(y,1.0/h);CHKERRQ(ierr);
371 
372   if ((ctx->vshift != 0.0) || (ctx->vscale != 1.0)) {
373     ierr = VecAXPBY(y,ctx->vshift,ctx->vscale,a);CHKERRQ(ierr);
374   }
375   if (ctx->dlscale) {
376     ierr = VecPointwiseMult(y,ctx->dlscale,y);CHKERRQ(ierr);
377   }
378   if (ctx->dshift) {
379     ierr = VecPointwiseMult(ctx->dshift,a,U);CHKERRQ(ierr);
380     ierr = VecAXPY(y,1.0,U);CHKERRQ(ierr);
381   }
382 
383   if (mat->nullsp) {ierr = MatNullSpaceRemove(mat->nullsp,y);CHKERRQ(ierr);}
384 
385   ierr = PetscLogEventEnd(MATMFFD_Mult,a,y,0,0);CHKERRQ(ierr);
386   PetscFunctionReturn(0);
387 }
388 
389 /*
390   MatGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix
391 
392         y ~= (F(u + ha) - F(u))/h,
393   where F = nonlinear function, as set by SNESSetFunction()
394         u = current iterate
395         h = difference interval
396 */
397 static PetscErrorCode MatGetDiagonal_MFFD(Mat mat,Vec a)
398 {
399   MatMFFD        ctx = (MatMFFD)mat->data;
400   PetscScalar    h,*aa,*ww,v;
401   PetscReal      epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON;
402   Vec            w,U;
403   PetscErrorCode ierr;
404   PetscInt       i,rstart,rend;
405 
406   PetscFunctionBegin;
407   if (!ctx->funci) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Requires calling MatMFFDSetFunctioni() first");
408 
409   w    = ctx->w;
410   U    = ctx->current_u;
411   ierr = (*ctx->func)(ctx->funcctx,U,a);CHKERRQ(ierr);
412   ierr = (*ctx->funcisetbase)(ctx->funcctx,U);CHKERRQ(ierr);
413   ierr = VecCopy(U,w);CHKERRQ(ierr);
414 
415   ierr = VecGetOwnershipRange(a,&rstart,&rend);CHKERRQ(ierr);
416   ierr = VecGetArray(a,&aa);CHKERRQ(ierr);
417   for (i=rstart; i<rend; i++) {
418     ierr = VecGetArray(w,&ww);CHKERRQ(ierr);
419     h    = ww[i-rstart];
420     if (h == 0.0) h = 1.0;
421     if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0)     h = umin;
422     else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin;
423     h *= epsilon;
424 
425     ww[i-rstart] += h;
426     ierr          = VecRestoreArray(w,&ww);CHKERRQ(ierr);
427     ierr          = (*ctx->funci)(ctx->funcctx,i,w,&v);CHKERRQ(ierr);
428     aa[i-rstart]  = (v - aa[i-rstart])/h;
429 
430     /* possibly shift and scale result */
431     if ((ctx->vshift != 0.0) || (ctx->vscale != 1.0)) {
432       aa[i - rstart] = ctx->vshift + ctx->vscale*aa[i-rstart];
433     }
434 
435     ierr          = VecGetArray(w,&ww);CHKERRQ(ierr);
436     ww[i-rstart] -= h;
437     ierr          = VecRestoreArray(w,&ww);CHKERRQ(ierr);
438   }
439   ierr = VecRestoreArray(a,&aa);CHKERRQ(ierr);
440   PetscFunctionReturn(0);
441 }
442 
443 static PetscErrorCode MatDiagonalScale_MFFD(Mat mat,Vec ll,Vec rr)
444 {
445   MatMFFD        aij = (MatMFFD)mat->data;
446   PetscErrorCode ierr;
447 
448   PetscFunctionBegin;
449   if (ll && !aij->dlscale) {
450     ierr = VecDuplicate(ll,&aij->dlscale);CHKERRQ(ierr);
451   }
452   if (rr && !aij->drscale) {
453     ierr = VecDuplicate(rr,&aij->drscale);CHKERRQ(ierr);
454   }
455   if (ll) {
456     ierr = VecCopy(ll,aij->dlscale);CHKERRQ(ierr);
457   }
458   if (rr) {
459     ierr = VecCopy(rr,aij->drscale);CHKERRQ(ierr);
460   }
461   PetscFunctionReturn(0);
462 }
463 
464 static PetscErrorCode MatDiagonalSet_MFFD(Mat mat,Vec ll,InsertMode mode)
465 {
466   MatMFFD        aij = (MatMFFD)mat->data;
467   PetscErrorCode ierr;
468 
469   PetscFunctionBegin;
470   if (mode == INSERT_VALUES) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"No diagonal set with INSERT_VALUES");
471   if (!aij->dshift) {
472     ierr = VecDuplicate(ll,&aij->dshift);CHKERRQ(ierr);
473   }
474   ierr = VecAXPY(aij->dshift,1.0,ll);CHKERRQ(ierr);
475   PetscFunctionReturn(0);
476 }
477 
478 static PetscErrorCode MatShift_MFFD(Mat Y,PetscScalar a)
479 {
480   MatMFFD shell = (MatMFFD)Y->data;
481 
482   PetscFunctionBegin;
483   shell->vshift += a;
484   PetscFunctionReturn(0);
485 }
486 
487 static PetscErrorCode MatScale_MFFD(Mat Y,PetscScalar a)
488 {
489   MatMFFD shell = (MatMFFD)Y->data;
490 
491   PetscFunctionBegin;
492   shell->vscale *= a;
493   PetscFunctionReturn(0);
494 }
495 
496 PETSC_EXTERN PetscErrorCode MatMFFDSetBase_MFFD(Mat J,Vec U,Vec F)
497 {
498   PetscErrorCode ierr;
499   MatMFFD        ctx = (MatMFFD)J->data;
500 
501   PetscFunctionBegin;
502   ierr = MatMFFDResetHHistory(J);CHKERRQ(ierr);
503 
504   ctx->current_u = U;
505   if (F) {
506     if (ctx->current_f_allocated) {ierr = VecDestroy(&ctx->current_f);CHKERRQ(ierr);}
507     ctx->current_f           = F;
508     ctx->current_f_allocated = PETSC_FALSE;
509   } else if (!ctx->current_f_allocated) {
510     ierr = MatCreateVecs(J,NULL,&ctx->current_f);CHKERRQ(ierr);
511 
512     ctx->current_f_allocated = PETSC_TRUE;
513   }
514   if (!ctx->w) {
515     ierr = VecDuplicate(ctx->current_u, &ctx->w);CHKERRQ(ierr);
516   }
517   J->assembled = PETSC_TRUE;
518   PetscFunctionReturn(0);
519 }
520 
521 typedef PetscErrorCode (*FCN3)(void*,Vec,Vec,PetscScalar*); /* force argument to next function to not be extern C*/
522 
523 static PetscErrorCode  MatMFFDSetCheckh_MFFD(Mat J,FCN3 fun,void *ectx)
524 {
525   MatMFFD ctx = (MatMFFD)J->data;
526 
527   PetscFunctionBegin;
528   ctx->checkh    = fun;
529   ctx->checkhctx = ectx;
530   PetscFunctionReturn(0);
531 }
532 
533 /*@C
534    MatMFFDSetOptionsPrefix - Sets the prefix used for searching for all
535    MatMFFD options in the database.
536 
537    Collective on Mat
538 
539    Input Parameter:
540 +  A - the Mat context
541 -  prefix - the prefix to prepend to all option names
542 
543    Notes:
544    A hyphen (-) must NOT be given at the beginning of the prefix name.
545    The first character of all runtime options is AUTOMATICALLY the hyphen.
546 
547    Level: advanced
548 
549 .keywords: SNES, matrix-free, parameters
550 
551 .seealso: MatSetFromOptions(), MatCreateSNESMF(), MatCreateMFFD()
552 @*/
553 PetscErrorCode  MatMFFDSetOptionsPrefix(Mat mat,const char prefix[])
554 
555 {
556   MatMFFD        mfctx = mat ? (MatMFFD)mat->data : (MatMFFD)NULL;
557   PetscErrorCode ierr;
558 
559   PetscFunctionBegin;
560   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
561   PetscValidHeaderSpecific(mfctx,MATMFFD_CLASSID,1);
562   ierr = PetscObjectSetOptionsPrefix((PetscObject)mfctx,prefix);CHKERRQ(ierr);
563   PetscFunctionReturn(0);
564 }
565 
566 static PetscErrorCode  MatSetFromOptions_MFFD(PetscOptionItems *PetscOptionsObject,Mat mat)
567 {
568   MatMFFD        mfctx = (MatMFFD)mat->data;
569   PetscErrorCode ierr;
570   PetscBool      flg;
571   char           ftype[256];
572 
573   PetscFunctionBegin;
574   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
575   PetscValidHeaderSpecific(mfctx,MATMFFD_CLASSID,1);
576   ierr = PetscObjectOptionsBegin((PetscObject)mfctx);CHKERRQ(ierr);
577   ierr = PetscOptionsFList("-mat_mffd_type","Matrix free type","MatMFFDSetType",MatMFFDList,((PetscObject)mfctx)->type_name,ftype,256,&flg);CHKERRQ(ierr);
578   if (flg) {
579     ierr = MatMFFDSetType(mat,ftype);CHKERRQ(ierr);
580   }
581 
582   ierr = PetscOptionsReal("-mat_mffd_err","set sqrt relative error in function","MatMFFDSetFunctionError",mfctx->error_rel,&mfctx->error_rel,0);CHKERRQ(ierr);
583   ierr = PetscOptionsInt("-mat_mffd_period","how often h is recomputed","MatMFFDSetPeriod",mfctx->recomputeperiod,&mfctx->recomputeperiod,0);CHKERRQ(ierr);
584 
585   flg  = PETSC_FALSE;
586   ierr = PetscOptionsBool("-mat_mffd_check_positivity","Insure that U + h*a is nonnegative","MatMFFDSetCheckh",flg,&flg,NULL);CHKERRQ(ierr);
587   if (flg) {
588     ierr = MatMFFDSetCheckh(mat,MatMFFDCheckPositivity,0);CHKERRQ(ierr);
589   }
590   if (mfctx->ops->setfromoptions) {
591     ierr = (*mfctx->ops->setfromoptions)(PetscOptionsObject,mfctx);CHKERRQ(ierr);
592   }
593   ierr = PetscOptionsEnd();CHKERRQ(ierr);
594   PetscFunctionReturn(0);
595 }
596 
597 static PetscErrorCode  MatMFFDSetPeriod_MFFD(Mat mat,PetscInt period)
598 {
599   MatMFFD ctx = (MatMFFD)mat->data;
600 
601   PetscFunctionBegin;
602   ctx->recomputeperiod = period;
603   PetscFunctionReturn(0);
604 }
605 
606 static PetscErrorCode  MatMFFDSetFunction_MFFD(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx)
607 {
608   MatMFFD ctx = (MatMFFD)mat->data;
609 
610   PetscFunctionBegin;
611   ctx->func    = func;
612   ctx->funcctx = funcctx;
613   PetscFunctionReturn(0);
614 }
615 
616 static PetscErrorCode  MatMFFDSetFunctionError_MFFD(Mat mat,PetscReal error)
617 {
618   MatMFFD ctx = (MatMFFD)mat->data;
619 
620   PetscFunctionBegin;
621   if (error != PETSC_DEFAULT) ctx->error_rel = error;
622   PetscFunctionReturn(0);
623 }
624 
625 static PetscErrorCode MatMissingDiagonal_MFFD(Mat A,PetscBool  *missing,PetscInt *d)
626 {
627   PetscFunctionBegin;
628   *missing = PETSC_FALSE;
629   PetscFunctionReturn(0);
630 }
631 
632 /*MC
633   MATMFFD - MATMFFD = "mffd" - A matrix free matrix type.
634 
635   Level: advanced
636 
637 .seealso: MatCreateMFFD(), MatCreateSNESMF(), MatMFFDSetFunction(), MatMFFDSetType(),
638           MatMFFDSetFunctionError(), MatMFFDDSSetUmin(), MatMFFDSetFunction()
639           MatMFFDSetHHistory(), MatMFFDResetHHistory(), MatCreateSNESMF(),
640           MatMFFDGetH(),
641 M*/
642 PETSC_EXTERN PetscErrorCode MatCreate_MFFD(Mat A)
643 {
644   MatMFFD        mfctx;
645   PetscErrorCode ierr;
646 
647   PetscFunctionBegin;
648   ierr = MatMFFDInitializePackage();CHKERRQ(ierr);
649 
650   ierr = PetscHeaderCreate(mfctx,MATMFFD_CLASSID,"MatMFFD","Matrix-free Finite Differencing","Mat",PetscObjectComm((PetscObject)A),MatDestroy_MFFD,MatView_MFFD);CHKERRQ(ierr);
651 
652   mfctx->error_rel                = PETSC_SQRT_MACHINE_EPSILON;
653   mfctx->recomputeperiod          = 1;
654   mfctx->count                    = 0;
655   mfctx->currenth                 = 0.0;
656   mfctx->historyh                 = NULL;
657   mfctx->ncurrenth                = 0;
658   mfctx->maxcurrenth              = 0;
659   ((PetscObject)mfctx)->type_name = 0;
660 
661   mfctx->vshift = 0.0;
662   mfctx->vscale = 1.0;
663 
664   /*
665      Create the empty data structure to contain compute-h routines.
666      These will be filled in below from the command line options or
667      a later call with MatMFFDSetType() or if that is not called
668      then it will default in the first use of MatMult_MFFD()
669   */
670   mfctx->ops->compute        = 0;
671   mfctx->ops->destroy        = 0;
672   mfctx->ops->view           = 0;
673   mfctx->ops->setfromoptions = 0;
674   mfctx->hctx                = 0;
675 
676   mfctx->func    = 0;
677   mfctx->funcctx = 0;
678   mfctx->w       = NULL;
679 
680   A->data = mfctx;
681 
682   A->ops->mult            = MatMult_MFFD;
683   A->ops->destroy         = MatDestroy_MFFD;
684   A->ops->setup           = MatSetUp_MFFD;
685   A->ops->view            = MatView_MFFD;
686   A->ops->assemblyend     = MatAssemblyEnd_MFFD;
687   A->ops->getdiagonal     = MatGetDiagonal_MFFD;
688   A->ops->scale           = MatScale_MFFD;
689   A->ops->shift           = MatShift_MFFD;
690   A->ops->diagonalscale   = MatDiagonalScale_MFFD;
691   A->ops->diagonalset     = MatDiagonalSet_MFFD;
692   A->ops->setfromoptions  = MatSetFromOptions_MFFD;
693   A->ops->missingdiagonal = MatMissingDiagonal_MFFD;
694   A->assembled            = PETSC_TRUE;
695 
696   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetBase_C",MatMFFDSetBase_MFFD);CHKERRQ(ierr);
697   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctioniBase_C",MatMFFDSetFunctioniBase_MFFD);CHKERRQ(ierr);
698   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctioni_C",MatMFFDSetFunctioni_MFFD);CHKERRQ(ierr);
699   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunction_C",MatMFFDSetFunction_MFFD);CHKERRQ(ierr);
700   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetCheckh_C",MatMFFDSetCheckh_MFFD);CHKERRQ(ierr);
701   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetPeriod_C",MatMFFDSetPeriod_MFFD);CHKERRQ(ierr);
702   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctionError_C",MatMFFDSetFunctionError_MFFD);CHKERRQ(ierr);
703   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDResetHHistory_C",MatMFFDResetHHistory_MFFD);CHKERRQ(ierr);
704 
705   mfctx->mat = A;
706 
707   ierr = PetscObjectChangeTypeName((PetscObject)A,MATMFFD);CHKERRQ(ierr);
708   PetscFunctionReturn(0);
709 }
710 
711 /*@
712    MatCreateMFFD - Creates a matrix-free matrix. See also MatCreateSNESMF()
713 
714    Collective on Vec
715 
716    Input Parameters:
717 +  comm - MPI communicator
718 .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
719            This value should be the same as the local size used in creating the
720            y vector for the matrix-vector product y = Ax.
721 .  n - This value should be the same as the local size used in creating the
722        x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have
723        calculated if N is given) For square matrices n is almost always m.
724 .  M - number of global rows (or PETSC_DETERMINE to have calculated if m is given)
725 -  N - number of global columns (or PETSC_DETERMINE to have calculated if n is given)
726 
727 
728    Output Parameter:
729 .  J - the matrix-free matrix
730 
731    Options Database Keys: call MatSetFromOptions() to trigger these
732 +  -mat_mffd_type - wp or ds (see MATMFFD_WP or MATMFFD_DS)
733 -  -mat_mffd_err - square root of estimated relative error in function evaluation
734 -  -mat_mffd_period - how often h is recomputed, defaults to 1, everytime
735 
736 
737    Level: advanced
738 
739    Notes:
740    The matrix-free matrix context merely contains the function pointers
741    and work space for performing finite difference approximations of
742    Jacobian-vector products, F'(u)*a,
743 
744    The default code uses the following approach to compute h
745 
746 .vb
747      F'(u)*a = [F(u+h*a) - F(u)]/h where
748      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
749        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
750  where
751      error_rel = square root of relative error in function evaluation
752      umin = minimum iterate parameter
753 .ve
754 
755    You can call SNESSetJacobian() with MatMFFDComputeJacobian() if you are using matrix and not a different
756    preconditioner matrix
757 
758    The user can set the error_rel via MatMFFDSetFunctionError() and
759    umin via MatMFFDDSSetUmin(); see Users-Manual: ch_snes for details.
760 
761    The user should call MatDestroy() when finished with the matrix-free
762    matrix context.
763 
764    Options Database Keys:
765 +  -mat_mffd_err <error_rel> - Sets error_rel
766 .  -mat_mffd_unim <umin> - Sets umin (for default PETSc routine that computes h only)
767 -  -mat_mffd_check_positivity
768 
769 .keywords: default, matrix-free, create, matrix
770 
771 .seealso: MatDestroy(), MatMFFDSetFunctionError(), MatMFFDDSSetUmin(), MatMFFDSetFunction()
772           MatMFFDSetHHistory(), MatMFFDResetHHistory(), MatCreateSNESMF(),
773           MatMFFDGetH(), MatMFFDRegister(), MatMFFDComputeJacobian()
774 
775 @*/
776 PetscErrorCode  MatCreateMFFD(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *J)
777 {
778   PetscErrorCode ierr;
779 
780   PetscFunctionBegin;
781   ierr = MatCreate(comm,J);CHKERRQ(ierr);
782   ierr = MatSetSizes(*J,m,n,M,N);CHKERRQ(ierr);
783   ierr = MatSetType(*J,MATMFFD);CHKERRQ(ierr);
784   ierr = MatSetUp(*J);CHKERRQ(ierr);
785   PetscFunctionReturn(0);
786 }
787 
788 /*@
789    MatMFFDGetH - Gets the last value that was used as the differencing
790    parameter.
791 
792    Not Collective
793 
794    Input Parameters:
795 .  mat - the matrix obtained with MatCreateSNESMF()
796 
797    Output Paramter:
798 .  h - the differencing step size
799 
800    Level: advanced
801 
802 .keywords: SNES, matrix-free, parameters
803 
804 .seealso: MatCreateSNESMF(),MatMFFDSetHHistory(), MatCreateMFFD(), MATMFFD, MatMFFDResetHHistory()
805 @*/
806 PetscErrorCode  MatMFFDGetH(Mat mat,PetscScalar *h)
807 {
808   MatMFFD        ctx = (MatMFFD)mat->data;
809   PetscErrorCode ierr;
810   PetscBool      match;
811 
812   PetscFunctionBegin;
813   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
814   PetscValidPointer(h,2);
815   ierr = PetscObjectTypeCompare((PetscObject)mat,MATMFFD,&match);CHKERRQ(ierr);
816   if (!match) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONG,"Not a MFFD matrix");
817 
818   *h = ctx->currenth;
819   PetscFunctionReturn(0);
820 }
821 
822 /*@C
823    MatMFFDSetFunction - Sets the function used in applying the matrix free.
824 
825    Logically Collective on Mat
826 
827    Input Parameters:
828 +  mat - the matrix free matrix created via MatCreateSNESMF() or MatCreateMFFD()
829 .  func - the function to use
830 -  funcctx - optional function context passed to function
831 
832    Calling Sequence of func:
833 $     func (void *funcctx, Vec x, Vec f)
834 
835 +  funcctx - user provided context
836 .  x - input vector
837 -  f - computed output function
838 
839    Level: advanced
840 
841    Notes:
842     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
843     matrix inside your compute Jacobian routine
844 
845     If this is not set then it will use the function set with SNESSetFunction() if MatCreateSNESMF() was used.
846 
847 .keywords: SNES, matrix-free, function
848 
849 .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD,
850           MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction()
851 @*/
852 PetscErrorCode  MatMFFDSetFunction(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx)
853 {
854   PetscErrorCode ierr;
855 
856   PetscFunctionBegin;
857   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
858   ierr = PetscTryMethod(mat,"MatMFFDSetFunction_C",(Mat,PetscErrorCode (*)(void*,Vec,Vec),void*),(mat,func,funcctx));CHKERRQ(ierr);
859   PetscFunctionReturn(0);
860 }
861 
862 /*@C
863    MatMFFDSetFunctioni - Sets the function for a single component
864 
865    Logically Collective on Mat
866 
867    Input Parameters:
868 +  mat - the matrix free matrix created via MatCreateSNESMF()
869 -  funci - the function to use
870 
871    Level: advanced
872 
873    Notes:
874     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
875     matrix inside your compute Jacobian routine
876 
877 
878 .keywords: SNES, matrix-free, function
879 
880 .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction()
881 
882 @*/
883 PetscErrorCode  MatMFFDSetFunctioni(Mat mat,PetscErrorCode (*funci)(void*,PetscInt,Vec,PetscScalar*))
884 {
885   PetscErrorCode ierr;
886 
887   PetscFunctionBegin;
888   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
889   ierr = PetscTryMethod(mat,"MatMFFDSetFunctioni_C",(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)),(mat,funci));CHKERRQ(ierr);
890   PetscFunctionReturn(0);
891 }
892 
893 /*@C
894    MatMFFDSetFunctioniBase - Sets the base vector for a single component function evaluation
895 
896    Logically Collective on Mat
897 
898    Input Parameters:
899 +  mat - the matrix free matrix created via MatCreateSNESMF()
900 -  func - the function to use
901 
902    Level: advanced
903 
904    Notes:
905     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
906     matrix inside your compute Jacobian routine
907 
908 
909 .keywords: SNES, matrix-free, function
910 
911 .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD
912           MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction()
913 @*/
914 PetscErrorCode  MatMFFDSetFunctioniBase(Mat mat,PetscErrorCode (*func)(void*,Vec))
915 {
916   PetscErrorCode ierr;
917 
918   PetscFunctionBegin;
919   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
920   ierr = PetscTryMethod(mat,"MatMFFDSetFunctioniBase_C",(Mat,PetscErrorCode (*)(void*,Vec)),(mat,func));CHKERRQ(ierr);
921   PetscFunctionReturn(0);
922 }
923 
924 /*@
925    MatMFFDSetPeriod - Sets how often h is recomputed, by default it is everytime
926 
927    Logically Collective on Mat
928 
929    Input Parameters:
930 +  mat - the matrix free matrix created via MatCreateSNESMF()
931 -  period - 1 for everytime, 2 for every second etc
932 
933    Options Database Keys:
934 +  -mat_mffd_period <period>
935 
936    Level: advanced
937 
938 
939 .keywords: SNES, matrix-free, parameters
940 
941 .seealso: MatCreateSNESMF(),MatMFFDGetH(),
942           MatMFFDSetHHistory(), MatMFFDResetHHistory()
943 @*/
944 PetscErrorCode  MatMFFDSetPeriod(Mat mat,PetscInt period)
945 {
946   PetscErrorCode ierr;
947 
948   PetscFunctionBegin;
949   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
950   PetscValidLogicalCollectiveInt(mat,period,2);
951   ierr = PetscTryMethod(mat,"MatMFFDSetPeriod_C",(Mat,PetscInt),(mat,period));CHKERRQ(ierr);
952   PetscFunctionReturn(0);
953 }
954 
955 /*@
956    MatMFFDSetFunctionError - Sets the error_rel for the approximation of
957    matrix-vector products using finite differences.
958 
959    Logically Collective on Mat
960 
961    Input Parameters:
962 +  mat - the matrix free matrix created via MatCreateMFFD() or MatCreateSNESMF()
963 -  error_rel - relative error (should be set to the square root of
964                the relative error in the function evaluations)
965 
966    Options Database Keys:
967 +  -mat_mffd_err <error_rel> - Sets error_rel
968 
969    Level: advanced
970 
971    Notes:
972    The default matrix-free matrix-vector product routine computes
973 .vb
974      F'(u)*a = [F(u+h*a) - F(u)]/h where
975      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
976        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   else
977 .ve
978 
979 .keywords: SNES, matrix-free, parameters
980 
981 .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD
982           MatMFFDSetHHistory(), MatMFFDResetHHistory()
983 @*/
984 PetscErrorCode  MatMFFDSetFunctionError(Mat mat,PetscReal error)
985 {
986   PetscErrorCode ierr;
987 
988   PetscFunctionBegin;
989   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
990   PetscValidLogicalCollectiveReal(mat,error,2);
991   ierr = PetscTryMethod(mat,"MatMFFDSetFunctionError_C",(Mat,PetscReal),(mat,error));CHKERRQ(ierr);
992   PetscFunctionReturn(0);
993 }
994 
995 /*@
996    MatMFFDSetHHistory - Sets an array to collect a history of the
997    differencing values (h) computed for the matrix-free product.
998 
999    Logically Collective on Mat
1000 
1001    Input Parameters:
1002 +  J - the matrix-free matrix context
1003 .  histroy - space to hold the history
1004 -  nhistory - number of entries in history, if more entries are generated than
1005               nhistory, then the later ones are discarded
1006 
1007    Level: advanced
1008 
1009    Notes:
1010    Use MatMFFDResetHHistory() to reset the history counter and collect
1011    a new batch of differencing parameters, h.
1012 
1013 .keywords: SNES, matrix-free, h history, differencing history
1014 
1015 .seealso: MatMFFDGetH(), MatCreateSNESMF(),
1016           MatMFFDResetHHistory(), MatMFFDSetFunctionError()
1017 
1018 @*/
1019 PetscErrorCode  MatMFFDSetHHistory(Mat J,PetscScalar history[],PetscInt nhistory)
1020 {
1021   MatMFFD        ctx = (MatMFFD)J->data;
1022   PetscErrorCode ierr;
1023   PetscBool      match;
1024 
1025   PetscFunctionBegin;
1026   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
1027   if (history) PetscValidPointer(history,2);
1028   PetscValidLogicalCollectiveInt(J,nhistory,3);
1029   ierr = PetscObjectTypeCompare((PetscObject)J,MATMFFD,&match);CHKERRQ(ierr);
1030   if (!match) SETERRQ(PetscObjectComm((PetscObject)J),PETSC_ERR_ARG_WRONG,"Not a MFFD matrix");
1031   ctx->historyh    = history;
1032   ctx->maxcurrenth = nhistory;
1033   ctx->currenth    = 0.;
1034   PetscFunctionReturn(0);
1035 }
1036 
1037 /*@
1038    MatMFFDResetHHistory - Resets the counter to zero to begin
1039    collecting a new set of differencing histories.
1040 
1041    Logically Collective on Mat
1042 
1043    Input Parameters:
1044 .  J - the matrix-free matrix context
1045 
1046    Level: advanced
1047 
1048    Notes:
1049    Use MatMFFDSetHHistory() to create the original history counter.
1050 
1051 .keywords: SNES, matrix-free, h history, differencing history
1052 
1053 .seealso: MatMFFDGetH(), MatCreateSNESMF(),
1054           MatMFFDSetHHistory(), MatMFFDSetFunctionError()
1055 
1056 @*/
1057 PetscErrorCode  MatMFFDResetHHistory(Mat J)
1058 {
1059   PetscErrorCode ierr;
1060 
1061   PetscFunctionBegin;
1062   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
1063   ierr = PetscTryMethod(J,"MatMFFDResetHHistory_C",(Mat),(J));CHKERRQ(ierr);
1064   PetscFunctionReturn(0);
1065 }
1066 
1067 /*@
1068     MatMFFDSetBase - Sets the vector U at which matrix vector products of the
1069         Jacobian are computed
1070 
1071     Logically Collective on Mat
1072 
1073     Input Parameters:
1074 +   J - the MatMFFD matrix
1075 .   U - the vector
1076 -   F - (optional) vector that contains F(u) if it has been already computed
1077 
1078     Notes: This is rarely used directly
1079 
1080     If F is provided then it is not recomputed. Otherwise the function is evaluated at the base
1081     point during the first MatMult() after each call to MatMFFDSetBase().
1082 
1083     Level: advanced
1084 
1085 @*/
1086 PetscErrorCode  MatMFFDSetBase(Mat J,Vec U,Vec F)
1087 {
1088   PetscErrorCode ierr;
1089 
1090   PetscFunctionBegin;
1091   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
1092   PetscValidHeaderSpecific(U,VEC_CLASSID,2);
1093   if (F) PetscValidHeaderSpecific(F,VEC_CLASSID,3);
1094   ierr = PetscTryMethod(J,"MatMFFDSetBase_C",(Mat,Vec,Vec),(J,U,F));CHKERRQ(ierr);
1095   PetscFunctionReturn(0);
1096 }
1097 
1098 /*@C
1099     MatMFFDSetCheckh - Sets a function that checks the computed h and adjusts
1100         it to satisfy some criteria
1101 
1102     Logically Collective on Mat
1103 
1104     Input Parameters:
1105 +   J - the MatMFFD matrix
1106 .   fun - the function that checks h
1107 -   ctx - any context needed by the function
1108 
1109     Options Database Keys:
1110 .   -mat_mffd_check_positivity
1111 
1112     Level: advanced
1113 
1114     Notes: For example, MatMFFDCheckPositivity() insures that all entries
1115        of U + h*a are non-negative
1116 
1117      The function you provide is called after the default h has been computed and allows you to
1118      modify it.
1119 
1120 .seealso:  MatMFFDCheckPositivity()
1121 @*/
1122 PetscErrorCode  MatMFFDSetCheckh(Mat J,PetscErrorCode (*fun)(void*,Vec,Vec,PetscScalar*),void *ctx)
1123 {
1124   PetscErrorCode ierr;
1125 
1126   PetscFunctionBegin;
1127   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
1128   ierr = PetscTryMethod(J,"MatMFFDSetCheckh_C",(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*),(J,fun,ctx));CHKERRQ(ierr);
1129   PetscFunctionReturn(0);
1130 }
1131 
1132 /*@
1133     MatMFFDCheckPositivity - Checks that all entries in U + h*a are positive or
1134         zero, decreases h until this is satisfied.
1135 
1136     Logically Collective on Vec
1137 
1138     Input Parameters:
1139 +   U - base vector that is added to
1140 .   a - vector that is added
1141 .   h - scaling factor on a
1142 -   dummy - context variable (unused)
1143 
1144     Options Database Keys:
1145 .   -mat_mffd_check_positivity
1146 
1147     Level: advanced
1148 
1149     Notes: This is rarely used directly, rather it is passed as an argument to
1150            MatMFFDSetCheckh()
1151 
1152 .seealso:  MatMFFDSetCheckh()
1153 @*/
1154 PetscErrorCode  MatMFFDCheckPositivity(void *dummy,Vec U,Vec a,PetscScalar *h)
1155 {
1156   PetscReal      val, minval;
1157   PetscScalar    *u_vec, *a_vec;
1158   PetscErrorCode ierr;
1159   PetscInt       i,n;
1160   MPI_Comm       comm;
1161 
1162   PetscFunctionBegin;
1163   PetscValidHeaderSpecific(U,VEC_CLASSID,2);
1164   PetscValidHeaderSpecific(a,VEC_CLASSID,3);
1165   PetscValidPointer(h,4);
1166   ierr   = PetscObjectGetComm((PetscObject)U,&comm);CHKERRQ(ierr);
1167   ierr   = VecGetArray(U,&u_vec);CHKERRQ(ierr);
1168   ierr   = VecGetArray(a,&a_vec);CHKERRQ(ierr);
1169   ierr   = VecGetLocalSize(U,&n);CHKERRQ(ierr);
1170   minval = PetscAbsScalar(*h*1.01);
1171   for (i=0; i<n; i++) {
1172     if (PetscRealPart(u_vec[i] + *h*a_vec[i]) <= 0.0) {
1173       val = PetscAbsScalar(u_vec[i]/a_vec[i]);
1174       if (val < minval) minval = val;
1175     }
1176   }
1177   ierr = VecRestoreArray(U,&u_vec);CHKERRQ(ierr);
1178   ierr = VecRestoreArray(a,&a_vec);CHKERRQ(ierr);
1179   ierr = MPIU_Allreduce(&minval,&val,1,MPIU_REAL,MPIU_MIN,comm);CHKERRQ(ierr);
1180   if (val <= PetscAbsScalar(*h)) {
1181     ierr = PetscInfo2(U,"Scaling back h from %g to %g\n",(double)PetscRealPart(*h),(double)(.99*val));CHKERRQ(ierr);
1182     if (PetscRealPart(*h) > 0.0) *h =  0.99*val;
1183     else                         *h = -0.99*val;
1184   }
1185   PetscFunctionReturn(0);
1186 }
1187