Actual source code: snesmfj2.c
1: /*$Id: snesmfj2.c,v 1.35 2001/08/21 21:03:55 bsmith Exp $*/
3: #include src/snes/snesimpl.h
5: EXTERN int DiffParameterCreate_More(SNES,Vec,void**);
6: EXTERN int DiffParameterCompute_More(SNES,void*,Vec,Vec,PetscReal*,PetscReal*);
7: EXTERN int DiffParameterDestroy_More(void*);
9: typedef struct { /* default context for matrix-free SNES */
10: SNES snes; /* SNES context */
11: Vec w; /* work vector */
12: MatNullSpace sp; /* null space context */
13: PetscReal error_rel; /* square root of relative error in computing function */
14: PetscReal umin; /* minimum allowable u'a value relative to |u|_1 */
15: int jorge; /* flag indicating use of Jorge's method for determining
16: the differencing parameter */
17: PetscReal h; /* differencing parameter */
18: int need_h; /* flag indicating whether we must compute h */
19: int need_err; /* flag indicating whether we must currently compute error_rel */
20: int compute_err; /* flag indicating whether we must ever compute error_rel */
21: int compute_err_iter; /* last iter where we've computer error_rel */
22: int compute_err_freq; /* frequency of computing error_rel */
23: void *data; /* implementation-specific data */
24: } MFCtx_Private;
26: int SNESMatrixFreeDestroy2_Private(Mat mat)
27: {
28: int ierr;
29: MFCtx_Private *ctx;
32: MatShellGetContext(mat,(void **)&ctx);
33: VecDestroy(ctx->w);
34: if (ctx->sp) {MatNullSpaceDestroy(ctx->sp);}
35: if (ctx->jorge || ctx->compute_err) {DiffParameterDestroy_More(ctx->data);}
36: PetscFree(ctx);
37: return(0);
38: }
40: /*
41: SNESMatrixFreeView2_Private - Views matrix-free parameters.
42: */
43: int SNESMatrixFreeView2_Private(Mat J,PetscViewer viewer)
44: {
45: int ierr;
46: MFCtx_Private *ctx;
47: PetscTruth isascii;
50: MatShellGetContext(J,(void **)&ctx);
51: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
52: if (isascii) {
53: PetscViewerASCIIPrintf(viewer," SNES matrix-free approximation:n");
54: if (ctx->jorge) {
55: PetscViewerASCIIPrintf(viewer," using Jorge's method of determining differencing parametern");
56: }
57: PetscViewerASCIIPrintf(viewer," err=%g (relative error in function evaluation)n",ctx->error_rel);
58: PetscViewerASCIIPrintf(viewer," umin=%g (minimum iterate parameter)n",ctx->umin);
59: if (ctx->compute_err) {
60: PetscViewerASCIIPrintf(viewer," freq_err=%d (frequency for computing err)n",ctx->compute_err_freq);
61: }
62: } else {
63: SETERRQ1(1,"Viewer type %s not supported by SNES matrix free Jorge",((PetscObject)viewer)->type_name);
64: }
65: return(0);
66: }
68: /*
69: SNESMatrixFreeMult2_Private - Default matrix-free form for Jacobian-vector
70: product, y = F'(u)*a:
71: y = (F(u + ha) - F(u)) /h,
72: where F = nonlinear function, as set by SNESSetFunction()
73: u = current iterate
74: h = difference interval
75: */
76: int SNESMatrixFreeMult2_Private(Mat mat,Vec a,Vec y)
77: {
78: MFCtx_Private *ctx;
79: SNES snes;
80: PetscReal h,norm,sum,umin,noise;
81: PetscScalar hs,dot,mone = -1.0;
82: Vec w,U,F;
83: int ierr,iter,(*eval_fct)(SNES,Vec,Vec);
84: MPI_Comm comm;
88: /* We log matrix-free matrix-vector products separately, so that we can
89: separate the performance monitoring from the cases that use conventional
90: storage. We may eventually modify event logging to associate events
91: with particular objects, hence alleviating the more general problem. */
92: PetscLogEventBegin(MAT_MultMatrixFree,a,y,0,0);
94: PetscObjectGetComm((PetscObject)mat,&comm);
95: MatShellGetContext(mat,(void **)&ctx);
96: snes = ctx->snes;
97: w = ctx->w;
98: umin = ctx->umin;
100: SNESGetSolution(snes,&U);
101: eval_fct = SNESComputeFunction;
102: SNESGetFunction(snes,&F,PETSC_NULL,PETSC_NULL);
104: /* Determine a "good" step size, h */
105: if (ctx->need_h) {
107: /* Use Jorge's method to compute h */
108: if (ctx->jorge) {
109: DiffParameterCompute_More(snes,ctx->data,U,a,&noise,&h);
111: /* Use the Brown/Saad method to compute h */
112: } else {
113: /* Compute error if desired */
114: SNESGetIterationNumber(snes,&iter);
115: if ((ctx->need_err) || ((ctx->compute_err_freq) && (ctx->compute_err_iter != iter) && (!((iter-1)%ctx->compute_err_freq)))) {
116: /* Use Jorge's method to compute noise */
117: DiffParameterCompute_More(snes,ctx->data,U,a,&noise,&h);
118: ctx->error_rel = sqrt(noise);
119: PetscLogInfo(snes,"SNESMatrixFreeMult2_Private: Using Jorge's noise: noise=%g, sqrt(noise)=%g, h_more=%gn",noise,ctx->error_rel,h);
120: ctx->compute_err_iter = iter;
121: ctx->need_err = 0;
122: }
124: VecDotBegin(U,a,&dot);
125: VecNormBegin(a,NORM_1,&sum);
126: VecNormBegin(a,NORM_2,&norm);
127: VecDotEnd(U,a,&dot);
128: VecNormEnd(a,NORM_1,&sum);
129: VecNormEnd(a,NORM_2,&norm);
132: /* Safeguard for step sizes too small */
133: if (sum == 0.0) {dot = 1.0; norm = 1.0;}
134: #if defined(PETSC_USE_COMPLEX)
135: else if (PetscAbsScalar(dot) < umin*sum && PetscRealPart(dot) >= 0.0) dot = umin*sum;
136: else if (PetscAbsScalar(dot) < 0.0 && PetscRealPart(dot) > -umin*sum) dot = -umin*sum;
137: #else
138: else if (dot < umin*sum && dot >= 0.0) dot = umin*sum;
139: else if (dot < 0.0 && dot > -umin*sum) dot = -umin*sum;
140: #endif
141: h = PetscRealPart(ctx->error_rel*dot/(norm*norm));
142: }
143: } else {
144: h = ctx->h;
145: }
146: if (!ctx->jorge || !ctx->need_h) PetscLogInfo(snes,"SNESMatrixFreeMult2_Private: h = %gn",h);
148: /* Evaluate function at F(u + ha) */
149: hs = h;
150: VecWAXPY(&hs,a,U,w);
151: eval_fct(snes,w,y);
152: VecAXPY(&mone,F,y);
153: hs = 1.0/hs;
154: VecScale(&hs,y);
155: if (ctx->sp) {MatNullSpaceRemove(ctx->sp,y,PETSC_NULL);}
157: PetscLogEventEnd(MAT_MultMatrixFree,a,y,0,0);
158: return(0);
159: }
161: /*@C
162: SNESMatrixFreeMatCreate2 - Creates a matrix-free matrix
163: context for use with a SNES solver. This matrix can be used as
164: the Jacobian argument for the routine SNESSetJacobian().
166: Input Parameters:
167: . snes - the SNES context
168: . x - vector where SNES solution is to be stored.
170: Output Parameter:
171: . J - the matrix-free matrix
173: Level: advanced
175: Notes:
176: The matrix-free matrix context merely contains the function pointers
177: and work space for performing finite difference approximations of
178: Jacobian-vector products, J(u)*a, via
180: $ J(u)*a = [J(u+h*a) - J(u)]/h,
181: $ where by default
182: $ h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1}
183: $ = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 otherwise
184: $ where
185: $ error_rel = square root of relative error in
186: $ function evaluation
187: $ umin = minimum iterate parameter
188: $ Alternatively, the differencing parameter, h, can be set using
189: $ Jorge's nifty new strategy if one specifies the option
190: $ -snes_mf_jorge
192: The user can set these parameters via MatSNESMFSetFunctionError().
193: See the nonlinear solvers chapter of the users manual for details.
195: The user should call MatDestroy() when finished with the matrix-free
196: matrix context.
198: Options Database Keys:
199: $ -snes_mf_err <error_rel>
200: $ -snes_mf_unim <umin>
201: $ -snes_mf_compute_err
202: $ -snes_mf_freq_err <freq>
203: $ -snes_mf_jorge
205: .keywords: SNES, default, matrix-free, create, matrix
207: .seealso: MatDestroy(), MatSNESMFSetFunctionError()
208: @*/
209: int SNESDefaultMatrixFreeCreate2(SNES snes,Vec x,Mat *J)
210: {
211: MPI_Comm comm;
212: MFCtx_Private *mfctx;
213: int n,nloc,ierr;
214: PetscTruth flg;
215: char p[64];
218: PetscNew(MFCtx_Private,&mfctx);
219: ierr = PetscMemzero(mfctx,sizeof(MFCtx_Private));
220: PetscLogObjectMemory(snes,sizeof(MFCtx_Private));
221: mfctx->sp = 0;
222: mfctx->snes = snes;
223: mfctx->error_rel = PETSC_SQRT_MACHINE_EPSILON;
224: mfctx->umin = 1.e-6;
225: mfctx->h = 0.0;
226: mfctx->need_h = 1;
227: mfctx->need_err = 0;
228: mfctx->compute_err = 0;
229: mfctx->compute_err_freq = 0;
230: mfctx->compute_err_iter = -1;
231: PetscOptionsGetReal(snes->prefix,"-snes_mf_err",&mfctx->error_rel,PETSC_NULL);
232: PetscOptionsGetReal(snes->prefix,"-snes_mf_umin",&mfctx->umin,PETSC_NULL);
233: PetscOptionsHasName(snes->prefix,"-snes_mf_jorge",(PetscTruth*)&mfctx->jorge);
234: PetscOptionsHasName(snes->prefix,"-snes_mf_compute_err",(PetscTruth*)&mfctx->compute_err);
235: PetscOptionsGetInt(snes->prefix,"-snes_mf_freq_err",&mfctx->compute_err_freq,&flg);
236: if (flg) {
237: if (mfctx->compute_err_freq < 0) mfctx->compute_err_freq = 0;
238: mfctx->compute_err = 1;
239: }
240: if (mfctx->compute_err == 1) mfctx->need_err = 1;
241: if (mfctx->jorge || mfctx->compute_err) {
242: DiffParameterCreate_More(snes,x,&mfctx->data);
243: } else mfctx->data = 0;
245: PetscOptionsHasName(PETSC_NULL,"-help",&flg);
246: PetscStrcpy(p,"-");
247: if (snes->prefix) PetscStrcat(p,snes->prefix);
248: if (flg) {
249: PetscPrintf(snes->comm," Matrix-free Options (via SNES):n");
250: PetscPrintf(snes->comm," %ssnes_mf_err <err>: set sqrt of relative error in function (default %g)n",p,mfctx->error_rel);
251: PetscPrintf(snes->comm," %ssnes_mf_umin <umin>: see users manual (default %g)n",p,mfctx->umin);
252: PetscPrintf(snes->comm," %ssnes_mf_jorge: use Jorge More's methodn",p);
253: PetscPrintf(snes->comm," %ssnes_mf_compute_err: compute sqrt or relative error in functionn",p);
254: PetscPrintf(snes->comm," %ssnes_mf_freq_err <freq>: frequency to recompute this (default only once)n",p);
255: PetscPrintf(snes->comm," %ssnes_mf_noise_file <file>: set file for printing noise infon",p);
256: }
257: VecDuplicate(x,&mfctx->w);
258: PetscObjectGetComm((PetscObject)x,&comm);
259: VecGetSize(x,&n);
260: VecGetLocalSize(x,&nloc);
261: MatCreateShell(comm,nloc,n,n,n,mfctx,J);
262: MatShellSetOperation(*J,MATOP_MULT,(void(*)(void))SNESMatrixFreeMult2_Private);
263: MatShellSetOperation(*J,MATOP_DESTROY,(void(*)(void))SNESMatrixFreeDestroy2_Private);
264: MatShellSetOperation(*J,MATOP_VIEW,(void(*)(void))SNESMatrixFreeView2_Private);
266: PetscLogObjectParent(*J,mfctx->w);
267: PetscLogObjectParent(snes,*J);
268: return(0);
269: }
271: /*@C
272: SNESDefaultMatrixFreeSetParameters2 - Sets the parameters for the approximation of
273: matrix-vector products using finite differences.
275: $ J(u)*a = [J(u+h*a) - J(u)]/h where
277: either the user sets h directly here, or this parameter is computed via
279: $ h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1}
280: $ = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 else
281: $
283: Input Parameters:
284: + mat - the matrix
285: . error_rel - relative error (should be set to the square root of
286: the relative error in the function evaluations)
287: . umin - minimum allowable u-value
288: - h - differencing parameter
290: Level: advanced
292: Notes:
293: If the user sets the parameter h directly, then this value will be used
294: instead of the default computation indicated above.
296: .keywords: SNES, matrix-free, parameters
298: .seealso: MatCreateSNESMF()
299: @*/
300: int SNESDefaultMatrixFreeSetParameters2(Mat mat,PetscReal error,PetscReal umin,PetscReal h)
301: {
302: MFCtx_Private *ctx;
303: int ierr;
306: MatShellGetContext(mat,(void **)&ctx);
307: if (ctx) {
308: if (error != PETSC_DEFAULT) ctx->error_rel = error;
309: if (umin != PETSC_DEFAULT) ctx->umin = umin;
310: if (h != PETSC_DEFAULT) {
311: ctx->h = h;
312: ctx->need_h = 0;
313: }
314: }
315: return(0);
316: }
318: int SNESUnSetMatrixFreeParameter(SNES snes)
319: {
320: MFCtx_Private *ctx;
321: int ierr;
322: Mat mat;
325: SNESGetJacobian(snes,&mat,PETSC_NULL,PETSC_NULL,PETSC_NULL);
326: MatShellGetContext(mat,(void **)&ctx);
327: if (ctx) ctx->need_h = 1;
328: return(0);
329: }
330: