Actual source code: snesmfj.c

  1: /*$Id: snesmfj.c,v 1.131 2001/09/05 18:45:40 bsmith Exp $*/

 3:  #include src/mat/matimpl.h
 4:  #include src/snes/mf/snesmfj.h

  6: PetscFList      MatSNESMPetscFList              = 0;
  7: PetscTruth MatSNESMFRegisterAllCalled = PETSC_FALSE;

 11: /*@C
 12:     MatSNESMFSetType - Sets the method that is used to compute the 
 13:     differencing parameter for finite differene matrix-free formulations. 

 15:     Input Parameters:
 16: +   mat - the "matrix-free" matrix created via MatCreateSNESMF(), or MatCreateMF()
 17:           or MatSetType(mat,MATMFFD);
 18: -   ftype - the type requested

 20:     Level: advanced

 22:     Notes:
 23:     For example, such routines can compute h for use in
 24:     Jacobian-vector products of the form

 26:                         F(x+ha) - F(x)
 27:           F'(u)a  ~=  ----------------
 28:                               h

 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: 

 42:   /* already set, so just return */
 43:   PetscTypeCompare((PetscObject)ctx,ftype,&match);
 44:   if (match) return(0);

 46:   /* destroy the old one if it exists */
 47:   if (ctx->ops->destroy) {
 48:     (*ctx->ops->destroy)(ctx);
 49:   }

 51:   /* Get the function pointers for the requrested method */
 52:   if (!MatSNESMFRegisterAllCalled) {MatSNESMFRegisterAll(PETSC_NULL);}

 54:    PetscFListFind(ctx->comm,MatSNESMPetscFList,ftype,(void (**)(void)) &r);

 56:   if (!r) SETERRQ(1,"Unknown MatSNESMF type given");

 58:   (*r)(ctx);

 60:   PetscObjectChangeTypeName((PetscObject)ctx,ftype);

 62:   return(0);
 63: }

 65: EXTERN_C_BEGIN
 68: int MatSNESMFSetFunctioniBase_FD(Mat mat,int (*func)(Vec,void *))
 69: {
 70:   MatSNESMFCtx ctx = (MatSNESMFCtx)mat->data;

 73:   ctx->funcisetbase = func;
 74:   return(0);
 75: }
 76: EXTERN_C_END

 78: EXTERN_C_BEGIN
 81: int MatSNESMFSetFunctioni_FD(Mat mat,int (*funci)(int,Vec,PetscScalar*,void *))
 82: {
 83:   MatSNESMFCtx ctx = (MatSNESMFCtx)mat->data;

 86:   ctx->funci = funci;
 87:   return(0);
 88: }
 89: EXTERN_C_END


 94: int MatSNESMFRegister(const char sname[],const char path[],const char name[],int (*function)(MatSNESMFCtx))
 95: {
 97:   char fullname[256];

100:   PetscFListConcat(path,name,fullname);
101:   PetscFListAdd(&MatSNESMPetscFList,sname,fullname,(void (*)(void))function);
102:   return(0);
103: }


108: /*@C
109:    MatSNESMFRegisterDestroy - Frees the list of MatSNESMF methods that were
110:    registered by MatSNESMFRegisterDynamic).

112:    Not Collective

114:    Level: developer

116: .keywords: MatSNESMF, register, destroy

118: .seealso: MatSNESMFRegisterDynamic), MatSNESMFRegisterAll()
119: @*/
120: int MatSNESMFRegisterDestroy(void)
121: {

125:   if (MatSNESMPetscFList) {
126:     PetscFListDestroy(&MatSNESMPetscFList);
127:     MatSNESMPetscFList = 0;
128:   }
129:   MatSNESMFRegisterAllCalled = PETSC_FALSE;
130:   return(0);
131: }

133: /* ----------------------------------------------------------------------------------------*/
136: int MatDestroy_MFFD(Mat mat)
137: {
138:   int          ierr;
139:   MatSNESMFCtx ctx = (MatSNESMFCtx)mat->data;

142:   if (ctx->w != PETSC_NULL) {
143:     VecDestroy(ctx->w);
144:   }
145:   if (ctx->ops->destroy) {(*ctx->ops->destroy)(ctx);}
146:   if (ctx->sp) {MatNullSpaceDestroy(ctx->sp);}
147:   PetscHeaderDestroy(ctx);
148:   return(0);
149: }

153: /*
154:    MatSNESMFView_MFFD - Views matrix-free parameters.

156: */
157: int MatView_MFFD(Mat J,PetscViewer viewer)
158: {
159:   int          ierr;
160:   MatSNESMFCtx ctx = (MatSNESMFCtx)J->data;
161:   PetscTruth   isascii;

164:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
165:   if (isascii) {
166:      PetscViewerASCIIPrintf(viewer,"  SNES matrix-free approximation:\n");
167:      PetscViewerASCIIPrintf(viewer,"    err=%g (relative error in function evaluation)\n",ctx->error_rel);
168:      if (!ctx->type_name) {
169:        PetscViewerASCIIPrintf(viewer,"    The compute h routine has not yet been set\n");
170:      } else {
171:        PetscViewerASCIIPrintf(viewer,"    Using %s compute h routine\n",ctx->type_name);
172:      }
173:      if (ctx->ops->view) {
174:        (*ctx->ops->view)(ctx,viewer);
175:      }
176:   } else {
177:     SETERRQ1(1,"Viewer type %s not supported for SNES matrix free matrix",((PetscObject)viewer)->type_name);
178:   }
179:   return(0);
180: }

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;

197:   MatSNESMFResetHHistory(J);
198:   if (j->usesnes) {
199:     SNESGetSolution(j->snes,&j->current_u);
200:     SNESGetFunction(j->snes,&j->current_f,PETSC_NULL,PETSC_NULL);
201:     if (j->w == PETSC_NULL) {
202:       VecDuplicate(j->current_u, &j->w);
203:     }
204:   }
205:   j->vshift = 0.0;
206:   j->vscale = 1.0;
207:   return(0);
208: }

212: /*
213:   MatMult_MFFD - Default matrix-free form for Jacobian-vector product, y = F'(u)*a:

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;

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:   PetscLogEventBegin(MAT_MultMatrixFree,a,y,0,0);

235:   snes = ctx->snes;
236:   w    = ctx->w;
237:   U    = ctx->current_u;

239:   /* 
240:       Compute differencing parameter 
241:   */
242:   if (!ctx->ops->compute) {
243:     MatSNESMFSetType(mat,MATSNESMF_WP);
244:     MatSNESMFSetFromOptions(mat);
245:   }
246:   (*ctx->ops->compute)(ctx,U,a,&h);

248:   if (ctx->checkh) {
249:     (*ctx->checkh)(U,a,&h,ctx->checkhctx);
250:   }

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++;

264:   /* w = u + ha */
265:   VecWAXPY(&h,a,U,w);

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:     (*eval_fct)(snes,w,y);
272:   } else {
273:     F = ctx->funcvec;
274:     /* compute func(U) as base for differencing */
275:     if (ctx->ncurrenth == 1) {
276:       (*ctx->func)(snes,U,F,ctx->funcctx);
277:     }
278:     (*ctx->func)(snes,w,y,ctx->funcctx);
279:   }

281:   VecAXPY(&mone,F,y);
282:   h    = 1.0/h;
283:   VecScale(&h,y);


286:   if (ctx->vshift != 0.0 && ctx->vscale != 1.0) {
287:     VecAXPBY(&ctx->vshift,&ctx->vscale,a,y);
288:   } else if (ctx->vscale != 1.0) {
289:     VecScale(&ctx->vscale,y);
290:   } else if (ctx->vshift != 0.0) {
291:     VecAXPY(&ctx->vshift,a,y);
292:   }

294:   if (ctx->sp) {MatNullSpaceRemove(ctx->sp,y,PETSC_NULL);}

296:   PetscLogEventEnd(MAT_MultMatrixFree,a,y,0,0);
297:   return(0);
298: }

302: /*
303:   MatGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix

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;

319:   if (!ctx->funci) {
320:     SETERRQ(1,"Requirers calling MatSNESMFSetFunctioni() first");
321:   }

323:   w    = ctx->w;
324:   U    = ctx->current_u;
325:   (*ctx->func)(0,U,a,ctx->funcctx);
326:   (*ctx->funcisetbase)(U,ctx->funcctx);
327:   VecCopy(U,w);

329:   VecGetOwnershipRange(a,&rstart,&rend);
330:   VecGetArray(a,&aa);
331:   for (i=rstart; i<rend; i++) {
332:     VecGetArray(w,&ww);
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:     VecRestoreArray(w,&ww);
346:     (*ctx->funci)(i,w,&v,ctx->funcctx);
347:     aa[i-rstart]  = (v - aa[i-rstart])/h;

349:     /* possibly shift and scale result */
350:     aa[i - rstart] = ctx->vshift + ctx->vscale*aa[i-rstart];

