Actual source code: wp.c
1: /*MC
2: MATSNESMF_WP - Implements an alternative approach for computing the differencing parameter
3: h used with the finite difference based matrix-free Jacobian. This code
4: implements the strategy of M. Pernice and H. Walker:
6: h = error_rel * sqrt(1 + ||U||) / ||a||
8: Notes:
9: 1) || U || does not change between linear iterations so can be reused
10: 2) In GMRES || a || == 1 and so does not need to ever be computed if you never
11: have a restart. Unfortunately a RESTART computes a matrix vector product
12: with ||a|| != 0 which breaks this
14: Reference: M. Pernice and H. F. Walker, "NITSOL: A Newton Iterative
15: Solver for Nonlinear Systems", SIAM J. Sci. Stat. Comput.", 1998,
16: vol 19, pp. 302--318.
18: Options Database Keys:
19: + -snes_mf_compute_norma - compute the norm of a everytime see MatSNESMFWPSetComputeNormA()
20: - -snes_mf_compute_normu -Compute the norm of u everytime see MatSNESMFWPSetComputeNormU()
23: Level: intermediate
25: Notes: Requires no global collectives when used with GMRES
27: Formula used:
28: F'(u)*a = [F(u+h*a) - F(u)]/h where
29: h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1}
30: = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 otherwise
31: where
32: error_rel = square root of relative error in function evaluation
33: umin = minimum iterate parameter
35: .seealso: MATMFFD, MatCreateMF(), MatCreateSNESMF(), MATSNESMF_DEFAULT
37: M*/
39: /*
40: This include file defines the data structure MatSNESMF that
41: includes information about the computation of h. It is shared by
42: all implementations that people provide.
44: See snesmfjdef.c for a full set of comments on the routines below.
45: */
46: #include src/mat/matimpl.h
47: #include src/snes/mf/snesmfj.h
49: typedef struct {
50: PetscReal normUfact; /* previous sqrt(1.0 + || U ||) */
51: PetscTruth computenorma,computenormU;
52: } MatSNESMFWP;
56: /*
57: MatSNESMFCompute_WP - Standard PETSc code for
58: computing h with matrix-free finite differences.
60: Input Parameters:
61: + ctx - the matrix free context
62: . U - the location at which you want the Jacobian
63: - a - the direction you want the derivative
65: Output Parameter:
66: . h - the scale computed
68: */
69: static PetscErrorCode MatSNESMFCompute_WP(MatSNESMFCtx ctx,Vec U,Vec a,PetscScalar *h)
70: {
71: MatSNESMFWP *hctx = (MatSNESMFWP*)ctx->hctx;
72: PetscReal normU,norma = 1.0;
77: if (!(ctx->count % ctx->recomputeperiod)) {
78: if (hctx->computenorma && (hctx->computenormU || !ctx->ncurrenth)) {
79: VecNormBegin(U,NORM_2,&normU);
80: VecNormBegin(a,NORM_2,&norma);
81: VecNormEnd(U,NORM_2,&normU);
82: VecNormEnd(a,NORM_2,&norma);
83: hctx->normUfact = sqrt(1.0+normU);
84: } else if (hctx->computenormU || !ctx->ncurrenth) {
85: VecNorm(U,NORM_2,&normU);
86: hctx->normUfact = sqrt(1.0+normU);
87: } else if (hctx->computenorma) {
88: VecNorm(a,NORM_2,&norma);
89: }
90: *h = ctx->error_rel*hctx->normUfact/norma;
91: } else {
92: *h = ctx->currenth;
93: }
94: return(0);
95: }
99: /*
100: MatSNESMFView_WP - Prints information about this particular
101: method for computing h. Note that this does not print the general
102: information about the matrix free, that is printed by the calling
103: routine.
105: Input Parameters:
106: + ctx - the matrix free context
107: - viewer - the PETSc viewer
109: */
110: static PetscErrorCode MatSNESMFView_WP(MatSNESMFCtx ctx,PetscViewer viewer)
111: {
112: MatSNESMFWP *hctx = (MatSNESMFWP *)ctx->hctx;
114: PetscTruth iascii;
117: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
118: if (iascii) {
119: if (hctx->computenorma){PetscViewerASCIIPrintf(viewer," Computes normA\n");}
120: else { PetscViewerASCIIPrintf(viewer," Does not compute normA\n");}
121: if (hctx->computenormU){ PetscViewerASCIIPrintf(viewer," Computes normU\n");}
122: else { PetscViewerASCIIPrintf(viewer," Does not compute normU\n");}
123: } else {
124: SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for SNES matrix-free WP",((PetscObject)viewer)->type_name);
125: }
126: return(0);
127: }
131: /*
132: MatSNESMFSetFromOptions_WP - Looks in the options database for
133: any options appropriate for this method
135: Input Parameter:
136: . ctx - the matrix free context
138: */
139: static PetscErrorCode MatSNESMFSetFromOptions_WP(MatSNESMFCtx ctx)
140: {
142: MatSNESMFWP *hctx = (MatSNESMFWP*)ctx->hctx;
145: PetscOptionsHead("Walker-Pernice options");
146: PetscOptionsLogical("-snes_mf_compute_norma","Compute the norm of a","MatSNESMFWPSetComputeNormA",
147: hctx->computenorma,&hctx->computenorma,0);
148: PetscOptionsLogical("-snes_mf_compute_normu","Compute the norm of u","MatSNESMFWPSetComputeNormU",
149: hctx->computenorma,&hctx->computenorma,0);
150: PetscOptionsTail();
151: return(0);
152: }
156: /*
157: MatSNESMFDestroy_WP - Frees the space allocated by
158: MatSNESMFCreate_WP().
160: Input Parameter:
161: . ctx - the matrix free context
163: Notes: does not free the ctx, that is handled by the calling routine
165: */
166: static PetscErrorCode MatSNESMFDestroy_WP(MatSNESMFCtx ctx)
167: {
170: PetscFree(ctx->hctx);
171: return(0);
172: }
177: PetscErrorCode MatSNESMFWPSetComputeNormA_P(Mat mat,PetscTruth flag)
178: {
179: MatSNESMFCtx ctx = (MatSNESMFCtx)mat->data;
180: MatSNESMFWP *hctx;
183: hctx = (MatSNESMFWP*)ctx->hctx;
184: hctx->computenorma = flag;
185: return(0);
186: }
191: /*@
192: MatSNESMFWPSetComputeNormA - Sets whether it computes the ||a|| used by the WP
193: PETSc routine for computing h. With GMRES since the ||a|| is always
194: one, you can save communication by setting this to false.
196: Input Parameters:
197: + A - the matrix created with MatCreateSNESMF()
198: - flag - PETSC_TRUE causes it to compute ||a||, PETSC_FALSE assumes it is 1.
200: Options Database Key:
201: . -snes_mf_compute_norma <true,false> - can be set false for GMRES to speed up code
203: Level: advanced
205: Notes:
206: See the manual page for MatCreateSNESMF() for a complete description of the
207: algorithm used to compute h.
209: .seealso: MatSNESMFSetFunctionError(), MatCreateSNESMF()
211: @*/
212: PetscErrorCode MatSNESMFWPSetComputeNormA(Mat A,PetscTruth flag)
213: {
214: PetscErrorCode ierr,(*f)(Mat,PetscTruth);
218: PetscObjectQueryFunction((PetscObject)A,"MatSNESMFWPSetComputeNormA_C",(void (**)(void))&f);
219: if (f) {
220: (*f)(A,flag);
221: }
222: return(0);
223: }
228: PetscErrorCode MatSNESMFWPSetComputeNormU_P(Mat mat,PetscTruth flag)
229: {
230: MatSNESMFCtx ctx = (MatSNESMFCtx)mat->data;
231: MatSNESMFWP *hctx;
234: hctx = (MatSNESMFWP*)ctx->hctx;
235: hctx->computenormU = flag;
236: return(0);
237: }
242: /*@
243: MatSNESMFWPSetComputeNormU - Sets whether it computes the ||U|| used by the WP
244: PETSc routine for computing h. With any Krylov solver this need only
245: be computed during the first iteration and kept for later.
247: Input Parameters:
248: + A - the matrix created with MatCreateSNESMF()
249: - flag - PETSC_TRUE causes it to compute ||U||, PETSC_FALSE uses the previous value
251: Options Database Key:
252: . -snes_mf_compute_normu <true,false> - true by default, false causes extra calculations
254: Level: advanced
256: Notes:
257: See the manual page for MatCreateSNESMF() for a complete description of the
258: algorithm used to compute h.
260: .seealso: MatSNESMFSetFunctionError(), MatCreateSNESMF()
262: @*/
263: PetscErrorCode MatSNESMFWPSetComputeNormU(Mat A,PetscTruth flag)
264: {
265: PetscErrorCode ierr,(*f)(Mat,PetscTruth);
269: PetscObjectQueryFunction((PetscObject)A,"MatSNESMFWPSetComputeNormU_C",(void (**)(void))&f);
270: if (f) {
271: (*f)(A,flag);
272: }
273: return(0);
274: }
279: /*
280: MatSNESMFCreate_WP - Standard PETSc code for
281: computing h with matrix-free finite differences.
283: Input Parameter:
284: . ctx - the matrix free context created by MatSNESMFCreate()
286: */
287: PetscErrorCode MatSNESMFCreate_WP(MatSNESMFCtx ctx)
288: {
290: MatSNESMFWP *hctx;
294: /* allocate my own private data structure */
295: PetscNew(MatSNESMFWP,&hctx);
296: ctx->hctx = (void*)hctx;
297: hctx->computenormU = PETSC_FALSE;
298: hctx->computenorma = PETSC_TRUE;
300: /* set the functions I am providing */
301: ctx->ops->compute = MatSNESMFCompute_WP;
302: ctx->ops->destroy = MatSNESMFDestroy_WP;
303: ctx->ops->view = MatSNESMFView_WP;
304: ctx->ops->setfromoptions = MatSNESMFSetFromOptions_WP;
306: PetscObjectComposeFunctionDynamic((PetscObject)ctx->mat,"MatSNESMFWPSetComputeNormA_C",
307: "MatSNESMFWPSetComputeNormA_P",
308: MatSNESMFWPSetComputeNormA_P);
309: PetscObjectComposeFunctionDynamic((PetscObject)ctx->mat,"MatSNESMFWPSetComputeNormU_C",
310: "MatSNESMFWPSetComputeNormU_P",
311: MatSNESMFWPSetComputeNormU_P);
313: return(0);
314: }