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