352:     VecGetArray(w,&ww);
353:     ww[i-rstart] -= h;
354:     VecRestoreArray(w,&ww);
355:   }
356:   VecRestoreArray(a,&aa);
357:   return(0);
358: }

362: int MatShift_MFFD(const PetscScalar *a,Mat Y)
363: {
364:   MatSNESMFCtx shell = (MatSNESMFCtx)Y->data;
366:   shell->vshift += *a;
367:   return(0);
368: }

372: int MatScale_MFFD(const PetscScalar *a,Mat Y)
373: {
374:   MatSNESMFCtx shell = (MatSNESMFCtx)Y->data;
376:   shell->vscale *= *a;
377:   return(0);
378: }


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().

388:    Collective on SNES and Vec

390:    Input Parameters:
391: +  snes - the SNES context
392: -  x - vector where SNES solution is to be stored.

394:    Output Parameter:
395: .  J - the matrix-free matrix

397:    Level: advanced

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, 

404:    The default code uses the following approach to compute h

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

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.

419:    The user should call MatDestroy() when finished with the matrix-free
420:    matrix context.

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

427: .keywords: SNES, default, matrix-free, create, matrix

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;

440:   MatCreateMF(x,J);

442:   mfctx          = (MatSNESMFCtx)(*J)->data;
443:   mfctx->snes    = snes;
444:   mfctx->usesnes = PETSC_TRUE;
445:   PetscLogObjectParent(snes,*J);
446:   return(0);
447: }

449: EXTERN_C_BEGIN
452: int MatSNESMFSetBase_FD(Mat J,Vec U)
453: {
454:   int          ierr;
455:   MatSNESMFCtx ctx = (MatSNESMFCtx)J->data;

458:   MatSNESMFResetHHistory(J);
459:   ctx->current_u = U;
460:   ctx->usesnes   = PETSC_FALSE;
461:   if (ctx->w == PETSC_NULL) {
462:     VecDuplicate(ctx->current_u, &ctx->w);
463:   }
464:   return(0);
465: }
466: EXTERN_C_END

468: EXTERN_C_BEGIN
471: int MatSNESMFSetCheckh_FD(Mat J,int (*fun)(Vec,Vec,PetscScalar*,void*),void*ectx)
472: {
473:   MatSNESMFCtx ctx = (MatSNESMFCtx)J->data;

476:   ctx->checkh    = fun;
477:   ctx->checkhctx = ectx;
478:   return(0);
479: }
480: EXTERN_C_END

484: /*@
485:    MatSNESMFSetFromOptions - Sets the MatSNESMF options from the command line
486:    parameter.

488:    Collective on Mat

490:    Input Parameters:
491: .  mat - the matrix obtained with MatCreateSNESMF()

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

498:    Level: advanced

500: .keywords: SNES, matrix-free, parameters

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];

513:   if (!MatSNESMFRegisterAllCalled) {MatSNESMFRegisterAll(PETSC_NULL);}
514: 
515:   PetscOptionsBegin(mfctx->comm,mfctx->prefix,"Set matrix free computation parameters","MatSNESMF");
516:   PetscOptionsList("-snes_mf_type","Matrix free type","MatSNESMFSetType",MatSNESMPetscFList,mfctx->type_name,ftype,256,&flg);
517:   if (flg) {
518:     MatSNESMFSetType(mat,ftype);
519:   }

521:   PetscOptionsReal("-snes_mf_err","set sqrt relative error in function","MatSNESMFSetFunctionError",mfctx->error_rel,&mfctx->error_rel,0);
522:   PetscOptionsInt("-snes_mf_period","how often h is recomputed","MatSNESMFSetPeriod",mfctx->recomputeperiod,&mfctx->recomputeperiod,0);
523:   if (mfctx->snes) {
524:     PetscOptionsName("-snes_mf_ksp_monitor","Monitor matrix-free parameters","MatSNESMFKSPMonitor",&flg);
525:     if (flg) {
526:       SLES sles;
527:       KSP  ksp;
528:       SNESGetSLES(mfctx->snes,&sles);
529:       SLESGetKSP(sles,&ksp);
530:       KSPSetMonitor(ksp,MatSNESMFKSPMonitor,PETSC_NULL,0);
531:     }
532:   }
533:   PetscOptionsName("-snes_mf_check_positivity","Insure that U + h*a is nonnegative","MatSNESMFSetCheckh",&flg);
534:   if (flg) {
535:     MatSNESMFSetCheckh(mat,MatSNESMFCheckPositivity,0);
536:   }
537:   if (mfctx->ops->setfromoptions) {
538:     (*mfctx->ops->setfromoptions)(mfctx);
539:   }
540:   PetscOptionsEnd();
541:   return(0);
542: }

544: /*MC
545:   MATMFFD - MATMFFD = "mffd" - A matrix free matrix type.

547:   Level: advanced

549: .seealso: MatCreateMF, MatCreateSNESMF
550: M*/

