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