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;

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

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

 18:     Level: advanced

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

 24:                         F(x+ha) - F(x)
 25:           F'(u)a  ~=  ----------------
 26:                               h

 28: .seealso: MatCreateSNESMF(), MatSNESMFRegisterDynamic)
 29: @*/
 30: int MatSNESMFSetType(Mat mat,MatSNESMFType ftype)
 31: {
 32:   int          ierr,(*r)(MatSNESMFCtx);
 33:   MatSNESMFCtx ctx = (MatSNESMFCtx)mat->data;
 34:   PetscTruth   match;
 35: 

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

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

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

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

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

 56:   (*r)(ctx);

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

 60:   return(0);
 61: }

 63: EXTERN_C_BEGIN
 64: int MatSNESMFSetFunctioniBase_FD(Mat mat,int (*func)(Vec,void *))
 65: {
 66:   MatSNESMFCtx ctx = (MatSNESMFCtx)mat->data;

 69:   ctx->funcisetbase = func;
 70:   return(0);
 71: }
 72: EXTERN_C_END

 74: EXTERN_C_BEGIN
 75: int MatSNESMFSetFunctioni_FD(Mat mat,int (*funci)(int,Vec,PetscScalar*,void *))
 76: {
 77:   MatSNESMFCtx ctx = (MatSNESMFCtx)mat->data;

 80:   ctx->funci = funci;
 81:   return(0);
 82: }
 83: EXTERN_C_END

 85: /*MC
 86:    MatSNESMFRegisterDynamic - Adds a method to the MatSNESMF registry.

 88:    Synopsis:
 89:    int MatSNESMFRegisterDynamic(char *name_solver,char *path,char *name_create,int (*routine_create)(MatSNESMF))

 91:    Not Collective

 93:    Input Parameters:
 94: +  name_solver - name of a new user-defined compute-h module
 95: .  path - path (either absolute or relative) the library containing this solver
 96: .  name_create - name of routine to create method context
 97: -  routine_create - routine to create method context

 99:    Level: developer

101:    Notes:
102:    MatSNESMFRegisterDynamic) may be called multiple times to add several user-defined solvers.

104:    If dynamic libraries are used, then the fourth input argument (routine_create)
105:    is ignored.

107:    Sample usage:
108: .vb
109:    MatSNESMFRegisterDynamic"my_h",/home/username/my_lib/lib/libO/solaris/mylib.a,
110:                "MyHCreate",MyHCreate);
111: .ve

113:    Then, your solver can be chosen with the procedural interface via
114: $     MatSNESMFSetType(mfctx,"my_h")
115:    or at runtime via the option
116: $     -snes_mf_type my_h

118: .keywords: MatSNESMF, register

120: .seealso: MatSNESMFRegisterAll(), MatSNESMFRegisterDestroy()
121: M*/

123: int MatSNESMFRegister(char *sname,char *path,char *name,int (*function)(MatSNESMFCtx))
124: {
126:   char fullname[256];

129:   PetscFListConcat(path,name,fullname);
130:   PetscFListAdd(&MatSNESMPetscFList,sname,fullname,(void (*)(void))function);
131:   return(0);
132: }


135: /*@C
136:    MatSNESMFRegisterDestroy - Frees the list of MatSNESMF methods that were
137:    registered by MatSNESMFRegisterDynamic).

139:    Not Collective

141:    Level: developer

143: .keywords: MatSNESMF, register, destroy

145: .seealso: MatSNESMFRegisterDynamic), MatSNESMFRegisterAll()
146: @*/
147: int MatSNESMFRegisterDestroy(void)
148: {

152:   if (MatSNESMPetscFList) {
153:     PetscFListDestroy(&MatSNESMPetscFList);
154:     MatSNESMPetscFList = 0;
155:   }
156:   MatSNESMFRegisterAllCalled = PETSC_FALSE;
157:   return(0);
158: }