554: EXTERN_C_BEGIN
555: int MatCreate_MFFD(Mat A)
556: {
557:   MatSNESMFCtx mfctx;
558:   int          ierr;

561: #ifndef PETSC_USE_DYNAMIC_LIBRARIES
562:   SNESInitializePackage(PETSC_NULL);
563: #endif

565:   PetscHeaderCreate(mfctx,_p_MatSNESMFCtx,struct _MFOps,MATSNESMFCTX_COOKIE,0,"SNESMF",A->comm,MatDestroy_MFFD,MatView_MFFD);
566:   PetscLogObjectCreate(mfctx);
567:   mfctx->sp              = 0;
568:   mfctx->snes            = 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:   mfctx->type_name       = 0;
577:   mfctx->usesnes         = PETSC_FALSE;

579:   mfctx->vshift          = 0.0;
580:   mfctx->vscale          = 1.0;

582:   /* 
583:      Create the empty data structure to contain compute-h routines.
584:      These will be filled in below from the command line options or 
585:      a later call with MatSNESMFSetType() or if that is not called 
586:      then it will default in the first use of MatMult_MFFD()
587:   */
588:   mfctx->ops->compute        = 0;
589:   mfctx->ops->destroy        = 0;
590:   mfctx->ops->view           = 0;
591:   mfctx->ops->setfromoptions = 0;
592:   mfctx->hctx                = 0;

594:   mfctx->func                = 0;
595:   mfctx->funcctx             = 0;
596:   mfctx->funcvec             = 0;
597:   mfctx->w                   = PETSC_NULL;

599:   A->data                = mfctx;

601:   A->ops->mult           = MatMult_MFFD;
602:   A->ops->destroy        = MatDestroy_MFFD;
603:   A->ops->view           = MatView_MFFD;
604:   A->ops->assemblyend    = MatAssemblyEnd_MFFD;
605:   A->ops->getdiagonal    = MatGetDiagonal_MFFD;
606:   A->ops->scale          = MatScale_MFFD;
607:   A->ops->shift          = MatShift_MFFD;
608:   A->ops->setfromoptions = MatSNESMFSetFromOptions;

610:   PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSNESMFSetBase_C","MatSNESMFSetBase_FD",MatSNESMFSetBase_FD);
611:   PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSNESMFSetFunctioniBase_C","MatSNESMFSetFunctioniBase_FD",MatSNESMFSetFunctioniBase_FD);
612:   PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSNESMFSetFunctioni_C","MatSNESMFSetFunctioni_FD",MatSNESMFSetFunctioni_FD);
613:   PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSNESMFSetCheckh_C","MatSNESMFSetCheckh_FD",MatSNESMFSetCheckh_FD);
614:   mfctx->mat = A;

616:   return(0);
617: }

619: EXTERN_C_END

