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