xref: /petsc/src/snes/mf/snesmfj.c (revision d34347eca8ac9657cccc548b3ba1c8201ddea9cf)
1 /*$Id: snesmfj.c,v 1.131 2001/09/05 18:45:40 bsmith Exp $*/
2 
3 #include "src/mat/matimpl.h"
4 #include "src/snes/mf/snesmfj.h"   /*I  "petscsnes.h"   I*/
5 
6 PetscFList      MatSNESMPetscFList              = 0;
7 PetscTruth MatSNESMFRegisterAllCalled = PETSC_FALSE;
8 
9 #undef __FUNCT__
10 #define __FUNCT__ "MatSNESMFSetType"
11 /*@C
12     MatSNESMFSetType - Sets the method that is used to compute the
13     differencing parameter for finite differene matrix-free formulations.
14 
15     Input Parameters:
16 +   mat - the "matrix-free" matrix created via MatCreateSNESMF(), or MatCreateMF()
17           or MatSetType(mat,MATMFFD);
18 -   ftype - the type requested
19 
20     Level: advanced
21 
22     Notes:
23     For example, such routines can compute h for use in
24     Jacobian-vector products of the form
25 
26                         F(x+ha) - F(x)
27           F'(u)a  ~=  ----------------
28                               h
29 
30 .seealso: MatCreateSNESMF(), MatSNESMFRegisterDynamic)
31 @*/
32 int MatSNESMFSetType(Mat mat,MatSNESMFType ftype)
33 {
34   int          ierr,(*r)(MatSNESMFCtx);
35   MatSNESMFCtx ctx = (MatSNESMFCtx)mat->data;
36   PetscTruth   match;
37 
38   PetscFunctionBegin;
39   PetscValidHeaderSpecific(mat,MAT_COOKIE);
40   PetscValidCharPointer(ftype);
41 
42   /* already set, so just return */
43   ierr = PetscTypeCompare((PetscObject)ctx,ftype,&match);CHKERRQ(ierr);
44   if (match) PetscFunctionReturn(0);
45 
46   /* destroy the old one if it exists */
47   if (ctx->ops->destroy) {
48     ierr = (*ctx->ops->destroy)(ctx);CHKERRQ(ierr);
49   }
50 
51   /* Get the function pointers for the requrested method */
52   if (!MatSNESMFRegisterAllCalled) {ierr = MatSNESMFRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
53 
54   ierr =  PetscFListFind(ctx->comm,MatSNESMPetscFList,ftype,(void (**)(void)) &r);CHKERRQ(ierr);
55 
56   if (!r) SETERRQ(1,"Unknown MatSNESMF type given");
57 
58   ierr = (*r)(ctx);CHKERRQ(ierr);
59 
60   ierr = PetscObjectChangeTypeName((PetscObject)ctx,ftype);CHKERRQ(ierr);
61 
62   PetscFunctionReturn(0);
63 }
64 
65 EXTERN_C_BEGIN
66 #undef __FUNCT__
67 #define __FUNCT__ "MatSNESMFSetFunctioniBase_FD"
68 int MatSNESMFSetFunctioniBase_FD(Mat mat,int (*func)(Vec,void *))
69 {
70   MatSNESMFCtx ctx = (MatSNESMFCtx)mat->data;
71 
72   PetscFunctionBegin;
73   ctx->funcisetbase = func;
74   PetscFunctionReturn(0);
75 }
76 EXTERN_C_END
77 
78 EXTERN_C_BEGIN
79 #undef __FUNCT__
80 #define __FUNCT__ "MatSNESMFSetFunctioni_FD"
81 int MatSNESMFSetFunctioni_FD(Mat mat,int (*funci)(int,Vec,PetscScalar*,void *))
82 {
83   MatSNESMFCtx ctx = (MatSNESMFCtx)mat->data;
84 
85   PetscFunctionBegin;
86   ctx->funci = funci;
87   PetscFunctionReturn(0);
88 }
89 EXTERN_C_END
90 
91 
92 #undef __FUNCT__
93 #define __FUNCT__ "MatSNESMFRegister"
94 int MatSNESMFRegister(const char sname[],const char path[],const char name[],int (*function)(MatSNESMFCtx))
95 {
96   int ierr;
97   char fullname[256];
98 
99   PetscFunctionBegin;
100   ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
101   ierr = PetscFListAdd(&MatSNESMPetscFList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
102   PetscFunctionReturn(0);
103 }
104 
105 
106 #undef __FUNCT__
107 #define __FUNCT__ "MatSNESMFRegisterDestroy"
108 /*@C
109    MatSNESMFRegisterDestroy - Frees the list of MatSNESMF methods that were
110    registered by MatSNESMFRegisterDynamic).
111 
112    Not Collective
113 
114    Level: developer
115 
116 .keywords: MatSNESMF, register, destroy
117 
118 .seealso: MatSNESMFRegisterDynamic), MatSNESMFRegisterAll()
119 @*/
120 int MatSNESMFRegisterDestroy(void)
121 {
122   int ierr;
123 
124   PetscFunctionBegin;
125   if (MatSNESMPetscFList) {
126     ierr = PetscFListDestroy(&MatSNESMPetscFList);CHKERRQ(ierr);
127     MatSNESMPetscFList = 0;
128   }
129   MatSNESMFRegisterAllCalled = PETSC_FALSE;
130   PetscFunctionReturn(0);
131 }
132 
133 /* ----------------------------------------------------------------------------------------*/
134 #undef __FUNCT__
135 #define __FUNCT__ "MatDestroy_MFFD"
136 int MatDestroy_MFFD(Mat mat)
137 {
138   int          ierr;
139   MatSNESMFCtx ctx = (MatSNESMFCtx)mat->data;
140 
141   PetscFunctionBegin;
142   if (ctx->w != PETSC_NULL) {
143     ierr = VecDestroy(ctx->w);CHKERRQ(ierr);
144   }
145   if (ctx->ops->destroy) {ierr = (*ctx->ops->destroy)(ctx);CHKERRQ(ierr);}
146   if (ctx->sp) {ierr = MatNullSpaceDestroy(ctx->sp);CHKERRQ(ierr);}
147   PetscHeaderDestroy(ctx);
148   PetscFunctionReturn(0);
149 }
150 
151 #undef __FUNCT__
152 #define __FUNCT__ "MatView_MFFD"
153 /*
154    MatSNESMFView_MFFD - Views matrix-free parameters.
155 
156 */
157 int MatView_MFFD(Mat J,PetscViewer viewer)
158 {
159   int          ierr;
160   MatSNESMFCtx ctx = (MatSNESMFCtx)J->data;
161   PetscTruth   isascii;
162 
163   PetscFunctionBegin;
164   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
165   if (isascii) {
166      ierr = PetscViewerASCIIPrintf(viewer,"  SNES matrix-free approximation:\n");CHKERRQ(ierr);
167      ierr = PetscViewerASCIIPrintf(viewer,"    err=%g (relative error in function evaluation)\n",ctx->error_rel);CHKERRQ(ierr);
168      if (!ctx->type_name) {
169        ierr = PetscViewerASCIIPrintf(viewer,"    The compute h routine has not yet been set\n");CHKERRQ(ierr);
170      } else {
171        ierr = PetscViewerASCIIPrintf(viewer,"    Using %s compute h routine\n",ctx->type_name);CHKERRQ(ierr);
172      }
173      if (ctx->ops->view) {
174        ierr = (*ctx->ops->view)(ctx,viewer);CHKERRQ(ierr);
175      }
176   } else {
177     SETERRQ1(1,"Viewer type %s not supported for SNES matrix free matrix",((PetscObject)viewer)->type_name);
178   }
179   PetscFunctionReturn(0);
180 }
181 
182 #undef __FUNCT__
183 #define __FUNCT__ "MatAssemblyEnd_MFFD"
184 /*
185    MatSNESMFAssemblyEnd_Private - Resets the ctx->ncurrenth to zero. This
186    allows the user to indicate the beginning of a new linear solve by calling
187    MatAssemblyXXX() on the matrix free matrix. This then allows the
188    MatSNESMFCreate_WP() to properly compute ||U|| only the first time
189    in the linear solver rather than every time.
190 */
191 int MatAssemblyEnd_MFFD(Mat J,MatAssemblyType mt)
192 {
193   int             ierr;
194   MatSNESMFCtx    j = (MatSNESMFCtx)J->data;
195 
196   PetscFunctionBegin;
197   ierr = MatSNESMFResetHHistory(J);CHKERRQ(ierr);
198   if (j->usesnes) {
199     ierr = SNESGetSolution(j->snes,&j->current_u);CHKERRQ(ierr);
200     ierr = SNESGetFunction(j->snes,&j->current_f,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
201     if (j->w == PETSC_NULL) {
202       ierr = VecDuplicate(j->current_u, &j->w);CHKERRQ(ierr);
203     }
204   }
205   j->vshift = 0.0;
206   j->vscale = 1.0;
207   PetscFunctionReturn(0);
208 }
209 
210 #undef __FUNCT__
211 #define __FUNCT__ "MatMult_MFFD"
212 /*
213   MatMult_MFFD - Default matrix-free form for Jacobian-vector product, y = F'(u)*a:
214 
215         y ~= (F(u + ha) - F(u))/h,
216   where F = nonlinear function, as set by SNESSetFunction()
217         u = current iterate
218         h = difference interval
219 */
220 int MatMult_MFFD(Mat mat,Vec a,Vec y)
221 {
222   MatSNESMFCtx    ctx = (MatSNESMFCtx)mat->data;
223   SNES            snes;
224   PetscScalar     h,mone = -1.0;
225   Vec             w,U,F;
226   int             ierr,(*eval_fct)(SNES,Vec,Vec)=0;
227 
228   PetscFunctionBegin;
229   /* We log matrix-free matrix-vector products separately, so that we can
230      separate the performance monitoring from the cases that use conventional
231      storage.  We may eventually modify event logging to associate events
232      with particular objects, hence alleviating the more general problem. */
233   ierr = PetscLogEventBegin(MAT_MultMatrixFree,a,y,0,0);CHKERRQ(ierr);
234 
235   snes = ctx->snes;
236   w    = ctx->w;
237   U    = ctx->current_u;
238 
239   /*
240       Compute differencing parameter
241   */
242   if (!ctx->ops->compute) {
243     ierr = MatSNESMFSetType(mat,MATSNESMF_WP);CHKERRQ(ierr);
244     ierr = MatSNESMFSetFromOptions(mat);CHKERRQ(ierr);
245   }
246   ierr = (*ctx->ops->compute)(ctx,U,a,&h);CHKERRQ(ierr);
247 
248   if (ctx->checkh) {
249     ierr = (*ctx->checkh)(U,a,&h,ctx->checkhctx);CHKERRQ(ierr);
250   }
251 
252   /* keep a record of the current differencing parameter h */
253   ctx->currenth = h;
254 #if defined(PETSC_USE_COMPLEX)
255   PetscLogInfo(mat,"MatMult_MFFD:Current differencing parameter: %g + %g i\n",PetscRealPart(h),PetscImaginaryPart(h));
256 #else
257   PetscLogInfo(mat,"MatMult_MFFD:Current differencing parameter: %15.12e\n",h);
258 #endif
259   if (ctx->historyh && ctx->ncurrenth < ctx->maxcurrenth) {
260     ctx->historyh[ctx->ncurrenth] = h;
261   }
262   ctx->ncurrenth++;
263 
264   /* w = u + ha */
265   ierr = VecWAXPY(&h,a,U,w);CHKERRQ(ierr);
266 
267   if (ctx->usesnes) {
268     eval_fct = SNESComputeFunction;
269     F    = ctx->current_f;
270     if (!F) SETERRQ(1,"You must call MatAssembly() even on matrix-free matrices");
271     ierr = (*eval_fct)(snes,w,y);CHKERRQ(ierr);
272   } else {
273     F = ctx->funcvec;
274     /* compute func(U) as base for differencing */
275     if (ctx->ncurrenth == 1) {
276       ierr = (*ctx->func)(snes,U,F,ctx->funcctx);CHKERRQ(ierr);
277     }
278     ierr = (*ctx->func)(snes,w,y,ctx->funcctx);CHKERRQ(ierr);
279   }
280 
281   ierr = VecAXPY(&mone,F,y);CHKERRQ(ierr);
282   h    = 1.0/h;
283   ierr = VecScale(&h,y);CHKERRQ(ierr);
284 
285 
286   if (ctx->vshift != 0.0 && ctx->vscale != 1.0) {
287     ierr = VecAXPBY(&ctx->vshift,&ctx->vscale,a,y);CHKERRQ(ierr);
288   } else if (ctx->vscale != 1.0) {
289     ierr = VecScale(&ctx->vscale,y);CHKERRQ(ierr);
290   } else if (ctx->vshift != 0.0) {
291     ierr = VecAXPY(&ctx->vshift,a,y);CHKERRQ(ierr);
292   }
293 
294   if (ctx->sp) {ierr = MatNullSpaceRemove(ctx->sp,y,PETSC_NULL);CHKERRQ(ierr);}
295 
296   ierr = PetscLogEventEnd(MAT_MultMatrixFree,a,y,0,0);CHKERRQ(ierr);
297   PetscFunctionReturn(0);
298 }
299 
300 #undef __FUNCT__
301 #define __FUNCT__ "MatGetDiagonal_MFFD"
302 /*
303   MatGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix
304 
305         y ~= (F(u + ha) - F(u))/h,
306   where F = nonlinear function, as set by SNESSetFunction()
307         u = current iterate
308         h = difference interval
309 */
310 int MatGetDiagonal_MFFD(Mat mat,Vec a)
311 {
312   MatSNESMFCtx ctx = (MatSNESMFCtx)mat->data;
313   PetscScalar  h,*aa,*ww,v;
314   PetscReal    epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON;
315   Vec          w,U;
316   int          i,ierr,rstart,rend;
317 
318   PetscFunctionBegin;
319   if (!ctx->funci) {
320     SETERRQ(1,"Requirers calling MatSNESMFSetFunctioni() first");
321   }
322 
323   w    = ctx->w;
324   U    = ctx->current_u;
325   ierr = (*ctx->func)(0,U,a,ctx->funcctx);CHKERRQ(ierr);
326   ierr = (*ctx->funcisetbase)(U,ctx->funcctx);CHKERRQ(ierr);
327   ierr = VecCopy(U,w);CHKERRQ(ierr);
328 
329   ierr = VecGetOwnershipRange(a,&rstart,&rend);CHKERRQ(ierr);
330   ierr = VecGetArray(a,&aa);CHKERRQ(ierr);
331   for (i=rstart; i<rend; i++) {
332     ierr = VecGetArray(w,&ww);CHKERRQ(ierr);
333     h  = ww[i-rstart];
334     if (h == 0.0) h = 1.0;
335 #if !defined(PETSC_USE_COMPLEX)
336     if (h < umin && h >= 0.0)      h = umin;
337     else if (h < 0.0 && h > -umin) h = -umin;
338 #else
339     if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0)     h = umin;
340     else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin;
341 #endif
342     h     *= epsilon;
343 
344     ww[i-rstart] += h;
345     ierr = VecRestoreArray(w,&ww);CHKERRQ(ierr);
346     ierr          = (*ctx->funci)(i,w,&v,ctx->funcctx);CHKERRQ(ierr);
347     aa[i-rstart]  = (v - aa[i-rstart])/h;
348 
349     /* possibly shift and scale result */
350     aa[i - rstart] = ctx->vshift + ctx->vscale*aa[i-rstart];
351 
352     ierr = VecGetArray(w,&ww);CHKERRQ(ierr);
353     ww[i-rstart] -= h;
354     ierr = VecRestoreArray(w,&ww);CHKERRQ(ierr);
355   }
356   ierr = VecRestoreArray(a,&aa);CHKERRQ(ierr);
357   PetscFunctionReturn(0);
358 }
359 
360 #undef __FUNCT__
361 #define __FUNCT__ "MatShift_MFFD"
362 int MatShift_MFFD(const PetscScalar *a,Mat Y)
363 {
364   MatSNESMFCtx shell = (MatSNESMFCtx)Y->data;
365   PetscFunctionBegin;
366   shell->vshift += *a;
367   PetscFunctionReturn(0);
368 }
369 
370 #undef __FUNCT__
371 #define __FUNCT__ "MatScale_MFFD"
372 int MatScale_MFFD(const PetscScalar *a,Mat Y)
373 {
374   MatSNESMFCtx shell = (MatSNESMFCtx)Y->data;
375   PetscFunctionBegin;
376   shell->vscale *= *a;
377   PetscFunctionReturn(0);
378 }
379 
380 
381 #undef __FUNCT__
382 #define __FUNCT__ "MatCreateSNESMF"
383 /*@C
384    MatCreateSNESMF - Creates a matrix-free matrix context for use with
385    a SNES solver.  This matrix can be used as the Jacobian argument for
386    the routine SNESSetJacobian().
387 
388    Collective on SNES and Vec
389 
390    Input Parameters:
391 +  snes - the SNES context
392 -  x - vector where SNES solution is to be stored.
393 
394    Output Parameter:
395 .  J - the matrix-free matrix
396 
397    Level: advanced
398 
399    Notes:
400    The matrix-free matrix context merely contains the function pointers
401    and work space for performing finite difference approximations of
402    Jacobian-vector products, F'(u)*a,
403 
404    The default code uses the following approach to compute h
405 
406 .vb
407      F'(u)*a = [F(u+h*a) - F(u)]/h where
408      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
409        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
410  where
411      error_rel = square root of relative error in function evaluation
412      umin = minimum iterate parameter
413 .ve
414 
415    The user can set the error_rel via MatSNESMFSetFunctionError() and
416    umin via MatSNESMFDefaultSetUmin(); see the nonlinear solvers chapter
417    of the users manual for details.
418 
419    The user should call MatDestroy() when finished with the matrix-free
420    matrix context.
421 
422    Options Database Keys:
423 +  -snes_mf_err <error_rel> - Sets error_rel
424 .  -snes_mf_unim <umin> - Sets umin (for default PETSc routine that computes h only)
425 -  -snes_mf_ksp_monitor - KSP monitor routine that prints differencing h
426 
427 .keywords: SNES, default, matrix-free, create, matrix
428 
429 .seealso: MatDestroy(), MatSNESMFSetFunctionError(), MatSNESMFDefaultSetUmin()
430           MatSNESMFSetHHistory(), MatSNESMFResetHHistory(), MatCreateMF(),
431           MatSNESMFGetH(),MatSNESMFKSPMonitor(), MatSNESMFRegisterDynamic), MatSNESMFComputeJacobian()
432 
433 @*/
434 int MatCreateSNESMF(SNES snes,Vec x,Mat *J)
435 {
436   MatSNESMFCtx mfctx;
437   int          ierr;
438 
439   PetscFunctionBegin;
440   ierr = MatCreateMF(x,J);CHKERRQ(ierr);
441 
442   mfctx          = (MatSNESMFCtx)(*J)->data;
443   mfctx->snes    = snes;
444   mfctx->usesnes = PETSC_TRUE;
445   PetscLogObjectParent(snes,*J);
446   PetscFunctionReturn(0);
447 }
448 
449 EXTERN_C_BEGIN
450 #undef __FUNCT__
451 #define __FUNCT__ "MatSNESMFSetBase_FD"
452 int MatSNESMFSetBase_FD(Mat J,Vec U)
453 {
454   int          ierr;
455   MatSNESMFCtx ctx = (MatSNESMFCtx)J->data;
456 
457   PetscFunctionBegin;
458   ierr = MatSNESMFResetHHistory(J);CHKERRQ(ierr);
459   ctx->current_u = U;
460   ctx->usesnes   = PETSC_FALSE;
461   if (ctx->w == PETSC_NULL) {
462     ierr = VecDuplicate(ctx->current_u, &ctx->w);CHKERRQ(ierr);
463   }
464   PetscFunctionReturn(0);
465 }
466 EXTERN_C_END
467 
468 EXTERN_C_BEGIN
469 #undef __FUNCT__
470 #define __FUNCT__ "MatSNESMFSetCheckh_FD"
471 int MatSNESMFSetCheckh_FD(Mat J,int (*fun)(Vec,Vec,PetscScalar*,void*),void*ectx)
472 {
473   MatSNESMFCtx ctx = (MatSNESMFCtx)J->data;
474 
475   PetscFunctionBegin;
476   ctx->checkh    = fun;
477   ctx->checkhctx = ectx;
478   PetscFunctionReturn(0);
479 }
480 EXTERN_C_END
481 
482 #undef __FUNCT__
483 #define __FUNCT__ "MatSNESMFSetFromOptions"
484 /*@
485    MatSNESMFSetFromOptions - Sets the MatSNESMF options from the command line
486    parameter.
487 
488    Collective on Mat
489 
490    Input Parameters:
491 .  mat - the matrix obtained with MatCreateSNESMF()
492 
493    Options Database Keys:
494 +  -snes_mf_type - <default,wp>
495 -  -snes_mf_err - square root of estimated relative error in function evaluation
496 -  -snes_mf_period - how often h is recomputed, defaults to 1, everytime
497 
498    Level: advanced
499 
500 .keywords: SNES, matrix-free, parameters
501 
502 .seealso: MatCreateSNESMF(),MatSNESMFSetHHistory(),
503           MatSNESMFResetHHistory(), MatSNESMFKSPMonitor()
504 @*/
505 int MatSNESMFSetFromOptions(Mat mat)
506 {
507   MatSNESMFCtx mfctx = (MatSNESMFCtx)mat->data;
508   int          ierr;
509   PetscTruth   flg;
510   char         ftype[256];
511 
512   PetscFunctionBegin;
513   if (!MatSNESMFRegisterAllCalled) {ierr = MatSNESMFRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
514 
515   ierr = PetscOptionsBegin(mfctx->comm,mfctx->prefix,"Set matrix free computation parameters","MatSNESMF");CHKERRQ(ierr);
516   ierr = PetscOptionsList("-snes_mf_type","Matrix free type","MatSNESMFSetType",MatSNESMPetscFList,mfctx->type_name,ftype,256,&flg);CHKERRQ(ierr);
517   if (flg) {
518     ierr = MatSNESMFSetType(mat,ftype);CHKERRQ(ierr);
519   }
520 
521   ierr = PetscOptionsReal("-snes_mf_err","set sqrt relative error in function","MatSNESMFSetFunctionError",mfctx->error_rel,&mfctx->error_rel,0);CHKERRQ(ierr);
522   ierr = PetscOptionsInt("-snes_mf_period","how often h is recomputed","MatSNESMFSetPeriod",mfctx->recomputeperiod,&mfctx->recomputeperiod,0);CHKERRQ(ierr);
523   if (mfctx->snes) {
524     ierr = PetscOptionsName("-snes_mf_ksp_monitor","Monitor matrix-free parameters","MatSNESMFKSPMonitor",&flg);CHKERRQ(ierr);
525     if (flg) {
526       SLES sles;
527       KSP  ksp;
528       ierr = SNESGetSLES(mfctx->snes,&sles);CHKERRQ(ierr);
529       ierr = SLESGetKSP(sles,&ksp);CHKERRQ(ierr);
530       ierr = KSPSetMonitor(ksp,MatSNESMFKSPMonitor,PETSC_NULL,0);CHKERRQ(ierr);
531     }
532   }
533   ierr = PetscOptionsName("-snes_mf_check_positivity","Insure that U + h*a is nonnegative","MatSNESMFSetCheckh",&flg);CHKERRQ(ierr);
534   if (flg) {
535     ierr = MatSNESMFSetCheckh(mat,MatSNESMFCheckPositivity,0);CHKERRQ(ierr);
536   }
537   if (mfctx->ops->setfromoptions) {
538     ierr = (*mfctx->ops->setfromoptions)(mfctx);CHKERRQ(ierr);
539   }
540   ierr = PetscOptionsEnd();CHKERRQ(ierr);
541   PetscFunctionReturn(0);
542 }
543 
544 #undef __FUNCT__
545 #define __FUNCT__ "MatCreate_MFFD"
546 EXTERN_C_BEGIN
547 int MatCreate_MFFD(Mat A)
548 {
549   MatSNESMFCtx mfctx;
550   int          ierr;
551 
552   PetscFunctionBegin;
553 #ifndef PETSC_USE_DYNAMIC_LIBRARIES
554   ierr = SNESInitializePackage(PETSC_NULL);                                                               CHKERRQ(ierr);
555 #endif
556 
557   PetscHeaderCreate(mfctx,_p_MatSNESMFCtx,struct _MFOps,MATSNESMFCTX_COOKIE,0,"SNESMF",A->comm,MatDestroy_MFFD,MatView_MFFD);
558   PetscLogObjectCreate(mfctx);
559   mfctx->sp              = 0;
560   mfctx->snes            = 0;
561   mfctx->error_rel       = PETSC_SQRT_MACHINE_EPSILON;
562   mfctx->recomputeperiod = 1;
563   mfctx->count           = 0;
564   mfctx->currenth        = 0.0;
565   mfctx->historyh        = PETSC_NULL;
566   mfctx->ncurrenth       = 0;
567   mfctx->maxcurrenth     = 0;
568   mfctx->type_name       = 0;
569   mfctx->usesnes         = PETSC_FALSE;
570 
571   mfctx->vshift          = 0.0;
572   mfctx->vscale          = 1.0;
573 
574   /*
575      Create the empty data structure to contain compute-h routines.
576      These will be filled in below from the command line options or
577      a later call with MatSNESMFSetType() or if that is not called
578      then it will default in the first use of MatMult_MFFD()
579   */
580   mfctx->ops->compute        = 0;
581   mfctx->ops->destroy        = 0;
582   mfctx->ops->view           = 0;
583   mfctx->ops->setfromoptions = 0;
584   mfctx->hctx                = 0;
585 
586   mfctx->func                = 0;
587   mfctx->funcctx             = 0;
588   mfctx->funcvec             = 0;
589   mfctx->w                   = PETSC_NULL;
590 
591   A->data                = mfctx;
592 
593   A->ops->mult           = MatMult_MFFD;
594   A->ops->destroy        = MatDestroy_MFFD;
595   A->ops->view           = MatView_MFFD;
596   A->ops->assemblyend    = MatAssemblyEnd_MFFD;
597   A->ops->getdiagonal    = MatGetDiagonal_MFFD;
598   A->ops->scale          = MatScale_MFFD;
599   A->ops->shift          = MatShift_MFFD;
600   A->ops->setfromoptions = MatSNESMFSetFromOptions;
601 
602   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSNESMFSetBase_C","MatSNESMFSetBase_FD",MatSNESMFSetBase_FD);CHKERRQ(ierr);
603   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSNESMFSetFunctioniBase_C","MatSNESMFSetFunctioniBase_FD",MatSNESMFSetFunctioniBase_FD);CHKERRQ(ierr);
604   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSNESMFSetFunctioni_C","MatSNESMFSetFunctioni_FD",MatSNESMFSetFunctioni_FD);CHKERRQ(ierr);
605   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSNESMFSetCheckh_C","MatSNESMFSetCheckh_FD",MatSNESMFSetCheckh_FD);CHKERRQ(ierr);
606   mfctx->mat = A;
607 
608   PetscFunctionReturn(0);
609 }
610 
611 EXTERN_C_END
612 
613 #undef __FUNCT__
614 #define __FUNCT__ "MatCreateMF"
615 /*@C
616    MatCreateMF - Creates a matrix-free matrix. See also MatCreateSNESMF()
617 
618    Collective on Vec
619 
620    Input Parameters:
621 .  x - vector that defines layout of the vectors and matrices
622 
623    Output Parameter:
624 .  J - the matrix-free matrix
625 
626    Level: advanced
627 
628    Notes:
629    The matrix-free matrix context merely contains the function pointers
630    and work space for performing finite difference approximations of
631    Jacobian-vector products, F'(u)*a,
632 
633    The default code uses the following approach to compute h
634 
635 .vb
636      F'(u)*a = [F(u+h*a) - F(u)]/h where
637      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
638        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
639  where
640      error_rel = square root of relative error in function evaluation
641      umin = minimum iterate parameter
642 .ve
643 
644    The user can set the error_rel via MatSNESMFSetFunctionError() and
645    umin via MatSNESMFDefaultSetUmin(); see the nonlinear solvers chapter
646    of the users manual for details.
647 
648    The user should call MatDestroy() when finished with the matrix-free
649    matrix context.
650 
651    Options Database Keys:
652 +  -snes_mf_err <error_rel> - Sets error_rel
653 .  -snes_mf_unim <umin> - Sets umin (for default PETSc routine that computes h only)
654 .  -snes_mf_ksp_monitor - KSP monitor routine that prints differencing h
655 -  -snes_mf_check_positivity
656 
657 .keywords: default, matrix-free, create, matrix
658 
659 .seealso: MatDestroy(), MatSNESMFSetFunctionError(), MatSNESMFDefaultSetUmin()
660           MatSNESMFSetHHistory(), MatSNESMFResetHHistory(), MatCreateSNESMF(),
661           MatSNESMFGetH(),MatSNESMFKSPMonitor(), MatSNESMFRegisterDynamic),, MatSNESMFComputeJacobian()
662 
663 @*/
664 int MatCreateMF(Vec x,Mat *J)
665 {
666   MPI_Comm     comm;
667   int          n,nloc,ierr;
668 
669   PetscFunctionBegin;
670   ierr = PetscObjectGetComm((PetscObject)x,&comm);CHKERRQ(ierr);
671   ierr = VecGetSize(x,&n);CHKERRQ(ierr);
672   ierr = VecGetLocalSize(x,&nloc);CHKERRQ(ierr);
673   ierr = MatCreate(comm,nloc,nloc,n,n,J);CHKERRQ(ierr);
674   ierr = MatRegister(MATMFFD,0,"MatCreate_MFFD",MatCreate_MFFD);CHKERRQ(ierr);
675   ierr = MatSetType(*J,MATMFFD);CHKERRQ(ierr);
676   PetscFunctionReturn(0);
677 }
678 
679 
680 #undef __FUNCT__
681 #define __FUNCT__ "MatSNESMFGetH"
682 /*@
683    MatSNESMFGetH - Gets the last value that was used as the differencing
684    parameter.
685 
686    Not Collective
687 
688    Input Parameters:
689 .  mat - the matrix obtained with MatCreateSNESMF()
690 
691    Output Paramter:
692 .  h - the differencing step size
693 
694    Level: advanced
695 
696 .keywords: SNES, matrix-free, parameters
697 
698 .seealso: MatCreateSNESMF(),MatSNESMFSetHHistory(),
699           MatSNESMFResetHHistory(),MatSNESMFKSPMonitor()
700 @*/
701 int MatSNESMFGetH(Mat mat,PetscScalar *h)
702 {
703   MatSNESMFCtx ctx = (MatSNESMFCtx)mat->data;
704 
705   PetscFunctionBegin;
706   *h = ctx->currenth;
707   PetscFunctionReturn(0);
708 }
709 
710 #undef __FUNCT__
711 #define __FUNCT__ "MatSNESMFKSPMonitor"
712 /*
713    MatSNESMFKSPMonitor - A KSP monitor for use with the default PETSc
714    SNES matrix free routines. Prints the differencing parameter used at
715    each step.
716 */
717 int MatSNESMFKSPMonitor(KSP ksp,int n,PetscReal rnorm,void *dummy)
718 {
719   PC             pc;
720   MatSNESMFCtx   ctx;
721   int            ierr;
722   Mat            mat;
723   MPI_Comm       comm;
724   PetscTruth     nonzeroinitialguess;
725 
726   PetscFunctionBegin;
727   ierr = PetscObjectGetComm((PetscObject)ksp,&comm);CHKERRQ(ierr);
728   ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
729   ierr = KSPGetInitialGuessNonzero(ksp,&nonzeroinitialguess);CHKERRQ(ierr);
730   ierr = PCGetOperators(pc,&mat,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
731   ctx  = (MatSNESMFCtx)mat->data;
732 
733   if (n > 0 || nonzeroinitialguess) {
734 #if defined(PETSC_USE_COMPLEX)
735     ierr = PetscPrintf(comm,"%d KSP Residual norm %14.12e h %g + %g i\n",n,rnorm,
736                 PetscRealPart(ctx->currenth),PetscImaginaryPart(ctx->currenth));CHKERRQ(ierr);
737 #else
738     ierr = PetscPrintf(comm,"%d KSP Residual norm %14.12e h %g \n",n,rnorm,ctx->currenth);CHKERRQ(ierr);
739 #endif
740   } else {
741     ierr = PetscPrintf(comm,"%d KSP Residual norm %14.12e\n",n,rnorm);CHKERRQ(ierr);
742   }
743   PetscFunctionReturn(0);
744 }
745 
746 #undef __FUNCT__
747 #define __FUNCT__ "MatSNESMFSetFunction"
748 /*@C
749    MatSNESMFSetFunction - Sets the function used in applying the matrix free.
750 
751    Collective on Mat
752 
753    Input Parameters:
754 +  mat - the matrix free matrix created via MatCreateSNESMF()
755 .  v   - workspace vector
756 .  func - the function to use
757 -  funcctx - optional function context passed to function
758 
759    Level: advanced
760 
761    Notes:
762     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
763     matrix inside your compute Jacobian routine
764 
765     If this is not set then it will use the function set with SNESSetFunction()
766 
767 .keywords: SNES, matrix-free, function
768 
769 .seealso: MatCreateSNESMF(),MatSNESMFGetH(),
770           MatSNESMFSetHHistory(), MatSNESMFResetHHistory(),
771           MatSNESMFKSPMonitor(), SNESetFunction()
772 @*/
773 int MatSNESMFSetFunction(Mat mat,Vec v,int (*func)(SNES,Vec,Vec,void *),void *funcctx)
774 {
775   MatSNESMFCtx ctx = (MatSNESMFCtx)mat->data;
776 
777   PetscFunctionBegin;
778   ctx->func    = func;
779   ctx->funcctx = funcctx;
780   ctx->funcvec = v;
781   PetscFunctionReturn(0);
782 }
783 
784 #undef __FUNCT__
785 #define __FUNCT__ "MatSNESMFSetFunctioni"
786 /*@C
787    MatSNESMFSetFunctioni - 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(),MatSNESMFGetH(),
805           MatSNESMFSetHHistory(), MatSNESMFResetHHistory(),
806           MatSNESMFKSPMonitor(), SNESetFunction()
807 @*/
808 int MatSNESMFSetFunctioni(Mat mat,int (*funci)(int,Vec,PetscScalar*,void *))
809 {
810   int  ierr,(*f)(Mat,int (*)(int,Vec,PetscScalar*,void *));
811 
812   PetscFunctionBegin;
813   PetscValidHeaderSpecific(mat,MAT_COOKIE);
814   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSNESMFSetFunctioni_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__ "MatSNESMFSetFunctioniBase"
824 /*@C
825    MatSNESMFSetFunctioniBase - 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(),MatSNESMFGetH(),
843           MatSNESMFSetHHistory(), MatSNESMFResetHHistory(),
844           MatSNESMFKSPMonitor(), SNESetFunction()
845 @*/
846 int MatSNESMFSetFunctioniBase(Mat mat,int (*func)(Vec,void *))
847 {
848   int  ierr,(*f)(Mat,int (*)(Vec,void *));
849 
850   PetscFunctionBegin;
851   PetscValidHeaderSpecific(mat,MAT_COOKIE);
852   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSNESMFSetFunctioniBase_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__ "MatSNESMFSetPeriod"
862 /*@
863    MatSNESMFSetPeriod - 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 +  -snes_mf_period <period>
873 
874    Level: advanced
875 
876 
877 .keywords: SNES, matrix-free, parameters
878 
879 .seealso: MatCreateSNESMF(),MatSNESMFGetH(),
880           MatSNESMFSetHHistory(), MatSNESMFResetHHistory(),
881           MatSNESMFKSPMonitor()
882 @*/
883 int MatSNESMFSetPeriod(Mat mat,int period)
884 {
885   MatSNESMFCtx ctx = (MatSNESMFCtx)mat->data;
886 
887   PetscFunctionBegin;
888   ctx->recomputeperiod = period;
889   PetscFunctionReturn(0);
890 }
891 
892 #undef __FUNCT__
893 #define __FUNCT__ "MatSNESMFSetFunctionError"
894 /*@
895    MatSNESMFSetFunctionError - 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 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 +  -snes_mf_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(),MatSNESMFGetH(),
921           MatSNESMFSetHHistory(), MatSNESMFResetHHistory(),
922           MatSNESMFKSPMonitor()
923 @*/
924 int MatSNESMFSetFunctionError(Mat mat,PetscReal error)
925 {
926   MatSNESMFCtx ctx = (MatSNESMFCtx)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__ "MatSNESMFAddNullSpace"
935 /*@
936    MatSNESMFAddNullSpace - 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(), MatSNESMFGetH(), MatCreateSNESMF(),
952           MatSNESMFSetHHistory(), MatSNESMFResetHHistory(),
953           MatSNESMFKSPMonitor(), MatSNESMFErrorRel()
954 @*/
955 int MatSNESMFAddNullSpace(Mat J,MatNullSpace nullsp)
956 {
957   int          ierr;
958   MatSNESMFCtx ctx = (MatSNESMFCtx)J->data;
959   MPI_Comm     comm;
960 
961   PetscFunctionBegin;
962   ierr = PetscObjectGetComm((PetscObject)J,&comm);CHKERRQ(ierr);
963 
964   ctx->sp = nullsp;
965   ierr    = PetscObjectReference((PetscObject)nullsp);CHKERRQ(ierr);
966   PetscFunctionReturn(0);
967 }
968 
969 #undef __FUNCT__
970 #define __FUNCT__ "MatSNESMFSetHHistory"
971 /*@
972    MatSNESMFSetHHistory - Sets an array to collect a history of the
973    differencing values (h) computed for the matrix-free product.
974 
975    Collective on Mat
976 
977    Input Parameters:
978 +  J - the matrix-free matrix context
979 .  histroy - space to hold the history
980 -  nhistory - number of entries in history, if more entries are generated than
981               nhistory, then the later ones are discarded
982 
983    Level: advanced
984 
985    Notes:
986    Use MatSNESMFResetHHistory() to reset the history counter and collect
987    a new batch of differencing parameters, h.
988 
989 .keywords: SNES, matrix-free, h history, differencing history
990 
991 .seealso: MatSNESMFGetH(), MatCreateSNESMF(),
992           MatSNESMFResetHHistory(),
993           MatSNESMFKSPMonitor(), MatSNESMFSetFunctionError()
994 
995 @*/
996 int MatSNESMFSetHHistory(Mat J,PetscScalar history[],int nhistory)
997 {
998   MatSNESMFCtx ctx = (MatSNESMFCtx)J->data;
999 
1000   PetscFunctionBegin;
1001   ctx->historyh    = history;
1002   ctx->maxcurrenth = nhistory;
1003   ctx->currenth    = 0;
1004   PetscFunctionReturn(0);
1005 }
1006 
1007 #undef __FUNCT__
1008 #define __FUNCT__ "MatSNESMFResetHHistory"
1009 /*@
1010    MatSNESMFResetHHistory - Resets the counter to zero to begin
1011    collecting a new set of differencing histories.
1012 
1013    Collective on Mat
1014 
1015    Input Parameters:
1016 .  J - the matrix-free matrix context
1017 
1018    Level: advanced
1019 
1020    Notes:
1021    Use MatSNESMFSetHHistory() to create the original history counter.
1022 
1023 .keywords: SNES, matrix-free, h history, differencing history
1024 
1025 .seealso: MatSNESMFGetH(), MatCreateSNESMF(),
1026           MatSNESMFSetHHistory(),
1027           MatSNESMFKSPMonitor(), MatSNESMFSetFunctionError()
1028 
1029 @*/
1030 int MatSNESMFResetHHistory(Mat J)
1031 {
1032   MatSNESMFCtx ctx = (MatSNESMFCtx)J->data;
1033 
1034   PetscFunctionBegin;
1035   ctx->ncurrenth    = 0;
1036   PetscFunctionReturn(0);
1037 }
1038 
1039 #undef __FUNCT__
1040 #define __FUNCT__ "MatSNESMFComputeJacobian"
1041 int MatSNESMFComputeJacobian(SNES snes,Vec x,Mat *jac,Mat *B,MatStructure *flag,void *dummy)
1042 {
1043   int ierr;
1044   PetscFunctionBegin;
1045   ierr = MatAssemblyBegin(*jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1046   ierr = MatAssemblyEnd(*jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1047   PetscFunctionReturn(0);
1048 }
1049 
1050 #undef __FUNCT__
1051 #define __FUNCT__ "MatSNESMFSetBase"
1052 /*@
1053     MatSNESMFSetBase - Sets the vector U at which matrix vector products of the
1054         Jacobian are computed
1055 
1056     Collective on Mat
1057 
1058     Input Parameters:
1059 +   J - the MatSNESMF matrix
1060 -   U - the vector
1061 
1062     Notes: This is rarely used directly
1063 
1064     Level: advanced
1065 
1066 @*/
1067 int MatSNESMFSetBase(Mat J,Vec U)
1068 {
1069   int  ierr,(*f)(Mat,Vec);
1070 
1071   PetscFunctionBegin;
1072   PetscValidHeaderSpecific(J,MAT_COOKIE);
1073   PetscValidHeaderSpecific(U,VEC_COOKIE);
1074   ierr = PetscObjectQueryFunction((PetscObject)J,"MatSNESMFSetBase_C",(void (**)(void))&f);CHKERRQ(ierr);
1075   if (f) {
1076     ierr = (*f)(J,U);CHKERRQ(ierr);
1077   }
1078   PetscFunctionReturn(0);
1079 }
1080 
1081 #undef __FUNCT__
1082 #define __FUNCT__ "MatSNESMFSetCheckh"
1083 /*@C
1084     MatSNESMFSetCheckh - Sets a function that checks the computed h and adjusts
1085         it to satisfy some criteria
1086 
1087     Collective on Mat
1088 
1089     Input Parameters:
1090 +   J - the MatSNESMF matrix
1091 .   fun - the function that checks h
1092 -   ctx - any context needed by the function
1093 
1094     Options Database Keys:
1095 .   -snes_mf_check_positivity
1096 
1097     Level: advanced
1098 
1099     Notes: For example, MatSNESMFSetCheckPositivity() insures that all entries
1100        of U + h*a are non-negative
1101 
1102 .seealso:  MatSNESMFSetCheckPositivity()
1103 @*/
1104 int MatSNESMFSetCheckh(Mat J,int (*fun)(Vec,Vec,PetscScalar*,void*),void* ctx)
1105 {
1106   int  ierr,(*f)(Mat,int (*)(Vec,Vec,PetscScalar*,void*),void*);
1107 
1108   PetscFunctionBegin;
1109   PetscValidHeaderSpecific(J,MAT_COOKIE);
1110   ierr = PetscObjectQueryFunction((PetscObject)J,"MatSNESMFSetCheckh_C",(void (**)(void))&f);CHKERRQ(ierr);
1111   if (f) {
1112     ierr = (*f)(J,fun,ctx);CHKERRQ(ierr);
1113   }
1114   PetscFunctionReturn(0);
1115 }
1116 
1117 #undef __FUNCT__
1118 #define __FUNCT__ "MatSNESMFSetCheckPositivity"
1119 /*@
1120     MatSNESMFCheckPositivity - Checks that all entries in U + h*a are positive or
1121         zero, decreases h until this is satisfied.
1122 
1123     Collective on Vec
1124 
1125     Input Parameters:
1126 +   U - base vector that is added to
1127 .   a - vector that is added
1128 .   h - scaling factor on a
1129 -   dummy - context variable (unused)
1130 
1131     Options Database Keys:
1132 .   -snes_mf_check_positivity
1133 
1134     Level: advanced
1135 
1136     Notes: This is rarely used directly, rather it is passed as an argument to
1137            MatSNESMFSetCheckh()
1138 
1139 .seealso:  MatSNESMFSetCheckh()
1140 @*/
1141 int MatSNESMFCheckPositivity(Vec U,Vec a,PetscScalar *h,void *dummy)
1142 {
1143   PetscReal     val, minval;
1144   PetscScalar   *u_vec, *a_vec;
1145   int           ierr, i, size;
1146   MPI_Comm      comm;
1147 
1148   PetscFunctionBegin;
1149   ierr = PetscObjectGetComm((PetscObject)U,&comm);CHKERRQ(ierr);
1150   ierr = VecGetArray(U,&u_vec);CHKERRQ(ierr);
1151   ierr = VecGetArray(a,&a_vec);CHKERRQ(ierr);
1152   ierr = VecGetLocalSize(U,&size);CHKERRQ(ierr);
1153   minval = PetscAbsScalar(*h*1.01);
1154   for(i=0;i<size;i++) {
1155     if (PetscRealPart(u_vec[i] + *h*a_vec[i]) <= 0.0) {
1156       val = PetscAbsScalar(u_vec[i]/a_vec[i]);
1157       if (val < minval) minval = val;
1158     }
1159   }
1160   ierr = VecRestoreArray(U,&u_vec);CHKERRQ(ierr);
1161   ierr = VecRestoreArray(a,&a_vec);CHKERRQ(ierr);
1162   ierr = PetscGlobalMin(&minval,&val,comm);CHKERRQ(ierr);
1163   if (val <= PetscAbsScalar(*h)) {
1164     PetscLogInfo(U,"MatSNESMFCheckPositivity: Scaling back h from %g to %g\n",PetscRealPart(*h),.99*val);
1165     if (PetscRealPart(*h) > 0.0) *h =  0.99*val;
1166     else                         *h = -0.99*val;
1167   }
1168   PetscFunctionReturn(0);
1169 }
1170 
1171 
1172 
1173 
1174 
1175 
1176 
1177