623: /*@C
624:    MatCreateMF - Creates a matrix-free matrix. See also MatCreateSNESMF() 

626:    Collective on Vec

628:    Input Parameters:
629: .  x - vector that defines layout of the vectors and matrices

631:    Output Parameter:
632: .  J - the matrix-free matrix

634:    Level: advanced

636:    Notes:
637:    The matrix-free matrix context merely contains the function pointers
638:    and work space for performing finite difference approximations of
639:    Jacobian-vector products, F'(u)*a, 

641:    The default code uses the following approach to compute h

643: .vb
644:      F'(u)*a = [F(u+h*a) - F(u)]/h where
645:      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
646:        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
647:  where
648:      error_rel = square root of relative error in function evaluation
649:      umin = minimum iterate parameter
650: .ve

652:    The user can set the error_rel via MatSNESMFSetFunctionError() and 
653:    umin via MatSNESMFDefaultSetUmin(); see the nonlinear solvers chapter
654:    of the users manual for details.

656:    The user should call MatDestroy() when finished with the matrix-free
657:    matrix context.

659:    Options Database Keys:
660: +  -snes_mf_err <error_rel> - Sets error_rel
661: .  -snes_mf_unim <umin> - Sets umin (for default PETSc routine that computes h only)
662: .  -snes_mf_ksp_monitor - KSP monitor routine that prints differencing h
663: -  -snes_mf_check_positivity

665: .keywords: default, matrix-free, create, matrix

667: .seealso: MatDestroy(), MatSNESMFSetFunctionError(), MatSNESMFDefaultSetUmin()
668:           MatSNESMFSetHHistory(), MatSNESMFResetHHistory(), MatCreateSNESMF(),
669:           MatSNESMFGetH(),MatSNESMFKSPMonitor(), MatSNESMFRegisterDynamic),, MatSNESMFComputeJacobian()
670:  
671: @*/
672: int MatCreateMF(Vec x,Mat *J)
673: {
674:   MPI_Comm     comm;
675:   int          n,nloc,ierr;

678:   PetscObjectGetComm((PetscObject)x,&comm);
679:   VecGetSize(x,&n);
680:   VecGetLocalSize(x,&nloc);
681:   MatCreate(comm,nloc,nloc,n,n,J);
682:   MatRegister(MATMFFD,0,"MatCreate_MFFD",MatCreate_MFFD);
683:   MatSetType(*J,MATMFFD);
684:   return(0);
685: }


690: /*@
691:    MatSNESMFGetH - Gets the last value that was used as the differencing 
692:    parameter.

694:    Not Collective

696:    Input Parameters:
697: .  mat - the matrix obtained with MatCreateSNESMF()

699:    Output Paramter:
700: .  h - the differencing step size

702:    Level: advanced

704: .keywords: SNES, matrix-free, parameters

706: .seealso: MatCreateSNESMF(),MatSNESMFSetHHistory(), 
707:           MatSNESMFResetHHistory(),MatSNESMFKSPMonitor()
708: @*/
709: int MatSNESMFGetH(Mat mat,PetscScalar *h)
710: {
711:   MatSNESMFCtx ctx = (MatSNESMFCtx)mat->data;

714:   *h = ctx->currenth;
715:   return(0);
716: }

720: /*
721:    MatSNESMFKSPMonitor - A KSP monitor for use with the default PETSc
722:    SNES matrix free routines. Prints the differencing parameter used at 
723:    each step.
724: */
725: int MatSNESMFKSPMonitor(KSP ksp,int n,PetscReal rnorm,void *dummy)
726: {
727:   PC             pc;
728:   MatSNESMFCtx   ctx;
729:   int            ierr;
730:   Mat            mat;
731:   MPI_Comm       comm;
732:   PetscTruth     nonzeroinitialguess;

735:   PetscObjectGetComm((PetscObject)ksp,&comm);
736:   KSPGetPC(ksp,&pc);
737:   KSPGetInitialGuessNonzero(ksp,&nonzeroinitialguess);
738:   PCGetOperators(pc,&mat,PETSC_NULL,PETSC_NULL);
739:   ctx  = (MatSNESMFCtx)mat->data;

741:   if (n > 0 || nonzeroinitialguess) {
742: #if defined(PETSC_USE_COMPLEX)
743:     PetscPrintf(comm,"%d KSP Residual norm %14.12e h %g + %g i\n",n,rnorm,
744:                 PetscRealPart(ctx->currenth),PetscImaginaryPart(ctx->currenth));
745: #else
746:     PetscPrintf(comm,"%d KSP Residual norm %14.12e h %g \n",n,rnorm,ctx->currenth);
747: #endif
748:   } else {
749:     PetscPrintf(comm,"%d KSP Residual norm %14.12e\n",n,rnorm);
750:   }
751:   return(0);
752: }

