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