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