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