Actual source code: snesmfj2.c

  1: #define PETSCSNES_DLL

 3:  #include private/snesimpl.h
  4: /* matimpl.h is needed only for logging of matrix operation */
 5:  #include private/matimpl.h

  7: EXTERN PetscErrorCode DiffParameterCreate_More(SNES,Vec,void**);
  8: EXTERN PetscErrorCode DiffParameterCompute_More(SNES,void*,Vec,Vec,PetscReal*,PetscReal*);
  9: EXTERN PetscErrorCode DiffParameterDestroy_More(void*);

 11: typedef struct {  /* default context for matrix-free SNES */
 12:   SNES         snes;             /* SNES context */
 13:   Vec          w;                /* work vector */
 14:   MatNullSpace sp;               /* null space context */
 15:   PetscReal    error_rel;        /* square root of relative error in computing function */
 16:   PetscReal    umin;             /* minimum allowable u'a value relative to |u|_1 */
 17:   PetscTruth   jorge;            /* flag indicating use of Jorge's method for determining
 18:                                    the differencing parameter */
 19:   PetscReal    h;                /* differencing parameter */
 20:   PetscTruth   need_h;           /* flag indicating whether we must compute h */
 21:   PetscTruth   need_err;         /* flag indicating whether we must currently compute error_rel */
 22:   PetscTruth   compute_err;      /* flag indicating whether we must ever compute error_rel */
 23:   PetscInt     compute_err_iter; /* last iter where we've computer error_rel */
 24:   PetscInt     compute_err_freq; /* frequency of computing error_rel */
 25:   void         *data;            /* implementation-specific data */
 26: } MFCtx_Private;

 30: PetscErrorCode SNESMatrixFreeDestroy2_Private(Mat mat)
 31: {
 33:   MFCtx_Private  *ctx;

 36:   MatShellGetContext(mat,(void **)&ctx);
 37:   VecDestroy(ctx->w);
 38:   if (ctx->sp) {MatNullSpaceDestroy(ctx->sp);}
 39:   if (ctx->jorge || ctx->compute_err) {DiffParameterDestroy_More(ctx->data);}
 40:   PetscFree(ctx);
 41:   return(0);
 42: }

 46: /*
 47:    SNESMatrixFreeView2_Private - Views matrix-free parameters.
 48:  */
 49: PetscErrorCode SNESMatrixFreeView2_Private(Mat J,PetscViewer viewer)
 50: {
 52:   MFCtx_Private  *ctx;
 53:   PetscTruth     iascii;

 56:   MatShellGetContext(J,(void **)&ctx);
 57:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
 58:   if (iascii) {
 59:      PetscViewerASCIIPrintf(viewer,"  SNES matrix-free approximation:\n");
 60:      if (ctx->jorge) {
 61:        PetscViewerASCIIPrintf(viewer,"    using Jorge's method of determining differencing parameter\n");
 62:      }
 63:      PetscViewerASCIIPrintf(viewer,"    err=%G (relative error in function evaluation)\n",ctx->error_rel);
 64:      PetscViewerASCIIPrintf(viewer,"    umin=%G (minimum iterate parameter)\n",ctx->umin);
 65:      if (ctx->compute_err) {
 66:        PetscViewerASCIIPrintf(viewer,"    freq_err=%D (frequency for computing err)\n",ctx->compute_err_freq);
 67:      }
 68:   } else {
 69:     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by SNES matrix free Jorge",((PetscObject)viewer)->type_name);
 70:   }
 71:   return(0);
 72: }

 76: /*
 77:   SNESMatrixFreeMult2_Private - Default matrix-free form for Jacobian-vector
 78:   product, y = F'(u)*a:
 79:         y = (F(u + ha) - F(u)) /h, 
 80:   where F = nonlinear function, as set by SNESSetFunction()
 81:         u = current iterate
 82:         h = difference interval
 83: */
 84: PetscErrorCode SNESMatrixFreeMult2_Private(Mat mat,Vec a,Vec y)
 85: {
 86:   MFCtx_Private  *ctx;
 87:   SNES           snes;
 88:   PetscReal      h,norm,sum,umin,noise;
 89:   PetscScalar    hs,dot;
 90:   Vec            w,U,F;
 91:   PetscErrorCode ierr,(*eval_fct)(SNES,Vec,Vec);
 92:   MPI_Comm       comm;
 93:   PetscInt       iter;


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

103:   PetscObjectGetComm((PetscObject)mat,&comm);
104:   MatShellGetContext(mat,(void **)&ctx);
105:   snes = ctx->snes;
106:   w    = ctx->w;
107:   umin = ctx->umin;

109:   SNESGetSolution(snes,&U);
110:   eval_fct = SNESComputeFunction;
111:   SNESGetFunction(snes,&F,PETSC_NULL,PETSC_NULL);

113:   /* Determine a "good" step size, h */
114:   if (ctx->need_h) {

116:     /* Use Jorge's method to compute h */
117:     if (ctx->jorge) {
118:       DiffParameterCompute_More(snes,ctx->data,U,a,&noise,&h);

120:     /* Use the Brown/Saad method to compute h */
121:     } else {
122:       /* Compute error if desired */
123:       SNESGetIterationNumber(snes,&iter);
124:       if ((ctx->need_err) || ((ctx->compute_err_freq) && (ctx->compute_err_iter != iter) && (!((iter-1)%ctx->compute_err_freq)))) {
125:         /* Use Jorge's method to compute noise */
126:         DiffParameterCompute_More(snes,ctx->data,U,a,&noise,&h);
127:         ctx->error_rel = sqrt(noise);
128:         PetscInfo3(snes,"Using Jorge's noise: noise=%G, sqrt(noise)=%G, h_more=%G\n",noise,ctx->error_rel,h);
129:         ctx->compute_err_iter = iter;
130:         ctx->need_err = PETSC_FALSE;
131:       }

133:       VecDotBegin(U,a,&dot);
134:       VecNormBegin(a,NORM_1,&sum);
135:       VecNormBegin(a,NORM_2,&norm);
136:       VecDotEnd(U,a,&dot);
137:       VecNormEnd(a,NORM_1,&sum);
138:       VecNormEnd(a,NORM_2,&norm);


141:       /* Safeguard for step sizes too small */
142:       if (sum == 0.0) {dot = 1.0; norm = 1.0;}
143: #if defined(PETSC_USE_COMPLEX)
144:       else if (PetscAbsScalar(dot) < umin*sum && PetscRealPart(dot) >= 0.0) dot = umin*sum;
145:       else if (PetscAbsScalar(dot) < 0.0 && PetscRealPart(dot) > -umin*sum) dot = -umin*sum;
146: #else
147:       else if (dot < umin*sum && dot >= 0.0) dot = umin*sum;
148:       else if (dot < 0.0 && dot > -umin*sum) dot = -umin*sum;
149: #endif
150:       h = PetscRealPart(ctx->error_rel*dot/(norm*norm));
151:     }
152:   } else {
153:     h = ctx->h;
154:   }
155:   if (!ctx->jorge || !ctx->need_h) {PetscInfo1(snes,"h = %G\n",h);}

157:   /* Evaluate function at F(u + ha) */
158:   hs = h;
159:   VecWAXPY(w,hs,a,U);
160:   eval_fct(snes,w,y);
161:   VecAXPY(y,-1.0,F);
162:   VecScale(y,1.0/hs);
163:   if (ctx->sp) {MatNullSpaceRemove(ctx->sp,y,PETSC_NULL);}

165:   PetscLogEventEnd(MATMFFD_Mult,a,y,0,0);
166:   return(0);
167: }

