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