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