Actual source code: snesmfjdef.c

  1: /*
  2:   Implements the default PETSc approach for computing the h 
  3:   parameter used with the finite difference based matrix-free 
  4:   Jacobian-vector products.

  6:   To make your own: clone this file and modify for your needs.

  8:   Mandatory functions:
  9:   -------------------
 10:       MatSNESMFCompute_  - for a given point and direction computes h

 12:       MatSNESMFCreate_ - fills in the MatSNESMF data structure
 13:                            for this particular implementation

 15:       
 16:    Optional functions:
 17:    -------------------
 18:       MatSNESMFView_ - prints information about the parameters being used.
 19:                        This is called when SNESView() or -snes_view is used.

 21:       MatSNESMFSetFromOptions_ - checks the options database for options that 
 22:                                apply to this method.

 24:       MatSNESMFDestroy_ - frees any space allocated by the routines above

 26: */

 28: /*
 29:     This include file defines the data structure  MatSNESMF that 
 30:    includes information about the computation of h. It is shared by 
 31:    all implementations that people provide
 32: */
 33:  #include src/mat/matimpl.h
 34:  #include src/snes/mf/snesmfj.h

 36: /*
 37:       The default method has one parameter that is used to 
 38:    "cutoff" very small values. This is stored in a data structure
 39:    that is only visible to this file. If your method has no parameters
 40:    it can omit this, if it has several simply reorganize the data structure.
 41:    The data structure is "hung-off" the MatSNESMF data structure in
 42:    the void *hctx; field.
 43: */
 44: typedef struct {
 45:   PetscReal umin;          /* minimum allowable u'a value relative to |u|_1 */
 46: } MatSNESMFDefault;

 50: /*
 51:    MatSNESMFCompute_Default - Standard PETSc code for computing the
 52:    differencing paramter (h) for use with matrix-free finite differences.

 54:    Input Parameters:
 55: +  ctx - the matrix free context
 56: .  U - the location at which you want the Jacobian
 57: -  a - the direction you want the derivative

 59:   
 60:    Output Parameter:
 61: .  h - the scale computed

 63: */
 64: static PetscErrorCode MatSNESMFCompute_Default(MatSNESMFCtx ctx,Vec U,Vec a,PetscScalar *h)
 65: {
 66:   MatSNESMFDefault *hctx = (MatSNESMFDefault*)ctx->hctx;
 67:   PetscReal        nrm,sum,umin = hctx->umin;
 68:   PetscScalar      dot;
 69:   PetscErrorCode   ierr;

 72:   if (!(ctx->count % ctx->recomputeperiod)) {
 73:     /*
 74:      This algorithm requires 2 norms and 1 inner product. Rather than
 75:      use directly the VecNorm() and VecDot() routines (and thus have 
 76:      three separate collective operations, we use the VecxxxBegin/End() routines
 77:     */
 78:     VecDotBegin(U,a,&dot);
 79:     VecNormBegin(a,NORM_1,&sum);
 80:     VecNormBegin(a,NORM_2,&nrm);
 81:     VecDotEnd(U,a,&dot);
 82:     VecNormEnd(a,NORM_1,&sum);
 83:     VecNormEnd(a,NORM_2,&nrm);

 85:     /* 
 86:       Safeguard for step sizes that are "too small"
 87:     */
 88:     if (!sum) {dot = 1.0; nrm = 1.0;}
 89: #if defined(PETSC_USE_COMPLEX)
 90:     else if (PetscAbsScalar(dot) < umin*sum && PetscRealPart(dot) >= 0.0) dot = umin*sum;
 91:     else if (PetscAbsScalar(dot) < 0.0 && PetscRealPart(dot) > -umin*sum) dot = -umin*sum;
 92: #else
 93:     else if (dot < umin*sum && dot >= 0.0) dot = umin*sum;
 94:     else if (dot < 0.0 && dot > -umin*sum) dot = -umin*sum;
 95: #endif
 96:     *h = ctx->error_rel*dot/(nrm*nrm);
 97:   } else {
 98:     *h = ctx->currenth;
 99:   }
100:   if (*h != *h) SETERRQ3(PETSC_ERR_PLIB,"Differencing parameter is not a number sum = %g dot = %g norm = %g",sum,PetscRealPart(dot),nrm);
101:   ctx->count++;
102:   return(0);
103: }

107: /*
108:    MatSNESMFView_Default - Prints information about this particular 
109:    method for computing h. Note that this does not print the general
110:    information about the matrix-free method, as such info is printed
111:    by the calling routine.

113:    Input Parameters:
114: +  ctx - the matrix free context
115: -  viewer - the PETSc viewer
116: */
117: static PetscErrorCode MatSNESMFView_Default(MatSNESMFCtx ctx,PetscViewer viewer)
118: {
119:   MatSNESMFDefault *hctx = (MatSNESMFDefault *)ctx->hctx;
120:   PetscErrorCode   ierr;
121:   PetscTruth       iascii;

124:   /*
125:      Currently this only handles the ascii file viewers, others
126:      could be added, but for this type of object other viewers
127:      make less sense
128:   */
129:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
130:   if (iascii) {
131:     PetscViewerASCIIPrintf(viewer,"    umin=%g (minimum iterate parameter)\n",hctx->umin);
132:   } else {
133:     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for this SNES matrix free matrix",((PetscObject)viewer)->type_name);
134:   }
135:   return(0);
136: }

