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