Actual source code: snesmfjdef.c
1: /*$Id: snesmfjdef.c,v 1.25 2001/04/10 19:36:54 bsmith Exp $*/
2: /*
3: Implements the default PETSc approach for computing the h
4: parameter used with the finite difference based matrix-free
5: Jacobian-vector products.
7: To make your own: clone this file and modify for your needs.
9: Mandatory functions:
10: -------------------
11: MatSNESMFCompute_ - for a given point and direction computes h
13: MatSNESMFCreate_ - fills in the MatSNESMF data structure
14: for this particular implementation
16:
17: Optional functions:
18: -------------------
19: MatSNESMFView_ - prints information about the parameters being used.
20: This is called when SNESView() or -snes_view is used.
22: MatSNESMFSetFromOptions_ - checks the options database for options that
23: apply to this method.
25: MatSNESMFDestroy_ - frees any space allocated by the routines above
27: */
29: /*
30: This include file defines the data structure MatSNESMF that
31: includes information about the computation of h. It is shared by
32: all implementations that people provide
33: */
34: #include "src/snes/mf/snesmfj.h" /*I "petscsnes.h" I*/
36: /*
37: The default method has one parameter that is used to
38: "cutoff" very small values. This is stored in a data structure
39: that is only visible to this file. If your method has no parameters
40: it can omit this, if it has several simply reorganize the data structure.
41: The data structure is "hung-off" the MatSNESMF data structure in
42: the void *hctx; field.
43: */
44: typedef struct {
45: PetscReal umin; /* minimum allowable u'a value relative to |u|_1 */
46: } MatSNESMFDefault;
48: /*
49: MatSNESMFCompute_Default - Standard PETSc code for computing the
50: differencing paramter (h) for use with matrix-free finite differences.
52: Input Parameters:
53: + ctx - the matrix free context
54: . U - the location at which you want the Jacobian
55: - a - the direction you want the derivative
57: Output Parameter:
58: . h - the scale computed
60: */
61: static int MatSNESMFCompute_Default(MatSNESMFCtx ctx,Vec U,Vec a,Scalar *h)
62: {
63: MatSNESMFDefault *hctx = (MatSNESMFDefault*)ctx->hctx;
64: PetscReal norm,sum,umin = hctx->umin;
65: Scalar dot;
66: int ierr;
69: if (!(ctx->count % ctx->recomputeperiod)) {
70: /*
71: This algorithm requires 2 norms and 1 inner product. Rather than
72: use directly the VecNorm() and VecDot() routines (and thus have
73: three separate collective operations, we use the VecxxxBegin/End() routines
74: */
75: VecDotBegin(U,a,&dot);
76: VecNormBegin(a,NORM_1,&sum);
77: VecNormBegin(a,NORM_2,&norm);
78: VecDotEnd(U,a,&dot);
79: VecNormEnd(a,NORM_1,&sum);
80: VecNormEnd(a,NORM_2,&norm);
82: /*
83: Safeguard for step sizes that are "too small"
84: */
85: if (!sum) {dot = 1.0; norm = 1.0;}
86: #if defined(PETSC_USE_COMPLEX)
87: else if (PetscAbsScalar(dot) < umin*sum && PetscRealPart(dot) >= 0.0) dot = umin*sum;
88: else if (PetscAbsScalar(dot) < 0.0 && PetscRealPart(dot) > -umin*sum) dot = -umin*sum;
89: #else
90: else if (dot < umin*sum && dot >= 0.0) dot = umin*sum;
91: else if (dot < 0.0 && dot > -umin*sum) dot = -umin*sum;
92: #endif
93: *h = ctx->error_rel*dot/(norm*norm);
94: } else {
95: *h = ctx->currenth;
96: }
97: if (*h != *h) SETERRQ(1,"Differencing parameter is not a number");
98: ctx->count++;
99: return(0);
100: }
102: /*
103: MatSNESMFView_Default - Prints information about this particular
104: method for computing h. Note that this does not print the general
105: information about the matrix-free method, as such info is printed
106: by the calling routine.
108: Input Parameters:
109: + ctx - the matrix free context
110: - viewer - the PETSc viewer
111: */
112: static int MatSNESMFView_Default(MatSNESMFCtx ctx,PetscViewer viewer)
113: {
114: MatSNESMFDefault *hctx = (MatSNESMFDefault *)ctx->hctx;
115: int ierr;
116: PetscTruth isascii;
119: /*
120: Currently this only handles the ascii file viewers, others
121: could be added, but for this type of object other viewers
122: make less sense
123: */
124: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
125: if (isascii) {
126: PetscViewerASCIIPrintf(viewer," umin=%g (minimum iterate parameter)n",hctx->umin);
127: } else {
128: SETERRQ1(1,"Viewer type %s not supported for this SNES matrix free matrix",((PetscObject)viewer)->type_name);
129: }
130: return(0);
131: }
133: /*
134: MatSNESMFSetFromOptions_Default - Looks in the options database for
135: any options appropriate for this method.
137: Input Parameter:
138: . ctx - the matrix free context
140: */
141: static int MatSNESMFSetFromOptions_Default(MatSNESMFCtx ctx)
142: {
143: int ierr;
144: MatSNESMFDefault *hctx = (MatSNESMFDefault*)ctx->hctx;
147: PetscOptionsHead("Default matrix free parameters");
148: PetscOptionsDouble("-snes_mf_umin","umin","MatSNESMFDefaultSetUmin",hctx->umin,&hctx->umin,0);
149: PetscOptionsTail();
150: return(0);
151: }
153: /*
154: MatSNESMFDestroy_Default - Frees the space allocated by
155: MatSNESMFCreate_Default().
157: Input Parameter:
158: . ctx - the matrix free context
160: Notes:
161: Does not free the ctx, that is handled by the calling routine
162: */
163: static int MatSNESMFDestroy_Default(MatSNESMFCtx ctx)
164: {
168: PetscFree(ctx->hctx);
169: return(0);
170: }
172: EXTERN_C_BEGIN
173: /*
174: The following two routines use the PetscObjectCompose() and PetscObjectQuery()
175: mechanism to allow the user to change the Umin parameter used in this method.
176: */
177: int MatSNESMFDefaultSetUmin_Private(Mat mat,PetscReal umin)
178: {
179: MatSNESMFCtx ctx;
180: MatSNESMFDefault *hctx;
181: int ierr;
184: MatShellGetContext(mat,(void **)&ctx);
185: if (!ctx) {
186: SETERRQ(1,"MatSNESMFDefaultSetUmin() attached to non-shell matrix");
187: }
188: hctx = (MatSNESMFDefault*)ctx->hctx;
189: hctx->umin = umin;
191: return(0);
192: }
193: EXTERN_C_END
195: /*@
196: MatSNESMFDefaultSetUmin - Sets the "umin" parameter used by the default
197: PETSc routine for computing the differencing parameter, h, which is used
198: for matrix-free Jacobian-vector products.
200: Input Parameters:
201: + A - the matrix created with MatCreateSNESMF()
202: - umin - the parameter
204: Level: advanced
206: Notes:
207: See the manual page for MatCreateSNESMF() for a complete description of the
208: algorithm used to compute h.
210: .seealso: MatSNESMFSetFunctionError(), MatCreateSNESMF()
212: @*/
213: int MatSNESMFDefaultSetUmin(Mat A,PetscReal umin)
214: {
215: int ierr,(*f)(Mat,PetscReal);
219: PetscObjectQueryFunction((PetscObject)A,"MatSNESMFDefaultSetUmin_C",(void (**)())&f);
220: if (f) {
221: (*f)(A,umin);
222: }
223: return(0);
224: }
226: EXTERN_C_BEGIN
227: /*
228: MatSNESMFCreate_Default - Standard PETSc code for
229: computing h with matrix-free finite differences.
231: Input Parameter:
232: . ctx - the matrix free context created by MatSNESMFCreate()
234: */
235: int MatSNESMFCreate_Default(MatSNESMFCtx ctx)
236: {
237: MatSNESMFDefault *hctx;
238: int ierr;
242: /* allocate my own private data structure */
243: ierr = PetscNew(MatSNESMFDefault,&hctx);
244: ctx->hctx = (void*)hctx;
245: /* set a default for my parameter */
246: hctx->umin = 1.e-6;
248: /* set the functions I am providing */
249: ctx->ops->compute = MatSNESMFCompute_Default;
250: ctx->ops->destroy = MatSNESMFDestroy_Default;
251: ctx->ops->view = MatSNESMFView_Default;
252: ctx->ops->setfromoptions = MatSNESMFSetFromOptions_Default;
254: PetscObjectComposeFunctionDynamic((PetscObject)ctx->mat,"MatSNESMFDefaultSetUmin_C",
255: "MatSNESMFDefaultSetUmin_Private",
256: MatSNESMFDefaultSetUmin_Private);
257: return(0);
258: }
259: EXTERN_C_END