30#include <gsl/gsl_linalg.h>
31#include <gsl/gsl_matrix.h>
32#include <gsl/gsl_bspline.h>
33#include <gsl/gsl_sort_vector.h>
34#include <gsl/gsl_multifit.h>
35#include <gsl/gsl_multimin.h>
36#include <gsl/gsl_statistics.h>
47 const size_t N = a->size;
48 const size_t stride = a->stride;
52 for (i = 0; i < N; i++)
53 sum += a->data[i * stride];
104 FREE(gsl_matrix_int_free, self->
mask);
106 FREE(gsl_matrix_free, self->
Vt);
107 FREE(gsl_vector_free, self->
S);
138static int svd_wrapper(gsl_matrix **A, gsl_matrix **V, gsl_vector **S)
140 int m = (*A)->size1, n = (*A)->size2, ret;
144 *V = gsl_matrix_alloc(n, n);
145 *S = gsl_vector_alloc(n);
146 w = gsl_vector_alloc(n);
147 ret = gsl_linalg_SV_decomp(*A , *V, *S, w);
150 gsl_matrix *At = gsl_matrix_alloc(n, m);
151 gsl_matrix_transpose_memcpy(At, *A);
153 *V = gsl_matrix_alloc(m, m);
154 *S = gsl_vector_alloc(m);
155 w = gsl_vector_alloc(m);
156 ret = gsl_linalg_SV_decomp(At , *V, *S, w);
178 gsl_vector *signs = gsl_vector_alloc(U->size2);
180 for (
size_t j = 0; j < U->size2; j++) {
181 gsl_vector_const_view Ucol = gsl_matrix_const_column(U, j);
182 double max_abs_val = -INFINITY;
183 size_t max_abs_ind = -1;
184 for(
size_t j = 0; j < U->size2; j++) {
185 double abs_val = fabs(gsl_vector_get(&Ucol.vector, j));
186 if (abs_val > max_abs_val) max_abs_ind = j;
189 double max_val = gsl_vector_get(&Ucol.vector, max_abs_ind);
190 gsl_vector_set(signs, j, max_val > 0 ? 1 : -1);
194 for (
size_t i = 0; i < U->size1; i++) {
195 gsl_vector_view Urow = gsl_matrix_row(U, i);
196 gsl_vector_mul(&Urow.vector, signs);
200 for (
size_t i = 0; i < V->size1; i++) {
201 gsl_vector_view Vrow = gsl_matrix_row(V, i);
202 gsl_vector_mul(&Vrow.vector, signs);
205 gsl_vector_free(signs);
221 cpl_ensure(data != NULL, CPL_ERROR_NULL_INPUT, NULL);
223 cpl_size n_obs = cpl_matrix_get_nrow(data);
224 cpl_size n_wave = cpl_matrix_get_ncol(data);
225 cpl_size nvalid = CPL_MIN(n_obs, n_wave);
228 gsl_matrix_const_view data_gsl_view = gsl_matrix_const_view_array(
229 cpl_matrix_get_data_const(data), n_obs, n_wave);
231 gsl_matrix_int *use_mask = gsl_matrix_int_alloc(n_obs, n_wave);
233 for(cpl_size i = 0; i < n_obs; i++) {
234 for (cpl_size j = 0; j < n_wave; j++)
235 gsl_matrix_int_set(use_mask, i, j,
236 cpl_matrix_get(mask, i, j) > 0 ? 0 : 1);
239 gsl_matrix_int_set_all(use_mask, 1);
243 gsl_vector* column_mean = gsl_vector_alloc(n_wave);
244 for(
int j = 0; j < n_wave; j++) {
245 gsl_vector_const_view col = gsl_matrix_const_column(&data_gsl_view.matrix, j);
247 gsl_vector_set(column_mean, j, mu);
251 gsl_matrix* mean_subtracted_data = gsl_matrix_alloc(n_obs, n_wave);
252 gsl_matrix_memcpy(mean_subtracted_data, &data_gsl_view.matrix);
253 for(cpl_size i = 0; i < n_obs; i++) {
254 gsl_vector_view mean_subtracted_row_view = gsl_matrix_row(mean_subtracted_data, i);
255 gsl_vector_sub(&mean_subtracted_row_view.vector, column_mean);
259 result->
n_obs = n_obs;
262 result->
data = mean_subtracted_data;
263 result->
mask = use_mask;
264 result->
mean = column_mean;
272 result->
signs = NULL;
285 cpl_ensure_code(self, CPL_ERROR_NULL_INPUT);
287 unsigned int n_obs = self->
n_obs;
288 unsigned int n_wave = self->
n_wave;
289 unsigned int n_valid = self->
n_valid;
295 gsl_matrix *U = gsl_matrix_alloc(n_obs, n_wave);
296 gsl_matrix_memcpy(U, self->
data);
300 return CPL_ERROR_ILLEGAL_INPUT;
322 self->
Vt = gsl_matrix_alloc(n_valid, n_wave);
323 gsl_matrix_transpose_memcpy(self->
Vt, V);
326 return CPL_ERROR_NONE;
340 const double ranges[2][2] = {{2.02, 2.04}, {2.34, 2.36}};
341 int range_indices[2][2];
344 int nwave = wave->size, irange = 0;
345 for (
int w = 0; w < nwave; w++) {
346 if (gsl_vector_get(wave, w) > ranges[irange / 2][irange % 2]) {
347 range_indices[irange / 2][irange % 2] = w;
354 int range_length[2] = {
355 range_indices[0][1] - range_indices[0][0],
356 range_indices[1][1] - range_indices[1][0]
358 int total_length = range_length[0] + range_length[1];
362 gsl_vector *visphi_range = gsl_vector_alloc(total_length);
363 for (
int r = 0; r < 2; r++) {
364 gsl_vector_const_view range = gsl_vector_const_subvector(
365 visphi, range_indices[r][0], range_length[r]);
366 gsl_vector_view range_dest = gsl_vector_subvector(
367 visphi_range, offset, range_length[r]);
368 gsl_vector_memcpy(&range_dest.vector, &range.vector);
369 offset += range_length[r];
374 gsl_sort_vector(visphi_range);
375 double result = gsl_stats_median_from_sorted_data(
376 visphi_range->data, 1, visphi_range->size
379 gsl_vector_free(visphi_range);
395 cpl_ensure_code(wave, CPL_ERROR_NULL_INPUT);
396 cpl_ensure_code(ncomp > 0 && ncomp <= self->n_valid, CPL_ERROR_ILLEGAL_INPUT);
398 const int n_wave = self->
n_wave;
399 gsl_vector_const_view wave_vec = gsl_vector_const_view_array(
400 cpl_array_get_data_double_const(wave), n_wave);
401 gsl_vector *signs = gsl_vector_alloc(ncomp);
403 for (
int c = 0; c < ncomp; c++) {
404 gsl_vector_const_view components = gsl_matrix_const_row(self->
Vt, c);
406 gsl_vector_set(signs, c, sign);
410 return CPL_ERROR_NONE;
425 cpl_ensure(results, CPL_ERROR_NULL_INPUT, NULL);
426 cpl_ensure(ncomp > 0 && ncomp <= results[0]->n_valid, CPL_ERROR_ILLEGAL_INPUT, NULL);
427 for (
int r = 0; r < nres; r++)
428 cpl_ensure(results[r]->signs, CPL_ERROR_ILLEGAL_INPUT, NULL);
430 const int n_wave = results[0]->
n_wave;
432 gsl_matrix *all_comp_vals = gsl_matrix_alloc(nres, n_wave);
433 gsl_matrix *median_comp_vals = gsl_matrix_alloc(ncomp, n_wave);
434 gsl_matrix *flip_vals = gsl_matrix_alloc(nres, n_wave);
435 gsl_vector *flip_vec = gsl_vector_alloc(n_wave);
437 for (
int c = 0; c < ncomp; c++) {
439 for (
int r = 0; r < nres; r++) {
440 gsl_vector_const_view rth_components = gsl_matrix_const_row(results[r]->Vt, c);
441 gsl_matrix_set_row(all_comp_vals, r, &rth_components.vector);
444 double sign = gsl_vector_get(results[r]->signs, c);
445 gsl_vector_set_all(flip_vec, sign);
446 gsl_matrix_set_row(flip_vals, r, flip_vec);
450 gsl_matrix_mul_elements(all_comp_vals, flip_vals);
453 for (
int j = 0; j < n_wave; j++) {
454 gsl_vector_view column_comp_vals = gsl_matrix_column(all_comp_vals, j);
455 gsl_sort_vector(&column_comp_vals.vector);
456 double median_val = gsl_stats_median_from_sorted_data(
457 column_comp_vals.vector.data, n_wave, nres
459 gsl_matrix_set(median_comp_vals, c, j, median_val);
462 gsl_matrix_free(all_comp_vals);
463 gsl_matrix_free(flip_vals);
464 gsl_vector_free(flip_vec);
485 int n_comp = cpl_matrix_get_nrow(components);
486 int n_wave = cpl_matrix_get_ncol(components);
489 gsl_matrix_const_view comps_gsl_view = gsl_matrix_const_view_array(
490 cpl_matrix_get_data_const(components), n_comp, n_wave);
491 gsl_matrix *my_comps = gsl_matrix_alloc(n_comp, n_wave);
492 gsl_matrix_memcpy(my_comps, &comps_gsl_view.matrix);
513 cpl_ensure_code(self, CPL_ERROR_NULL_INPUT);
516 cpl_ensure_code(wave, CPL_ERROR_NULL_INPUT);
517 cpl_ensure_code(degree > 0, CPL_ERROR_ILLEGAL_INPUT);
518 cpl_ensure_code(ncomp > 0 && ncomp <= self->n_comp, CPL_ERROR_ILLEGAL_INPUT);
520 const int n_coeffs = degree;
522 const int nbreak = n_coeffs - k + 2;
523 const int n_wave = self->
n_wave;
525 gsl_bspline_workspace *bw = gsl_bspline_alloc(k, nbreak);
526 const gsl_multifit_robust_type *T = gsl_multifit_robust_default;
527 gsl_multifit_robust_workspace *mw = gsl_multifit_robust_alloc(T, n_wave, n_coeffs);
528 gsl_multifit_robust_maxiter(1000, mw);
530 gsl_vector *B = gsl_vector_alloc(n_coeffs);
531 gsl_vector *coeffs = gsl_vector_alloc(n_coeffs);
532 gsl_matrix *cov = gsl_matrix_alloc(n_coeffs, n_coeffs);
533 double wvi, Bj, fit_val, fit_err;
535 gsl_bspline_knots_uniform(cpl_array_get(wave, 0, NULL), cpl_array_get(wave, n_wave-1, NULL), bw);
537 gsl_matrix *result = gsl_matrix_alloc(ncomp, n_wave);
540 for (
int c = 0; c < ncomp; c++) {
544 gsl_matrix *X = gsl_matrix_alloc(n_wave, n_coeffs);
545 for (
int i = 0; i < n_wave; i++) {
546 wvi = cpl_array_get(wave, i, NULL);
547 gsl_bspline_eval(wvi, B, bw);
549 for (
int j = 0; j < n_coeffs; j++) {
550 Bj = gsl_vector_get(B, j);
551 gsl_matrix_set(X, i, j, Bj);
556 gsl_multifit_robust(X, &cv.vector, coeffs, cov, mw);
559 for (
int j = 0; j < n_wave; j++) {
560 wvi = cpl_array_get(wave, j, NULL);
561 gsl_bspline_eval(wvi, B, bw);
562 gsl_multifit_robust_est(B, coeffs, cov, &fit_val, &fit_err);
563 gsl_matrix_set(result, c, j, fit_val);
572 gsl_vector_free(coeffs);
573 gsl_matrix_free(cov);
574 gsl_bspline_free(bw);
575 gsl_multifit_robust_free(mw);
576 return CPL_ERROR_NONE;
589 cpl_ensure_code(self, CPL_ERROR_NULL_INPUT);
592 cpl_ensure_code(wave, CPL_ERROR_NULL_INPUT);
593 cpl_ensure_code(degree > 0, CPL_ERROR_ILLEGAL_INPUT);
594 cpl_ensure_code(ncomp > 0 && ncomp <= self->n_comp, CPL_ERROR_ILLEGAL_INPUT);
596 const int n_wave = self->
n_wave;
597 const int n_coeffs = degree;
599 const gsl_multifit_robust_type *T = gsl_multifit_robust_default;
600 gsl_multifit_robust_workspace *mw = gsl_multifit_robust_alloc(T, n_wave, n_coeffs);
601 gsl_multifit_robust_maxiter(1000, mw);
603 gsl_vector *coeffs = gsl_vector_alloc(n_coeffs);
604 gsl_matrix *cov = gsl_matrix_alloc(n_coeffs, n_coeffs);
605 double wvi, Xij, fit_val, fit_err;
607 gsl_matrix *result = gsl_matrix_alloc(ncomp, n_wave);
610 for (
int c = 0; c < ncomp; c++) {
614 gsl_matrix *X = gsl_matrix_alloc(n_wave, n_coeffs);
615 for (
int i = 0; i < n_wave; i++) {
616 wvi = cpl_array_get(wave, i, NULL);
618 for (
int j = 0; j < n_coeffs; j++) {
620 gsl_matrix_set(X, i, j, Xij);
625 gsl_multifit_robust(X, &cv.vector, coeffs, cov, mw);
628 for (
int j = 0; j < n_wave; j++) {
629 gsl_vector_const_view Xi = gsl_matrix_const_row(X, j);
630 gsl_multifit_robust_est(&Xi.vector, coeffs, cov, &fit_val, &fit_err);
631 gsl_matrix_set(result, c, j, fit_val);
639 gsl_vector_free(coeffs);
640 gsl_matrix_free(cov);
641 gsl_multifit_robust_free(mw);
642 return CPL_ERROR_NONE;
654 cpl_ensure_code(self, CPL_ERROR_NULL_INPUT);
655 cpl_ensure_code(residual, CPL_ERROR_NULL_INPUT);
657 int n_obs = self->
n_obs;
658 int n_wave = self->
n_wave;
661 gsl_matrix_const_view resid_gsl_view = gsl_matrix_const_view_array(
662 cpl_matrix_get_data_const(residual), n_obs, n_wave);
665 gsl_vector* refined_mean = gsl_vector_alloc(n_wave);
666 for(
int j = 0; j < n_wave; j++) {
667 gsl_vector_const_view col = gsl_matrix_const_column(&resid_gsl_view.matrix, j);
669 gsl_vector_set(refined_mean, j, mu);
673 for (
int i = 0; i < n_obs; i++) {
674 gsl_vector_view row = gsl_matrix_row(self->
data, i);
675 gsl_vector_add(&row.vector, self->
mean);
676 gsl_vector_sub(&row.vector, refined_mean);
680 gsl_vector_swap(self->
mean, refined_mean);
681 gsl_vector_free(refined_mean);
683 return CPL_ERROR_NONE;
696 cpl_ensure(self, CPL_ERROR_NULL_INPUT, NULL);
698 cpl_ensure(component > 0 && component <= self->n_comp, CPL_ERROR_ILLEGAL_INPUT, NULL);
700 cpl_vector *result = cpl_vector_new(self->
n_wave);
701 gsl_vector_const_view row = gsl_matrix_const_row(self->
median_components, component-1);
702 for(
int i = 0; i < self->
n_wave; i++)
703 cpl_vector_set(result, i, gsl_vector_get(&row.vector, i));
718 cpl_ensure(self, CPL_ERROR_NULL_INPUT, NULL);
720 cpl_ensure(component > 0 && component <= self->n_comp, CPL_ERROR_ILLEGAL_INPUT, NULL);
722 cpl_vector *result = cpl_vector_new(self->
n_wave);
723 gsl_vector_const_view row = gsl_matrix_const_row(self->
component_fits, component-1);
724 for(
int i = 0; i < self->
n_wave; i++)
725 cpl_vector_set(result, i, gsl_vector_get(&row.vector, i));
737 gsl_vector_int_const_view
mask;
748 gsl_vector *model = gsl_vector_calloc(p->
n_data);
749 gsl_vector *scratch = gsl_vector_alloc(p->
n_data);
750 gsl_vector *chi2 = gsl_vector_alloc(p->
n_data);
754 gsl_vector_const_view cv = gsl_matrix_const_row(p->
components, i);
755 gsl_vector_memcpy(scratch, &cv.vector);
756 gsl_vector_scale(scratch, gsl_vector_get(x, i));
757 gsl_vector_add(model, scratch);
761 gsl_vector_memcpy(chi2, &(p->
data.vector));
762 gsl_vector_sub(chi2, model);
763 gsl_vector_mul(chi2, chi2);
766 double chi2_sum = 0.0;
767 for (
int i = 0; i < p->
n_data; i++) {
768 if(gsl_vector_int_get(&(p->
mask.vector), i) > 0) {
769 chi2_sum += gsl_vector_get(chi2, i);
775 if (n_valid <= p->n_components)
780 gsl_vector_free(chi2);
781 gsl_vector_free(scratch);
782 gsl_vector_free(model);
794 cpl_ensure_code(self, CPL_ERROR_NULL_INPUT);
795 if (fit_mean_subtracted)
796 cpl_ensure_code(self->
mean, CPL_ERROR_ILLEGAL_INPUT);
798 cpl_ensure_code(model, CPL_ERROR_NULL_INPUT);
801 const int MAX_ITERATIONS = 1000;
802 const int n_obs = self->
n_obs;
803 const int n_wave = model->
n_wave;
804 const int n_comp = model->
n_comp;
806 const gsl_multimin_fminimizer_type *T = gsl_multimin_fminimizer_nmsimplex2rand;
807 gsl_multimin_fminimizer *mini = gsl_multimin_fminimizer_alloc(T, n_comp);
810 gsl_vector *coeffs = gsl_vector_alloc(n_comp);
811 gsl_vector *step_size = gsl_vector_alloc(n_comp);
814 gsl_matrix *data = gsl_matrix_calloc(n_obs, n_wave);
815 if (!fit_mean_subtracted) {
816 for (
int i = 0; i < n_obs; i++)
817 gsl_matrix_set_row(data, i, self->
mean);
819 gsl_matrix_add(data, self->
data);
821 self->
data_fit = gsl_matrix_calloc(n_obs, n_wave);
825 for (
int i = 0; i < n_obs; i++) {
826 gsl_vector_set_all(coeffs, 1.0);
827 gsl_vector_set_all(step_size, 0.1);
831 .
n_data = n_wave, .n_components = n_comp,
833 .data = gsl_matrix_const_row(data, i),
834 .mask = gsl_matrix_int_const_row(self->
mask, i)
838 gsl_multimin_function func = {
844 gsl_multimin_fminimizer_set(mini, &func, coeffs, step_size);
849 printf(
"iterations c_0 c_1 size chi2\n");
852 status = gsl_multimin_fminimizer_iterate(mini);
859 if (mini->fval < 0) {
860 cpl_error_set_message(cpl_func, CPL_ERROR_ILLEGAL_INPUT,
"No valid visphi data for fitting\n");
865 size = gsl_multimin_fminimizer_size(mini);
866 status = gsl_multimin_test_size(size, 1e-3);
869 printf(
"%d %f %f %f %f\n", iterations, mini->x->data[0], mini->x->data[1], size, mini->fval);
870 }
while (status == GSL_CONTINUE && iterations < MAX_ITERATIONS);
873 if (status == GSL_CONTINUE) {
874 cpl_msg_warning(cpl_func,
"model-fitting did not converge");
875 }
else if (status != GSL_SUCCESS) {
876 cpl_msg_error(cpl_func,
"model-fitting failed with error code %d (%s)", status, gsl_strerror(status));
879 gsl_vector_view fit_vals = gsl_matrix_row(self->
data_fit, i);
880 for (
int c = 0; c < n_comp; c++) {
881 gsl_vector *term = gsl_vector_alloc_row_from_matrix(model->
component_fits, c);
882 gsl_vector_scale(term, gsl_vector_get(mini->x, c));
883 gsl_vector_add(&fit_vals.vector, term);
884 gsl_vector_free(term);
887 if (fit_mean_subtracted) {
888 gsl_vector_add(&fit_vals.vector, self->
mean);
892 gsl_matrix_free(data);
893 gsl_vector_free(coeffs);
894 gsl_vector_free(step_size);
895 gsl_multimin_fminimizer_free(mini);
896 return status ? CPL_ERROR_ILLEGAL_OUTPUT : CPL_ERROR_NONE;
909 cpl_ensure(self, CPL_ERROR_NULL_INPUT, NULL);
910 cpl_ensure(component >= 0 && component <= self->n_valid, CPL_ERROR_ILLEGAL_INPUT, NULL);
912 cpl_vector *result = cpl_vector_new(self->
n_wave);
913 if (component == 0) {
915 for(
int i = 0; i < self->
n_wave; i++)
916 cpl_vector_set(result, i, gsl_vector_get(self->
mean, i));
919 gsl_vector_const_view row = gsl_matrix_const_row(self->
Vt, component-1);
920 for(
int i = 0; i < self->
n_wave; i++) {
921 cpl_vector_set(result, i, gsl_vector_get(&row.vector, i));
923 cpl_vector_multiply_scalar(result, gsl_vector_get(self->
signs, component-1));
937 cpl_ensure(self, CPL_ERROR_NULL_INPUT, NULL);
938 cpl_ensure(self->
data_fit, CPL_ERROR_NULL_INPUT, NULL);
940 cpl_matrix *result = cpl_matrix_new(self->
n_obs, self->
n_wave);
941 for(
int i = 0; i < self->
n_obs; i++) {
942 for(
int j = 0; j < self->
n_wave; j++) {
943 cpl_matrix_set(result, i, j, gsl_matrix_get(self->
data_fit, i, j));
958 cpl_ensure(self, CPL_ERROR_NULL_INPUT, NULL);
959 cpl_ensure(self->
data, CPL_ERROR_NULL_INPUT, NULL);
960 cpl_ensure(self->
mean, CPL_ERROR_NULL_INPUT, NULL);
961 cpl_ensure(self->
data_fit, CPL_ERROR_NULL_INPUT, NULL);
964 gsl_matrix *resid = gsl_matrix_alloc(self->
n_obs, self->
n_wave);
965 gsl_matrix_memcpy(resid, self->
data);
967 for(
int i = 0; i < self->
n_obs; i++) {
968 gsl_vector_view row_view = gsl_matrix_row(resid, i);
969 gsl_vector_add(&row_view.vector, self->
mean);
972 gsl_matrix_sub(resid, self->
data_fit);
974 cpl_matrix *result = cpl_matrix_new(self->
n_obs, self->
n_wave);
975 for(
int i = 0; i < self->
n_obs; i++) {
976 for(
int j = 0; j < self->
n_wave; j++) {
977 cpl_matrix_set(result, i, j, gsl_matrix_get(resid, i, j));
cpl_error_code gravi_pca_fit_components_polynomial(gravi_pca_model *self, const cpl_array *wave, int degree, int ncomp)
Fit polynomial model to each of a set of median-averaged PCA components.
cpl_error_code gravi_pca_refine_mean(gravi_pca_result *self, const cpl_matrix *residual)
Override the mean component calculated from the PCA decomposition. In conjunction with the fit_mean_s...
cpl_matrix * gravi_pca_get_data_residual(const gravi_pca_result *self)
Get residual (data - model).
static int svd_wrapper(gsl_matrix **A, gsl_matrix **V, gsl_vector **S)
Handle computing SVD for MxN matrix, M<N.
cpl_error_code gravi_pca_set_component_signs(gravi_pca_result *self, const cpl_array *wave, const int ncomp)
Impose sign convention on PCA components.
gravi_pca_model * gravi_pca_create_model(const gravi_pca_result **results, int nres, int ncomp)
Compute median values of PCA components over a set of decomposition results.
static double gravi_pca_model_chi2(const gsl_vector *x, void *params)
Evaluate chi^2 statistic for fit.
void gravi_pca_result_delete(gravi_pca_result *self)
Deallocate a gravi_pca_result object.
gravi_pca_result * gravi_pca_create_result(const cpl_matrix *data, const cpl_matrix *mask)
Construct a new gravi_pca_result object from a matrix of visphi data.
cpl_vector * gravi_pca_get_component(const gravi_pca_result *self, int component)
Get components from PCA decomposition.
cpl_matrix * gravi_pca_get_data_fit(const gravi_pca_result *self)
Get noise-free model.
cpl_vector * gravi_pca_get_component_median(const gravi_pca_model *self, int component)
Get median-averaged component from a set of PCA decompositions.
struct _gravi_pca_model_params_ gravi_pca_model_params
Type to hold parameters for fitting PCA components to data.
cpl_error_code gravi_pca_decomp_matrix_svd(gravi_pca_result *self)
Perform PCA decomposition by calculating singular value decomposition.
static double _gsl_vector_sum(const gsl_vector *a)
void gravi_pca_model_delete(gravi_pca_model *self)
Deallocate a gravi_pca_model object.
cpl_error_code gravi_pca_fit_components_bspline(gravi_pca_model *self, const cpl_array *wave, int degree, int ncomp)
Fit B-spline model to each of a set of median-averaged PCA components.
static double gravi_pca_get_component_sign(const gsl_vector *visphi, const gsl_vector *wave)
Calculate sign convention on visphi components.
cpl_vector * gravi_pca_get_component_fit(const gravi_pca_model *self, int component)
Get fit to components from PCA decomposition.
static void svd_flip(gsl_matrix *U, gsl_matrix *V)
Adjust U and V so that loadings for most significant components are positive.
gravi_pca_model * gravi_pca_load_model(const cpl_matrix *components)
Create PCA model from existing set of components.
cpl_error_code gravi_pca_fit_model(gravi_pca_result *self, const gravi_pca_model *model, cpl_boolean fit_mean_subtracted, cpl_boolean verbose)
Fit model formed from linear combination of PCA components to data.
#define FREE(function, variable)
Type to hold average (median) components obtained from a set of PCA decompositions and/or best-fit mo...
gsl_matrix * median_components
gsl_matrix * component_fits
Type to hold parameters for fitting PCA components to data.
gsl_vector_const_view data
gsl_vector_int_const_view mask
const gsl_matrix * components
Type to hold results of a PCA decomposition.