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