140: /*
141:    MatSNESMFSetFromOptions_Default - Looks in the options database for 
142:    any options appropriate for this method.

144:    Input Parameter:
145: .  ctx - the matrix free context

147: */
148: static PetscErrorCode MatSNESMFSetFromOptions_Default(MatSNESMFCtx ctx)
149: {
150:   PetscErrorCode   ierr;
151:   MatSNESMFDefault *hctx = (MatSNESMFDefault*)ctx->hctx;

154:   PetscOptionsHead("Default matrix free parameters");
155:     PetscOptionsReal("-snes_mf_umin","umin","MatSNESMFDefaultSetUmin",hctx->umin,&hctx->umin,0);
156:   PetscOptionsTail();
157:   return(0);
158: }

162: /*
163:    MatSNESMFDestroy_Default - Frees the space allocated by 
164:    MatSNESMFCreate_Default(). 

166:    Input Parameter:
167: .  ctx - the matrix free context

169:    Notes: 
170:    Does not free the ctx, that is handled by the calling routine
171: */
172: static PetscErrorCode MatSNESMFDestroy_Default(MatSNESMFCtx ctx)
173: {

177:   PetscFree(ctx->hctx);
178:   return(0);
179: }

184: /*
185:    The following two routines use the PetscObjectCompose() and PetscObjectQuery()
186:    mechanism to allow the user to change the Umin parameter used in this method.
187: */
188: PetscErrorCode MatSNESMFDefaultSetUmin_Private(Mat mat,PetscReal umin)
189: {
190:   MatSNESMFCtx     ctx = (MatSNESMFCtx)mat->data;
191:   MatSNESMFDefault *hctx;

194:   if (!ctx) {
195:     SETERRQ(PETSC_ERR_ARG_WRONG,"MatSNESMFDefaultSetUmin() attached to non-shell matrix");
196:   }
197:   hctx = (MatSNESMFDefault*)ctx->hctx;
198:   hctx->umin = umin;
199:   return(0);
200: }

205: /*@
206:     MatSNESMFDefaultSetUmin - Sets the "umin" parameter used by the default
207:     PETSc routine for computing the differencing parameter, h, which is used
208:     for matrix-free Jacobian-vector products.

210:    Input Parameters:
211: +  A - the matrix created with MatCreateSNESMF()
212: -  umin - the parameter

214:    Level: advanced

216:    Notes:
217:    See the manual page for MatCreateSNESMF() for a complete description of the
218:    algorithm used to compute h.

220: .seealso: MatSNESMFSetFunctionError(), MatCreateSNESMF()

222: @*/
223: PetscErrorCode MatSNESMFDefaultSetUmin(Mat A,PetscReal umin)
224: {
225:   PetscErrorCode ierr,(*f)(Mat,PetscReal);

229:   PetscObjectQueryFunction((PetscObject)A,"MatSNESMFDefaultSetUmin_C",(void (**)(void))&f);
230:   if (f) {
231:     (*f)(A,umin);
232:   }
233:   return(0);
234: }

236: /*MC
237:      MATSNESMF_DEFAULT - the default code for compute the "h" used in the finite difference
238:             matrix-free matrix vector product

240:    Options Database Keys:
241: .  -snes_mf_umin <umin> see MatSNESMFDefaultSetUmin()


244:    Level: intermediate

246:    Notes: Requires 2 norms and 1 inner product, but they are computed together
247:        so only one parallel collective operation is needed. See MATSNESMF_WP for a method
248:        (with GMRES) that requires NO collective operations.

250:    Formula used:
251:      F'(u)*a = [F(u+h*a) - F(u)]/h where
252:      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
253:        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
254:  where
255:      error_rel = square root of relative error in function evaluation
256:      umin = minimum iterate parameter

258: .seealso: MATMFFD, MatCreateMF(), MatCreateSNESMF(), MATSNESMF_WP, MatSNESMFDefaultSetUmin()

260: M*/
264: PetscErrorCode MatSNESMFCreate_Default(MatSNESMFCtx ctx)
265: {
266:   MatSNESMFDefault *hctx;
267:   PetscErrorCode   ierr;


271:   /* allocate my own private data structure */
272:   PetscNew(MatSNESMFDefault,&hctx);
273:   ctx->hctx  = (void*)hctx;
274:   /* set a default for my parameter */
275:   hctx->umin = 1.e-6;

277:   /* set the functions I am providing */
278:   ctx->ops->compute        = MatSNESMFCompute_Default;
279:   ctx->ops->destroy        = MatSNESMFDestroy_Default;
280:   ctx->ops->view           = MatSNESMFView_Default;
281:   ctx->ops->setfromoptions = MatSNESMFSetFromOptions_Default;

283:   PetscObjectComposeFunctionDynamic((PetscObject)ctx->mat,"MatSNESMFDefaultSetUmin_C",
284:                             "MatSNESMFDefaultSetUmin_Private",
285:                              MatSNESMFDefaultSetUmin_Private);
286:   return(0);
287: }