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