xref: /petsc/src/snes/mf/snesmfj.c (revision a2ec6df8280ae283e8b0f45d20ebe04f5d3dfc7b)
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       SLES sles;
527       KSP  ksp;
528       ierr = SNESGetSLES(mfctx->snes,&sles);CHKERRQ(ierr);
529       ierr = SLESGetKSP(sles,&ksp);CHKERRQ(ierr);
530       ierr = KSPSetMonitor(ksp,MatSNESMFKSPMonitor,PETSC_NULL,0);CHKERRQ(ierr);
531     }
532   }
533   ierr = PetscOptionsName("-snes_mf_check_positivity","Insure that U + h*a is nonnegative","MatSNESMFSetCheckh",&flg);CHKERRQ(ierr);
534   if (flg) {
535     ierr = MatSNESMFSetCheckh(mat,MatSNESMFCheckPositivity,0);CHKERRQ(ierr);
536   }
537   if (mfctx->ops->setfromoptions) {
538     ierr = (*mfctx->ops->setfromoptions)(mfctx);CHKERRQ(ierr);
539   }
540   ierr = PetscOptionsEnd();CHKERRQ(ierr);
541   PetscFunctionReturn(0);
542 }
543 
544 /*MC
545   MATMFFD = "mffd" - A matrix free matrix type.
546 
547   Level: advanced
548 
549 .seealso: MatCreateMF, MatCreateSNESMF
550 M*/
551 
552 #undef __FUNCT__
553 #define __FUNCT__ "MatCreate_MFFD"
554 EXTERN_C_BEGIN
555 int MatCreate_MFFD(Mat A)
556 {
557   MatSNESMFCtx mfctx;
558   int          ierr;
559 
560   PetscFunctionBegin;
561 #ifndef PETSC_USE_DYNAMIC_LIBRARIES
562   ierr = SNESInitializePackage(PETSC_NULL);                                                               CHKERRQ(ierr);
563 #endif
564 
565   PetscHeaderCreate(mfctx,_p_MatSNESMFCtx,struct _MFOps,MATSNESMFCTX_COOKIE,0,"SNESMF",A->comm,MatDestroy_MFFD,MatView_MFFD);
566   PetscLogObjectCreate(mfctx);
567   mfctx->sp              = 0;
568   mfctx->snes            = 0;
569   mfctx->error_rel       = PETSC_SQRT_MACHINE_EPSILON;
570   mfctx->recomputeperiod = 1;
571   mfctx->count           = 0;
572   mfctx->currenth        = 0.0;
573   mfctx->historyh        = PETSC_NULL;
574   mfctx->ncurrenth       = 0;
575   mfctx->maxcurrenth     = 0;
576   mfctx->type_name       = 0;
577   mfctx->usesnes         = PETSC_FALSE;
578 
579   mfctx->vshift          = 0.0;
580   mfctx->vscale          = 1.0;
581 
582   /*
583      Create the empty data structure to contain compute-h routines.
584      These will be filled in below from the command line options or
585      a later call with MatSNESMFSetType() or if that is not called
586      then it will default in the first use of MatMult_MFFD()
587   */
588   mfctx->ops->compute        = 0;
589   mfctx->ops->destroy        = 0;
590   mfctx->ops->view           = 0;
591   mfctx->ops->setfromoptions = 0;
592   mfctx->hctx                = 0;
593 
594   mfctx->func                = 0;
595   mfctx->funcctx             = 0;
596   mfctx->funcvec             = 0;
597   mfctx->w                   = PETSC_NULL;
598 
599   A->data                = mfctx;
600 
601   A->ops->mult           = MatMult_MFFD;
602   A->ops->destroy        = MatDestroy_MFFD;
603   A->ops->view           = MatView_MFFD;
604   A->ops->assemblyend    = MatAssemblyEnd_MFFD;
605   A->ops->getdiagonal    = MatGetDiagonal_MFFD;
606   A->ops->scale          = MatScale_MFFD;
607   A->ops->shift          = MatShift_MFFD;
608   A->ops->setfromoptions = MatSNESMFSetFromOptions;
609 
610   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSNESMFSetBase_C","MatSNESMFSetBase_FD",MatSNESMFSetBase_FD);CHKERRQ(ierr);
611   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSNESMFSetFunctioniBase_C","MatSNESMFSetFunctioniBase_FD",MatSNESMFSetFunctioniBase_FD);CHKERRQ(ierr);
612   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSNESMFSetFunctioni_C","MatSNESMFSetFunctioni_FD",MatSNESMFSetFunctioni_FD);CHKERRQ(ierr);
613   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSNESMFSetCheckh_C","MatSNESMFSetCheckh_FD",MatSNESMFSetCheckh_FD);CHKERRQ(ierr);
614   mfctx->mat = A;
615 
616   PetscFunctionReturn(0);
617 }
618 
619 EXTERN_C_END
620 
621 #undef __FUNCT__
622 #define __FUNCT__ "MatCreateMF"
623 /*@C
624    MatCreateMF - Creates a matrix-free matrix. See also MatCreateSNESMF()
625 
626    Collective on Vec
627 
628    Input Parameters:
629 .  x - vector that defines layout of the vectors and matrices
630 
631    Output Parameter:
632 .  J - the matrix-free matrix
633 
634    Level: advanced
635 
636    Notes:
637    The matrix-free matrix context merely contains the function pointers
638    and work space for performing finite difference approximations of
639    Jacobian-vector products, F'(u)*a,
640 
641    The default code uses the following approach to compute h
642 
643 .vb
644      F'(u)*a = [F(u+h*a) - F(u)]/h where
645      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
646        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
647  where
648      error_rel = square root of relative error in function evaluation
649      umin = minimum iterate parameter
650 .ve
651 
652    The user can set the error_rel via MatSNESMFSetFunctionError() and
653    umin via MatSNESMFDefaultSetUmin(); see the nonlinear solvers chapter
654    of the users manual for details.
655 
656    The user should call MatDestroy() when finished with the matrix-free
657    matrix context.
658 
659    Options Database Keys:
660 +  -snes_mf_err <error_rel> - Sets error_rel
661 .  -snes_mf_unim <umin> - Sets umin (for default PETSc routine that computes h only)
662 .  -snes_mf_ksp_monitor - KSP monitor routine that prints differencing h
663 -  -snes_mf_check_positivity
664 
665 .keywords: default, matrix-free, create, matrix
666 
667 .seealso: MatDestroy(), MatSNESMFSetFunctionError(), MatSNESMFDefaultSetUmin()
668           MatSNESMFSetHHistory(), MatSNESMFResetHHistory(), MatCreateSNESMF(),
669           MatSNESMFGetH(),MatSNESMFKSPMonitor(), MatSNESMFRegisterDynamic),, MatSNESMFComputeJacobian()
670 
671 @*/
672 int MatCreateMF(Vec x,Mat *J)
673 {
674   MPI_Comm     comm;
675   int          n,nloc,ierr;
676 
677   PetscFunctionBegin;
678   ierr = PetscObjectGetComm((PetscObject)x,&comm);CHKERRQ(ierr);
679   ierr = VecGetSize(x,&n);CHKERRQ(ierr);
680   ierr = VecGetLocalSize(x,&nloc);CHKERRQ(ierr);
681   ierr = MatCreate(comm,nloc,nloc,n,n,J);CHKERRQ(ierr);
682   ierr = MatRegister(MATMFFD,0,"MatCreate_MFFD",MatCreate_MFFD);CHKERRQ(ierr);
683   ierr = MatSetType(*J,MATMFFD);CHKERRQ(ierr);
684   PetscFunctionReturn(0);
685 }
686 
687 
688 #undef __FUNCT__
689 #define __FUNCT__ "MatSNESMFGetH"
690 /*@
691    MatSNESMFGetH - Gets the last value that was used as the differencing
692    parameter.
693 
694    Not Collective
695 
696    Input Parameters:
697 .  mat - the matrix obtained with MatCreateSNESMF()
698 
699    Output Paramter:
700 .  h - the differencing step size
701 
702    Level: advanced
703 
704 .keywords: SNES, matrix-free, parameters
705 
706 .seealso: MatCreateSNESMF(),MatSNESMFSetHHistory(),
707           MatSNESMFResetHHistory(),MatSNESMFKSPMonitor()
708 @*/
709 int MatSNESMFGetH(Mat mat,PetscScalar *h)
710 {
711   MatSNESMFCtx ctx = (MatSNESMFCtx)mat->data;
712 
713   PetscFunctionBegin;
714   *h = ctx->currenth;
715   PetscFunctionReturn(0);
716 }
717 
718 #undef __FUNCT__
719 #define __FUNCT__ "MatSNESMFKSPMonitor"
720 /*
721    MatSNESMFKSPMonitor - A KSP monitor for use with the default PETSc
722    SNES matrix free routines. Prints the differencing parameter used at
723    each step.
724 */
725 int MatSNESMFKSPMonitor(KSP ksp,int n,PetscReal rnorm,void *dummy)
726 {
727   PC             pc;
728   MatSNESMFCtx   ctx;
729   int            ierr;
730   Mat            mat;
731   MPI_Comm       comm;
732   PetscTruth     nonzeroinitialguess;
733 
734   PetscFunctionBegin;
735   ierr = PetscObjectGetComm((PetscObject)ksp,&comm);CHKERRQ(ierr);
736   ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
737   ierr = KSPGetInitialGuessNonzero(ksp,&nonzeroinitialguess);CHKERRQ(ierr);
738   ierr = PCGetOperators(pc,&mat,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
739   ctx  = (MatSNESMFCtx)mat->data;
740 
741   if (n > 0 || nonzeroinitialguess) {
742 #if defined(PETSC_USE_COMPLEX)
743     ierr = PetscPrintf(comm,"%d KSP Residual norm %14.12e h %g + %g i\n",n,rnorm,
744                 PetscRealPart(ctx->currenth),PetscImaginaryPart(ctx->currenth));CHKERRQ(ierr);
745 #else
746     ierr = PetscPrintf(comm,"%d KSP Residual norm %14.12e h %g \n",n,rnorm,ctx->currenth);CHKERRQ(ierr);
747 #endif
748   } else {
749     ierr = PetscPrintf(comm,"%d KSP Residual norm %14.12e\n",n,rnorm);CHKERRQ(ierr);
750   }
751   PetscFunctionReturn(0);
752 }
753 
754 #undef __FUNCT__
755 #define __FUNCT__ "MatSNESMFSetFunction"
756 /*@C
757    MatSNESMFSetFunction - Sets the function used in applying the matrix free.
758 
759    Collective on Mat
760 
761    Input Parameters:
762 +  mat - the matrix free matrix created via MatCreateSNESMF()
763 .  v   - workspace vector
764 .  func - the function to use
765 -  funcctx - optional function context passed to function
766 
767    Level: advanced
768 
769    Notes:
770     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
771     matrix inside your compute Jacobian routine
772 
773     If this is not set then it will use the function set with SNESSetFunction()
774 
775 .keywords: SNES, matrix-free, function
776 
777 .seealso: MatCreateSNESMF(),MatSNESMFGetH(),
778           MatSNESMFSetHHistory(), MatSNESMFResetHHistory(),
779           MatSNESMFKSPMonitor(), SNESetFunction()
780 @*/
781 int MatSNESMFSetFunction(Mat mat,Vec v,int (*func)(SNES,Vec,Vec,void *),void *funcctx)
782 {
783   MatSNESMFCtx ctx = (MatSNESMFCtx)mat->data;
784 
785   PetscFunctionBegin;
786   ctx->func    = func;
787   ctx->funcctx = funcctx;
788   ctx->funcvec = v;
789   PetscFunctionReturn(0);
790 }
791 
792 #undef __FUNCT__
793 #define __FUNCT__ "MatSNESMFSetFunctioni"
794 /*@C
795    MatSNESMFSetFunctioni - Sets the function for a single component
796 
797    Collective on Mat
798 
799    Input Parameters:
800 +  mat - the matrix free matrix created via MatCreateSNESMF()
801 -  funci - the function to use
802 
803    Level: advanced
804 
805    Notes:
806     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
807     matrix inside your compute Jacobian routine
808 
809 
810 .keywords: SNES, matrix-free, function
811 
812 .seealso: MatCreateSNESMF(),MatSNESMFGetH(),
813           MatSNESMFSetHHistory(), MatSNESMFResetHHistory(),
814           MatSNESMFKSPMonitor(), SNESetFunction()
815 @*/
816 int MatSNESMFSetFunctioni(Mat mat,int (*funci)(int,Vec,PetscScalar*,void *))
817 {
818   int  ierr,(*f)(Mat,int (*)(int,Vec,PetscScalar*,void *));
819 
820   PetscFunctionBegin;
821   PetscValidHeaderSpecific(mat,MAT_COOKIE);
822   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSNESMFSetFunctioni_C",(void (**)(void))&f);CHKERRQ(ierr);
823   if (f) {
824     ierr = (*f)(mat,funci);CHKERRQ(ierr);
825   }
826   PetscFunctionReturn(0);
827 }
828 
829 
830 #undef __FUNCT__
831 #define __FUNCT__ "MatSNESMFSetFunctioniBase"
832 /*@C
833    MatSNESMFSetFunctioniBase - Sets the base vector for a single component function evaluation
834 
835    Collective on Mat
836 
837    Input Parameters:
838 +  mat - the matrix free matrix created via MatCreateSNESMF()
839 -  func - the function to use
840 
841    Level: advanced
842 
843    Notes:
844     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
845     matrix inside your compute Jacobian routine
846 
847 
848 .keywords: SNES, matrix-free, function
849 
850 .seealso: MatCreateSNESMF(),MatSNESMFGetH(),
851           MatSNESMFSetHHistory(), MatSNESMFResetHHistory(),
852           MatSNESMFKSPMonitor(), SNESetFunction()
853 @*/
854 int MatSNESMFSetFunctioniBase(Mat mat,int (*func)(Vec,void *))
855 {
856   int  ierr,(*f)(Mat,int (*)(Vec,void *));
857 
858   PetscFunctionBegin;
859   PetscValidHeaderSpecific(mat,MAT_COOKIE);
860   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSNESMFSetFunctioniBase_C",(void (**)(void))&f);CHKERRQ(ierr);
861   if (f) {
862     ierr = (*f)(mat,func);CHKERRQ(ierr);
863   }
864   PetscFunctionReturn(0);
865 }
866 
867 
868 #undef __FUNCT__
869 #define __FUNCT__ "MatSNESMFSetPeriod"
870 /*@
871    MatSNESMFSetPeriod - Sets how often h is recomputed, by default it is everytime
872 
873    Collective on Mat
874 
875    Input Parameters:
876 +  mat - the matrix free matrix created via MatCreateSNESMF()
877 -  period - 1 for everytime, 2 for every second etc
878 
879    Options Database Keys:
880 +  -snes_mf_period <period>
881 
882    Level: advanced
883 
884 
885 .keywords: SNES, matrix-free, parameters
886 
887 .seealso: MatCreateSNESMF(),MatSNESMFGetH(),
888           MatSNESMFSetHHistory(), MatSNESMFResetHHistory(),
889           MatSNESMFKSPMonitor()
890 @*/
891 int MatSNESMFSetPeriod(Mat mat,int period)
892 {
893   MatSNESMFCtx ctx = (MatSNESMFCtx)mat->data;
894 
895   PetscFunctionBegin;
896   ctx->recomputeperiod = period;
897   PetscFunctionReturn(0);
898 }
899 
900 #undef __FUNCT__
901 #define __FUNCT__ "MatSNESMFSetFunctionError"
902 /*@
903    MatSNESMFSetFunctionError - Sets the error_rel for the approximation of
904    matrix-vector products using finite differences.
905 
906    Collective on Mat
907 
908    Input Parameters:
909 +  mat - the matrix free matrix created via MatCreateSNESMF()
910 -  error_rel - relative error (should be set to the square root of
911                the relative error in the function evaluations)
912 
913    Options Database Keys:
914 +  -snes_mf_err <error_rel> - Sets error_rel
915 
916    Level: advanced
917 
918    Notes:
919    The default matrix-free matrix-vector product routine computes
920 .vb
921      F'(u)*a = [F(u+h*a) - F(u)]/h where
922      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
923        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   else
924 .ve
925 
926 .keywords: SNES, matrix-free, parameters
927 
928 .seealso: MatCreateSNESMF(),MatSNESMFGetH(),
929           MatSNESMFSetHHistory(), MatSNESMFResetHHistory(),
930           MatSNESMFKSPMonitor()
931 @*/
932 int MatSNESMFSetFunctionError(Mat mat,PetscReal error)
933 {
934   MatSNESMFCtx ctx = (MatSNESMFCtx)mat->data;
935 
936   PetscFunctionBegin;
937   if (error != PETSC_DEFAULT) ctx->error_rel = error;
938   PetscFunctionReturn(0);
939 }
940 
941 #undef __FUNCT__
942 #define __FUNCT__ "MatSNESMFAddNullSpace"
943 /*@
944    MatSNESMFAddNullSpace - Provides a null space that an operator is
945    supposed to have.  Since roundoff will create a small component in
946    the null space, if you know the null space you may have it
947    automatically removed.
948 
949    Collective on Mat
950 
951    Input Parameters:
952 +  J - the matrix-free matrix context
953 -  nullsp - object created with MatNullSpaceCreate()
954 
955    Level: advanced
956 
957 .keywords: SNES, matrix-free, null space
958 
959 .seealso: MatNullSpaceCreate(), MatSNESMFGetH(), MatCreateSNESMF(),
960           MatSNESMFSetHHistory(), MatSNESMFResetHHistory(),
961           MatSNESMFKSPMonitor(), MatSNESMFErrorRel()
962 @*/
963 int MatSNESMFAddNullSpace(Mat J,MatNullSpace nullsp)
964 {
965   int          ierr;
966   MatSNESMFCtx ctx = (MatSNESMFCtx)J->data;
967   MPI_Comm     comm;
968 
969   PetscFunctionBegin;
970   ierr = PetscObjectGetComm((PetscObject)J,&comm);CHKERRQ(ierr);
971 
972   ctx->sp = nullsp;
973   ierr    = PetscObjectReference((PetscObject)nullsp);CHKERRQ(ierr);
974   PetscFunctionReturn(0);
975 }
976 
977 #undef __FUNCT__
978 #define __FUNCT__ "MatSNESMFSetHHistory"
979 /*@
980    MatSNESMFSetHHistory - Sets an array to collect a history of the
981    differencing values (h) computed for the matrix-free product.
982 
983    Collective on Mat
984 
985    Input Parameters:
986 +  J - the matrix-free matrix context
987 .  histroy - space to hold the history
988 -  nhistory - number of entries in history, if more entries are generated than
989               nhistory, then the later ones are discarded
990 
991    Level: advanced
992 
993    Notes:
994    Use MatSNESMFResetHHistory() to reset the history counter and collect
995    a new batch of differencing parameters, h.
996 
997 .keywords: SNES, matrix-free, h history, differencing history
998 
999 .seealso: MatSNESMFGetH(), MatCreateSNESMF(),
1000           MatSNESMFResetHHistory(),
1001           MatSNESMFKSPMonitor(), MatSNESMFSetFunctionError()
1002 
1003 @*/
1004 int MatSNESMFSetHHistory(Mat J,PetscScalar history[],int nhistory)
1005 {
1006   MatSNESMFCtx ctx = (MatSNESMFCtx)J->data;
1007 
1008   PetscFunctionBegin;
1009   ctx->historyh    = history;
1010   ctx->maxcurrenth = nhistory;
1011   ctx->currenth    = 0;
1012   PetscFunctionReturn(0);
1013 }
1014 
1015 #undef __FUNCT__
1016 #define __FUNCT__ "MatSNESMFResetHHistory"
1017 /*@
1018    MatSNESMFResetHHistory - Resets the counter to zero to begin
1019    collecting a new set of differencing histories.
1020 
1021    Collective on Mat
1022 
1023    Input Parameters:
1024 .  J - the matrix-free matrix context
1025 
1026    Level: advanced
1027 
1028    Notes:
1029    Use MatSNESMFSetHHistory() to create the original history counter.
1030 
1031 .keywords: SNES, matrix-free, h history, differencing history
1032 
1033 .seealso: MatSNESMFGetH(), MatCreateSNESMF(),
1034           MatSNESMFSetHHistory(),
1035           MatSNESMFKSPMonitor(), MatSNESMFSetFunctionError()
1036 
1037 @*/
1038 int MatSNESMFResetHHistory(Mat J)
1039 {
1040   MatSNESMFCtx ctx = (MatSNESMFCtx)J->data;
1041 
1042   PetscFunctionBegin;
1043   ctx->ncurrenth    = 0;
1044   PetscFunctionReturn(0);
1045 }
1046 
1047 #undef __FUNCT__
1048 #define __FUNCT__ "MatSNESMFComputeJacobian"
1049 int MatSNESMFComputeJacobian(SNES snes,Vec x,Mat *jac,Mat *B,MatStructure *flag,void *dummy)
1050 {
1051   int ierr;
1052   PetscFunctionBegin;
1053   ierr = MatAssemblyBegin(*jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1054   ierr = MatAssemblyEnd(*jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1055   PetscFunctionReturn(0);
1056 }
1057 
1058 #undef __FUNCT__
1059 #define __FUNCT__ "MatSNESMFSetBase"
1060 /*@
1061     MatSNESMFSetBase - Sets the vector U at which matrix vector products of the
1062         Jacobian are computed
1063 
1064     Collective on Mat
1065 
1066     Input Parameters:
1067 +   J - the MatSNESMF matrix
1068 -   U - the vector
1069 
1070     Notes: This is rarely used directly
1071 
1072     Level: advanced
1073 
1074 @*/
1075 int MatSNESMFSetBase(Mat J,Vec U)
1076 {
1077   int  ierr,(*f)(Mat,Vec);
1078 
1079   PetscFunctionBegin;
1080   PetscValidHeaderSpecific(J,MAT_COOKIE);
1081   PetscValidHeaderSpecific(U,VEC_COOKIE);
1082   ierr = PetscObjectQueryFunction((PetscObject)J,"MatSNESMFSetBase_C",(void (**)(void))&f);CHKERRQ(ierr);
1083   if (f) {
1084     ierr = (*f)(J,U);CHKERRQ(ierr);
1085   }
1086   PetscFunctionReturn(0);
1087 }
1088 
1089 #undef __FUNCT__
1090 #define __FUNCT__ "MatSNESMFSetCheckh"
1091 /*@C
1092     MatSNESMFSetCheckh - Sets a function that checks the computed h and adjusts
1093         it to satisfy some criteria
1094 
1095     Collective on Mat
1096 
1097     Input Parameters:
1098 +   J - the MatSNESMF matrix
1099 .   fun - the function that checks h
1100 -   ctx - any context needed by the function
1101 
1102     Options Database Keys:
1103 .   -snes_mf_check_positivity
1104 
1105     Level: advanced
1106 
1107     Notes: For example, MatSNESMFSetCheckPositivity() insures that all entries
1108        of U + h*a are non-negative
1109 
1110 .seealso:  MatSNESMFSetCheckPositivity()
1111 @*/
1112 int MatSNESMFSetCheckh(Mat J,int (*fun)(Vec,Vec,PetscScalar*,void*),void* ctx)
1113 {
1114   int  ierr,(*f)(Mat,int (*)(Vec,Vec,PetscScalar*,void*),void*);
1115 
1116   PetscFunctionBegin;
1117   PetscValidHeaderSpecific(J,MAT_COOKIE);
1118   ierr = PetscObjectQueryFunction((PetscObject)J,"MatSNESMFSetCheckh_C",(void (**)(void))&f);CHKERRQ(ierr);
1119   if (f) {
1120     ierr = (*f)(J,fun,ctx);CHKERRQ(ierr);
1121   }
1122   PetscFunctionReturn(0);
1123 }
1124 
1125 #undef __FUNCT__
1126 #define __FUNCT__ "MatSNESMFSetCheckPositivity"
1127 /*@
1128     MatSNESMFCheckPositivity - Checks that all entries in U + h*a are positive or
1129         zero, decreases h until this is satisfied.
1130 
1131     Collective on Vec
1132 
1133     Input Parameters:
1134 +   U - base vector that is added to
1135 .   a - vector that is added
1136 .   h - scaling factor on a
1137 -   dummy - context variable (unused)
1138 
1139     Options Database Keys:
1140 .   -snes_mf_check_positivity
1141 
1142     Level: advanced
1143 
1144     Notes: This is rarely used directly, rather it is passed as an argument to
1145            MatSNESMFSetCheckh()
1146 
1147 .seealso:  MatSNESMFSetCheckh()
1148 @*/
1149 int MatSNESMFCheckPositivity(Vec U,Vec a,PetscScalar *h,void *dummy)
1150 {
1151   PetscReal     val, minval;
1152   PetscScalar   *u_vec, *a_vec;
1153   int           ierr, i, size;
1154   MPI_Comm      comm;
1155 
1156   PetscFunctionBegin;
1157   ierr = PetscObjectGetComm((PetscObject)U,&comm);CHKERRQ(ierr);
1158   ierr = VecGetArray(U,&u_vec);CHKERRQ(ierr);
1159   ierr = VecGetArray(a,&a_vec);CHKERRQ(ierr);
1160   ierr = VecGetLocalSize(U,&size);CHKERRQ(ierr);
1161   minval = PetscAbsScalar(*h*1.01);
1162   for(i=0;i<size;i++) {
1163     if (PetscRealPart(u_vec[i] + *h*a_vec[i]) <= 0.0) {
1164       val = PetscAbsScalar(u_vec[i]/a_vec[i]);
1165       if (val < minval) minval = val;
1166     }
1167   }
1168   ierr = VecRestoreArray(U,&u_vec);CHKERRQ(ierr);
1169   ierr = VecRestoreArray(a,&a_vec);CHKERRQ(ierr);
1170   ierr = PetscGlobalMin(&minval,&val,comm);CHKERRQ(ierr);
1171   if (val <= PetscAbsScalar(*h)) {
1172     PetscLogInfo(U,"MatSNESMFCheckPositivity: Scaling back h from %g to %g\n",PetscRealPart(*h),.99*val);
1173     if (PetscRealPart(*h) > 0.0) *h =  0.99*val;
1174     else                         *h = -0.99*val;
1175   }
1176   PetscFunctionReturn(0);
1177 }
1178 
1179 
1180 
1181 
1182 
1183 
1184 
1185