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