xref: /petsc/src/mat/impls/mffd/mffd.c (revision 2205254efee3a00a594e5e2a3a70f74dcb40bc03)
1 
2 #include <petsc-private/matimpl.h>
3 #include <../src/mat/impls/mffd/mffdimpl.h>   /*I  "petscmat.h"   I*/
4 
5 PetscFunctionList MatMFFDList              = 0;
6 PetscBool         MatMFFDRegisterAllCalled = PETSC_FALSE;
7 
8 PetscClassId  MATMFFD_CLASSID;
9 PetscLogEvent MATMFFD_Mult;
10 
11 static PetscBool MatMFFDPackageInitialized = PETSC_FALSE;
12 #undef __FUNCT__
13 #define __FUNCT__ "MatMFFDFinalizePackage"
14 /*@C
15   MatMFFDFinalizePackage - This function destroys everything in the MatMFFD package. It is
16   called from PetscFinalize().
17 
18   Level: developer
19 
20 .keywords: Petsc, destroy, package
21 .seealso: PetscFinalize()
22 @*/
23 PetscErrorCode  MatMFFDFinalizePackage(void)
24 {
25   PetscFunctionBegin;
26   MatMFFDPackageInitialized = PETSC_FALSE;
27   MatMFFDRegisterAllCalled  = PETSC_FALSE;
28   MatMFFDList               = PETSC_NULL;
29   PetscFunctionReturn(0);
30 }
31 
32 #undef __FUNCT__
33 #define __FUNCT__ "MatMFFDInitializePackage"
34 /*@C
35   MatMFFDInitializePackage - This function initializes everything in the MatMFFD package. It is called
36   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to MatCreate_MFFD()
37   when using static libraries.
38 
39   Input Parameter:
40 . path - The dynamic library path, or PETSC_NULL
41 
42   Level: developer
43 
44 .keywords: Vec, initialize, package
45 .seealso: PetscInitialize()
46 @*/
47 PetscErrorCode  MatMFFDInitializePackage(const char path[])
48 {
49   char           logList[256];
50   char           *className;
51   PetscBool      opt;
52   PetscErrorCode ierr;
53 
54   PetscFunctionBegin;
55   if (MatMFFDPackageInitialized) PetscFunctionReturn(0);
56   MatMFFDPackageInitialized = PETSC_TRUE;
57   /* Register Classes */
58   ierr = PetscClassIdRegister("MatMFFD",&MATMFFD_CLASSID);CHKERRQ(ierr);
59   /* Register Constructors */
60   ierr = MatMFFDRegisterAll(path);CHKERRQ(ierr);
61   /* Register Events */
62   ierr = PetscLogEventRegister("MatMult MF",          MATMFFD_CLASSID,&MATMFFD_Mult);CHKERRQ(ierr);
63 
64   /* Process info exclusions */
65   ierr = PetscOptionsGetString(PETSC_NULL, "-info_exclude", logList, 256, &opt);CHKERRQ(ierr);
66   if (opt) {
67     ierr = PetscStrstr(logList, "matmffd", &className);CHKERRQ(ierr);
68     if (className) {
69       ierr = PetscInfoDeactivateClass(MATMFFD_CLASSID);CHKERRQ(ierr);
70     }
71   }
72   /* Process summary exclusions */
73   ierr = PetscOptionsGetString(PETSC_NULL, "-log_summary_exclude", logList, 256, &opt);CHKERRQ(ierr);
74   if (opt) {
75     ierr = PetscStrstr(logList, "matmffd", &className);CHKERRQ(ierr);
76     if (className) {
77       ierr = PetscLogEventDeactivateClass(MATMFFD_CLASSID);CHKERRQ(ierr);
78     }
79   }
80   ierr = PetscRegisterFinalize(MatMFFDFinalizePackage);CHKERRQ(ierr);
81   PetscFunctionReturn(0);
82 }
83 
84 #undef __FUNCT__
85 #define __FUNCT__ "MatMFFDSetType"
86 /*@C
87     MatMFFDSetType - Sets the method that is used to compute the
88     differencing parameter for finite differene matrix-free formulations.
89 
90     Input Parameters:
91 +   mat - the "matrix-free" matrix created via MatCreateSNESMF(), or MatCreateMFFD()
92           or MatSetType(mat,MATMFFD);
93 -   ftype - the type requested, either MATMFFD_WP or MATMFFD_DS
94 
95     Level: advanced
96 
97     Notes:
98     For example, such routines can compute h for use in
99     Jacobian-vector products of the form
100 
101                         F(x+ha) - F(x)
102           F'(u)a  ~=  ----------------
103                               h
104 
105 .seealso: MatCreateSNESMF(), MatMFFDRegisterDynamic(), MatMFFDSetFunction()
106 @*/
107 PetscErrorCode  MatMFFDSetType(Mat mat,MatMFFDType ftype)
108 {
109   PetscErrorCode ierr,(*r)(MatMFFD);
110   MatMFFD        ctx = (MatMFFD)mat->data;
111   PetscBool      match;
112 
113   PetscFunctionBegin;
114   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
115   PetscValidCharPointer(ftype,2);
116 
117   ierr = PetscObjectTypeCompare((PetscObject)mat,MATMFFD,&match);CHKERRQ(ierr);
118   if (!match) PetscFunctionReturn(0);
119 
120   /* already set, so just return */
121   ierr = PetscObjectTypeCompare((PetscObject)ctx,ftype,&match);CHKERRQ(ierr);
122   if (match) PetscFunctionReturn(0);
123 
124   /* destroy the old one if it exists */
125   if (ctx->ops->destroy) {
126     ierr = (*ctx->ops->destroy)(ctx);CHKERRQ(ierr);
127   }
128 
129   ierr =  PetscFunctionListFind(((PetscObject)ctx)->comm,MatMFFDList,ftype,PETSC_TRUE,(void (**)(void)) &r);CHKERRQ(ierr);
130   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown MatMFFD type %s given",ftype);
131   ierr = (*r)(ctx);CHKERRQ(ierr);
132   ierr = PetscObjectChangeTypeName((PetscObject)ctx,ftype);CHKERRQ(ierr);
133   PetscFunctionReturn(0);
134 }
135 
136 typedef PetscErrorCode (*FCN1)(void*,Vec); /* force argument to next function to not be extern C*/
137 EXTERN_C_BEGIN
138 #undef __FUNCT__
139 #define __FUNCT__ "MatMFFDSetFunctioniBase_MFFD"
140 PetscErrorCode  MatMFFDSetFunctioniBase_MFFD(Mat mat,FCN1 func)
141 {
142   MatMFFD ctx = (MatMFFD)mat->data;
143 
144   PetscFunctionBegin;
145   ctx->funcisetbase = func;
146   PetscFunctionReturn(0);
147 }
148 EXTERN_C_END
149 
150 typedef PetscErrorCode (*FCN2)(void*,PetscInt,Vec,PetscScalar*); /* force argument to next function to not be extern C*/
151 EXTERN_C_BEGIN
152 #undef __FUNCT__
153 #define __FUNCT__ "MatMFFDSetFunctioni_MFFD"
154 PetscErrorCode  MatMFFDSetFunctioni_MFFD(Mat mat,FCN2 funci)
155 {
156   MatMFFD ctx = (MatMFFD)mat->data;
157 
158   PetscFunctionBegin;
159   ctx->funci = funci;
160   PetscFunctionReturn(0);
161 }
162 EXTERN_C_END
163 
164 EXTERN_C_BEGIN
165 #undef __FUNCT__
166 #define __FUNCT__ "MatMFFDResetHHistory_MFFD"
167 PetscErrorCode  MatMFFDResetHHistory_MFFD(Mat J)
168 {
169   MatMFFD ctx = (MatMFFD)J->data;
170 
171   PetscFunctionBegin;
172   ctx->ncurrenth = 0;
173   PetscFunctionReturn(0);
174 }
175 EXTERN_C_END
176 
177 #undef __FUNCT__
178 #define __FUNCT__ "MatMFFDRegister"
179 PetscErrorCode  MatMFFDRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(MatMFFD))
180 {
181   PetscErrorCode ierr;
182   char           fullname[PETSC_MAX_PATH_LEN];
183 
184   PetscFunctionBegin;
185   ierr = PetscFunctionListConcat(path,name,fullname);CHKERRQ(ierr);
186   ierr = PetscFunctionListAdd(PETSC_COMM_WORLD,&MatMFFDList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
187   PetscFunctionReturn(0);
188 }
189 
190 
191 #undef __FUNCT__
192 #define __FUNCT__ "MatMFFDRegisterDestroy"
193 /*@C
194    MatMFFDRegisterDestroy - Frees the list of MatMFFD methods that were
195    registered by MatMFFDRegisterDynamic).
196 
197    Not Collective
198 
199    Level: developer
200 
201 .keywords: MatMFFD, register, destroy
202 
203 .seealso: MatMFFDRegisterDynamic), MatMFFDRegisterAll()
204 @*/
205 PetscErrorCode  MatMFFDRegisterDestroy(void)
206 {
207   PetscErrorCode ierr;
208 
209   PetscFunctionBegin;
210   ierr = PetscFunctionListDestroy(&MatMFFDList);CHKERRQ(ierr);
211 
212   MatMFFDRegisterAllCalled = PETSC_FALSE;
213   PetscFunctionReturn(0);
214 }
215 
216 EXTERN_C_BEGIN
217 #undef __FUNCT__
218 #define __FUNCT__ "MatMFFDAddNullSpace_MFFD"
219 PetscErrorCode  MatMFFDAddNullSpace_MFFD(Mat J,MatNullSpace nullsp)
220 {
221   PetscErrorCode ierr;
222   MatMFFD        ctx = (MatMFFD)J->data;
223 
224   PetscFunctionBegin;
225   ierr = PetscObjectReference((PetscObject)nullsp);CHKERRQ(ierr);
226   if (ctx->sp) { ierr = MatNullSpaceDestroy(&ctx->sp);CHKERRQ(ierr); }
227   ctx->sp = nullsp;
228   PetscFunctionReturn(0);
229 }
230 EXTERN_C_END
231 
232 /* ----------------------------------------------------------------------------------------*/
233 #undef __FUNCT__
234 #define __FUNCT__ "MatDestroy_MFFD"
235 PetscErrorCode MatDestroy_MFFD(Mat mat)
236 {
237   PetscErrorCode ierr;
238   MatMFFD        ctx = (MatMFFD)mat->data;
239 
240   PetscFunctionBegin;
241   ierr = VecDestroy(&ctx->w);CHKERRQ(ierr);
242   ierr = VecDestroy(&ctx->drscale);CHKERRQ(ierr);
243   ierr = VecDestroy(&ctx->dlscale);CHKERRQ(ierr);
244   ierr = VecDestroy(&ctx->dshift);CHKERRQ(ierr);
245   if (ctx->current_f_allocated) {
246     ierr = VecDestroy(&ctx->current_f);CHKERRQ(ierr);
247   }
248   if (ctx->ops->destroy) {ierr = (*ctx->ops->destroy)(ctx);CHKERRQ(ierr);}
249   ierr      = MatNullSpaceDestroy(&ctx->sp);CHKERRQ(ierr);
250   ierr      = PetscHeaderDestroy(&ctx);CHKERRQ(ierr);
251   mat->data = 0;
252 
253   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetBase_C","",PETSC_NULL);CHKERRQ(ierr);
254   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctioniBase_C","",PETSC_NULL);CHKERRQ(ierr);
255   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctioni_C","",PETSC_NULL);CHKERRQ(ierr);
256   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunction_C","",PETSC_NULL);CHKERRQ(ierr);
257   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctionError_C","",PETSC_NULL);CHKERRQ(ierr);
258   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetCheckh_C","",PETSC_NULL);CHKERRQ(ierr);
259   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetPeriod_C","",PETSC_NULL);CHKERRQ(ierr);
260   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDResetHHistory_C","",PETSC_NULL);CHKERRQ(ierr);
261   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDAddNullSpace_C","",PETSC_NULL);CHKERRQ(ierr);
262   PetscFunctionReturn(0);
263 }
264 
265 #undef __FUNCT__
266 #define __FUNCT__ "MatView_MFFD"
267 /*
268    MatMFFDView_MFFD - Views matrix-free parameters.
269 
270 */
271 PetscErrorCode MatView_MFFD(Mat J,PetscViewer viewer)
272 {
273   PetscErrorCode ierr;
274   MatMFFD        ctx = (MatMFFD)J->data;
275   PetscBool      iascii, viewbase, viewfunction;
276   const char     *prefix;
277 
278   PetscFunctionBegin;
279   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
280   if (iascii) {
281     ierr = PetscViewerASCIIPrintf(viewer,"Matrix-free approximation:\n");CHKERRQ(ierr);
282     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
283     ierr = PetscViewerASCIIPrintf(viewer,"err=%G (relative error in function evaluation)\n",ctx->error_rel);CHKERRQ(ierr);
284     if (!((PetscObject)ctx)->type_name) {
285       ierr = PetscViewerASCIIPrintf(viewer,"The compute h routine has not yet been set\n");CHKERRQ(ierr);
286     } else {
287       ierr = PetscViewerASCIIPrintf(viewer,"Using %s compute h routine\n",((PetscObject)ctx)->type_name);CHKERRQ(ierr);
288     }
289     if (ctx->ops->view) {
290       ierr = (*ctx->ops->view)(ctx,viewer);CHKERRQ(ierr);
291     }
292     ierr = PetscObjectGetOptionsPrefix((PetscObject)J, &prefix);CHKERRQ(ierr);
293 
294     ierr = PetscOptionsHasName(prefix, "-mat_mffd_view_base", &viewbase);CHKERRQ(ierr);
295     if (viewbase) {
296       ierr = PetscViewerASCIIPrintf(viewer, "Base:\n");CHKERRQ(ierr);
297       ierr = VecView(ctx->current_u, viewer);CHKERRQ(ierr);
298     }
299     ierr = PetscOptionsHasName(prefix, "-mat_mffd_view_function", &viewfunction);CHKERRQ(ierr);
300     if (viewfunction) {
301       ierr = PetscViewerASCIIPrintf(viewer, "Function:\n");CHKERRQ(ierr);
302       ierr = VecView(ctx->current_f, viewer);CHKERRQ(ierr);
303     }
304     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
305   }
306   PetscFunctionReturn(0);
307 }
308 
309 #undef __FUNCT__
310 #define __FUNCT__ "MatAssemblyEnd_MFFD"
311 /*
312    MatAssemblyEnd_MFFD - Resets the ctx->ncurrenth to zero. This
313    allows the user to indicate the beginning of a new linear solve by calling
314    MatAssemblyXXX() on the matrix free matrix. This then allows the
315    MatCreateMFFD_WP() to properly compute ||U|| only the first time
316    in the linear solver rather than every time.
317 */
318 PetscErrorCode MatAssemblyEnd_MFFD(Mat J,MatAssemblyType mt)
319 {
320   PetscErrorCode ierr;
321   MatMFFD        j = (MatMFFD)J->data;
322 
323   PetscFunctionBegin;
324   ierr      = MatMFFDResetHHistory(J);CHKERRQ(ierr);
325   j->vshift = 0.0;
326   j->vscale = 1.0;
327   PetscFunctionReturn(0);
328 }
329 
330 #undef __FUNCT__
331 #define __FUNCT__ "MatMult_MFFD"
332 /*
333   MatMult_MFFD - Default matrix-free form for Jacobian-vector product, y = F'(u)*a:
334 
335         y ~= (F(u + ha) - F(u))/h,
336   where F = nonlinear function, as set by SNESSetFunction()
337         u = current iterate
338         h = difference interval
339 */
340 PetscErrorCode MatMult_MFFD(Mat mat,Vec a,Vec y)
341 {
342   MatMFFD        ctx = (MatMFFD)mat->data;
343   PetscScalar    h;
344   Vec            w,U,F;
345   PetscErrorCode ierr;
346   PetscBool      zeroa;
347 
348   PetscFunctionBegin;
349   if (!ctx->current_u) SETERRQ(((PetscObject)mat)->comm,PETSC_ERR_ARG_WRONGSTATE,"MatMFFDSetBase() has not been called, this is often caused by forgetting to call \n\t\tMatAssemblyBegin/End on the first Mat in the SNES compute function");
350   /* We log matrix-free matrix-vector products separately, so that we can
351      separate the performance monitoring from the cases that use conventional
352      storage.  We may eventually modify event logging to associate events
353      with particular objects, hence alleviating the more general problem. */
354   ierr = PetscLogEventBegin(MATMFFD_Mult,a,y,0,0);CHKERRQ(ierr);
355 
356   w = ctx->w;
357   U = ctx->current_u;
358   F = ctx->current_f;
359   /*
360       Compute differencing parameter
361   */
362   if (!ctx->ops->compute) {
363     ierr = MatMFFDSetType(mat,MATMFFD_WP);CHKERRQ(ierr);
364     ierr = MatSetFromOptions(mat);CHKERRQ(ierr);
365   }
366   ierr = (*ctx->ops->compute)(ctx,U,a,&h,&zeroa);CHKERRQ(ierr);
367   if (zeroa) {
368     ierr = VecSet(y,0.0);CHKERRQ(ierr);
369     PetscFunctionReturn(0);
370   }
371 
372   if (PetscIsInfOrNanScalar(h)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Computed Nan differencing parameter h");
373   if (ctx->checkh) {
374     ierr = (*ctx->checkh)(ctx->checkhctx,U,a,&h);CHKERRQ(ierr);
375   }
376 
377   /* keep a record of the current differencing parameter h */
378   ctx->currenth = h;
379 #if defined(PETSC_USE_COMPLEX)
380   ierr = PetscInfo2(mat,"Current differencing parameter: %G + %G i\n",PetscRealPart(h),PetscImaginaryPart(h));CHKERRQ(ierr);
381 #else
382   ierr = PetscInfo1(mat,"Current differencing parameter: %15.12e\n",h);CHKERRQ(ierr);
383 #endif
384   if (ctx->historyh && ctx->ncurrenth < ctx->maxcurrenth) {
385     ctx->historyh[ctx->ncurrenth] = h;
386   }
387   ctx->ncurrenth++;
388 
389   /* w = u + ha */
390   if (ctx->drscale) {
391     ierr = VecPointwiseMult(ctx->drscale,a,U);CHKERRQ(ierr);
392     ierr = VecAYPX(U,h,w);CHKERRQ(ierr);
393   } else {
394     ierr = VecWAXPY(w,h,a,U);CHKERRQ(ierr);
395   }
396 
397   /* compute func(U) as base for differencing; only needed first time in and not when provided by user */
398   if (ctx->ncurrenth == 1 && ctx->current_f_allocated) {
399     ierr = (*ctx->func)(ctx->funcctx,U,F);CHKERRQ(ierr);
400   }
401   ierr = (*ctx->func)(ctx->funcctx,w,y);CHKERRQ(ierr);
402 
403   ierr = VecAXPY(y,-1.0,F);CHKERRQ(ierr);
404   ierr = VecScale(y,1.0/h);CHKERRQ(ierr);
405 
406   if ((ctx->vshift != 0.0) || (ctx->vscale != 1.0)) {
407     ierr = VecAXPBY(y,ctx->vshift,ctx->vscale,a);CHKERRQ(ierr);
408   }
409   if (ctx->dlscale) {
410     ierr = VecPointwiseMult(y,ctx->dlscale,y);CHKERRQ(ierr);
411   }
412   if (ctx->dshift) {
413     ierr = VecPointwiseMult(ctx->dshift,a,U);CHKERRQ(ierr);
414     ierr = VecAXPY(y,1.0,U);CHKERRQ(ierr);
415   }
416 
417   if (ctx->sp) {ierr = MatNullSpaceRemove(ctx->sp,y,PETSC_NULL);CHKERRQ(ierr);}
418 
419   ierr = PetscLogEventEnd(MATMFFD_Mult,a,y,0,0);CHKERRQ(ierr);
420   PetscFunctionReturn(0);
421 }
422 
423 #undef __FUNCT__
424 #define __FUNCT__ "MatGetDiagonal_MFFD"
425 /*
426   MatGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix
427 
428         y ~= (F(u + ha) - F(u))/h,
429   where F = nonlinear function, as set by SNESSetFunction()
430         u = current iterate
431         h = difference interval
432 */
433 PetscErrorCode MatGetDiagonal_MFFD(Mat mat,Vec a)
434 {
435   MatMFFD        ctx = (MatMFFD)mat->data;
436   PetscScalar    h,*aa,*ww,v;
437   PetscReal      epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON;
438   Vec            w,U;
439   PetscErrorCode ierr;
440   PetscInt       i,rstart,rend;
441 
442   PetscFunctionBegin;
443   if (!ctx->funci) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Requires calling MatMFFDSetFunctioni() first");
444 
445   w    = ctx->w;
446   U    = ctx->current_u;
447   ierr = (*ctx->func)(ctx->funcctx,U,a);CHKERRQ(ierr);
448   ierr = (*ctx->funcisetbase)(ctx->funcctx,U);CHKERRQ(ierr);
449   ierr = VecCopy(U,w);CHKERRQ(ierr);
450 
451   ierr = VecGetOwnershipRange(a,&rstart,&rend);CHKERRQ(ierr);
452   ierr = VecGetArray(a,&aa);CHKERRQ(ierr);
453   for (i=rstart; i<rend; i++) {
454     ierr = VecGetArray(w,&ww);CHKERRQ(ierr);
455     h    = ww[i-rstart];
456     if (h == 0.0) h = 1.0;
457     if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0)     h = umin;
458     else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin;
459     h *= epsilon;
460 
461     ww[i-rstart] += h;
462     ierr = VecRestoreArray(w,&ww);CHKERRQ(ierr);
463     ierr          = (*ctx->funci)(ctx->funcctx,i,w,&v);CHKERRQ(ierr);
464     aa[i-rstart]  = (v - aa[i-rstart])/h;
465 
466     /* possibly shift and scale result */
467     if ((ctx->vshift != 0.0) || (ctx->vscale != 1.0)) {
468       aa[i - rstart] = ctx->vshift + ctx->vscale*aa[i-rstart];
469     }
470 
471     ierr          = VecGetArray(w,&ww);CHKERRQ(ierr);
472     ww[i-rstart] -= h;
473     ierr          = VecRestoreArray(w,&ww);CHKERRQ(ierr);
474   }
475   ierr = VecRestoreArray(a,&aa);CHKERRQ(ierr);
476   PetscFunctionReturn(0);
477 }
478 
479 #undef __FUNCT__
480 #define __FUNCT__ "MatDiagonalScale_MFFD"
481 PetscErrorCode MatDiagonalScale_MFFD(Mat mat,Vec ll,Vec rr)
482 {
483   MatMFFD        aij = (MatMFFD)mat->data;
484   PetscErrorCode ierr;
485 
486   PetscFunctionBegin;
487   if (ll && !aij->dlscale) {
488     ierr = VecDuplicate(ll,&aij->dlscale);CHKERRQ(ierr);
489   }
490   if (rr && !aij->drscale) {
491     ierr = VecDuplicate(rr,&aij->drscale);CHKERRQ(ierr);
492   }
493   if (ll) {
494     ierr = VecCopy(ll,aij->dlscale);CHKERRQ(ierr);
495   }
496   if (rr) {
497     ierr = VecCopy(rr,aij->drscale);CHKERRQ(ierr);
498   }
499   PetscFunctionReturn(0);
500 }
501 
502 #undef __FUNCT__
503 #define __FUNCT__ "MatDiagonalSet_MFFD"
504 PetscErrorCode MatDiagonalSet_MFFD(Mat mat,Vec ll,InsertMode mode)
505 {
506   MatMFFD        aij = (MatMFFD)mat->data;
507   PetscErrorCode ierr;
508 
509   PetscFunctionBegin;
510   if (mode == INSERT_VALUES) SETERRQ(((PetscObject)mat)->comm,PETSC_ERR_SUP,"No diagonal set with INSERT_VALUES");
511   if (!aij->dshift) {
512     ierr = VecDuplicate(ll,&aij->dshift);CHKERRQ(ierr);
513   }
514   ierr = VecAXPY(aij->dshift,1.0,ll);CHKERRQ(ierr);
515   PetscFunctionReturn(0);
516 }
517 
518 #undef __FUNCT__
519 #define __FUNCT__ "MatShift_MFFD"
520 PetscErrorCode MatShift_MFFD(Mat Y,PetscScalar a)
521 {
522   MatMFFD shell = (MatMFFD)Y->data;
523 
524   PetscFunctionBegin;
525   shell->vshift += a;
526   PetscFunctionReturn(0);
527 }
528 
529 #undef __FUNCT__
530 #define __FUNCT__ "MatScale_MFFD"
531 PetscErrorCode MatScale_MFFD(Mat Y,PetscScalar a)
532 {
533   MatMFFD shell = (MatMFFD)Y->data;
534 
535   PetscFunctionBegin;
536   shell->vscale *= a;
537   PetscFunctionReturn(0);
538 }
539 
540 EXTERN_C_BEGIN
541 #undef __FUNCT__
542 #define __FUNCT__ "MatMFFDSetBase_MFFD"
543 PetscErrorCode  MatMFFDSetBase_MFFD(Mat J,Vec U,Vec F)
544 {
545   PetscErrorCode ierr;
546   MatMFFD        ctx = (MatMFFD)J->data;
547 
548   PetscFunctionBegin;
549   ierr = MatMFFDResetHHistory(J);CHKERRQ(ierr);
550 
551   ctx->current_u = U;
552   if (F) {
553     if (ctx->current_f_allocated) {ierr = VecDestroy(&ctx->current_f);CHKERRQ(ierr);}
554     ctx->current_f           = F;
555     ctx->current_f_allocated = PETSC_FALSE;
556   } else if (!ctx->current_f_allocated) {
557     ierr = VecDuplicate(ctx->current_u, &ctx->current_f);CHKERRQ(ierr);
558 
559     ctx->current_f_allocated = PETSC_TRUE;
560   }
561   if (!ctx->w) {
562     ierr = VecDuplicate(ctx->current_u, &ctx->w);CHKERRQ(ierr);
563   }
564   J->assembled = PETSC_TRUE;
565   PetscFunctionReturn(0);
566 }
567 EXTERN_C_END
568 typedef PetscErrorCode (*FCN3)(void*,Vec,Vec,PetscScalar*); /* force argument to next function to not be extern C*/
569 EXTERN_C_BEGIN
570 #undef __FUNCT__
571 #define __FUNCT__ "MatMFFDSetCheckh_MFFD"
572 PetscErrorCode  MatMFFDSetCheckh_MFFD(Mat J,FCN3 fun,void *ectx)
573 {
574   MatMFFD ctx = (MatMFFD)J->data;
575 
576   PetscFunctionBegin;
577   ctx->checkh    = fun;
578   ctx->checkhctx = ectx;
579   PetscFunctionReturn(0);
580 }
581 EXTERN_C_END
582 
583 #undef __FUNCT__
584 #define __FUNCT__ "MatMFFDSetOptionsPrefix"
585 /*@C
586    MatMFFDSetOptionsPrefix - Sets the prefix used for searching for all
587    MatMFFD options in the database.
588 
589    Collective on Mat
590 
591    Input Parameter:
592 +  A - the Mat context
593 -  prefix - the prefix to prepend to all option names
594 
595    Notes:
596    A hyphen (-) must NOT be given at the beginning of the prefix name.
597    The first character of all runtime options is AUTOMATICALLY the hyphen.
598 
599    Level: advanced
600 
601 .keywords: SNES, matrix-free, parameters
602 
603 .seealso: MatSetFromOptions(), MatCreateSNESMF()
604 @*/
605 PetscErrorCode  MatMFFDSetOptionsPrefix(Mat mat,const char prefix[])
606 
607 {
608   MatMFFD        mfctx = mat ? (MatMFFD)mat->data : (MatMFFD)PETSC_NULL;
609   PetscErrorCode ierr;
610 
611   PetscFunctionBegin;
612   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
613   PetscValidHeaderSpecific(mfctx,MATMFFD_CLASSID,1);
614   ierr = PetscObjectSetOptionsPrefix((PetscObject)mfctx,prefix);CHKERRQ(ierr);
615   PetscFunctionReturn(0);
616 }
617 
618 #undef __FUNCT__
619 #define __FUNCT__ "MatSetFromOptions_MFFD"
620 PetscErrorCode  MatSetFromOptions_MFFD(Mat mat)
621 {
622   MatMFFD        mfctx = (MatMFFD)mat->data;
623   PetscErrorCode ierr;
624   PetscBool      flg;
625   char           ftype[256];
626 
627   PetscFunctionBegin;
628   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
629   PetscValidHeaderSpecific(mfctx,MATMFFD_CLASSID,1);
630   ierr = PetscObjectOptionsBegin((PetscObject)mfctx);CHKERRQ(ierr);
631   ierr = PetscOptionsList("-mat_mffd_type","Matrix free type","MatMFFDSetType",MatMFFDList,((PetscObject)mfctx)->type_name,ftype,256,&flg);CHKERRQ(ierr);
632   if (flg) {
633     ierr = MatMFFDSetType(mat,ftype);CHKERRQ(ierr);
634   }
635 
636   ierr = PetscOptionsReal("-mat_mffd_err","set sqrt relative error in function","MatMFFDSetFunctionError",mfctx->error_rel,&mfctx->error_rel,0);CHKERRQ(ierr);
637   ierr = PetscOptionsInt("-mat_mffd_period","how often h is recomputed","MatMFFDSetPeriod",mfctx->recomputeperiod,&mfctx->recomputeperiod,0);CHKERRQ(ierr);
638 
639   flg  = PETSC_FALSE;
640   ierr = PetscOptionsBool("-mat_mffd_check_positivity","Insure that U + h*a is nonnegative","MatMFFDSetCheckh",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
641   if (flg) {
642     ierr = MatMFFDSetCheckh(mat,MatMFFDCheckPositivity,0);CHKERRQ(ierr);
643   }
644   if (mfctx->ops->setfromoptions) {
645     ierr = (*mfctx->ops->setfromoptions)(mfctx);CHKERRQ(ierr);
646   }
647   ierr = PetscOptionsEnd();CHKERRQ(ierr);
648   PetscFunctionReturn(0);
649 }
650 
651 EXTERN_C_BEGIN
652 #undef __FUNCT__
653 #define __FUNCT__ "MatMFFDSetPeriod_MFFD"
654 PetscErrorCode  MatMFFDSetPeriod_MFFD(Mat mat,PetscInt period)
655 {
656   MatMFFD ctx = (MatMFFD)mat->data;
657 
658   PetscFunctionBegin;
659   PetscValidLogicalCollectiveInt(mat,period,2);
660   ctx->recomputeperiod = period;
661   PetscFunctionReturn(0);
662 }
663 EXTERN_C_END
664 
665 EXTERN_C_BEGIN
666 #undef __FUNCT__
667 #define __FUNCT__ "MatMFFDSetFunction_MFFD"
668 PetscErrorCode  MatMFFDSetFunction_MFFD(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx)
669 {
670   MatMFFD ctx = (MatMFFD)mat->data;
671 
672   PetscFunctionBegin;
673   ctx->func    = func;
674   ctx->funcctx = funcctx;
675   PetscFunctionReturn(0);
676 }
677 EXTERN_C_END
678 
679 EXTERN_C_BEGIN
680 #undef __FUNCT__
681 #define __FUNCT__ "MatMFFDSetFunctionError_MFFD"
682 PetscErrorCode  MatMFFDSetFunctionError_MFFD(Mat mat,PetscReal error)
683 {
684   MatMFFD ctx = (MatMFFD)mat->data;
685 
686   PetscFunctionBegin;
687   PetscValidLogicalCollectiveReal(mat,error,2);
688   if (error != PETSC_DEFAULT) ctx->error_rel = error;
689   PetscFunctionReturn(0);
690 }
691 EXTERN_C_END
692 
693 /*MC
694   MATMFFD - MATMFFD = "mffd" - A matrix free matrix type.
695 
696   Level: advanced
697 
698 .seealso: MatCreateMFFD(), MatCreateSNESMF(), MatMFFDSetFunction()
699 M*/
700 EXTERN_C_BEGIN
701 #undef __FUNCT__
702 #define __FUNCT__ "MatCreate_MFFD"
703 PetscErrorCode  MatCreate_MFFD(Mat A)
704 {
705   MatMFFD        mfctx;
706   PetscErrorCode ierr;
707 
708   PetscFunctionBegin;
709 #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
710   ierr = MatMFFDInitializePackage(PETSC_NULL);CHKERRQ(ierr);
711 #endif
712 
713   ierr = PetscHeaderCreate(mfctx,_p_MatMFFD,struct _MFOps,MATMFFD_CLASSID,0,"MatMFFD","Matrix-free Finite Differencing","Mat",((PetscObject)A)->comm,MatDestroy_MFFD,MatView_MFFD);CHKERRQ(ierr);
714 
715   mfctx->sp                       = 0;
716   mfctx->error_rel                = PETSC_SQRT_MACHINE_EPSILON;
717   mfctx->recomputeperiod          = 1;
718   mfctx->count                    = 0;
719   mfctx->currenth                 = 0.0;
720   mfctx->historyh                 = PETSC_NULL;
721   mfctx->ncurrenth                = 0;
722   mfctx->maxcurrenth              = 0;
723   ((PetscObject)mfctx)->type_name = 0;
724 
725   mfctx->vshift = 0.0;
726   mfctx->vscale = 1.0;
727 
728   /*
729      Create the empty data structure to contain compute-h routines.
730      These will be filled in below from the command line options or
731      a later call with MatMFFDSetType() or if that is not called
732      then it will default in the first use of MatMult_MFFD()
733   */
734   mfctx->ops->compute        = 0;
735   mfctx->ops->destroy        = 0;
736   mfctx->ops->view           = 0;
737   mfctx->ops->setfromoptions = 0;
738   mfctx->hctx                = 0;
739 
740   mfctx->func    = 0;
741   mfctx->funcctx = 0;
742   mfctx->w       = PETSC_NULL;
743 
744   A->data = mfctx;
745 
746   A->ops->mult           = MatMult_MFFD;
747   A->ops->destroy        = MatDestroy_MFFD;
748   A->ops->view           = MatView_MFFD;
749   A->ops->assemblyend    = MatAssemblyEnd_MFFD;
750   A->ops->getdiagonal    = MatGetDiagonal_MFFD;
751   A->ops->scale          = MatScale_MFFD;
752   A->ops->shift          = MatShift_MFFD;
753   A->ops->diagonalscale  = MatDiagonalScale_MFFD;
754   A->ops->diagonalset    = MatDiagonalSet_MFFD;
755   A->ops->setfromoptions = MatSetFromOptions_MFFD;
756   A->assembled           = PETSC_TRUE;
757 
758   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
759   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
760 
761   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDSetBase_C","MatMFFDSetBase_MFFD",MatMFFDSetBase_MFFD);CHKERRQ(ierr);
762   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDSetFunctioniBase_C","MatMFFDSetFunctioniBase_MFFD",MatMFFDSetFunctioniBase_MFFD);CHKERRQ(ierr);
763   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDSetFunctioni_C","MatMFFDSetFunctioni_MFFD",MatMFFDSetFunctioni_MFFD);CHKERRQ(ierr);
764   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDSetFunction_C","MatMFFDSetFunction_MFFD",MatMFFDSetFunction_MFFD);CHKERRQ(ierr);
765   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDSetCheckh_C","MatMFFDSetCheckh_MFFD",MatMFFDSetCheckh_MFFD);CHKERRQ(ierr);
766   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDSetPeriod_C","MatMFFDSetPeriod_MFFD",MatMFFDSetPeriod_MFFD);CHKERRQ(ierr);
767   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDSetFunctionError_C","MatMFFDSetFunctionError_MFFD",MatMFFDSetFunctionError_MFFD);CHKERRQ(ierr);
768   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDResetHHistory_C","MatMFFDResetHHistory_MFFD",MatMFFDResetHHistory_MFFD);CHKERRQ(ierr);
769   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDAddNullSpace_C","MatMFFDAddNullSpace_MFFD",MatMFFDAddNullSpace_MFFD);CHKERRQ(ierr);
770 
771   mfctx->mat = A;
772 
773   ierr = PetscObjectChangeTypeName((PetscObject)A,MATMFFD);CHKERRQ(ierr);
774   PetscFunctionReturn(0);
775 }
776 EXTERN_C_END
777 
778 #undef __FUNCT__
779 #define __FUNCT__ "MatCreateMFFD"
780 /*@
781    MatCreateMFFD - Creates a matrix-free matrix. See also MatCreateSNESMF()
782 
783    Collective on Vec
784 
785    Input Parameters:
786 +  comm - MPI communicator
787 .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
788            This value should be the same as the local size used in creating the
789            y vector for the matrix-vector product y = Ax.
790 .  n - This value should be the same as the local size used in creating the
791        x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have
792        calculated if N is given) For square matrices n is almost always m.
793 .  M - number of global rows (or PETSC_DETERMINE to have calculated if m is given)
794 -  N - number of global columns (or PETSC_DETERMINE to have calculated if n is given)
795 
796 
797    Output Parameter:
798 .  J - the matrix-free matrix
799 
800    Options Database Keys: call MatSetFromOptions() to trigger these
801 +  -mat_mffd_type - wp or ds (see MATMFFD_WP or MATMFFD_DS)
802 -  -mat_mffd_err - square root of estimated relative error in function evaluation
803 -  -mat_mffd_period - how often h is recomputed, defaults to 1, everytime
804 
805 
806    Level: advanced
807 
808    Notes:
809    The matrix-free matrix context merely contains the function pointers
810    and work space for performing finite difference approximations of
811    Jacobian-vector products, F'(u)*a,
812 
813    The default code uses the following approach to compute h
814 
815 .vb
816      F'(u)*a = [F(u+h*a) - F(u)]/h where
817      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
818        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
819  where
820      error_rel = square root of relative error in function evaluation
821      umin = minimum iterate parameter
822 .ve
823 
824    The user can set the error_rel via MatMFFDSetFunctionError() and
825    umin via MatMFFDDSSetUmin(); see the <A href="../../docs/manual.pdf#nameddest=ch_snes">SNES chapter of the users manual</A> for details.
826 
827    The user should call MatDestroy() when finished with the matrix-free
828    matrix context.
829 
830    Options Database Keys:
831 +  -mat_mffd_err <error_rel> - Sets error_rel
832 .  -mat_mffd_unim <umin> - Sets umin (for default PETSc routine that computes h only)
833 -  -mat_mffd_check_positivity
834 
835 .keywords: default, matrix-free, create, matrix
836 
837 .seealso: MatDestroy(), MatMFFDSetFunctionError(), MatMFFDDSSetUmin(), MatMFFDSetFunction()
838           MatMFFDSetHHistory(), MatMFFDResetHHistory(), MatCreateSNESMF(),
839           MatMFFDGetH(), MatMFFDRegisterDynamic), MatMFFDComputeJacobian()
840 
841 @*/
842 PetscErrorCode  MatCreateMFFD(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *J)
843 {
844   PetscErrorCode ierr;
845 
846   PetscFunctionBegin;
847   ierr = MatCreate(comm,J);CHKERRQ(ierr);
848   ierr = MatSetSizes(*J,m,n,M,N);CHKERRQ(ierr);
849   ierr = MatSetType(*J,MATMFFD);CHKERRQ(ierr);
850   ierr = MatSetUp(*J);CHKERRQ(ierr);
851   PetscFunctionReturn(0);
852 }
853 
854 
855 #undef __FUNCT__
856 #define __FUNCT__ "MatMFFDGetH"
857 /*@
858    MatMFFDGetH - Gets the last value that was used as the differencing
859    parameter.
860 
861    Not Collective
862 
863    Input Parameters:
864 .  mat - the matrix obtained with MatCreateSNESMF()
865 
866    Output Paramter:
867 .  h - the differencing step size
868 
869    Level: advanced
870 
871 .keywords: SNES, matrix-free, parameters
872 
873 .seealso: MatCreateSNESMF(),MatMFFDSetHHistory(), MatCreateMFFD(), MATMFFD, MatMFFDResetHHistory()
874 @*/
875 PetscErrorCode  MatMFFDGetH(Mat mat,PetscScalar *h)
876 {
877   MatMFFD        ctx = (MatMFFD)mat->data;
878   PetscErrorCode ierr;
879   PetscBool      match;
880 
881   PetscFunctionBegin;
882   ierr = PetscObjectTypeCompare((PetscObject)mat,MATMFFD,&match);CHKERRQ(ierr);
883   if (!match) SETERRQ(((PetscObject)mat)->comm,PETSC_ERR_ARG_WRONG,"Not a MFFD matrix");
884 
885   *h = ctx->currenth;
886   PetscFunctionReturn(0);
887 }
888 
889 #undef __FUNCT__
890 #define __FUNCT__ "MatMFFDSetFunction"
891 /*@C
892    MatMFFDSetFunction - Sets the function used in applying the matrix free.
893 
894    Logically Collective on Mat
895 
896    Input Parameters:
897 +  mat - the matrix free matrix created via MatCreateSNESMF()
898 .  func - the function to use
899 -  funcctx - optional function context passed to function
900 
901    Level: advanced
902 
903    Notes:
904     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
905     matrix inside your compute Jacobian routine
906 
907     If this is not set then it will use the function set with SNESSetFunction() if MatCreateSNESMF() was used.
908 
909 .keywords: SNES, matrix-free, function
910 
911 .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD,
912           MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction()
913 @*/
914 PetscErrorCode  MatMFFDSetFunction(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx)
915 {
916   PetscErrorCode ierr;
917 
918   PetscFunctionBegin;
919   ierr = PetscTryMethod(mat,"MatMFFDSetFunction_C",(Mat,PetscErrorCode (*)(void*,Vec,Vec),void*),(mat,func,funcctx));CHKERRQ(ierr);
920   PetscFunctionReturn(0);
921 }
922 
923 #undef __FUNCT__
924 #define __FUNCT__ "MatMFFDSetFunctioni"
925 /*@C
926    MatMFFDSetFunctioni - Sets the function for a single component
927 
928    Logically Collective on Mat
929 
930    Input Parameters:
931 +  mat - the matrix free matrix created via MatCreateSNESMF()
932 -  funci - the function to use
933 
934    Level: advanced
935 
936    Notes:
937     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
938     matrix inside your compute Jacobian routine
939 
940 
941 .keywords: SNES, matrix-free, function
942 
943 .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction()
944 
945 @*/
946 PetscErrorCode  MatMFFDSetFunctioni(Mat mat,PetscErrorCode (*funci)(void*,PetscInt,Vec,PetscScalar*))
947 {
948   PetscErrorCode ierr;
949 
950   PetscFunctionBegin;
951   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
952   ierr = PetscTryMethod(mat,"MatMFFDSetFunctioni_C",(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)),(mat,funci));CHKERRQ(ierr);
953   PetscFunctionReturn(0);
954 }
955 
956 
957 #undef __FUNCT__
958 #define __FUNCT__ "MatMFFDSetFunctioniBase"
959 /*@C
960    MatMFFDSetFunctioniBase - Sets the base vector for a single component function evaluation
961 
962    Logically Collective on Mat
963 
964    Input Parameters:
965 +  mat - the matrix free matrix created via MatCreateSNESMF()
966 -  func - the function to use
967 
968    Level: advanced
969 
970    Notes:
971     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
972     matrix inside your compute Jacobian routine
973 
974 
975 .keywords: SNES, matrix-free, function
976 
977 .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD
978           MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction()
979 @*/
980 PetscErrorCode  MatMFFDSetFunctioniBase(Mat mat,PetscErrorCode (*func)(void*,Vec))
981 {
982   PetscErrorCode ierr;
983 
984   PetscFunctionBegin;
985   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
986   ierr = PetscTryMethod(mat,"MatMFFDSetFunctioniBase_C",(Mat,PetscErrorCode (*)(void*,Vec)),(mat,func));CHKERRQ(ierr);
987   PetscFunctionReturn(0);
988 }
989 
990 #undef __FUNCT__
991 #define __FUNCT__ "MatMFFDSetPeriod"
992 /*@
993    MatMFFDSetPeriod - Sets how often h is recomputed, by default it is everytime
994 
995    Logically Collective on Mat
996 
997    Input Parameters:
998 +  mat - the matrix free matrix created via MatCreateSNESMF()
999 -  period - 1 for everytime, 2 for every second etc
1000 
1001    Options Database Keys:
1002 +  -mat_mffd_period <period>
1003 
1004    Level: advanced
1005 
1006 
1007 .keywords: SNES, matrix-free, parameters
1008 
1009 .seealso: MatCreateSNESMF(),MatMFFDGetH(),
1010           MatMFFDSetHHistory(), MatMFFDResetHHistory()
1011 @*/
1012 PetscErrorCode  MatMFFDSetPeriod(Mat mat,PetscInt period)
1013 {
1014   PetscErrorCode ierr;
1015 
1016   PetscFunctionBegin;
1017   ierr = PetscTryMethod(mat,"MatMFFDSetPeriod_C",(Mat,PetscInt),(mat,period));CHKERRQ(ierr);
1018   PetscFunctionReturn(0);
1019 }
1020 
1021 #undef __FUNCT__
1022 #define __FUNCT__ "MatMFFDSetFunctionError"
1023 /*@
1024    MatMFFDSetFunctionError - Sets the error_rel for the approximation of
1025    matrix-vector products using finite differences.
1026 
1027    Logically Collective on Mat
1028 
1029    Input Parameters:
1030 +  mat - the matrix free matrix created via MatCreateMFFD() or MatCreateSNESMF()
1031 -  error_rel - relative error (should be set to the square root of
1032                the relative error in the function evaluations)
1033 
1034    Options Database Keys:
1035 +  -mat_mffd_err <error_rel> - Sets error_rel
1036 
1037    Level: advanced
1038 
1039    Notes:
1040    The default matrix-free matrix-vector product routine computes
1041 .vb
1042      F'(u)*a = [F(u+h*a) - F(u)]/h where
1043      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
1044        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   else
1045 .ve
1046 
1047 .keywords: SNES, matrix-free, parameters
1048 
1049 .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD
1050           MatMFFDSetHHistory(), MatMFFDResetHHistory()
1051 @*/
1052 PetscErrorCode  MatMFFDSetFunctionError(Mat mat,PetscReal error)
1053 {
1054   PetscErrorCode ierr;
1055 
1056   PetscFunctionBegin;
1057   ierr = PetscTryMethod(mat,"MatMFFDSetFunctionError_C",(Mat,PetscReal),(mat,error));CHKERRQ(ierr);
1058   PetscFunctionReturn(0);
1059 }
1060 
1061 #undef __FUNCT__
1062 #define __FUNCT__ "MatMFFDAddNullSpace"
1063 /*@
1064    MatMFFDAddNullSpace - Provides a null space that an operator is
1065    supposed to have.  Since roundoff will create a small component in
1066    the null space, if you know the null space you may have it
1067    automatically removed.
1068 
1069    Logically Collective on Mat
1070 
1071    Input Parameters:
1072 +  J - the matrix-free matrix context
1073 -  nullsp - object created with MatNullSpaceCreate()
1074 
1075    Level: advanced
1076 
1077 .keywords: SNES, matrix-free, null space
1078 
1079 .seealso: MatNullSpaceCreate(), MatMFFDGetH(), MatCreateSNESMF(), MatCreateMFFD(), MATMFFD
1080           MatMFFDSetHHistory(), MatMFFDResetHHistory()
1081 @*/
1082 PetscErrorCode  MatMFFDAddNullSpace(Mat J,MatNullSpace nullsp)
1083 {
1084   PetscErrorCode ierr;
1085 
1086   PetscFunctionBegin;
1087   ierr = PetscTryMethod(J,"MatMFFDAddNullSpace_C",(Mat,MatNullSpace),(J,nullsp));CHKERRQ(ierr);
1088   PetscFunctionReturn(0);
1089 }
1090 
1091 #undef __FUNCT__
1092 #define __FUNCT__ "MatMFFDSetHHistory"
1093 /*@
1094    MatMFFDSetHHistory - Sets an array to collect a history of the
1095    differencing values (h) computed for the matrix-free product.
1096 
1097    Logically Collective on Mat
1098 
1099    Input Parameters:
1100 +  J - the matrix-free matrix context
1101 .  histroy - space to hold the history
1102 -  nhistory - number of entries in history, if more entries are generated than
1103               nhistory, then the later ones are discarded
1104 
1105    Level: advanced
1106 
1107    Notes:
1108    Use MatMFFDResetHHistory() to reset the history counter and collect
1109    a new batch of differencing parameters, h.
1110 
1111 .keywords: SNES, matrix-free, h history, differencing history
1112 
1113 .seealso: MatMFFDGetH(), MatCreateSNESMF(),
1114           MatMFFDResetHHistory(), MatMFFDSetFunctionError()
1115 
1116 @*/
1117 PetscErrorCode  MatMFFDSetHHistory(Mat J,PetscScalar history[],PetscInt nhistory)
1118 {
1119   MatMFFD        ctx = (MatMFFD)J->data;
1120   PetscErrorCode ierr;
1121   PetscBool      match;
1122 
1123   PetscFunctionBegin;
1124   ierr = PetscObjectTypeCompare((PetscObject)J,MATMFFD,&match);CHKERRQ(ierr);
1125   if (!match) SETERRQ(((PetscObject)J)->comm,PETSC_ERR_ARG_WRONG,"Not a MFFD matrix");
1126   ctx->historyh    = history;
1127   ctx->maxcurrenth = nhistory;
1128   ctx->currenth    = 0.;
1129   PetscFunctionReturn(0);
1130 }
1131 
1132 
1133 #undef __FUNCT__
1134 #define __FUNCT__ "MatMFFDResetHHistory"
1135 /*@
1136    MatMFFDResetHHistory - Resets the counter to zero to begin
1137    collecting a new set of differencing histories.
1138 
1139    Logically Collective on Mat
1140 
1141    Input Parameters:
1142 .  J - the matrix-free matrix context
1143 
1144    Level: advanced
1145 
1146    Notes:
1147    Use MatMFFDSetHHistory() to create the original history counter.
1148 
1149 .keywords: SNES, matrix-free, h history, differencing history
1150 
1151 .seealso: MatMFFDGetH(), MatCreateSNESMF(),
1152           MatMFFDSetHHistory(), MatMFFDSetFunctionError()
1153 
1154 @*/
1155 PetscErrorCode  MatMFFDResetHHistory(Mat J)
1156 {
1157   PetscErrorCode ierr;
1158 
1159   PetscFunctionBegin;
1160   ierr = PetscTryMethod(J,"MatMFFDResetHHistory_C",(Mat),(J));CHKERRQ(ierr);
1161   PetscFunctionReturn(0);
1162 }
1163 
1164 
1165 #undef __FUNCT__
1166 #define __FUNCT__ "MatMFFDSetBase"
1167 /*@
1168     MatMFFDSetBase - Sets the vector U at which matrix vector products of the
1169         Jacobian are computed
1170 
1171     Logically Collective on Mat
1172 
1173     Input Parameters:
1174 +   J - the MatMFFD matrix
1175 .   U - the vector
1176 -   F - (optional) vector that contains F(u) if it has been already computed
1177 
1178     Notes: This is rarely used directly
1179 
1180     If F is provided then it is not recomputed. Otherwise the function is evaluated at the base
1181     point during the first MatMult() after each call to MatMFFDSetBase().
1182 
1183     Level: advanced
1184 
1185 @*/
1186 PetscErrorCode  MatMFFDSetBase(Mat J,Vec U,Vec F)
1187 {
1188   PetscErrorCode ierr;
1189 
1190   PetscFunctionBegin;
1191   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
1192   PetscValidHeaderSpecific(U,VEC_CLASSID,2);
1193   if (F) PetscValidHeaderSpecific(F,VEC_CLASSID,3);
1194   ierr = PetscTryMethod(J,"MatMFFDSetBase_C",(Mat,Vec,Vec),(J,U,F));CHKERRQ(ierr);
1195   PetscFunctionReturn(0);
1196 }
1197 
1198 #undef __FUNCT__
1199 #define __FUNCT__ "MatMFFDSetCheckh"
1200 /*@C
1201     MatMFFDSetCheckh - Sets a function that checks the computed h and adjusts
1202         it to satisfy some criteria
1203 
1204     Logically Collective on Mat
1205 
1206     Input Parameters:
1207 +   J - the MatMFFD matrix
1208 .   fun - the function that checks h
1209 -   ctx - any context needed by the function
1210 
1211     Options Database Keys:
1212 .   -mat_mffd_check_positivity
1213 
1214     Level: advanced
1215 
1216     Notes: For example, MatMFFDSetCheckPositivity() insures that all entries
1217        of U + h*a are non-negative
1218 
1219 .seealso:  MatMFFDSetCheckPositivity()
1220 @*/
1221 PetscErrorCode  MatMFFDSetCheckh(Mat J,PetscErrorCode (*fun)(void*,Vec,Vec,PetscScalar*),void *ctx)
1222 {
1223   PetscErrorCode ierr;
1224 
1225   PetscFunctionBegin;
1226   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
1227   ierr = PetscTryMethod(J,"MatMFFDSetCheckh_C",(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*),(J,fun,ctx));CHKERRQ(ierr);
1228   PetscFunctionReturn(0);
1229 }
1230 
1231 #undef __FUNCT__
1232 #define __FUNCT__ "MatMFFDSetCheckPositivity"
1233 /*@
1234     MatMFFDCheckPositivity - Checks that all entries in U + h*a are positive or
1235         zero, decreases h until this is satisfied.
1236 
1237     Logically Collective on Vec
1238 
1239     Input Parameters:
1240 +   U - base vector that is added to
1241 .   a - vector that is added
1242 .   h - scaling factor on a
1243 -   dummy - context variable (unused)
1244 
1245     Options Database Keys:
1246 .   -mat_mffd_check_positivity
1247 
1248     Level: advanced
1249 
1250     Notes: This is rarely used directly, rather it is passed as an argument to
1251            MatMFFDSetCheckh()
1252 
1253 .seealso:  MatMFFDSetCheckh()
1254 @*/
1255 PetscErrorCode  MatMFFDCheckPositivity(void *dummy,Vec U,Vec a,PetscScalar *h)
1256 {
1257   PetscReal      val, minval;
1258   PetscScalar    *u_vec, *a_vec;
1259   PetscErrorCode ierr;
1260   PetscInt       i,n;
1261   MPI_Comm       comm;
1262 
1263   PetscFunctionBegin;
1264   ierr   = PetscObjectGetComm((PetscObject)U,&comm);CHKERRQ(ierr);
1265   ierr   = VecGetArray(U,&u_vec);CHKERRQ(ierr);
1266   ierr   = VecGetArray(a,&a_vec);CHKERRQ(ierr);
1267   ierr   = VecGetLocalSize(U,&n);CHKERRQ(ierr);
1268   minval = PetscAbsScalar(*h*1.01);
1269   for (i=0; i<n; i++) {
1270     if (PetscRealPart(u_vec[i] + *h*a_vec[i]) <= 0.0) {
1271       val = PetscAbsScalar(u_vec[i]/a_vec[i]);
1272       if (val < minval) minval = val;
1273     }
1274   }
1275   ierr = VecRestoreArray(U,&u_vec);CHKERRQ(ierr);
1276   ierr = VecRestoreArray(a,&a_vec);CHKERRQ(ierr);
1277   ierr = MPI_Allreduce(&minval,&val,1,MPIU_REAL,MPIU_MIN,comm);CHKERRQ(ierr);
1278   if (val <= PetscAbsScalar(*h)) {
1279     ierr = PetscInfo2(U,"Scaling back h from %G to %G\n",PetscRealPart(*h),.99*val);CHKERRQ(ierr);
1280     if (PetscRealPart(*h) > 0.0) *h =  0.99*val;
1281     else                         *h = -0.99*val;
1282   }
1283   PetscFunctionReturn(0);
1284 }
1285 
1286 
1287 
1288 
1289 
1290 
1291 
1292