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