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