756: /*@C
757:    MatSNESMFSetFunction - Sets the function used in applying the matrix free.

759:    Collective on Mat

761:    Input Parameters:
762: +  mat - the matrix free matrix created via MatCreateSNESMF()
763: .  v   - workspace vector
764: .  func - the function to use
765: -  funcctx - optional function context passed to function

767:    Level: advanced

769:    Notes:
770:     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
771:     matrix inside your compute Jacobian routine

773:     If this is not set then it will use the function set with SNESSetFunction()

775: .keywords: SNES, matrix-free, function

777: .seealso: MatCreateSNESMF(),MatSNESMFGetH(),
778:           MatSNESMFSetHHistory(), MatSNESMFResetHHistory(),
779:           MatSNESMFKSPMonitor(), SNESetFunction()
780: @*/
781: int MatSNESMFSetFunction(Mat mat,Vec v,int (*func)(SNES,Vec,Vec,void *),void *funcctx)
782: {
783:   MatSNESMFCtx ctx = (MatSNESMFCtx)mat->data;

786:   ctx->func    = func;
787:   ctx->funcctx = funcctx;
788:   ctx->funcvec = v;
789:   return(0);
790: }

794: /*@C
795:    MatSNESMFSetFunctioni - Sets the function for a single component

797:    Collective on Mat

799:    Input Parameters:
800: +  mat - the matrix free matrix created via MatCreateSNESMF()
801: -  funci - the function to use

803:    Level: advanced

805:    Notes:
806:     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
807:     matrix inside your compute Jacobian routine


810: .keywords: SNES, matrix-free, function

812: .seealso: MatCreateSNESMF(),MatSNESMFGetH(),
813:           MatSNESMFSetHHistory(), MatSNESMFResetHHistory(),
814:           MatSNESMFKSPMonitor(), SNESetFunction()
815: @*/
816: int MatSNESMFSetFunctioni(Mat mat,int (*funci)(int,Vec,PetscScalar*,void *))
817: {
818:   int  ierr,(*f)(Mat,int (*)(int,Vec,PetscScalar*,void *));

822:   PetscObjectQueryFunction((PetscObject)mat,"MatSNESMFSetFunctioni_C",(void (**)(void))&f);
823:   if (f) {
824:     (*f)(mat,funci);
825:   }
826:   return(0);
827: }


832: /*@C
833:    MatSNESMFSetFunctioniBase - Sets the base vector for a single component function evaluation

835:    Collective on Mat

837:    Input Parameters:
838: +  mat - the matrix free matrix created via MatCreateSNESMF()
839: -  func - the function to use

841:    Level: advanced

843:    Notes:
844:     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
845:     matrix inside your compute Jacobian routine


848: .keywords: SNES, matrix-free, function

850: .seealso: MatCreateSNESMF(),MatSNESMFGetH(),
851:           MatSNESMFSetHHistory(), MatSNESMFResetHHistory(),
852:           MatSNESMFKSPMonitor(), SNESetFunction()
853: @*/
854: int MatSNESMFSetFunctioniBase(Mat mat,int (*func)(Vec,void *))
855: {
856:   int  ierr,(*f)(Mat,int (*)(Vec,void *));

860:   PetscObjectQueryFunction((PetscObject)mat,"MatSNESMFSetFunctioniBase_C",(void (**)(void))&f);
861:   if (f) {
862:     (*f)(mat,func);
863:   }
864:   return(0);
865: }


870: /*@
871:    MatSNESMFSetPeriod - Sets how often h is recomputed, by default it is everytime

873:    Collective on Mat

875:    Input Parameters:
876: +  mat - the matrix free matrix created via MatCreateSNESMF()
877: -  period - 1 for everytime, 2 for every second etc

879:    Options Database Keys:
880: +  -snes_mf_period <period>

882:    Level: advanced


885: .keywords: SNES, matrix-free, parameters

887: .seealso: MatCreateSNESMF(),MatSNESMFGetH(),
888:           MatSNESMFSetHHistory(), MatSNESMFResetHHistory(),
889:           MatSNESMFKSPMonitor()
890: @*/
891: int MatSNESMFSetPeriod(Mat mat,int period)
892: {
893:   MatSNESMFCtx ctx = (MatSNESMFCtx)mat->data;

896:   ctx->recomputeperiod = period;
897:   return(0);
898: }

