xref: /petsc/src/mat/impls/mffd/mffd.c (revision 95cf157145500822f889726955674eed6376a2ab)
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   /* This "if" prevents PETSc from erroring when the mat is rectangular */
390   if ((ctx->vshift != 0.0) || (ctx->vscale != 1.0)) {
391     ierr = VecAXPBY(y,ctx->vshift,ctx->vscale,a);CHKERRQ(ierr);
392   }
393 
394   if (ctx->dlscale) {
395     ierr = VecPointwiseMult(y,ctx->dlscale,y);CHKERRQ(ierr);
396   }
397   if (ctx->dshift) {
398     if (!ctx->dshiftw) {
399       ierr = VecDuplicate(y,&ctx->dshiftw);CHKERRQ(ierr);
400     }
401     ierr = VecPointwiseMult(ctx->dshift,a,ctx->dshiftw);CHKERRQ(ierr);
402     ierr = VecAXPY(y,1.0,ctx->dshiftw);CHKERRQ(ierr);
403   }
404 
405   if (mat->nullsp) {ierr = MatNullSpaceRemove(mat->nullsp,y);CHKERRQ(ierr);}
406 
407   ierr = PetscLogEventEnd(MATMFFD_Mult,a,y,0,0);CHKERRQ(ierr);
408   PetscFunctionReturn(0);
409 }
410 
411 /*
412   MatGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix
413 
414         y ~= (F(u + ha) - F(u))/h,
415   where F = nonlinear function, as set by SNESSetFunction()
416         u = current iterate
417         h = difference interval
418 */
419 PetscErrorCode MatGetDiagonal_MFFD(Mat mat,Vec a)
420 {
421   MatMFFD        ctx = (MatMFFD)mat->data;
422   PetscScalar    h,*aa,*ww,v;
423   PetscReal      epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON;
424   Vec            w,U;
425   PetscErrorCode ierr;
426   PetscInt       i,rstart,rend;
427 
428   PetscFunctionBegin;
429   if (!ctx->func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Requires calling MatMFFDSetFunction() first");
430   if (!ctx->funci) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Requires calling MatMFFDSetFunctioni() first");
431   if (!ctx->funcisetbase) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Requires calling MatMFFDSetFunctioniBase() first");
432   w    = ctx->w;
433   U    = ctx->current_u;
434   ierr = (*ctx->func)(ctx->funcctx,U,a);CHKERRQ(ierr);
435   ierr = (*ctx->funcisetbase)(ctx->funcctx,U);CHKERRQ(ierr);
436   ierr = VecCopy(U,w);CHKERRQ(ierr);
437 
438   ierr = VecGetOwnershipRange(a,&rstart,&rend);CHKERRQ(ierr);
439   ierr = VecGetArray(a,&aa);CHKERRQ(ierr);
440   for (i=rstart; i<rend; i++) {
441     ierr = VecGetArray(w,&ww);CHKERRQ(ierr);
442     h    = ww[i-rstart];
443     if (h == 0.0) h = 1.0;
444     if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0)     h = umin;
445     else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin;
446     h *= epsilon;
447 
448     ww[i-rstart] += h;
449     ierr          = VecRestoreArray(w,&ww);CHKERRQ(ierr);
450     ierr          = (*ctx->funci)(ctx->funcctx,i,w,&v);CHKERRQ(ierr);
451     aa[i-rstart]  = (v - aa[i-rstart])/h;
452 
453     /* possibly shift and scale result */
454     if ((ctx->vshift != 0.0) || (ctx->vscale != 1.0)) {
455       aa[i - rstart] = ctx->vshift + ctx->vscale*aa[i-rstart];
456     }
457 
458     ierr          = VecGetArray(w,&ww);CHKERRQ(ierr);
459     ww[i-rstart] -= h;
460     ierr          = VecRestoreArray(w,&ww);CHKERRQ(ierr);
461   }
462   ierr = VecRestoreArray(a,&aa);CHKERRQ(ierr);
463   PetscFunctionReturn(0);
464 }
465 
466 static PetscErrorCode MatDiagonalScale_MFFD(Mat mat,Vec ll,Vec rr)
467 {
468   MatMFFD        aij = (MatMFFD)mat->data;
469   PetscErrorCode ierr;
470 
471   PetscFunctionBegin;
472   if (ll && !aij->dlscale) {
473     ierr = VecDuplicate(ll,&aij->dlscale);CHKERRQ(ierr);
474   }
475   if (rr && !aij->drscale) {
476     ierr = VecDuplicate(rr,&aij->drscale);CHKERRQ(ierr);
477   }
478   if (ll) {
479     ierr = VecCopy(ll,aij->dlscale);CHKERRQ(ierr);
480   }
481   if (rr) {
482     ierr = VecCopy(rr,aij->drscale);CHKERRQ(ierr);
483   }
484   PetscFunctionReturn(0);
485 }
486 
487 static PetscErrorCode MatDiagonalSet_MFFD(Mat mat,Vec ll,InsertMode mode)
488 {
489   MatMFFD        aij = (MatMFFD)mat->data;
490   PetscErrorCode ierr;
491 
492   PetscFunctionBegin;
493   if (mode == INSERT_VALUES) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"No diagonal set with INSERT_VALUES");
494   if (!aij->dshift) {
495     ierr = VecDuplicate(ll,&aij->dshift);CHKERRQ(ierr);
496   }
497   ierr = VecAXPY(aij->dshift,1.0,ll);CHKERRQ(ierr);
498   PetscFunctionReturn(0);
499 }
500 
501 static PetscErrorCode MatShift_MFFD(Mat Y,PetscScalar a)
502 {
503   MatMFFD shell = (MatMFFD)Y->data;
504 
505   PetscFunctionBegin;
506   shell->vshift += a;
507   PetscFunctionReturn(0);
508 }
509 
510 static PetscErrorCode MatScale_MFFD(Mat Y,PetscScalar a)
511 {
512   MatMFFD shell = (MatMFFD)Y->data;
513 
514   PetscFunctionBegin;
515   shell->vscale *= a;
516   PetscFunctionReturn(0);
517 }
518 
519 PETSC_EXTERN PetscErrorCode MatMFFDSetBase_MFFD(Mat J,Vec U,Vec F)
520 {
521   PetscErrorCode ierr;
522   MatMFFD        ctx = (MatMFFD)J->data;
523 
524   PetscFunctionBegin;
525   ierr = MatMFFDResetHHistory(J);CHKERRQ(ierr);
526   if (!ctx->current_u) {
527     ierr = VecDuplicate(U,&ctx->current_u);CHKERRQ(ierr);
528     ierr = VecLockReadPush(ctx->current_u);CHKERRQ(ierr);
529   }
530   ierr = VecLockReadPop(ctx->current_u);CHKERRQ(ierr);
531   ierr = VecCopy(U,ctx->current_u);CHKERRQ(ierr);
532   ierr = VecLockReadPush(ctx->current_u);CHKERRQ(ierr);
533   if (F) {
534     if (ctx->current_f_allocated) {ierr = VecDestroy(&ctx->current_f);CHKERRQ(ierr);}
535     ctx->current_f           = F;
536     ctx->current_f_allocated = PETSC_FALSE;
537   } else if (!ctx->current_f_allocated) {
538     ierr = MatCreateVecs(J,NULL,&ctx->current_f);CHKERRQ(ierr);
539 
540     ctx->current_f_allocated = PETSC_TRUE;
541   }
542   if (!ctx->w) {
543     ierr = VecDuplicate(ctx->current_u,&ctx->w);CHKERRQ(ierr);
544   }
545   J->assembled = PETSC_TRUE;
546   PetscFunctionReturn(0);
547 }
548 
549 typedef PetscErrorCode (*FCN3)(void*,Vec,Vec,PetscScalar*); /* force argument to next function to not be extern C*/
550 
551 static PetscErrorCode  MatMFFDSetCheckh_MFFD(Mat J,FCN3 fun,void *ectx)
552 {
553   MatMFFD ctx = (MatMFFD)J->data;
554 
555   PetscFunctionBegin;
556   ctx->checkh    = fun;
557   ctx->checkhctx = ectx;
558   PetscFunctionReturn(0);
559 }
560 
561 /*@C
562    MatMFFDSetOptionsPrefix - Sets the prefix used for searching for all
563    MatMFFD options in the database.
564 
565    Collective on Mat
566 
567    Input Parameter:
568 +  A - the Mat context
569 -  prefix - the prefix to prepend to all option names
570 
571    Notes:
572    A hyphen (-) must NOT be given at the beginning of the prefix name.
573    The first character of all runtime options is AUTOMATICALLY the hyphen.
574 
575    Level: advanced
576 
577 .keywords: SNES, matrix-free, parameters
578 
579 .seealso: MatSetFromOptions(), MatCreateSNESMF(), MatCreateMFFD()
580 @*/
581 PetscErrorCode  MatMFFDSetOptionsPrefix(Mat mat,const char prefix[])
582 
583 {
584   MatMFFD        mfctx = mat ? (MatMFFD)mat->data : (MatMFFD)NULL;
585   PetscErrorCode ierr;
586 
587   PetscFunctionBegin;
588   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
589   PetscValidHeaderSpecific(mfctx,MATMFFD_CLASSID,1);
590   ierr = PetscObjectSetOptionsPrefix((PetscObject)mfctx,prefix);CHKERRQ(ierr);
591   PetscFunctionReturn(0);
592 }
593 
594 static PetscErrorCode  MatSetFromOptions_MFFD(PetscOptionItems *PetscOptionsObject,Mat mat)
595 {
596   MatMFFD        mfctx = (MatMFFD)mat->data;
597   PetscErrorCode ierr;
598   PetscBool      flg;
599   char           ftype[256];
600 
601   PetscFunctionBegin;
602   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
603   PetscValidHeaderSpecific(mfctx,MATMFFD_CLASSID,1);
604   ierr = PetscObjectOptionsBegin((PetscObject)mfctx);CHKERRQ(ierr);
605   ierr = PetscOptionsFList("-mat_mffd_type","Matrix free type","MatMFFDSetType",MatMFFDList,((PetscObject)mfctx)->type_name,ftype,256,&flg);CHKERRQ(ierr);
606   if (flg) {
607     ierr = MatMFFDSetType(mat,ftype);CHKERRQ(ierr);
608   }
609 
610   ierr = PetscOptionsReal("-mat_mffd_err","set sqrt relative error in function","MatMFFDSetFunctionError",mfctx->error_rel,&mfctx->error_rel,0);CHKERRQ(ierr);
611   ierr = PetscOptionsInt("-mat_mffd_period","how often h is recomputed","MatMFFDSetPeriod",mfctx->recomputeperiod,&mfctx->recomputeperiod,0);CHKERRQ(ierr);
612 
613   flg  = PETSC_FALSE;
614   ierr = PetscOptionsBool("-mat_mffd_check_positivity","Insure that U + h*a is nonnegative","MatMFFDSetCheckh",flg,&flg,NULL);CHKERRQ(ierr);
615   if (flg) {
616     ierr = MatMFFDSetCheckh(mat,MatMFFDCheckPositivity,0);CHKERRQ(ierr);
617   }
618 #if defined(PETSC_USE_COMPLEX)
619   ierr = PetscOptionsBool("-mat_mffd_complex","Use Lyness complex number trick to compute the matrix-vector product","None",mfctx->usecomplex,&mfctx->usecomplex,NULL);CHKERRQ(ierr);
620 #endif
621   if (mfctx->ops->setfromoptions) {
622     ierr = (*mfctx->ops->setfromoptions)(PetscOptionsObject,mfctx);CHKERRQ(ierr);
623   }
624   ierr = PetscOptionsEnd();CHKERRQ(ierr);
625   PetscFunctionReturn(0);
626 }
627 
628 static PetscErrorCode  MatMFFDSetPeriod_MFFD(Mat mat,PetscInt period)
629 {
630   MatMFFD ctx = (MatMFFD)mat->data;
631 
632   PetscFunctionBegin;
633   ctx->recomputeperiod = period;
634   PetscFunctionReturn(0);
635 }
636 
637 static PetscErrorCode  MatMFFDSetFunction_MFFD(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx)
638 {
639   MatMFFD ctx = (MatMFFD)mat->data;
640 
641   PetscFunctionBegin;
642   ctx->func    = func;
643   ctx->funcctx = funcctx;
644   PetscFunctionReturn(0);
645 }
646 
647 static PetscErrorCode  MatMFFDSetFunctionError_MFFD(Mat mat,PetscReal error)
648 {
649   MatMFFD ctx = (MatMFFD)mat->data;
650 
651   PetscFunctionBegin;
652   if (error != PETSC_DEFAULT) ctx->error_rel = error;
653   PetscFunctionReturn(0);
654 }
655 
656 static PetscErrorCode MatMissingDiagonal_MFFD(Mat A,PetscBool  *missing,PetscInt *d)
657 {
658   PetscFunctionBegin;
659   *missing = PETSC_FALSE;
660   PetscFunctionReturn(0);
661 }
662 
663 /*MC
664   MATMFFD - MATMFFD = "mffd" - A matrix free matrix type.
665 
666   Level: advanced
667 
668 .seealso: MatCreateMFFD(), MatCreateSNESMF(), MatMFFDSetFunction(), MatMFFDSetType(),
669           MatMFFDSetFunctionError(), MatMFFDDSSetUmin(), MatMFFDSetFunction()
670           MatMFFDSetHHistory(), MatMFFDResetHHistory(), MatCreateSNESMF(),
671           MatMFFDGetH(),
672 M*/
673 PETSC_EXTERN PetscErrorCode MatCreate_MFFD(Mat A)
674 {
675   MatMFFD        mfctx;
676   PetscErrorCode ierr;
677 
678   PetscFunctionBegin;
679   ierr = MatMFFDInitializePackage();CHKERRQ(ierr);
680 
681   ierr = PetscHeaderCreate(mfctx,MATMFFD_CLASSID,"MatMFFD","Matrix-free Finite Differencing","Mat",PetscObjectComm((PetscObject)A),MatDestroy_MFFD,MatView_MFFD);CHKERRQ(ierr);
682 
683   mfctx->error_rel                = PETSC_SQRT_MACHINE_EPSILON;
684   mfctx->recomputeperiod          = 1;
685   mfctx->count                    = 0;
686   mfctx->currenth                 = 0.0;
687   mfctx->historyh                 = NULL;
688   mfctx->ncurrenth                = 0;
689   mfctx->maxcurrenth              = 0;
690   ((PetscObject)mfctx)->type_name = 0;
691 
692   mfctx->vshift = 0.0;
693   mfctx->vscale = 1.0;
694 
695   /*
696      Create the empty data structure to contain compute-h routines.
697      These will be filled in below from the command line options or
698      a later call with MatMFFDSetType() or if that is not called
699      then it will default in the first use of MatMult_MFFD()
700   */
701   mfctx->ops->compute        = 0;
702   mfctx->ops->destroy        = 0;
703   mfctx->ops->view           = 0;
704   mfctx->ops->setfromoptions = 0;
705   mfctx->hctx                = 0;
706 
707   mfctx->func    = 0;
708   mfctx->funcctx = 0;
709   mfctx->w       = NULL;
710 
711   A->data = mfctx;
712 
713   A->ops->mult            = MatMult_MFFD;
714   A->ops->destroy         = MatDestroy_MFFD;
715   A->ops->view            = MatView_MFFD;
716   A->ops->assemblyend     = MatAssemblyEnd_MFFD;
717   A->ops->scale           = MatScale_MFFD;
718   A->ops->shift           = MatShift_MFFD;
719   A->ops->diagonalscale   = MatDiagonalScale_MFFD;
720   A->ops->diagonalset     = MatDiagonalSet_MFFD;
721   A->ops->setfromoptions  = MatSetFromOptions_MFFD;
722   A->ops->missingdiagonal = MatMissingDiagonal_MFFD;
723   A->assembled            = PETSC_TRUE;
724 
725   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetBase_C",MatMFFDSetBase_MFFD);CHKERRQ(ierr);
726   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctioniBase_C",MatMFFDSetFunctioniBase_MFFD);CHKERRQ(ierr);
727   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctioni_C",MatMFFDSetFunctioni_MFFD);CHKERRQ(ierr);
728   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunction_C",MatMFFDSetFunction_MFFD);CHKERRQ(ierr);
729   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetCheckh_C",MatMFFDSetCheckh_MFFD);CHKERRQ(ierr);
730   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetPeriod_C",MatMFFDSetPeriod_MFFD);CHKERRQ(ierr);
731   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctionError_C",MatMFFDSetFunctionError_MFFD);CHKERRQ(ierr);
732   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDResetHHistory_C",MatMFFDResetHHistory_MFFD);CHKERRQ(ierr);
733 
734   mfctx->mat = A;
735 
736   ierr = PetscObjectChangeTypeName((PetscObject)A,MATMFFD);CHKERRQ(ierr);
737   PetscFunctionReturn(0);
738 }
739 
740 /*@
741    MatCreateMFFD - Creates a matrix-free matrix. See also MatCreateSNESMF()
742 
743    Collective on Vec
744 
745    Input Parameters:
746 +  comm - MPI communicator
747 .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
748            This value should be the same as the local size used in creating the
749            y vector for the matrix-vector product y = Ax.
750 .  n - This value should be the same as the local size used in creating the
751        x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have
752        calculated if N is given) For square matrices n is almost always m.
753 .  M - number of global rows (or PETSC_DETERMINE to have calculated if m is given)
754 -  N - number of global columns (or PETSC_DETERMINE to have calculated if n is given)
755 
756 
757    Output Parameter:
758 .  J - the matrix-free matrix
759 
760    Options Database Keys: call MatSetFromOptions() to trigger these
761 +  -mat_mffd_type - wp or ds (see MATMFFD_WP or MATMFFD_DS)
762 .  -mat_mffd_err - square root of estimated relative error in function evaluation
763 .  -mat_mffd_period - how often h is recomputed, defaults to 1, everytime
764 .  -mat_mffd_check_positivity - possibly decrease h until U + h*a has only positive values
765 -  -mat_mffd_complex - use the Lyness trick with complex numbers to compute the matrix-vector product instead of differencing
766                        (requires real valued functions but that PETSc be configured for complex numbers)
767 
768 
769    Level: advanced
770 
771    Notes:
772    The matrix-free matrix context merely contains the function pointers
773    and work space for performing finite difference approximations of
774    Jacobian-vector products, F'(u)*a,
775 
776    The default code uses the following approach to compute h
777 
778 .vb
779      F'(u)*a = [F(u+h*a) - F(u)]/h where
780      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
781        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
782  where
783      error_rel = square root of relative error in function evaluation
784      umin = minimum iterate parameter
785 .ve
786 
787    You can call SNESSetJacobian() with MatMFFDComputeJacobian() if you are using matrix and not a different
788    preconditioner matrix
789 
790    The user can set the error_rel via MatMFFDSetFunctionError() and
791    umin via MatMFFDDSSetUmin(); see Users-Manual: ch_snes for details.
792 
793    The user should call MatDestroy() when finished with the matrix-free
794    matrix context.
795 
796    Options Database Keys:
797 +  -mat_mffd_err <error_rel> - Sets error_rel
798 .  -mat_mffd_unim <umin> - Sets umin (for default PETSc routine that computes h only)
799 -  -mat_mffd_check_positivity
800 
801 .keywords: default, matrix-free, create, matrix
802 
803 .seealso: MatDestroy(), MatMFFDSetFunctionError(), MatMFFDDSSetUmin(), MatMFFDSetFunction()
804           MatMFFDSetHHistory(), MatMFFDResetHHistory(), MatCreateSNESMF(),
805           MatMFFDGetH(), MatMFFDRegister(), MatMFFDComputeJacobian()
806 
807 @*/
808 PetscErrorCode  MatCreateMFFD(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *J)
809 {
810   PetscErrorCode ierr;
811 
812   PetscFunctionBegin;
813   ierr = MatCreate(comm,J);CHKERRQ(ierr);
814   ierr = MatSetSizes(*J,m,n,M,N);CHKERRQ(ierr);
815   ierr = MatSetType(*J,MATMFFD);CHKERRQ(ierr);
816   ierr = MatSetUp(*J);CHKERRQ(ierr);
817   PetscFunctionReturn(0);
818 }
819 
820 /*@
821    MatMFFDGetH - Gets the last value that was used as the differencing
822    parameter.
823 
824    Not Collective
825 
826    Input Parameters:
827 .  mat - the matrix obtained with MatCreateSNESMF()
828 
829    Output Paramter:
830 .  h - the differencing step size
831 
832    Level: advanced
833 
834 .keywords: SNES, matrix-free, parameters
835 
836 .seealso: MatCreateSNESMF(),MatMFFDSetHHistory(), MatCreateMFFD(), MATMFFD, MatMFFDResetHHistory()
837 @*/
838 PetscErrorCode  MatMFFDGetH(Mat mat,PetscScalar *h)
839 {
840   MatMFFD        ctx = (MatMFFD)mat->data;
841   PetscErrorCode ierr;
842   PetscBool      match;
843 
844   PetscFunctionBegin;
845   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
846   PetscValidPointer(h,2);
847   ierr = PetscObjectTypeCompare((PetscObject)mat,MATMFFD,&match);CHKERRQ(ierr);
848   if (!match) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONG,"Not a MFFD matrix");
849 
850   *h = ctx->currenth;
851   PetscFunctionReturn(0);
852 }
853 
854 /*@C
855    MatMFFDSetFunction - Sets the function used in applying the matrix free.
856 
857    Logically Collective on Mat
858 
859    Input Parameters:
860 +  mat - the matrix free matrix created via MatCreateSNESMF() or MatCreateMFFD()
861 .  func - the function to use
862 -  funcctx - optional function context passed to function
863 
864    Calling Sequence of func:
865 $     func (void *funcctx, Vec x, Vec f)
866 
867 +  funcctx - user provided context
868 .  x - input vector
869 -  f - computed output function
870 
871    Level: advanced
872 
873    Notes:
874     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
875     matrix inside your compute Jacobian routine
876 
877     If this is not set then it will use the function set with SNESSetFunction() if MatCreateSNESMF() was used.
878 
879 .keywords: SNES, matrix-free, function
880 
881 .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD,
882           MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction()
883 @*/
884 PetscErrorCode  MatMFFDSetFunction(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx)
885 {
886   PetscErrorCode ierr;
887 
888   PetscFunctionBegin;
889   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
890   ierr = PetscTryMethod(mat,"MatMFFDSetFunction_C",(Mat,PetscErrorCode (*)(void*,Vec,Vec),void*),(mat,func,funcctx));CHKERRQ(ierr);
891   PetscFunctionReturn(0);
892 }
893 
894 /*@C
895    MatMFFDSetFunctioni - Sets the function for a single component
896 
897    Logically Collective on Mat
898 
899    Input Parameters:
900 +  mat - the matrix free matrix created via MatCreateSNESMF()
901 -  funci - the function to use
902 
903    Level: advanced
904 
905    Notes:
906     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
907     matrix inside your compute Jacobian routine.
908     This function is necessary to compute the diagonal of the matrix.
909     funci must not contain any MPI call as it is called inside a loop on the local portion of the vector.
910 
911 .keywords: SNES, matrix-free, function
912 
913 .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction(), MatGetDiagonal()
914 
915 @*/
916 PetscErrorCode  MatMFFDSetFunctioni(Mat mat,PetscErrorCode (*funci)(void*,PetscInt,Vec,PetscScalar*))
917 {
918   PetscErrorCode ierr;
919 
920   PetscFunctionBegin;
921   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
922   ierr = PetscTryMethod(mat,"MatMFFDSetFunctioni_C",(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)),(mat,funci));CHKERRQ(ierr);
923   PetscFunctionReturn(0);
924 }
925 
926 /*@C
927    MatMFFDSetFunctioniBase - Sets the base vector for a single component function evaluation
928 
929    Logically Collective on Mat
930 
931    Input Parameters:
932 +  mat - the matrix free matrix created via MatCreateSNESMF()
933 -  func - the function to use
934 
935    Level: advanced
936 
937    Notes:
938     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
939     matrix inside your compute Jacobian routine.
940     This function is necessary to compute the diagonal of the matrix.
941 
942 
943 .keywords: SNES, matrix-free, function
944 
945 .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD
946           MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction(), MatGetDiagonal()
947 @*/
948 PetscErrorCode  MatMFFDSetFunctioniBase(Mat mat,PetscErrorCode (*func)(void*,Vec))
949 {
950   PetscErrorCode ierr;
951 
952   PetscFunctionBegin;
953   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
954   ierr = PetscTryMethod(mat,"MatMFFDSetFunctioniBase_C",(Mat,PetscErrorCode (*)(void*,Vec)),(mat,func));CHKERRQ(ierr);
955   PetscFunctionReturn(0);
956 }
957 
958 /*@
959    MatMFFDSetPeriod - Sets how often h is recomputed, by default it is everytime
960 
961    Logically Collective on Mat
962 
963    Input Parameters:
964 +  mat - the matrix free matrix created via MatCreateSNESMF()
965 -  period - 1 for everytime, 2 for every second etc
966 
967    Options Database Keys:
968 +  -mat_mffd_period <period>
969 
970    Level: advanced
971 
972 
973 .keywords: SNES, matrix-free, parameters
974 
975 .seealso: MatCreateSNESMF(),MatMFFDGetH(),
976           MatMFFDSetHHistory(), MatMFFDResetHHistory()
977 @*/
978 PetscErrorCode  MatMFFDSetPeriod(Mat mat,PetscInt period)
979 {
980   PetscErrorCode ierr;
981 
982   PetscFunctionBegin;
983   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
984   PetscValidLogicalCollectiveInt(mat,period,2);
985   ierr = PetscTryMethod(mat,"MatMFFDSetPeriod_C",(Mat,PetscInt),(mat,period));CHKERRQ(ierr);
986   PetscFunctionReturn(0);
987 }
988 
989 /*@
990    MatMFFDSetFunctionError - Sets the error_rel for the approximation of
991    matrix-vector products using finite differences.
992 
993    Logically Collective on Mat
994 
995    Input Parameters:
996 +  mat - the matrix free matrix created via MatCreateMFFD() or MatCreateSNESMF()
997 -  error_rel - relative error (should be set to the square root of
998                the relative error in the function evaluations)
999 
1000    Options Database Keys:
1001 +  -mat_mffd_err <error_rel> - Sets error_rel
1002 
1003    Level: advanced
1004 
1005    Notes:
1006    The default matrix-free matrix-vector product routine computes
1007 .vb
1008      F'(u)*a = [F(u+h*a) - F(u)]/h where
1009      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
1010        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   else
1011 .ve
1012 
1013 .keywords: SNES, matrix-free, parameters
1014 
1015 .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD
1016           MatMFFDSetHHistory(), MatMFFDResetHHistory()
1017 @*/
1018 PetscErrorCode  MatMFFDSetFunctionError(Mat mat,PetscReal error)
1019 {
1020   PetscErrorCode ierr;
1021 
1022   PetscFunctionBegin;
1023   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
1024   PetscValidLogicalCollectiveReal(mat,error,2);
1025   ierr = PetscTryMethod(mat,"MatMFFDSetFunctionError_C",(Mat,PetscReal),(mat,error));CHKERRQ(ierr);
1026   PetscFunctionReturn(0);
1027 }
1028 
1029 /*@
1030    MatMFFDSetHHistory - Sets an array to collect a history of the
1031    differencing values (h) computed for the matrix-free product.
1032 
1033    Logically Collective on Mat
1034 
1035    Input Parameters:
1036 +  J - the matrix-free matrix context
1037 .  histroy - space to hold the history
1038 -  nhistory - number of entries in history, if more entries are generated than
1039               nhistory, then the later ones are discarded
1040 
1041    Level: advanced
1042 
1043    Notes:
1044    Use MatMFFDResetHHistory() to reset the history counter and collect
1045    a new batch of differencing parameters, h.
1046 
1047 .keywords: SNES, matrix-free, h history, differencing history
1048 
1049 .seealso: MatMFFDGetH(), MatCreateSNESMF(),
1050           MatMFFDResetHHistory(), MatMFFDSetFunctionError()
1051 
1052 @*/
1053 PetscErrorCode  MatMFFDSetHHistory(Mat J,PetscScalar history[],PetscInt nhistory)
1054 {
1055   MatMFFD        ctx = (MatMFFD)J->data;
1056   PetscErrorCode ierr;
1057   PetscBool      match;
1058 
1059   PetscFunctionBegin;
1060   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
1061   if (history) PetscValidPointer(history,2);
1062   PetscValidLogicalCollectiveInt(J,nhistory,3);
1063   ierr = PetscObjectTypeCompare((PetscObject)J,MATMFFD,&match);CHKERRQ(ierr);
1064   if (!match) SETERRQ(PetscObjectComm((PetscObject)J),PETSC_ERR_ARG_WRONG,"Not a MFFD matrix");
1065   ctx->historyh    = history;
1066   ctx->maxcurrenth = nhistory;
1067   ctx->currenth    = 0.;
1068   PetscFunctionReturn(0);
1069 }
1070 
1071 /*@
1072    MatMFFDResetHHistory - Resets the counter to zero to begin
1073    collecting a new set of differencing histories.
1074 
1075    Logically Collective on Mat
1076 
1077    Input Parameters:
1078 .  J - the matrix-free matrix context
1079 
1080    Level: advanced
1081 
1082    Notes:
1083    Use MatMFFDSetHHistory() to create the original history counter.
1084 
1085 .keywords: SNES, matrix-free, h history, differencing history
1086 
1087 .seealso: MatMFFDGetH(), MatCreateSNESMF(),
1088           MatMFFDSetHHistory(), MatMFFDSetFunctionError()
1089 
1090 @*/
1091 PetscErrorCode  MatMFFDResetHHistory(Mat J)
1092 {
1093   PetscErrorCode ierr;
1094 
1095   PetscFunctionBegin;
1096   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
1097   ierr = PetscTryMethod(J,"MatMFFDResetHHistory_C",(Mat),(J));CHKERRQ(ierr);
1098   PetscFunctionReturn(0);
1099 }
1100 
1101 /*@
1102     MatMFFDSetBase - Sets the vector U at which matrix vector products of the
1103         Jacobian are computed
1104 
1105     Logically Collective on Mat
1106 
1107     Input Parameters:
1108 +   J - the MatMFFD matrix
1109 .   U - the vector
1110 -   F - (optional) vector that contains F(u) if it has been already computed
1111 
1112     Notes:
1113     This is rarely used directly
1114 
1115     If F is provided then it is not recomputed. Otherwise the function is evaluated at the base
1116     point during the first MatMult() after each call to MatMFFDSetBase().
1117 
1118     Level: advanced
1119 
1120 @*/
1121 PetscErrorCode  MatMFFDSetBase(Mat J,Vec U,Vec F)
1122 {
1123   PetscErrorCode ierr;
1124 
1125   PetscFunctionBegin;
1126   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
1127   PetscValidHeaderSpecific(U,VEC_CLASSID,2);
1128   if (F) PetscValidHeaderSpecific(F,VEC_CLASSID,3);
1129   ierr = PetscTryMethod(J,"MatMFFDSetBase_C",(Mat,Vec,Vec),(J,U,F));CHKERRQ(ierr);
1130   PetscFunctionReturn(0);
1131 }
1132 
1133 /*@C
1134     MatMFFDSetCheckh - Sets a function that checks the computed h and adjusts
1135         it to satisfy some criteria
1136 
1137     Logically Collective on Mat
1138 
1139     Input Parameters:
1140 +   J - the MatMFFD matrix
1141 .   fun - the function that checks h
1142 -   ctx - any context needed by the function
1143 
1144     Options Database Keys:
1145 .   -mat_mffd_check_positivity
1146 
1147     Level: advanced
1148 
1149     Notes:
1150     For example, MatMFFDCheckPositivity() insures that all entries
1151        of U + h*a are non-negative
1152 
1153      The function you provide is called after the default h has been computed and allows you to
1154      modify it.
1155 
1156 .seealso:  MatMFFDCheckPositivity()
1157 @*/
1158 PetscErrorCode  MatMFFDSetCheckh(Mat J,PetscErrorCode (*fun)(void*,Vec,Vec,PetscScalar*),void *ctx)
1159 {
1160   PetscErrorCode ierr;
1161 
1162   PetscFunctionBegin;
1163   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
1164   ierr = PetscTryMethod(J,"MatMFFDSetCheckh_C",(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*),(J,fun,ctx));CHKERRQ(ierr);
1165   PetscFunctionReturn(0);
1166 }
1167 
1168 /*@
1169     MatMFFDCheckPositivity - Checks that all entries in U + h*a are positive or
1170         zero, decreases h until this is satisfied.
1171 
1172     Logically Collective on Vec
1173 
1174     Input Parameters:
1175 +   U - base vector that is added to
1176 .   a - vector that is added
1177 .   h - scaling factor on a
1178 -   dummy - context variable (unused)
1179 
1180     Options Database Keys:
1181 .   -mat_mffd_check_positivity
1182 
1183     Level: advanced
1184 
1185     Notes:
1186     This is rarely used directly, rather it is passed as an argument to
1187            MatMFFDSetCheckh()
1188 
1189 .seealso:  MatMFFDSetCheckh()
1190 @*/
1191 PetscErrorCode  MatMFFDCheckPositivity(void *dummy,Vec U,Vec a,PetscScalar *h)
1192 {
1193   PetscReal      val, minval;
1194   PetscScalar    *u_vec, *a_vec;
1195   PetscErrorCode ierr;
1196   PetscInt       i,n;
1197   MPI_Comm       comm;
1198 
1199   PetscFunctionBegin;
1200   PetscValidHeaderSpecific(U,VEC_CLASSID,2);
1201   PetscValidHeaderSpecific(a,VEC_CLASSID,3);
1202   PetscValidPointer(h,4);
1203   ierr   = PetscObjectGetComm((PetscObject)U,&comm);CHKERRQ(ierr);
1204   ierr   = VecGetArray(U,&u_vec);CHKERRQ(ierr);
1205   ierr   = VecGetArray(a,&a_vec);CHKERRQ(ierr);
1206   ierr   = VecGetLocalSize(U,&n);CHKERRQ(ierr);
1207   minval = PetscAbsScalar(*h)*PetscRealConstant(1.01);
1208   for (i=0; i<n; i++) {
1209     if (PetscRealPart(u_vec[i] + *h*a_vec[i]) <= 0.0) {
1210       val = PetscAbsScalar(u_vec[i]/a_vec[i]);
1211       if (val < minval) minval = val;
1212     }
1213   }
1214   ierr = VecRestoreArray(U,&u_vec);CHKERRQ(ierr);
1215   ierr = VecRestoreArray(a,&a_vec);CHKERRQ(ierr);
1216   ierr = MPIU_Allreduce(&minval,&val,1,MPIU_REAL,MPIU_MIN,comm);CHKERRQ(ierr);
1217   if (val <= PetscAbsScalar(*h)) {
1218     ierr = PetscInfo2(U,"Scaling back h from %g to %g\n",(double)PetscRealPart(*h),(double)(.99*val));CHKERRQ(ierr);
1219     if (PetscRealPart(*h) > 0.0) *h =  0.99*val;
1220     else                         *h = -0.99*val;
1221   }
1222   PetscFunctionReturn(0);
1223 }
1224