Actual source code: snesmfjdef.c
1: /*$Id: snesmfjdef.c,v 1.31 2001/08/30 21:37:18 balay 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/mat/matimpl.h
35: #include src/snes/mf/snesmfj.h
37: /*
38: The default method has one parameter that is used to
39: "cutoff" very small values. This is stored in a data structure
40: that is only visible to this file. If your method has no parameters
41: it can omit this, if it has several simply reorganize the data structure.
42: The data structure is "hung-off" the MatSNESMF data structure in
43: the void *hctx; field.
44: */
45: typedef struct {
46: PetscReal umin; /* minimum allowable u'a value relative to |u|_1 */
47: } MatSNESMFDefault;
49: /*
50: MatSNESMFCompute_Default - Standard PETSc code for computing the
51: differencing paramter (h) for use with matrix-free finite differences.
53: Input Parameters:
54: + ctx - the matrix free context
55: . U - the location at which you want the Jacobian
56: - a - the direction you want the derivative
58: Output Parameter:
59: . h - the scale computed
61: */
62: static int MatSNESMFCompute_Default(MatSNESMFCtx ctx,Vec U,Vec a,PetscScalar *h)
63: {
64: MatSNESMFDefault *hctx = (MatSNESMFDefault*)ctx->hctx;
65: PetscReal nrm,sum,umin = hctx->umin;
66: PetscScalar dot;
67: int ierr;
70: if (!(ctx->count % ctx->recomputeperiod)) {
71: /*
72: This algorithm requires 2 norms and 1 inner product. Rather than
73: use directly the VecNorm() and VecDot() routines (and thus have
74: three separate collective operations, we use the VecxxxBegin/End() routines
75: */
76: VecDotBegin(U,a,&dot);
77: VecNormBegin(a,NORM_1,&sum);
78: VecNormBegin(a,NORM_2,&nrm);
79: VecDotEnd(U,a,&dot);
80: VecNormEnd(a,NORM_1,&sum);
81: VecNormEnd(a,NORM_2,&nrm);
83: /*
84: Safeguard for step sizes that are "too small"
85: */
86: if (!sum) {dot = 1.0; nrm = 1.0;}
87: #if defined(PETSC_USE_COMPLEX)
88: else if (PetscAbsScalar(dot) < umin*sum && PetscRealPart(dot) >= 0.0) dot = umin*sum;
89: else if (PetscAbsScalar(dot) < 0.0 && PetscRealPart(dot) > -umin*sum) dot = -umin*sum;
90: #else
91: else if (dot < umin*sum && dot >= 0.0) dot = umin*sum;
92: else if (dot < 0.0 && dot > -umin*sum) dot = -umin*sum;
93: #endif
94: *h = ctx->error_rel*dot/(nrm*nrm);
95: } else {
96: *h = ctx->currenth;
97: }
98: if (*h != *h) SETERRQ3(1,"Differencing parameter is not a number sum = %g dot = %g norm = %g",sum,PetscRealPart(dot),nrm);
99: ctx->count++;
100: return(0);
101: }
103: /*
104: MatSNESMFView_Default - Prints information about this particular
105: method for computing h. Note that this does not print the general
106: information about the matrix-free method, as such info is printed
107: by the calling routine.
109: Input Parameters:
110: + ctx - the matrix free context
111: - viewer - the PETSc viewer
112: */
113: static int MatSNESMFView_Default(MatSNESMFCtx ctx,PetscViewer viewer)
114: {
115: MatSNESMFDefault *hctx = (MatSNESMFDefault *)ctx->hctx;
116: int ierr;
117: PetscTruth isascii;
120: /*
121: Currently this only handles the ascii file viewers, others
122: could be added, but for this type of object other viewers
123: make less sense
124: */
125: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
126: if (isascii) {
127: PetscViewerASCIIPrintf(viewer," umin=%g (minimum iterate parameter)n",hctx->umin);
128: } else {
129: SETERRQ1(1,"Viewer type %s not supported for this SNES matrix free matrix",((PetscObject)viewer)->type_name);
130: }
131: return(0);
132: }
134: /*
135: MatSNESMFSetFromOptions_Default - Looks in the options database for
136: any options appropriate for this method.
138: Input Parameter:
139: . ctx - the matrix free context
141: */
142: static int MatSNESMFSetFromOptions_Default(MatSNESMFCtx ctx)
143: {
144: int ierr;
145: MatSNESMFDefault *hctx = (MatSNESMFDefault*)ctx->hctx;
148: PetscOptionsHead("Default matrix free parameters");
149: PetscOptionsReal("-snes_mf_umin","umin","MatSNESMFDefaultSetUmin",hctx->umin,&hctx->umin,0);
150: PetscOptionsTail();
151: return(0);
152: }
154: /*
155: MatSNESMFDestroy_Default - Frees the space allocated by
156: MatSNESMFCreate_Default().
158: Input Parameter:
159: . ctx - the matrix free context
161: Notes:
162: Does not free the ctx, that is handled by the calling routine
163: */
164: static int MatSNESMFDestroy_Default(MatSNESMFCtx ctx)
165: {
169: PetscFree(ctx->hctx);
170: return(0);
171: }
173: EXTERN_C_BEGIN
174: /*
175: The following two routines use the PetscObjectCompose() and PetscObjectQuery()
176: mechanism to allow the user to change the Umin parameter used in this method.
177: */
178: int MatSNESMFDefaultSetUmin_Private(Mat mat,PetscReal umin)
179: {
180: MatSNESMFCtx ctx = (MatSNESMFCtx)mat->data;
181: MatSNESMFDefault *hctx;
184: if (!ctx) {
185: SETERRQ(1,"MatSNESMFDefaultSetUmin() attached to non-shell matrix");
186: }
187: hctx = (MatSNESMFDefault*)ctx->hctx;
188: hctx->umin = umin;
189: return(0);
190: }
191: EXTERN_C_END
193: /*@
194: MatSNESMFDefaultSetUmin - Sets the "umin" parameter used by the default
195: PETSc routine for computing the differencing parameter, h, which is used
196: for matrix-free Jacobian-vector products.
198: Input Parameters:
199: + A - the matrix created with MatCreateSNESMF()
200: - umin - the parameter
202: Level: advanced
204: Notes:
205: See the manual page for MatCreateSNESMF() for a complete description of the
206: algorithm used to compute h.
208: .seealso: MatSNESMFSetFunctionError(), MatCreateSNESMF()
210: @*/
211: int MatSNESMFDefaultSetUmin(Mat A,PetscReal umin)
212: {
213: int ierr,(*f)(Mat,PetscReal);
217: PetscObjectQueryFunction((PetscObject)A,"MatSNESMFDefaultSetUmin_C",(void (**)(void))&f);
218: if (f) {
219: (*f)(A,umin);
220: }
221: return(0);
222: }
224: EXTERN_C_BEGIN
225: /*
226: MatSNESMFCreate_Default - Standard PETSc code for
227: computing h with matrix-free finite differences.
229: Input Parameter:
230: . ctx - the matrix free context created by MatSNESMFCreate()
232: */
233: int MatSNESMFCreate_Default(MatSNESMFCtx ctx)
234: {
235: MatSNESMFDefault *hctx;
236: int ierr;
240: /* allocate my own private data structure */
241: ierr = PetscNew(MatSNESMFDefault,&hctx);
242: ctx->hctx = (void*)hctx;
243: /* set a default for my parameter */
244: hctx->umin = 1.e-6;
246: /* set the functions I am providing */
247: ctx->ops->compute = MatSNESMFCompute_Default;
248: ctx->ops->destroy = MatSNESMFDestroy_Default;
249: ctx->ops->view = MatSNESMFView_Default;
250: ctx->ops->setfromoptions = MatSNESMFSetFromOptions_Default;
252: PetscObjectComposeFunctionDynamic((PetscObject)ctx->mat,"MatSNESMFDefaultSetUmin_C",
253: "MatSNESMFDefaultSetUmin_Private",
254: MatSNESMFDefaultSetUmin_Private);
255: return(0);
256: }
257: EXTERN_C_END