902: /*@
903:    MatSNESMFSetFunctionError - Sets the error_rel for the approximation of
904:    matrix-vector products using finite differences.

906:    Collective on Mat

908:    Input Parameters:
909: +  mat - the matrix free matrix created via MatCreateSNESMF()
910: -  error_rel - relative error (should be set to the square root of
911:                the relative error in the function evaluations)

913:    Options Database Keys:
914: +  -snes_mf_err <error_rel> - Sets error_rel

916:    Level: advanced

918:    Notes:
919:    The default matrix-free matrix-vector product routine computes
920: .vb
921:      F'(u)*a = [F(u+h*a) - F(u)]/h where
922:      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
923:        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   else
924: .ve

926: .keywords: SNES, matrix-free, parameters

928: .seealso: MatCreateSNESMF(),MatSNESMFGetH(),
929:           MatSNESMFSetHHistory(), MatSNESMFResetHHistory(),
930:           MatSNESMFKSPMonitor()
931: @*/
932: int MatSNESMFSetFunctionError(Mat mat,PetscReal error)
933: {
934:   MatSNESMFCtx ctx = (MatSNESMFCtx)mat->data;

937:   if (error != PETSC_DEFAULT) ctx->error_rel = error;
938:   return(0);
939: }

943: /*@
944:    MatSNESMFAddNullSpace - Provides a null space that an operator is
945:    supposed to have.  Since roundoff will create a small component in
946:    the null space, if you know the null space you may have it
947:    automatically removed.

949:    Collective on Mat 

951:    Input Parameters:
952: +  J - the matrix-free matrix context
953: -  nullsp - object created with MatNullSpaceCreate()

955:    Level: advanced

957: .keywords: SNES, matrix-free, null space

959: .seealso: MatNullSpaceCreate(), MatSNESMFGetH(), MatCreateSNESMF(),
960:           MatSNESMFSetHHistory(), MatSNESMFResetHHistory(),
961:           MatSNESMFKSPMonitor(), MatSNESMFErrorRel()
962: @*/
963: int MatSNESMFAddNullSpace(Mat J,MatNullSpace nullsp)
964: {
965:   int          ierr;
966:   MatSNESMFCtx ctx = (MatSNESMFCtx)J->data;
967:   MPI_Comm     comm;

970:   PetscObjectGetComm((PetscObject)J,&comm);

972:   ctx->sp = nullsp;
973:   PetscObjectReference((PetscObject)nullsp);
974:   return(0);
975: }

979: /*@
980:    MatSNESMFSetHHistory - Sets an array to collect a history of the
981:    differencing values (h) computed for the matrix-free product.

983:    Collective on Mat 

985:    Input Parameters:
986: +  J - the matrix-free matrix context
987: .  histroy - space to hold the history
988: -  nhistory - number of entries in history, if more entries are generated than
989:               nhistory, then the later ones are discarded

991:    Level: advanced

993:    Notes:
994:    Use MatSNESMFResetHHistory() to reset the history counter and collect
995:    a new batch of differencing parameters, h.

997: .keywords: SNES, matrix-free, h history, differencing history

999: .seealso: MatSNESMFGetH(), MatCreateSNESMF(),
1000:           MatSNESMFResetHHistory(),
1001:           MatSNESMFKSPMonitor(), MatSNESMFSetFunctionError()

1003: @*/
1004: int MatSNESMFSetHHistory(Mat J,PetscScalar history[],int nhistory)
1005: {
1006:   MatSNESMFCtx ctx = (MatSNESMFCtx)J->data;

1009:   ctx->historyh    = history;
1010:   ctx->maxcurrenth = nhistory;
1011:   ctx->currenth    = 0;
1012:   return(0);
1013: }

1017: /*@
1018:    MatSNESMFResetHHistory - Resets the counter to zero to begin 
1019:    collecting a new set of differencing histories.

1021:    Collective on Mat 

1023:    Input Parameters:
1024: .  J - the matrix-free matrix context

1026:    Level: advanced

1028:    Notes:
1029:    Use MatSNESMFSetHHistory() to create the original history counter.

1031: .keywords: SNES, matrix-free, h history, differencing history

1033: .seealso: MatSNESMFGetH(), MatCreateSNESMF(),
1034:           MatSNESMFSetHHistory(),
1035:           MatSNESMFKSPMonitor(), MatSNESMFSetFunctionError()

1037: @*/
1038: int MatSNESMFResetHHistory(Mat J)
1039: {
1040:   MatSNESMFCtx ctx = (MatSNESMFCtx)J->data;

1043:   ctx->ncurrenth    = 0;
1044:   return(0);
1045: }

