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