xref: /petsc/src/snes/mf/snesmfj.c (revision 4cfa41b92b75d5220e24ab129c849dd026e4e576)
1 
2 #ifdef PETSC_RCS_HEADER
3 static char vcid[] = "$Id: snesmfj.c,v 1.75 1998/12/03 04:05:28 bsmith Exp balay $";
4 #endif
5 
6 #include "src/snes/snesimpl.h"
7 #include "src/snes/mf/snesmfj.h"   /*I  "snes.h"   I*/
8 
9 FList MatSNESFDMFList              = 0;
10 int   MatSNESFDMFRegisterAllCalled = 0;
11 
12 
13 #undef __FUNC__
14 #define __FUNC__ "MatSNESFDMFSetType"
15 /*@
16       MatSNESFDMFSetType - Sets the method that is used to compute the h in the
17             finite difference matrix free formulation.
18 
19    Input Parameters:
20 +     mat - the matrix free matrix created via MatCreateSNESFDMF()
21 -     ftype - the type requested
22 
23 .seealso: MatCreateSNESFDMF(), MatSNESFDMFRegister()
24 @*/
25 int MatSNESFDMFSetType(Mat mat,char *ftype)
26 {
27   int            ierr, (*r)(MatSNESFDMFCtx);
28   MatSNESFDMFCtx ctx;
29 
30   PetscFunctionBegin;
31   ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr);
32 
33   /* already set, so just return */
34   if (!PetscStrcmp(ctx->type_name,ftype)) PetscFunctionReturn(0);
35 
36   /* destroy the old one if it exists */
37   if (ctx->ops->destroy) {
38     ierr = (*ctx->ops->destroy)(ctx);CHKERRQ(ierr);
39   }
40 
41   /* Get the function pointers for the iterative method requested */
42   if (!MatSNESFDMFRegisterAllCalled) {ierr = MatSNESFDMFRegisterAll(PETSC_NULL); CHKERRQ(ierr);}
43 
44   ierr =  FListFind(ctx->comm, MatSNESFDMFList, ftype,(int (**)(void *)) &r );CHKERRQ(ierr);
45 
46   if (!r) SETERRQ(1,1,"Unknown MatSNESFDMF type given");
47 
48   ierr = (*r)(ctx); CHKERRQ(ierr);
49   ierr = PetscStrncpy(ctx->type_name,ftype,256);CHKERRQ(ierr);
50 
51   PetscFunctionReturn(0);
52 }
53 
54 /*MC
55    MatSNESFDMFRegister - Adds a method to the MatSNESFDMF registry
56 
57    Synopsis:
58    MatSNESFDMFRegister(char *name_solver,char *path,char *name_create,int (*routine_create)(MatSNESFDMF))
59 
60    Not Collective
61 
62    Input Parameters:
63 +  name_solver - name of a new user-defined compute-h module
64 .  path - path (either absolute or relative) the library containing this solver
65 .  name_create - name of routine to create method context
66 -  routine_create - routine to create method context
67 
68    Notes:
69    MatSNESFDMFRegister() may be called multiple times to add several user-defined solvers.
70 
71    If dynamic libraries are used, then the fourth input argument (routine_create)
72    is ignored.
73 
74    Sample usage:
75 .vb
76    MatSNESFDMFRegister("my_h",/home/username/my_lib/lib/libO/solaris/mylib.a,
77                "MyHCreate",MyHCreate);
78 .ve
79 
80    Then, your solver can be chosen with the procedural interface via
81 $     MatSNESFDMFSetType(mfctx,"my_h")
82    or at runtime via the option
83 $     -snes_mf_type my_h
84 
85 .keywords: MatSNESFDMF, register
86 
87 .seealso: MatSNESFDMFRegisterAll(), MatSNESFDMFRegisterDestroy()
88 M*/
89 
90 #undef __FUNC__
91 #define __FUNC__ "MatSNESFDMFRegister_Private"
92 int MatSNESFDMFRegister_Private(char *sname,char *path,char *name,int (*function)(MatSNESFDMFCtx))
93 {
94   int ierr;
95   char fullname[256];
96 
97   PetscFunctionBegin;
98   PetscStrcpy(fullname,path); PetscStrcat(fullname,":");PetscStrcat(fullname,name);
99   ierr = FListAdd_Private(&MatSNESFDMFList,sname,fullname,(int (*)(void*))function);CHKERRQ(ierr);
100   PetscFunctionReturn(0);
101 }
102 
103 
104 #undef __FUNC__
105 #define __FUNC__ "MatSNESFDMFRegisterDestroy"
106 /*@C
107    MatSNESFDMFRegisterDestroy - Frees the list of MatSNESFDMF methods that were
108    registered by MatSNESFDMFRegister().
109 
110    Not Collective
111 
112 .keywords: MatSNESFDMF, register, destroy
113 
114 .seealso: MatSNESFDMFRegister(), MatSNESFDMFRegisterAll()
115 @*/
116 int MatSNESFDMFRegisterDestroy(void)
117 {
118   int ierr;
119 
120   PetscFunctionBegin;
121   if (MatSNESFDMFList) {
122     ierr = FListDestroy( MatSNESFDMFList );CHKERRQ(ierr);
123     MatSNESFDMFList = 0;
124   }
125   MatSNESFDMFRegisterAllCalled = 0;
126   PetscFunctionReturn(0);
127 }
128 
129 /* ----------------------------------------------------------------------------------------*/
130 #undef __FUNC__
131 #define __FUNC__ "MatSNESFDMFDestroy_Private"
132 int MatSNESFDMFDestroy_Private(Mat mat)
133 {
134   int            ierr;
135   MatSNESFDMFCtx ctx;
136 
137   PetscFunctionBegin;
138   ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr);
139   ierr = VecDestroy(ctx->w); CHKERRQ(ierr);
140   if (ctx->ops->destroy) {ierr = (*ctx->ops->destroy)(ctx);CHKERRQ(ierr);}
141   if (ctx->sp) {ierr = PCNullSpaceDestroy(ctx->sp);CHKERRQ(ierr);}
142   PetscFree(ctx->ops);
143   PetscFree(ctx);
144   PetscFunctionReturn(0);
145 }
146 
147 #undef __FUNC__
148 #define __FUNC__ "MatSNESFDMFView_Private"
149 /*
150    MatSNESFDMFView_Private - Views matrix-free parameters.
151 
152 */
153 int MatSNESFDMFView_Private(Mat J,Viewer viewer)
154 {
155   int            ierr;
156   MatSNESFDMFCtx ctx;
157   MPI_Comm       comm;
158   FILE           *fd;
159   ViewerType     vtype;
160 
161   PetscFunctionBegin;
162   PetscObjectGetComm((PetscObject)J,&comm);
163   ierr = MatShellGetContext(J,(void **)&ctx); CHKERRQ(ierr);
164   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
165   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
166   if (!PetscStrcmp(vtype,ASCII_VIEWER)) {
167      PetscFPrintf(comm,fd,"  SNES matrix-free approximation:\n");
168      PetscFPrintf(comm,fd,"    err=%g (relative error in function evaluation)\n",ctx->error_rel);
169      PetscFPrintf(ctx->comm,fd,"    Using %s compute h routine\n",ctx->type_name);
170      if (ctx->ops->view) {
171        ierr = (*ctx->ops->view)(ctx,viewer);CHKERRQ(ierr);
172      }
173   } else {
174     SETERRQ(1,1,"Viewer type not supported for this object");
175   }
176   PetscFunctionReturn(0);
177 }
178 
179 #undef __FUNC__
180 #define __FUNC__ "MatSNESFDMFAssemblyEnd_Private"
181 /*
182    MatSNESFDMFAssemblyEnd_Private - Resets the ctx->ncurrenth to zero. This
183      allows the user to indicate the beginning of a new linear solve by call
184      MatAssemblyXXX() on the matrix free matrix. This then allows the
185      MatSNESFDMFCreate_WP() to properly compute the || U|| only the first
186      time in the linear solver rather than every time
187 
188 */
189 int MatSNESFDMFAssemblyEnd_Private(Mat J)
190 {
191   int            ierr;
192 
193   PetscFunctionBegin;
194   ierr = MatSNESFDMFResetHHistory(J);CHKERRQ(ierr);
195   PetscFunctionReturn(0);
196 }
197 
198 
199 #undef __FUNC__
200 #define __FUNC__ "MatSNESFDMFMult_Private"
201 /*
202   MatSNESFDMFMult_Private - Default matrix-free form for Jacobian-vector
203   product, y = F'(u)*a:
204 
205         y ~= ( F(u + ha) - F(u) )/h,
206   where F = nonlinear function, as set by SNESSetFunction()
207         u = current iterate
208         h = difference interval
209 */
210 int MatSNESFDMFMult_Private(Mat mat,Vec a,Vec y)
211 {
212   MatSNESFDMFCtx ctx;
213   SNES           snes;
214   Scalar         h, mone = -1.0;
215   Vec            w,U,F;
216   int            ierr, (*eval_fct)(SNES,Vec,Vec)=0;
217 
218   PetscFunctionBegin;
219   /* We log matrix-free matrix-vector products separately, so that we can
220      separate the performance monitoring from the cases that use conventional
221      storage.  We may eventually modify event logging to associate events
222      with particular objects, hence alleviating the more general problem. */
223   PLogEventBegin(MAT_MatrixFreeMult,a,y,0,0);
224 
225   ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr);
226   snes = ctx->snes;
227   w    = ctx->w;
228 
229   ierr = SNESGetSolution(snes,&U); CHKERRQ(ierr);
230   if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
231     eval_fct = SNESComputeFunction;
232     ierr = SNESGetFunction(snes,&F); CHKERRQ(ierr);
233   } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
234     eval_fct = SNESComputeGradient;
235     ierr = SNESGetGradient(snes,&F); CHKERRQ(ierr);
236   } else SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Invalid method class");
237 
238   if (!ctx->ops->compute) {
239     ierr = MatSNESFDMFSetType(mat,"default");CHKERRQ(ierr);
240     ierr = MatSNESFDMFSetFromOptions(mat);CHKERRQ(ierr);
241   }
242   ierr = (*ctx->ops->compute)(ctx,U,a,&h); CHKERRQ(ierr);
243 
244   /* keep a record of the current differencing parameter h */
245   ctx->currenth = h;
246 #if defined(USE_PETSC_COMPLEX)
247   PLogInfo(mat,"Current differencing parameter: %g + %g i\n",PetscReal(h),PetscImaginary(h));
248 #else
249   PLogInfo(mat,"Current differencing parameter: %g\n",h);
250 #endif
251   if (ctx->historyh && ctx->ncurrenth < ctx->maxcurrenth) {
252     ctx->historyh[ctx->ncurrenth++] = h;
253   } else {
254     ctx->ncurrenth++;
255   }
256 
257   /* Evaluate function at F(u + ha) */
258   ierr = VecWAXPY(&h,a,U,w); CHKERRQ(ierr);
259   ierr = eval_fct(snes,w,y); CHKERRQ(ierr);
260 
261   ierr = VecAXPY(&mone,F,y); CHKERRQ(ierr);
262   h = 1.0/h;
263   ierr = VecScale(&h,y); CHKERRQ(ierr);
264   if (ctx->sp) {ierr = PCNullSpaceRemove(ctx->sp,y); CHKERRQ(ierr);}
265 
266   PLogEventEnd(MAT_MatrixFreeMult,a,y,0,0);
267   PetscFunctionReturn(0);
268 }
269 
270 #undef __FUNC__
271 #define __FUNC__ "MatCreateSNESFDMF"
272 /*@C
273    MatCreateSNESFDMF - Creates a matrix-free matrix
274    context for use with a SNES solver.  This matrix can be used as
275    the Jacobian argument for the routine SNESSetJacobian().
276 
277    Collective on SNES and Vec
278 
279    Input Parameters:
280 +  snes - the SNES context
281 -  x - vector where SNES solution is to be stored.
282 
283    Output Parameter:
284 .  J - the matrix-free matrix
285 
286    Notes:
287    The matrix-free matrix context merely contains the function pointers
288    and work space for performing finite difference approximations of
289    Jacobian-vector products, J(u)*a,
290 
291    The default code uses the following approach to compute h
292 
293 .vb
294      J(u)*a = [J(u+h*a) - J(u)]/h where
295      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
296        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
297  where
298      error_rel = square root of relative error in function evaluation
299      umin = minimum iterate parameter
300 .ve
301 
302    The user can set the error_rel via MatSNESFDMFSetFunctionError() and
303    umin via MatSNESFDMFDefaultSetUmin()
304    See the nonlinear solvers chapter of the users manual for details.
305 
306    The user should call MatDestroy() when finished with the matrix-free
307    matrix context.
308 
309    Options Database Keys:
310 +  -snes_mf_err <error_rel> - Sets error_rel
311 .  -snes_mf_unim <umin> - Sets umin (for default PETSc routine that computes h only)
312 -  -snes_mf_ksp_monitor - KSP monitor routine that prints differencing h
313 
314 .keywords: SNES, default, matrix-free, create, matrix
315 
316 .seealso: MatDestroy(), MatSNESFDMFSetFunctionError(), MatSNESFDMFDefaultSetUmin()
317           MatSNESFDMFSetHHistory(), MatSNESFDMFResetHHistory(),
318           MatSNESFDMFGetH(),MatSNESFDMFKSPMonitor(), MatSNESFDMFRegister()
319 
320 @*/
321 int MatCreateSNESFDMF(SNES snes,Vec x, Mat *J)
322 {
323   MPI_Comm       comm;
324   MatSNESFDMFCtx mfctx;
325   int            n, nloc, ierr;
326 
327   PetscFunctionBegin;
328   mfctx = (MatSNESFDMFCtx) PetscMalloc(sizeof(struct _p_MatSNESFDMFCtx)); CHKPTRQ(mfctx);
329   PLogObjectMemory(snes,sizeof(MatSNESFDMFCtx));
330   mfctx->comm         = snes->comm;
331   mfctx->sp           = 0;
332   mfctx->snes         = snes;
333   mfctx->error_rel    = 1.e-8; /* assumes double precision */
334   mfctx->currenth     = 0.0;
335   mfctx->historyh     = PETSC_NULL;
336   mfctx->ncurrenth    = 0;
337   mfctx->maxcurrenth  = 0;
338 
339   /*
340      Create the empty data structure to contain compute-h routines.
341      These will be filled in below from the command line options or
342      a later call with MatSNESFDMFSetType() or if that is not called
343      then it will default in the first use of MatSNESFDMFMult_private()
344   */
345   mfctx->ops                 = (MFOps *)PetscMalloc(sizeof(MFOps)); CHKPTRQ(mfctx->ops);
346   mfctx->ops->compute        = 0;
347   mfctx->ops->destroy        = 0;
348   mfctx->ops->view           = 0;
349   mfctx->ops->printhelp      = 0;
350   mfctx->ops->setfromoptions = 0;
351   mfctx->hctx                = 0;
352 
353   ierr = VecDuplicate(x,&mfctx->w); CHKERRQ(ierr);
354   ierr = PetscObjectGetComm((PetscObject)x,&comm); CHKERRQ(ierr);
355   ierr = VecGetSize(x,&n); CHKERRQ(ierr);
356   ierr = VecGetLocalSize(x,&nloc); CHKERRQ(ierr);
357   ierr = MatCreateShell(comm,nloc,n,n,n,mfctx,J); CHKERRQ(ierr);
358   ierr = MatShellSetOperation(*J,MATOP_MULT,(void*)MatSNESFDMFMult_Private);CHKERRQ(ierr);
359   ierr = MatShellSetOperation(*J,MATOP_DESTROY,(void *)MatSNESFDMFDestroy_Private);CHKERRQ(ierr);
360   ierr = MatShellSetOperation(*J,MATOP_VIEW,(void *)MatSNESFDMFView_Private); CHKERRQ(ierr);
361   ierr = MatShellSetOperation(*J,MATOP_ASSEMBLY_END,(void *)MatSNESFDMFAssemblyEnd_Private);CHKERRQ(ierr);
362   PLogObjectParent(*J,mfctx->w);
363   PLogObjectParent(snes,*J);
364 
365   mfctx->mat = *J;
366 
367 
368   PetscFunctionReturn(0);
369 }
370 
371 #undef __FUNC__
372 #define __FUNC__ "MatSNESFDMFGetH"
373 /*@
374    MatSNESFDMFSetFromOptions - Sets the MatSNESFDMF options from the command line
375      parameter.
376 
377    Collective on Mat
378 
379    Input Parameters:
380 .   mat - the matrix obtained with MatCreateSNESFDMF()
381 
382 .keywords: SNES, matrix-free, parameters
383 
384 .seealso: MatCreateSNESFDMF(),MatSNESFDMFSetHHistory(),
385           MatSNESFDMFResetHHistory(),MatSNESFDMFKSPMonitor()
386 @*/
387 int MatSNESFDMFSetFromOptions(Mat mat)
388 {
389   MatSNESFDMFCtx mfctx;
390   int            ierr,flg;
391   char           ftype[256],p[64];
392 
393   PetscFunctionBegin;
394   ierr = MatShellGetContext(mat,(void **)&mfctx); CHKERRQ(ierr);
395   if (mfctx) {
396     /* allow user to set the type */
397     ierr = OptionsGetString(mfctx->snes->prefix,"-snes_mf_type",ftype,256,&flg);CHKERRQ(ierr);
398     if (flg) {
399       ierr = MatSNESFDMFSetType(mat,ftype);CHKERRQ(ierr);
400     }
401 
402     ierr = OptionsGetDouble(mfctx->snes->prefix,"-snes_mf_err",&mfctx->error_rel,&flg);CHKERRQ(ierr);
403     if (mfctx->ops->setfromoptions) {
404       ierr = (*mfctx->ops->setfromoptions)(mfctx);CHKERRQ(ierr);
405     }
406 
407     ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr);
408     PetscStrcpy(p,"-");
409     if (mfctx->snes->prefix) PetscStrcat(p,mfctx->snes->prefix);
410     if (flg) {
411       (*PetscHelpPrintf)(mfctx->snes->comm,"   %ssnes_mf_err <err>: set sqrt rel error in function (default %g)\n",p,mfctx->error_rel);
412       if (mfctx->ops->printhelp) {
413         (*mfctx->ops->printhelp)(mfctx);CHKERRQ(ierr);
414       }
415     }
416   }
417   PetscFunctionReturn(0);
418 }
419 
420 #undef __FUNC__
421 #define __FUNC__ "MatSNESFDMFGetH"
422 /*@
423    MatSNESFDMFGetH - Gets the last h that was used as the differencing
424      parameter.
425 
426    Not Collective
427 
428    Input Parameters:
429 .   mat - the matrix obtained with MatCreateSNESFDMF()
430 
431    Output Paramter:
432 .  h - the differencing step size
433 
434 .keywords: SNES, matrix-free, parameters
435 
436 .seealso: MatCreateSNESFDMF(),MatSNESFDMFSetHHistory(),
437           MatSNESFDMFResetHHistory(),MatSNESFDMFKSPMonitor()
438 @*/
439 int MatSNESFDMFGetH(Mat mat,Scalar *h)
440 {
441   MatSNESFDMFCtx ctx;
442   int            ierr;
443 
444   PetscFunctionBegin;
445   ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr);
446   if (ctx) {
447     *h = ctx->currenth;
448   }
449   PetscFunctionReturn(0);
450 }
451 
452 #undef __FUNC__
453 #define __FUNC__ "MatSNESFDMFKSPMonitor"
454 /*
455    MatSNESFDMFKSPMonitor - A KSP monitor for use with the default PETSc
456       SNES matrix free routines. Prints the h differencing parameter used at each
457       timestep.
458 
459 */
460 int MatSNESFDMFKSPMonitor(KSP ksp,int n,double rnorm,void *dummy)
461 {
462   PC             pc;
463   MatSNESFDMFCtx ctx;
464   int            ierr;
465   Mat            mat;
466   MPI_Comm       comm;
467   PetscTruth     nonzeroinitialguess;
468 
469   PetscFunctionBegin;
470   ierr = PetscObjectGetComm((PetscObject)ksp,&comm);CHKERRQ(ierr);
471   ierr = KSPGetPC(ksp,&pc); CHKERRQ(ierr);
472   ierr = KSPGetInitialGuessNonzero(ksp,&nonzeroinitialguess);CHKERRQ(ierr);
473   ierr = PCGetOperators(pc,&mat,PETSC_NULL,PETSC_NULL); CHKERRQ(ierr);
474   ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr);
475   if (!ctx) {
476     SETERRQ(1,1,"Matrix is not a matrix free shell matrix");
477   }
478   if (n > 0 || nonzeroinitialguess) {
479 #if defined(USE_PETSC_COMPLEX)
480     PetscPrintf(comm,"%d KSP Residual norm %14.12e h %g + %g i\n",n,rnorm,
481                 PetscReal(ctx->currenth),PetscImaginary(ctx->currenth));
482 #else
483     PetscPrintf(comm,"%d KSP Residual norm %14.12e h %g \n",n,rnorm,ctx->currenth);
484 #endif
485   } else {
486     PetscPrintf(comm,"%d KSP Residual norm %14.12e\n",n,rnorm);
487   }
488   PetscFunctionReturn(0);
489 }
490 
491 #undef __FUNC__
492 #define __FUNC__ "MatSNESFDMFSetFunctionError"
493 /*@
494    MatSNESFDMFSetFunctionError - Sets the error_rel for the approximation of
495    matrix-vector products using finite differences.
496 
497    Collective on Mat
498 
499    Input Parameters:
500 +  mat - the matrix free matrix created via MatCreateSNESFDMF()
501 -  error_rel - relative error (should be set to the square root of
502                the relative error in the function evaluations)
503 
504    Notes:
505    The default matrix-free matrix-vector product routine computes
506 .vb
507      J(u)*a = [J(u+h*a) - J(u)]/h where
508      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
509        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   else
510 .ve
511 
512    Options Database Keys:
513 +  -snes_mf_err <error_rel> - Sets error_rel
514 
515 .keywords: SNES, matrix-free, parameters
516 
517 .seealso: MatCreateSNESFDMF(),MatSNESFDMFGetH(),
518           MatSNESFDMFSetHHistory(), MatSNESFDMFResetHHistory(),
519           MatSNESFDMFKSPMonitor()
520 @*/
521 int MatSNESFDMFSetFunctionError(Mat mat,double error)
522 {
523   MatSNESFDMFCtx ctx;
524   int            ierr;
525 
526   PetscFunctionBegin;
527   ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr);
528   if (ctx) {
529     if (error != PETSC_DEFAULT) ctx->error_rel = error;
530   }
531   PetscFunctionReturn(0);
532 }
533 
534 #undef __FUNC__
535 #define __FUNC__ "MatSNESFDMFAddNullSpace"
536 /*@
537    MatSNESFDMFAddNullSpace - Provides a null space that
538    an operator is supposed to have.  Since roundoff will create a
539    small component in the null space, if you know the null space
540    you may have it automatically removed.
541 
542    Collective on Mat
543 
544    Input Parameters:
545 +  J - the matrix-free matrix context
546 .  has_cnst - PETSC_TRUE or PETSC_FALSE, indicating if null space has constants
547 .  n - number of vectors (excluding constant vector) in null space
548 -  vecs - the vectors that span the null space (excluding the constant vector);
549           these vectors must be orthonormal
550 
551 .keywords: SNES, matrix-free, null space
552 
553 .seealso: MatSNESFDMFGetH(), MatCreateSNESFDMF(),
554           MatSNESFDMFSetHHistory(), MatSNESFDMFResetHHistory(),
555           MatSNESFDMFKSPMonitor(), MatSNESFDMFErrorRel()
556 
557 @*/
558 int MatSNESFDMFAddNullSpace(Mat J,int has_cnst,int n,Vec *vecs)
559 {
560   int            ierr;
561   MatSNESFDMFCtx ctx;
562   MPI_Comm       comm;
563 
564   PetscFunctionBegin;
565   PetscObjectGetComm((PetscObject)J,&comm);
566 
567   ierr = MatShellGetContext(J,(void **)&ctx); CHKERRQ(ierr);
568   /* no context indicates that it is not the "matrix free" matrix type */
569   if (!ctx) PetscFunctionReturn(0);
570   ierr = PCNullSpaceCreate(comm,has_cnst,n,vecs,&ctx->sp); CHKERRQ(ierr);
571   PetscFunctionReturn(0);
572 }
573 
574 #undef __FUNC__
575 #define __FUNC__ "MatSNESFDMFSetHHistory"
576 /*@
577    MatSNESFDMFSetHHistory - Sets an array to collect a history
578       of the differencing values h computed for the matrix free product
579 
580    Collective on Mat
581 
582    Input Parameters:
583 +  J - the matrix-free matrix context
584 .  histroy - space to hold the h history
585 -  nhistory - number of entries in history, if more h are generated than
586               nhistory the later ones are discarded
587 
588    Notes:
589     Use MatSNESFDMFResetHHistory() to reset the history counter
590     and collect a new batch of h.
591 
592 .keywords: SNES, matrix-free, h history, differencing history
593 
594 .seealso: MatSNESFDMFGetH(), MatCreateSNESFDMF(),
595           MatSNESFDMFResetHHistory(),
596           MatSNESFDMFKSPMonitor(), MatSNESFDMFSetFunctionError()
597 
598 @*/
599 int MatSNESFDMFSetHHistory(Mat J,Scalar *history,int nhistory)
600 {
601   int            ierr;
602   MatSNESFDMFCtx ctx;
603 
604   PetscFunctionBegin;
605 
606   ierr = MatShellGetContext(J,(void **)&ctx); CHKERRQ(ierr);
607   /* no context indicates that it is not the "matrix free" matrix type */
608   if (!ctx) PetscFunctionReturn(0);
609   ctx->historyh    = history;
610   ctx->maxcurrenth = nhistory;
611   ctx->currenth    = 0;
612 
613   PetscFunctionReturn(0);
614 }
615 
616 #undef __FUNC__
617 #define __FUNC__ "MatSNESFDMFResetHHistory"
618 /*@
619    MatSNESFDMFResetHHistory - Resets the counter to zero to begin
620       collecting a new set of differencing histories.
621 
622    Collective on Mat
623 
624    Input Parameters:
625 .  J - the matrix-free matrix context
626 
627    Notes:
628     Use MatSNESFDMFSetHHistory() to create the original history counter
629 
630 .keywords: SNES, matrix-free, h history, differencing history
631 
632 .seealso: MatSNESFDMFGetH(), MatCreateSNESFDMF(),
633           MatSNESFDMFSetHHistory(),
634           MatSNESFDMFKSPMonitor(), MatSNESFDMFSetFunctionError()
635 
636 @*/
637 int MatSNESFDMFResetHHistory(Mat J)
638 {
639   int            ierr;
640   MatSNESFDMFCtx ctx;
641 
642   PetscFunctionBegin;
643 
644   ierr = MatShellGetContext(J,(void **)&ctx); CHKERRQ(ierr);
645   /* no context indicates that it is not the "matrix free" matrix type */
646   if (!ctx) PetscFunctionReturn(0);
647   ctx->ncurrenth    = 0;
648 
649   PetscFunctionReturn(0);
650 }
651 
652