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