160: /* ----------------------------------------------------------------------------------------*/
161: int MatDestroy_MFFD(Mat mat)
162: {
163:   int          ierr;
164:   MatSNESMFCtx ctx = (MatSNESMFCtx)mat->data;

167:   VecDestroy(ctx->w);
168:   if (ctx->ops->destroy) {(*ctx->ops->destroy)(ctx);}
169:   if (ctx->sp) {MatNullSpaceDestroy(ctx->sp);}
170:   PetscHeaderDestroy(ctx);
171:   return(0);
172: }

174: /*
175:    MatSNESMFView_MFFD - Views matrix-free parameters.

177: */
178: int MatView_MFFD(Mat J,PetscViewer viewer)
179: {
180:   int          ierr;
181:   MatSNESMFCtx ctx = (MatSNESMFCtx)J->data;
182:   PetscTruth   isascii;

185:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
186:   if (isascii) {
187:      PetscViewerASCIIPrintf(viewer,"  SNES matrix-free approximation:n");
188:      PetscViewerASCIIPrintf(viewer,"    err=%g (relative error in function evaluation)n",ctx->error_rel);
189:      if (!ctx->type_name) {
190:        PetscViewerASCIIPrintf(viewer,"    The compute h routine has not yet been setn");
191:      } else {
192:        PetscViewerASCIIPrintf(viewer,"    Using %s compute h routinen",ctx->type_name);
193:      }
194:      if (ctx->ops->view) {
195:        (*ctx->ops->view)(ctx,viewer);
196:      }
197:   } else {
198:     SETERRQ1(1,"Viewer type %s not supported for SNES matrix free matrix",((PetscObject)viewer)->type_name);
199:   }
200:   return(0);
201: }

203: /*
204:    MatSNESMFAssemblyEnd_Private - Resets the ctx->ncurrenth to zero. This 
205:    allows the user to indicate the beginning of a new linear solve by calling
206:    MatAssemblyXXX() on the matrix free matrix. This then allows the 
207:    MatSNESMFCreate_WP() to properly compute ||U|| only the first time
208:    in the linear solver rather than every time.
209: */
210: int MatAssemblyEnd_MFFD(Mat J,MatAssemblyType mt)
211: {
212:   int             ierr;
213:   MatSNESMFCtx    j = (MatSNESMFCtx)J->data;

216:   MatSNESMFResetHHistory(J);
217:   if (j->usesnes) {
218:     SNESGetSolution(j->snes,&j->current_u);
219:     SNESGetFunction(j->snes,&j->current_f,PETSC_NULL,PETSC_NULL);
220:   }
221:   j->vshift = 0.0;
222:   j->vscale = 1.0;
223:   return(0);
224: }

