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