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