226: /*
227:   MatSNESMFMult_Private - Default matrix-free form for Jacobian-vector
228:   product, y = F'(u)*a:

230:         y ~= (F(u + ha) - F(u))/h, 
231:   where F = nonlinear function, as set by SNESSetFunction()
232:         u = current iterate
233:         h = difference interval
234: */
235: int MatMult_MFFD(Mat mat,Vec a,Vec y)
236: {
237:   MatSNESMFCtx    ctx = (MatSNESMFCtx)mat->data;
238:   SNES            snes;
239:   PetscScalar     h,mone = -1.0;
240:   Vec             w,U,F;
241:   int             ierr,(*eval_fct)(SNES,Vec,Vec)=0;

244:   /* We log matrix-free matrix-vector products separately, so that we can
245:      separate the performance monitoring from the cases that use conventional
246:      storage.  We may eventually modify event logging to associate events
247:      with particular objects, hence alleviating the more general problem. */
248:   PetscLogEventBegin(MAT_MultMatrixFree,a,y,0,0);

250:   snes = ctx->snes;
251:   w    = ctx->w;
252:   U    = ctx->current_u;

254:   /* 
255:       Compute differencing parameter 
256:   */
257:   if (!ctx->ops->compute) {
258:     MatSNESMFSetType(mat,MATSNESMF_DEFAULT);
259:     MatSNESMFSetFromOptions(mat);
260:   }
261:   (*ctx->ops->compute)(ctx,U,a,&h);

263:   /* keep a record of the current differencing parameter h */
264:   ctx->currenth = h;
265: #if defined(PETSC_USE_COMPLEX)
266:   PetscLogInfo(mat,"MatMult_MFFD:Current differencing parameter: %g + %g in",PetscRealPart(h),PetscImaginaryPart(h));
267: #else
268:   PetscLogInfo(mat,"MatMult_MFFD:Current differencing parameter: %15.12en",h);
269: #endif
270:   if (ctx->historyh && ctx->ncurrenth < ctx->maxcurrenth) {
271:     ctx->historyh[ctx->ncurrenth] = h;
272:   }
273:   ctx->ncurrenth++;

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

278:   if (ctx->usesnes) {
279:     eval_fct = SNESComputeFunction;
280:     F    = ctx->current_f;
281:     if (!F) SETERRQ(1,"You must call MatAssembly() even on matrix-free matrices");
282:     (*eval_fct)(snes,w,y);
283:   } else {
284:     F = ctx->funcvec;
285:     /* compute func(U) as base for differencing */
286:     if (ctx->ncurrenth == 1) {
287:       (*ctx->func)(snes,U,F,ctx->funcctx);
288:     }
289:     (*ctx->func)(snes,w,y,ctx->funcctx);
290:   }

292:   VecAXPY(&mone,F,y);
293:   h    = 1.0/h;
294:   VecScale(&h,y);


297:   if (ctx->vshift != 0.0 && ctx->vscale != 1.0) {
298:     VecAXPBY(&ctx->vshift,&ctx->vscale,a,y);
299:   } else if (ctx->vscale != 1.0) {
300:     VecScale(&ctx->vscale,y);
301:   } else if (ctx->vshift != 0.0) {
302:     VecAXPY(&ctx->vshift,a,y);
303:   }

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

307:   PetscLogEventEnd(MAT_MultMatrixFree,a,y,0,0);
308:   return(0);
309: }

311: /*
312:   MatGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix

314:         y ~= (F(u + ha) - F(u))/h, 
315:   where F = nonlinear function, as set by SNESSetFunction()
316:         u = current iterate
317:         h = difference interval
318: */
319: int MatGetDiagonal_MFFD(Mat mat,Vec a)
320: {
321:   MatSNESMFCtx ctx = (MatSNESMFCtx)mat->data;
322:   PetscScalar  h,*aa,*ww,v;
323:   PetscReal    epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON;
324:   Vec          w,U;
325:   int          i,ierr,rstart,rend;

328:   if (!ctx->funci) {
329:     SETERRQ(1,"Requirers calling MatSNESMFSetFunctioni() first");
330:   }

332:   w    = ctx->w;
333:   U    = ctx->current_u;
334:   (*ctx->func)(0,U,a,ctx->funcctx);
335:   (*ctx->funcisetbase)(U,ctx->funcctx);
336:   VecCopy(U,w);

338:   VecGetOwnershipRange(a,&rstart,&rend);
339:   VecGetArray(a,&aa);
340:   for (i=rstart; i<rend; i++) {
341:     VecGetArray(w,&ww);
342:     h  = ww[i-rstart];
343:     if (h == 0.0) h = 1.0;
344: #if !defined(PETSC_USE_COMPLEX)
345:     if (h < umin && h >= 0.0)      h = umin;
346:     else if (h < 0.0 && h > -umin) h = -umin;
347: #else
348:     if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0)     h = umin;
349:     else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin;
350: #endif
351:     h     *= epsilon;
352: 
353:     ww[i-rstart] += h;
354:     VecRestoreArray(w,&ww);
355:     ierr          = (*ctx->funci)(i,w,&v,ctx->funcctx);
356:     aa[i-rstart]  = (v - aa[i-rstart])/h;

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

361:     VecGetArray(w,&ww);
362:     ww[i-rstart] -= h;
363:     VecRestoreArray(w,&ww);
364:   }
365:   VecRestoreArray(a,&aa);
366:   return(0);
367: }

369: int MatShift_MFFD(PetscScalar *a,Mat Y)
370: {
371:   MatSNESMFCtx shell = (MatSNESMFCtx)Y->data;
373:   shell->vshift += *a;
374:   return(0);
375: }

377: int MatScale_MFFD(PetscScalar *a,Mat Y)
378: {
379:   MatSNESMFCtx shell = (MatSNESMFCtx)Y->data;
381:   shell->vscale *= *a;
382:   return(0);
383: }


