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