1 /* PSPP - a program for statistical analysis.
2 Copyright (C) 2009, 2010 Free Software Foundation, Inc.
4 This program is free software: you can redistribute it and/or modify
5 it under the terms of the GNU General Public License as published by
6 the Free Software Foundation, either version 3 of the License, or
7 (at your option) any later version.
9 This program is distributed in the hope that it will be useful,
10 but WITHOUT ANY WARRANTY; without even the implied warranty of
11 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12 GNU General Public License for more details.
14 You should have received a copy of the GNU General Public License
15 along with this program. If not, see <http://www.gnu.org/licenses/>. */
20 #include <gsl/gsl_vector.h>
21 #include <gsl/gsl_linalg.h>
22 #include <gsl/gsl_matrix.h>
23 #include <gsl/gsl_eigen.h>
24 #include <gsl/gsl_blas.h>
25 #include <gsl/gsl_sort_vector.h>
27 #include <math/covariance.h>
29 #include <math/correlation.h>
30 #include <math/moments.h>
31 #include <data/procedure.h>
32 #include <language/lexer/variable-parser.h>
33 #include <language/lexer/value-parser.h>
34 #include <language/command.h>
35 #include <language/lexer/lexer.h>
37 #include <data/casegrouper.h>
38 #include <data/casereader.h>
39 #include <data/casewriter.h>
40 #include <data/dictionary.h>
41 #include <data/format.h>
42 #include <data/subcase.h>
44 #include <libpspp/misc.h>
45 #include <libpspp/message.h>
47 #include <output/tab.h>
49 #include <output/charts/scree.h>
50 #include <output/chart-item.h>
53 #define _(msgid) gettext (msgid)
54 #define N_(msgid) msgid
69 enum extraction_method
78 PLOT_ROTATION = 0x0002
83 PRINT_UNIVARIATE = 0x0001,
84 PRINT_DETERMINANT = 0x0002,
88 PRINT_COVARIANCE = 0x0020,
89 PRINT_CORRELATION = 0x0040,
90 PRINT_ROTATION = 0x0080,
91 PRINT_EXTRACTION = 0x0100,
92 PRINT_INITIAL = 0x0200,
106 typedef void (*rotation_coefficients) (double *x, double *y,
107 double a, double b, double c, double d,
108 const gsl_matrix *loadings );
112 varimax_coefficients (double *x, double *y,
113 double a, double b, double c, double d,
114 const gsl_matrix *loadings )
116 *x = d - 2 * a * b / loadings->size1;
117 *y = c - (a * a - b * b) / loadings->size1;
121 equamax_coefficients (double *x, double *y,
122 double a, double b, double c, double d,
123 const gsl_matrix *loadings )
125 *x = d - loadings->size2 * a * b / loadings->size1;
126 *y = c - loadings->size2 * (a * a - b * b) / (2 * loadings->size1);
130 quartimax_coefficients (double *x, double *y,
131 double a UNUSED, double b UNUSED, double c, double d,
132 const gsl_matrix *loadings UNUSED)
138 static const rotation_coefficients rotation_coeff[3] = {
139 varimax_coefficients,
140 equamax_coefficients,
141 quartimax_coefficients
148 const struct variable **vars;
150 const struct variable *wv;
153 enum missing_type missing_type;
154 enum mv_class exclude;
155 enum print_opts print;
156 enum extraction_method extraction;
158 enum rotation_type rotation;
160 /* Extraction Criteria */
173 /* Intermediate values used in calculation */
175 const gsl_matrix *corr ; /* The correlation matrix */
176 const gsl_matrix *cov ; /* The covariance matrix */
177 const gsl_matrix *n ; /* Matrix of number of samples */
179 gsl_vector *eval ; /* The eigenvalues */
180 gsl_matrix *evec ; /* The eigenvectors */
184 gsl_vector *msr ; /* Multiple Squared Regressions */
187 static struct idata *
188 idata_alloc (size_t n_vars)
190 struct idata *id = xzalloc (sizeof (*id));
192 id->n_extractions = 0;
193 id->msr = gsl_vector_alloc (n_vars);
195 id->eval = gsl_vector_alloc (n_vars);
196 id->evec = gsl_matrix_alloc (n_vars, n_vars);
202 idata_free (struct idata *id)
204 gsl_vector_free (id->msr);
205 gsl_vector_free (id->eval);
206 gsl_matrix_free (id->evec);
213 dump_matrix (const gsl_matrix *m)
217 for (i = 0 ; i < m->size1; ++i)
219 for (j = 0 ; j < m->size2; ++j)
220 printf ("%02f ", gsl_matrix_get (m, i, j));
227 dump_matrix_permute (const gsl_matrix *m, const gsl_permutation *p)
231 for (i = 0 ; i < m->size1; ++i)
233 for (j = 0 ; j < m->size2; ++j)
234 printf ("%02f ", gsl_matrix_get (m, gsl_permutation_get (p, i), j));
241 dump_vector (const gsl_vector *v)
244 for (i = 0 ; i < v->size; ++i)
246 printf ("%02f\n", gsl_vector_get (v, i));
253 n_extracted_factors (const struct cmd_factor *factor, struct idata *idata)
257 /* If there is a cached value, then return that. */
258 if ( idata->n_extractions != 0)
259 return idata->n_extractions;
261 /* Otherwise, if the number of factors has been explicitly requested,
263 if (factor->n_factors > 0)
265 idata->n_extractions = factor->n_factors;
269 /* Use the MIN_EIGEN setting. */
270 for (i = 0 ; i < idata->eval->size; ++i)
272 double evali = fabs (gsl_vector_get (idata->eval, i));
274 idata->n_extractions = i;
276 if (evali < factor->min_eigen)
281 return idata->n_extractions;
285 /* Returns a newly allocated matrix identical to M.
286 It it the callers responsibility to free the returned value.
289 matrix_dup (const gsl_matrix *m)
291 gsl_matrix *n = gsl_matrix_alloc (m->size1, m->size2);
293 gsl_matrix_memcpy (n, m);
301 /* Copy of the subject */
306 gsl_permutation *perm;
313 static struct smr_workspace *ws_create (const gsl_matrix *input)
315 struct smr_workspace *ws = xmalloc (sizeof (*ws));
317 ws->m = gsl_matrix_alloc (input->size1, input->size2);
318 ws->inverse = gsl_matrix_calloc (input->size1 - 1, input->size2 - 1);
319 ws->perm = gsl_permutation_alloc (input->size1 - 1);
320 ws->result1 = gsl_matrix_calloc (input->size1 - 1, 1);
321 ws->result2 = gsl_matrix_calloc (1, 1);
327 ws_destroy (struct smr_workspace *ws)
329 gsl_matrix_free (ws->result2);
330 gsl_matrix_free (ws->result1);
331 gsl_permutation_free (ws->perm);
332 gsl_matrix_free (ws->inverse);
333 gsl_matrix_free (ws->m);
340 Return the square of the regression coefficient for VAR regressed against all other variables.
343 squared_multiple_correlation (const gsl_matrix *corr, int var, struct smr_workspace *ws)
345 /* For an explanation of what this is doing, see
346 http://www.visualstatistics.net/Visual%20Statistics%20Multimedia/multiple_regression_analysis.htm
352 gsl_matrix_memcpy (ws->m, corr);
354 gsl_matrix_swap_rows (ws->m, 0, var);
355 gsl_matrix_swap_columns (ws->m, 0, var);
357 rxx = gsl_matrix_submatrix (ws->m, 1, 1, ws->m->size1 - 1, ws->m->size1 - 1);
359 gsl_linalg_LU_decomp (&rxx.matrix, ws->perm, &signum);
361 gsl_linalg_LU_invert (&rxx.matrix, ws->perm, ws->inverse);
364 gsl_matrix_const_view rxy = gsl_matrix_const_submatrix (ws->m, 1, 0, ws->m->size1 - 1, 1);
365 gsl_matrix_const_view ryx = gsl_matrix_const_submatrix (ws->m, 0, 1, 1, ws->m->size1 - 1);
367 gsl_blas_dgemm (CblasNoTrans, CblasNoTrans,
368 1.0, ws->inverse, &rxy.matrix, 0.0, ws->result1);
370 gsl_blas_dgemm (CblasNoTrans, CblasNoTrans,
371 1.0, &ryx.matrix, ws->result1, 0.0, ws->result2);
374 return gsl_matrix_get (ws->result2, 0, 0);
379 static double the_communality (const gsl_matrix *evec, const gsl_vector *eval, int n, int n_factors);
382 struct factor_matrix_workspace
385 gsl_eigen_symmv_workspace *eigen_ws;
395 static struct factor_matrix_workspace *
396 factor_matrix_workspace_alloc (size_t n, size_t nf)
398 struct factor_matrix_workspace *ws = xmalloc (sizeof (*ws));
401 ws->gamma = gsl_matrix_calloc (nf, nf);
402 ws->eigen_ws = gsl_eigen_symmv_alloc (n);
403 ws->eval = gsl_vector_alloc (n);
404 ws->evec = gsl_matrix_alloc (n, n);
405 ws->r = gsl_matrix_alloc (n, n);
411 factor_matrix_workspace_free (struct factor_matrix_workspace *ws)
413 gsl_eigen_symmv_free (ws->eigen_ws);
414 gsl_vector_free (ws->eval);
415 gsl_matrix_free (ws->evec);
416 gsl_matrix_free (ws->gamma);
417 gsl_matrix_free (ws->r);
422 Shift P left by OFFSET places, and overwrite TARGET
423 with the shifted result.
424 Positions in TARGET less than OFFSET are unchanged.
427 perm_shift_apply (gsl_permutation *target, const gsl_permutation *p,
431 assert (target->size == p->size);
432 assert (offset <= target->size);
434 for (i = 0; i < target->size - offset; ++i)
436 target->data[i] = p->data [i + offset];
442 Indirectly sort the rows of matrix INPUT, storing the sort order in PERM.
443 The sort criteria are as follows:
445 Rows are sorted on the first column, until the absolute value of an
446 element in a subsequent column is greater than that of the first
447 column. Thereafter, rows will be sorted on the second column,
448 until the absolute value of an element in a subsequent column
449 exceeds that of the second column ...
452 sort_matrix_indirect (const gsl_matrix *input, gsl_permutation *perm)
454 const size_t n = perm->size;
455 const size_t m = input->size2;
462 assert (perm->size == input->size1);
464 p = gsl_permutation_alloc (n);
466 /* Copy INPUT into MAT, discarding the sign */
467 mat = gsl_matrix_alloc (n, m);
468 for (i = 0 ; i < mat->size1; ++i)
470 for (j = 0 ; j < mat->size2; ++j)
472 double x = gsl_matrix_get (input, i, j);
473 gsl_matrix_set (mat, i, j, fabs (x));
477 while (column_n < m && row_n < n)
479 gsl_vector_const_view columni = gsl_matrix_const_column (mat, column_n);
480 gsl_sort_vector_index (p, &columni.vector);
482 for (i = 0 ; i < n; ++i)
484 gsl_vector_view row = gsl_matrix_row (mat, p->data[n - 1 - i]);
485 size_t maxindex = gsl_vector_max_index (&row.vector);
487 if ( maxindex > column_n )
490 /* All subsequent elements of this row, are of no interest.
491 So set them all to a highly negative value */
492 for (j = column_n + 1; j < row.vector.size ; ++j)
493 gsl_vector_set (&row.vector, j, -DBL_MAX);
496 perm_shift_apply (perm, p, row_n);
502 gsl_permutation_free (p);
503 gsl_matrix_free (mat);
505 assert ( 0 == gsl_permutation_valid (perm));
507 /* We want the biggest value to be first */
508 gsl_permutation_reverse (perm);
513 drot_go (double phi, double *l0, double *l1)
515 double r0 = cos (phi) * *l0 + sin (phi) * *l1;
516 double r1 = - sin (phi) * *l0 + cos (phi) * *l1;
524 clone_matrix (const gsl_matrix *m)
527 gsl_matrix *c = gsl_matrix_calloc (m->size1, m->size2);
529 for (j = 0 ; j < c->size1; ++j)
531 for (k = 0 ; k < c->size2; ++k)
533 const double *v = gsl_matrix_const_ptr (m, j, k);
534 gsl_matrix_set (c, j, k, *v);
543 rotate (const gsl_matrix *unrot, const gsl_vector *communalities, enum rotation_type rot_type, gsl_matrix *result)
548 /* First get a normalised version of UNROT */
549 gsl_matrix *normalised = gsl_matrix_calloc (unrot->size1, unrot->size2);
550 gsl_matrix *h_sqrt = gsl_matrix_calloc (communalities->size, communalities->size);
551 gsl_matrix *h_sqrt_inv ;
553 /* H is the diagonal matrix containing the absolute values of the communalities */
554 for (i = 0 ; i < communalities->size ; ++i)
556 double *ptr = gsl_matrix_ptr (h_sqrt, i, i);
557 *ptr = fabs (gsl_vector_get (communalities, i));
560 /* Take the square root of the communalities */
561 gsl_linalg_cholesky_decomp (h_sqrt);
564 /* Save a copy of h_sqrt and invert it */
565 h_sqrt_inv = clone_matrix (h_sqrt);
566 gsl_linalg_cholesky_decomp (h_sqrt_inv);
567 gsl_linalg_cholesky_invert (h_sqrt_inv);
569 /* normalised vertion is H^{1/2} x UNROT */
570 gsl_blas_dgemm (CblasNoTrans, CblasNoTrans, 1.0, h_sqrt_inv, unrot, 0.0, normalised);
572 gsl_matrix_free (h_sqrt_inv);
574 /* Now perform the rotation iterations */
575 for (i = 0 ; i < 25 ; ++i)
577 for (j = 0 ; j < normalised->size2; ++j)
579 for (k = j + 1 ; k < normalised->size2; ++k)
588 for (p = 0; p < normalised->size1; ++p)
590 double jv = gsl_matrix_get (normalised, p, j);
591 double kv = gsl_matrix_get (normalised, p, k);
593 double u = jv * jv - kv * kv;
594 double v = 2 * jv * kv;
601 rotation_coeff [rot_type] (&x, &y, a, b, c, d, normalised);
603 phi = atan2 (x, y) / 4.0 ;
605 for (p = 0; p < normalised->size1; ++p)
607 double *lambda0 = gsl_matrix_ptr (normalised, p, j);
608 double *lambda1 = gsl_matrix_ptr (normalised, p, k);
609 drot_go (phi, lambda0, lambda1);
615 gsl_blas_dgemm (CblasNoTrans, CblasNoTrans, 1.0,
616 h_sqrt, normalised, 0.0, result);
618 gsl_matrix_free (h_sqrt);
621 /* reflect negative sums */
622 for (i = 0 ; i < result->size2; ++i)
625 for (j = 0 ; j < result->size1; ++j)
627 sum += gsl_matrix_get (result, j, i);
631 for (j = 0 ; j < result->size1; ++j)
633 double *lambda = gsl_matrix_ptr (result, j, i);
638 // dump_matrix (result);
643 Get an approximation for the factor matrix into FACTORS, and the communalities into COMMUNALITIES.
644 R is the matrix to be analysed.
645 WS is a pointer to a structure which must have been initialised with factor_matrix_workspace_init.
648 iterate_factor_matrix (const gsl_matrix *r, gsl_vector *communalities, gsl_matrix *factors,
649 struct factor_matrix_workspace *ws)
654 assert (r->size1 == r->size2);
655 assert (r->size1 == communalities->size);
657 assert (factors->size1 == r->size1);
658 assert (factors->size2 == ws->n_factors);
660 gsl_matrix_memcpy (ws->r, r);
662 /* Apply Communalities to diagonal of correlation matrix */
663 for (i = 0 ; i < communalities->size ; ++i)
665 double *x = gsl_matrix_ptr (ws->r, i, i);
666 *x = gsl_vector_get (communalities, i);
669 gsl_eigen_symmv (ws->r, ws->eval, ws->evec, ws->eigen_ws);
671 mv = gsl_matrix_submatrix (ws->evec, 0, 0, ws->evec->size1, ws->n_factors);
673 /* Gamma is the diagonal matrix containing the absolute values of the eigenvalues */
674 for (i = 0 ; i < ws->n_factors ; ++i)
676 double *ptr = gsl_matrix_ptr (ws->gamma, i, i);
677 *ptr = fabs (gsl_vector_get (ws->eval, i));
680 /* Take the square root of gamma */
681 gsl_linalg_cholesky_decomp (ws->gamma);
683 gsl_blas_dgemm (CblasNoTrans, CblasNoTrans, 1.0, &mv.matrix, ws->gamma, 0.0, factors);
685 for (i = 0 ; i < r->size1 ; ++i)
687 double h = the_communality (ws->evec, ws->eval, i, ws->n_factors);
688 gsl_vector_set (communalities, i, h);
694 static bool run_factor (struct dataset *ds, const struct cmd_factor *factor);
698 cmd_factor (struct lexer *lexer, struct dataset *ds)
700 bool extraction_seen = false;
701 const struct dictionary *dict = dataset_dict (ds);
703 struct cmd_factor factor;
706 factor.method = METHOD_CORR;
707 factor.missing_type = MISS_LISTWISE;
708 factor.exclude = MV_ANY;
709 factor.print = PRINT_INITIAL | PRINT_EXTRACTION | PRINT_ROTATION;
710 factor.extraction = EXTRACTION_PC;
711 factor.n_factors = 0;
712 factor.min_eigen = SYSMIS;
713 factor.iterations = 25;
714 factor.econverge = 0.001;
718 factor.rotation = ROT_VARIMAX;
720 factor.wv = dict_get_weight (dict);
722 lex_match (lexer, '/');
724 if (!lex_force_match_id (lexer, "VARIABLES"))
729 lex_match (lexer, '=');
731 if (!parse_variables_const (lexer, dict, &factor.vars, &factor.n_vars,
732 PV_NO_DUPLICATE | PV_NUMERIC))
735 if (factor.n_vars < 2)
736 msg (MW, _("Factor analysis on a single variable is not useful."));
738 while (lex_token (lexer) != '.')
740 lex_match (lexer, '/');
742 if (lex_match_id (lexer, "PLOT"))
744 lex_match (lexer, '=');
745 while (lex_token (lexer) != '.' && lex_token (lexer) != '/')
747 if (lex_match_id (lexer, "EIGEN"))
749 factor.plot |= PLOT_SCREE;
751 #if FACTOR_FULLY_IMPLEMENTED
752 else if (lex_match_id (lexer, "ROTATION"))
758 lex_error (lexer, NULL);
763 else if (lex_match_id (lexer, "METHOD"))
765 lex_match (lexer, '=');
766 while (lex_token (lexer) != '.' && lex_token (lexer) != '/')
768 if (lex_match_id (lexer, "COVARIANCE"))
770 factor.method = METHOD_COV;
772 else if (lex_match_id (lexer, "CORRELATION"))
774 factor.method = METHOD_CORR;
778 lex_error (lexer, NULL);
783 else if (lex_match_id (lexer, "ROTATION"))
785 lex_match (lexer, '=');
786 while (lex_token (lexer) != '.' && lex_token (lexer) != '/')
788 /* VARIMAX and DEFAULT are defaults */
789 if (lex_match_id (lexer, "VARIMAX") || lex_match_id (lexer, "DEFAULT"))
791 factor.rotation = ROT_VARIMAX;
793 else if (lex_match_id (lexer, "EQUAMAX"))
795 factor.rotation = ROT_EQUAMAX;
797 else if (lex_match_id (lexer, "QUARTIMAX"))
799 factor.rotation = ROT_QUARTIMAX;
801 else if (lex_match_id (lexer, "NOROTATE"))
803 factor.rotation = ROT_NONE;
807 lex_error (lexer, NULL);
812 else if (lex_match_id (lexer, "CRITERIA"))
814 lex_match (lexer, '=');
815 while (lex_token (lexer) != '.' && lex_token (lexer) != '/')
817 if (lex_match_id (lexer, "FACTORS"))
819 if ( lex_force_match (lexer, '('))
821 lex_force_int (lexer);
822 factor.n_factors = lex_integer (lexer);
824 lex_force_match (lexer, ')');
827 else if (lex_match_id (lexer, "MINEIGEN"))
829 if ( lex_force_match (lexer, '('))
831 lex_force_num (lexer);
832 factor.min_eigen = lex_number (lexer);
834 lex_force_match (lexer, ')');
837 else if (lex_match_id (lexer, "ECONVERGE"))
839 if ( lex_force_match (lexer, '('))
841 lex_force_num (lexer);
842 factor.econverge = lex_number (lexer);
844 lex_force_match (lexer, ')');
847 else if (lex_match_id (lexer, "ITERATE"))
849 if ( lex_force_match (lexer, '('))
851 lex_force_int (lexer);
852 factor.iterations = lex_integer (lexer);
854 lex_force_match (lexer, ')');
857 else if (lex_match_id (lexer, "DEFAULT"))
859 factor.n_factors = 0;
860 factor.min_eigen = 1;
861 factor.iterations = 25;
865 lex_error (lexer, NULL);
870 else if (lex_match_id (lexer, "EXTRACTION"))
872 extraction_seen = true;
873 lex_match (lexer, '=');
874 while (lex_token (lexer) != '.' && lex_token (lexer) != '/')
876 if (lex_match_id (lexer, "PAF"))
878 factor.extraction = EXTRACTION_PAF;
880 else if (lex_match_id (lexer, "PC"))
882 factor.extraction = EXTRACTION_PC;
884 else if (lex_match_id (lexer, "PA1"))
886 factor.extraction = EXTRACTION_PC;
888 else if (lex_match_id (lexer, "DEFAULT"))
890 factor.extraction = EXTRACTION_PC;
894 lex_error (lexer, NULL);
899 else if (lex_match_id (lexer, "FORMAT"))
901 lex_match (lexer, '=');
902 while (lex_token (lexer) != '.' && lex_token (lexer) != '/')
904 if (lex_match_id (lexer, "SORT"))
908 else if (lex_match_id (lexer, "BLANK"))
910 if ( lex_force_match (lexer, '('))
912 lex_force_num (lexer);
913 factor.blank = lex_number (lexer);
915 lex_force_match (lexer, ')');
918 else if (lex_match_id (lexer, "DEFAULT"))
925 lex_error (lexer, NULL);
930 else if (lex_match_id (lexer, "PRINT"))
933 lex_match (lexer, '=');
934 while (lex_token (lexer) != '.' && lex_token (lexer) != '/')
936 if (lex_match_id (lexer, "UNIVARIATE"))
938 factor.print |= PRINT_UNIVARIATE;
940 else if (lex_match_id (lexer, "DET"))
942 factor.print |= PRINT_DETERMINANT;
944 #if FACTOR_FULLY_IMPLEMENTED
945 else if (lex_match_id (lexer, "INV"))
948 else if (lex_match_id (lexer, "AIC"))
952 else if (lex_match_id (lexer, "SIG"))
954 factor.print |= PRINT_SIG;
956 else if (lex_match_id (lexer, "CORRELATION"))
958 factor.print |= PRINT_CORRELATION;
960 #if FACTOR_FULLY_IMPLEMENTED
961 else if (lex_match_id (lexer, "COVARIANCE"))
965 else if (lex_match_id (lexer, "ROTATION"))
967 factor.print |= PRINT_ROTATION;
969 else if (lex_match_id (lexer, "EXTRACTION"))
971 factor.print |= PRINT_EXTRACTION;
973 else if (lex_match_id (lexer, "INITIAL"))
975 factor.print |= PRINT_INITIAL;
977 #if FACTOR_FULLY_IMPLEMENTED
978 else if (lex_match_id (lexer, "KMO"))
981 else if (lex_match_id (lexer, "REPR"))
984 else if (lex_match_id (lexer, "FSCORE"))
988 else if (lex_match (lexer, T_ALL))
990 factor.print = 0xFFFF;
992 else if (lex_match_id (lexer, "DEFAULT"))
994 factor.print |= PRINT_INITIAL ;
995 factor.print |= PRINT_EXTRACTION ;
996 factor.print |= PRINT_ROTATION ;
1000 lex_error (lexer, NULL);
1005 else if (lex_match_id (lexer, "MISSING"))
1007 lex_match (lexer, '=');
1008 while (lex_token (lexer) != '.' && lex_token (lexer) != '/')
1010 if (lex_match_id (lexer, "INCLUDE"))
1012 factor.exclude = MV_SYSTEM;
1014 else if (lex_match_id (lexer, "EXCLUDE"))
1016 factor.exclude = MV_ANY;
1018 else if (lex_match_id (lexer, "LISTWISE"))
1020 factor.missing_type = MISS_LISTWISE;
1022 else if (lex_match_id (lexer, "PAIRWISE"))
1024 factor.missing_type = MISS_PAIRWISE;
1026 else if (lex_match_id (lexer, "MEANSUB"))
1028 factor.missing_type = MISS_MEANSUB;
1032 lex_error (lexer, NULL);
1039 lex_error (lexer, NULL);
1044 if ( ! run_factor (ds, &factor))
1055 static void do_factor (const struct cmd_factor *factor, struct casereader *group);
1059 run_factor (struct dataset *ds, const struct cmd_factor *factor)
1061 struct dictionary *dict = dataset_dict (ds);
1063 struct casereader *group;
1065 struct casegrouper *grouper = casegrouper_create_splits (proc_open (ds), dict);
1067 while (casegrouper_get_next_group (grouper, &group))
1069 if ( factor->missing_type == MISS_LISTWISE )
1070 group = casereader_create_filter_missing (group, factor->vars, factor->n_vars,
1073 do_factor (factor, group);
1076 ok = casegrouper_destroy (grouper);
1077 ok = proc_commit (ds) && ok;
1083 /* Return the communality of variable N, calculated to N_FACTORS */
1085 the_communality (const gsl_matrix *evec, const gsl_vector *eval, int n, int n_factors)
1092 assert (n < eval->size);
1093 assert (n < evec->size1);
1094 assert (n_factors <= eval->size);
1096 for (i = 0 ; i < n_factors; ++i)
1098 double evali = fabs (gsl_vector_get (eval, i));
1100 double eveci = gsl_matrix_get (evec, n, i);
1102 comm += pow2 (eveci) * evali;
1108 /* Return the communality of variable N, calculated to N_FACTORS */
1110 communality (struct idata *idata, int n, int n_factors)
1112 return the_communality (idata->evec, idata->eval, n, n_factors);
1117 show_scree (const struct cmd_factor *f, struct idata *idata)
1122 if ( !(f->plot & PLOT_SCREE) )
1126 label = f->extraction == EXTRACTION_PC ? _("Component Number") : _("Factor Number");
1128 s = scree_create (idata->eval, label);
1134 show_communalities (const struct cmd_factor * factor,
1135 const gsl_vector *initial, const gsl_vector *extracted)
1139 const int heading_columns = 1;
1140 int nc = heading_columns;
1141 const int heading_rows = 1;
1142 const int nr = heading_rows + factor->n_vars;
1143 struct tab_table *t;
1145 if (factor->print & PRINT_EXTRACTION)
1148 if (factor->print & PRINT_INITIAL)
1151 /* No point having a table with only headings */
1155 t = tab_create (nc, nr);
1157 tab_title (t, _("Communalities"));
1159 tab_headers (t, heading_columns, 0, heading_rows, 0);
1162 if (factor->print & PRINT_INITIAL)
1163 tab_text (t, c++, 0, TAB_CENTER | TAT_TITLE, _("Initial"));
1165 if (factor->print & PRINT_EXTRACTION)
1166 tab_text (t, c++, 0, TAB_CENTER | TAT_TITLE, _("Extraction"));
1168 /* Outline the box */
1175 /* Vertical lines */
1182 tab_hline (t, TAL_1, 0, nc - 1, heading_rows);
1183 tab_vline (t, TAL_2, heading_columns, 0, nr - 1);
1185 for (i = 0 ; i < factor->n_vars; ++i)
1188 tab_text (t, c++, i + heading_rows, TAT_TITLE, var_to_string (factor->vars[i]));
1190 if (factor->print & PRINT_INITIAL)
1191 tab_double (t, c++, i + heading_rows, 0, gsl_vector_get (initial, i), NULL);
1193 if (factor->print & PRINT_EXTRACTION)
1194 tab_double (t, c++, i + heading_rows, 0, gsl_vector_get (extracted, i), NULL);
1202 show_factor_matrix (const struct cmd_factor *factor, struct idata *idata, const char *title, const gsl_matrix *fm)
1205 const int n_factors = idata->n_extractions;
1207 const int heading_columns = 1;
1208 const int heading_rows = 2;
1209 const int nr = heading_rows + factor->n_vars;
1210 const int nc = heading_columns + n_factors;
1211 gsl_permutation *perm;
1213 struct tab_table *t = tab_create (nc, nr);
1216 if ( factor->extraction == EXTRACTION_PC )
1217 tab_title (t, _("Component Matrix"));
1219 tab_title (t, _("Factor Matrix"));
1222 tab_title (t, title);
1224 tab_headers (t, heading_columns, 0, heading_rows, 0);
1226 if ( factor->extraction == EXTRACTION_PC )
1230 TAB_CENTER | TAT_TITLE, _("Component"));
1235 TAB_CENTER | TAT_TITLE, _("Factor"));
1238 tab_hline (t, TAL_1, heading_columns, nc - 1, 1);
1241 /* Outline the box */
1248 /* Vertical lines */
1255 tab_hline (t, TAL_1, 0, nc - 1, heading_rows);
1256 tab_vline (t, TAL_2, heading_columns, 0, nr - 1);
1259 /* Initialise to the identity permutation */
1260 perm = gsl_permutation_calloc (factor->n_vars);
1263 sort_matrix_indirect (fm, perm);
1265 for (i = 0 ; i < n_factors; ++i)
1267 tab_text_format (t, heading_columns + i, 1, TAB_CENTER | TAT_TITLE, _("%d"), i + 1);
1270 for (i = 0 ; i < factor->n_vars; ++i)
1273 const int matrix_row = perm->data[i];
1274 tab_text (t, 0, i + heading_rows, TAT_TITLE, var_to_string (factor->vars[matrix_row]));
1276 for (j = 0 ; j < n_factors; ++j)
1278 double x = gsl_matrix_get (fm, matrix_row, j);
1280 if ( fabs (x) < factor->blank)
1283 tab_double (t, heading_columns + j, heading_rows + i, 0, x, NULL);
1287 gsl_permutation_free (perm);
1294 show_explained_variance (const struct cmd_factor * factor, struct idata *idata,
1295 const gsl_vector *initial_eigenvalues,
1296 const gsl_vector *extracted_eigenvalues)
1300 const int heading_columns = 1;
1301 const int heading_rows = 2;
1302 const int nr = heading_rows + factor->n_vars;
1304 struct tab_table *t ;
1306 double i_total = 0.0;
1309 double e_total = 0.0;
1312 int nc = heading_columns;
1314 if (factor->print & PRINT_EXTRACTION)
1317 if (factor->print & PRINT_INITIAL)
1320 if (factor->print & PRINT_ROTATION)
1323 /* No point having a table with only headings */
1324 if ( nc <= heading_columns)
1327 t = tab_create (nc, nr);
1329 tab_title (t, _("Total Variance Explained"));
1331 tab_headers (t, heading_columns, 0, heading_rows, 0);
1333 /* Outline the box */
1340 /* Vertical lines */
1347 tab_hline (t, TAL_1, 0, nc - 1, heading_rows);
1348 tab_hline (t, TAL_1, 1, nc - 1, 1);
1350 tab_vline (t, TAL_2, heading_columns, 0, nr - 1);
1353 if ( factor->extraction == EXTRACTION_PC)
1354 tab_text (t, 0, 1, TAB_LEFT | TAT_TITLE, _("Component"));
1356 tab_text (t, 0, 1, TAB_LEFT | TAT_TITLE, _("Factor"));
1359 if (factor->print & PRINT_INITIAL)
1361 tab_joint_text (t, c, 0, c + 2, 0, TAB_CENTER | TAT_TITLE, _("Initial Eigenvalues"));
1365 if (factor->print & PRINT_EXTRACTION)
1367 tab_joint_text (t, c, 0, c + 2, 0, TAB_CENTER | TAT_TITLE, _("Extraction Sums of Squared Loadings"));
1371 if (factor->print & PRINT_ROTATION)
1373 tab_joint_text (t, c, 0, c + 2, 0, TAB_CENTER | TAT_TITLE, _("Rotation Sums of Squared Loadings"));
1377 for (i = 0; i < (nc - heading_columns) / 3 ; ++i)
1379 tab_text (t, i * 3 + 1, 1, TAB_CENTER | TAT_TITLE, _("Total"));
1380 /* xgettext:no-c-format */
1381 tab_text (t, i * 3 + 2, 1, TAB_CENTER | TAT_TITLE, _("% of Variance"));
1382 tab_text (t, i * 3 + 3, 1, TAB_CENTER | TAT_TITLE, _("Cumulative %"));
1384 tab_vline (t, TAL_2, heading_columns + i * 3, 0, nr - 1);
1387 for (i = 0 ; i < initial_eigenvalues->size; ++i)
1388 i_total += gsl_vector_get (initial_eigenvalues, i);
1390 if ( factor->extraction == EXTRACTION_PAF)
1392 e_total = factor->n_vars;
1400 for (i = 0 ; i < factor->n_vars; ++i)
1402 const double i_lambda = gsl_vector_get (initial_eigenvalues, i);
1403 double i_percent = 100.0 * i_lambda / i_total ;
1405 const double e_lambda = gsl_vector_get (extracted_eigenvalues, i);
1406 double e_percent = 100.0 * e_lambda / e_total ;
1410 tab_text_format (t, c++, i + heading_rows, TAB_LEFT | TAT_TITLE, _("%d"), i + 1);
1415 /* Initial Eigenvalues */
1416 if (factor->print & PRINT_INITIAL)
1418 tab_double (t, c++, i + heading_rows, 0, i_lambda, NULL);
1419 tab_double (t, c++, i + heading_rows, 0, i_percent, NULL);
1420 tab_double (t, c++, i + heading_rows, 0, i_cum, NULL);
1423 if (factor->print & PRINT_EXTRACTION)
1425 if (i < idata->n_extractions)
1427 /* Sums of squared loadings */
1428 tab_double (t, c++, i + heading_rows, 0, e_lambda, NULL);
1429 tab_double (t, c++, i + heading_rows, 0, e_percent, NULL);
1430 tab_double (t, c++, i + heading_rows, 0, e_cum, NULL);
1440 show_correlation_matrix (const struct cmd_factor *factor, const struct idata *idata)
1442 struct tab_table *t ;
1444 int y_pos_corr = -1;
1446 int suffix_rows = 0;
1448 const int heading_rows = 1;
1449 const int heading_columns = 2;
1451 int nc = heading_columns ;
1452 int nr = heading_rows ;
1453 int n_data_sets = 0;
1455 if (factor->print & PRINT_CORRELATION)
1457 y_pos_corr = n_data_sets;
1459 nc = heading_columns + factor->n_vars;
1462 if (factor->print & PRINT_SIG)
1464 y_pos_sig = n_data_sets;
1466 nc = heading_columns + factor->n_vars;
1469 nr += n_data_sets * factor->n_vars;
1471 if (factor->print & PRINT_DETERMINANT)
1474 /* If the table would contain only headings, don't bother rendering it */
1475 if (nr <= heading_rows && suffix_rows == 0)
1478 t = tab_create (nc, nr + suffix_rows);
1480 tab_title (t, _("Correlation Matrix"));
1482 tab_hline (t, TAL_1, 0, nc - 1, heading_rows);
1484 if (nr > heading_rows)
1486 tab_headers (t, heading_columns, 0, heading_rows, 0);
1488 tab_vline (t, TAL_2, 2, 0, nr - 1);
1490 /* Outline the box */
1497 /* Vertical lines */
1505 for (i = 0; i < factor->n_vars; ++i)
1506 tab_text (t, heading_columns + i, 0, TAT_TITLE, var_to_string (factor->vars[i]));
1509 for (i = 0 ; i < n_data_sets; ++i)
1511 int y = heading_rows + i * factor->n_vars;
1513 for (v = 0; v < factor->n_vars; ++v)
1514 tab_text (t, 1, y + v, TAT_TITLE, var_to_string (factor->vars[v]));
1516 tab_hline (t, TAL_1, 0, nc - 1, y);
1519 if (factor->print & PRINT_CORRELATION)
1521 const double y = heading_rows + y_pos_corr;
1522 tab_text (t, 0, y, TAT_TITLE, _("Correlations"));
1524 for (i = 0; i < factor->n_vars; ++i)
1526 for (j = 0; j < factor->n_vars; ++j)
1527 tab_double (t, heading_columns + i, y + j, 0, gsl_matrix_get (idata->corr, i, j), NULL);
1531 if (factor->print & PRINT_SIG)
1533 const double y = heading_rows + y_pos_sig * factor->n_vars;
1534 tab_text (t, 0, y, TAT_TITLE, _("Sig. 1-tailed"));
1536 for (i = 0; i < factor->n_vars; ++i)
1538 for (j = 0; j < factor->n_vars; ++j)
1540 double rho = gsl_matrix_get (idata->corr, i, j);
1541 double w = gsl_matrix_get (idata->n, i, j);
1546 tab_double (t, heading_columns + i, y + j, 0, significance_of_correlation (rho, w), NULL);
1552 if (factor->print & PRINT_DETERMINANT)
1557 const int size = idata->corr->size1;
1558 gsl_permutation *p = gsl_permutation_calloc (size);
1559 gsl_matrix *tmp = gsl_matrix_calloc (size, size);
1560 gsl_matrix_memcpy (tmp, idata->corr);
1562 gsl_linalg_LU_decomp (tmp, p, &sign);
1563 det = gsl_linalg_LU_det (tmp, sign);
1564 gsl_permutation_free (p);
1565 gsl_matrix_free (tmp);
1568 tab_text (t, 0, nr, TAB_LEFT | TAT_TITLE, _("Determinant"));
1569 tab_double (t, 1, nr, 0, det, NULL);
1578 do_factor (const struct cmd_factor *factor, struct casereader *r)
1581 const gsl_matrix *var_matrix;
1582 const gsl_matrix *mean_matrix;
1584 const gsl_matrix *analysis_matrix;
1585 struct idata *idata = idata_alloc (factor->n_vars);
1587 struct covariance *cov = covariance_create (factor->n_vars, factor->vars,
1588 factor->wv, factor->exclude);
1590 for ( ; (c = casereader_read (r) ); case_unref (c))
1592 covariance_accumulate (cov, c);
1595 idata->cov = covariance_calculate (cov);
1597 var_matrix = covariance_moments (cov, MOMENT_VARIANCE);
1598 mean_matrix = covariance_moments (cov, MOMENT_MEAN);
1599 idata->n = covariance_moments (cov, MOMENT_NONE);
1601 if ( factor->method == METHOD_CORR)
1603 idata->corr = correlation_from_covariance (idata->cov, var_matrix);
1604 analysis_matrix = idata->corr;
1607 analysis_matrix = idata->cov;
1609 if ( factor->print & PRINT_UNIVARIATE)
1613 const struct fmt_spec *wfmt = factor->wv ? var_get_print_format (factor->wv) : & F_8_0;
1616 const int heading_columns = 1;
1617 const int heading_rows = 1;
1619 const int nr = heading_rows + factor->n_vars;
1621 struct tab_table *t = tab_create (nc, nr);
1622 tab_title (t, _("Descriptive Statistics"));
1624 tab_headers (t, heading_columns, 0, heading_rows, 0);
1626 /* Outline the box */
1633 /* Vertical lines */
1640 tab_hline (t, TAL_1, 0, nc - 1, heading_rows);
1641 tab_vline (t, TAL_2, heading_columns, 0, nr - 1);
1643 tab_text (t, 1, 0, TAB_CENTER | TAT_TITLE, _("Mean"));
1644 tab_text (t, 2, 0, TAB_CENTER | TAT_TITLE, _("Std. Deviation"));
1645 tab_text (t, 3, 0, TAB_CENTER | TAT_TITLE, _("Analysis N"));
1647 for (i = 0 ; i < factor->n_vars; ++i)
1649 const struct variable *v = factor->vars[i];
1650 tab_text (t, 0, i + heading_rows, TAB_LEFT | TAT_TITLE, var_to_string (v));
1652 tab_double (t, 1, i + heading_rows, 0, gsl_matrix_get (mean_matrix, i, i), NULL);
1653 tab_double (t, 2, i + heading_rows, 0, sqrt (gsl_matrix_get (var_matrix, i, i)), NULL);
1654 tab_double (t, 3, i + heading_rows, 0, gsl_matrix_get (idata->n, i, i), wfmt);
1660 show_correlation_matrix (factor, idata);
1664 gsl_eigen_symmv_workspace *workspace = gsl_eigen_symmv_alloc (factor->n_vars);
1666 gsl_eigen_symmv (matrix_dup (analysis_matrix), idata->eval, idata->evec, workspace);
1668 gsl_eigen_symmv_free (workspace);
1671 gsl_eigen_symmv_sort (idata->eval, idata->evec, GSL_EIGEN_SORT_ABS_DESC);
1674 idata->n_extractions = n_extracted_factors (factor, idata);
1676 if (idata->n_extractions == 0)
1678 msg (MW, _("The FACTOR criteria result in zero factors extracted. Therefore no analysis will be performed."));
1682 if (idata->n_extractions > factor->n_vars)
1684 msg (MW, _("The FACTOR criteria result in more factors than variables, which is not meaningful. No analysis will be performed."));
1689 const gsl_vector *extracted_eigenvalues = NULL;
1690 gsl_vector *initial_communalities = gsl_vector_alloc (factor->n_vars);
1691 gsl_vector *extracted_communalities = gsl_vector_alloc (factor->n_vars);
1693 struct factor_matrix_workspace *fmw = factor_matrix_workspace_alloc (idata->msr->size, idata->n_extractions);
1694 gsl_matrix *factor_matrix = gsl_matrix_calloc (factor->n_vars, fmw->n_factors);
1696 if ( factor->extraction == EXTRACTION_PAF)
1698 gsl_vector *diff = gsl_vector_alloc (idata->msr->size);
1699 struct smr_workspace *ws = ws_create (analysis_matrix);
1701 for (i = 0 ; i < factor->n_vars ; ++i)
1703 double r2 = squared_multiple_correlation (analysis_matrix, i, ws);
1705 gsl_vector_set (idata->msr, i, r2);
1709 gsl_vector_memcpy (initial_communalities, idata->msr);
1711 for (i = 0; i < factor->iterations; ++i)
1714 gsl_vector_memcpy (diff, idata->msr);
1716 iterate_factor_matrix (analysis_matrix, idata->msr, factor_matrix, fmw);
1718 gsl_vector_sub (diff, idata->msr);
1720 gsl_vector_minmax (diff, &min, &max);
1722 if ( fabs (min) < factor->econverge && fabs (max) < factor->econverge)
1725 gsl_vector_free (diff);
1729 gsl_vector_memcpy (extracted_communalities, idata->msr);
1730 extracted_eigenvalues = fmw->eval;
1732 else if (factor->extraction == EXTRACTION_PC)
1734 for (i = 0; i < factor->n_vars; ++i)
1735 gsl_vector_set (initial_communalities, i, communality (idata, i, factor->n_vars));
1737 gsl_vector_memcpy (extracted_communalities, initial_communalities);
1739 iterate_factor_matrix (analysis_matrix, extracted_communalities, factor_matrix, fmw);
1742 extracted_eigenvalues = idata->eval;
1746 show_communalities (factor, initial_communalities, extracted_communalities);
1748 show_explained_variance (factor, idata, idata->eval, extracted_eigenvalues);
1750 factor_matrix_workspace_free (fmw);
1752 show_scree (factor, idata);
1754 show_factor_matrix (factor, idata,
1755 factor->extraction == EXTRACTION_PC ? _("Component Matrix") : _("Factor Matrix"),
1758 if ( factor->rotation != ROT_NONE)
1760 gsl_matrix *rotated_factors = gsl_matrix_calloc (factor_matrix->size1, factor_matrix->size2);
1762 rotate (factor_matrix, extracted_communalities, factor->rotation, rotated_factors);
1764 show_factor_matrix (factor, idata,
1765 factor->extraction == EXTRACTION_PC ? _("Rotated Component Matrix") : _("Rotated Factor Matrix"),
1768 gsl_matrix_free (rotated_factors);
1772 gsl_vector_free (initial_communalities);
1773 gsl_vector_free (extracted_communalities);
1780 casereader_destroy (r);