386: /*@C
387:    MatCreateSNESMF - Creates a matrix-free matrix context for use with
388:    a SNES solver.  This matrix can be used as the Jacobian argument for
389:    the routine SNESSetJacobian().

391:    Collective on SNES and Vec

393:    Input Parameters:
394: +  snes - the SNES context
395: -  x - vector where SNES solution is to be stored.

397:    Output Parameter:
398: .  J - the matrix-free matrix

400:    Level: advanced

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

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

409: .vb
410:      F'(u)*a = [F(u+h*a) - F(u)]/h where
411:      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
412:        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
413:  where
414:      error_rel = square root of relative error in function evaluation
415:      umin = minimum iterate parameter
416: .ve

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

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

425:    Options Database Keys:
426: +  -snes_mf_err <error_rel> - Sets error_rel
427: .  -snes_mf_unim <umin> - Sets umin (for default PETSc routine that computes h only)
428: -  -snes_mf_ksp_monitor - KSP monitor routine that prints differencing h

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

432: .seealso: MatDestroy(), MatSNESMFSetFunctionError(), MatSNESMFDefaultSetUmin()
433:           MatSNESMFSetHHistory(), MatSNESMFResetHHistory(), MatCreateMF(),
434:           MatSNESMFGetH(),MatSNESMFKSPMonitor(), MatSNESMFRegisterDynamic), MatSNESMFComputeJacobian()
435:  
436: @*/
437: int MatCreateSNESMF(SNES snes,Vec x,Mat *J)
438: {
439:   MatSNESMFCtx mfctx;
440:   int          ierr;

443:   MatCreateMF(x,J);

445:   mfctx          = (MatSNESMFCtx)(*J)->data;
446:   mfctx->snes    = snes;
447:   mfctx->usesnes = PETSC_TRUE;
448:   PetscLogObjectParent(snes,*J);
449:   return(0);
450: }

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

459:   MatSNESMFResetHHistory(J);
460:   ctx->current_u = U;
461:   ctx->usesnes   = PETSC_FALSE;
462:   return(0);
463: }
464: EXTERN_C_END

466: /*@
467:    MatSNESMFSetFromOptions - Sets the MatSNESMF options from the command line
468:    parameter.

470:    Collective on Mat

472:    Input Parameters:
473: .  mat - the matrix obtained with MatCreateSNESMF()

475:    Options Database Keys:
476: +  -snes_mf_type - <default,wp>
477: -  -snes_mf_err - square root of estimated relative error in function evaluation
478: -  -snes_mf_period - how often h is recomputed, defaults to 1, everytime

480:    Level: advanced

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

484: .seealso: MatCreateSNESMF(),MatSNESMFSetHHistory(), 
485:           MatSNESMFResetHHistory(), MatSNESMFKSPMonitor()
486: @*/
487: int MatSNESMFSetFromOptions(Mat mat)
488: {
489:   MatSNESMFCtx mfctx = (MatSNESMFCtx)mat->data;
490:   int          ierr;
491:   PetscTruth   flg;
492:   char         ftype[256];

495:   if (!MatSNESMFRegisterAllCalled) {MatSNESMFRegisterAll(PETSC_NULL);}
496: 
497:   PetscOptionsBegin(mfctx->comm,mfctx->prefix,"Set matrix free computation parameters","MatSNESMF");
498:   PetscOptionsList("-snes_mf_type","Matrix free type","MatSNESMFSetType",MatSNESMPetscFList,mfctx->type_name,ftype,256,&flg);
499:   if (flg) {
500:     MatSNESMFSetType(mat,ftype);
501:   }

503:   PetscOptionsReal("-snes_mf_err","set sqrt relative error in function","MatSNESMFSetFunctionError",mfctx->error_rel,&mfctx->error_rel,0);
504:   PetscOptionsInt("-snes_mf_period","how often h is recomputed","MatSNESMFSetPeriod",mfctx->recomputeperiod,&mfctx->recomputeperiod,0);
505:   if (mfctx->snes) {
506:     PetscOptionsName("-snes_mf_ksp_monitor","Monitor matrix-free parameters","MatSNESMFKSPMonitor",&flg);
507:     if (flg) {
508:       SLES sles;
509:       KSP  ksp;
510:       SNESGetSLES(mfctx->snes,&sles);
511:       SLESGetKSP(sles,&ksp);
512:       KSPSetMonitor(ksp,MatSNESMFKSPMonitor,PETSC_NULL,0);
513:     }
514:   }
515:   if (mfctx->ops->setfromoptions) {
516:     (*mfctx->ops->setfromoptions)(mfctx);
517:   }
518:   PetscOptionsEnd();
519:   return(0);
520: }

