Actual source code: wp.c
1: #define PETSCMAT_DLL
3: /*MC
4: MATMFFD_WP - Implements an alternative approach for computing the differencing parameter
5: h used with the finite difference based matrix-free Jacobian. This code
6: implements the strategy of M. Pernice and H. Walker:
8: h = error_rel * sqrt(1 + ||U||) / ||a||
10: Notes:
11: 1) || U || does not change between linear iterations so is reused
12: 2) In GMRES || a || == 1 and so does not need to ever be computed except at restart
13: when it is recomputed.
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: Options Database Keys:
20: . -mat_mffd_compute_normu -Compute the norm of u everytime see MatMFFDWPSetComputeNormU()
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
30: .seealso: MATMFFD, MatCreateMFFD(), MatCreateSNESMF(), MATMFFD_DS
32: M*/
34: /*
35: This include file defines the data structure MatMFFD that
36: includes information about the computation of h. It is shared by
37: all implementations that people provide.
39: See snesmfjdef.c for a full set of comments on the routines below.
40: */
41: #include private/matimpl.h
42: #include ../src/mat/impls/mffd/mffdimpl.h
44: typedef struct {
45: PetscReal normUfact; /* previous sqrt(1.0 + || U ||) */
46: PetscTruth computenorma,computenormU;
47: } MatMFFD_WP;
51: /*
52: MatMFFDCompute_WP - Standard PETSc code for
53: computing h with matrix-free finite differences.
55: Input Parameters:
56: + ctx - the matrix free context
57: . U - the location at which you want the Jacobian
58: - a - the direction you want the derivative
60: Output Parameter:
61: . h - the scale computed
63: */
64: static PetscErrorCode MatMFFDCompute_WP(MatMFFD ctx,Vec U,Vec a,PetscScalar *h,PetscTruth *zeroa)
65: {
66: MatMFFD_WP *hctx = (MatMFFD_WP*)ctx->hctx;
67: PetscReal normU,norma = 1.0;
71: if (!(ctx->count % ctx->recomputeperiod)) {
72: if (hctx->computenorma && (hctx->computenormU || !ctx->ncurrenth)) {
73: VecNormBegin(U,NORM_2,&normU);
74: VecNormBegin(a,NORM_2,&norma);
75: VecNormEnd(U,NORM_2,&normU);
76: VecNormEnd(a,NORM_2,&norma);
77: hctx->normUfact = sqrt(1.0+normU);
78: } else if (hctx->computenormU || !ctx->ncurrenth) {
79: VecNorm(U,NORM_2,&normU);
80: hctx->normUfact = sqrt(1.0+normU);
81: } else if (hctx->computenorma) {
82: VecNorm(a,NORM_2,&norma);
83: }
84: if (norma == 0.0) {
85: *zeroa = PETSC_TRUE;
86: return(0);
87: }
88: *zeroa = PETSC_FALSE;
89: *h = ctx->error_rel*hctx->normUfact/norma;
90: } else {
91: *h = ctx->currenth;
92: }
93: return(0);
94: }
98: /*
99: MatMFFDView_WP - Prints information about this particular
100: method for computing h. Note that this does not print the general
101: information about the matrix free, that is printed by the calling
102: routine.
104: Input Parameters:
105: + ctx - the matrix free context
106: - viewer - the PETSc viewer
108: */
109: static PetscErrorCode MatMFFDView_WP(MatMFFD ctx,PetscViewer viewer)
110: {
111: MatMFFD_WP *hctx = (MatMFFD_WP *)ctx->hctx;
113: PetscTruth iascii;
116: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
117: if (iascii) {
118: if (hctx->computenorma){PetscViewerASCIIPrintf(viewer," Computes normA\n");}
119: else { PetscViewerASCIIPrintf(viewer," Does not compute normA\n");}
120: if (hctx->computenormU){ PetscViewerASCIIPrintf(viewer," Computes normU\n");}
121: else { PetscViewerASCIIPrintf(viewer," Does not compute normU\n");}
122: } else {
123: SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for SNES matrix-free WP",((PetscObject)viewer)->type_name);
124: }
125: return(0);
126: }
130: /*
131: MatMFFDSetFromOptions_WP - Looks in the options database for
132: any options appropriate for this method
134: Input Parameter:
135: . ctx - the matrix free context
137: */
138: static PetscErrorCode MatMFFDSetFromOptions_WP(MatMFFD ctx)
139: {
141: MatMFFD_WP *hctx = (MatMFFD_WP*)ctx->hctx;
144: PetscOptionsHead("Walker-Pernice options");
145: PetscOptionsTruth("-mat_mffd_compute_normu","Compute the norm of u","MatMFFDWPSetComputeNormU",
146: hctx->computenorma,&hctx->computenorma,0);
147: PetscOptionsTail();
148: return(0);
149: }
153: /*
154: MatMFFDDestroy_WP - Frees the space allocated by
155: MatMFFDCreate_WP().
157: Input Parameter:
158: . ctx - the matrix free context
160: Notes: does not free the ctx, that is handled by the calling routine
162: */
163: static PetscErrorCode MatMFFDDestroy_WP(MatMFFD ctx)
164: {
167: PetscFree(ctx->hctx);
168: return(0);
169: }
174: PetscErrorCode MatMFFDWPSetComputeNormU_P(Mat mat,PetscTruth flag)
175: {
176: MatMFFD ctx = (MatMFFD)mat->data;
177: MatMFFD_WP *hctx = (MatMFFD_WP*)ctx->hctx;
180: hctx->computenormU = flag;
181: return(0);
182: }
187: /*@
188: MatMFFDWPSetComputeNormU - Sets whether it computes the ||U|| used by the WP
189: PETSc routine for computing h. With any Krylov solver this need only
190: be computed during the first iteration and kept for later.
192: Input Parameters:
193: + A - the matrix created with MatCreateSNESMF()
194: - flag - PETSC_TRUE causes it to compute ||U||, PETSC_FALSE uses the previous value
196: Options Database Key:
197: . -mat_mffd_compute_normu <true,false> - true by default, false can save calculations but you
198: must be sure that ||U|| has not changed in the mean time.
200: Level: advanced
202: Notes:
203: See the manual page for MATMFFD_WP for a complete description of the
204: algorithm used to compute h.
206: .seealso: MatMFFDSetFunctionError(), MatCreateSNESMF()
208: @*/
209: PetscErrorCode MatMFFDWPSetComputeNormU(Mat A,PetscTruth flag)
210: {
211: PetscErrorCode ierr,(*f)(Mat,PetscTruth);
215: PetscObjectQueryFunction((PetscObject)A,"MatMFFDWPSetComputeNormU_C",(void (**)(void))&f);
216: if (f) {
217: (*f)(A,flag);
218: }
219: return(0);
220: }
225: /*
226: MatMFFDCreate_WP - Standard PETSc code for
227: computing h with matrix-free finite differences.
229: Input Parameter:
230: . ctx - the matrix free context created by MatMFFDCreate()
232: */
233: PetscErrorCode MatMFFDCreate_WP(MatMFFD ctx)
234: {
236: MatMFFD_WP *hctx;
240: /* allocate my own private data structure */
241: PetscNewLog(ctx,MatMFFD_WP,&hctx);
242: ctx->hctx = (void*)hctx;
243: hctx->computenormU = PETSC_FALSE;
244: hctx->computenorma = PETSC_TRUE;
246: /* set the functions I am providing */
247: ctx->ops->compute = MatMFFDCompute_WP;
248: ctx->ops->destroy = MatMFFDDestroy_WP;
249: ctx->ops->view = MatMFFDView_WP;
250: ctx->ops->setfromoptions = MatMFFDSetFromOptions_WP;
252: PetscObjectComposeFunctionDynamic((PetscObject)ctx->mat,"MatMFFDWPSetComputeNormU_C",
253: "MatMFFDWPSetComputeNormU_P",
254: MatMFFDWPSetComputeNormU_P);
256: return(0);
257: }