171: /*@C
172:    SNESMatrixFreeCreate2 - Creates a matrix-free matrix
173:    context for use with a SNES solver.  This matrix can be used as
174:    the Jacobian argument for the routine SNESSetJacobian().

176:    Input Parameters:
177: .  snes - the SNES context
178: .  x - vector where SNES solution is to be stored.

180:    Output Parameter:
181: .  J - the matrix-free matrix

183:    Level: advanced

185:    Notes:
186:    The matrix-free matrix context merely contains the function pointers
187:    and work space for performing finite difference approximations of
188:    Jacobian-vector products, J(u)*a, via

190: $       J(u)*a = [J(u+h*a) - J(u)]/h,
191: $   where by default
192: $        h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
193: $          = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
194: $   where
195: $        error_rel = square root of relative error in
196: $                    function evaluation
197: $        umin = minimum iterate parameter
198: $   Alternatively, the differencing parameter, h, can be set using
199: $   Jorge's nifty new strategy if one specifies the option 
200: $          -snes_mf_jorge

202:    The user can set these parameters via MatMFFDSetFunctionError().
203:    See the nonlinear solvers chapter of the users manual for details.

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

208:    Options Database Keys:
209: $  -snes_mf_err <error_rel>
210: $  -snes_mf_unim <umin>
211: $  -snes_mf_compute_err
212: $  -snes_mf_freq_err <freq>
213: $  -snes_mf_jorge

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

217: .seealso: MatDestroy(), MatMFFDSetFunctionError()
218: @*/
219: PetscErrorCode  SNESDefaultMatrixFreeCreate2(SNES snes,Vec x,Mat *J)
220: {
221:   MPI_Comm       comm;
222:   MFCtx_Private  *mfctx;
224:   PetscInt       n,nloc;
225:   PetscTruth     flg;
226:   char           p[64];

229:   PetscNewLog(snes,MFCtx_Private,&mfctx);
230:   mfctx->sp   = 0;
231:   mfctx->snes = snes;
232:   mfctx->error_rel        = PETSC_SQRT_MACHINE_EPSILON;
233:   mfctx->umin             = 1.e-6;
234:   mfctx->h                = 0.0;
235:   mfctx->need_h           = PETSC_TRUE;
236:   mfctx->need_err         = PETSC_FALSE;
237:   mfctx->compute_err      = PETSC_FALSE;
238:   mfctx->compute_err_freq = 0;
239:   mfctx->compute_err_iter = -1;
240:   mfctx->compute_err      = PETSC_FALSE;
241:   mfctx->jorge            = PETSC_FALSE;
242:   PetscOptionsGetReal(((PetscObject)snes)->prefix,"-snes_mf_err",&mfctx->error_rel,PETSC_NULL);
243:   PetscOptionsGetReal(((PetscObject)snes)->prefix,"-snes_mf_umin",&mfctx->umin,PETSC_NULL);
244:   PetscOptionsGetTruth(((PetscObject)snes)->prefix,"-snes_mf_jorge",&mfctx->jorge,PETSC_NULL);
245:   PetscOptionsGetTruth(((PetscObject)snes)->prefix,"-snes_mf_compute_err",&mfctx->compute_err,PETSC_NULL);
246:   PetscOptionsGetInt(((PetscObject)snes)->prefix,"-snes_mf_freq_err",&mfctx->compute_err_freq,&flg);
247:   if (flg) {
248:     if (mfctx->compute_err_freq < 0) mfctx->compute_err_freq = 0;
249:     mfctx->compute_err = PETSC_TRUE;
250:   }
251:   if (mfctx->compute_err) mfctx->need_err = PETSC_TRUE;
252:   if (mfctx->jorge || mfctx->compute_err) {
253:     DiffParameterCreate_More(snes,x,&mfctx->data);
254:   } else mfctx->data = 0;

256:   PetscOptionsHasName(PETSC_NULL,"-help",&flg);
257:   PetscStrcpy(p,"-");
258:   if (((PetscObject)snes)->prefix) PetscStrcat(p,((PetscObject)snes)->prefix);
259:   if (flg) {
260:     PetscPrintf(((PetscObject)snes)->comm," Matrix-free Options (via SNES):\n");
261:     PetscPrintf(((PetscObject)snes)->comm,"   %ssnes_mf_err <err>: set sqrt of relative error in function (default %G)\n",p,mfctx->error_rel);
262:     PetscPrintf(((PetscObject)snes)->comm,"   %ssnes_mf_umin <umin>: see users manual (default %G)\n",p,mfctx->umin);
263:     PetscPrintf(((PetscObject)snes)->comm,"   %ssnes_mf_jorge: use Jorge More's method\n",p);
264:     PetscPrintf(((PetscObject)snes)->comm,"   %ssnes_mf_compute_err: compute sqrt or relative error in function\n",p);
265:     PetscPrintf(((PetscObject)snes)->comm,"   %ssnes_mf_freq_err <freq>: frequency to recompute this (default only once)\n",p);
266:     PetscPrintf(((PetscObject)snes)->comm,"   %ssnes_mf_noise_file <file>: set file for printing noise info\n",p);
267:   }
268:   VecDuplicate(x,&mfctx->w);
269:   PetscObjectGetComm((PetscObject)x,&comm);
270:   VecGetSize(x,&n);
271:   VecGetLocalSize(x,&nloc);
272:   MatCreate(comm,J);
273:   MatSetSizes(*J,nloc,n,n,n);
274:   MatSetType(*J,MATSHELL);
275:   MatShellSetContext(*J,mfctx);
276:   MatShellSetOperation(*J,MATOP_MULT,(void(*)(void))SNESMatrixFreeMult2_Private);
277:   MatShellSetOperation(*J,MATOP_DESTROY,(void(*)(void))SNESMatrixFreeDestroy2_Private);
278:   MatShellSetOperation(*J,MATOP_VIEW,(void(*)(void))SNESMatrixFreeView2_Private);

280:   PetscLogObjectParent(*J,mfctx->w);
281:   PetscLogObjectParent(snes,*J);
282:   return(0);
283: }

