Actual source code: wp.c
1: /*$Id: wp.c,v 1.33 2001/04/10 19:36:54 bsmith Exp $*/
2: /*
3: Implements an alternative approach for computing the differencing parameter
4: h used with the finite difference based matrix-free Jacobian. This code
5: implements the strategy of M. Pernice and H. Walker:
7: h = error_rel * sqrt(1 + ||U||) / ||a||
9: Notes:
10: 1) || U || does not change between linear iterations so can be reused
11: 2) In GMRES || a || == 1 and so does not need to ever be computed if you never
12: have a restart. Unfortunately a RESTART computes a matrix vector product
13: with ||a|| != 0 which breaks this
15: Reference: M. Pernice and H. F. Walker, "NITSOL: A Newton Iterative
16: Solver for Nonlinear Systems", SIAM J. Sci. Stat. Comput.", 1998,
17: vol 19, pp. 302--318.
19: See snesmfjdef.c for a full set of comments on the routines below.
20: */
22: /*
23: This include file defines the data structure MatSNESMF that
24: includes information about the computation of h. It is shared by
25: all implementations that people provide.
26: */
27: #include src/snes/mf/snesmfj.h
29: typedef struct {
30: PetscReal normUfact; /* previous sqrt(1.0 + || U ||) */
31: PetscTruth computenorma,computenormU;
32: } MatSNESMFWP;
34: /*
35: MatSNESMFCompute_WP - Standard PETSc code for
36: computing h with matrix-free finite differences.
38: Input Parameters:
39: + ctx - the matrix free context
40: . U - the location at which you want the Jacobian
41: - a - the direction you want the derivative
43: Output Parameter:
44: . h - the scale computed
46: */
47: static int MatSNESMFCompute_WP(MatSNESMFCtx ctx,Vec U,Vec a,Scalar *h)
48: {
49: MatSNESMFWP *hctx = (MatSNESMFWP*)ctx->hctx;
50: PetscReal normU,norma = 1.0;
51: int ierr;
55: if (!(ctx->count % ctx->recomputeperiod)) {
56: if (hctx->computenorma && (hctx->computenormU || !ctx->ncurrenth)) {
57: VecNormBegin(U,NORM_2,&normU);
58: VecNormBegin(a,NORM_2,&norma);
59: VecNormEnd(U,NORM_2,&normU);
60: VecNormEnd(a,NORM_2,&norma);
61: hctx->normUfact = sqrt(1.0+normU);
62: } else if (hctx->computenormU || !ctx->ncurrenth) {
63: VecNorm(U,NORM_2,&normU);
64: hctx->normUfact = sqrt(1.0+normU);
65: } else if (hctx->computenorma) {
66: VecNorm(a,NORM_2,&norma);
67: }
68: *h = ctx->error_rel*hctx->normUfact/norma;
69: } else {
70: *h = ctx->currenth;
71: }
72: return(0);
73: }
75: /*
76: MatSNESMFView_WP - Prints information about this particular
77: method for computing h. Note that this does not print the general
78: information about the matrix free, that is printed by the calling
79: routine.
81: Input Parameters:
82: + ctx - the matrix free context
83: - viewer - the PETSc viewer
85: */
86: static int MatSNESMFView_WP(MatSNESMFCtx ctx,PetscViewer viewer)
87: {
88: MatSNESMFWP *hctx = (MatSNESMFWP *)ctx->hctx;
89: int ierr;
90: PetscTruth isascii;
93: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
94: if (isascii) {
95: if (hctx->computenorma){PetscViewerASCIIPrintf(viewer," Computes normAn");}
96: else { PetscViewerASCIIPrintf(viewer," Does not compute normAn");}
97: if (hctx->computenormU){ PetscViewerASCIIPrintf(viewer," Computes normUn");}
98: else { PetscViewerASCIIPrintf(viewer," Does not compute normUn");}
99: } else {
100: SETERRQ1(1,"Viewer type %s not supported for SNES matrix-free WP",((PetscObject)viewer)->type_name);
101: }
102: return(0);
103: }
105: /*
106: MatSNESMFSetFromOptions_WP - Looks in the options database for
107: any options appropriate for this method
109: Input Parameter:
110: . ctx - the matrix free context
112: */
113: static int MatSNESMFSetFromOptions_WP(MatSNESMFCtx ctx)
114: {
115: int ierr;
116: MatSNESMFWP *hctx = (MatSNESMFWP*)ctx->hctx;
119: PetscOptionsHead("Walker-Pernice options");
120: PetscOptionsLogical("-snes_mf_compute_norma","Compute the norm of a","MatSNESMFWPSetComputeNormA",
121: hctx->computenorma,&hctx->computenorma,0);
122: PetscOptionsLogical("-snes_mf_compute_normu","Compute the norm of u","MatSNESMFWPSetComputeNormU",
123: hctx->computenorma,&hctx->computenorma,0);
124: PetscOptionsTail();
125: return(0);
126: }
128: /*
129: MatSNESMFDestroy_WP - Frees the space allocated by
130: MatSNESMFCreate_WP().
132: Input Parameter:
133: . ctx - the matrix free context
135: Notes: does not free the ctx, that is handled by the calling routine
137: */
138: static int MatSNESMFDestroy_WP(MatSNESMFCtx ctx)
139: {
142: PetscFree(ctx->hctx);
143: return(0);
144: }
146: EXTERN_C_BEGIN
147: int MatSNESMFWPSetComputeNormA_P(Mat mat,PetscTruth flag)
148: {
149: MatSNESMFCtx ctx;
150: MatSNESMFWP *hctx;
151: int ierr;
154: MatShellGetContext(mat,(void **)&ctx);
155: if (!ctx) {
156: SETERRQ(1,"MatSNESMFWPSetComputeNormA() attached to non-shell matrix");
157: }
158: hctx = (MatSNESMFWP*)ctx->hctx;
159: hctx->computenorma = flag;
161: return(0);
162: }
163: EXTERN_C_END
165: /*@
166: MatSNESMFWPSetComputeNormA - Sets whether it computes the ||a|| used by the WP
167: PETSc routine for computing h. With GMRES since the ||a|| is always
168: one, you can save communication by setting this to false.
170: Input Parameters:
171: + A - the matrix created with MatCreateSNESMF()
172: - flag - PETSC_TRUE causes it to compute ||a||, PETSC_FALSE assumes it is 1.
174: Level: advanced
176: Notes:
177: See the manual page for MatCreateSNESMF() for a complete description of the
178: algorithm used to compute h.
180: .seealso: MatSNESMFSetFunctionError(), MatCreateSNESMF()
182: @*/
183: int MatSNESMFWPSetComputeNormA(Mat A,PetscTruth flag)
184: {
185: int ierr,(*f)(Mat,PetscTruth);
189: PetscObjectQueryFunction((PetscObject)A,"MatSNESMFWPSetComputeNormA_C",(void (**)())&f);
190: if (f) {
191: (*f)(A,flag);
192: }
193: return(0);
194: }
196: EXTERN_C_BEGIN
197: int MatSNESMFWPSetComputeNormU_P(Mat mat,PetscTruth flag)
198: {
199: MatSNESMFCtx ctx;
200: MatSNESMFWP *hctx;
201: int ierr;
204: MatShellGetContext(mat,(void **)&ctx);
205: if (!ctx) {
206: SETERRQ(1,"MatSNESMFWPSetComputeNormU() attached to non-shell matrix");
207: }
208: hctx = (MatSNESMFWP*)ctx->hctx;
209: hctx->computenormU = flag;
211: return(0);
212: }
213: EXTERN_C_END
215: /*@
216: MatSNESMFWPSetComputeNormU - Sets whether it computes the ||U|| used by the WP
217: PETSc routine for computing h. With any Krylov solver this need only
218: be computed during the first iteration and kept for later.
220: Input Parameters:
221: + A - the matrix created with MatCreateSNESMF()
222: - flag - PETSC_TRUE causes it to compute ||U||, PETSC_FALSE uses the previous value
224: Level: advanced
226: Notes:
227: See the manual page for MatCreateSNESMF() for a complete description of the
228: algorithm used to compute h.
230: .seealso: MatSNESMFSetFunctionError(), MatCreateSNESMF()
232: @*/
233: int MatSNESMFWPSetComputeNormU(Mat A,PetscTruth flag)
234: {
235: int ierr,(*f)(Mat,PetscTruth);
239: PetscObjectQueryFunction((PetscObject)A,"MatSNESMFWPSetComputeNormU_C",(void (**)())&f);
240: if (f) {
241: (*f)(A,flag);
242: }
243: return(0);
244: }
246: EXTERN_C_BEGIN
247: /*
248: MatSNESMFCreate_WP - Standard PETSc code for
249: computing h with matrix-free finite differences.
251: Input Parameter:
252: . ctx - the matrix free context created by MatSNESMFCreate()
254: */
255: int MatSNESMFCreate_WP(MatSNESMFCtx ctx)
256: {
257: int ierr;
258: MatSNESMFWP *hctx;
262: /* allocate my own private data structure */
263: ierr = PetscNew(MatSNESMFWP,&hctx);
264: ctx->hctx = (void*)hctx;
265: hctx->computenormU = PETSC_FALSE;
266: hctx->computenorma = PETSC_TRUE;
268: /* set the functions I am providing */
269: ctx->ops->compute = MatSNESMFCompute_WP;
270: ctx->ops->destroy = MatSNESMFDestroy_WP;
271: ctx->ops->view = MatSNESMFView_WP;
272: ctx->ops->setfromoptions = MatSNESMFSetFromOptions_WP;
274: PetscObjectComposeFunctionDynamic((PetscObject)ctx->mat,"MatSNESMFWPSetComputeNormA_C",
275: "MatSNESMFWPSetComputeNormA_P",
276: MatSNESMFWPSetComputeNormA_P);
277: PetscObjectComposeFunctionDynamic((PetscObject)ctx->mat,"MatSNESMFWPSetComputeNormU_C",
278: "MatSNESMFWPSetComputeNormU_P",
279: MatSNESMFWPSetComputeNormU_P);
281: return(0);
282: }
283: EXTERN_C_END