Actual source code: snesdnest.c

petsc-3.10.5 2019-03-28
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: PETSC_INTERN PetscErrorCode SNESNoise_dnest_(PetscInt*,PetscScalar*,PetscScalar*,PetscScalar*,PetscScalar*,PetscScalar*,PetscInt*,PetscScalar*);

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

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

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


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

 34: /*     ********** */

 36: /*     Subroutine dnest */

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

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

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

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

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

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

 62: /*     The subroutine statement is */

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

 66: /*     where */

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

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

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

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

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

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

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

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

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

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

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

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

111: /*     ********** */
112:   /* Parameter adjustments */
113:   --eps;
114:   --fval;

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

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

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

167:       /* Computing MAX */
168:       d__1 = fval[i__];
169:       d__2 = scale;
170:       d__3 = PetscAbsScalar(d__1);
171:       scale = PetscMax(d__2,d__3);
172:     }

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

212:   /*     Determine the noise level. */
213:   /* Computing MIN */
214:   d__1 = PetscMin(eps[4],eps[5]);
215:   emin = PetscMin(d__1,eps[6]);

217:   /* Computing MAX */
218:   d__1 = PetscMax(eps[4],eps[5]);
219:   emax = PetscMax(d__1,eps[6]);

221:   if (emax <= emin * 4 && dnoise) {
222:     *fnoise = (eps[4] + eps[5] + eps[6]) / 3;
223:     if (*fder2 != 0.) {
224:       *info = 1;
225:       *hopt = PetscSqrtScalar(*fnoise / PetscAbsScalar(*fder2)) * 1.68;
226:     } else {
227:       *info = 4;
228:       *hopt = *h__ * 10;
229:     }
230:     return(0);
231:   }

233:   /* Computing MIN */
234:   d__1 = PetscMin(eps[3],eps[4]);
235:   emin = PetscMin(d__1,eps[5]);

237:   /* Computing MAX */
238:   d__1 = PetscMax(eps[3],eps[4]);
239:   emax = PetscMax(d__1,eps[5]);

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