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