void gsl_complex_arccos (complex_t const *a, complex_t *res) { /* z = arccos(a) */ gnm_float R = GSL_REAL (a), I = GSL_IMAG (a); if (I == 0) { gsl_complex_arccos_real (R, res); } else { gnm_float x = gnm_abs (R); gnm_float y = gnm_abs (I); gnm_float r = gnm_hypot (x + 1, y); gnm_float s = gnm_hypot (x - 1, y); gnm_float A = 0.5 * (r + s); gnm_float B = x / A; gnm_float y2 = y * y; gnm_float real, imag; const gnm_float A_crossover = 1.5; const gnm_float B_crossover = 0.6417; if (B <= B_crossover) { real = gnm_acos (B); } else { if (x <= 1) { gnm_float D = 0.5 * (A + x) * (y2 / (r + x + 1) + (s + (1 - x))); real = gnm_atan (gnm_sqrt (D) / x); } else { gnm_float Apx = A + x; gnm_float D = 0.5 * (Apx / (r + x + 1) + Apx / (s + (x - 1))); real = gnm_atan ((y * gnm_sqrt (D)) / x); } } if (A <= A_crossover) { gnm_float Am1; if (x < 1) { Am1 = 0.5 * (y2 / (r + (x + 1)) + y2 / (s + (1 - x))); } else { Am1 = 0.5 * (y2 / (r + (x + 1)) + (s + (x - 1))); } imag = gnm_log1p (Am1 + gnm_sqrt (Am1 * (A + 1))); } else { imag = gnm_log (A + gnm_sqrt (A * A - 1)); } complex_init (res, (R >= 0) ? real : M_PIgnum - real, (I >= 0) ? -imag : imag); } }
static char * try_auto_date (GnmValue *value, const GOFormat *format, GODateConventions const *date_conv) { gnm_float v, vr, vs; GOFormat *actual; char *res; gboolean needs_date, needs_time, needs_frac_sec; gboolean is_date; int is_time; GString *xlfmt; GDate date; format = gnm_format_specialize (format, value); is_date = go_format_is_date (format) > 0; is_time = go_format_is_time (format); if (!is_date && is_time <= 0) return NULL; /* We don't want to coerce strings. */ if (!VALUE_IS_FLOAT (value)) return NULL; /* Verify that the date is valid. */ if (!datetime_value_to_g (&date, value, date_conv)) return NULL; v = value_get_as_float (value); vr = gnm_fake_round (v); vs = (24 * 60 * 60) * gnm_abs (v - vr); needs_date = is_time < 2 && (is_date || gnm_abs (v) >= 1); needs_time = is_time > 0 || gnm_abs (v - vr) > 1e-9; needs_frac_sec = needs_time && gnm_abs (vs - gnm_fake_round (vs)) >= 0.5e-3; xlfmt = g_string_new (NULL); if (needs_date) g_string_append (xlfmt, "yyyy/mm/dd"); if (needs_time) { if (needs_date) g_string_append_c (xlfmt, ' '); if (is_time == 2) g_string_append (xlfmt, "[h]:mm:ss"); else g_string_append (xlfmt, "hh:mm:ss"); if (needs_frac_sec) g_string_append (xlfmt, ".000"); } actual = go_format_new_from_XL (xlfmt->str); g_string_free (xlfmt, TRUE); res = format_value (actual, value, -1, date_conv); go_format_unref (actual); return res; }
static gnm_float * compute_gradient (GnmNlsolve *nl, const gnm_float *xs) { gnm_float *g; gnm_float y0; const int n = nl->vars->len; int i; set_vector (nl, xs); y0 = get_value (nl); g = g_new (gnm_float, n); for (i = 0; i < n; i++) { gnm_float x0 = xs[i]; gnm_float dx; gnm_float y1; gnm_float eps = gnm_pow2 (-25); if (x0 == 0) dx = eps; else dx = gnm_abs (x0) * eps; set_value (nl, i, x0 + dx); y1 = get_value (nl); g[i] = (y1 - y0) / dx; set_value (nl, i, x0); } return g; }
static void pochhammer_small_n (gnm_float x, gnm_float n, GnmQuad *res) { GnmQuad qx, qn, qr, qs, f1, f2, f3, f4, f5; gnm_float r; gboolean debug = FALSE; g_return_if_fail (x >= 20); g_return_if_fail (gnm_abs (n) <= 1); /* * G(x) = c * x^(x-1/2) * exp(-x) * E(x) * G(x+n) = c * (x+n)^(x+n-1/2) * exp(-(x+n)) * E(x+n) * = c * (x+n)^(x-1/2) * (x+n)^n * exp(-x) * exp(-n) * E(x+n) * * G(x+n)/G(x) * = (1+n/x)^(x-1/2) * (x+n)^n * exp(-n) * E(x+n)/E(x) * = (1+n/x)^x / sqrt(1+n/x) * (x+n)^n * exp(-n) * E(x+n)/E(x) * = exp(x*log(1+n/x) - n) / sqrt(1+n/x) * (x+n)^n * E(x+n)/E(x) * = exp(x*log1p(n/x) - n) / sqrt(1+n/x) * (x+n)^n * E(x+n)/E(x) * = exp(x*(log1pmx(n/x)+n/x) - n) / sqrt(1+n/x) * (x+n)^n * E(x+n)/E(x) * = exp(x*log1pmx(n/x) + n - n) / sqrt(1+n/x) * (x+n)^n * E(x+n)/E(x) * = exp(x*log1pmx(n/x)) / sqrt(1+n/x) * (x+n)^n * E(x+n)/E(x) */ gnm_quad_init (&qx, x); gnm_quad_init (&qn, n); gnm_quad_div (&qr, &qn, &qx); r = gnm_quad_value (&qr); gnm_quad_add (&qs, &qx, &qn); /* exp(x*log1pmx(n/x)) */ gnm_quad_mul12 (&f1, log1pmx (r), x); /* sub-opt */ gnm_quad_exp (&f1, NULL, &f1); if (debug) g_printerr ("f1=%.20g\n", gnm_quad_value (&f1)); /* sqrt(1+n/x) */ gnm_quad_add (&f2, &gnm_quad_one, &qr); gnm_quad_sqrt (&f2, &f2); if (debug) g_printerr ("f2=%.20g\n", gnm_quad_value (&f2)); /* (x+n)^n */ gnm_quad_pow (&f3, NULL, &qs, &qn); if (debug) g_printerr ("f3=%.20g\n", gnm_quad_value (&f3)); /* E(x+n) */ gamma_error_factor (&f4, &qs); if (debug) g_printerr ("f4=%.20g\n", gnm_quad_value (&f4)); /* E(x) */ gamma_error_factor (&f5, &qx); if (debug) g_printerr ("f5=%.20g\n", gnm_quad_value (&f5)); gnm_quad_div (res, &f1, &f2); gnm_quad_mul (res, res, &f3); gnm_quad_mul (res, res, &f4); gnm_quad_div (res, res, &f5); }
double lgamma_r (double x, int *signp) { *signp = +1; if (gnm_isnan (x)) return gnm_nan; if (x > 0) { gnm_float f = 1; while (x < 10) { f *= x; x++; } return (M_LN_SQRT_2PI + (x - 0.5) * gnm_log(x) - x + lgammacor(x)) - gnm_log (f); } else { gnm_float axm2 = gnm_fmod (-x, 2.0); gnm_float y = gnm_sinpi (axm2) / M_PIgnum; *signp = axm2 > 1.0 ? +1 : -1; return y == 0 ? gnm_nan : - gnm_log (gnm_abs (y)) - lgamma1p (-x); } }
static void gsl_complex_arcsin_real (gnm_float a, complex_t *res) { /* z = arcsin(a) */ if (gnm_abs (a) <= 1.0) { complex_init (res, gnm_asin (a), 0.0); } else { if (a < 0.0) { complex_init (res, -M_PI_2gnum, gnm_acosh (-a)); } else { complex_init (res, M_PI_2gnum, -gnm_acosh (a)); } } }
gnm_float psnorm (gnm_float x, gnm_float shape, gnm_float location, gnm_float scale, gboolean lower_tail, gboolean log_p) { gnm_float result, h; if (gnm_isnan (x) || gnm_isnan (shape) || gnm_isnan (location) || gnm_isnan (scale)) return gnm_nan; if (shape == 0.) return pnorm (x, location, scale, lower_tail, log_p); /* Normalize */ h = (x - location) / scale; /* Flip to a lower-tail problem. */ if (!lower_tail) { h = -h; shape = -shape; lower_tail = !lower_tail; } if (gnm_abs (shape) < 10) { gnm_float s = pnorm (h, 0, 1, lower_tail, FALSE); gnm_float t = 2 * gnm_owent (h, shape); result = s - t; } else { /* * Make use of this result for Owen's T: * * T(h,a) = .5N(h) + .5N(ha) - N(h)N(ha) - T(ha,1/a) */ gnm_float s = pnorm (h * shape, 0, 1, TRUE, FALSE); gnm_float u = gnm_erf (h / M_SQRT2gnum); gnm_float t = 2 * gnm_owent (h * shape, 1 / shape); result = s * u + t; } /* * Negatives can occur due to rounding errors and hopefully for no * other reason. */ result= CLAMP (result, 0.0, 1.0); if (log_p) return gnm_log (result); else return result; }
static gboolean polish_iter (GnmNlsolve *nl) { GnmSolver *sol = nl->parent; const int n = nl->vars->len; gnm_float *x; gnm_float step; gboolean any_at_all = FALSE; x = g_new (gnm_float, n); for (step = gnm_pow2 (-10); step > GNM_EPSILON; step *= 0.75) { int c, s; for (c = 0; c < n; c++) { for (s = 0; s <= 1; s++) { gnm_float y; gnm_float dx = step * gnm_abs (nl->xk[c]); if (dx == 0) dx = step; if (s) dx = -dx; memcpy (x, nl->xk, n * sizeof (gnm_float)); x[c] += dx; set_vector (nl, x); y = get_value (nl); if (y < nl->yk && gnm_solver_check_constraints (sol)) { nl->yk = y; memcpy (nl->xk, x, n * sizeof (gnm_float)); any_at_all = TRUE; if (nl->debug) g_printerr ("Polish step %.15" GNM_FORMAT_g " in direction %d\n", dx, c); break; } } } } g_free (x); if (any_at_all) gnm_nlsolve_set_solution (nl); return any_at_all; }
void gsl_complex_arctan (complex_t const *a, complex_t *res) { /* z = arctan(a) */ gnm_float R = GSL_REAL (a), I = GSL_IMAG (a); if (I == 0) { complex_init (res, gnm_atan (R), 0); } else { /* FIXME: This is a naive implementation which does not fully * take into account cancellation errors, overflow, underflow * etc. It would benefit from the Hull et al treatment. */ gnm_float r = gnm_hypot (R, I); gnm_float imag; gnm_float u = 2 * I / (1 + r * r); /* FIXME: the following cross-over should be optimized but 0.1 * seems to work ok */ if (gnm_abs (u) < 0.1) { imag = 0.25 * (gnm_log1p (u) - gnm_log1p (-u)); } else { gnm_float A = gnm_hypot (R, I + 1); gnm_float B = gnm_hypot (R, I - 1); imag = 0.5 * gnm_log (A / B); } if (R == 0) { if (I > 1) { complex_init (res, M_PI_2gnum, imag); } else if (I < -1) { complex_init (res, -M_PI_2gnum, imag); } else { complex_init (res, 0, imag); } } else { complex_init (res, 0.5 * gnm_atan2 (2 * R, ((1 + r) * (1 - r))), imag); } } }
void gsl_complex_tanh (complex_t const *a, complex_t *res) { /* z = tanh(a) */ gnm_float R = GSL_REAL (a), I = GSL_IMAG (a); if (gnm_abs (R) < 1.0) { gnm_float D = gnm_pow (gnm_cos (I), 2.0) + gnm_pow (gnm_sinh (R), 2.0); complex_init (res, gnm_sinh (R) * cosh (R) / D, 0.5 * gnm_sin (2 * I) / D); } else { gnm_float D = gnm_pow (gnm_cos (I), 2.0) + gnm_pow (gnm_sinh (R), 2.0); gnm_float F = 1 + gnm_pow (gnm_cos (I) / gnm_sinh (R), 2.0); complex_init (res, 1.0 / (gnm_tanh (R) * F), 0.5 * gnm_sin (2 * I) / D); } }
/** * gnm_lbeta3: * @a: a number * @b: a number * @sign: (out): the sign * * Returns: the logarithm of the absolute value of the Beta function * evaluated at @a and @b. The sign will be stored in @sign as -1 or * +1. This function is useful because the result of the beta * function can be too large for doubles. */ gnm_float gnm_lbeta3 (gnm_float a, gnm_float b, int *sign) { int sign_a, sign_b, sign_ab; gnm_float ab = a + b; gnm_float res_a, res_b, res_ab; GnmQuad r; int e; switch (qbetaf (a, b, &r, &e)) { case 0: { gnm_float m = gnm_quad_value (&r); *sign = (m >= 0 ? +1 : -1); return gnm_log (gnm_abs (m)) + e * M_LN2gnum; } case 1: /* Overflow */ break; default: *sign = 1; return gnm_nan; } if (a > 0 && b > 0) { *sign = 1; return gnm_lbeta (a, b); } /* This is awful */ res_a = gnm_lgamma_r (a, &sign_a); res_b = gnm_lgamma_r (b, &sign_b); res_ab = gnm_lgamma_r (ab, &sign_ab); *sign = sign_a * sign_b * sign_ab; return res_a + res_b - res_ab; }
static GnmValue * gnumeric_periodogram (GnmFuncEvalInfo *ei, GnmValue const * const *argv) { gnm_float *ord, *absc; int filter, interp; int n0, n1, nb; GnmValue *error = NULL; GnmValue *res; CollectFlags flags; GnmEvalPos const * const ep = ei->pos; GnmValue const * const Pt = argv[0]; int i; GSList *missing0 = NULL, *missing1 = NULL; gnm_complex *in, *out = NULL; int const cols = value_area_get_width (Pt, ep); int const rows = value_area_get_height (Pt, ep); if (cols == 1) nb=rows; else { if (rows == 1) nb=cols; else nb=0; } if (nb == 0) { res = value_new_error_std (ei->pos, GNM_ERROR_VALUE); return res; } flags=COLLECT_IGNORE_BLANKS | COLLECT_IGNORE_STRINGS | COLLECT_IGNORE_BOOLS; ord = collect_floats_value_with_info (argv[0], ei->pos, flags, &n0, &missing0, &error); if (error) { g_slist_free (missing0); return error; } if (n0 == 0) { res = value_new_error_std (ei->pos, GNM_ERROR_VALUE); return res; } if (argv[1]) { filter = (int) gnm_floor (value_get_as_float (argv[1])); if (filter < 0 || filter > FILTER_WELCH) { g_slist_free (missing0); g_free (ord); res = value_new_error_std (ei->pos, GNM_ERROR_VALUE); return res; } } else filter = FILTER_NONE; if (argv[2]) { gnm_float *interpolated, *new_ord, start, incr; int n2; INTERPPROC(interpproc) = NULL; absc = collect_floats_value_with_info (argv[2], ei->pos, flags, &n1, &missing1, &error); if (n1 == 1) { g_slist_free (missing1); g_free (absc); goto no_absc; } if (error) { g_slist_free (missing0); g_slist_free (missing1); g_free (absc); return error; } if (n1 == 0) { g_slist_free (missing0); g_slist_free (missing1); g_free (absc); g_free (ord); return value_new_error_std (ei->pos, GNM_ERROR_VALUE); } if (argv[3]) { interp = (int) gnm_floor (value_get_as_float (argv[3])); if (interp < 0 || interp > INTERPOLATION_SPLINE_AVG) { g_slist_free (missing0); g_slist_free (missing1); g_free (absc); g_free (ord); return error; } } else interp = INTERPOLATION_LINEAR; if (missing0 || missing1) { GSList *missing = gnm_slist_sort_merge (missing0, missing1); gnm_strip_missing (ord, &n0, missing); gnm_strip_missing (absc, &n1, missing); g_slist_free (missing); if (n0 != n1) g_warning ("This should not happen. n0=%d n1=%d\n", n0, n1); } n0 = n1 = MIN (n0, n1); /* here we test if there is abscissas are always increasing, if not, an error is returned */ if (n0 < 2 || !gnm_range_increasing (absc, n0) || n0 == 0) { g_free (absc); g_free (ord); return value_new_error_std (ei->pos, GNM_ERROR_VALUE); } if (argv[4]) { n1 = (int) gnm_floor (value_get_as_float (argv[4])); if (n1 < n0) { g_free (absc); g_free (ord); return value_new_error_std (ei->pos, GNM_ERROR_VALUE); } nb = 1; while (nb < n1) nb *= 2; } else { n1 = 1; while (n1 < n0) n1 *= 2; nb = n1; } incr = (absc[n0 - 1] - absc[0]) / n1; switch (interp) { case INTERPOLATION_LINEAR: interpproc = linear_interpolation; start = absc[0]; n2 = n1; break; case INTERPOLATION_LINEAR_AVG: interpproc = linear_averaging; start = absc[0] - incr / 2.; n2 = n1 + 1; break; case INTERPOLATION_STAIRCASE: interpproc = staircase_interpolation; start = absc[0]; n2 = n1; break; case INTERPOLATION_STAIRCASE_AVG: interpproc = staircase_averaging; start = absc[0] - incr / 2.; n2 = n1 + 1; break; case INTERPOLATION_SPLINE: interpproc = spline_interpolation; start = absc[0]; n2 = n1; break; case INTERPOLATION_SPLINE_AVG: interpproc = spline_averaging; start = absc[0] - incr / 2.; n2 = n1 + 1; break; default: g_free (absc); g_free (ord); return value_new_error_std (ei->pos, GNM_ERROR_NA); } interpolated = g_new (gnm_float, n1); for (i = 0; i < n2; i++) interpolated[i] = start + i * incr; new_ord = interpproc (absc, ord, n0, interpolated, n1); g_free (ord); ord = new_ord; if (ord == NULL) { g_free (absc); g_free (interpolated); return value_new_error_std (ei->pos, GNM_ERROR_NA); } n0 = n1; } else { no_absc: /* we have no interpolation to apply, so just take the values */ if (missing0) { gnm_strip_missing (ord, &n0, missing0); g_slist_free (missing0); } nb = 1; while (nb < n0) nb *= 2; } /* Now apply the filter if any */ if (filter != FILTER_NONE) { gnm_float factor; switch (filter) { case FILTER_BARTLETT: factor = n0 / 2.; for (i = 0; i < n0; i++) ord[i] *= 1. - gnm_abs ((i / factor - 1)); break; case FILTER_HANN: factor = 2. * M_PIgnum / n0; for (i = 0; i < n0; i++) ord[i] *= 0.5 * (1 - gnm_cos (factor * i)); break; case FILTER_WELCH: factor = n0 / 2.; for (i = 0; i < n0; i++) ord[i] *= 1. - (i / factor - 1.) * (i / factor - 1.); break; } } /* Transform and return the result */ in = g_new0 (gnm_complex, nb); for (i = 0; i < n0; i++){ in[i].re = ord[i]; } g_free (ord); gnm_fourier_fft (in, nb, 1, &out, FALSE); g_free (in); nb /= 2; if (out && nb > 0) { res = value_new_array_non_init (1 , nb); res->v_array.vals[0] = g_new (GnmValue *, nb); for (i = 0; i < nb; i++) res->v_array.vals[0][i] = value_new_float (gnm_sqrt ( out[i].re * out[i].re + out[i].im * out[i].im)); g_free (out); } else
gnm_float pochhammer (gnm_float x, gnm_float n) { gnm_float rn, rx, lr; GnmQuad m1, m2; int e1, e2; if (gnm_isnan (x) || gnm_isnan (n)) return gnm_nan; if (n == 0) return 1; rx = gnm_floor (x); rn = gnm_floor (n); /* * Use naive multiplication when n is a small integer. * We don't want to use this if x is also an integer * (but we might do so below if x is insanely large). */ if (n == rn && x != rx && n >= 0 && n < 40) return pochhammer_naive (x, (int)n); if (!qfactf (x + n - 1, &m1, &e1) && !qfactf (x - 1, &m2, &e2)) { void *state = gnm_quad_start (); int de = e1 - e2; GnmQuad qr; gnm_float r; gnm_quad_div (&qr, &m1, &m2); r = gnm_quad_value (&qr); gnm_quad_end (state); return gnm_ldexp (r, de); } if (x == rx && x <= 0) { if (n != rn) return 0; if (x == 0) return (n > 0) ? 0 : ((gnm_fmod (-n, 2) == 0 ? +1 : -1) / gnm_fact (-n)); if (n > -x) return gnm_nan; } /* * We have left the common cases. One of x+n and x is * insanely big, possibly both. */ if (gnm_abs (x) < 1) return gnm_pinf; if (n < 0) return 1 / pochhammer (x + n, -n); if (n == rn && n >= 0 && n < 100) return pochhammer_naive (x, (int)n); if (gnm_abs (n) < 1) { /* x is big. */ void *state = gnm_quad_start (); GnmQuad qr; gnm_float r; pochhammer_small_n (x, n, &qr); r = gnm_quad_value (&qr); gnm_quad_end (state); return r; } /* Panic mode. */ g_printerr ("x=%.20g n=%.20g\n", x, n); lr = ((x - 0.5) * gnm_log1p (n / x) + n * gnm_log (x + n) - n + (lgammacor (x + n) - lgammacor (x))); return gnm_exp (lr); }
/* Compute gnm_log(gamma(a+1)) accurately also for small a (0 < a < 0.5). */ gnm_float lgamma1p (gnm_float a) { const gnm_float eulers_const = GNM_const(0.5772156649015328606065120900824024); /* coeffs[i] holds (zeta(i+2)-1)/(i+2) , i = 1:N, N = 40 : */ const int N = 40; static const gnm_float coeffs[40] = { GNM_const(0.3224670334241132182362075833230126e-0), GNM_const(0.6735230105319809513324605383715000e-1), GNM_const(0.2058080842778454787900092413529198e-1), GNM_const(0.7385551028673985266273097291406834e-2), GNM_const(0.2890510330741523285752988298486755e-2), GNM_const(0.1192753911703260977113935692828109e-2), GNM_const(0.5096695247430424223356548135815582e-3), GNM_const(0.2231547584535793797614188036013401e-3), GNM_const(0.9945751278180853371459589003190170e-4), GNM_const(0.4492623673813314170020750240635786e-4), GNM_const(0.2050721277567069155316650397830591e-4), GNM_const(0.9439488275268395903987425104415055e-5), GNM_const(0.4374866789907487804181793223952411e-5), GNM_const(0.2039215753801366236781900709670839e-5), GNM_const(0.9551412130407419832857179772951265e-6), GNM_const(0.4492469198764566043294290331193655e-6), GNM_const(0.2120718480555466586923135901077628e-6), GNM_const(0.1004322482396809960872083050053344e-6), GNM_const(0.4769810169363980565760193417246730e-7), GNM_const(0.2271109460894316491031998116062124e-7), GNM_const(0.1083865921489695409107491757968159e-7), GNM_const(0.5183475041970046655121248647057669e-8), GNM_const(0.2483674543802478317185008663991718e-8), GNM_const(0.1192140140586091207442548202774640e-8), GNM_const(0.5731367241678862013330194857961011e-9), GNM_const(0.2759522885124233145178149692816341e-9), GNM_const(0.1330476437424448948149715720858008e-9), GNM_const(0.6422964563838100022082448087644648e-10), GNM_const(0.3104424774732227276239215783404066e-10), GNM_const(0.1502138408075414217093301048780668e-10), GNM_const(0.7275974480239079662504549924814047e-11), GNM_const(0.3527742476575915083615072228655483e-11), GNM_const(0.1711991790559617908601084114443031e-11), GNM_const(0.8315385841420284819798357793954418e-12), GNM_const(0.4042200525289440065536008957032895e-12), GNM_const(0.1966475631096616490411045679010286e-12), GNM_const(0.9573630387838555763782200936508615e-13), GNM_const(0.4664076026428374224576492565974577e-13), GNM_const(0.2273736960065972320633279596737272e-13), GNM_const(0.1109139947083452201658320007192334e-13) }; const gnm_float c = GNM_const(0.2273736845824652515226821577978691e-12);/* zeta(N+2)-1 */ gnm_float lgam; int i; if (gnm_abs (a) >= 0.5) return gnm_lgamma (a + 1); /* Abramowitz & Stegun 6.1.33, * also http://functions.wolfram.com/06.11.06.0008.01 */ lgam = c * gnm_logcf (-a / 2, N + 2, 1); for (i = N - 1; i >= 0; i--) lgam = coeffs[i] - a * lgam; return (a * lgam - eulers_const) * a - log1pmx (a); } /* lgamma1p */
static gboolean rosenbrock_iter (GnmNlsolve *nl) { GnmSolver *sol = nl->parent; const int n = nl->vars->len; int i, j; const gnm_float alpha = 3; const gnm_float beta = 0.5; gboolean any_at_all = FALSE; gnm_float *d, **A, *x, *dx, *t; char *state; int dones = 0; gnm_float ykm1 = nl->yk, *xkm1; gnm_float eps = gnm_pow2 (-16); int safety = 0; if (nl->tentative) { nl->tentative--; if (nl->tentative == 0) { if (nl->debug) g_printerr ("Tentative move rejected\n"); rosenbrock_tentative_end (nl, FALSE); } } if (nl->k % 20 == 0) { for (i = 0; i < n; i++) for (j = 0; j < n; j++) nl->xi[i][j] = (i == j); } A = g_new (gnm_float *, n); for (i = 0; i < n; i++) A[i] = g_new (gnm_float, n); dx = g_new (gnm_float, n); for (i = 0; i < n; i++) dx[i] = 0; x = g_new (gnm_float, n); t = g_new (gnm_float, n); d = g_new (gnm_float, n); for (i = 0; i < n; i++) { d[i] = (nl->xk[i] == 0) ? eps : gnm_abs (nl->xk[i]) * eps; } xkm1 = g_memdup (nl->xk, n * sizeof (gnm_float)); state = g_new0 (char, n); while (dones < n) { /* * A safety that shouldn't get hit, but might if the function * being optimized is non-deterministic. */ if (safety++ > n * GNM_MANT_DIG) break; for (i = 0; i < n; i++) { gnm_float y; if (state[i] == 2) continue; /* x = xk + (d[i] * xi[i]) */ for (j = 0; j < n; j++) x[j] = nl->xk[j] + d[i] * nl->xi[i][j]; set_vector (nl, x); y = get_value (nl); if (y <= nl->yk && gnm_solver_check_constraints (sol)) { if (y < nl->yk) { nl->yk = y; memcpy (nl->xk, x, n * sizeof (gnm_float)); dx[i] += d[i]; any_at_all = TRUE; } switch (state[i]) { case 0: state[i] = 1; /* Fall through */ case 1: d[i] *= alpha; break; default: case 2: break; } } else { switch (state[i]) { case 1: state[i] = 2; dones++; /* Fall through */ case 0: d[i] *= -beta; break; default: case 2: /* No sign change. */ d[i] *= 0.5; break; } } } } if (any_at_all) { gnm_float div, sum; for (j = n - 1; j >= 0; j--) for (i = 0; i < n; i++) A[j][i] = (j == n - 1 ? 0 : A[j + 1][i]) + dx[j] * nl->xi[j][i]; sum = 0; for (i = n - 1; i >= 0; i--) { sum += dx[i] * dx[i]; t[i] = sum; } for (i = n - 1; i > 0; i--) { div = gnm_sqrt (t[i - 1] * t[i]); if (div != 0) for (j = 0; j < n; j++) { nl->xi[i][j] = (dx[i - 1] * A[i][j] - nl->xi[i - 1][j] * t[i]) / div; g_assert (gnm_finite (nl->xi[i][j])); } } gnm_range_hypot (dx, n, &div); if (div != 0) { for (i = 0; i < n; i++) { nl->xi[0][i] = A[0][i] / div; if (!gnm_finite (nl->xi[0][i])) { g_printerr ("%g %g %g\n", div, A[0][i], nl->xi[0][i]); g_assert (gnm_finite (nl->xi[0][i])); } } } /* ---------------------------------------- */ if (!nl->tentative) { set_vector (nl, nl->xk); gnm_nlsolve_set_solution (nl); } if (nl->tentative) { if (nl->yk < nl->tentative_yk) { if (nl->debug) g_printerr ("Tentative move accepted!\n"); rosenbrock_tentative_end (nl, TRUE); } } else if (gnm_abs (nl->yk - ykm1) > gnm_abs (ykm1) * 0.01) { /* A big step. */ nl->smallsteps = 0; } else { nl->smallsteps++; } if (!nl->tentative && nl->smallsteps > 50) { gnm_float yk = nl->yk; nl->tentative = 10; nl->tentative_xk = g_memdup (nl->xk, n * sizeof (gnm_float)); nl->tentative_yk = yk; for (i = 0; i < 4; i++) { gnm_float ymax = yk + gnm_abs (yk) * (0.10 / (i + 1)); if (i > 0) ymax = MIN (ymax, nl->yk); if (!newton_improve (nl, nl->xk, &nl->yk, ymax)) break; } if (nl->debug) print_vector ("Tentative move to", nl->xk, n); } } g_free (x); g_free (xkm1); g_free (dx); g_free (t); g_free (d); free_matrix (A, n); g_free (state); return any_at_all; }
static gboolean glpk_affine_func (GString *dst, GnmCell *target, GnmSubSolver *ssol, gnm_float const *x1, gnm_float const *x2, gboolean zero_too, gnm_float cst, GError **err) { GnmSolver *sol = GNM_SOLVER (ssol); unsigned ui; gboolean any = FALSE; gnm_float y; gboolean ok = TRUE; GPtrArray *input_cells = sol->input_cells; gnm_float *cs; if (!target) { gnm_string_add_number (dst, cst); return TRUE; } gnm_solver_set_vars (sol, x1); gnm_cell_eval (target); y = cst + value_get_as_float (target->value); cs = gnm_solver_get_lp_coeffs (sol, target, x1, x2, err); if (!cs) goto fail; /* Adjust constant for choice of x1. */ for (ui = 0; ui < input_cells->len; ui++) y -= x1[ui] * cs[ui]; for (ui = 0; ui < input_cells->len; ui++) { GnmCell *cell = g_ptr_array_index (input_cells, ui); gnm_float x = cs[ui]; if (x == 0 && !zero_too) continue; if (any) { if (x < 0) g_string_append (dst, " - "); else g_string_append (dst, " + "); } else { if (x < 0) g_string_append_c (dst, '-'); } x = gnm_abs (x); if (x != 1) { gnm_string_add_number (dst, x); g_string_append_c (dst, ' '); } g_string_append (dst, glpk_var_name (ssol, cell)); any = TRUE; } if (!any || y) { if (any) { g_string_append_c (dst, ' '); if (y > 0) g_string_append_c (dst, '+'); } gnm_string_add_number (dst, y); } fail: g_free (cs); return ok; }
static GnmValue * gnumeric_randdiscrete (GnmFuncEvalInfo *ei, GnmValue const * const *argv) { GnmValue *res = NULL; gnm_float *values = NULL; gnm_float *probs = NULL; int nv, np, i; gnm_float p; values = collect_floats_value (argv[0], ei->pos, COLLECT_IGNORE_STRINGS | COLLECT_IGNORE_BOOLS | COLLECT_IGNORE_BLANKS, &nv, &res); if (res) goto out; if (argv[1]) { probs = collect_floats_value (argv[1], ei->pos, COLLECT_IGNORE_STRINGS | COLLECT_IGNORE_BOOLS | COLLECT_IGNORE_BLANKS, &np, &res); if (res) goto out; } else np = nv; if (nv < 1 || nv != np) goto error; if (probs) { gnm_float pmin, psum; gnm_range_min (probs, np, &pmin); if (pmin < 0) goto error; gnm_range_sum (probs, np, &psum); if (gnm_abs (psum - 1) > 1e-10) goto error; } p = random_01 (); if (probs) { for (i = 0; i < np; i++) { p -= probs[i]; if (p < 0) break; } } else { /* Uniform. */ i = (int)gnm_floor (p * nv); } /* MIN is needed because of the sum grace. */ res = value_new_float (values[MIN (i, nv - 1)]); out: g_free (values); g_free (probs); return res; error: res = value_new_error_NUM (ei->pos); goto out; }