522: EXTERN_C_BEGIN
523: int MatCreate_MFFD(Mat A)
524: {
525:   MatSNESMFCtx mfctx;
526:   int          ierr;

529: #ifndef PETSC_USE_DYNAMIC_LIBRARIES
530:   SNESInitializePackage(PETSC_NULL);
531: #endif

533:   PetscHeaderCreate(mfctx,_p_MatSNESMFCtx,struct _MFOps,MATSNESMFCTX_COOKIE,0,"SNESMF",A->comm,MatDestroy_MFFD,MatView_MFFD);
534:   PetscLogObjectCreate(mfctx);
535:   mfctx->sp              = 0;
536:   mfctx->snes            = 0;
537:   mfctx->error_rel       = PETSC_SQRT_MACHINE_EPSILON;
538:   mfctx->recomputeperiod = 1;
539:   mfctx->count           = 0;
540:   mfctx->currenth        = 0.0;
541:   mfctx->historyh        = PETSC_NULL;
542:   mfctx->ncurrenth       = 0;
543:   mfctx->maxcurrenth     = 0;
544:   mfctx->type_name       = 0;
545:   mfctx->usesnes         = PETSC_FALSE;

547:   mfctx->vshift          = 0.0;
548:   mfctx->vscale          = 1.0;

550:   /* 
551:      Create the empty data structure to contain compute-h routines.
552:      These will be filled in below from the command line options or 
553:      a later call with MatSNESMFSetType() or if that is not called 
554:      then it will default in the first use of MatMult_MFFD()
555:   */
556:   mfctx->ops->compute        = 0;
557:   mfctx->ops->destroy        = 0;
558:   mfctx->ops->view           = 0;
559:   mfctx->ops->setfromoptions = 0;
560:   mfctx->hctx                = 0;

562:   mfctx->func                = 0;
563:   mfctx->funcctx             = 0;
564:   mfctx->funcvec             = 0;

566:   A->data                = mfctx;

568:   A->ops->mult           = MatMult_MFFD;
569:   A->ops->destroy        = MatDestroy_MFFD;
570:   A->ops->view           = MatView_MFFD;
571:   A->ops->assemblyend    = MatAssemblyEnd_MFFD;
572:   A->ops->getdiagonal    = MatGetDiagonal_MFFD;
573:   A->ops->scale          = MatScale_MFFD;
574:   A->ops->shift          = MatShift_MFFD;
575:   A->ops->setfromoptions = MatSNESMFSetFromOptions;

577:   PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSNESMFSetBase_C","MatSNESMFSetBase_FD",MatSNESMFSetBase_FD);
578:   PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSNESMFSetFunctioniBase_C","MatSNESMFSetFunctioniBase_FD",MatSNESMFSetFunctioniBase_FD);
579:   PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSNESMFSetFunctioni_C","MatSNESMFSetFunctioni_FD",MatSNESMFSetFunctioni_FD);
580:   mfctx->mat = A;
581:   VecCreateMPI(A->comm,A->n,A->N,&mfctx->w);

583:   return(0);
584: }

586: EXTERN_C_END

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

591:    Collective on Vec

593:    Input Parameters:
594: .  x - vector that defines layout of the vectors and matrices

596:    Output Parameter:
597: .  J - the matrix-free matrix

599:    Level: advanced

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

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

608: .vb
609:      F'(u)*a = [F(u+h*a) - F(u)]/h where
610:      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
611:        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
612:  where
613:      error_rel = square root of relative error in function evaluation
614:      umin = minimum iterate parameter
615: .ve

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

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

