Actual source code: snesdnest.c

petsc-3.6.1 2015-08-06
Report Typos and Errors
  2: /* fnoise/snesdnest.F -- translated by f2c (version 20020314).
  3: */
  4: #include <petscsys.h>
  5: #define FALSE_ 0
  6: #define TRUE_ 1

  8: /*  Noise estimation routine, written by Jorge More'.  Details are below. */

 10: /* Subroutine */ PetscErrorCode SNESNoise_dnest_(PetscInt *nf, double *fval,double *h__,double *fnoise, double *fder2, double *hopt, PetscInt *info, double *eps)
 11: {
 12:   /* Initialized data */

 14:   static double const__[15] = { .71,.41,.23,.12,.063,.033,.018,.0089,
 15:                                 .0046,.0024,.0012,6.1e-4,3.1e-4,1.6e-4,8e-5 };

 17:   /* System generated locals */
 18:   PetscInt i__1;
 19:   double   d__1, d__2, d__3, d__4;


 22:   /* Local variables */
 23:   static double   emin, emax;
 24:   static PetscInt dsgn[6];
 25:   static double   f_max, f_min, stdv;
 26:   static PetscInt i__, j;
 27:   static double   scale;
 28:   static PetscInt mh;
 29:   static PetscInt cancel[6], dnoise;
 30:   static double   err2, est1, est2, est3, est4;

 32: /*     ********** */

 34: /*     Subroutine dnest */

 36: /*     This subroutine estimates the noise in a function */
 37: /*     and provides estimates of the optimal difference parameter */
 38: /*     for a forward-difference approximation. */

 40: /*     The user must provide a difference parameter h, and the */
 41: /*     function value at nf points centered around the current point. */
 42: /*     For example, if nf = 7, the user must provide */

 44: /*        f(x-2*h), f(x-h), f(x), f(x+h),  f(x+2*h), */

 46: /*     in the array fval. The use of nf = 7 function evaluations is */
 47: /*     recommended. */

 49: /*     The noise in the function is roughly defined as the variance in */
 50: /*     the computed value of the function. The noise in the function */
 51: /*     provides valuable information. For example, function values */
 52: /*     smaller than the noise should be considered to be zero. */

 54: /*     This subroutine requires an initial estimate for h. Under estimates */
 55: /*     are usually preferred. If noise is not detected, the user should */
 56: /*     increase or decrease h according to the ouput value of info. */
 57: /*     In most cases, the subroutine detects noise with the initial */
 58: /*     value of h. */

 60: /*     The subroutine statement is */

 62: /*       subroutine dnest(nf,fval,h,hopt,fnoise,info,eps) */

 64: /*     where */

 66: /*       nf is an int variable. */
 67: /*         On entry nf is the number of function values. */
 68: /*         On exit nf is unchanged. */

 70: /*       f is a double precision array of dimension nf. */
 71: /*         On entry f contains the function values. */
 72: /*         On exit f is overwritten. */

 74: /*       h is a double precision variable. */
 75: /*         On entry h is an estimate of the optimal difference parameter. */
 76: /*         On exit h is unchanged. */

 78: /*       fnoise is a double precision variable. */
 79: /*         On entry fnoise need not be specified. */
 80: /*         On exit fnoise is set to an estimate of the function noise */
 81: /*            if noise is detected; otherwise fnoise is set to zero. */

 83: /*       hopt is a double precision variable. */
 84: /*         On entry hopt need not be specified. */
 85: /*         On exit hopt is set to an estimate of the optimal difference */
 86: /*            parameter if noise is detected; otherwise hopt is set to zero. */

 88: /*       info is an int variable. */
 89: /*         On entry info need not be specified. */
 90: /*         On exit info is set as follows: */

 92: /*            info = 1  Noise has been detected. */

 94: /*            info = 2  Noise has not been detected; h is too small. */
 95: /*                      Try 100*h for the next value of h. */

 97: /*            info = 3  Noise has not been detected; h is too large. */
 98: /*                      Try h/100 for the next value of h. */

100: /*            info = 4  Noise has been detected but the estimate of hopt */
101: /*                      is not reliable; h is too small. */

103: /*       eps is a double precision work array of dimension nf. */

105: /*     MINPACK-2 Project. April 1997. */
106: /*     Argonne National Laboratory. */
107: /*     Jorge J. More'. */

109: /*     ********** */
110:   /* Parameter adjustments */
111:   --eps;
112:   --fval;

114:   /* Function Body */
115:   *fnoise = 0.;
116:   *fder2  = 0.;
117:   *hopt   = 0.;
118: /*     Compute an estimate of the second derivative and */
119: /*     determine a bound on the error. */
120:   mh   = (*nf + 1) / 2;
121:   est1 = (fval[mh + 1] - fval[mh] * 2 + fval[mh - 1]) / *h__ / *h__;
122:   est2 = (fval[mh + 2] - fval[mh] * 2 + fval[mh - 2]) / (*h__ * 2) / (*h__ * 2);
123:   est3 = (fval[mh + 3] - fval[mh] * 2 + fval[mh - 3]) / (*h__ * 3) / (*h__ * 3);
124:   est4 = (est1 + est2 + est3) / 3;
125: /* Computing MAX */
126: /* Computing PETSCMAX */
127:   d__3 = PetscMax(est1,est2);
128: /* Computing MIN */
129:   d__4 = PetscMin(est1,est2);
130:   d__1 = PetscMax(d__3,est3) - est4, d__2 = est4 - PetscMin(d__4,est3);
131:   err2 = PetscMax(d__1,d__2);
132: /*      write (2,123) est1, est2, est3 */
133: /* 123  format ('Second derivative estimates', 3d12.2) */
134:   if (err2 <= PetscAbsScalar(est4) * .1) *fder2 = est4;
135:   else if (err2 < PetscAbsScalar(est4))  *fder2 = est3;
136:   else *fder2 = 0.;

138: /*     Compute the range of function values. */
139:   f_min = fval[1];
140:   f_max = fval[1];
141:   i__1  = *nf;
142:   for (i__ = 2; i__ <= i__1; ++i__) {
143:     /* Computing MIN */
144:     d__1 = f_min, d__2 = fval[i__];
145:     f_min = PetscMin(d__1,d__2);

147:     /* Computing MAX */
148:     d__1 = f_max, d__2 = fval[i__];
149:     f_max = PetscMax(d__1,d__2);
150:   }
151: /*     Construct the difference table. */
152:   dnoise = FALSE_;
153:   for (j = 1; j <= 6; ++j) {
154:     dsgn[j - 1]   = FALSE_;
155:     cancel[j - 1] = FALSE_;
156:     scale         = 0.;
157:     i__1          = *nf - j;
158:     for (i__ = 1; i__ <= i__1; ++i__) {
159:       fval[i__] = fval[i__ + 1] - fval[i__];
160:       if (fval[i__] == 0.) cancel[j - 1] = TRUE_;

162:       /* Computing MAX */
163:       d__2 = scale, d__3 = (d__1 = fval[i__], PetscAbsScalar(d__1));
164:       scale = PetscMax(d__2,d__3);
165:     }

167:     /*        Compute the estimates for the noise level. */
168:     if (scale == 0.) stdv = 0.;
169:     else {
170:       stdv = 0.;
171:       i__1 = *nf - j;
172:       for (i__ = 1; i__ <= i__1; ++i__) {
173:         /* Computing 2nd power */
174:         d__1 = fval[i__] / scale;
175:         stdv += d__1 * d__1;
176:       }
177:       stdv = scale * PetscSqrtScalar(stdv / (*nf - j));
178:     }
179:     eps[j] = const__[j - 1] * stdv;
180: /*        Determine differences in sign. */
181:     i__1 = *nf - j - 1;
182:     for (i__ = 1; i__ <= i__1; ++i__) {
183:       /* Computing MIN */
184:       d__1 = fval[i__], d__2 = fval[i__ + 1];
185:       /* Computing MAX */
186:       d__3 = fval[i__], d__4 = fval[i__ + 1];
187:       if (PetscMin(d__1,d__2) < 0. && PetscMax(d__3,d__4) > 0.) dsgn[j - 1] = TRUE_;
188:     }
189:   }
190:   /*     First requirement for detection of noise. */
191:   dnoise = dsgn[3];
192:   /*     Check for h too small or too large. */
193:   *info = 0;
194:   if (f_max == f_min) *info = 2;
195:   else /* if (complicated condition) */ {
196:     /* Computing MIN */
197:     d__1 = PetscAbsScalar(f_max), d__2 = PetscAbsScalar(f_min);
198:     if (f_max - f_min > PetscMin(d__1,d__2) * .1) *info = 3;
199:   }
200:   if (*info != 0) return(0);

202:   /*     Determine the noise level. */
203:   /* Computing MIN */
204:   d__1 = PetscMin(eps[4],eps[5]);
205:   emin = PetscMin(d__1,eps[6]);

207:   /* Computing MAX */
208:   d__1 = PetscMax(eps[4],eps[5]);
209:   emax = PetscMax(d__1,eps[6]);

211:   if (emax <= emin * 4 && dnoise) {
212:     *fnoise = (eps[4] + eps[5] + eps[6]) / 3;
213:     if (*fder2 != 0.) {
214:       *info = 1;
215:       *hopt = PetscSqrtScalar(*fnoise / PetscAbsScalar(*fder2)) * 1.68;
216:     } else {
217:       *info = 4;
218:       *hopt = *h__ * 10;
219:     }
220:     return(0);
221:   }

223:   /* Computing MIN */
224:   d__1 = PetscMin(eps[3],eps[4]);
225:   emin = PetscMin(d__1,eps[5]);

227:   /* Computing MAX */
228:   d__1 = PetscMax(eps[3],eps[4]);
229:   emax = PetscMax(d__1,eps[5]);

231:   if (emax <= emin * 4 && dnoise) {
232:     *fnoise = (eps[3] + eps[4] + eps[5]) / 3;
233:     if (*fder2 != 0.) {
234:       *info = 1;
235:       *hopt = PetscSqrtScalar(*fnoise / PetscAbsScalar(*fder2)) * 1.68;
236:     } else {
237:       *info = 4;
238:       *hopt = *h__ * 10;
239:     }
240:     return(0);
241:   }
242: /*     Noise not detected; decide if h is too small or too large. */
243:   if (!cancel[3]) {
244:     if (dsgn[3]) *info = 2;
245:     else *info = 3;
246:     return(0);
247:   }
248:   if (!cancel[2]) {
249:     if (dsgn[2]) *info = 2;
250:     else *info = 3;
251:     return(0);
252:   }
253: /*     If there is cancelllation on the third and fourth column */
254: /*     then h is too small */
255:   *info = 2;
256:   return(0);
257: /*      if (cancel .or. dsgn(3)) then */
258: /*         info = 2 */
259: /*      else */
260: /*         info = 3 */
261: /*      end if */
262: } /* dnest_ */