void magcal_free(magcal_workspace *w) { if (w->Ex) free(w->Ex); if (w->Ey) free(w->Ey); if (w->Ez) free(w->Ez); if (w->F) free(w->F); if (w->fdf_s) gsl_multifit_fdfsolver_free(w->fdf_s); if (w->fdf_ridge) gsl_multifit_fdfridge_free(w->fdf_ridge); if (w->covar) gsl_matrix_free(w->covar); free(w); }
static int magcal_init(const satdata_mag *data, magcal_workspace *w) { int s = 0; size_t i; size_t n = 0; for (i = 0; i < data->n; ++i) { /* don't store flagged data */ if (data->flags[i]) continue; /* don't process high latitude data */ if (fabs(data->latitude[i]) > MAGCAL_MAX_LATITUDE) continue; w->Ex[n] = SATDATA_VEC_X(data->B_VFM, i); w->Ey[n] = SATDATA_VEC_Y(data->B_VFM, i); w->Ez[n] = SATDATA_VEC_Z(data->B_VFM, i); w->F[n] = data->F[i]; ++n; } if (n < 200) { fprintf(stderr, "magcal_init: insufficient data points for calibration: %zu\n", n); return -1; } if (n != w->n) { gsl_multifit_fdfsolver_free(w->fdf_s); gsl_multifit_fdfridge_free(w->fdf_ridge); w->fdf_s = gsl_multifit_fdfsolver_alloc(w->fdf_type, n, w->p); w->fdf_ridge = gsl_multifit_fdfridge_alloc(w->fdf_type, n, w->p); w->n = n; } #if MAGCAL_SCALE w->B_s = GSL_MAX(gsl_stats_sd(w->Ex, 1, n), GSL_MAX(gsl_stats_sd(w->Ey, 1, n), gsl_stats_sd(w->Ez, 1, n))); #endif /* center and scale data arrays */ for (i = 0; i < n; ++i) { w->Ex[i] /= w->B_s; w->Ey[i] /= w->B_s; w->Ez[i] /= w->B_s; w->F[i] /= w->B_s; } return s; } /* magcal_init() */
gsl_multifit_fdfridge * gsl_multifit_fdfridge_alloc (const gsl_multifit_fdfsolver_type * T, const size_t n, const size_t p) { gsl_multifit_fdfridge * work; work = calloc(1, sizeof(gsl_multifit_fdfridge)); if (work == NULL) { GSL_ERROR_VAL("failed to allocate workspace", GSL_ENOMEM, 0); } work->s = gsl_multifit_fdfsolver_alloc (T, n + p, p); if (work->s == NULL) { gsl_multifit_fdfridge_free(work); GSL_ERROR_VAL("failed to allocate space for fdfsolver", GSL_ENOMEM, 0); } work->wts = gsl_vector_alloc(n + p); if (work->wts == NULL) { gsl_multifit_fdfridge_free(work); GSL_ERROR_VAL("failed to allocate space for weight vector", GSL_ENOMEM, 0); } work->n = n; work->p = p; work->lambda = 0.0; /* initialize weights to 1 (for augmented portion of vector) */ gsl_vector_set_all(work->wts, 1.0); return work; } /* gsl_multifit_fdfridge_alloc() */
static void test_fdfridge(const gsl_multifit_fdfsolver_type * T, const double xtol, const double gtol, const double ftol, const double epsrel, const double x0_scale, test_fdf_problem *problem, const double *wts) { gsl_multifit_function_fdf *fdf = problem->fdf; const size_t n = fdf->n; const size_t p = fdf->p; const size_t max_iter = 1500; gsl_vector *x0 = gsl_vector_alloc(p); gsl_vector_view x0v = gsl_vector_view_array(problem->x0, p); gsl_multifit_fdfridge *w = gsl_multifit_fdfridge_alloc (T, n, p); const char *pname = problem->name; char sname[2048]; int status, info; double lambda = 0.0; sprintf(sname, "ridge/%s", gsl_multifit_fdfridge_name(w)); /* scale starting point x0 */ gsl_vector_memcpy(x0, &x0v.vector); test_scale_x0(x0, x0_scale); /* test undamped case with lambda = 0 */ if (wts) { gsl_vector_const_view wv = gsl_vector_const_view_array(wts, n); gsl_multifit_fdfridge_wset(w, fdf, x0, lambda, &wv.vector); } else gsl_multifit_fdfridge_set(w, fdf, x0, lambda); status = gsl_multifit_fdfridge_driver(w, max_iter, xtol, gtol, ftol, &info); gsl_test(status, "%s/%s did not converge, status=%s", sname, pname, gsl_strerror(status)); /* check solution */ test_fdf_checksol(sname, pname, epsrel, w->s, problem); /* test for self consisent solution with L = \lambda I */ { const double eps = 1.0e-10; gsl_matrix *L = gsl_matrix_calloc(p, p); gsl_vector_view diag = gsl_matrix_diagonal(L); gsl_multifit_fdfridge *w2 = gsl_multifit_fdfridge_alloc (T, n, p); gsl_vector *y0 = gsl_vector_alloc(p); size_t i; /* pick some value for lambda and set L = \lambda I */ lambda = 5.0; gsl_vector_set_all(&diag.vector, lambda); /* scale initial vector */ gsl_vector_memcpy(x0, &x0v.vector); test_scale_x0(x0, x0_scale); gsl_vector_memcpy(y0, x0); if (wts) { gsl_vector_const_view wv = gsl_vector_const_view_array(wts, n); gsl_multifit_fdfridge_wset(w, fdf, x0, lambda, &wv.vector); gsl_multifit_fdfridge_wset3(w2, fdf, y0, L, &wv.vector); } else { gsl_multifit_fdfridge_set(w, fdf, x0, lambda); gsl_multifit_fdfridge_set3(w2, fdf, y0, L); } /* solve with scalar lambda routine */ status = gsl_multifit_fdfridge_driver(w, max_iter, xtol, gtol, ftol, &info); gsl_test(status, "%s/lambda/%s did not converge, status=%s", sname, pname, gsl_strerror(status)); /* solve with general matrix routine */ status = gsl_multifit_fdfridge_driver(w2, max_iter, xtol, gtol, ftol, &info); gsl_test(status, "%s/L/%s did not converge, status=%s", sname, pname, gsl_strerror(status)); /* test x = y */ for (i = 0; i < p; ++i) { double xi = gsl_vector_get(w->s->x, i); double yi = gsl_vector_get(w2->s->x, i); if (fabs(xi) < eps) { gsl_test_abs(yi, xi, eps, "%s/%s ridge lambda=%g i="F_ZU, sname, pname, lambda, i); } else { gsl_test_rel(yi, xi, eps, "%s/%s ridge lambda=%g i="F_ZU, sname, pname, lambda, i); } } gsl_matrix_free(L); gsl_vector_free(y0); gsl_multifit_fdfridge_free(w2); } gsl_multifit_fdfridge_free(w); gsl_vector_free(x0); }