xref: /petsc/src/snes/mf/snesmfj.c (revision 46129b977243d4edc5b889e5cc50ffbb698cfa6e)
1 
2 #include "src/mat/matimpl.h"
3 #include "src/snes/mf/snesmfj.h"   /*I  "petscsnes.h"   I*/
4 
5 PetscFList MatSNESMPetscFList         = 0;
6 PetscTruth MatSNESMFRegisterAllCalled = PETSC_FALSE;
7 
8 PetscCookie MATSNESMFCTX_COOKIE = 0;
9 PetscEvent  MATSNESMF_Mult = 0;
10 
11 #undef __FUNCT__
12 #define __FUNCT__ "MatSNESMFSetType"
13 /*@C
14     MatSNESMFSetType - Sets the method that is used to compute the
15     differencing parameter for finite differene matrix-free formulations.
16 
17     Input Parameters:
18 +   mat - the "matrix-free" matrix created via MatCreateSNESMF(), or MatCreateMF()
19           or MatSetType(mat,MATMFFD);
20 -   ftype - the type requested
21 
22     Level: advanced
23 
24     Notes:
25     For example, such routines can compute h for use in
26     Jacobian-vector products of the form
27 
28                         F(x+ha) - F(x)
29           F'(u)a  ~=  ----------------
30                               h
31 
32 .seealso: MatCreateSNESMF(), MatSNESMFRegisterDynamic)
33 @*/
34 PetscErrorCode MatSNESMFSetType(Mat mat,const MatSNESMFType ftype)
35 {
36   PetscErrorCode ierr,(*r)(MatSNESMFCtx);
37   MatSNESMFCtx   ctx = (MatSNESMFCtx)mat->data;
38   PetscTruth     match;
39 
40   PetscFunctionBegin;
41   PetscValidHeaderSpecific(mat,MAT_COOKIE,1);
42   PetscValidCharPointer(ftype,2);
43 
44   /* already set, so just return */
45   ierr = PetscTypeCompare((PetscObject)ctx,ftype,&match);CHKERRQ(ierr);
46   if (match) PetscFunctionReturn(0);
47 
48   /* destroy the old one if it exists */
49   if (ctx->ops->destroy) {
50     ierr = (*ctx->ops->destroy)(ctx);CHKERRQ(ierr);
51   }
52 
53   /* Get the function pointers for the requrested method */
54   if (!MatSNESMFRegisterAllCalled) {ierr = MatSNESMFRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
55   ierr =  PetscFListFind(ctx->comm,MatSNESMPetscFList,ftype,(void (**)(void)) &r);CHKERRQ(ierr);
56   if (!r) SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown MatSNESMF type %s given",ftype);
57   ierr = (*r)(ctx);CHKERRQ(ierr);
58   ierr = PetscObjectChangeTypeName((PetscObject)ctx,ftype);CHKERRQ(ierr);
59   PetscFunctionReturn(0);
60 }
61 
62 typedef PetscErrorCode (*FCN1)(Vec,void*); /* force argument to next function to not be extern C*/
63 EXTERN_C_BEGIN
64 #undef __FUNCT__
65 #define __FUNCT__ "MatSNESMFSetFunctioniBase_FD"
66 PetscErrorCode MatSNESMFSetFunctioniBase_FD(Mat mat,FCN1 func)
67 {
68   MatSNESMFCtx ctx = (MatSNESMFCtx)mat->data;
69 
70   PetscFunctionBegin;
71   ctx->funcisetbase = func;
72   PetscFunctionReturn(0);
73 }
74 EXTERN_C_END
75 
76 typedef PetscErrorCode (*FCN2)(PetscInt,Vec,PetscScalar*,void*); /* force argument to next function to not be extern C*/
77 EXTERN_C_BEGIN
78 #undef __FUNCT__
79 #define __FUNCT__ "MatSNESMFSetFunctioni_FD"
80 PetscErrorCode MatSNESMFSetFunctioni_FD(Mat mat,FCN2 funci)
81 {
82   MatSNESMFCtx ctx = (MatSNESMFCtx)mat->data;
83 
84   PetscFunctionBegin;
85   ctx->funci = funci;
86   PetscFunctionReturn(0);
87 }
88 EXTERN_C_END
89 
90 
91 #undef __FUNCT__
92 #define __FUNCT__ "MatSNESMFRegister"
93 PetscErrorCode MatSNESMFRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(MatSNESMFCtx))
94 {
95   PetscErrorCode ierr;
96   char           fullname[PETSC_MAX_PATH_LEN];
97 
98   PetscFunctionBegin;
99   ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
100   ierr = PetscFListAdd(&MatSNESMPetscFList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
101   PetscFunctionReturn(0);
102 }
103 
104 
105 #undef __FUNCT__
106 #define __FUNCT__ "MatSNESMFRegisterDestroy"
107 /*@C
108    MatSNESMFRegisterDestroy - Frees the list of MatSNESMF methods that were
109    registered by MatSNESMFRegisterDynamic).
110 
111    Not Collective
112 
113    Level: developer
114 
115 .keywords: MatSNESMF, register, destroy
116 
117 .seealso: MatSNESMFRegisterDynamic), MatSNESMFRegisterAll()
118 @*/
119 PetscErrorCode MatSNESMFRegisterDestroy(void)
120 {
121   PetscErrorCode ierr;
122 
123   PetscFunctionBegin;
124   if (MatSNESMPetscFList) {
125     ierr = PetscFListDestroy(&MatSNESMPetscFList);CHKERRQ(ierr);
126     MatSNESMPetscFList = 0;
127   }
128   MatSNESMFRegisterAllCalled = PETSC_FALSE;
129   PetscFunctionReturn(0);
130 }
131 
132 /* ----------------------------------------------------------------------------------------*/
133 #undef __FUNCT__
134 #define __FUNCT__ "MatDestroy_MFFD"
135 PetscErrorCode MatDestroy_MFFD(Mat mat)
136 {
137   PetscErrorCode ierr;
138   MatSNESMFCtx   ctx = (MatSNESMFCtx)mat->data;
139 
140   PetscFunctionBegin;
141   if (ctx->w) {
142     ierr = VecDestroy(ctx->w);CHKERRQ(ierr);
143   }
144   if (ctx->ops->destroy) {ierr = (*ctx->ops->destroy)(ctx);CHKERRQ(ierr);}
145   if (ctx->sp) {ierr = MatNullSpaceDestroy(ctx->sp);CHKERRQ(ierr);}
146   ierr = PetscHeaderDestroy(ctx);CHKERRQ(ierr);
147 
148   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatSNESMFSetBase_C","",PETSC_NULL);CHKERRQ(ierr);
149   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatSNESMFSetFunctioniBase_C","",PETSC_NULL);CHKERRQ(ierr);
150   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatSNESMFSetFunctioni_C","",PETSC_NULL);CHKERRQ(ierr);
151   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatSNESMFSetCheckh_C","",PETSC_NULL);CHKERRQ(ierr);
152 
153   PetscFunctionReturn(0);
154 }
155 
156 #undef __FUNCT__
157 #define __FUNCT__ "MatView_MFFD"
158 /*
159    MatSNESMFView_MFFD - Views matrix-free parameters.
160 
161 */
162 PetscErrorCode MatView_MFFD(Mat J,PetscViewer viewer)
163 {
164   PetscErrorCode ierr;
165   MatSNESMFCtx   ctx = (MatSNESMFCtx)J->data;
166   PetscTruth     iascii;
167 
168   PetscFunctionBegin;
169   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
170   if (iascii) {
171      ierr = PetscViewerASCIIPrintf(viewer,"  SNES matrix-free approximation:\n");CHKERRQ(ierr);
172      ierr = PetscViewerASCIIPrintf(viewer,"    err=%g (relative error in function evaluation)\n",ctx->error_rel);CHKERRQ(ierr);
173      if (!ctx->type_name) {
174        ierr = PetscViewerASCIIPrintf(viewer,"    The compute h routine has not yet been set\n");CHKERRQ(ierr);
175      } else {
176        ierr = PetscViewerASCIIPrintf(viewer,"    Using %s compute h routine\n",ctx->type_name);CHKERRQ(ierr);
177      }
178      if (ctx->ops->view) {
179        ierr = (*ctx->ops->view)(ctx,viewer);CHKERRQ(ierr);
180      }
181   } else {
182     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for SNES matrix free matrix",((PetscObject)viewer)->type_name);
183   }
184   PetscFunctionReturn(0);
185 }
186 
187 #undef __FUNCT__
188 #define __FUNCT__ "MatAssemblyEnd_MFFD"
189 /*
190    MatAssemblyEnd_MFFD - Resets the ctx->ncurrenth to zero. This
191    allows the user to indicate the beginning of a new linear solve by calling
192    MatAssemblyXXX() on the matrix free matrix. This then allows the
193    MatSNESMFCreate_WP() to properly compute ||U|| only the first time
194    in the linear solver rather than every time.
195 */
196 PetscErrorCode MatAssemblyEnd_MFFD(Mat J,MatAssemblyType mt)
197 {
198   PetscErrorCode ierr;
199   MatSNESMFCtx   j = (MatSNESMFCtx)J->data;
200 
201   PetscFunctionBegin;
202   ierr = MatSNESMFResetHHistory(J);CHKERRQ(ierr);
203   if (j->usesnes) {
204     ierr = SNESGetSolution(j->snes,&j->current_u);CHKERRQ(ierr);
205     ierr = SNESGetFunction(j->snes,&j->current_f,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
206     if (!j->w) {
207       ierr = VecDuplicate(j->current_u, &j->w);CHKERRQ(ierr);
208     }
209   }
210   j->vshift = 0.0;
211   j->vscale = 1.0;
212   PetscFunctionReturn(0);
213 }
214 
215 #undef __FUNCT__
216 #define __FUNCT__ "MatMult_MFFD"
217 /*
218   MatMult_MFFD - Default matrix-free form for Jacobian-vector product, y = F'(u)*a:
219 
220         y ~= (F(u + ha) - F(u))/h,
221   where F = nonlinear function, as set by SNESSetFunction()
222         u = current iterate
223         h = difference interval
224 */
225 PetscErrorCode MatMult_MFFD(Mat mat,Vec a,Vec y)
226 {
227   MatSNESMFCtx    ctx = (MatSNESMFCtx)mat->data;
228   SNES            snes;
229   PetscScalar     h,mone = -1.0;
230   Vec             w,U,F;
231   PetscErrorCode ierr,(*eval_fct)(SNES,Vec,Vec)=0;
232 
233   PetscFunctionBegin;
234   /* We log matrix-free matrix-vector products separately, so that we can
235      separate the performance monitoring from the cases that use conventional
236      storage.  We may eventually modify event logging to associate events
237      with particular objects, hence alleviating the more general problem. */
238   ierr = PetscLogEventBegin(MATSNESMF_Mult,a,y,0,0);CHKERRQ(ierr);
239 
240   snes = ctx->snes;
241   w    = ctx->w;
242   U    = ctx->current_u;
243 
244   /*
245       Compute differencing parameter
246   */
247   if (!ctx->ops->compute) {
248     ierr = MatSNESMFSetType(mat,MATSNESMF_WP);CHKERRQ(ierr);
249     ierr = MatSNESMFSetFromOptions(mat);CHKERRQ(ierr);
250   }
251   ierr = (*ctx->ops->compute)(ctx,U,a,&h);CHKERRQ(ierr);
252 
253   if (ctx->checkh) {
254     ierr = (*ctx->checkh)(U,a,&h,ctx->checkhctx);CHKERRQ(ierr);
255   }
256 
257   /* keep a record of the current differencing parameter h */
258   ctx->currenth = h;
259 #if defined(PETSC_USE_COMPLEX)
260   ierr = PetscLogInfo((mat,"MatMult_MFFD:Current differencing parameter: %g + %g i\n",PetscRealPart(h),PetscImaginaryPart(h)));CHKERRQ(ierr);
261 #else
262   ierr = PetscLogInfo((mat,"MatMult_MFFD:Current differencing parameter: %15.12e\n",h));CHKERRQ(ierr);
263 #endif
264   if (ctx->historyh && ctx->ncurrenth < ctx->maxcurrenth) {
265     ctx->historyh[ctx->ncurrenth] = h;
266   }
267   ctx->ncurrenth++;
268 
269   /* w = u + ha */
270   ierr = VecWAXPY(&h,a,U,w);CHKERRQ(ierr);
271 
272   if (ctx->usesnes) {
273     eval_fct = SNESComputeFunction;
274     F    = ctx->current_f;
275     if (!F) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"You must call MatAssembly() even on matrix-free matrices");
276     ierr = (*eval_fct)(snes,w,y);CHKERRQ(ierr);
277   } else {
278     F = ctx->funcvec;
279     /* compute func(U) as base for differencing */
280     if (ctx->ncurrenth == 1) {
281       ierr = (*ctx->func)(snes,U,F,ctx->funcctx);CHKERRQ(ierr);
282     }
283     ierr = (*ctx->func)(snes,w,y,ctx->funcctx);CHKERRQ(ierr);
284   }
285 
286   ierr = VecAXPY(&mone,F,y);CHKERRQ(ierr);
287   h    = 1.0/h;
288   ierr = VecScale(&h,y);CHKERRQ(ierr);
289 
290   ierr = VecAXPBY(&ctx->vshift,&ctx->vscale,a,y);CHKERRQ(ierr);
291 
292   if (ctx->sp) {ierr = MatNullSpaceRemove(ctx->sp,y,PETSC_NULL);CHKERRQ(ierr);}
293 
294   ierr = PetscLogEventEnd(MATSNESMF_Mult,a,y,0,0);CHKERRQ(ierr);
295   PetscFunctionReturn(0);
296 }
297 
298 #undef __FUNCT__
299 #define __FUNCT__ "MatGetDiagonal_MFFD"
300 /*
301   MatGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix
302 
303         y ~= (F(u + ha) - F(u))/h,
304   where F = nonlinear function, as set by SNESSetFunction()
305         u = current iterate
306         h = difference interval
307 */
308 PetscErrorCode MatGetDiagonal_MFFD(Mat mat,Vec a)
309 {
310   MatSNESMFCtx   ctx = (MatSNESMFCtx)mat->data;
311   PetscScalar    h,*aa,*ww,v;
312   PetscReal      epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON;
313   Vec            w,U;
314   PetscErrorCode ierr;
315   PetscInt       i,rstart,rend;
316 
317   PetscFunctionBegin;
318   if (!ctx->funci) {
319     SETERRQ(PETSC_ERR_ORDER,"Requires calling MatSNESMFSetFunctioni() first");
320   }
321 
322   w    = ctx->w;
323   U    = ctx->current_u;
324   ierr = (*ctx->func)(0,U,a,ctx->funcctx);CHKERRQ(ierr);
325   ierr = (*ctx->funcisetbase)(U,ctx->funcctx);CHKERRQ(ierr);
326   ierr = VecCopy(U,w);CHKERRQ(ierr);
327 
328   ierr = VecGetOwnershipRange(a,&rstart,&rend);CHKERRQ(ierr);
329   ierr = VecGetArray(a,&aa);CHKERRQ(ierr);
330   for (i=rstart; i<rend; i++) {
331     ierr = VecGetArray(w,&ww);CHKERRQ(ierr);
332     h  = ww[i-rstart];
333     if (h == 0.0) h = 1.0;
334 #if !defined(PETSC_USE_COMPLEX)
335     if (h < umin && h >= 0.0)      h = umin;
336     else if (h < 0.0 && h > -umin) h = -umin;
337 #else
338     if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0)     h = umin;
339     else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin;
340 #endif
341     h     *= epsilon;
342 
343     ww[i-rstart] += h;
344     ierr = VecRestoreArray(w,&ww);CHKERRQ(ierr);
345     ierr          = (*ctx->funci)(i,w,&v,ctx->funcctx);CHKERRQ(ierr);
346     aa[i-rstart]  = (v - aa[i-rstart])/h;
347 
348     /* possibly shift and scale result */
349     aa[i - rstart] = ctx->vshift + ctx->vscale*aa[i-rstart];
350 
351     ierr = VecGetArray(w,&ww);CHKERRQ(ierr);
352     ww[i-rstart] -= h;
353     ierr = VecRestoreArray(w,&ww);CHKERRQ(ierr);
354   }
355   ierr = VecRestoreArray(a,&aa);CHKERRQ(ierr);
356   PetscFunctionReturn(0);
357 }
358 
359 #undef __FUNCT__
360 #define __FUNCT__ "MatShift_MFFD"
361 PetscErrorCode MatShift_MFFD(const PetscScalar *a,Mat Y)
362 {
363   MatSNESMFCtx shell = (MatSNESMFCtx)Y->data;
364   PetscFunctionBegin;
365   shell->vshift += *a;
366   PetscFunctionReturn(0);
367 }
368 
369 #undef __FUNCT__
370 #define __FUNCT__ "MatScale_MFFD"
371 PetscErrorCode MatScale_MFFD(const PetscScalar *a,Mat Y)
372 {
373   MatSNESMFCtx shell = (MatSNESMFCtx)Y->data;
374   PetscFunctionBegin;
375   shell->vscale *= *a;
376   PetscFunctionReturn(0);
377 }
378 
379 
380 #undef __FUNCT__
381 #define __FUNCT__ "MatCreateSNESMF"
382 /*@C
383    MatCreateSNESMF - Creates a matrix-free matrix context for use with
384    a SNES solver.  This matrix can be used as the Jacobian argument for
385    the routine SNESSetJacobian().
386 
387    Collective on SNES and Vec
388 
389    Input Parameters:
390 +  snes - the SNES context
391 -  x - vector where SNES solution is to be stored.
392 
393    Output Parameter:
394 .  J - the matrix-free matrix
395 
396    Level: advanced
397 
398    Notes:
399    The matrix-free matrix context merely contains the function pointers
400    and work space for performing finite difference approximations of
401    Jacobian-vector products, F'(u)*a,
402 
403    The default code uses the following approach to compute h
404 
405 .vb
406      F'(u)*a = [F(u+h*a) - F(u)]/h where
407      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
408        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
409  where
410      error_rel = square root of relative error in function evaluation
411      umin = minimum iterate parameter
412 .ve
413 
414    The user can set the error_rel via MatSNESMFSetFunctionError() and
415    umin via MatSNESMFDefaultSetUmin(); see the nonlinear solvers chapter
416    of the users manual for details.
417 
418    The user should call MatDestroy() when finished with the matrix-free
419    matrix context.
420 
421    Options Database Keys:
422 +  -snes_mf_err <error_rel> - Sets error_rel
423 .  -snes_mf_unim <umin> - Sets umin (for default PETSc routine that computes h only)
424 -  -snes_mf_ksp_monitor - KSP monitor routine that prints differencing h
425 
426 .keywords: SNES, default, matrix-free, create, matrix
427 
428 .seealso: MatDestroy(), MatSNESMFSetFunctionError(), MatSNESMFDefaultSetUmin()
429           MatSNESMFSetHHistory(), MatSNESMFResetHHistory(), MatCreateMF(),
430           MatSNESMFGetH(),MatSNESMFKSPMonitor(), MatSNESMFRegisterDynamic), MatSNESMFComputeJacobian()
431 
432 @*/
433 PetscErrorCode MatCreateSNESMF(SNES snes,Vec x,Mat *J)
434 {
435   MatSNESMFCtx   mfctx;
436   PetscErrorCode ierr;
437 
438   PetscFunctionBegin;
439   ierr = MatCreateMF(x,J);CHKERRQ(ierr);
440 
441   mfctx          = (MatSNESMFCtx)(*J)->data;
442   mfctx->snes    = snes;
443   mfctx->usesnes = PETSC_TRUE;
444   ierr = PetscLogObjectParent(snes,*J);CHKERRQ(ierr);
445   PetscFunctionReturn(0);
446 }
447 
448 EXTERN_C_BEGIN
449 #undef __FUNCT__
450 #define __FUNCT__ "MatSNESMFSetBase_FD"
451 PetscErrorCode MatSNESMFSetBase_FD(Mat J,Vec U)
452 {
453   PetscErrorCode ierr;
454   MatSNESMFCtx   ctx = (MatSNESMFCtx)J->data;
455 
456   PetscFunctionBegin;
457   ierr = MatSNESMFResetHHistory(J);CHKERRQ(ierr);
458   ctx->current_u = U;
459   ctx->usesnes   = PETSC_FALSE;
460   if (!ctx->w) {
461     ierr = VecDuplicate(ctx->current_u, &ctx->w);CHKERRQ(ierr);
462   }
463   J->assembled = PETSC_TRUE;
464   PetscFunctionReturn(0);
465 }
466 EXTERN_C_END
467 
468 typedef PetscErrorCode (*FCN3)(Vec,Vec,PetscScalar*,void*); /* force argument to next function to not be extern C*/
469 EXTERN_C_BEGIN
470 #undef __FUNCT__
471 #define __FUNCT__ "MatSNESMFSetCheckh_FD"
472 PetscErrorCode MatSNESMFSetCheckh_FD(Mat J,FCN3 fun,void*ectx)
473 {
474   MatSNESMFCtx ctx = (MatSNESMFCtx)J->data;
475 
476   PetscFunctionBegin;
477   ctx->checkh    = fun;
478   ctx->checkhctx = ectx;
479   PetscFunctionReturn(0);
480 }
481 EXTERN_C_END
482 
483 #undef __FUNCT__
484 #define __FUNCT__ "MatSNESMFSetFromOptions"
485 /*@
486    MatSNESMFSetFromOptions - Sets the MatSNESMF options from the command line
487    parameter.
488 
489    Collective on Mat
490 
491    Input Parameters:
492 .  mat - the matrix obtained with MatCreateSNESMF()
493 
494    Options Database Keys:
495 +  -snes_mf_type - <default,wp>
496 -  -snes_mf_err - square root of estimated relative error in function evaluation
497 -  -snes_mf_period - how often h is recomputed, defaults to 1, everytime
498 
499    Level: advanced
500 
501 .keywords: SNES, matrix-free, parameters
502 
503 .seealso: MatCreateSNESMF(),MatSNESMFSetHHistory(),
504           MatSNESMFResetHHistory(), MatSNESMFKSPMonitor()
505 @*/
506 PetscErrorCode MatSNESMFSetFromOptions(Mat mat)
507 {
508   MatSNESMFCtx   mfctx = (MatSNESMFCtx)mat->data;
509   PetscErrorCode ierr;
510   PetscTruth     flg;
511   char           ftype[256];
512 
513   PetscFunctionBegin;
514   if (!MatSNESMFRegisterAllCalled) {ierr = MatSNESMFRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
515 
516   ierr = PetscOptionsBegin(mfctx->comm,mfctx->prefix,"Set matrix free computation parameters","MatSNESMF");CHKERRQ(ierr);
517   ierr = PetscOptionsList("-snes_mf_type","Matrix free type","MatSNESMFSetType",MatSNESMPetscFList,mfctx->type_name,ftype,256,&flg);CHKERRQ(ierr);
518   if (flg) {
519     ierr = MatSNESMFSetType(mat,ftype);CHKERRQ(ierr);
520   }
521 
522   ierr = PetscOptionsReal("-snes_mf_err","set sqrt relative error in function","MatSNESMFSetFunctionError",mfctx->error_rel,&mfctx->error_rel,0);CHKERRQ(ierr);
523   ierr = PetscOptionsInt("-snes_mf_period","how often h is recomputed","MatSNESMFSetPeriod",mfctx->recomputeperiod,&mfctx->recomputeperiod,0);CHKERRQ(ierr);
524   if (mfctx->snes) {
525     ierr = PetscOptionsName("-snes_mf_ksp_monitor","Monitor matrix-free parameters","MatSNESMFKSPMonitor",&flg);CHKERRQ(ierr);
526     if (flg) {
527       KSP ksp;
528       ierr = SNESGetKSP(mfctx->snes,&ksp);CHKERRQ(ierr);
529       ierr = KSPSetMonitor(ksp,MatSNESMFKSPMonitor,PETSC_NULL,0);CHKERRQ(ierr);
530     }
531   }
532   ierr = PetscOptionsName("-snes_mf_check_positivity","Insure that U + h*a is nonnegative","MatSNESMFSetCheckh",&flg);CHKERRQ(ierr);
533   if (flg) {
534     ierr = MatSNESMFSetCheckh(mat,MatSNESMFCheckPositivity,0);CHKERRQ(ierr);
535   }
536   if (mfctx->ops->setfromoptions) {
537     ierr = (*mfctx->ops->setfromoptions)(mfctx);CHKERRQ(ierr);
538   }
539   ierr = PetscOptionsEnd();CHKERRQ(ierr);
540   PetscFunctionReturn(0);
541 }
542 
543 /*MC
544   MATMFFD - MATMFFD = "mffd" - A matrix free matrix type.
545 
546   Level: advanced
547 
548 .seealso: MatCreateMF(), MatCreateSNESMF()
549 M*/
550 EXTERN_C_BEGIN
551 #undef __FUNCT__
552 #define __FUNCT__ "MatCreate_MFFD"
553 PetscErrorCode MatCreate_MFFD(Mat A)
554 {
555   MatSNESMFCtx mfctx;
556   PetscErrorCode ierr;
557 
558   PetscFunctionBegin;
559 #ifndef PETSC_USE_DYNAMIC_LIBRARIES
560   ierr = SNESInitializePackage(PETSC_NULL);CHKERRQ(ierr);
561 #endif
562 
563   ierr = PetscHeaderCreate(mfctx,_p_MatSNESMFCtx,struct _MFOps,MATSNESMFCTX_COOKIE,0,"SNESMF",A->comm,MatDestroy_MFFD,MatView_MFFD);CHKERRQ(ierr);
564   mfctx->sp              = 0;
565   mfctx->snes            = 0;
566   mfctx->error_rel       = PETSC_SQRT_MACHINE_EPSILON;
567   mfctx->recomputeperiod = 1;
568   mfctx->count           = 0;
569   mfctx->currenth        = 0.0;
570   mfctx->historyh        = PETSC_NULL;
571   mfctx->ncurrenth       = 0;
572   mfctx->maxcurrenth     = 0;
573   mfctx->type_name       = 0;
574   mfctx->usesnes         = PETSC_FALSE;
575 
576   mfctx->vshift          = 0.0;
577   mfctx->vscale          = 1.0;
578 
579   /*
580      Create the empty data structure to contain compute-h routines.
581      These will be filled in below from the command line options or
582      a later call with MatSNESMFSetType() or if that is not called
583      then it will default in the first use of MatMult_MFFD()
584   */
585   mfctx->ops->compute        = 0;
586   mfctx->ops->destroy        = 0;
587   mfctx->ops->view           = 0;
588   mfctx->ops->setfromoptions = 0;
589   mfctx->hctx                = 0;
590 
591   mfctx->func                = 0;
592   mfctx->funcctx             = 0;
593   mfctx->funcvec             = 0;
594   mfctx->w                   = PETSC_NULL;
595 
596   A->data                = mfctx;
597 
598   A->ops->mult           = MatMult_MFFD;
599   A->ops->destroy        = MatDestroy_MFFD;
600   A->ops->view           = MatView_MFFD;
601   A->ops->assemblyend    = MatAssemblyEnd_MFFD;
602   A->ops->getdiagonal    = MatGetDiagonal_MFFD;
603   A->ops->scale          = MatScale_MFFD;
604   A->ops->shift          = MatShift_MFFD;
605   A->ops->setfromoptions = MatSNESMFSetFromOptions;
606   A->assembled = PETSC_TRUE;
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 EXTERN_C_END
617 
618 #undef __FUNCT__
619 #define __FUNCT__ "MatCreateMF"
620 /*@C
621    MatCreateMF - Creates a matrix-free matrix. See also MatCreateSNESMF()
622 
623    Collective on Vec
624 
625    Input Parameters:
626 .  x - vector that defines layout of the vectors and matrices
627 
628    Output Parameter:
629 .  J - the matrix-free matrix
630 
631    Level: advanced
632 
633    Notes:
634    The matrix-free matrix context merely contains the function pointers
635    and work space for performing finite difference approximations of
636    Jacobian-vector products, F'(u)*a,
637 
638    The default code uses the following approach to compute h
639 
640 .vb
641      F'(u)*a = [F(u+h*a) - F(u)]/h where
642      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
643        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
644  where
645      error_rel = square root of relative error in function evaluation
646      umin = minimum iterate parameter
647 .ve
648 
649    The user can set the error_rel via MatSNESMFSetFunctionError() and
650    umin via MatSNESMFDefaultSetUmin(); see the nonlinear solvers chapter
651    of the users manual for details.
652 
653    The user should call MatDestroy() when finished with the matrix-free
654    matrix context.
655 
656    Options Database Keys:
657 +  -snes_mf_err <error_rel> - Sets error_rel
658 .  -snes_mf_unim <umin> - Sets umin (for default PETSc routine that computes h only)
659 .  -snes_mf_ksp_monitor - KSP monitor routine that prints differencing h
660 -  -snes_mf_check_positivity
661 
662 .keywords: default, matrix-free, create, matrix
663 
664 .seealso: MatDestroy(), MatSNESMFSetFunctionError(), MatSNESMFDefaultSetUmin()
665           MatSNESMFSetHHistory(), MatSNESMFResetHHistory(), MatCreateSNESMF(),
666           MatSNESMFGetH(),MatSNESMFKSPMonitor(), MatSNESMFRegisterDynamic),, MatSNESMFComputeJacobian()
667 
668 @*/
669 PetscErrorCode MatCreateMF(Vec x,Mat *J)
670 {
671   MPI_Comm       comm;
672   PetscErrorCode ierr;
673   PetscInt       n,nloc;
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 = MatRegisterDynamic(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 PetscErrorCode 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 PetscErrorCode MatSNESMFKSPMonitor(KSP ksp,PetscInt n,PetscReal rnorm,void *dummy)
724 {
725   PC             pc;
726   MatSNESMFCtx   ctx;
727   PetscErrorCode 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 PetscErrorCode MatSNESMFSetFunction(Mat mat,Vec v,PetscErrorCode (*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 PetscErrorCode MatSNESMFSetFunctioni(Mat mat,PetscErrorCode (*funci)(PetscInt,Vec,PetscScalar*,void *))
815 {
816   PetscErrorCode ierr,(*f)(Mat,PetscErrorCode (*)(PetscInt,Vec,PetscScalar*,void *));
817 
818   PetscFunctionBegin;
819   PetscValidHeaderSpecific(mat,MAT_COOKIE,1);
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 PetscErrorCode MatSNESMFSetFunctioniBase(Mat mat,PetscErrorCode (*func)(Vec,void *))
853 {
854   PetscErrorCode ierr,(*f)(Mat,PetscErrorCode (*)(Vec,void *));
855 
856   PetscFunctionBegin;
857   PetscValidHeaderSpecific(mat,MAT_COOKIE,1);
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 PetscErrorCode MatSNESMFSetPeriod(Mat mat,PetscInt 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 PetscErrorCode 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 PetscErrorCode MatSNESMFAddNullSpace(Mat J,MatNullSpace nullsp)
962 {
963   PetscErrorCode 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 PetscErrorCode MatSNESMFSetHHistory(Mat J,PetscScalar history[],PetscInt 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 PetscErrorCode 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 PetscErrorCode MatSNESMFComputeJacobian(SNES snes,Vec x,Mat *jac,Mat *B,MatStructure *flag,void *dummy)
1048 {
1049   PetscErrorCode 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 PetscErrorCode MatSNESMFSetBase(Mat J,Vec U)
1074 {
1075   PetscErrorCode ierr,(*f)(Mat,Vec);
1076 
1077   PetscFunctionBegin;
1078   PetscValidHeaderSpecific(J,MAT_COOKIE,1);
1079   PetscValidHeaderSpecific(U,VEC_COOKIE,2);
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 PetscErrorCode MatSNESMFSetCheckh(Mat J,PetscErrorCode (*fun)(Vec,Vec,PetscScalar*,void*),void* ctx)
1111 {
1112   PetscErrorCode ierr,(*f)(Mat,PetscErrorCode (*)(Vec,Vec,PetscScalar*,void*),void*);
1113 
1114   PetscFunctionBegin;
1115   PetscValidHeaderSpecific(J,MAT_COOKIE,1);
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 PetscErrorCode MatSNESMFCheckPositivity(Vec U,Vec a,PetscScalar *h,void *dummy)
1148 {
1149   PetscReal      val, minval;
1150   PetscScalar    *u_vec, *a_vec;
1151   PetscErrorCode ierr;
1152   PetscInt       i,n;
1153   MPI_Comm       comm;
1154 
1155   PetscFunctionBegin;
1156   ierr = PetscObjectGetComm((PetscObject)U,&comm);CHKERRQ(ierr);
1157   ierr = VecGetArray(U,&u_vec);CHKERRQ(ierr);
1158   ierr = VecGetArray(a,&a_vec);CHKERRQ(ierr);
1159   ierr = VecGetLocalSize(U,&n);CHKERRQ(ierr);
1160   minval = PetscAbsScalar(*h*1.01);
1161   for(i=0;i<n;i++) {
1162     if (PetscRealPart(u_vec[i] + *h*a_vec[i]) <= 0.0) {
1163       val = PetscAbsScalar(u_vec[i]/a_vec[i]);
1164       if (val < minval) minval = val;
1165     }
1166   }
1167   ierr = VecRestoreArray(U,&u_vec);CHKERRQ(ierr);
1168   ierr = VecRestoreArray(a,&a_vec);CHKERRQ(ierr);
1169   ierr = PetscGlobalMin(&minval,&val,comm);CHKERRQ(ierr);
1170   if (val <= PetscAbsScalar(*h)) {
1171     ierr = PetscLogInfo((U,"MatSNESMFCheckPositivity: Scaling back h from %g to %g\n",PetscRealPart(*h),.99*val));CHKERRQ(ierr);
1172     if (PetscRealPart(*h) > 0.0) *h =  0.99*val;
1173     else                         *h = -0.99*val;
1174   }
1175   PetscFunctionReturn(0);
1176 }
1177 
1178 
1179 
1180 
1181 
1182 
1183 
1184