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