xref: /petsc/src/mat/impls/mffd/mffd.c (revision fd705b320d8d44969be9ca25a36dbdd35fbe8e12)
1 #define PETSCMAT_DLL
2 
3 #include "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 = 0;
10 PetscEvent  MATMFFD_Mult = 0;
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 = PetscLogClassRegister(&MATMFFD_COOKIE,     "MatMFFD");CHKERRQ(ierr);
40   /* Register Constructors */
41   ierr = MatMFFDRegisterAll(path);CHKERRQ(ierr);
42   /* Register Events */
43   ierr = PetscLogEventRegister(&MATMFFD_Mult, "MatMult MF",          MATMFFD_COOKIE);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,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 */
321   if (ctx->ncurrenth == 1) {
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   ierr = PetscOptionsName("-mat_mffd_check_positivity","Insure that U + h*a is nonnegative","MatMFFDSetCheckh",&flg);CHKERRQ(ierr);
500   if (flg) {
501     ierr = MatMFFDSetCheckh(mat,MatMFFDCheckPositivity,0);CHKERRQ(ierr);
502   }
503   if (mfctx->ops->setfromoptions) {
504     ierr = (*mfctx->ops->setfromoptions)(mfctx);CHKERRQ(ierr);
505   }
506   ierr = PetscOptionsEnd();CHKERRQ(ierr);
507   PetscFunctionReturn(0);
508 }
509 
510 /*MC
511   MATMFFD - MATMFFD = "mffd" - A matrix free matrix type.
512 
513   Level: advanced
514 
515 .seealso: MatCreateMFFD(), MatCreateSNESMF()
516 M*/
517 EXTERN_C_BEGIN
518 #undef __FUNCT__
519 #define __FUNCT__ "MatCreate_MFFD"
520 PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_MFFD(Mat A)
521 {
522   MatMFFD         mfctx;
523   PetscErrorCode  ierr;
524 
525   PetscFunctionBegin;
526 #ifndef PETSC_USE_DYNAMIC_LIBRARIES
527   ierr = MatMFFDInitializePackage(PETSC_NULL);CHKERRQ(ierr);
528 #endif
529 
530   ierr = PetscHeaderCreate(mfctx,_p_MatMFFD,struct _MFOps,MATMFFD_COOKIE,0,"MatMFFD",((PetscObject)A)->comm,MatDestroy_MFFD,MatView_MFFD);CHKERRQ(ierr);
531   mfctx->sp              = 0;
532   mfctx->error_rel       = PETSC_SQRT_MACHINE_EPSILON;
533   mfctx->recomputeperiod = 1;
534   mfctx->count           = 0;
535   mfctx->currenth        = 0.0;
536   mfctx->historyh        = PETSC_NULL;
537   mfctx->ncurrenth       = 0;
538   mfctx->maxcurrenth     = 0;
539   ((PetscObject)mfctx)->type_name       = 0;
540 
541   mfctx->vshift          = 0.0;
542   mfctx->vscale          = 1.0;
543 
544   /*
545      Create the empty data structure to contain compute-h routines.
546      These will be filled in below from the command line options or
547      a later call with MatMFFDSetType() or if that is not called
548      then it will default in the first use of MatMult_MFFD()
549   */
550   mfctx->ops->compute        = 0;
551   mfctx->ops->destroy        = 0;
552   mfctx->ops->view           = 0;
553   mfctx->ops->setfromoptions = 0;
554   mfctx->hctx                = 0;
555 
556   mfctx->func                = 0;
557   mfctx->funcctx             = 0;
558   mfctx->w                   = PETSC_NULL;
559 
560   A->data                = mfctx;
561 
562   A->ops->mult           = MatMult_MFFD;
563   A->ops->destroy        = MatDestroy_MFFD;
564   A->ops->view           = MatView_MFFD;
565   A->ops->assemblyend    = MatAssemblyEnd_MFFD;
566   A->ops->getdiagonal    = MatGetDiagonal_MFFD;
567   A->ops->scale          = MatScale_MFFD;
568   A->ops->shift          = MatShift_MFFD;
569   A->ops->setfromoptions = MatMFFDSetFromOptions;
570   A->assembled = PETSC_TRUE;
571 
572   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDSetBase_C","MatMFFDSetBase_MFFD",MatMFFDSetBase_MFFD);CHKERRQ(ierr);
573   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDSetFunctioniBase_C","MatMFFDSetFunctioniBase_MFFD",MatMFFDSetFunctioniBase_MFFD);CHKERRQ(ierr);
574   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDSetFunctioni_C","MatMFFDSetFunctioni_MFFD",MatMFFDSetFunctioni_MFFD);CHKERRQ(ierr);
575   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDSetCheckh_C","MatMFFDSetCheckh_MFFD",MatMFFDSetCheckh_MFFD);CHKERRQ(ierr);
576   mfctx->mat = A;
577   ierr = PetscObjectChangeTypeName((PetscObject)A,MATMFFD);CHKERRQ(ierr);
578   PetscFunctionReturn(0);
579 }
580 EXTERN_C_END
581 
582 #undef __FUNCT__
583 #define __FUNCT__ "MatCreateMFFD"
584 /*@
585    MatCreateMFFD - Creates a matrix-free matrix. See also MatCreateSNESMF()
586 
587    Collective on Vec
588 
589    Input Parameters:
590 +  comm - MPI communicator
591 .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
592            This value should be the same as the local size used in creating the
593            y vector for the matrix-vector product y = Ax.
594 .  n - This value should be the same as the local size used in creating the
595        x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have
596        calculated if N is given) For square matrices n is almost always m.
597 .  M - number of global rows (or PETSC_DETERMINE to have calculated if m is given)
598 -  N - number of global columns (or PETSC_DETERMINE to have calculated if n is given)
599 
600 
601    Output Parameter:
602 .  J - the matrix-free matrix
603 
604    Level: advanced
605 
606    Notes:
607    The matrix-free matrix context merely contains the function pointers
608    and work space for performing finite difference approximations of
609    Jacobian-vector products, F'(u)*a,
610 
611    The default code uses the following approach to compute h
612 
613 .vb
614      F'(u)*a = [F(u+h*a) - F(u)]/h where
615      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
616        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
617  where
618      error_rel = square root of relative error in function evaluation
619      umin = minimum iterate parameter
620 .ve
621 
622    The user can set the error_rel via MatMFFDSetFunctionError() and
623    umin via MatMFFDDefaultSetUmin(); see the nonlinear solvers chapter
624    of the users manual for details.
625 
626    The user should call MatDestroy() when finished with the matrix-free
627    matrix context.
628 
629    Options Database Keys:
630 +  -mat_mffd_err <error_rel> - Sets error_rel
631 .  -mat_mffd_unim <umin> - Sets umin (for default PETSc routine that computes h only)
632 .  -mat_mffd_ksp_monitor - KSP monitor routine that prints differencing h
633 -  -mat_mffd_check_positivity
634 
635 .keywords: default, matrix-free, create, matrix
636 
637 .seealso: MatDestroy(), MatMFFDSetFunctionError(), MatMFFDDefaultSetUmin(), MatMFFDSetFunction()
638           MatMFFDSetHHistory(), MatMFFDResetHHistory(), MatCreateSNESMF(),
639           MatMFFDGetH(),MatMFFDKSPMonitor(), MatMFFDRegisterDynamic),, MatMFFDComputeJacobian()
640 
641 @*/
642 PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMFFD(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *J)
643 {
644   PetscErrorCode ierr;
645 
646   PetscFunctionBegin;
647   ierr = MatCreate(comm,J);CHKERRQ(ierr);
648   ierr = MatSetSizes(*J,m,n,M,N);CHKERRQ(ierr);
649   ierr = MatSetType(*J,MATMFFD);CHKERRQ(ierr);
650   PetscFunctionReturn(0);
651 }
652 
653 
654 #undef __FUNCT__
655 #define __FUNCT__ "MatMFFDGetH"
656 /*@
657    MatMFFDGetH - Gets the last value that was used as the differencing
658    parameter.
659 
660    Not Collective
661 
662    Input Parameters:
663 .  mat - the matrix obtained with MatCreateSNESMF()
664 
665    Output Paramter:
666 .  h - the differencing step size
667 
668    Level: advanced
669 
670 .keywords: SNES, matrix-free, parameters
671 
672 .seealso: MatCreateSNESMF(),MatMFFDSetHHistory(), MatCreateMFFD(), MATMFFD
673           MatMFFDResetHHistory(),MatMFFDKSPMonitor()
674 @*/
675 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDGetH(Mat mat,PetscScalar *h)
676 {
677   MatMFFD ctx = (MatMFFD)mat->data;
678 
679   PetscFunctionBegin;
680   *h = ctx->currenth;
681   PetscFunctionReturn(0);
682 }
683 
684 
685 #undef __FUNCT__
686 #define __FUNCT__ "MatMFFDSetFunction"
687 /*@C
688    MatMFFDSetFunction - Sets the function used in applying the matrix free.
689 
690    Collective on Mat
691 
692    Input Parameters:
693 +  mat - the matrix free matrix created via MatCreateSNESMF()
694 .  func - the function to use
695 -  funcctx - optional function context passed to function
696 
697    Level: advanced
698 
699    Notes:
700     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
701     matrix inside your compute Jacobian routine
702 
703     If this is not set then it will use the function set with SNESSetFunction() if MatCreateSNESMF() was used.
704 
705 .keywords: SNES, matrix-free, function
706 
707 .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD
708           MatMFFDSetHHistory(), MatMFFDResetHHistory(),
709           MatMFFDKSPMonitor(), SNESetFunction()
710 @*/
711 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunction(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx)
712 {
713   MatMFFD ctx = (MatMFFD)mat->data;
714 
715   PetscFunctionBegin;
716   ctx->func    = func;
717   ctx->funcctx = funcctx;
718   PetscFunctionReturn(0);
719 }
720 
721 #undef __FUNCT__
722 #define __FUNCT__ "MatMFFDSetFunctioni"
723 /*@C
724    MatMFFDSetFunctioni - Sets the function for a single component
725 
726    Collective on Mat
727 
728    Input Parameters:
729 +  mat - the matrix free matrix created via MatCreateSNESMF()
730 -  funci - the function to use
731 
732    Level: advanced
733 
734    Notes:
735     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
736     matrix inside your compute Jacobian routine
737 
738 
739 .keywords: SNES, matrix-free, function
740 
741 .seealso: MatCreateSNESMF(),MatMFFDGetH(),
742           MatMFFDSetHHistory(), MatMFFDResetHHistory(),
743           MatMFFDKSPMonitor(), SNESetFunction()
744 @*/
745 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctioni(Mat mat,PetscErrorCode (*funci)(void*,PetscInt,Vec,PetscScalar*))
746 {
747   PetscErrorCode ierr,(*f)(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*));
748 
749   PetscFunctionBegin;
750   PetscValidHeaderSpecific(mat,MAT_COOKIE,1);
751   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatMFFDSetFunctioni_C",(void (**)(void))&f);CHKERRQ(ierr);
752   if (f) {
753     ierr = (*f)(mat,funci);CHKERRQ(ierr);
754   }
755   PetscFunctionReturn(0);
756 }
757 
758 
759 #undef __FUNCT__
760 #define __FUNCT__ "MatMFFDSetFunctioniBase"
761 /*@C
762    MatMFFDSetFunctioniBase - Sets the base vector for a single component function evaluation
763 
764    Collective on Mat
765 
766    Input Parameters:
767 +  mat - the matrix free matrix created via MatCreateSNESMF()
768 -  func - the function to use
769 
770    Level: advanced
771 
772    Notes:
773     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
774     matrix inside your compute Jacobian routine
775 
776 
777 .keywords: SNES, matrix-free, function
778 
779 .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD
780           MatMFFDSetHHistory(), MatMFFDResetHHistory(),
781           MatMFFDKSPMonitor(), SNESetFunction()
782 @*/
783 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctioniBase(Mat mat,PetscErrorCode (*func)(void*,Vec))
784 {
785   PetscErrorCode ierr,(*f)(Mat,PetscErrorCode (*)(void*,Vec));
786 
787   PetscFunctionBegin;
788   PetscValidHeaderSpecific(mat,MAT_COOKIE,1);
789   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatMFFDSetFunctioniBase_C",(void (**)(void))&f);CHKERRQ(ierr);
790   if (f) {
791     ierr = (*f)(mat,func);CHKERRQ(ierr);
792   }
793   PetscFunctionReturn(0);
794 }
795 
796 
797 #undef __FUNCT__
798 #define __FUNCT__ "MatMFFDSetPeriod"
799 /*@
800    MatMFFDSetPeriod - Sets how often h is recomputed, by default it is everytime
801 
802    Collective on Mat
803 
804    Input Parameters:
805 +  mat - the matrix free matrix created via MatCreateSNESMF()
806 -  period - 1 for everytime, 2 for every second etc
807 
808    Options Database Keys:
809 +  -mat_mffd_period <period>
810 
811    Level: advanced
812 
813 
814 .keywords: SNES, matrix-free, parameters
815 
816 .seealso: MatCreateSNESMF(),MatMFFDGetH(),
817           MatMFFDSetHHistory(), MatMFFDResetHHistory(),
818           MatMFFDKSPMonitor()
819 @*/
820 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetPeriod(Mat mat,PetscInt period)
821 {
822   MatMFFD ctx = (MatMFFD)mat->data;
823 
824   PetscFunctionBegin;
825   ctx->recomputeperiod = period;
826   PetscFunctionReturn(0);
827 }
828 
829 #undef __FUNCT__
830 #define __FUNCT__ "MatMFFDSetFunctionError"
831 /*@
832    MatMFFDSetFunctionError - Sets the error_rel for the approximation of
833    matrix-vector products using finite differences.
834 
835    Collective on Mat
836 
837    Input Parameters:
838 +  mat - the matrix free matrix created via MatCreateMFFD() or MatCreateSNESMF()
839 -  error_rel - relative error (should be set to the square root of
840                the relative error in the function evaluations)
841 
842    Options Database Keys:
843 +  -mat_mffd_err <error_rel> - Sets error_rel
844 
845    Level: advanced
846 
847    Notes:
848    The default matrix-free matrix-vector product routine computes
849 .vb
850      F'(u)*a = [F(u+h*a) - F(u)]/h where
851      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
852        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   else
853 .ve
854 
855 .keywords: SNES, matrix-free, parameters
856 
857 .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD
858           MatMFFDSetHHistory(), MatMFFDResetHHistory(),
859           MatMFFDKSPMonitor()
860 @*/
861 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctionError(Mat mat,PetscReal error)
862 {
863   MatMFFD ctx = (MatMFFD)mat->data;
864 
865   PetscFunctionBegin;
866   if (error != PETSC_DEFAULT) ctx->error_rel = error;
867   PetscFunctionReturn(0);
868 }
869 
870 #undef __FUNCT__
871 #define __FUNCT__ "MatMFFDAddNullSpace"
872 /*@
873    MatMFFDAddNullSpace - Provides a null space that an operator is
874    supposed to have.  Since roundoff will create a small component in
875    the null space, if you know the null space you may have it
876    automatically removed.
877 
878    Collective on Mat
879 
880    Input Parameters:
881 +  J - the matrix-free matrix context
882 -  nullsp - object created with MatNullSpaceCreate()
883 
884    Level: advanced
885 
886 .keywords: SNES, matrix-free, null space
887 
888 .seealso: MatNullSpaceCreate(), MatMFFDGetH(), MatCreateSNESMF(), MatCreateMFFD(), MATMFFD
889           MatMFFDSetHHistory(), MatMFFDResetHHistory(),
890           MatMFFDKSPMonitor(), MatMFFDErrorRel()
891 @*/
892 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDAddNullSpace(Mat J,MatNullSpace nullsp)
893 {
894   PetscErrorCode ierr;
895   MatMFFD      ctx = (MatMFFD)J->data;
896 
897   PetscFunctionBegin;
898   ierr = PetscObjectReference((PetscObject)nullsp);CHKERRQ(ierr);
899   if (ctx->sp) { ierr = MatNullSpaceDestroy(ctx->sp);CHKERRQ(ierr); }
900   ctx->sp = nullsp;
901   PetscFunctionReturn(0);
902 }
903 
904 #undef __FUNCT__
905 #define __FUNCT__ "MatMFFDSetHHistory"
906 /*@
907    MatMFFDSetHHistory - Sets an array to collect a history of the
908    differencing values (h) computed for the matrix-free product.
909 
910    Collective on Mat
911 
912    Input Parameters:
913 +  J - the matrix-free matrix context
914 .  histroy - space to hold the history
915 -  nhistory - number of entries in history, if more entries are generated than
916               nhistory, then the later ones are discarded
917 
918    Level: advanced
919 
920    Notes:
921    Use MatMFFDResetHHistory() to reset the history counter and collect
922    a new batch of differencing parameters, h.
923 
924 .keywords: SNES, matrix-free, h history, differencing history
925 
926 .seealso: MatMFFDGetH(), MatCreateSNESMF(),
927           MatMFFDResetHHistory(),
928           MatMFFDKSPMonitor(), MatMFFDSetFunctionError()
929 
930 @*/
931 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetHHistory(Mat J,PetscScalar history[],PetscInt nhistory)
932 {
933   MatMFFD ctx = (MatMFFD)J->data;
934 
935   PetscFunctionBegin;
936   ctx->historyh    = history;
937   ctx->maxcurrenth = nhistory;
938   ctx->currenth    = 0;
939   PetscFunctionReturn(0);
940 }
941 
942 #undef __FUNCT__
943 #define __FUNCT__ "MatMFFDResetHHistory"
944 /*@
945    MatMFFDResetHHistory - Resets the counter to zero to begin
946    collecting a new set of differencing histories.
947 
948    Collective on Mat
949 
950    Input Parameters:
951 .  J - the matrix-free matrix context
952 
953    Level: advanced
954 
955    Notes:
956    Use MatMFFDSetHHistory() to create the original history counter.
957 
958 .keywords: SNES, matrix-free, h history, differencing history
959 
960 .seealso: MatMFFDGetH(), MatCreateSNESMF(),
961           MatMFFDSetHHistory(),
962           MatMFFDKSPMonitor(), MatMFFDSetFunctionError()
963 
964 @*/
965 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDResetHHistory(Mat J)
966 {
967   MatMFFD ctx = (MatMFFD)J->data;
968 
969   PetscFunctionBegin;
970   ctx->ncurrenth    = 0;
971   PetscFunctionReturn(0);
972 }
973 
974 
975 #undef __FUNCT__
976 #define __FUNCT__ "MatMFFDSetBase"
977 /*@
978     MatMFFDSetBase - Sets the vector U at which matrix vector products of the
979         Jacobian are computed
980 
981     Collective on Mat
982 
983     Input Parameters:
984 +   J - the MatMFFD matrix
985 .   U - the vector
986 -   F - vector that contains F(u) if it has been already computed
987 
988     Notes: This is rarely used directly
989 
990     Warning: for a given Mat, this must be called either ALWAYS with an F or never
991 
992     Level: advanced
993 
994 @*/
995 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetBase(Mat J,Vec U,Vec F)
996 {
997   PetscErrorCode ierr,(*f)(Mat,Vec,Vec);
998 
999   PetscFunctionBegin;
1000   PetscValidHeaderSpecific(J,MAT_COOKIE,1);
1001   PetscValidHeaderSpecific(U,VEC_COOKIE,2);
1002   if (F) PetscValidHeaderSpecific(F,VEC_COOKIE,3);
1003   ierr = PetscObjectQueryFunction((PetscObject)J,"MatMFFDSetBase_C",(void (**)(void))&f);CHKERRQ(ierr);
1004   if (f) {
1005     ierr = (*f)(J,U,F);CHKERRQ(ierr);
1006   }
1007   PetscFunctionReturn(0);
1008 }
1009 
1010 #undef __FUNCT__
1011 #define __FUNCT__ "MatMFFDSetCheckh"
1012 /*@C
1013     MatMFFDSetCheckh - Sets a function that checks the computed h and adjusts
1014         it to satisfy some criteria
1015 
1016     Collective on Mat
1017 
1018     Input Parameters:
1019 +   J - the MatMFFD matrix
1020 .   fun - the function that checks h
1021 -   ctx - any context needed by the function
1022 
1023     Options Database Keys:
1024 .   -mat_mffd_check_positivity
1025 
1026     Level: advanced
1027 
1028     Notes: For example, MatMFFDSetCheckPositivity() insures that all entries
1029        of U + h*a are non-negative
1030 
1031 .seealso:  MatMFFDSetCheckPositivity()
1032 @*/
1033 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetCheckh(Mat J,PetscErrorCode (*fun)(void*,Vec,Vec,PetscScalar*),void* ctx)
1034 {
1035   PetscErrorCode ierr,(*f)(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*);
1036 
1037   PetscFunctionBegin;
1038   PetscValidHeaderSpecific(J,MAT_COOKIE,1);
1039   ierr = PetscObjectQueryFunction((PetscObject)J,"MatMFFDSetCheckh_C",(void (**)(void))&f);CHKERRQ(ierr);
1040   if (f) {
1041     ierr = (*f)(J,fun,ctx);CHKERRQ(ierr);
1042   }
1043   PetscFunctionReturn(0);
1044 }
1045 
1046 #undef __FUNCT__
1047 #define __FUNCT__ "MatMFFDSetCheckPositivity"
1048 /*@
1049     MatMFFDCheckPositivity - Checks that all entries in U + h*a are positive or
1050         zero, decreases h until this is satisfied.
1051 
1052     Collective on Vec
1053 
1054     Input Parameters:
1055 +   U - base vector that is added to
1056 .   a - vector that is added
1057 .   h - scaling factor on a
1058 -   dummy - context variable (unused)
1059 
1060     Options Database Keys:
1061 .   -mat_mffd_check_positivity
1062 
1063     Level: advanced
1064 
1065     Notes: This is rarely used directly, rather it is passed as an argument to
1066            MatMFFDSetCheckh()
1067 
1068 .seealso:  MatMFFDSetCheckh()
1069 @*/
1070 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDCheckPositivity(void* dummy,Vec U,Vec a,PetscScalar *h)
1071 {
1072   PetscReal      val, minval;
1073   PetscScalar    *u_vec, *a_vec;
1074   PetscErrorCode ierr;
1075   PetscInt       i,n;
1076   MPI_Comm       comm;
1077 
1078   PetscFunctionBegin;
1079   ierr = PetscObjectGetComm((PetscObject)U,&comm);CHKERRQ(ierr);
1080   ierr = VecGetArray(U,&u_vec);CHKERRQ(ierr);
1081   ierr = VecGetArray(a,&a_vec);CHKERRQ(ierr);
1082   ierr = VecGetLocalSize(U,&n);CHKERRQ(ierr);
1083   minval = PetscAbsScalar(*h*1.01);
1084   for(i=0;i<n;i++) {
1085     if (PetscRealPart(u_vec[i] + *h*a_vec[i]) <= 0.0) {
1086       val = PetscAbsScalar(u_vec[i]/a_vec[i]);
1087       if (val < minval) minval = val;
1088     }
1089   }
1090   ierr = VecRestoreArray(U,&u_vec);CHKERRQ(ierr);
1091   ierr = VecRestoreArray(a,&a_vec);CHKERRQ(ierr);
1092   ierr = PetscGlobalMin(&minval,&val,comm);CHKERRQ(ierr);
1093   if (val <= PetscAbsScalar(*h)) {
1094     ierr = PetscInfo2(U,"Scaling back h from %G to %G\n",PetscRealPart(*h),.99*val);CHKERRQ(ierr);
1095     if (PetscRealPart(*h) > 0.0) *h =  0.99*val;
1096     else                         *h = -0.99*val;
1097   }
1098   PetscFunctionReturn(0);
1099 }
1100 
1101 
1102 
1103 
1104 
1105 
1106 
1107