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