Actual source code: mffddef.c

petsc-3.4.5 2014-06-29
  2: /*
  3:   Implements the DS PETSc approach for computing the h
  4:   parameter used with the finite difference based matrix-free
  5:   Jacobian-vector products.

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

  9:   Mandatory functions:
 10:   -------------------
 11:       MatMFFDCompute_  - for a given point and direction computes h

 13:       MatCreateMFFD _ - fills in the MatMFFD data structure
 14:                            for this particular implementation


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

 22:       MatMFFDSetFromOptions_ - checks the options database for options that
 23:                                apply to this method.

 25:       MatMFFDDestroy_ - frees any space allocated by the routines above

 27: */

 29: /*
 30:     This include file defines the data structure  MatMFFD that
 31:    includes information about the computation of h. It is shared by
 32:    all implementations that people provide
 33: */
 34: #include <petsc-private/matimpl.h>
 35: #include <../src/mat/impls/mffd/mffdimpl.h>   /*I  "petscmat.h"   I*/

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

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

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


 61:    Output Parameter:
 62: .  h - the scale computed

 64: */
 65: static PetscErrorCode MatMFFDCompute_DS(MatMFFD ctx,Vec U,Vec a,PetscScalar *h,PetscBool  *zeroa)
 66: {
 67:   MatMFFD_DS     *hctx = (MatMFFD_DS*)ctx->hctx;
 68:   PetscReal      nrm,sum,umin = hctx->umin;
 69:   PetscScalar    dot;

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

 86:     if (nrm == 0.0) {
 87:       *zeroa = PETSC_TRUE;
 88:       return(0);
 89:     }
 90:     *zeroa = PETSC_FALSE;

 92:     /*
 93:       Safeguard for step sizes that are "too small"
 94:     */
 95:     if (PetscAbsScalar(dot) < umin*sum && PetscRealPart(dot) >= 0.0) dot = umin*sum;
 96:     else if (PetscAbsScalar(dot) < 0.0 && PetscRealPart(dot) > -umin*sum) dot = -umin*sum;
 97:     *h = ctx->error_rel*dot/(nrm*nrm);
 98:   } else {
 99:     *h = ctx->currenth;
100:   }
101:   if (*h != *h) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Differencing parameter is not a number sum = %G dot = %G norm = %G",sum,PetscRealPart(dot),nrm);
102:   ctx->count++;
103:   return(0);
104: }

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

114:    Input Parameters:
115: +  ctx - the matrix free context
116: -  viewer - the PETSc viewer
117: */
118: static PetscErrorCode MatMFFDView_DS(MatMFFD ctx,PetscViewer viewer)
119: {
120:   MatMFFD_DS     *hctx = (MatMFFD_DS*)ctx->hctx;
122:   PetscBool      iascii;

125:   /*
126:      Currently this only handles the ascii file viewers, others
127:      could be added, but for this type of object other viewers
128:      make less sense
129:   */
130:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
131:   if (iascii) {
132:     PetscViewerASCIIPrintf(viewer,"    umin=%G (minimum iterate parameter)\n",hctx->umin);
133:   }
134:   return(0);
135: }

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

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

146: */
147: static PetscErrorCode MatMFFDSetFromOptions_DS(MatMFFD ctx)
148: {
150:   MatMFFD_DS     *hctx = (MatMFFD_DS*)ctx->hctx;

153:   PetscOptionsHead("Finite difference matrix free parameters");
154:   PetscOptionsReal("-mat_mffd_umin","umin","MatMFFDDSSetUmin",hctx->umin,&hctx->umin,0);
155:   PetscOptionsTail();
156:   return(0);
157: }

161: /*
162:    MatMFFDDestroy_DS - Frees the space allocated by
163:    MatCreateMFFD_DS().

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

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

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

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

192:   if (!ctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"MatMFFDDSSetUmin() attached to non-shell matrix");
193:   hctx       = (MatMFFD_DS*)ctx->hctx;
194:   hctx->umin = umin;
195:   return(0);
196: }

200: /*@
201:     MatMFFDDSSetUmin - Sets the "umin" parameter used by the
202:     PETSc routine for computing the differencing parameter, h, which is used
203:     for matrix-free Jacobian-vector products.

205:    Input Parameters:
206: +  A - the matrix created with MatCreateSNESMF()
207: -  umin - the parameter

209:    Level: advanced

211:    Notes:
212:    See the manual page for MatCreateSNESMF() for a complete description of the
213:    algorithm used to compute h.

215: .seealso: MatMFFDSetFunctionError(), MatCreateSNESMF()

217: @*/
218: PetscErrorCode  MatMFFDDSSetUmin(Mat A,PetscReal umin)
219: {

224:   PetscTryMethod(A,"MatMFFDDSSetUmin_C",(Mat,PetscReal),(A,umin));
225:   return(0);
226: }

228: /*MC
229:      MATMFFD_DS - the code for compute the "h" used in the finite difference
230:             matrix-free matrix vector product.  This code
231:         implements the strategy in Dennis and Schnabel, "Numerical Methods for Unconstrained
232:         Optimization and Nonlinear Equations".

234:    Options Database Keys:
235: .  -mat_mffd_umin <umin> see MatMFFDDSSetUmin()

237:    Level: intermediate

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

243:    Formula used:
244:      F'(u)*a = [F(u+h*a) - F(u)]/h where
245:      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
246:        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
247:  where
248:      error_rel = square root of relative error in function evaluation
249:      umin = minimum iterate parameter

251: .seealso: MATMFFD, MatCreateMFFD(), MatCreateSNESMF(), MATMFFD_WP, MatMFFDDSSetUmin()

253: M*/
256: PETSC_EXTERN PetscErrorCode MatCreateMFFD_DS(MatMFFD ctx)
257: {
258:   MatMFFD_DS     *hctx;

262:   /* allocate my own private data structure */
263:   PetscNewLog(ctx,MatMFFD_DS,&hctx);
264:   ctx->hctx = (void*)hctx;
265:   /* set a default for my parameter */
266:   hctx->umin = 1.e-6;

268:   /* set the functions I am providing */
269:   ctx->ops->compute        = MatMFFDCompute_DS;
270:   ctx->ops->destroy        = MatMFFDDestroy_DS;
271:   ctx->ops->view           = MatMFFDView_DS;
272:   ctx->ops->setfromoptions = MatMFFDSetFromOptions_DS;

274:   PetscObjectComposeFunction((PetscObject)ctx->mat,"MatMFFDDSSetUmin_C",MatMFFDDSSetUmin_DS);
275:   return(0);
276: }