xref: /petsc/src/snes/interface/noise/snesdnest.c (revision 63dd3a1af947c5d07342e150a17d503f381a847e)
1*63dd3a1aSKris Buschelman #define PETSCSNES_DLL
2*63dd3a1aSKris Buschelman 
354f251a9SBarry Smith /* fnoise/snesdnest.F -- translated by f2c (version 20020314).
454f251a9SBarry Smith */
554f251a9SBarry Smith #include "petsc.h"
654f251a9SBarry Smith #define FALSE_ 0
754f251a9SBarry Smith #define TRUE_ 1
854f251a9SBarry Smith 
954f251a9SBarry Smith /*  Noise estimation routine, written by Jorge More'.  Details are below. */
1054f251a9SBarry Smith 
11a7cc72afSBarry Smith /* Subroutine */ PetscErrorCode dnest_(PetscInt *nf, double *fval,double *h__,double *fnoise, double *fder2, double *hopt, PetscInt *info, double *eps)
1254f251a9SBarry Smith {
1354f251a9SBarry Smith     /* Initialized data */
1454f251a9SBarry Smith 
1554f251a9SBarry Smith     static double const__[15] = { .71,.41,.23,.12,.063,.033,.018,.0089,
1654f251a9SBarry Smith 	    .0046,.0024,.0012,6.1e-4,3.1e-4,1.6e-4,8e-5 };
1754f251a9SBarry Smith 
1854f251a9SBarry Smith     /* System generated locals */
19a7cc72afSBarry Smith     PetscInt i__1;
2054f251a9SBarry Smith     double d__1, d__2, d__3, d__4;
2154f251a9SBarry Smith 
2254f251a9SBarry Smith 
2354f251a9SBarry Smith     /* Local variables */
2454f251a9SBarry Smith     static double emin, emax;
25a7cc72afSBarry Smith     static PetscInt dsgn[6];
268e653690SSatish Balay     static double f_max, f_min, stdv;
27a7cc72afSBarry Smith     static PetscInt i__, j;
2854f251a9SBarry Smith     static double scale;
29a7cc72afSBarry Smith     static PetscInt mh;
30a7cc72afSBarry Smith     static PetscInt cancel[6], dnoise;
3154f251a9SBarry Smith     static double err2, est1, est2, est3, est4;
3254f251a9SBarry Smith 
3354f251a9SBarry Smith /*     ********** */
3454f251a9SBarry Smith 
3554f251a9SBarry Smith /*     Subroutine dnest */
3654f251a9SBarry Smith 
3754f251a9SBarry Smith /*     This subroutine estimates the noise in a function */
3854f251a9SBarry Smith /*     and provides estimates of the optimal difference parameter */
3954f251a9SBarry Smith /*     for a forward-difference approximation. */
4054f251a9SBarry Smith 
4154f251a9SBarry Smith /*     The user must provide a difference parameter h, and the */
4254f251a9SBarry Smith /*     function value at nf points centered around the current point. */
4354f251a9SBarry Smith /*     For example, if nf = 7, the user must provide */
4454f251a9SBarry Smith 
4554f251a9SBarry Smith /*        f(x-2*h), f(x-h), f(x), f(x+h),  f(x+2*h), */
4654f251a9SBarry Smith 
4754f251a9SBarry Smith /*     in the array fval. The use of nf = 7 function evaluations is */
4854f251a9SBarry Smith /*     recommended. */
4954f251a9SBarry Smith 
5054f251a9SBarry Smith /*     The noise in the function is roughly defined as the variance in */
5154f251a9SBarry Smith /*     the computed value of the function. The noise in the function */
5254f251a9SBarry Smith /*     provides valuable information. For example, function values */
5354f251a9SBarry Smith /*     smaller than the noise should be considered to be zero. */
5454f251a9SBarry Smith 
5554f251a9SBarry Smith /*     This subroutine requires an initial estimate for h. Under estimates */
5654f251a9SBarry Smith /*     are usually preferred. If noise is not detected, the user should */
5754f251a9SBarry Smith /*     increase or decrease h according to the ouput value of info. */
5854f251a9SBarry Smith /*     In most cases, the subroutine detects noise with the initial */
5954f251a9SBarry Smith /*     value of h. */
6054f251a9SBarry Smith 
6154f251a9SBarry Smith /*     The subroutine statement is */
6254f251a9SBarry Smith 
6354f251a9SBarry Smith /*       subroutine dnest(nf,fval,h,hopt,fnoise,info,eps) */
6454f251a9SBarry Smith 
6554f251a9SBarry Smith /*     where */
6654f251a9SBarry Smith 
6754f251a9SBarry Smith /*       nf is an int variable. */
6854f251a9SBarry Smith /*         On entry nf is the number of function values. */
6954f251a9SBarry Smith /*         On exit nf is unchanged. */
7054f251a9SBarry Smith 
7154f251a9SBarry Smith /*       f is a double precision array of dimension nf. */
7254f251a9SBarry Smith /*         On entry f contains the function values. */
7354f251a9SBarry Smith /*         On exit f is overwritten. */
7454f251a9SBarry Smith 
7554f251a9SBarry Smith /*       h is a double precision variable. */
7654f251a9SBarry Smith /*         On entry h is an estimate of the optimal difference parameter. */
7754f251a9SBarry Smith /*         On exit h is unchanged. */
7854f251a9SBarry Smith 
7954f251a9SBarry Smith /*       fnoise is a double precision variable. */
8054f251a9SBarry Smith /*         On entry fnoise need not be specified. */
8154f251a9SBarry Smith /*         On exit fnoise is set to an estimate of the function noise */
8254f251a9SBarry Smith /*            if noise is detected; otherwise fnoise is set to zero. */
8354f251a9SBarry Smith 
8454f251a9SBarry Smith /*       hopt is a double precision variable. */
8554f251a9SBarry Smith /*         On entry hopt need not be specified. */
8654f251a9SBarry Smith /*         On exit hopt is set to an estimate of the optimal difference */
8754f251a9SBarry Smith /*            parameter if noise is detected; otherwise hopt is set to zero. */
8854f251a9SBarry Smith 
8954f251a9SBarry Smith /*       info is an int variable. */
9054f251a9SBarry Smith /*         On entry info need not be specified. */
9154f251a9SBarry Smith /*         On exit info is set as follows: */
9254f251a9SBarry Smith 
9354f251a9SBarry Smith /*            info = 1  Noise has been detected. */
9454f251a9SBarry Smith 
9554f251a9SBarry Smith /*            info = 2  Noise has not been detected; h is too small. */
9654f251a9SBarry Smith /*                      Try 100*h for the next value of h. */
9754f251a9SBarry Smith 
9854f251a9SBarry Smith /*            info = 3  Noise has not been detected; h is too large. */
9954f251a9SBarry Smith /*                      Try h/100 for the next value of h. */
10054f251a9SBarry Smith 
10154f251a9SBarry Smith /*            info = 4  Noise has been detected but the estimate of hopt */
10254f251a9SBarry Smith /*                      is not reliable; h is too small. */
10354f251a9SBarry Smith 
10454f251a9SBarry Smith /*       eps is a double precision work array of dimension nf. */
10554f251a9SBarry Smith 
10654f251a9SBarry Smith /*     MINPACK-2 Project. April 1997. */
10754f251a9SBarry Smith /*     Argonne National Laboratory. */
10854f251a9SBarry Smith /*     Jorge J. More'. */
10954f251a9SBarry Smith 
11054f251a9SBarry Smith /*     ********** */
11154f251a9SBarry Smith     /* Parameter adjustments */
11254f251a9SBarry Smith     --eps;
11354f251a9SBarry Smith     --fval;
11454f251a9SBarry Smith 
11554f251a9SBarry Smith     /* Function Body */
11654f251a9SBarry Smith     *fnoise = 0.;
11754f251a9SBarry Smith     *fder2 = 0.;
11854f251a9SBarry Smith     *hopt = 0.;
11954f251a9SBarry Smith /*     Compute an estimate of the second derivative and */
12054f251a9SBarry Smith /*     determine a bound on the error. */
12154f251a9SBarry Smith     mh = (*nf + 1) / 2;
12254f251a9SBarry Smith     est1 = (fval[mh + 1] - fval[mh] * 2 + fval[mh - 1]) / *h__ / *h__;
12354f251a9SBarry Smith     est2 = (fval[mh + 2] - fval[mh] * 2 + fval[mh - 2]) / (*h__ * 2) / (*h__ *
12454f251a9SBarry Smith 	     2);
12554f251a9SBarry Smith     est3 = (fval[mh + 3] - fval[mh] * 2 + fval[mh - 3]) / (*h__ * 3) / (*h__ *
12654f251a9SBarry Smith 	     3);
12754f251a9SBarry Smith     est4 = (est1 + est2 + est3) / 3;
12854f251a9SBarry Smith /* Computing MAX */
12954f251a9SBarry Smith /* Computing PETSCMAX */
13054f251a9SBarry Smith     d__3 = PetscMax(est1,est2);
13154f251a9SBarry Smith /* Computing MIN */
13254f251a9SBarry Smith     d__4 = PetscMin(est1,est2);
13354f251a9SBarry Smith     d__1 = PetscMax(d__3,est3) - est4, d__2 = est4 - PetscMin(d__4,est3);
13454f251a9SBarry Smith     err2 = PetscMax(d__1,d__2);
13554f251a9SBarry Smith /*      write (2,123) est1, est2, est3 */
13654f251a9SBarry Smith /* 123  format ('Second derivative estimates', 3d12.2) */
13754f251a9SBarry Smith     if (err2 <= PetscAbsScalar(est4) * .1) {
13854f251a9SBarry Smith 	*fder2 = est4;
13954f251a9SBarry Smith     } else if (err2 < PetscAbsScalar(est4)) {
14054f251a9SBarry Smith 	*fder2 = est3;
14154f251a9SBarry Smith     } else {
14254f251a9SBarry Smith 	*fder2 = 0.;
14354f251a9SBarry Smith     }
14454f251a9SBarry Smith /*     Compute the range of function values. */
1458e653690SSatish Balay     f_min = fval[1];
1468e653690SSatish Balay     f_max = fval[1];
14754f251a9SBarry Smith     i__1 = *nf;
14854f251a9SBarry Smith     for (i__ = 2; i__ <= i__1; ++i__) {
14954f251a9SBarry Smith /* Computing MIN */
1508e653690SSatish Balay 	d__1 = f_min, d__2 = fval[i__];
1518e653690SSatish Balay 	f_min = PetscMin(d__1,d__2);
15254f251a9SBarry Smith /* Computing MAX */
1538e653690SSatish Balay 	d__1 = f_max, d__2 = fval[i__];
1548e653690SSatish Balay 	f_max = PetscMax(d__1,d__2);
15554f251a9SBarry Smith     }
15654f251a9SBarry Smith /*     Construct the difference table. */
15754f251a9SBarry Smith     dnoise = FALSE_;
15854f251a9SBarry Smith     for (j = 1; j <= 6; ++j) {
15954f251a9SBarry Smith 	dsgn[j - 1] = FALSE_;
16054f251a9SBarry Smith 	cancel[j - 1] = FALSE_;
16154f251a9SBarry Smith 	scale = 0.;
16254f251a9SBarry Smith 	i__1 = *nf - j;
16354f251a9SBarry Smith 	for (i__ = 1; i__ <= i__1; ++i__) {
16454f251a9SBarry Smith 	    fval[i__] = fval[i__ + 1] - fval[i__];
16554f251a9SBarry Smith 	    if (fval[i__] == 0.) {
16654f251a9SBarry Smith 		cancel[j - 1] = TRUE_;
16754f251a9SBarry Smith 	    }
16854f251a9SBarry Smith /* Computing MAX */
16954f251a9SBarry Smith 	    d__2 = scale, d__3 = (d__1 = fval[i__], PetscAbsScalar(d__1));
17054f251a9SBarry Smith 	    scale = PetscMax(d__2,d__3);
17154f251a9SBarry Smith 	}
17254f251a9SBarry Smith /*        Compute the estimates for the noise level. */
17354f251a9SBarry Smith 	if (scale == 0.) {
17454f251a9SBarry Smith 	    stdv = 0.;
17554f251a9SBarry Smith 	} else {
17654f251a9SBarry Smith 	    stdv = 0.;
17754f251a9SBarry Smith 	    i__1 = *nf - j;
17854f251a9SBarry Smith 	    for (i__ = 1; i__ <= i__1; ++i__) {
17954f251a9SBarry Smith /* Computing 2nd power */
18054f251a9SBarry Smith 		d__1 = fval[i__] / scale;
18154f251a9SBarry Smith 		stdv += d__1 * d__1;
18254f251a9SBarry Smith 	    }
18354f251a9SBarry Smith 	    stdv = scale * PetscSqrtScalar(stdv / (*nf - j));
18454f251a9SBarry Smith 	}
18554f251a9SBarry Smith 	eps[j] = const__[j - 1] * stdv;
18654f251a9SBarry Smith /*        Determine differences in sign. */
18754f251a9SBarry Smith 	i__1 = *nf - j - 1;
18854f251a9SBarry Smith 	for (i__ = 1; i__ <= i__1; ++i__) {
18954f251a9SBarry Smith /* Computing MIN */
19054f251a9SBarry Smith 	    d__1 = fval[i__], d__2 = fval[i__ + 1];
19154f251a9SBarry Smith /* Computing MAX */
19254f251a9SBarry Smith 	    d__3 = fval[i__], d__4 = fval[i__ + 1];
19354f251a9SBarry Smith 	    if (PetscMin(d__1,d__2) < 0. && PetscMax(d__3,d__4) > 0.) {
19454f251a9SBarry Smith 		dsgn[j - 1] = TRUE_;
19554f251a9SBarry Smith 	    }
19654f251a9SBarry Smith 	}
19754f251a9SBarry Smith     }
19854f251a9SBarry Smith /*     First requirement for detection of noise. */
19954f251a9SBarry Smith     dnoise = dsgn[3];
20054f251a9SBarry Smith /*     Check for h too small or too large. */
20154f251a9SBarry Smith     *info = 0;
2028e653690SSatish Balay     if (f_max == f_min) {
20354f251a9SBarry Smith 	*info = 2;
20454f251a9SBarry Smith     } else /* if(complicated condition) */ {
20554f251a9SBarry Smith /* Computing MIN */
2068e653690SSatish Balay 	d__1 = PetscAbsScalar(f_max), d__2 = PetscAbsScalar(f_min);
2078e653690SSatish Balay 	if (f_max - f_min > PetscMin(d__1,d__2) * .1) {
20854f251a9SBarry Smith 	    *info = 3;
20954f251a9SBarry Smith 	}
21054f251a9SBarry Smith     }
21154f251a9SBarry Smith     if (*info != 0) {
212a7cc72afSBarry Smith 	PetscFunctionReturn(0);
21354f251a9SBarry Smith     }
21454f251a9SBarry Smith /*     Determine the noise level. */
21554f251a9SBarry Smith /* Computing MIN */
21654f251a9SBarry Smith     d__1 = PetscMin(eps[4],eps[5]);
21754f251a9SBarry Smith     emin = PetscMin(d__1,eps[6]);
21854f251a9SBarry Smith /* Computing MAX */
21954f251a9SBarry Smith     d__1 = PetscMax(eps[4],eps[5]);
22054f251a9SBarry Smith     emax = PetscMax(d__1,eps[6]);
22154f251a9SBarry Smith     if (emax <= emin * 4 && dnoise) {
22254f251a9SBarry Smith 	*fnoise = (eps[4] + eps[5] + eps[6]) / 3;
22354f251a9SBarry Smith 	if (*fder2 != 0.) {
22454f251a9SBarry Smith 	    *info = 1;
22554f251a9SBarry Smith 	    *hopt = PetscSqrtScalar(*fnoise / PetscAbsScalar(*fder2)) * 1.68;
22654f251a9SBarry Smith 	} else {
22754f251a9SBarry Smith 	    *info = 4;
22854f251a9SBarry Smith 	    *hopt = *h__ * 10;
22954f251a9SBarry Smith 	}
230a7cc72afSBarry Smith 	PetscFunctionReturn(0);
23154f251a9SBarry Smith     }
23254f251a9SBarry Smith /* Computing MIN */
23354f251a9SBarry Smith     d__1 = PetscMin(eps[3],eps[4]);
23454f251a9SBarry Smith     emin = PetscMin(d__1,eps[5]);
23554f251a9SBarry Smith /* Computing MAX */
23654f251a9SBarry Smith     d__1 = PetscMax(eps[3],eps[4]);
23754f251a9SBarry Smith     emax = PetscMax(d__1,eps[5]);
23854f251a9SBarry Smith     if (emax <= emin * 4 && dnoise) {
23954f251a9SBarry Smith 	*fnoise = (eps[3] + eps[4] + eps[5]) / 3;
24054f251a9SBarry Smith 	if (*fder2 != 0.) {
24154f251a9SBarry Smith 	    *info = 1;
24254f251a9SBarry Smith 	    *hopt = PetscSqrtScalar(*fnoise / PetscAbsScalar(*fder2)) * 1.68;
24354f251a9SBarry Smith 	} else {
24454f251a9SBarry Smith 	    *info = 4;
24554f251a9SBarry Smith 	    *hopt = *h__ * 10;
24654f251a9SBarry Smith 	}
247a7cc72afSBarry Smith 	PetscFunctionReturn(0);
24854f251a9SBarry Smith     }
24954f251a9SBarry Smith /*     Noise not detected; decide if h is too small or too large. */
25054f251a9SBarry Smith     if (! cancel[3]) {
25154f251a9SBarry Smith 	if (dsgn[3]) {
25254f251a9SBarry Smith 	    *info = 2;
25354f251a9SBarry Smith 	} else {
25454f251a9SBarry Smith 	    *info = 3;
25554f251a9SBarry Smith 	}
256a7cc72afSBarry Smith 	PetscFunctionReturn(0);
25754f251a9SBarry Smith     }
25854f251a9SBarry Smith     if (! cancel[2]) {
25954f251a9SBarry Smith 	if (dsgn[2]) {
26054f251a9SBarry Smith 	    *info = 2;
26154f251a9SBarry Smith 	} else {
26254f251a9SBarry Smith 	    *info = 3;
26354f251a9SBarry Smith 	}
264a7cc72afSBarry Smith 	PetscFunctionReturn(0);
26554f251a9SBarry Smith     }
26654f251a9SBarry Smith /*     If there is cancelllation on the third and fourth column */
26754f251a9SBarry Smith /*     then h is too small */
26854f251a9SBarry Smith     *info = 2;
269a7cc72afSBarry Smith     PetscFunctionReturn(0);
27054f251a9SBarry Smith /*      if (cancel .or. dsgn(3)) then */
27154f251a9SBarry Smith /*         info = 2 */
27254f251a9SBarry Smith /*      else */
27354f251a9SBarry Smith /*         info = 3 */
27454f251a9SBarry Smith /*      end if */
27554f251a9SBarry Smith } /* dnest_ */
27654f251a9SBarry Smith 
277