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