287: /*@C
288:    SNESDefaultMatrixFreeSetParameters2 - Sets the parameters for the approximation of
289:    matrix-vector products using finite differences.

291: $       J(u)*a = [J(u+h*a) - J(u)]/h where

293:    either the user sets h directly here, or this parameter is computed via

295: $        h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
296: $          = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   else
297: $

299:    Input Parameters:
300: +  mat - the matrix
301: .  error_rel - relative error (should be set to the square root of
302:      the relative error in the function evaluations)
303: .  umin - minimum allowable u-value
304: -  h - differencing parameter

306:    Level: advanced

308:    Notes:
309:    If the user sets the parameter h directly, then this value will be used
310:    instead of the default computation indicated above.

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

314: .seealso: MatCreateSNESMF()
315: @*/
316: PetscErrorCode  SNESDefaultMatrixFreeSetParameters2(Mat mat,PetscReal error,PetscReal umin,PetscReal h)
317: {
318:   MFCtx_Private  *ctx;

322:   MatShellGetContext(mat,(void **)&ctx);
323:   if (ctx) {
324:     if (error != PETSC_DEFAULT) ctx->error_rel = error;
325:     if (umin  != PETSC_DEFAULT) ctx->umin = umin;
326:     if (h     != PETSC_DEFAULT) {
327:       ctx->h = h;
328:       ctx->need_h = PETSC_FALSE;
329:     }
330:   }
331:   return(0);
332: }

334: PetscErrorCode  SNESUnSetMatrixFreeParameter(SNES snes)
335: {
336:   MFCtx_Private  *ctx;
338:   Mat            mat;

341:   SNESGetJacobian(snes,&mat,PETSC_NULL,PETSC_NULL,PETSC_NULL);
342:   MatShellGetContext(mat,(void **)&ctx);
343:   if (ctx) ctx->need_h = PETSC_TRUE;
344:   return(0);
345: }
346: