Blame multifit_nlinear/test_brown1.c

Packit 67cb25
#define brown1_N         20
Packit 67cb25
#define brown1_P         4
Packit 67cb25
Packit 67cb25
static double brown1_x0[brown1_P] = { 25, 5, -5, -1 };
Packit 67cb25
static double brown1_epsrel = 1.0e-5;
Packit 67cb25
Packit 67cb25
static void
Packit 67cb25
brown1_checksol(const double x[], const double sumsq,
Packit 67cb25
                const double epsrel, const char *sname,
Packit 67cb25
                const char *pname)
Packit 67cb25
{
Packit 67cb25
  size_t i;
Packit 67cb25
  const double sumsq_exact = 8.582220162635628e+04;
Packit 67cb25
  const double brown1_x[brown1_P] = {
Packit 67cb25
    -1.159443990239263e+01, 1.320363005221244e+01,
Packit 67cb25
    -4.034395456782477e-01, 2.367789088597534e-01 };
Packit 67cb25
Packit 67cb25
  gsl_test_rel(sumsq, sumsq_exact, epsrel, "%s/%s sumsq",
Packit 67cb25
               sname, pname);
Packit 67cb25
Packit 67cb25
  for (i = 0; i < brown1_P; ++i)
Packit 67cb25
    {
Packit 67cb25
      gsl_test_rel(x[i], brown1_x[i], epsrel, "%s/%s i=%zu",
Packit 67cb25
                   sname, pname, i);
Packit 67cb25
    }
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
static int
Packit 67cb25
brown1_f (const gsl_vector * x, void *params, gsl_vector * f)
Packit 67cb25
{
Packit 67cb25
  double x0 = gsl_vector_get (x, 0);
Packit 67cb25
  double x1 = gsl_vector_get (x, 1);
Packit 67cb25
  double x2 = gsl_vector_get (x, 2);
Packit 67cb25
  double x3 = gsl_vector_get (x, 3);
Packit 67cb25
  size_t i;
Packit 67cb25
Packit 67cb25
  for (i = 0; i < brown1_N; i++)
Packit 67cb25
    {
Packit 67cb25
      double ti = 0.2 * (i + 1);
Packit 67cb25
      double ui = x0 + x1 * ti - exp (ti);
Packit 67cb25
      double vi = x2 + x3 * sin (ti) - cos (ti);
Packit 67cb25
Packit 67cb25
      gsl_vector_set (f, i, ui * ui + vi * vi);
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
  (void)params; /* avoid unused parameter warning */
Packit 67cb25
Packit 67cb25
  return GSL_SUCCESS;
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
static int
Packit 67cb25
brown1_df (const gsl_vector * x, void *params, gsl_matrix * df)
Packit 67cb25
{
Packit 67cb25
  double x0 = gsl_vector_get (x, 0);
Packit 67cb25
  double x1 = gsl_vector_get (x, 1);
Packit 67cb25
  double x2 = gsl_vector_get (x, 2);
Packit 67cb25
  double x3 = gsl_vector_get (x, 3);
Packit 67cb25
  size_t i;
Packit 67cb25
Packit 67cb25
  for (i = 0; i < brown1_N; i++)
Packit 67cb25
    {
Packit 67cb25
      double ti = 0.2 * (i + 1);
Packit 67cb25
      double ui = x0 + x1 * ti - exp (ti);
Packit 67cb25
      double vi = x2 + x3 * sin (ti) - cos (ti);
Packit 67cb25
Packit 67cb25
      gsl_matrix_set (df, i, 0, 2 * ui);
Packit 67cb25
      gsl_matrix_set (df, i, 1, 2 * ui * ti);
Packit 67cb25
      gsl_matrix_set (df, i, 2, 2 * vi);
Packit 67cb25
      gsl_matrix_set (df, i, 3, 2 * vi * sin (ti));
Packit 67cb25
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
  (void)params; /* avoid unused parameter warning */
Packit 67cb25
Packit 67cb25
  return GSL_SUCCESS;
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
static int
Packit 67cb25
brown1_fvv (const gsl_vector * x, const gsl_vector * v,
Packit 67cb25
            void *params, gsl_vector * fvv)
Packit 67cb25
{
Packit 67cb25
  double v0 = gsl_vector_get (v, 0);
Packit 67cb25
  double v1 = gsl_vector_get (v, 1);
Packit 67cb25
  double v2 = gsl_vector_get (v, 2);
Packit 67cb25
  double v3 = gsl_vector_get (v, 3);
Packit 67cb25
  size_t i;
Packit 67cb25
Packit 67cb25
  for (i = 0; i < brown1_N; i++)
Packit 67cb25
    {
Packit 67cb25
      double ti = 0.2 * (i + 1);
Packit 67cb25
      double term1 = v0 + ti*v1;
Packit 67cb25
      double term2 = v3*sin(ti);
Packit 67cb25
Packit 67cb25
      gsl_vector_set (fvv, i, 2.0 * (term1*term1 + v2*v2 +
Packit 67cb25
                                     term2 * (2*v2 + term2)));
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
  (void)x;      /* avoid unused parameter warning */
Packit 67cb25
  (void)params; /* avoid unused parameter warning */
Packit 67cb25
Packit 67cb25
  return GSL_SUCCESS;
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
static gsl_multifit_nlinear_fdf brown1_func =
Packit 67cb25
{
Packit 67cb25
  brown1_f,
Packit 67cb25
  brown1_df,
Packit 67cb25
  brown1_fvv,
Packit 67cb25
  brown1_N,
Packit 67cb25
  brown1_P,
Packit 67cb25
  NULL,
Packit 67cb25
  0,
Packit 67cb25
  0,
Packit 67cb25
  0
Packit 67cb25
};
Packit 67cb25
Packit 67cb25
static test_fdf_problem brown1_problem =
Packit 67cb25
{
Packit 67cb25
  "brown_dennis",
Packit 67cb25
  brown1_x0,
Packit 67cb25
  NULL,
Packit 67cb25
  NULL,
Packit 67cb25
  &brown1_epsrel,
Packit 67cb25
  &brown1_checksol,
Packit 67cb25
  &brown1_func
Packit 67cb25
};