624:    Options Database Keys:
625: +  -snes_mf_err <error_rel> - Sets error_rel
626: .  -snes_mf_unim <umin> - Sets umin (for default PETSc routine that computes h only)
627: -  -snes_mf_ksp_monitor - KSP monitor routine that prints differencing h

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

631: .seealso: MatDestroy(), MatSNESMFSetFunctionError(), MatSNESMFDefaultSetUmin()
632:           MatSNESMFSetHHistory(), MatSNESMFResetHHistory(), MatCreateSNESMF(),
633:           MatSNESMFGetH(),MatSNESMFKSPMonitor(), MatSNESMFRegisterDynamic),, MatSNESMFComputeJacobian()
634:  
635: @*/
636: int MatCreateMF(Vec x,Mat *J)
637: {
638:   MPI_Comm     comm;
639:   int          n,nloc,ierr;

642:   PetscObjectGetComm((PetscObject)x,&comm);
643:   VecGetSize(x,&n);
644:   VecGetLocalSize(x,&nloc);
645:   MatCreate(comm,nloc,nloc,n,n,J);
646:   MatRegister(MATMFFD,0,"MatCreate_MFFD",MatCreate_MFFD);
647:   MatSetType(*J,MATMFFD);
648:   return(0);
649: }


652: /*@
653:    MatSNESMFGetH - Gets the last value that was used as the differencing 
654:    parameter.

656:    Not Collective

658:    Input Parameters:
659: .  mat - the matrix obtained with MatCreateSNESMF()

661:    Output Paramter:
662: .  h - the differencing step size

664:    Level: advanced

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

668: .seealso: MatCreateSNESMF(),MatSNESMFSetHHistory(), 
669:           MatSNESMFResetHHistory(),MatSNESMFKSPMonitor()
670: @*/
671: int MatSNESMFGetH(Mat mat,PetscScalar *h)
672: {
673:   MatSNESMFCtx ctx = (MatSNESMFCtx)mat->data;

676:   *h = ctx->currenth;
677:   return(0);
678: }

680: /*
681:    MatSNESMFKSPMonitor - A KSP monitor for use with the default PETSc
682:    SNES matrix free routines. Prints the differencing parameter used at 
683:    each step.
684: */
685: int MatSNESMFKSPMonitor(KSP ksp,int n,PetscReal rnorm,void *dummy)
686: {
687:   PC             pc;
688:   MatSNESMFCtx   ctx;
689:   int            ierr;
690:   Mat            mat;
691:   MPI_Comm       comm;
692:   PetscTruth     nonzeroinitialguess;

695:   PetscObjectGetComm((PetscObject)ksp,&comm);
696:   KSPGetPC(ksp,&pc);
697:   KSPGetInitialGuessNonzero(ksp,&nonzeroinitialguess);
698:   PCGetOperators(pc,&mat,PETSC_NULL,PETSC_NULL);
699:   ctx  = (MatSNESMFCtx)mat->data;

701:   if (n > 0 || nonzeroinitialguess) {
702: #if defined(PETSC_USE_COMPLEX)
703:     PetscPrintf(comm,"%d KSP Residual norm %14.12e h %g + %g in",n,rnorm,
704:                 PetscRealPart(ctx->currenth),PetscImaginaryPart(ctx->currenth));
705: #else
706:     PetscPrintf(comm,"%d KSP Residual norm %14.12e h %g n",n,rnorm,ctx->currenth);
707: #endif
708:   } else {
709:     PetscPrintf(comm,"%d KSP Residual norm %14.12en",n,rnorm);
710:   }
711:   return(0);
712: }

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

717:    Collective on Mat

719:    Input Parameters:
720: +  mat - the matrix free matrix created via MatCreateSNESMF()
721: .  v   - workspace vector
722: .  func - the function to use
723: -  funcctx - optional function context passed to function

725:    Level: advanced

727:    Notes:
728:     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
729:     matrix inside your compute Jacobian routine

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

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

