Actual source code: wp.c
1: /*$Id: wp.c,v 1.36 2001/08/06 21:17:13 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/mat/matimpl.h
28: #include src/snes/mf/snesmfj.h
30: typedef struct {
31: PetscReal normUfact; /* previous sqrt(1.0 + || U ||) */
32: PetscTruth computenorma,computenormU;
33: } MatSNESMFWP;
35: /*
36: MatSNESMFCompute_WP - Standard PETSc code for
37: computing h with matrix-free finite differences.
39: Input Parameters:
40: + ctx - the matrix free context
41: . U - the location at which you want the Jacobian
42: - a - the direction you want the derivative
44: Output Parameter:
45: . h - the scale computed
47: */
48: static int MatSNESMFCompute_WP(MatSNESMFCtx ctx,Vec U,Vec a,PetscScalar *h)
49: {
50: MatSNESMFWP *hctx = (MatSNESMFWP*)ctx->hctx;
51: PetscReal normU,norma = 1.0;
52: int ierr;
56: if (!(ctx->count % ctx->recomputeperiod)) {
57: if (hctx->computenorma && (hctx->computenormU || !ctx->ncurrenth)) {
58: VecNormBegin(U,NORM_2,&normU);
59: VecNormBegin(a,NORM_2,&norma);
60: VecNormEnd(U,NORM_2,&normU);
61: VecNormEnd(a,NORM_2,&norma);
62: hctx->normUfact = sqrt(1.0+normU);
63: } else if (hctx->computenormU || !ctx->ncurrenth) {
64: VecNorm(U,NORM_2,&normU);
65: hctx->normUfact = sqrt(1.0+normU);
66: } else if (hctx->computenorma) {
67: VecNorm(a,NORM_2,&norma);
68: }
69: *h = ctx->error_rel*hctx->normUfact/norma;
70: } else {
71: *h = ctx->currenth;
72: }
73: return(0);
74: }
76: /*
77: MatSNESMFView_WP - Prints information about this particular
78: method for computing h. Note that this does not print the general
79: information about the matrix free, that is printed by the calling
80: routine.
82: Input Parameters:
83: + ctx - the matrix free context
84: - viewer - the PETSc viewer
86: */
87: static int MatSNESMFView_WP(MatSNESMFCtx ctx,PetscViewer viewer)
88: {
89: MatSNESMFWP *hctx = (MatSNESMFWP *)ctx->hctx;
90: int ierr;
91: PetscTruth isascii;
94: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
95: if (isascii) {
96: if (hctx->computenorma){PetscViewerASCIIPrintf(viewer," Computes normAn");}
97: else { PetscViewerASCIIPrintf(viewer," Does not compute normAn");}
98: if (hctx->computenormU){ PetscViewerASCIIPrintf(viewer," Computes normUn");}
99: else { PetscViewerASCIIPrintf(viewer," Does not compute normUn");}
100: } else {
101: SETERRQ1(1,"Viewer type %s not supported for SNES matrix-free WP",((PetscObject)viewer)->type_name);
102: }
103: return(0);
104: }
106: /*
107: MatSNESMFSetFromOptions_WP - Looks in the options database for
108: any options appropriate for this method
110: Input Parameter:
111: . ctx - the matrix free context
113: */
114: static int MatSNESMFSetFromOptions_WP(MatSNESMFCtx ctx)
115: {
116: int ierr;
117: MatSNESMFWP *hctx = (MatSNESMFWP*)ctx->hctx;
120: PetscOptionsHead("Walker-Pernice options");
121: PetscOptionsLogical("-snes_mf_compute_norma","Compute the norm of a","MatSNESMFWPSetComputeNormA",
122: hctx->computenorma,&hctx->computenorma,0);
123: PetscOptionsLogical("-snes_mf_compute_normu","Compute the norm of u","MatSNESMFWPSetComputeNormU",
124: hctx->computenorma,&hctx->computenorma,0);
125: PetscOptionsTail();
126: return(0);
127: }
129: /*
130: MatSNESMFDestroy_WP - Frees the space allocated by
131: MatSNESMFCreate_WP().
133: Input Parameter:
134: . ctx - the matrix free context
136: Notes: does not free the ctx, that is handled by the calling routine
138: */
139: static int MatSNESMFDestroy_WP(MatSNESMFCtx ctx)
140: {
143: PetscFree(ctx->hctx);
144: return(0);
145: }
147: EXTERN_C_BEGIN
148: int MatSNESMFWPSetComputeNormA_P(Mat mat,PetscTruth flag)
149: {
150: MatSNESMFCtx ctx = (MatSNESMFCtx)mat->data;
151: MatSNESMFWP *hctx;
154: hctx = (MatSNESMFWP*)ctx->hctx;
155: hctx->computenorma = flag;
156: return(0);
157: }
158: EXTERN_C_END
160: /*@
161: MatSNESMFWPSetComputeNormA - Sets whether it computes the ||a|| used by the WP
162: PETSc routine for computing h. With GMRES since the ||a|| is always
163: one, you can save communication by setting this to false.
165: Input Parameters:
166: + A - the matrix created with MatCreateSNESMF()
167: - flag - PETSC_TRUE causes it to compute ||a||, PETSC_FALSE assumes it is 1.
169: Level: advanced
171: Notes:
172: See the manual page for MatCreateSNESMF() for a complete description of the
173: algorithm used to compute h.
175: .seealso: MatSNESMFSetFunctionError(), MatCreateSNESMF()
177: @*/
178: int MatSNESMFWPSetComputeNormA(Mat A,PetscTruth flag)
179: {
180: int ierr,(*f)(Mat,PetscTruth);
184: PetscObjectQueryFunction((PetscObject)A,"MatSNESMFWPSetComputeNormA_C",(void (**)(void))&f);
185: if (f) {
186: (*f)(A,flag);
187: }
188: return(0);
189: }
191: EXTERN_C_BEGIN
192: int MatSNESMFWPSetComputeNormU_P(Mat mat,PetscTruth flag)
193: {
194: MatSNESMFCtx ctx = (MatSNESMFCtx)mat->data;
195: MatSNESMFWP *hctx;
198: hctx = (MatSNESMFWP*)ctx->hctx;
199: hctx->computenormU = flag;
200: return(0);
201: }
202: EXTERN_C_END
204: /*@
205: MatSNESMFWPSetComputeNormU - Sets whether it computes the ||U|| used by the WP
206: PETSc routine for computing h. With any Krylov solver this need only
207: be computed during the first iteration and kept for later.
209: Input Parameters:
210: + A - the matrix created with MatCreateSNESMF()
211: - flag - PETSC_TRUE causes it to compute ||U||, PETSC_FALSE uses the previous value
213: Level: advanced
215: Notes:
216: See the manual page for MatCreateSNESMF() for a complete description of the
217: algorithm used to compute h.
219: .seealso: MatSNESMFSetFunctionError(), MatCreateSNESMF()
221: @*/
222: int MatSNESMFWPSetComputeNormU(Mat A,PetscTruth flag)
223: {
224: int ierr,(*f)(Mat,PetscTruth);
228: PetscObjectQueryFunction((PetscObject)A,"MatSNESMFWPSetComputeNormU_C",(void (**)(void))&f);
229: if (f) {
230: (*f)(A,flag);
231: }
232: return(0);
233: }
235: EXTERN_C_BEGIN
236: /*
237: MatSNESMFCreate_WP - Standard PETSc code for
238: computing h with matrix-free finite differences.
240: Input Parameter:
241: . ctx - the matrix free context created by MatSNESMFCreate()
243: */
244: int MatSNESMFCreate_WP(MatSNESMFCtx ctx)
245: {
246: int ierr;
247: MatSNESMFWP *hctx;
251: /* allocate my own private data structure */
252: ierr = PetscNew(MatSNESMFWP,&hctx);
253: ctx->hctx = (void*)hctx;
254: hctx->computenormU = PETSC_FALSE;
255: hctx->computenorma = PETSC_TRUE;
257: /* set the functions I am providing */
258: ctx->ops->compute = MatSNESMFCompute_WP;
259: ctx->ops->destroy = MatSNESMFDestroy_WP;
260: ctx->ops->view = MatSNESMFView_WP;
261: ctx->ops->setfromoptions = MatSNESMFSetFromOptions_WP;
263: PetscObjectComposeFunctionDynamic((PetscObject)ctx->mat,"MatSNESMFWPSetComputeNormA_C",
264: "MatSNESMFWPSetComputeNormA_P",
265: MatSNESMFWPSetComputeNormA_P);
266: PetscObjectComposeFunctionDynamic((PetscObject)ctx->mat,"MatSNESMFWPSetComputeNormU_C",
267: "MatSNESMFWPSetComputeNormU_P",
268: MatSNESMFWPSetComputeNormU_P);
270: return(0);
271: }
272: EXTERN_C_END