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