xref: /petsc/src/snes/mf/snesmfj.c (revision aea2a34e36487be060e5e9669bb39cba3c3f8e3b)
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   ierr = PetscFListDestroy(&MatSNESMPetscFList);CHKERRQ(ierr);
124   MatSNESMFRegisterAllCalled = PETSC_FALSE;
125   PetscFunctionReturn(0);
126 }
127 
128 /* ----------------------------------------------------------------------------------------*/
129 #undef __FUNCT__
130 #define __FUNCT__ "MatDestroy_MFFD"
131 PetscErrorCode MatDestroy_MFFD(Mat mat)
132 {
133   PetscErrorCode ierr;
134   MatSNESMF      ctx = (MatSNESMF)mat->data;
135 
136   PetscFunctionBegin;
137   if (ctx->w) {
138     ierr = VecDestroy(ctx->w);CHKERRQ(ierr);
139   }
140   if (ctx->ops->destroy) {ierr = (*ctx->ops->destroy)(ctx);CHKERRQ(ierr);}
141   if (ctx->sp) {ierr = MatNullSpaceDestroy(ctx->sp);CHKERRQ(ierr);}
142   ierr = PetscHeaderDestroy(ctx);CHKERRQ(ierr);
143 
144   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatSNESMFSetBase_C","",PETSC_NULL);CHKERRQ(ierr);
145   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatSNESMFSetFunctioniBase_C","",PETSC_NULL);CHKERRQ(ierr);
146   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatSNESMFSetFunctioni_C","",PETSC_NULL);CHKERRQ(ierr);
147   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatSNESMFSetCheckh_C","",PETSC_NULL);CHKERRQ(ierr);
148 
149   PetscFunctionReturn(0);
150 }
151 
152 #undef __FUNCT__
153 #define __FUNCT__ "MatView_MFFD"
154 /*
155    MatSNESMFView_MFFD - Views matrix-free parameters.
156 
157 */
158 PetscErrorCode MatView_MFFD(Mat J,PetscViewer viewer)
159 {
160   PetscErrorCode ierr;
161   MatSNESMF      ctx = (MatSNESMF)J->data;
162   PetscTruth     iascii;
163 
164   PetscFunctionBegin;
165   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
166   if (iascii) {
167      ierr = PetscViewerASCIIPrintf(viewer,"  SNES matrix-free approximation:\n");CHKERRQ(ierr);
168      ierr = PetscViewerASCIIPrintf(viewer,"    err=%G (relative error in function evaluation)\n",ctx->error_rel);CHKERRQ(ierr);
169      if (!ctx->type_name) {
170        ierr = PetscViewerASCIIPrintf(viewer,"    The compute h routine has not yet been set\n");CHKERRQ(ierr);
171      } else {
172        ierr = PetscViewerASCIIPrintf(viewer,"    Using %s compute h routine\n",ctx->type_name);CHKERRQ(ierr);
173      }
174      if (ctx->ops->view) {
175        ierr = (*ctx->ops->view)(ctx,viewer);CHKERRQ(ierr);
176      }
177   } else {
178     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for SNES matrix free matrix",((PetscObject)viewer)->type_name);
179   }
180   PetscFunctionReturn(0);
181 }
182 
183 #undef __FUNCT__
184 #define __FUNCT__ "MatAssemblyEnd_MFFD"
185 /*
186    MatAssemblyEnd_MFFD - Resets the ctx->ncurrenth to zero. This
187    allows the user to indicate the beginning of a new linear solve by calling
188    MatAssemblyXXX() on the matrix free matrix. This then allows the
189    MatSNESMFCreate_WP() to properly compute ||U|| only the first time
190    in the linear solver rather than every time.
191 */
192 PetscErrorCode MatAssemblyEnd_MFFD(Mat J,MatAssemblyType mt)
193 {
194   PetscErrorCode ierr;
195   MatSNESMF      j = (MatSNESMF)J->data;
196 
197   PetscFunctionBegin;
198   ierr = MatSNESMFResetHHistory(J);CHKERRQ(ierr);
199   if (j->usesnes) {
200     ierr = SNESGetSolution(j->snes,&j->current_u);CHKERRQ(ierr);
201     ierr = SNESGetFunction(j->snes,&j->current_f,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
202     if (!j->w) {
203       ierr = VecDuplicate(j->current_u, &j->w);CHKERRQ(ierr);
204     }
205   }
206   j->vshift = 0.0;
207   j->vscale = 1.0;
208   PetscFunctionReturn(0);
209 }
210 
211 #undef __FUNCT__
212 #define __FUNCT__ "MatMult_MFFD"
213 /*
214   MatMult_MFFD - Default matrix-free form for Jacobian-vector product, y = F'(u)*a:
215 
216         y ~= (F(u + ha) - F(u))/h,
217   where F = nonlinear function, as set by SNESSetFunction()
218         u = current iterate
219         h = difference interval
220 */
221 PetscErrorCode MatMult_MFFD(Mat mat,Vec a,Vec y)
222 {
223   MatSNESMF      ctx = (MatSNESMF)mat->data;
224   SNES           snes;
225   PetscScalar    h;
226   Vec            w,U,F;
227   PetscErrorCode ierr,(*eval_fct)(SNES,Vec,Vec)=0;
228   PetscTruth     zeroa;
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(MATSNESMF_Mult,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,&zeroa);CHKERRQ(ierr);
249   if (zeroa) {
250     ierr = VecSet(y,0.0);CHKERRQ(ierr);
251     PetscFunctionReturn(0);
252   }
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 = PetscInfo2(mat,"Current differencing parameter: %G + %G i\n",PetscRealPart(h),PetscImaginaryPart(h));CHKERRQ(ierr);
262 #else
263   ierr = PetscInfo1(mat,"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,-1.0,F);CHKERRQ(ierr);
288   ierr = VecScale(y,1.0/h);CHKERRQ(ierr);
289 
290   ierr = VecAXPBY(y,ctx->vshift,ctx->vscale,a);CHKERRQ(ierr);
291 
292   if (ctx->sp) {ierr = MatNullSpaceRemove(ctx->sp,y,PETSC_NULL);CHKERRQ(ierr);}
293 
294   ierr = PetscLogEventEnd(MATSNESMF_Mult,a,y,0,0);CHKERRQ(ierr);
295   PetscFunctionReturn(0);
296 }
297 
298 #undef __FUNCT__
299 #define __FUNCT__ "MatGetDiagonal_MFFD"
300 /*
301   MatGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix
302 
303         y ~= (F(u + ha) - F(u))/h,
304   where F = nonlinear function, as set by SNESSetFunction()
305         u = current iterate
306         h = difference interval
307 */
308 PetscErrorCode MatGetDiagonal_MFFD(Mat mat,Vec a)
309 {
310   MatSNESMF      ctx = (MatSNESMF)mat->data;
311   PetscScalar    h,*aa,*ww,v;
312   PetscReal      epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON;
313   Vec            w,U;
314   PetscErrorCode ierr;
315   PetscInt       i,rstart,rend;
316 
317   PetscFunctionBegin;
318   if (!ctx->funci) {
319     SETERRQ(PETSC_ERR_ORDER,"Requires calling MatSNESMFSetFunctioni() first");
320   }
321 
322   w    = ctx->w;
323   U    = ctx->current_u;
324   ierr = (*ctx->func)(0,U,a,ctx->funcctx);CHKERRQ(ierr);
325   ierr = (*ctx->funcisetbase)(U,ctx->funcctx);CHKERRQ(ierr);
326   ierr = VecCopy(U,w);CHKERRQ(ierr);
327 
328   ierr = VecGetOwnershipRange(a,&rstart,&rend);CHKERRQ(ierr);
329   ierr = VecGetArray(a,&aa);CHKERRQ(ierr);
330   for (i=rstart; i<rend; i++) {
331     ierr = VecGetArray(w,&ww);CHKERRQ(ierr);
332     h  = ww[i-rstart];
333     if (h == 0.0) h = 1.0;
334 #if !defined(PETSC_USE_COMPLEX)
335     if (h < umin && h >= 0.0)      h = umin;
336     else if (h < 0.0 && h > -umin) h = -umin;
337 #else
338     if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0)     h = umin;
339     else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin;
340 #endif
341     h     *= epsilon;
342 
343     ww[i-rstart] += h;
344     ierr = VecRestoreArray(w,&ww);CHKERRQ(ierr);
345     ierr          = (*ctx->funci)(i,w,&v,ctx->funcctx);CHKERRQ(ierr);
346     aa[i-rstart]  = (v - aa[i-rstart])/h;
347 
348     /* possibly shift and scale result */
349     aa[i - rstart] = ctx->vshift + ctx->vscale*aa[i-rstart];
350 
351     ierr = VecGetArray(w,&ww);CHKERRQ(ierr);
352     ww[i-rstart] -= h;
353     ierr = VecRestoreArray(w,&ww);CHKERRQ(ierr);
354   }
355   ierr = VecRestoreArray(a,&aa);CHKERRQ(ierr);
356   PetscFunctionReturn(0);
357 }
358 
359 #undef __FUNCT__
360 #define __FUNCT__ "MatShift_MFFD"
361 PetscErrorCode MatShift_MFFD(Mat Y,PetscScalar a)
362 {
363   MatSNESMF shell = (MatSNESMF)Y->data;
364   PetscFunctionBegin;
365   shell->vshift += a;
366   PetscFunctionReturn(0);
367 }
368 
369 #undef __FUNCT__
370 #define __FUNCT__ "MatScale_MFFD"
371 PetscErrorCode MatScale_MFFD(Mat Y,PetscScalar a)
372 {
373   MatSNESMF shell = (MatSNESMF)Y->data;
374   PetscFunctionBegin;
375   shell->vscale *= a;
376   PetscFunctionReturn(0);
377 }
378 
379 
380 #undef __FUNCT__
381 #define __FUNCT__ "MatCreateSNESMF"
382 /*@
383    MatCreateSNESMF - Creates a matrix-free matrix context for use with
384    a SNES solver.  This matrix can be used as the Jacobian argument for
385    the routine SNESSetJacobian().
386 
387    Collective on SNES and Vec
388 
389    Input Parameters:
390 +  snes - the SNES context
391 -  x - vector where SNES solution is to be stored.
392 
393    Output Parameter:
394 .  J - the matrix-free matrix
395 
396    Level: advanced
397 
398    Notes:
399    The matrix-free matrix context merely contains the function pointers
400    and work space for performing finite difference approximations of
401    Jacobian-vector products, F'(u)*a,
402 
403    The default code uses the following approach to compute h
404 
405 .vb
406      F'(u)*a = [F(u+h*a) - F(u)]/h where
407      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
408        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
409  where
410      error_rel = square root of relative error in function evaluation
411      umin = minimum iterate parameter
412 .ve
413    (see MATSNESMF_WP or MATSNESMF_DS)
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_type - wp or ds (see MATSNESMF_WP or MATSNESMF_DS)
425 .  -snes_mf_unim <umin> - Sets umin (for default PETSc routine that computes h only)
426 -  -snes_mf_ksp_monitor - KSP monitor routine that prints differencing h
427 
428 .keywords: SNES, default, matrix-free, create, matrix
429 
430 .seealso: MatDestroy(), MatSNESMFSetFunctionError(), MatSNESMFDefaultSetUmin()
431           MatSNESMFSetHHistory(), MatSNESMFResetHHistory(), MatCreateMF(),
432           MatSNESMFGetH(),MatSNESMFKSPMonitor(), MatSNESMFRegisterDynamic), MatSNESMFComputeJacobian()
433 
434 @*/
435 PetscErrorCode PETSCSNES_DLLEXPORT MatCreateSNESMF(SNES snes,Vec x,Mat *J)
436 {
437   MatSNESMF      mfctx;
438   PetscErrorCode ierr;
439 
440   PetscFunctionBegin;
441   ierr = MatCreateMF(x,J);CHKERRQ(ierr);
442 
443   mfctx          = (MatSNESMF)(*J)->data;
444   mfctx->snes    = snes;
445   mfctx->usesnes = PETSC_TRUE;
446   ierr = PetscLogObjectParent(snes,*J);CHKERRQ(ierr);
447   PetscFunctionReturn(0);
448 }
449 
450 EXTERN_C_BEGIN
451 #undef __FUNCT__
452 #define __FUNCT__ "MatSNESMFSetBase_FD"
453 PetscErrorCode PETSCSNES_DLLEXPORT MatSNESMFSetBase_FD(Mat J,Vec U)
454 {
455   PetscErrorCode ierr;
456   MatSNESMF      ctx = (MatSNESMF)J->data;
457 
458   PetscFunctionBegin;
459   ierr = MatSNESMFResetHHistory(J);CHKERRQ(ierr);
460   ctx->current_u = U;
461   ctx->usesnes   = PETSC_FALSE;
462   if (!ctx->w) {
463     ierr = VecDuplicate(ctx->current_u, &ctx->w);CHKERRQ(ierr);
464   }
465   J->assembled = PETSC_TRUE;
466   PetscFunctionReturn(0);
467 }
468 EXTERN_C_END
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   MatSNESMF ctx = (MatSNESMF)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 - wp or ds (see MATSNESMF_WP or MATSNESMF_DS)
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   MatSNESMF      mfctx = (MatSNESMF)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 = KSPMonitorSet(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   MatSNESMF      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_MatSNESMF,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   ierr = PetscObjectChangeTypeName((PetscObject)A,MATMFFD);CHKERRQ(ierr);
615   PetscFunctionReturn(0);
616 }
617 EXTERN_C_END
618 
619 #undef __FUNCT__
620 #define __FUNCT__ "MatCreateMF"
621 /*@
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(), MatSNESMFSetFunction()
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   MatSNESMF ctx = (MatSNESMF)mat->data;
712 
713   PetscFunctionBegin;
714   *h = ctx->currenth;
715   PetscFunctionReturn(0);
716 }
717 
718 #undef __FUNCT__
719 #define __FUNCT__ "MatSNESMFKSPMonitor"
720 /*@C
721    MatSNESMFKSPMonitor - A KSP monitor for use with the  PETSc
722                 SNES matrix free routines. Prints the differencing parameter used at
723                 each step.
724 
725    Collective on KSP
726 
727    Input Parameters:
728 +    ksp - the Krylov solver object
729 .    n  - the current iteration
730 .    rnorm  - the current residual norm (may be preconditioned or not depending on solver and options
731 _    dummy  - unused argument
732 
733   Options Database:
734 .   -snes_mf_ksp_monitor - turn this monitor on
735 
736    Notes: Use KSPMonitorSet(ksp,MatSNESMFKSPMonitor,PETSC_NULL,PETSC_NULL);
737 
738 .seealso:  KSPMonitorSet()
739 
740 @*/
741 PetscErrorCode PETSCSNES_DLLEXPORT MatSNESMFKSPMonitor(KSP ksp,PetscInt n,PetscReal rnorm,void *dummy)
742 {
743   MatSNESMF      ctx;
744   PetscErrorCode ierr;
745   Mat            mat;
746   MPI_Comm       comm;
747   PetscTruth     nonzeroinitialguess;
748 
749   PetscFunctionBegin;
750   ierr = PetscObjectGetComm((PetscObject)ksp,&comm);CHKERRQ(ierr);
751   ierr = KSPGetInitialGuessNonzero(ksp,&nonzeroinitialguess);CHKERRQ(ierr);
752   ierr = KSPGetOperators(ksp,&mat,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
753   ctx  = (MatSNESMF)mat->data;
754 
755   if (n > 0 || nonzeroinitialguess) {
756 #if defined(PETSC_USE_COMPLEX)
757     ierr = PetscPrintf(comm,"%D KSP Residual norm %14.12e h %G + %G i\n",n,rnorm,
758                 PetscRealPart(ctx->currenth),PetscImaginaryPart(ctx->currenth));CHKERRQ(ierr);
759 #else
760     ierr = PetscPrintf(comm,"%D KSP Residual norm %14.12e h %G \n",n,rnorm,ctx->currenth);CHKERRQ(ierr);
761 #endif
762   } else {
763     ierr = PetscPrintf(comm,"%D KSP Residual norm %14.12e\n",n,rnorm);CHKERRQ(ierr);
764   }
765   PetscFunctionReturn(0);
766 }
767 
768 #undef __FUNCT__
769 #define __FUNCT__ "MatSNESMFSetFunction"
770 /*@C
771    MatSNESMFSetFunction - Sets the function used in applying the matrix free.
772 
773    Collective on Mat
774 
775    Input Parameters:
776 +  mat - the matrix free matrix created via MatCreateSNESMF()
777 .  v   - workspace vector
778 .  func - the function to use
779 -  funcctx - optional function context passed to function
780 
781    Level: advanced
782 
783    Notes:
784     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
785     matrix inside your compute Jacobian routine
786 
787     If this is not set then it will use the function set with SNESSetFunction()
788 
789 .keywords: SNES, matrix-free, function
790 
791 .seealso: MatCreateSNESMF(),MatSNESMFGetH(),
792           MatSNESMFSetHHistory(), MatSNESMFResetHHistory(),
793           MatSNESMFKSPMonitor(), SNESetFunction()
794 @*/
795 PetscErrorCode PETSCSNES_DLLEXPORT MatSNESMFSetFunction(Mat mat,Vec v,PetscErrorCode (*func)(SNES,Vec,Vec,void *),void *funcctx)
796 {
797   MatSNESMF ctx = (MatSNESMF)mat->data;
798 
799   PetscFunctionBegin;
800   ctx->func    = func;
801   ctx->funcctx = funcctx;
802   ctx->funcvec = v;
803   PetscFunctionReturn(0);
804 }
805 
806 #undef __FUNCT__
807 #define __FUNCT__ "MatSNESMFSetFunctioni"
808 /*@C
809    MatSNESMFSetFunctioni - Sets the function for a single component
810 
811    Collective on Mat
812 
813    Input Parameters:
814 +  mat - the matrix free matrix created via MatCreateSNESMF()
815 -  funci - the function to use
816 
817    Level: advanced
818 
819    Notes:
820     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
821     matrix inside your compute Jacobian routine
822 
823 
824 .keywords: SNES, matrix-free, function
825 
826 .seealso: MatCreateSNESMF(),MatSNESMFGetH(),
827           MatSNESMFSetHHistory(), MatSNESMFResetHHistory(),
828           MatSNESMFKSPMonitor(), SNESetFunction()
829 @*/
830 PetscErrorCode PETSCSNES_DLLEXPORT MatSNESMFSetFunctioni(Mat mat,PetscErrorCode (*funci)(PetscInt,Vec,PetscScalar*,void *))
831 {
832   PetscErrorCode ierr,(*f)(Mat,PetscErrorCode (*)(PetscInt,Vec,PetscScalar*,void *));
833 
834   PetscFunctionBegin;
835   PetscValidHeaderSpecific(mat,MAT_COOKIE,1);
836   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSNESMFSetFunctioni_C",(void (**)(void))&f);CHKERRQ(ierr);
837   if (f) {
838     ierr = (*f)(mat,funci);CHKERRQ(ierr);
839   }
840   PetscFunctionReturn(0);
841 }
842 
843 
844 #undef __FUNCT__
845 #define __FUNCT__ "MatSNESMFSetFunctioniBase"
846 /*@C
847    MatSNESMFSetFunctioniBase - Sets the base vector for a single component function evaluation
848 
849    Collective on Mat
850 
851    Input Parameters:
852 +  mat - the matrix free matrix created via MatCreateSNESMF()
853 -  func - the function to use
854 
855    Level: advanced
856 
857    Notes:
858     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
859     matrix inside your compute Jacobian routine
860 
861 
862 .keywords: SNES, matrix-free, function
863 
864 .seealso: MatCreateSNESMF(),MatSNESMFGetH(),
865           MatSNESMFSetHHistory(), MatSNESMFResetHHistory(),
866           MatSNESMFKSPMonitor(), SNESetFunction()
867 @*/
868 PetscErrorCode PETSCSNES_DLLEXPORT MatSNESMFSetFunctioniBase(Mat mat,PetscErrorCode (*func)(Vec,void *))
869 {
870   PetscErrorCode ierr,(*f)(Mat,PetscErrorCode (*)(Vec,void *));
871 
872   PetscFunctionBegin;
873   PetscValidHeaderSpecific(mat,MAT_COOKIE,1);
874   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSNESMFSetFunctioniBase_C",(void (**)(void))&f);CHKERRQ(ierr);
875   if (f) {
876     ierr = (*f)(mat,func);CHKERRQ(ierr);
877   }
878   PetscFunctionReturn(0);
879 }
880 
881 
882 #undef __FUNCT__
883 #define __FUNCT__ "MatSNESMFSetPeriod"
884 /*@
885    MatSNESMFSetPeriod - Sets how often h is recomputed, by default it is everytime
886 
887    Collective on Mat
888 
889    Input Parameters:
890 +  mat - the matrix free matrix created via MatCreateSNESMF()
891 -  period - 1 for everytime, 2 for every second etc
892 
893    Options Database Keys:
894 +  -snes_mf_period <period>
895 
896    Level: advanced
897 
898 
899 .keywords: SNES, matrix-free, parameters
900 
901 .seealso: MatCreateSNESMF(),MatSNESMFGetH(),
902           MatSNESMFSetHHistory(), MatSNESMFResetHHistory(),
903           MatSNESMFKSPMonitor()
904 @*/
905 PetscErrorCode PETSCSNES_DLLEXPORT MatSNESMFSetPeriod(Mat mat,PetscInt period)
906 {
907   MatSNESMF ctx = (MatSNESMF)mat->data;
908 
909   PetscFunctionBegin;
910   ctx->recomputeperiod = period;
911   PetscFunctionReturn(0);
912 }
913 
914 #undef __FUNCT__
915 #define __FUNCT__ "MatSNESMFSetFunctionError"
916 /*@
917    MatSNESMFSetFunctionError - Sets the error_rel for the approximation of
918    matrix-vector products using finite differences.
919 
920    Collective on Mat
921 
922    Input Parameters:
923 +  mat - the matrix free matrix created via MatCreateSNESMF()
924 -  error_rel - relative error (should be set to the square root of
925                the relative error in the function evaluations)
926 
927    Options Database Keys:
928 +  -snes_mf_err <error_rel> - Sets error_rel
929 
930    Level: advanced
931 
932    Notes:
933    The default matrix-free matrix-vector product routine computes
934 .vb
935      F'(u)*a = [F(u+h*a) - F(u)]/h where
936      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
937        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   else
938 .ve
939 
940 .keywords: SNES, matrix-free, parameters
941 
942 .seealso: MatCreateSNESMF(),MatSNESMFGetH(),
943           MatSNESMFSetHHistory(), MatSNESMFResetHHistory(),
944           MatSNESMFKSPMonitor()
945 @*/
946 PetscErrorCode PETSCSNES_DLLEXPORT MatSNESMFSetFunctionError(Mat mat,PetscReal error)
947 {
948   MatSNESMF ctx = (MatSNESMF)mat->data;
949 
950   PetscFunctionBegin;
951   if (error != PETSC_DEFAULT) ctx->error_rel = error;
952   PetscFunctionReturn(0);
953 }
954 
955 #undef __FUNCT__
956 #define __FUNCT__ "MatSNESMFAddNullSpace"
957 /*@
958    MatSNESMFAddNullSpace - Provides a null space that an operator is
959    supposed to have.  Since roundoff will create a small component in
960    the null space, if you know the null space you may have it
961    automatically removed.
962 
963    Collective on Mat
964 
965    Input Parameters:
966 +  J - the matrix-free matrix context
967 -  nullsp - object created with MatNullSpaceCreate()
968 
969    Level: advanced
970 
971 .keywords: SNES, matrix-free, null space
972 
973 .seealso: MatNullSpaceCreate(), MatSNESMFGetH(), MatCreateSNESMF(),
974           MatSNESMFSetHHistory(), MatSNESMFResetHHistory(),
975           MatSNESMFKSPMonitor(), MatSNESMFErrorRel()
976 @*/
977 PetscErrorCode PETSCSNES_DLLEXPORT MatSNESMFAddNullSpace(Mat J,MatNullSpace nullsp)
978 {
979   PetscErrorCode ierr;
980   MatSNESMF      ctx = (MatSNESMF)J->data;
981 
982   PetscFunctionBegin;
983   ierr = PetscObjectReference((PetscObject)nullsp);CHKERRQ(ierr);
984   if (ctx->sp) { ierr = MatNullSpaceDestroy(ctx->sp);CHKERRQ(ierr); }
985   ctx->sp = nullsp;
986   PetscFunctionReturn(0);
987 }
988 
989 #undef __FUNCT__
990 #define __FUNCT__ "MatSNESMFSetHHistory"
991 /*@
992    MatSNESMFSetHHistory - Sets an array to collect a history of the
993    differencing values (h) computed for the matrix-free product.
994 
995    Collective on Mat
996 
997    Input Parameters:
998 +  J - the matrix-free matrix context
999 .  histroy - space to hold the history
1000 -  nhistory - number of entries in history, if more entries are generated than
1001               nhistory, then the later ones are discarded
1002 
1003    Level: advanced
1004 
1005    Notes:
1006    Use MatSNESMFResetHHistory() to reset the history counter and collect
1007    a new batch of differencing parameters, h.
1008 
1009 .keywords: SNES, matrix-free, h history, differencing history
1010 
1011 .seealso: MatSNESMFGetH(), MatCreateSNESMF(),
1012           MatSNESMFResetHHistory(),
1013           MatSNESMFKSPMonitor(), MatSNESMFSetFunctionError()
1014 
1015 @*/
1016 PetscErrorCode PETSCSNES_DLLEXPORT MatSNESMFSetHHistory(Mat J,PetscScalar history[],PetscInt nhistory)
1017 {
1018   MatSNESMF ctx = (MatSNESMF)J->data;
1019 
1020   PetscFunctionBegin;
1021   ctx->historyh    = history;
1022   ctx->maxcurrenth = nhistory;
1023   ctx->currenth    = 0;
1024   PetscFunctionReturn(0);
1025 }
1026 
1027 #undef __FUNCT__
1028 #define __FUNCT__ "MatSNESMFResetHHistory"
1029 /*@
1030    MatSNESMFResetHHistory - Resets the counter to zero to begin
1031    collecting a new set of differencing histories.
1032 
1033    Collective on Mat
1034 
1035    Input Parameters:
1036 .  J - the matrix-free matrix context
1037 
1038    Level: advanced
1039 
1040    Notes:
1041    Use MatSNESMFSetHHistory() to create the original history counter.
1042 
1043 .keywords: SNES, matrix-free, h history, differencing history
1044 
1045 .seealso: MatSNESMFGetH(), MatCreateSNESMF(),
1046           MatSNESMFSetHHistory(),
1047           MatSNESMFKSPMonitor(), MatSNESMFSetFunctionError()
1048 
1049 @*/
1050 PetscErrorCode PETSCSNES_DLLEXPORT MatSNESMFResetHHistory(Mat J)
1051 {
1052   MatSNESMF ctx = (MatSNESMF)J->data;
1053 
1054   PetscFunctionBegin;
1055   ctx->ncurrenth    = 0;
1056   PetscFunctionReturn(0);
1057 }
1058 
1059 #undef __FUNCT__
1060 #define __FUNCT__ "MatSNESMFComputeJacobian"
1061 /*@
1062    MatSNESMFComputeJacobian - Tells the matrix-free Jacobian object the new location at which
1063        Jacobian matrix vector products will be computed at, i.e. J(x) * a.
1064 
1065    Collective on SNES
1066 
1067    Input Parameters:
1068 +   snes - the nonlinear solver context
1069 .   x - the point at which the Jacobian vector products will be performed
1070 .   jac - the matrix-free Jacobian object
1071 .   B - either the same as jac or another matrix type (ignored)
1072 .   flag - not relevent for matrix-free form
1073 -   dummy - the user context (ignored)
1074 
1075    Level: developer
1076 
1077    Notes:
1078      This can be passed into SNESSetJacobian() when using a completely matrix-free solver,
1079      that is the B matrix is also the same matrix operator. This is used when you select
1080      -snes_mf but rarely used directly by users.
1081 
1082 .seealso: MatSNESMFGetH(), MatCreateSNESMF(),
1083           MatSNESMFSetHHistory(),
1084           MatSNESMFKSPMonitor(), MatSNESMFSetFunctionError(), MatSNESMFCreate(), SNESSetJacobian()
1085 
1086 @*/
1087 PetscErrorCode PETSCSNES_DLLEXPORT MatSNESMFComputeJacobian(SNES snes,Vec x,Mat *jac,Mat *B,MatStructure *flag,void *dummy)
1088 {
1089   PetscErrorCode ierr;
1090   PetscFunctionBegin;
1091   ierr = MatAssemblyBegin(*jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1092   ierr = MatAssemblyEnd(*jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1093   PetscFunctionReturn(0);
1094 }
1095 
1096 #undef __FUNCT__
1097 #define __FUNCT__ "MatSNESMFSetBase"
1098 /*@
1099     MatSNESMFSetBase - Sets the vector U at which matrix vector products of the
1100         Jacobian are computed
1101 
1102     Collective on Mat
1103 
1104     Input Parameters:
1105 +   J - the MatSNESMF matrix
1106 -   U - the vector
1107 
1108     Notes: This is rarely used directly
1109 
1110     Level: advanced
1111 
1112 @*/
1113 PetscErrorCode PETSCSNES_DLLEXPORT MatSNESMFSetBase(Mat J,Vec U)
1114 {
1115   PetscErrorCode ierr,(*f)(Mat,Vec);
1116 
1117   PetscFunctionBegin;
1118   PetscValidHeaderSpecific(J,MAT_COOKIE,1);
1119   PetscValidHeaderSpecific(U,VEC_COOKIE,2);
1120   ierr = PetscObjectQueryFunction((PetscObject)J,"MatSNESMFSetBase_C",(void (**)(void))&f);CHKERRQ(ierr);
1121   if (f) {
1122     ierr = (*f)(J,U);CHKERRQ(ierr);
1123   }
1124   PetscFunctionReturn(0);
1125 }
1126 
1127 #undef __FUNCT__
1128 #define __FUNCT__ "MatSNESMFSetCheckh"
1129 /*@C
1130     MatSNESMFSetCheckh - Sets a function that checks the computed h and adjusts
1131         it to satisfy some criteria
1132 
1133     Collective on Mat
1134 
1135     Input Parameters:
1136 +   J - the MatSNESMF matrix
1137 .   fun - the function that checks h
1138 -   ctx - any context needed by the function
1139 
1140     Options Database Keys:
1141 .   -snes_mf_check_positivity
1142 
1143     Level: advanced
1144 
1145     Notes: For example, MatSNESMFSetCheckPositivity() insures that all entries
1146        of U + h*a are non-negative
1147 
1148 .seealso:  MatSNESMFSetCheckPositivity()
1149 @*/
1150 PetscErrorCode PETSCSNES_DLLEXPORT MatSNESMFSetCheckh(Mat J,PetscErrorCode (*fun)(Vec,Vec,PetscScalar*,void*),void* ctx)
1151 {
1152   PetscErrorCode ierr,(*f)(Mat,PetscErrorCode (*)(Vec,Vec,PetscScalar*,void*),void*);
1153 
1154   PetscFunctionBegin;
1155   PetscValidHeaderSpecific(J,MAT_COOKIE,1);
1156   ierr = PetscObjectQueryFunction((PetscObject)J,"MatSNESMFSetCheckh_C",(void (**)(void))&f);CHKERRQ(ierr);
1157   if (f) {
1158     ierr = (*f)(J,fun,ctx);CHKERRQ(ierr);
1159   }
1160   PetscFunctionReturn(0);
1161 }
1162 
1163 #undef __FUNCT__
1164 #define __FUNCT__ "MatSNESMFSetCheckPositivity"
1165 /*@
1166     MatSNESMFCheckPositivity - Checks that all entries in U + h*a are positive or
1167         zero, decreases h until this is satisfied.
1168 
1169     Collective on Vec
1170 
1171     Input Parameters:
1172 +   U - base vector that is added to
1173 .   a - vector that is added
1174 .   h - scaling factor on a
1175 -   dummy - context variable (unused)
1176 
1177     Options Database Keys:
1178 .   -snes_mf_check_positivity
1179 
1180     Level: advanced
1181 
1182     Notes: This is rarely used directly, rather it is passed as an argument to
1183            MatSNESMFSetCheckh()
1184 
1185 .seealso:  MatSNESMFSetCheckh()
1186 @*/
1187 PetscErrorCode PETSCSNES_DLLEXPORT MatSNESMFCheckPositivity(Vec U,Vec a,PetscScalar *h,void *dummy)
1188 {
1189   PetscReal      val, minval;
1190   PetscScalar    *u_vec, *a_vec;
1191   PetscErrorCode ierr;
1192   PetscInt       i,n;
1193   MPI_Comm       comm;
1194 
1195   PetscFunctionBegin;
1196   ierr = PetscObjectGetComm((PetscObject)U,&comm);CHKERRQ(ierr);
1197   ierr = VecGetArray(U,&u_vec);CHKERRQ(ierr);
1198   ierr = VecGetArray(a,&a_vec);CHKERRQ(ierr);
1199   ierr = VecGetLocalSize(U,&n);CHKERRQ(ierr);
1200   minval = PetscAbsScalar(*h*1.01);
1201   for(i=0;i<n;i++) {
1202     if (PetscRealPart(u_vec[i] + *h*a_vec[i]) <= 0.0) {
1203       val = PetscAbsScalar(u_vec[i]/a_vec[i]);
1204       if (val < minval) minval = val;
1205     }
1206   }
1207   ierr = VecRestoreArray(U,&u_vec);CHKERRQ(ierr);
1208   ierr = VecRestoreArray(a,&a_vec);CHKERRQ(ierr);
1209   ierr = PetscGlobalMin(&minval,&val,comm);CHKERRQ(ierr);
1210   if (val <= PetscAbsScalar(*h)) {
1211     ierr = PetscInfo2(U,"Scaling back h from %G to %G\n",PetscRealPart(*h),.99*val);CHKERRQ(ierr);
1212     if (PetscRealPart(*h) > 0.0) *h =  0.99*val;
1213     else                         *h = -0.99*val;
1214   }
1215   PetscFunctionReturn(0);
1216 }
1217 
1218 
1219 
1220 
1221 
1222 
1223 
1224