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