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