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