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