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