735: .seealso: MatCreateSNESMF(),MatSNESMFGetH(),
736:           MatSNESMFSetHHistory(), MatSNESMFResetHHistory(),
737:           MatSNESMFKSPMonitor(), SNESetFunction()
738: @*/
739: int MatSNESMFSetFunction(Mat mat,Vec v,int (*func)(SNES,Vec,Vec,void *),void *funcctx)
740: {
741:   MatSNESMFCtx ctx = (MatSNESMFCtx)mat->data;

744:   ctx->func    = func;
745:   ctx->funcctx = funcctx;
746:   ctx->funcvec = v;
747:   return(0);
748: }

750: /*@C
751:    MatSNESMFSetFunctioni - Sets the function for a single component

753:    Collective on Mat

755:    Input Parameters:
756: +  mat - the matrix free matrix created via MatCreateSNESMF()
757: -  funci - the function to use

759:    Level: advanced

761:    Notes:
762:     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
763:     matrix inside your compute Jacobian routine


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

768: .seealso: MatCreateSNESMF(),MatSNESMFGetH(),
769:           MatSNESMFSetHHistory(), MatSNESMFResetHHistory(),
770:           MatSNESMFKSPMonitor(), SNESetFunction()
771: @*/
772: int MatSNESMFSetFunctioni(Mat mat,int (*funci)(int,Vec,PetscScalar*,void *))
773: {
774:   int  ierr,(*f)(Mat,int (*)(int,Vec,PetscScalar*,void *));

778:   PetscObjectQueryFunction((PetscObject)mat,"MatSNESMFSetFunctioni_C",(void (**)(void))&f);
779:   if (f) {
780:     (*f)(mat,funci);
781:   }
782:   return(0);
783: }


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

789:    Collective on Mat

791:    Input Parameters:
792: +  mat - the matrix free matrix created via MatCreateSNESMF()
793: -  func - the function to use

795:    Level: advanced

797:    Notes:
798:     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
799:     matrix inside your compute Jacobian routine


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

804: .seealso: MatCreateSNESMF(),MatSNESMFGetH(),
805:           MatSNESMFSetHHistory(), MatSNESMFResetHHistory(),
806:           MatSNESMFKSPMonitor(), SNESetFunction()
807: @*/
808: int MatSNESMFSetFunctioniBase(Mat mat,int (*func)(Vec,void *))
809: {
810:   int  ierr,(*f)(Mat,int (*)(Vec,void *));

814:   PetscObjectQueryFunction((PetscObject)mat,"MatSNESMFSetFunctioniBase_C",(void (**)(void))&f);
815:   if (f) {
816:     (*f)(mat,func);
817:   }
818:   return(0);
819: }


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

825:    Collective on Mat

827:    Input Parameters:
828: +  mat - the matrix free matrix created via MatCreateSNESMF()
829: -  period - 1 for everytime, 2 for every second etc

831:    Options Database Keys:
832: +  -snes_mf_period <period>

834:    Level: advanced


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

839: .seealso: MatCreateSNESMF(),MatSNESMFGetH(),
840:           MatSNESMFSetHHistory(), MatSNESMFResetHHistory(),
841:           MatSNESMFKSPMonitor()
842: @*/
843: int MatSNESMFSetPeriod(Mat mat,int period)
844: {
845:   MatSNESMFCtx ctx = (MatSNESMFCtx)mat->data;

848:   ctx->recomputeperiod = period;
849:   return(0);
850: }

852: /*@
853:    MatSNESMFSetFunctionError - Sets the error_rel for the approximation of
854:    matrix-vector products using finite differences.

856:    Collective on Mat

858:    Input Parameters:
859: +  mat - the matrix free matrix created via MatCreateSNESMF()
860: -  error_rel - relative error (should be set to the square root of
861:                the relative error in the function evaluations)

863:    Options Database Keys:
864: +  -snes_mf_err <error_rel> - Sets error_rel

866:    Level: advanced

868:    Notes:
869:    The default matrix-free matrix-vector product routine computes
870: .vb
871:      F'(u)*a = [F(u+h*a) - F(u)]/h where
872:      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
873:        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   else
874: .ve

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

878: .seealso: MatCreateSNESMF(),MatSNESMFGetH(),
879:           MatSNESMFSetHHistory(), MatSNESMFResetHHistory(),
880:           MatSNESMFKSPMonitor()
881: @*/
882: int MatSNESMFSetFunctionError(Mat mat,PetscReal error)
883: {
884:   MatSNESMFCtx ctx = (MatSNESMFCtx)mat->data;

887:   if (error != PETSC_DEFAULT) ctx->error_rel = error;
888:   return(0);
889: }

