xref: /petsc/src/snes/interface/noise/snesdnest.c (revision 8bc8193efbc389280f83b3d41dffa9e2d23e2ace)
1 /* fnoise/snesdnest.F -- translated by f2c (version 20020314).
2 */
3 #include "petsc.h"
4 #define FALSE_ 0
5 #define TRUE_ 1
6 
7 /*  Noise estimation routine, written by Jorge More'.  Details are below. */
8 
9 /* Subroutine */ PetscErrorCode dnest_(PetscInt *nf, double *fval,double *h__,double *fnoise, double *fder2, double *hopt, PetscInt *info, double *eps)
10 {
11     /* Initialized data */
12 
13     static double const__[15] = { .71,.41,.23,.12,.063,.033,.018,.0089,
14 	    .0046,.0024,.0012,6.1e-4,3.1e-4,1.6e-4,8e-5 };
15 
16     /* System generated locals */
17     PetscInt i__1;
18     double d__1, d__2, d__3, d__4;
19 
20 
21     /* Local variables */
22     static double emin, emax;
23     static PetscInt dsgn[6];
24     static double f_max, f_min, stdv;
25     static PetscInt i__, j;
26     static double scale;
27     static PetscInt mh;
28     static PetscInt cancel[6], dnoise;
29     static double err2, est1, est2, est3, est4;
30 
31 /*     ********** */
32 
33 /*     Subroutine dnest */
34 
35 /*     This subroutine estimates the noise in a function */
36 /*     and provides estimates of the optimal difference parameter */
37 /*     for a forward-difference approximation. */
38 
39 /*     The user must provide a difference parameter h, and the */
40 /*     function value at nf points centered around the current point. */
41 /*     For example, if nf = 7, the user must provide */
42 
43 /*        f(x-2*h), f(x-h), f(x), f(x+h),  f(x+2*h), */
44 
45 /*     in the array fval. The use of nf = 7 function evaluations is */
46 /*     recommended. */
47 
48 /*     The noise in the function is roughly defined as the variance in */
49 /*     the computed value of the function. The noise in the function */
50 /*     provides valuable information. For example, function values */
51 /*     smaller than the noise should be considered to be zero. */
52 
53 /*     This subroutine requires an initial estimate for h. Under estimates */
54 /*     are usually preferred. If noise is not detected, the user should */
55 /*     increase or decrease h according to the ouput value of info. */
56 /*     In most cases, the subroutine detects noise with the initial */
57 /*     value of h. */
58 
59 /*     The subroutine statement is */
60 
61 /*       subroutine dnest(nf,fval,h,hopt,fnoise,info,eps) */
62 
63 /*     where */
64 
65 /*       nf is an int variable. */
66 /*         On entry nf is the number of function values. */
67 /*         On exit nf is unchanged. */
68 
69 /*       f is a double precision array of dimension nf. */
70 /*         On entry f contains the function values. */
71 /*         On exit f is overwritten. */
72 
73 /*       h is a double precision variable. */
74 /*         On entry h is an estimate of the optimal difference parameter. */
75 /*         On exit h is unchanged. */
76 
77 /*       fnoise is a double precision variable. */
78 /*         On entry fnoise need not be specified. */
79 /*         On exit fnoise is set to an estimate of the function noise */
80 /*            if noise is detected; otherwise fnoise is set to zero. */
81 
82 /*       hopt is a double precision variable. */
83 /*         On entry hopt need not be specified. */
84 /*         On exit hopt is set to an estimate of the optimal difference */
85 /*            parameter if noise is detected; otherwise hopt is set to zero. */
86 
87 /*       info is an int variable. */
88 /*         On entry info need not be specified. */
89 /*         On exit info is set as follows: */
90 
91 /*            info = 1  Noise has been detected. */
92 
93 /*            info = 2  Noise has not been detected; h is too small. */
94 /*                      Try 100*h for the next value of h. */
95 
96 /*            info = 3  Noise has not been detected; h is too large. */
97 /*                      Try h/100 for the next value of h. */
98 
99 /*            info = 4  Noise has been detected but the estimate of hopt */
100 /*                      is not reliable; h is too small. */
101 
102 /*       eps is a double precision work array of dimension nf. */
103 
104 /*     MINPACK-2 Project. April 1997. */
105 /*     Argonne National Laboratory. */
106 /*     Jorge J. More'. */
107 
108 /*     ********** */
109     /* Parameter adjustments */
110     --eps;
111     --fval;
112 
113     /* Function Body */
114     *fnoise = 0.;
115     *fder2 = 0.;
116     *hopt = 0.;
117 /*     Compute an estimate of the second derivative and */
118 /*     determine a bound on the error. */
119     mh = (*nf + 1) / 2;
120     est1 = (fval[mh + 1] - fval[mh] * 2 + fval[mh - 1]) / *h__ / *h__;
121     est2 = (fval[mh + 2] - fval[mh] * 2 + fval[mh - 2]) / (*h__ * 2) / (*h__ *
122 	     2);
123     est3 = (fval[mh + 3] - fval[mh] * 2 + fval[mh - 3]) / (*h__ * 3) / (*h__ *
124 	     3);
125     est4 = (est1 + est2 + est3) / 3;
126 /* Computing MAX */
127 /* Computing PETSCMAX */
128     d__3 = PetscMax(est1,est2);
129 /* Computing MIN */
130     d__4 = PetscMin(est1,est2);
131     d__1 = PetscMax(d__3,est3) - est4, d__2 = est4 - PetscMin(d__4,est3);
132     err2 = PetscMax(d__1,d__2);
133 /*      write (2,123) est1, est2, est3 */
134 /* 123  format ('Second derivative estimates', 3d12.2) */
135     if (err2 <= PetscAbsScalar(est4) * .1) {
136 	*fder2 = est4;
137     } else if (err2 < PetscAbsScalar(est4)) {
138 	*fder2 = est3;
139     } else {
140 	*fder2 = 0.;
141     }
142 /*     Compute the range of function values. */
143     f_min = fval[1];
144     f_max = fval[1];
145     i__1 = *nf;
146     for (i__ = 2; i__ <= i__1; ++i__) {
147 /* Computing MIN */
148 	d__1 = f_min, d__2 = fval[i__];
149 	f_min = PetscMin(d__1,d__2);
150 /* Computing MAX */
151 	d__1 = f_max, d__2 = fval[i__];
152 	f_max = PetscMax(d__1,d__2);
153     }
154 /*     Construct the difference table. */
155     dnoise = FALSE_;
156     for (j = 1; j <= 6; ++j) {
157 	dsgn[j - 1] = FALSE_;
158 	cancel[j - 1] = FALSE_;
159 	scale = 0.;
160 	i__1 = *nf - j;
161 	for (i__ = 1; i__ <= i__1; ++i__) {
162 	    fval[i__] = fval[i__ + 1] - fval[i__];
163 	    if (fval[i__] == 0.) {
164 		cancel[j - 1] = TRUE_;
165 	    }
166 /* Computing MAX */
167 	    d__2 = scale, d__3 = (d__1 = fval[i__], PetscAbsScalar(d__1));
168 	    scale = PetscMax(d__2,d__3);
169 	}
170 /*        Compute the estimates for the noise level. */
171 	if (scale == 0.) {
172 	    stdv = 0.;
173 	} else {
174 	    stdv = 0.;
175 	    i__1 = *nf - j;
176 	    for (i__ = 1; i__ <= i__1; ++i__) {
177 /* Computing 2nd power */
178 		d__1 = fval[i__] / scale;
179 		stdv += d__1 * d__1;
180 	    }
181 	    stdv = scale * PetscSqrtScalar(stdv / (*nf - j));
182 	}
183 	eps[j] = const__[j - 1] * stdv;
184 /*        Determine differences in sign. */
185 	i__1 = *nf - j - 1;
186 	for (i__ = 1; i__ <= i__1; ++i__) {
187 /* Computing MIN */
188 	    d__1 = fval[i__], d__2 = fval[i__ + 1];
189 /* Computing MAX */
190 	    d__3 = fval[i__], d__4 = fval[i__ + 1];
191 	    if (PetscMin(d__1,d__2) < 0. && PetscMax(d__3,d__4) > 0.) {
192 		dsgn[j - 1] = TRUE_;
193 	    }
194 	}
195     }
196 /*     First requirement for detection of noise. */
197     dnoise = dsgn[3];
198 /*     Check for h too small or too large. */
199     *info = 0;
200     if (f_max == f_min) {
201 	*info = 2;
202     } else /* if(complicated condition) */ {
203 /* Computing MIN */
204 	d__1 = PetscAbsScalar(f_max), d__2 = PetscAbsScalar(f_min);
205 	if (f_max - f_min > PetscMin(d__1,d__2) * .1) {
206 	    *info = 3;
207 	}
208     }
209     if (*info != 0) {
210 	PetscFunctionReturn(0);
211     }
212 /*     Determine the noise level. */
213 /* Computing MIN */
214     d__1 = PetscMin(eps[4],eps[5]);
215     emin = PetscMin(d__1,eps[6]);
216 /* Computing MAX */
217     d__1 = PetscMax(eps[4],eps[5]);
218     emax = PetscMax(d__1,eps[6]);
219     if (emax <= emin * 4 && dnoise) {
220 	*fnoise = (eps[4] + eps[5] + eps[6]) / 3;
221 	if (*fder2 != 0.) {
222 	    *info = 1;
223 	    *hopt = PetscSqrtScalar(*fnoise / PetscAbsScalar(*fder2)) * 1.68;
224 	} else {
225 	    *info = 4;
226 	    *hopt = *h__ * 10;
227 	}
228 	PetscFunctionReturn(0);
229     }
230 /* Computing MIN */
231     d__1 = PetscMin(eps[3],eps[4]);
232     emin = PetscMin(d__1,eps[5]);
233 /* Computing MAX */
234     d__1 = PetscMax(eps[3],eps[4]);
235     emax = PetscMax(d__1,eps[5]);
236     if (emax <= emin * 4 && dnoise) {
237 	*fnoise = (eps[3] + eps[4] + eps[5]) / 3;
238 	if (*fder2 != 0.) {
239 	    *info = 1;
240 	    *hopt = PetscSqrtScalar(*fnoise / PetscAbsScalar(*fder2)) * 1.68;
241 	} else {
242 	    *info = 4;
243 	    *hopt = *h__ * 10;
244 	}
245 	PetscFunctionReturn(0);
246     }
247 /*     Noise not detected; decide if h is too small or too large. */
248     if (! cancel[3]) {
249 	if (dsgn[3]) {
250 	    *info = 2;
251 	} else {
252 	    *info = 3;
253 	}
254 	PetscFunctionReturn(0);
255     }
256     if (! cancel[2]) {
257 	if (dsgn[2]) {
258 	    *info = 2;
259 	} else {
260 	    *info = 3;
261 	}
262 	PetscFunctionReturn(0);
263     }
264 /*     If there is cancelllation on the third and fourth column */
265 /*     then h is too small */
266     *info = 2;
267     PetscFunctionReturn(0);
268 /*      if (cancel .or. dsgn(3)) then */
269 /*         info = 2 */
270 /*      else */
271 /*         info = 3 */
272 /*      end if */
273 } /* dnest_ */
274 
275