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