1049: int MatSNESMFComputeJacobian(SNES snes,Vec x,Mat *jac,Mat *B,MatStructure *flag,void *dummy)
1050: {
1053:   MatAssemblyBegin(*jac,MAT_FINAL_ASSEMBLY);
1054:   MatAssemblyEnd(*jac,MAT_FINAL_ASSEMBLY);
1055:   return(0);
1056: }

1060: /*@
1061:     MatSNESMFSetBase - Sets the vector U at which matrix vector products of the 
1062:         Jacobian are computed

1064:     Collective on Mat

1066:     Input Parameters:
1067: +   J - the MatSNESMF matrix
1068: -   U - the vector

1070:     Notes: This is rarely used directly

1072:     Level: advanced

1074: @*/
1075: int MatSNESMFSetBase(Mat J,Vec U)
1076: {
1077:   int  ierr,(*f)(Mat,Vec);

1082:   PetscObjectQueryFunction((PetscObject)J,"MatSNESMFSetBase_C",(void (**)(void))&f);
1083:   if (f) {
1084:     (*f)(J,U);
1085:   }
1086:   return(0);
1087: }

1091: /*@C
1092:     MatSNESMFSetCheckh - Sets a function that checks the computed h and adjusts
1093:         it to satisfy some criteria

1095:     Collective on Mat

1097:     Input Parameters:
1098: +   J - the MatSNESMF matrix
1099: .   fun - the function that checks h
1100: -   ctx - any context needed by the function

1102:     Options Database Keys:
1103: .   -snes_mf_check_positivity

1105:     Level: advanced

1107:     Notes: For example, MatSNESMFSetCheckPositivity() insures that all entries
1108:        of U + h*a are non-negative

1110: .seealso:  MatSNESMFSetCheckPositivity()
1111: @*/
1112: int MatSNESMFSetCheckh(Mat J,int (*fun)(Vec,Vec,PetscScalar*,void*),void* ctx)
1113: {
1114:   int  ierr,(*f)(Mat,int (*)(Vec,Vec,PetscScalar*,void*),void*);

1118:   PetscObjectQueryFunction((PetscObject)J,"MatSNESMFSetCheckh_C",(void (**)(void))&f);
1119:   if (f) {
1120:     (*f)(J,fun,ctx);
1121:   }
1122:   return(0);
1123: }

1127: /*@
1128:     MatSNESMFCheckPositivity - Checks that all entries in U + h*a are positive or
1129:         zero, decreases h until this is satisfied.

1131:     Collective on Vec

1133:     Input Parameters:
1134: +   U - base vector that is added to
1135: .   a - vector that is added
1136: .   h - scaling factor on a
1137: -   dummy - context variable (unused)

1139:     Options Database Keys:
1140: .   -snes_mf_check_positivity

1142:     Level: advanced

1144:     Notes: This is rarely used directly, rather it is passed as an argument to 
1145:            MatSNESMFSetCheckh()

1147: .seealso:  MatSNESMFSetCheckh()
1148: @*/
1149: int MatSNESMFCheckPositivity(Vec U,Vec a,PetscScalar *h,void *dummy)
1150: {
1151:   PetscReal     val, minval;
1152:   PetscScalar   *u_vec, *a_vec;
1153:   int           ierr, i, size;
1154:   MPI_Comm      comm;

1157:   PetscObjectGetComm((PetscObject)U,&comm);
1158:   VecGetArray(U,&u_vec);
1159:   VecGetArray(a,&a_vec);
1160:   VecGetLocalSize(U,&size);
1161:   minval = PetscAbsScalar(*h*1.01);
1162:   for(i=0;i<size;i++) {
1163:     if (PetscRealPart(u_vec[i] + *h*a_vec[i]) <= 0.0) {
1164:       val = PetscAbsScalar(u_vec[i]/a_vec[i]);
1165:       if (val < minval) minval = val;
1166:     }
1167:   }
1168:   VecRestoreArray(U,&u_vec);
1169:   VecRestoreArray(a,&a_vec);
1170:   PetscGlobalMin(&minval,&val,comm);
1171:   if (val <= PetscAbsScalar(*h)) {
1172:     PetscLogInfo(U,"MatSNESMFCheckPositivity: Scaling back h from %g to %g\n",PetscRealPart(*h),.99*val);
1173:     if (PetscRealPart(*h) > 0.0) *h =  0.99*val;
1174:     else                         *h = -0.99*val;
1175:   }
1176:   return(0);
1177: }