891: /*@
892:    MatSNESMFAddNullSpace - Provides a null space that an operator is
893:    supposed to have.  Since roundoff will create a small component in
894:    the null space, if you know the null space you may have it
895:    automatically removed.

897:    Collective on Mat 

899:    Input Parameters:
900: +  J - the matrix-free matrix context
901: -  nullsp - object created with MatNullSpaceCreate()

903:    Level: advanced

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

907: .seealso: MatNullSpaceCreate(), MatSNESMFGetH(), MatCreateSNESMF(),
908:           MatSNESMFSetHHistory(), MatSNESMFResetHHistory(),
909:           MatSNESMFKSPMonitor(), MatSNESMFErrorRel()
910: @*/
911: int MatSNESMFAddNullSpace(Mat J,MatNullSpace nullsp)
912: {
913:   int          ierr;
914:   MatSNESMFCtx ctx = (MatSNESMFCtx)J->data;
915:   MPI_Comm     comm;

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

920:   ctx->sp = nullsp;
921:   ierr    = PetscObjectReference((PetscObject)nullsp);
922:   return(0);
923: }

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

929:    Collective on Mat 

931:    Input Parameters:
932: +  J - the matrix-free matrix context
933: .  histroy - space to hold the history
934: -  nhistory - number of entries in history, if more entries are generated than
935:               nhistory, then the later ones are discarded

937:    Level: advanced

939:    Notes:
940:    Use MatSNESMFResetHHistory() to reset the history counter and collect
941:    a new batch of differencing parameters, h.

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

945: .seealso: MatSNESMFGetH(), MatCreateSNESMF(),
946:           MatSNESMFResetHHistory(),
947:           MatSNESMFKSPMonitor(), MatSNESMFSetFunctionError()

949: @*/
950: int MatSNESMFSetHHistory(Mat J,PetscScalar *history,int nhistory)
951: {
952:   MatSNESMFCtx ctx = (MatSNESMFCtx)J->data;

955:   ctx->historyh    = history;
956:   ctx->maxcurrenth = nhistory;
957:   ctx->currenth    = 0;
958:   return(0);
959: }

961: /*@
962:    MatSNESMFResetHHistory - Resets the counter to zero to begin 
963:    collecting a new set of differencing histories.

965:    Collective on Mat 

967:    Input Parameters:
968: .  J - the matrix-free matrix context

970:    Level: advanced

972:    Notes:
973:    Use MatSNESMFSetHHistory() to create the original history counter.

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

977: .seealso: MatSNESMFGetH(), MatCreateSNESMF(),
978:           MatSNESMFSetHHistory(),
979:           MatSNESMFKSPMonitor(), MatSNESMFSetFunctionError()

981: @*/
982: int MatSNESMFResetHHistory(Mat J)
983: {
984:   MatSNESMFCtx ctx = (MatSNESMFCtx)J->data;

987:   ctx->ncurrenth    = 0;
988:   return(0);
989: }

991: int MatSNESMFComputeJacobian(SNES snes,Vec x,Mat *jac,Mat *B,MatStructure *flag,void *dummy)
992: {
995:   MatAssemblyBegin(*jac,MAT_FINAL_ASSEMBLY);
996:   MatAssemblyEnd(*jac,MAT_FINAL_ASSEMBLY);
997:   return(0);
998: }

1000: int MatSNESMFSetBase(Mat J,Vec U)
1001: {
1002:   int  ierr,(*f)(Mat,Vec);

1007:   PetscObjectQueryFunction((PetscObject)J,"MatSNESMFSetBase_C",(void (**)(void))&f);
1008:   if (f) {
1009:     (*f)(J,U);
1010:   }
1011:   return(0);
1012: }