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