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