163dd3a1aSKris Buschelman 254f251a9SBarry Smith /* fnoise/snesdnest.F -- translated by f2c (version 20020314). 354f251a9SBarry Smith */ 4c6db04a5SJed Brown #include <petscsys.h> 554f251a9SBarry Smith #define FALSE_ 0 654f251a9SBarry Smith #define TRUE_ 1 754f251a9SBarry Smith 854f251a9SBarry Smith /* Noise estimation routine, written by Jorge More'. Details are below. */ 954f251a9SBarry Smith 1025acbd8eSLisandro Dalcin PETSC_INTERN PetscErrorCode SNESNoise_dnest_(PetscInt*,PetscScalar*,PetscScalar*,PetscScalar*,PetscScalar*,PetscScalar*,PetscInt*,PetscScalar*); 1125acbd8eSLisandro Dalcin 1225acbd8eSLisandro Dalcin PetscErrorCode SNESNoise_dnest_(PetscInt *nf, double *fval,double *h__,double *fnoise, double *fder2, double *hopt, PetscInt *info, double *eps) 1354f251a9SBarry Smith { 1454f251a9SBarry Smith /* Initialized data */ 1554f251a9SBarry Smith 1654f251a9SBarry Smith static double const__[15] = { .71,.41,.23,.12,.063,.033,.018,.0089, 1754f251a9SBarry Smith .0046,.0024,.0012,6.1e-4,3.1e-4,1.6e-4,8e-5 }; 1854f251a9SBarry Smith 1954f251a9SBarry Smith /* System generated locals */ 20a7cc72afSBarry Smith PetscInt i__1; 2154f251a9SBarry Smith double d__1, d__2, d__3, d__4; 2254f251a9SBarry Smith 2354f251a9SBarry Smith 2454f251a9SBarry Smith /* Local variables */ 2554f251a9SBarry Smith static double emin, emax; 26a7cc72afSBarry Smith static PetscInt dsgn[6]; 278e653690SSatish Balay static double f_max, f_min, stdv; 28a7cc72afSBarry Smith static PetscInt i__, j; 2954f251a9SBarry Smith static double scale; 30a7cc72afSBarry Smith static PetscInt mh; 31a7cc72afSBarry Smith static PetscInt cancel[6], dnoise; 3254f251a9SBarry Smith static double err2, est1, est2, est3, est4; 3354f251a9SBarry Smith 3454f251a9SBarry Smith /* ********** */ 3554f251a9SBarry Smith 3654f251a9SBarry Smith /* Subroutine dnest */ 3754f251a9SBarry Smith 3854f251a9SBarry Smith /* This subroutine estimates the noise in a function */ 3954f251a9SBarry Smith /* and provides estimates of the optimal difference parameter */ 4054f251a9SBarry Smith /* for a forward-difference approximation. */ 4154f251a9SBarry Smith 4254f251a9SBarry Smith /* The user must provide a difference parameter h, and the */ 4354f251a9SBarry Smith /* function value at nf points centered around the current point. */ 4454f251a9SBarry Smith /* For example, if nf = 7, the user must provide */ 4554f251a9SBarry Smith 4654f251a9SBarry Smith /* f(x-2*h), f(x-h), f(x), f(x+h), f(x+2*h), */ 4754f251a9SBarry Smith 4854f251a9SBarry Smith /* in the array fval. The use of nf = 7 function evaluations is */ 4954f251a9SBarry Smith /* recommended. */ 5054f251a9SBarry Smith 5154f251a9SBarry Smith /* The noise in the function is roughly defined as the variance in */ 5254f251a9SBarry Smith /* the computed value of the function. The noise in the function */ 5354f251a9SBarry Smith /* provides valuable information. For example, function values */ 5454f251a9SBarry Smith /* smaller than the noise should be considered to be zero. */ 5554f251a9SBarry Smith 5654f251a9SBarry Smith /* This subroutine requires an initial estimate for h. Under estimates */ 5754f251a9SBarry Smith /* are usually preferred. If noise is not detected, the user should */ 5854f251a9SBarry Smith /* increase or decrease h according to the ouput value of info. */ 5954f251a9SBarry Smith /* In most cases, the subroutine detects noise with the initial */ 6054f251a9SBarry Smith /* value of h. */ 6154f251a9SBarry Smith 6254f251a9SBarry Smith /* The subroutine statement is */ 6354f251a9SBarry Smith 6454f251a9SBarry Smith /* subroutine dnest(nf,fval,h,hopt,fnoise,info,eps) */ 6554f251a9SBarry Smith 6654f251a9SBarry Smith /* where */ 6754f251a9SBarry Smith 68*5bd1e576SStefano Zampini /* nf is a PetscInt variable. */ 6954f251a9SBarry Smith /* On entry nf is the number of function values. */ 7054f251a9SBarry Smith /* On exit nf is unchanged. */ 7154f251a9SBarry Smith 7254f251a9SBarry Smith /* f is a double precision array of dimension nf. */ 7354f251a9SBarry Smith /* On entry f contains the function values. */ 7454f251a9SBarry Smith /* On exit f is overwritten. */ 7554f251a9SBarry Smith 7654f251a9SBarry Smith /* h is a double precision variable. */ 7754f251a9SBarry Smith /* On entry h is an estimate of the optimal difference parameter. */ 7854f251a9SBarry Smith /* On exit h is unchanged. */ 7954f251a9SBarry Smith 8054f251a9SBarry Smith /* fnoise is a double precision variable. */ 8154f251a9SBarry Smith /* On entry fnoise need not be specified. */ 8254f251a9SBarry Smith /* On exit fnoise is set to an estimate of the function noise */ 8354f251a9SBarry Smith /* if noise is detected; otherwise fnoise is set to zero. */ 8454f251a9SBarry Smith 8554f251a9SBarry Smith /* hopt is a double precision variable. */ 8654f251a9SBarry Smith /* On entry hopt need not be specified. */ 8754f251a9SBarry Smith /* On exit hopt is set to an estimate of the optimal difference */ 8854f251a9SBarry Smith /* parameter if noise is detected; otherwise hopt is set to zero. */ 8954f251a9SBarry Smith 90*5bd1e576SStefano Zampini /* info is a PetscInt variable. */ 9154f251a9SBarry Smith /* On entry info need not be specified. */ 9254f251a9SBarry Smith /* On exit info is set as follows: */ 9354f251a9SBarry Smith 9454f251a9SBarry Smith /* info = 1 Noise has been detected. */ 9554f251a9SBarry Smith 9654f251a9SBarry Smith /* info = 2 Noise has not been detected; h is too small. */ 9754f251a9SBarry Smith /* Try 100*h for the next value of h. */ 9854f251a9SBarry Smith 9954f251a9SBarry Smith /* info = 3 Noise has not been detected; h is too large. */ 10054f251a9SBarry Smith /* Try h/100 for the next value of h. */ 10154f251a9SBarry Smith 10254f251a9SBarry Smith /* info = 4 Noise has been detected but the estimate of hopt */ 10354f251a9SBarry Smith /* is not reliable; h is too small. */ 10454f251a9SBarry Smith 10554f251a9SBarry Smith /* eps is a double precision work array of dimension nf. */ 10654f251a9SBarry Smith 10754f251a9SBarry Smith /* MINPACK-2 Project. April 1997. */ 10854f251a9SBarry Smith /* Argonne National Laboratory. */ 10954f251a9SBarry Smith /* Jorge J. More'. */ 11054f251a9SBarry Smith 11154f251a9SBarry Smith /* ********** */ 11254f251a9SBarry Smith /* Parameter adjustments */ 11354f251a9SBarry Smith --eps; 11454f251a9SBarry Smith --fval; 11554f251a9SBarry Smith 11654f251a9SBarry Smith /* Function Body */ 11754f251a9SBarry Smith *fnoise = 0.; 11854f251a9SBarry Smith *fder2 = 0.; 11954f251a9SBarry Smith *hopt = 0.; 12054f251a9SBarry Smith /* Compute an estimate of the second derivative and */ 12154f251a9SBarry Smith /* determine a bound on the error. */ 12254f251a9SBarry Smith mh = (*nf + 1) / 2; 12354f251a9SBarry Smith est1 = (fval[mh + 1] - fval[mh] * 2 + fval[mh - 1]) / *h__ / *h__; 12454eb447fSKarl Rupp est2 = (fval[mh + 2] - fval[mh] * 2 + fval[mh - 2]) / (*h__ * 2) / (*h__ * 2); 12554eb447fSKarl Rupp est3 = (fval[mh + 3] - fval[mh] * 2 + fval[mh - 3]) / (*h__ * 3) / (*h__ * 3); 12654f251a9SBarry Smith est4 = (est1 + est2 + est3) / 3; 12754f251a9SBarry Smith /* Computing MAX */ 12854f251a9SBarry Smith /* Computing PETSCMAX */ 12954f251a9SBarry Smith d__3 = PetscMax(est1,est2); 13054f251a9SBarry Smith /* Computing MIN */ 13154f251a9SBarry Smith d__4 = PetscMin(est1,est2); 1322a274a78SSatish Balay d__1 = PetscMax(d__3,est3) - est4; 1332a274a78SSatish Balay 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) */ 137f5af7f23SKarl Rupp if (err2 <= PetscAbsScalar(est4) * .1) *fder2 = est4; 138f5af7f23SKarl Rupp else if (err2 < PetscAbsScalar(est4)) *fder2 = est3; 139f5af7f23SKarl Rupp else *fder2 = 0.; 140f5af7f23SKarl Rupp 14154f251a9SBarry Smith /* Compute the range of function values. */ 1428e653690SSatish Balay f_min = fval[1]; 1438e653690SSatish Balay f_max = fval[1]; 14454f251a9SBarry Smith i__1 = *nf; 14554f251a9SBarry Smith for (i__ = 2; i__ <= i__1; ++i__) { 14654f251a9SBarry Smith /* Computing MIN */ 1472a274a78SSatish Balay d__1 = f_min; 1482a274a78SSatish Balay d__2 = fval[i__]; 1498e653690SSatish Balay f_min = PetscMin(d__1,d__2); 150f5af7f23SKarl Rupp 15154f251a9SBarry Smith /* Computing MAX */ 1522a274a78SSatish Balay d__1 = f_max; 1532a274a78SSatish Balay 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__]; 165f5af7f23SKarl Rupp if (fval[i__] == 0.) cancel[j - 1] = TRUE_; 166f5af7f23SKarl Rupp 16754f251a9SBarry Smith /* Computing MAX */ 1682a274a78SSatish Balay d__1 = fval[i__]; 1692a274a78SSatish Balay d__2 = scale; 1702a274a78SSatish Balay d__3 = PetscAbsScalar(d__1); 17154f251a9SBarry Smith scale = PetscMax(d__2,d__3); 17254f251a9SBarry Smith } 173f5af7f23SKarl Rupp 17454f251a9SBarry Smith /* Compute the estimates for the noise level. */ 175f5af7f23SKarl Rupp if (scale == 0.) stdv = 0.; 176f5af7f23SKarl Rupp else { 17754f251a9SBarry Smith stdv = 0.; 17854f251a9SBarry Smith i__1 = *nf - j; 17954f251a9SBarry Smith for (i__ = 1; i__ <= i__1; ++i__) { 18054f251a9SBarry Smith /* Computing 2nd power */ 18154f251a9SBarry Smith d__1 = fval[i__] / scale; 18254f251a9SBarry Smith stdv += d__1 * d__1; 18354f251a9SBarry Smith } 18454f251a9SBarry Smith stdv = scale * PetscSqrtScalar(stdv / (*nf - j)); 18554f251a9SBarry Smith } 18654f251a9SBarry Smith eps[j] = const__[j - 1] * stdv; 18754f251a9SBarry Smith /* Determine differences in sign. */ 18854f251a9SBarry Smith i__1 = *nf - j - 1; 18954f251a9SBarry Smith for (i__ = 1; i__ <= i__1; ++i__) { 19054f251a9SBarry Smith /* Computing MIN */ 1912a274a78SSatish Balay d__1 = fval[i__]; 1922a274a78SSatish Balay d__2 = fval[i__ + 1]; 19354f251a9SBarry Smith /* Computing MAX */ 1942a274a78SSatish Balay d__3 = fval[i__]; 1952a274a78SSatish Balay d__4 = fval[i__ + 1]; 196f5af7f23SKarl Rupp if (PetscMin(d__1,d__2) < 0. && PetscMax(d__3,d__4) > 0.) dsgn[j - 1] = TRUE_; 19754f251a9SBarry Smith } 19854f251a9SBarry Smith } 19954f251a9SBarry Smith /* First requirement for detection of noise. */ 20054f251a9SBarry Smith dnoise = dsgn[3]; 20154f251a9SBarry Smith /* Check for h too small or too large. */ 20254f251a9SBarry Smith *info = 0; 203f5af7f23SKarl Rupp if (f_max == f_min) *info = 2; 204f5af7f23SKarl Rupp else /* if (complicated condition) */ { 20554f251a9SBarry Smith /* Computing MIN */ 2062a274a78SSatish Balay d__1 = PetscAbsScalar(f_max); 2072a274a78SSatish Balay d__2 = PetscAbsScalar(f_min); 208f5af7f23SKarl Rupp if (f_max - f_min > PetscMin(d__1,d__2) * .1) *info = 3; 20954f251a9SBarry Smith } 210f5af7f23SKarl Rupp if (*info != 0) PetscFunctionReturn(0); 211f5af7f23SKarl Rupp 21254f251a9SBarry Smith /* Determine the noise level. */ 21354f251a9SBarry Smith /* Computing MIN */ 21454f251a9SBarry Smith d__1 = PetscMin(eps[4],eps[5]); 21554f251a9SBarry Smith emin = PetscMin(d__1,eps[6]); 216f5af7f23SKarl Rupp 21754f251a9SBarry Smith /* Computing MAX */ 21854f251a9SBarry Smith d__1 = PetscMax(eps[4],eps[5]); 21954f251a9SBarry Smith emax = PetscMax(d__1,eps[6]); 220f5af7f23SKarl Rupp 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 } 232f5af7f23SKarl Rupp 23354f251a9SBarry Smith /* Computing MIN */ 23454f251a9SBarry Smith d__1 = PetscMin(eps[3],eps[4]); 23554f251a9SBarry Smith emin = PetscMin(d__1,eps[5]); 236f5af7f23SKarl Rupp 23754f251a9SBarry Smith /* Computing MAX */ 23854f251a9SBarry Smith d__1 = PetscMax(eps[3],eps[4]); 23954f251a9SBarry Smith emax = PetscMax(d__1,eps[5]); 240f5af7f23SKarl Rupp 24154f251a9SBarry Smith if (emax <= emin * 4 && dnoise) { 24254f251a9SBarry Smith *fnoise = (eps[3] + eps[4] + eps[5]) / 3; 24354f251a9SBarry Smith if (*fder2 != 0.) { 24454f251a9SBarry Smith *info = 1; 24554f251a9SBarry Smith *hopt = PetscSqrtScalar(*fnoise / PetscAbsScalar(*fder2)) * 1.68; 24654f251a9SBarry Smith } else { 24754f251a9SBarry Smith *info = 4; 24854f251a9SBarry Smith *hopt = *h__ * 10; 24954f251a9SBarry Smith } 250a7cc72afSBarry Smith PetscFunctionReturn(0); 25154f251a9SBarry Smith } 25254f251a9SBarry Smith /* Noise not detected; decide if h is too small or too large. */ 25354f251a9SBarry Smith if (!cancel[3]) { 254f5af7f23SKarl Rupp if (dsgn[3]) *info = 2; 255f5af7f23SKarl Rupp else *info = 3; 256a7cc72afSBarry Smith PetscFunctionReturn(0); 25754f251a9SBarry Smith } 25854f251a9SBarry Smith if (!cancel[2]) { 259f5af7f23SKarl Rupp if (dsgn[2]) *info = 2; 260f5af7f23SKarl Rupp else *info = 3; 261a7cc72afSBarry Smith PetscFunctionReturn(0); 26254f251a9SBarry Smith } 26354f251a9SBarry Smith /* If there is cancelllation on the third and fourth column */ 26454f251a9SBarry Smith /* then h is too small */ 26554f251a9SBarry Smith *info = 2; 266a7cc72afSBarry Smith PetscFunctionReturn(0); 26754f251a9SBarry Smith /* if (cancel .or. dsgn(3)) then */ 26854f251a9SBarry Smith /* info = 2 */ 26954f251a9SBarry Smith /* else */ 27054f251a9SBarry Smith /* info = 3 */ 27154f251a9SBarry Smith /* end if */ 27254f251a9SBarry Smith } /* dnest_ */ 27354f251a9SBarry Smith 274