1 /* PSPP - a program for statistical analysis.
2 Copyright (C) 2009, 2010, 2011, 2012, 2014, 2015,
3 2016, 2017 Free Software Foundation, Inc.
5 This program is free software: you can redistribute it and/or modify
6 it under the terms of the GNU General Public License as published by
7 the Free Software Foundation, either version 3 of the License, or
8 (at your option) any later version.
10 This program is distributed in the hope that it will be useful,
11 but WITHOUT ANY WARRANTY; without even the implied warranty of
12 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13 GNU General Public License for more details.
15 You should have received a copy of the GNU General Public License
16 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>
26 #include <gsl/gsl_cdf.h>
28 #include "data/any-reader.h"
29 #include "data/casegrouper.h"
30 #include "data/casereader.h"
31 #include "data/casewriter.h"
32 #include "data/dataset.h"
33 #include "data/dictionary.h"
34 #include "data/format.h"
35 #include "data/subcase.h"
36 #include "language/command.h"
37 #include "language/lexer/lexer.h"
38 #include "language/lexer/value-parser.h"
39 #include "language/lexer/variable-parser.h"
40 #include "language/data-io/file-handle.h"
41 #include "language/data-io/matrix-reader.h"
42 #include "libpspp/cast.h"
43 #include "libpspp/message.h"
44 #include "libpspp/misc.h"
45 #include "math/correlation.h"
46 #include "math/covariance.h"
47 #include "math/moments.h"
48 #include "output/charts/scree.h"
49 #include "output/pivot-table.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,
107 typedef void (*rotation_coefficients) (double *x, double *y,
108 double a, double b, double c, double d,
109 const gsl_matrix *loadings);
113 varimax_coefficients (double *x, double *y,
114 double a, double b, double c, double d,
115 const gsl_matrix *loadings)
117 *x = d - 2 * a * b / loadings->size1;
118 *y = c - (a * a - b * b) / loadings->size1;
122 equamax_coefficients (double *x, double *y,
123 double a, double b, double c, double d,
124 const gsl_matrix *loadings)
126 *x = d - loadings->size2 * a * b / loadings->size1;
127 *y = c - loadings->size2 * (a * a - b * b) / (2 * loadings->size1);
131 quartimax_coefficients (double *x, double *y,
132 double a UNUSED, double b UNUSED, double c, double d,
133 const gsl_matrix *loadings UNUSED)
139 static const rotation_coefficients rotation_coeff[] = {
140 varimax_coefficients,
141 equamax_coefficients,
142 quartimax_coefficients,
143 varimax_coefficients /* PROMAX is identical to VARIMAX */
147 /* return diag (C'C) ^ {-0.5} */
149 diag_rcp_sqrt (const gsl_matrix *C)
152 gsl_matrix *d = gsl_matrix_calloc (C->size1, C->size2);
153 gsl_matrix *r = gsl_matrix_calloc (C->size1, C->size2);
155 assert (C->size1 == C->size2);
157 gsl_linalg_matmult_mod (C, GSL_LINALG_MOD_TRANSPOSE,
158 C, GSL_LINALG_MOD_NONE,
161 for (j = 0 ; j < d->size2; ++j)
163 double e = gsl_matrix_get (d, j, j);
165 gsl_matrix_set (r, j, j, e);
175 /* return diag ((C'C)^-1) ^ {-0.5} */
177 diag_rcp_inv_sqrt (const gsl_matrix *CCinv)
180 gsl_matrix *r = gsl_matrix_calloc (CCinv->size1, CCinv->size2);
182 assert (CCinv->size1 == CCinv->size2);
184 for (j = 0 ; j < CCinv->size2; ++j)
186 double e = gsl_matrix_get (CCinv, j, j);
188 gsl_matrix_set (r, j, j, e);
201 const struct variable **vars;
203 const struct variable *wv;
206 enum missing_type missing_type;
207 enum mv_class exclude;
208 enum print_opts print;
209 enum extraction_method extraction;
211 enum rotation_type rotation;
212 int rotation_iterations;
215 /* Extraction Criteria */
219 int extraction_iterations;
231 /* Intermediate values used in calculation */
232 struct matrix_material mm;
234 gsl_matrix *analysis_matrix; /* A pointer to either mm.corr or mm.cov */
236 gsl_vector *eval ; /* The eigenvalues */
237 gsl_matrix *evec ; /* The eigenvectors */
241 gsl_vector *msr ; /* Multiple Squared Regressions */
243 double detR; /* The determinant of the correlation matrix */
245 gsl_matrix *ai_cov; /* The anti-image covariance matrix */
246 gsl_matrix *ai_cor; /* The anti-image correlation matrix */
247 struct covariance *cvm;
250 static struct idata *
251 idata_alloc (size_t n_vars)
253 struct idata *id = xzalloc (sizeof (*id));
255 id->n_extractions = 0;
256 id->msr = gsl_vector_alloc (n_vars);
258 id->eval = gsl_vector_alloc (n_vars);
259 id->evec = gsl_matrix_alloc (n_vars, n_vars);
265 idata_free (struct idata *id)
267 gsl_vector_free (id->msr);
268 gsl_vector_free (id->eval);
269 gsl_matrix_free (id->evec);
270 gsl_matrix_free (id->ai_cov);
271 gsl_matrix_free (id->ai_cor);
276 /* Return the sum of squares of all the elements in row J excluding column J */
278 ssq_row_od_n (const gsl_matrix *m, int j)
282 assert (m->size1 == m->size2);
284 assert (j < m->size1);
286 for (i = 0; i < m->size1; ++i)
288 if (i == j) continue;
289 ss += pow2 (gsl_matrix_get (m, i, j));
295 /* Return the sum of squares of all the elements excluding row N */
297 ssq_od_n (const gsl_matrix *m, int n)
301 assert (m->size1 == m->size2);
303 assert (n < m->size1);
305 for (i = 0; i < m->size1; ++i)
307 for (j = 0; j < m->size2; ++j)
309 if (i == j) continue;
310 ss += pow2 (gsl_matrix_get (m, i, j));
319 anti_image_corr (const gsl_matrix *m, const struct idata *idata)
323 assert (m->size1 == m->size2);
325 a = gsl_matrix_alloc (m->size1, m->size2);
327 for (i = 0; i < m->size1; ++i)
329 for (j = 0; j < m->size2; ++j)
331 double *p = gsl_matrix_ptr (a, i, j);
332 *p = gsl_matrix_get (m, i, j);
333 *p /= sqrt (gsl_matrix_get (m, i, i) *
334 gsl_matrix_get (m, j, j));
338 for (i = 0; i < m->size1; ++i)
340 double r = ssq_row_od_n (idata->mm.corr, i);
341 double u = ssq_row_od_n (a, i);
342 gsl_matrix_set (a, i, i, r / (r + u));
349 anti_image_cov (const gsl_matrix *m)
353 assert (m->size1 == m->size2);
355 a = gsl_matrix_alloc (m->size1, m->size2);
357 for (i = 0; i < m->size1; ++i)
359 for (j = 0; j < m->size2; ++j)
361 double *p = gsl_matrix_ptr (a, i, j);
362 *p = gsl_matrix_get (m, i, j);
363 *p /= gsl_matrix_get (m, i, i);
364 *p /= gsl_matrix_get (m, j, j);
373 dump_matrix (const gsl_matrix *m)
377 for (i = 0 ; i < m->size1; ++i)
379 for (j = 0 ; j < m->size2; ++j)
380 printf ("%02f ", gsl_matrix_get (m, i, j));
386 dump_matrix_permute (const gsl_matrix *m, const gsl_permutation *p)
390 for (i = 0 ; i < m->size1; ++i)
392 for (j = 0 ; j < m->size2; ++j)
393 printf ("%02f ", gsl_matrix_get (m, gsl_permutation_get (p, i), j));
400 dump_vector (const gsl_vector *v)
403 for (i = 0 ; i < v->size; ++i)
405 printf ("%02f\n", gsl_vector_get (v, i));
413 n_extracted_factors (const struct cmd_factor *factor, struct idata *idata)
417 /* If there is a cached value, then return that. */
418 if (idata->n_extractions != 0)
419 return idata->n_extractions;
421 /* Otherwise, if the number of factors has been explicitly requested,
423 if (factor->n_factors > 0)
425 idata->n_extractions = factor->n_factors;
429 /* Use the MIN_EIGEN setting. */
430 for (i = 0 ; i < idata->eval->size; ++i)
432 double evali = fabs (gsl_vector_get (idata->eval, i));
434 idata->n_extractions = i;
436 if (evali < factor->min_eigen)
441 return idata->n_extractions;
445 /* Returns a newly allocated matrix identical to M.
446 It it the callers responsibility to free the returned value.
449 matrix_dup (const gsl_matrix *m)
451 gsl_matrix *n = gsl_matrix_alloc (m->size1, m->size2);
453 gsl_matrix_memcpy (n, m);
461 /* Copy of the subject */
466 gsl_permutation *perm;
473 static struct smr_workspace *ws_create (const gsl_matrix *input)
475 struct smr_workspace *ws = xmalloc (sizeof (*ws));
477 ws->m = gsl_matrix_alloc (input->size1, input->size2);
478 ws->inverse = gsl_matrix_calloc (input->size1 - 1, input->size2 - 1);
479 ws->perm = gsl_permutation_alloc (input->size1 - 1);
480 ws->result1 = gsl_matrix_calloc (input->size1 - 1, 1);
481 ws->result2 = gsl_matrix_calloc (1, 1);
487 ws_destroy (struct smr_workspace *ws)
489 gsl_matrix_free (ws->result2);
490 gsl_matrix_free (ws->result1);
491 gsl_permutation_free (ws->perm);
492 gsl_matrix_free (ws->inverse);
493 gsl_matrix_free (ws->m);
500 Return the square of the regression coefficient for VAR regressed against all other variables.
503 squared_multiple_correlation (const gsl_matrix *corr, int var, struct smr_workspace *ws)
505 /* For an explanation of what this is doing, see
506 http://www.visualstatistics.net/Visual%20Statistics%20Multimedia/multiple_regression_analysis.htm
512 gsl_matrix_memcpy (ws->m, corr);
514 gsl_matrix_swap_rows (ws->m, 0, var);
515 gsl_matrix_swap_columns (ws->m, 0, var);
517 rxx = gsl_matrix_submatrix (ws->m, 1, 1, ws->m->size1 - 1, ws->m->size1 - 1);
519 gsl_linalg_LU_decomp (&rxx.matrix, ws->perm, &signum);
521 gsl_linalg_LU_invert (&rxx.matrix, ws->perm, ws->inverse);
524 gsl_matrix_const_view rxy = gsl_matrix_const_submatrix (ws->m, 1, 0, ws->m->size1 - 1, 1);
525 gsl_matrix_const_view ryx = gsl_matrix_const_submatrix (ws->m, 0, 1, 1, ws->m->size1 - 1);
527 gsl_blas_dgemm (CblasNoTrans, CblasNoTrans,
528 1.0, ws->inverse, &rxy.matrix, 0.0, ws->result1);
530 gsl_blas_dgemm (CblasNoTrans, CblasNoTrans,
531 1.0, &ryx.matrix, ws->result1, 0.0, ws->result2);
534 return gsl_matrix_get (ws->result2, 0, 0);
539 static double the_communality (const gsl_matrix *evec, const gsl_vector *eval, int n, int n_factors);
542 struct factor_matrix_workspace
545 gsl_eigen_symmv_workspace *eigen_ws;
555 static struct factor_matrix_workspace *
556 factor_matrix_workspace_alloc (size_t n, size_t nf)
558 struct factor_matrix_workspace *ws = xmalloc (sizeof (*ws));
561 ws->gamma = gsl_matrix_calloc (nf, nf);
562 ws->eigen_ws = gsl_eigen_symmv_alloc (n);
563 ws->eval = gsl_vector_alloc (n);
564 ws->evec = gsl_matrix_alloc (n, n);
565 ws->r = gsl_matrix_alloc (n, n);
571 factor_matrix_workspace_free (struct factor_matrix_workspace *ws)
573 gsl_eigen_symmv_free (ws->eigen_ws);
574 gsl_vector_free (ws->eval);
575 gsl_matrix_free (ws->evec);
576 gsl_matrix_free (ws->gamma);
577 gsl_matrix_free (ws->r);
582 Shift P left by OFFSET places, and overwrite TARGET
583 with the shifted result.
584 Positions in TARGET less than OFFSET are unchanged.
587 perm_shift_apply (gsl_permutation *target, const gsl_permutation *p,
591 assert (target->size == p->size);
592 assert (offset <= target->size);
594 for (i = 0; i < target->size - offset; ++i)
596 target->data[i] = p->data [i + offset];
602 Indirectly sort the rows of matrix INPUT, storing the sort order in PERM.
603 The sort criteria are as follows:
605 Rows are sorted on the first column, until the absolute value of an
606 element in a subsequent column is greater than that of the first
607 column. Thereafter, rows will be sorted on the second column,
608 until the absolute value of an element in a subsequent column
609 exceeds that of the second column ...
612 sort_matrix_indirect (const gsl_matrix *input, gsl_permutation *perm)
614 const size_t n = perm->size;
615 const size_t m = input->size2;
622 assert (perm->size == input->size1);
624 p = gsl_permutation_alloc (n);
626 /* Copy INPUT into MAT, discarding the sign */
627 mat = gsl_matrix_alloc (n, m);
628 for (i = 0 ; i < mat->size1; ++i)
630 for (j = 0 ; j < mat->size2; ++j)
632 double x = gsl_matrix_get (input, i, j);
633 gsl_matrix_set (mat, i, j, fabs (x));
637 while (column_n < m && row_n < n)
639 gsl_vector_const_view columni = gsl_matrix_const_column (mat, column_n);
640 gsl_sort_vector_index (p, &columni.vector);
642 for (i = 0 ; i < n; ++i)
644 gsl_vector_view row = gsl_matrix_row (mat, p->data[n - 1 - i]);
645 size_t maxindex = gsl_vector_max_index (&row.vector);
647 if (maxindex > column_n)
650 /* All subsequent elements of this row, are of no interest.
651 So set them all to a highly negative value */
652 for (j = column_n + 1; j < row.vector.size ; ++j)
653 gsl_vector_set (&row.vector, j, -DBL_MAX);
656 perm_shift_apply (perm, p, row_n);
662 gsl_permutation_free (p);
663 gsl_matrix_free (mat);
665 assert (0 == gsl_permutation_valid (perm));
667 /* We want the biggest value to be first */
668 gsl_permutation_reverse (perm);
673 drot_go (double phi, double *l0, double *l1)
675 double r0 = cos (phi) * *l0 + sin (phi) * *l1;
676 double r1 = - sin (phi) * *l0 + cos (phi) * *l1;
684 clone_matrix (const gsl_matrix *m)
687 gsl_matrix *c = gsl_matrix_calloc (m->size1, m->size2);
689 for (j = 0 ; j < c->size1; ++j)
691 for (k = 0 ; k < c->size2; ++k)
693 const double *v = gsl_matrix_const_ptr (m, j, k);
694 gsl_matrix_set (c, j, k, *v);
703 initial_sv (const gsl_matrix *fm)
708 for (j = 0 ; j < fm->size2; ++j)
713 for (k = j + 1 ; k < fm->size2; ++k)
715 double lambda = gsl_matrix_get (fm, k, j);
716 double lambda_sq = lambda * lambda;
717 double lambda_4 = lambda_sq * lambda_sq;
722 sv += (fm->size1 * l4s - (l2s * l2s)) / (fm->size1 * fm->size1);
728 rotate (const struct cmd_factor *cf, const gsl_matrix *unrot,
729 const gsl_vector *communalities,
731 gsl_vector *rotated_loadings,
732 gsl_matrix *pattern_matrix,
733 gsl_matrix *factor_correlation_matrix
740 /* First get a normalised version of UNROT */
741 gsl_matrix *normalised = gsl_matrix_calloc (unrot->size1, unrot->size2);
742 gsl_matrix *h_sqrt = gsl_matrix_calloc (communalities->size, communalities->size);
743 gsl_matrix *h_sqrt_inv ;
745 /* H is the diagonal matrix containing the absolute values of the communalities */
746 for (i = 0 ; i < communalities->size ; ++i)
748 double *ptr = gsl_matrix_ptr (h_sqrt, i, i);
749 *ptr = fabs (gsl_vector_get (communalities, i));
752 /* Take the square root of the communalities */
753 gsl_linalg_cholesky_decomp (h_sqrt);
756 /* Save a copy of h_sqrt and invert it */
757 h_sqrt_inv = clone_matrix (h_sqrt);
758 gsl_linalg_cholesky_decomp (h_sqrt_inv);
759 gsl_linalg_cholesky_invert (h_sqrt_inv);
761 /* normalised vertion is H^{1/2} x UNROT */
762 gsl_blas_dgemm (CblasNoTrans, CblasNoTrans, 1.0, h_sqrt_inv, unrot, 0.0, normalised);
764 gsl_matrix_free (h_sqrt_inv);
767 /* Now perform the rotation iterations */
769 prev_sv = initial_sv (normalised);
770 for (i = 0 ; i < cf->rotation_iterations ; ++i)
773 for (j = 0 ; j < normalised->size2; ++j)
775 /* These variables relate to the convergence criterium */
779 for (k = j + 1 ; k < normalised->size2; ++k)
789 for (p = 0; p < normalised->size1; ++p)
791 double jv = gsl_matrix_get (normalised, p, j);
792 double kv = gsl_matrix_get (normalised, p, k);
794 double u = jv * jv - kv * kv;
795 double v = 2 * jv * kv;
802 rotation_coeff [cf->rotation] (&x, &y, a, b, c, d, normalised);
804 phi = atan2 (x, y) / 4.0 ;
806 /* Don't bother rotating if the angle is small */
807 if (fabs (sin (phi)) <= pow (10.0, -15.0))
810 for (p = 0; p < normalised->size1; ++p)
812 double *lambda0 = gsl_matrix_ptr (normalised, p, j);
813 double *lambda1 = gsl_matrix_ptr (normalised, p, k);
814 drot_go (phi, lambda0, lambda1);
817 /* Calculate the convergence criterium */
819 double lambda = gsl_matrix_get (normalised, k, j);
820 double lambda_sq = lambda * lambda;
821 double lambda_4 = lambda_sq * lambda_sq;
827 sv += (normalised->size1 * l4s - (l2s * l2s)) / (normalised->size1 * normalised->size1);
830 if (fabs (sv - prev_sv) <= cf->rconverge)
836 gsl_blas_dgemm (CblasNoTrans, CblasNoTrans, 1.0,
837 h_sqrt, normalised, 0.0, result);
839 gsl_matrix_free (h_sqrt);
840 gsl_matrix_free (normalised);
842 if (cf->rotation == ROT_PROMAX)
844 /* general purpose m by m matrix, where m is the number of factors */
845 gsl_matrix *mm1 = gsl_matrix_calloc (unrot->size2, unrot->size2);
846 gsl_matrix *mm2 = gsl_matrix_calloc (unrot->size2, unrot->size2);
848 /* general purpose m by p matrix, where p is the number of variables */
849 gsl_matrix *mp1 = gsl_matrix_calloc (unrot->size2, unrot->size1);
851 gsl_matrix *pm1 = gsl_matrix_calloc (unrot->size1, unrot->size2);
853 gsl_permutation *perm = gsl_permutation_alloc (unrot->size2);
859 /* The following variables follow the notation by SPSS Statistical Algorithms
861 gsl_matrix *L = gsl_matrix_calloc (unrot->size2, unrot->size2);
862 gsl_matrix *P = clone_matrix (result);
867 /* Vector of length p containing (indexed by i)
868 \Sum^m_j {\lambda^2_{ij}} */
869 gsl_vector *rssq = gsl_vector_calloc (unrot->size1);
871 for (i = 0; i < P->size1; ++i)
874 for (j = 0; j < P->size2; ++j)
876 sum += gsl_matrix_get (result, i, j)
877 * gsl_matrix_get (result, i, j);
881 gsl_vector_set (rssq, i, sqrt (sum));
884 for (i = 0; i < P->size1; ++i)
886 for (j = 0; j < P->size2; ++j)
888 double l = gsl_matrix_get (result, i, j);
889 double r = gsl_vector_get (rssq, i);
890 gsl_matrix_set (P, i, j, pow (fabs (l / r), cf->promax_power + 1) * r / l);
894 gsl_vector_free (rssq);
896 gsl_linalg_matmult_mod (result,
897 GSL_LINALG_MOD_TRANSPOSE,
902 gsl_linalg_LU_decomp (mm1, perm, &signum);
903 gsl_linalg_LU_invert (mm1, perm, mm2);
905 gsl_linalg_matmult_mod (mm2, GSL_LINALG_MOD_NONE,
906 result, GSL_LINALG_MOD_TRANSPOSE,
909 gsl_linalg_matmult_mod (mp1, GSL_LINALG_MOD_NONE,
910 P, GSL_LINALG_MOD_NONE,
913 D = diag_rcp_sqrt (L);
914 Q = gsl_matrix_calloc (unrot->size2, unrot->size2);
916 gsl_linalg_matmult_mod (L, GSL_LINALG_MOD_NONE,
917 D, GSL_LINALG_MOD_NONE,
920 gsl_matrix *QQinv = gsl_matrix_calloc (unrot->size2, unrot->size2);
922 gsl_linalg_matmult_mod (Q, GSL_LINALG_MOD_TRANSPOSE,
923 Q, GSL_LINALG_MOD_NONE,
926 gsl_linalg_cholesky_decomp (QQinv);
927 gsl_linalg_cholesky_invert (QQinv);
930 gsl_matrix *C = diag_rcp_inv_sqrt (QQinv);
931 gsl_matrix *Cinv = clone_matrix (C);
933 gsl_linalg_cholesky_decomp (Cinv);
934 gsl_linalg_cholesky_invert (Cinv);
937 gsl_linalg_matmult_mod (result, GSL_LINALG_MOD_NONE,
938 Q, GSL_LINALG_MOD_NONE,
941 gsl_linalg_matmult_mod (pm1, GSL_LINALG_MOD_NONE,
942 Cinv, GSL_LINALG_MOD_NONE,
946 gsl_linalg_matmult_mod (C, GSL_LINALG_MOD_NONE,
947 QQinv, GSL_LINALG_MOD_NONE,
950 gsl_linalg_matmult_mod (mm1, GSL_LINALG_MOD_NONE,
951 C, GSL_LINALG_MOD_TRANSPOSE,
952 factor_correlation_matrix);
954 gsl_linalg_matmult_mod (pattern_matrix, GSL_LINALG_MOD_NONE,
955 factor_correlation_matrix, GSL_LINALG_MOD_NONE,
958 gsl_matrix_memcpy (result, pm1);
961 gsl_matrix_free (QQinv);
963 gsl_matrix_free (Cinv);
970 gsl_permutation_free (perm);
972 gsl_matrix_free (mm1);
973 gsl_matrix_free (mm2);
974 gsl_matrix_free (mp1);
975 gsl_matrix_free (pm1);
979 /* reflect negative sums and populate the rotated loadings vector*/
980 for (i = 0 ; i < result->size2; ++i)
984 for (j = 0 ; j < result->size1; ++j)
986 double s = gsl_matrix_get (result, j, i);
991 gsl_vector_set (rotated_loadings, i, ssq);
994 for (j = 0 ; j < result->size1; ++j)
996 double *lambda = gsl_matrix_ptr (result, j, i);
1004 Get an approximation for the factor matrix into FACTORS, and the communalities into COMMUNALITIES.
1005 R is the matrix to be analysed.
1006 WS is a pointer to a structure which must have been initialised with factor_matrix_workspace_init.
1009 iterate_factor_matrix (const gsl_matrix *r, gsl_vector *communalities, gsl_matrix *factors,
1010 struct factor_matrix_workspace *ws)
1013 gsl_matrix_view mv ;
1015 assert (r->size1 == r->size2);
1016 assert (r->size1 == communalities->size);
1018 assert (factors->size1 == r->size1);
1019 assert (factors->size2 == ws->n_factors);
1021 gsl_matrix_memcpy (ws->r, r);
1023 /* Apply Communalities to diagonal of correlation matrix */
1024 for (i = 0 ; i < communalities->size ; ++i)
1026 double *x = gsl_matrix_ptr (ws->r, i, i);
1027 *x = gsl_vector_get (communalities, i);
1030 gsl_eigen_symmv (ws->r, ws->eval, ws->evec, ws->eigen_ws);
1032 mv = gsl_matrix_submatrix (ws->evec, 0, 0, ws->evec->size1, ws->n_factors);
1034 /* Gamma is the diagonal matrix containing the absolute values of the eigenvalues */
1035 for (i = 0 ; i < ws->n_factors ; ++i)
1037 double *ptr = gsl_matrix_ptr (ws->gamma, i, i);
1038 *ptr = fabs (gsl_vector_get (ws->eval, i));
1041 /* Take the square root of gamma */
1042 gsl_linalg_cholesky_decomp (ws->gamma);
1044 gsl_blas_dgemm (CblasNoTrans, CblasNoTrans, 1.0, &mv.matrix, ws->gamma, 0.0, factors);
1046 for (i = 0 ; i < r->size1 ; ++i)
1048 double h = the_communality (ws->evec, ws->eval, i, ws->n_factors);
1049 gsl_vector_set (communalities, i, h);
1055 static bool run_factor (struct dataset *ds, const struct cmd_factor *factor);
1057 static void do_factor_by_matrix (const struct cmd_factor *factor, struct idata *idata);
1062 cmd_factor (struct lexer *lexer, struct dataset *ds)
1064 struct dictionary *dict = NULL;
1065 int n_iterations = 25;
1066 struct cmd_factor factor;
1069 factor.method = METHOD_CORR;
1070 factor.missing_type = MISS_LISTWISE;
1071 factor.exclude = MV_ANY;
1072 factor.print = PRINT_INITIAL | PRINT_EXTRACTION | PRINT_ROTATION;
1073 factor.extraction = EXTRACTION_PC;
1074 factor.n_factors = 0;
1075 factor.min_eigen = SYSMIS;
1076 factor.extraction_iterations = 25;
1077 factor.rotation_iterations = 25;
1078 factor.econverge = 0.001;
1081 factor.sort = false;
1083 factor.rotation = ROT_VARIMAX;
1086 factor.rconverge = 0.0001;
1088 lex_match (lexer, T_SLASH);
1090 struct matrix_reader *mr = NULL;
1091 struct casereader *matrix_reader = NULL;
1093 if (lex_match_id (lexer, "VARIABLES"))
1095 lex_match (lexer, T_EQUALS);
1096 dict = dataset_dict (ds);
1097 factor.wv = dict_get_weight (dict);
1099 if (!parse_variables_const (lexer, dict, &factor.vars, &factor.n_vars,
1100 PV_NO_DUPLICATE | PV_NUMERIC))
1103 else if (lex_match_id (lexer, "MATRIX"))
1105 lex_match (lexer, T_EQUALS);
1106 if (! lex_force_match_id (lexer, "IN"))
1108 if (!lex_force_match (lexer, T_LPAREN))
1112 if (lex_match_id (lexer, "CORR"))
1115 else if (lex_match_id (lexer, "COV"))
1120 lex_error (lexer, _("Matrix input for %s must be either COV or CORR"), "FACTOR");
1123 if (! lex_force_match (lexer, T_EQUALS))
1125 if (lex_match (lexer, T_ASTERISK))
1127 dict = dataset_dict (ds);
1128 matrix_reader = casereader_clone (dataset_source (ds));
1132 struct file_handle *fh = fh_parse (lexer, FH_REF_FILE, NULL);
1137 = any_reader_open_and_decode (fh, NULL, &dict, NULL);
1139 if (! (matrix_reader && dict))
1145 if (! lex_force_match (lexer, T_RPAREN))
1148 mr = create_matrix_reader_from_case_reader (dict, matrix_reader,
1149 &factor.vars, &factor.n_vars);
1156 while (lex_token (lexer) != T_ENDCMD)
1158 lex_match (lexer, T_SLASH);
1160 if (lex_match_id (lexer, "ANALYSIS"))
1162 struct const_var_set *vs;
1163 const struct variable **vars;
1167 lex_match (lexer, T_EQUALS);
1169 vs = const_var_set_create_from_array (factor.vars, factor.n_vars);
1170 ok = parse_const_var_set_vars (lexer, vs, &vars, &n_vars,
1171 PV_NO_DUPLICATE | PV_NUMERIC);
1172 const_var_set_destroy (vs);
1179 factor.n_vars = n_vars;
1181 else if (lex_match_id (lexer, "PLOT"))
1183 lex_match (lexer, T_EQUALS);
1184 while (lex_token (lexer) != T_ENDCMD && lex_token (lexer) != T_SLASH)
1186 if (lex_match_id (lexer, "EIGEN"))
1188 factor.plot |= PLOT_SCREE;
1190 #if FACTOR_FULLY_IMPLEMENTED
1191 else if (lex_match_id (lexer, "ROTATION"))
1197 lex_error (lexer, NULL);
1202 else if (lex_match_id (lexer, "METHOD"))
1204 lex_match (lexer, T_EQUALS);
1205 while (lex_token (lexer) != T_ENDCMD && lex_token (lexer) != T_SLASH)
1207 if (lex_match_id (lexer, "COVARIANCE"))
1209 factor.method = METHOD_COV;
1211 else if (lex_match_id (lexer, "CORRELATION"))
1213 factor.method = METHOD_CORR;
1217 lex_error (lexer, NULL);
1222 else if (lex_match_id (lexer, "ROTATION"))
1224 lex_match (lexer, T_EQUALS);
1225 while (lex_token (lexer) != T_ENDCMD && lex_token (lexer) != T_SLASH)
1227 /* VARIMAX and DEFAULT are defaults */
1228 if (lex_match_id (lexer, "VARIMAX") || lex_match_id (lexer, "DEFAULT"))
1230 factor.rotation = ROT_VARIMAX;
1232 else if (lex_match_id (lexer, "EQUAMAX"))
1234 factor.rotation = ROT_EQUAMAX;
1236 else if (lex_match_id (lexer, "QUARTIMAX"))
1238 factor.rotation = ROT_QUARTIMAX;
1240 else if (lex_match_id (lexer, "PROMAX"))
1242 factor.promax_power = 5;
1243 if (lex_match (lexer, T_LPAREN)
1244 && lex_force_int (lexer))
1246 factor.promax_power = lex_integer (lexer);
1248 if (! lex_force_match (lexer, T_RPAREN))
1251 factor.rotation = ROT_PROMAX;
1253 else if (lex_match_id (lexer, "NOROTATE"))
1255 factor.rotation = ROT_NONE;
1259 lex_error (lexer, NULL);
1263 factor.rotation_iterations = n_iterations;
1265 else if (lex_match_id (lexer, "CRITERIA"))
1267 lex_match (lexer, T_EQUALS);
1268 while (lex_token (lexer) != T_ENDCMD && lex_token (lexer) != T_SLASH)
1270 if (lex_match_id (lexer, "FACTORS"))
1272 if (lex_force_match (lexer, T_LPAREN)
1273 && lex_force_int (lexer))
1275 factor.n_factors = lex_integer (lexer);
1277 if (! lex_force_match (lexer, T_RPAREN))
1281 else if (lex_match_id (lexer, "MINEIGEN"))
1283 if (lex_force_match (lexer, T_LPAREN)
1284 && lex_force_num (lexer))
1286 factor.min_eigen = lex_number (lexer);
1288 if (! lex_force_match (lexer, T_RPAREN))
1292 else if (lex_match_id (lexer, "ECONVERGE"))
1294 if (lex_force_match (lexer, T_LPAREN)
1295 && lex_force_num (lexer))
1297 factor.econverge = lex_number (lexer);
1299 if (! lex_force_match (lexer, T_RPAREN))
1303 else if (lex_match_id (lexer, "RCONVERGE"))
1305 if (lex_force_match (lexer, T_LPAREN)
1306 && lex_force_num (lexer))
1308 factor.rconverge = lex_number (lexer);
1310 if (! lex_force_match (lexer, T_RPAREN))
1314 else if (lex_match_id (lexer, "ITERATE"))
1316 if (lex_force_match (lexer, T_LPAREN)
1317 && lex_force_int_range (lexer, "ITERATE", 0, INT_MAX))
1319 n_iterations = lex_integer (lexer);
1321 if (! lex_force_match (lexer, T_RPAREN))
1325 else if (lex_match_id (lexer, "DEFAULT"))
1327 factor.n_factors = 0;
1328 factor.min_eigen = 1;
1333 lex_error (lexer, NULL);
1338 else if (lex_match_id (lexer, "EXTRACTION"))
1340 lex_match (lexer, T_EQUALS);
1341 while (lex_token (lexer) != T_ENDCMD && lex_token (lexer) != T_SLASH)
1343 if (lex_match_id (lexer, "PAF"))
1345 factor.extraction = EXTRACTION_PAF;
1347 else if (lex_match_id (lexer, "PC"))
1349 factor.extraction = EXTRACTION_PC;
1351 else if (lex_match_id (lexer, "PA1"))
1353 factor.extraction = EXTRACTION_PC;
1355 else if (lex_match_id (lexer, "DEFAULT"))
1357 factor.extraction = EXTRACTION_PC;
1361 lex_error (lexer, NULL);
1365 factor.extraction_iterations = n_iterations;
1367 else if (lex_match_id (lexer, "FORMAT"))
1369 lex_match (lexer, T_EQUALS);
1370 while (lex_token (lexer) != T_ENDCMD && lex_token (lexer) != T_SLASH)
1372 if (lex_match_id (lexer, "SORT"))
1376 else if (lex_match_id (lexer, "BLANK"))
1378 if (lex_force_match (lexer, T_LPAREN)
1379 && lex_force_num (lexer))
1381 factor.blank = lex_number (lexer);
1383 if (! lex_force_match (lexer, T_RPAREN))
1387 else if (lex_match_id (lexer, "DEFAULT"))
1390 factor.sort = false;
1394 lex_error (lexer, NULL);
1399 else if (lex_match_id (lexer, "PRINT"))
1402 lex_match (lexer, T_EQUALS);
1403 while (lex_token (lexer) != T_ENDCMD && lex_token (lexer) != T_SLASH)
1405 if (lex_match_id (lexer, "UNIVARIATE"))
1407 factor.print |= PRINT_UNIVARIATE;
1409 else if (lex_match_id (lexer, "DET"))
1411 factor.print |= PRINT_DETERMINANT;
1413 #if FACTOR_FULLY_IMPLEMENTED
1414 else if (lex_match_id (lexer, "INV"))
1418 else if (lex_match_id (lexer, "AIC"))
1420 factor.print |= PRINT_AIC;
1422 else if (lex_match_id (lexer, "SIG"))
1424 factor.print |= PRINT_SIG;
1426 else if (lex_match_id (lexer, "CORRELATION"))
1428 factor.print |= PRINT_CORRELATION;
1430 else if (lex_match_id (lexer, "COVARIANCE"))
1432 factor.print |= PRINT_COVARIANCE;
1434 else if (lex_match_id (lexer, "ROTATION"))
1436 factor.print |= PRINT_ROTATION;
1438 else if (lex_match_id (lexer, "EXTRACTION"))
1440 factor.print |= PRINT_EXTRACTION;
1442 else if (lex_match_id (lexer, "INITIAL"))
1444 factor.print |= PRINT_INITIAL;
1446 else if (lex_match_id (lexer, "KMO"))
1448 factor.print |= PRINT_KMO;
1450 #if FACTOR_FULLY_IMPLEMENTED
1451 else if (lex_match_id (lexer, "REPR"))
1454 else if (lex_match_id (lexer, "FSCORE"))
1458 else if (lex_match (lexer, T_ALL))
1460 factor.print = 0xFFFF;
1462 else if (lex_match_id (lexer, "DEFAULT"))
1464 factor.print |= PRINT_INITIAL ;
1465 factor.print |= PRINT_EXTRACTION ;
1466 factor.print |= PRINT_ROTATION ;
1470 lex_error (lexer, NULL);
1475 else if (lex_match_id (lexer, "MISSING"))
1477 lex_match (lexer, T_EQUALS);
1478 while (lex_token (lexer) != T_ENDCMD && lex_token (lexer) != T_SLASH)
1480 if (lex_match_id (lexer, "INCLUDE"))
1482 factor.exclude = MV_SYSTEM;
1484 else if (lex_match_id (lexer, "EXCLUDE"))
1486 factor.exclude = MV_ANY;
1488 else if (lex_match_id (lexer, "LISTWISE"))
1490 factor.missing_type = MISS_LISTWISE;
1492 else if (lex_match_id (lexer, "PAIRWISE"))
1494 factor.missing_type = MISS_PAIRWISE;
1496 else if (lex_match_id (lexer, "MEANSUB"))
1498 factor.missing_type = MISS_MEANSUB;
1502 lex_error (lexer, NULL);
1509 lex_error (lexer, NULL);
1514 if (factor.rotation == ROT_NONE)
1515 factor.print &= ~PRINT_ROTATION;
1517 if (factor.n_vars < 2)
1518 msg (MW, _("Factor analysis on a single variable is not useful."));
1520 if (factor.n_vars < 1)
1522 msg (ME, _("Factor analysis without variables is not possible."));
1528 struct idata *id = idata_alloc (factor.n_vars);
1530 while (next_matrix_from_reader (&id->mm, mr,
1531 factor.vars, factor.n_vars))
1533 do_factor_by_matrix (&factor, id);
1535 gsl_matrix_free (id->ai_cov);
1537 gsl_matrix_free (id->ai_cor);
1539 gsl_matrix_free (id->mm.corr);
1541 gsl_matrix_free (id->mm.cov);
1548 if (! run_factor (ds, &factor))
1552 destroy_matrix_reader (mr);
1557 destroy_matrix_reader (mr);
1562 static void do_factor (const struct cmd_factor *factor, struct casereader *group);
1566 run_factor (struct dataset *ds, const struct cmd_factor *factor)
1568 struct dictionary *dict = dataset_dict (ds);
1570 struct casereader *group;
1572 struct casegrouper *grouper = casegrouper_create_splits (proc_open (ds), dict);
1574 while (casegrouper_get_next_group (grouper, &group))
1576 if (factor->missing_type == MISS_LISTWISE)
1577 group = casereader_create_filter_missing (group, factor->vars, factor->n_vars,
1580 do_factor (factor, group);
1583 ok = casegrouper_destroy (grouper);
1584 ok = proc_commit (ds) && ok;
1590 /* Return the communality of variable N, calculated to N_FACTORS */
1592 the_communality (const gsl_matrix *evec, const gsl_vector *eval, int n, int n_factors)
1599 assert (n < eval->size);
1600 assert (n < evec->size1);
1601 assert (n_factors <= eval->size);
1603 for (i = 0 ; i < n_factors; ++i)
1605 double evali = fabs (gsl_vector_get (eval, i));
1607 double eveci = gsl_matrix_get (evec, n, i);
1609 comm += pow2 (eveci) * evali;
1615 /* Return the communality of variable N, calculated to N_FACTORS */
1617 communality (const struct idata *idata, int n, int n_factors)
1619 return the_communality (idata->evec, idata->eval, n, n_factors);
1624 show_scree (const struct cmd_factor *f, const struct idata *idata)
1629 if (!(f->plot & PLOT_SCREE))
1633 label = f->extraction == EXTRACTION_PC ? _("Component Number") : _("Factor Number");
1635 s = scree_create (idata->eval, label);
1641 show_communalities (const struct cmd_factor * factor,
1642 const gsl_vector *initial, const gsl_vector *extracted)
1644 if (!(factor->print & (PRINT_INITIAL | PRINT_EXTRACTION)))
1647 struct pivot_table *table = pivot_table_create (N_("Communalities"));
1649 struct pivot_dimension *communalities = pivot_dimension_create (
1650 table, PIVOT_AXIS_COLUMN, N_("Communalities"));
1651 if (factor->print & PRINT_INITIAL)
1652 pivot_category_create_leaves (communalities->root, N_("Initial"));
1653 if (factor->print & PRINT_EXTRACTION)
1654 pivot_category_create_leaves (communalities->root, N_("Extraction"));
1656 struct pivot_dimension *variables = pivot_dimension_create (
1657 table, PIVOT_AXIS_ROW, N_("Variables"));
1659 for (size_t i = 0 ; i < factor->n_vars; ++i)
1661 int row = pivot_category_create_leaf (
1662 variables->root, pivot_value_new_variable (factor->vars[i]));
1665 if (factor->print & PRINT_INITIAL)
1666 pivot_table_put2 (table, col++, row, pivot_value_new_number (
1667 gsl_vector_get (initial, i)));
1668 if (factor->print & PRINT_EXTRACTION)
1669 pivot_table_put2 (table, col++, row, pivot_value_new_number (
1670 gsl_vector_get (extracted, i)));
1673 pivot_table_submit (table);
1676 static struct pivot_dimension *
1677 create_numeric_dimension (struct pivot_table *table,
1678 enum pivot_axis_type axis_type, const char *name,
1679 size_t n, bool show_label)
1681 struct pivot_dimension *d = pivot_dimension_create (table, axis_type, name);
1682 d->root->show_label = show_label;
1683 for (int i = 0 ; i < n; ++i)
1684 pivot_category_create_leaf (d->root, pivot_value_new_integer (i + 1));
1689 show_factor_matrix (const struct cmd_factor *factor, const struct idata *idata, const char *title, const gsl_matrix *fm)
1691 struct pivot_table *table = pivot_table_create (title);
1693 const int n_factors = idata->n_extractions;
1694 create_numeric_dimension (
1695 table, PIVOT_AXIS_COLUMN,
1696 factor->extraction == EXTRACTION_PC ? N_("Component") : N_("Factor"),
1699 struct pivot_dimension *variables = pivot_dimension_create (
1700 table, PIVOT_AXIS_ROW, N_("Variables"));
1702 /* Initialise to the identity permutation */
1703 gsl_permutation *perm = gsl_permutation_calloc (factor->n_vars);
1706 sort_matrix_indirect (fm, perm);
1708 for (size_t i = 0 ; i < factor->n_vars; ++i)
1710 const int matrix_row = perm->data[i];
1712 int var_idx = pivot_category_create_leaf (
1713 variables->root, pivot_value_new_variable (factor->vars[matrix_row]));
1715 for (size_t j = 0 ; j < n_factors; ++j)
1717 double x = gsl_matrix_get (fm, matrix_row, j);
1718 if (fabs (x) < factor->blank)
1721 pivot_table_put2 (table, j, var_idx, pivot_value_new_number (x));
1725 gsl_permutation_free (perm);
1727 pivot_table_submit (table);
1731 put_variance (struct pivot_table *table, int row, int phase_idx,
1732 double lambda, double percent, double cum)
1734 double entries[] = { lambda, percent, cum };
1735 for (size_t i = 0; i < sizeof entries / sizeof *entries; i++)
1736 pivot_table_put3 (table, i, phase_idx, row,
1737 pivot_value_new_number (entries[i]));
1741 show_explained_variance (const struct cmd_factor * factor,
1742 const struct idata *idata,
1743 const gsl_vector *initial_eigenvalues,
1744 const gsl_vector *extracted_eigenvalues,
1745 const gsl_vector *rotated_loadings)
1747 if (!(factor->print & (PRINT_INITIAL | PRINT_EXTRACTION | PRINT_ROTATION)))
1750 struct pivot_table *table = pivot_table_create (
1751 N_("Total Variance Explained"));
1753 pivot_dimension_create (table, PIVOT_AXIS_COLUMN, N_("Statistics"),
1754 N_("Total"), PIVOT_RC_OTHER,
1755 /* xgettext:no-c-format */
1756 N_("% of Variance"), PIVOT_RC_PERCENT,
1757 /* xgettext:no-c-format */
1758 N_("Cumulative %"), PIVOT_RC_PERCENT);
1760 struct pivot_dimension *phase = pivot_dimension_create (
1761 table, PIVOT_AXIS_COLUMN, N_("Phase"));
1762 if (factor->print & PRINT_INITIAL)
1763 pivot_category_create_leaves (phase->root, N_("Initial Eigenvalues"));
1765 if (factor->print & PRINT_EXTRACTION)
1766 pivot_category_create_leaves (phase->root,
1767 N_("Extraction Sums of Squared Loadings"));
1769 if (factor->print & PRINT_ROTATION)
1770 pivot_category_create_leaves (phase->root,
1771 N_("Rotation Sums of Squared Loadings"));
1773 struct pivot_dimension *components = pivot_dimension_create (
1774 table, PIVOT_AXIS_ROW,
1775 factor->extraction == EXTRACTION_PC ? N_("Component") : N_("Factor"));
1777 double i_total = 0.0;
1778 for (size_t i = 0 ; i < initial_eigenvalues->size; ++i)
1779 i_total += gsl_vector_get (initial_eigenvalues, i);
1781 double e_total = (factor->extraction == EXTRACTION_PAF
1788 for (size_t i = 0 ; i < factor->n_vars; ++i)
1790 const double i_lambda = gsl_vector_get (initial_eigenvalues, i);
1791 double i_percent = 100.0 * i_lambda / i_total ;
1794 const double e_lambda = gsl_vector_get (extracted_eigenvalues, i);
1795 double e_percent = 100.0 * e_lambda / e_total ;
1798 int row = pivot_category_create_leaf (
1799 components->root, pivot_value_new_integer (i + 1));
1803 /* Initial Eigenvalues */
1804 if (factor->print & PRINT_INITIAL)
1805 put_variance (table, row, phase_idx++, i_lambda, i_percent, i_cum);
1807 if (i < idata->n_extractions)
1809 if (factor->print & PRINT_EXTRACTION)
1810 put_variance (table, row, phase_idx++, e_lambda, e_percent, e_cum);
1812 if (rotated_loadings != NULL && factor->print & PRINT_ROTATION)
1814 double r_lambda = gsl_vector_get (rotated_loadings, i);
1815 double r_percent = 100.0 * r_lambda / e_total ;
1816 if (factor->rotation == ROT_PROMAX)
1817 r_lambda = r_percent = SYSMIS;
1820 put_variance (table, row, phase_idx++, r_lambda, r_percent,
1826 pivot_table_submit (table);
1830 show_factor_correlation (const struct cmd_factor * factor, const gsl_matrix *fcm)
1832 struct pivot_table *table = pivot_table_create (
1833 N_("Factor Correlation Matrix"));
1835 create_numeric_dimension (
1836 table, PIVOT_AXIS_ROW,
1837 factor->extraction == EXTRACTION_PC ? N_("Component") : N_("Factor"),
1840 create_numeric_dimension (table, PIVOT_AXIS_COLUMN, N_("Factor 2"),
1843 for (size_t i = 0 ; i < fcm->size1; ++i)
1844 for (size_t j = 0 ; j < fcm->size2; ++j)
1845 pivot_table_put2 (table, j, i,
1846 pivot_value_new_number (gsl_matrix_get (fcm, i, j)));
1848 pivot_table_submit (table);
1852 add_var_dims (struct pivot_table *table, const struct cmd_factor *factor)
1854 for (int i = 0; i < 2; i++)
1856 struct pivot_dimension *d = pivot_dimension_create (
1857 table, i ? PIVOT_AXIS_ROW : PIVOT_AXIS_COLUMN,
1860 for (size_t j = 0; j < factor->n_vars; j++)
1861 pivot_category_create_leaf (
1862 d->root, pivot_value_new_variable (factor->vars[j]));
1867 show_aic (const struct cmd_factor *factor, const struct idata *idata)
1869 if ((factor->print & PRINT_AIC) == 0)
1872 struct pivot_table *table = pivot_table_create (N_("Anti-Image Matrices"));
1874 add_var_dims (table, factor);
1876 pivot_dimension_create (table, PIVOT_AXIS_ROW, N_("Statistics"),
1877 N_("Anti-image Covariance"),
1878 N_("Anti-image Correlation"));
1880 for (size_t i = 0; i < factor->n_vars; ++i)
1881 for (size_t j = 0; j < factor->n_vars; ++j)
1883 double cov = gsl_matrix_get (idata->ai_cov, i, j);
1884 pivot_table_put3 (table, i, j, 0, pivot_value_new_number (cov));
1886 double corr = gsl_matrix_get (idata->ai_cor, i, j);
1887 pivot_table_put3 (table, i, j, 1, pivot_value_new_number (corr));
1890 pivot_table_submit (table);
1894 show_correlation_matrix (const struct cmd_factor *factor, const struct idata *idata)
1896 if (!(factor->print & (PRINT_CORRELATION | PRINT_SIG | PRINT_DETERMINANT)))
1899 struct pivot_table *table = pivot_table_create (N_("Correlation Matrix"));
1901 if (factor->print & (PRINT_CORRELATION | PRINT_SIG))
1903 add_var_dims (table, factor);
1905 struct pivot_dimension *statistics = pivot_dimension_create (
1906 table, PIVOT_AXIS_ROW, N_("Statistics"));
1907 if (factor->print & PRINT_CORRELATION)
1908 pivot_category_create_leaves (statistics->root, N_("Correlation"),
1909 PIVOT_RC_CORRELATION);
1910 if (factor->print & PRINT_SIG)
1911 pivot_category_create_leaves (statistics->root, N_("Sig. (1-tailed)"),
1912 PIVOT_RC_SIGNIFICANCE);
1915 if (factor->print & PRINT_CORRELATION)
1917 for (int i = 0; i < factor->n_vars; ++i)
1918 for (int j = 0; j < factor->n_vars; ++j)
1920 double corr = gsl_matrix_get (idata->mm.corr, i, j);
1921 pivot_table_put3 (table, j, i, stat_idx,
1922 pivot_value_new_number (corr));
1927 if (factor->print & PRINT_SIG)
1929 for (int i = 0; i < factor->n_vars; ++i)
1930 for (int j = 0; j < factor->n_vars; ++j)
1933 double rho = gsl_matrix_get (idata->mm.corr, i, j);
1934 double w = gsl_matrix_get (idata->mm.n, i, j);
1935 double sig = significance_of_correlation (rho, w);
1936 pivot_table_put3 (table, j, i, stat_idx,
1937 pivot_value_new_number (sig));
1943 if (factor->print & PRINT_DETERMINANT)
1944 table->caption = pivot_value_new_user_text_nocopy (
1945 xasprintf ("%s: %.2f", _("Determinant"), idata->detR));
1947 pivot_table_submit (table);
1951 show_covariance_matrix (const struct cmd_factor *factor, const struct idata *idata)
1953 if (!(factor->print & PRINT_COVARIANCE))
1956 struct pivot_table *table = pivot_table_create (N_("Covariance Matrix"));
1957 add_var_dims (table, factor);
1959 for (int i = 0; i < factor->n_vars; ++i)
1960 for (int j = 0; j < factor->n_vars; ++j)
1962 double cov = gsl_matrix_get (idata->mm.cov, i, j);
1963 pivot_table_put2 (table, j, i, pivot_value_new_number (cov));
1966 pivot_table_submit (table);
1971 do_factor (const struct cmd_factor *factor, struct casereader *r)
1974 struct idata *idata = idata_alloc (factor->n_vars);
1976 idata->cvm = covariance_1pass_create (factor->n_vars, factor->vars,
1977 factor->wv, factor->exclude, true);
1979 for (; (c = casereader_read (r)); case_unref (c))
1981 covariance_accumulate (idata->cvm, c);
1984 idata->mm.cov = covariance_calculate (idata->cvm);
1986 if (idata->mm.cov == NULL)
1988 msg (MW, _("The dataset contains no complete observations. No analysis will be performed."));
1989 covariance_destroy (idata->cvm);
1993 idata->mm.var_matrix = covariance_moments (idata->cvm, MOMENT_VARIANCE);
1994 idata->mm.mean_matrix = covariance_moments (idata->cvm, MOMENT_MEAN);
1995 idata->mm.n = covariance_moments (idata->cvm, MOMENT_NONE);
1997 do_factor_by_matrix (factor, idata);
2000 gsl_matrix_free (idata->mm.corr);
2001 gsl_matrix_free (idata->mm.cov);
2004 casereader_destroy (r);
2008 do_factor_by_matrix (const struct cmd_factor *factor, struct idata *idata)
2010 if (!idata->mm.cov && !idata->mm.corr)
2012 msg (ME, _("The dataset has no complete covariance or correlation matrix."));
2016 if (idata->mm.cov && !idata->mm.corr)
2017 idata->mm.corr = correlation_from_covariance (idata->mm.cov, idata->mm.var_matrix);
2018 if (idata->mm.corr && !idata->mm.cov)
2019 idata->mm.cov = covariance_from_correlation (idata->mm.corr, idata->mm.var_matrix);
2020 if (factor->method == METHOD_CORR)
2021 idata->analysis_matrix = idata->mm.corr;
2023 idata->analysis_matrix = idata->mm.cov;
2026 r_inv = clone_matrix (idata->mm.corr);
2027 gsl_linalg_cholesky_decomp (r_inv);
2028 gsl_linalg_cholesky_invert (r_inv);
2030 idata->ai_cov = anti_image_cov (r_inv);
2031 idata->ai_cor = anti_image_corr (r_inv, idata);
2034 double sum_ssq_r = 0;
2035 double sum_ssq_a = 0;
2036 for (i = 0; i < r_inv->size1; ++i)
2038 sum_ssq_r += ssq_od_n (idata->mm.corr, i);
2039 sum_ssq_a += ssq_od_n (idata->ai_cor, i);
2042 gsl_matrix_free (r_inv);
2044 if (factor->print & PRINT_DETERMINANT
2045 || factor->print & PRINT_KMO)
2049 const int size = idata->mm.corr->size1;
2050 gsl_permutation *p = gsl_permutation_calloc (size);
2051 gsl_matrix *tmp = gsl_matrix_calloc (size, size);
2052 gsl_matrix_memcpy (tmp, idata->mm.corr);
2054 gsl_linalg_LU_decomp (tmp, p, &sign);
2055 idata->detR = gsl_linalg_LU_det (tmp, sign);
2056 gsl_permutation_free (p);
2057 gsl_matrix_free (tmp);
2060 if (factor->print & PRINT_UNIVARIATE)
2062 struct pivot_table *table = pivot_table_create (
2063 N_("Descriptive Statistics"));
2064 pivot_table_set_weight_var (table, factor->wv);
2066 pivot_dimension_create (table, PIVOT_AXIS_COLUMN, N_("Statistics"),
2067 N_("Mean"), PIVOT_RC_OTHER,
2068 N_("Std. Deviation"), PIVOT_RC_OTHER,
2069 N_("Analysis N"), PIVOT_RC_COUNT);
2071 struct pivot_dimension *variables = pivot_dimension_create (
2072 table, PIVOT_AXIS_ROW, N_("Variables"));
2074 for (i = 0 ; i < factor->n_vars; ++i)
2076 const struct variable *v = factor->vars[i];
2078 int row = pivot_category_create_leaf (
2079 variables->root, pivot_value_new_variable (v));
2081 double entries[] = {
2082 gsl_matrix_get (idata->mm.mean_matrix, i, i),
2083 sqrt (gsl_matrix_get (idata->mm.var_matrix, i, i)),
2084 gsl_matrix_get (idata->mm.n, i, i),
2086 for (size_t j = 0; j < sizeof entries / sizeof *entries; j++)
2087 pivot_table_put2 (table, j, row,
2088 pivot_value_new_number (entries[j]));
2091 pivot_table_submit (table);
2094 if (factor->print & PRINT_KMO)
2096 struct pivot_table *table = pivot_table_create (
2097 N_("KMO and Bartlett's Test"));
2099 struct pivot_dimension *statistics = pivot_dimension_create (
2100 table, PIVOT_AXIS_ROW, N_("Statistics"),
2101 N_("Kaiser-Meyer-Olkin Measure of Sampling Adequacy"), PIVOT_RC_OTHER);
2102 pivot_category_create_group (
2103 statistics->root, N_("Bartlett's Test of Sphericity"),
2104 N_("Approx. Chi-Square"), PIVOT_RC_OTHER,
2105 N_("df"), PIVOT_RC_INTEGER,
2106 N_("Sig."), PIVOT_RC_SIGNIFICANCE);
2108 /* The literature doesn't say what to do for the value of W when
2109 missing values are involved. The best thing I can think of
2110 is to take the mean average. */
2112 for (i = 0; i < idata->mm.n->size1; ++i)
2113 w += gsl_matrix_get (idata->mm.n, i, i);
2114 w /= idata->mm.n->size1;
2116 double xsq = ((w - 1 - (2 * factor->n_vars + 5) / 6.0)
2117 * -log (idata->detR));
2118 double df = factor->n_vars * (factor->n_vars - 1) / 2;
2119 double entries[] = {
2120 sum_ssq_r / (sum_ssq_r + sum_ssq_a),
2123 gsl_cdf_chisq_Q (xsq, df)
2125 for (size_t i = 0; i < sizeof entries / sizeof *entries; i++)
2126 pivot_table_put1 (table, i, pivot_value_new_number (entries[i]));
2128 pivot_table_submit (table);
2131 show_correlation_matrix (factor, idata);
2132 show_covariance_matrix (factor, idata);
2134 covariance_destroy (idata->cvm);
2137 gsl_matrix *am = matrix_dup (idata->analysis_matrix);
2138 gsl_eigen_symmv_workspace *workspace = gsl_eigen_symmv_alloc (factor->n_vars);
2140 gsl_eigen_symmv (am, idata->eval, idata->evec, workspace);
2142 gsl_eigen_symmv_free (workspace);
2143 gsl_matrix_free (am);
2146 gsl_eigen_symmv_sort (idata->eval, idata->evec, GSL_EIGEN_SORT_ABS_DESC);
2148 idata->n_extractions = n_extracted_factors (factor, idata);
2150 if (idata->n_extractions == 0)
2152 msg (MW, _("The %s criteria result in zero factors extracted. Therefore no analysis will be performed."), "FACTOR");
2156 if (idata->n_extractions > factor->n_vars)
2159 _("The %s criteria result in more factors than variables, which is not meaningful. No analysis will be performed."),
2165 gsl_matrix *rotated_factors = NULL;
2166 gsl_matrix *pattern_matrix = NULL;
2167 gsl_matrix *fcm = NULL;
2168 gsl_vector *rotated_loadings = NULL;
2170 const gsl_vector *extracted_eigenvalues = NULL;
2171 gsl_vector *initial_communalities = gsl_vector_alloc (factor->n_vars);
2172 gsl_vector *extracted_communalities = gsl_vector_alloc (factor->n_vars);
2174 struct factor_matrix_workspace *fmw = factor_matrix_workspace_alloc (idata->msr->size, idata->n_extractions);
2175 gsl_matrix *factor_matrix = gsl_matrix_calloc (factor->n_vars, fmw->n_factors);
2177 if (factor->extraction == EXTRACTION_PAF)
2179 gsl_vector *diff = gsl_vector_alloc (idata->msr->size);
2180 struct smr_workspace *ws = ws_create (idata->analysis_matrix);
2182 for (i = 0 ; i < factor->n_vars ; ++i)
2184 double r2 = squared_multiple_correlation (idata->analysis_matrix, i, ws);
2186 gsl_vector_set (idata->msr, i, r2);
2190 gsl_vector_memcpy (initial_communalities, idata->msr);
2192 for (i = 0; i < factor->extraction_iterations; ++i)
2195 gsl_vector_memcpy (diff, idata->msr);
2197 iterate_factor_matrix (idata->analysis_matrix, idata->msr, factor_matrix, fmw);
2199 gsl_vector_sub (diff, idata->msr);
2201 gsl_vector_minmax (diff, &min, &max);
2203 if (fabs (min) < factor->econverge && fabs (max) < factor->econverge)
2206 gsl_vector_free (diff);
2210 gsl_vector_memcpy (extracted_communalities, idata->msr);
2211 extracted_eigenvalues = fmw->eval;
2213 else if (factor->extraction == EXTRACTION_PC)
2215 for (i = 0; i < factor->n_vars; ++i)
2216 gsl_vector_set (initial_communalities, i, communality (idata, i, factor->n_vars));
2218 gsl_vector_memcpy (extracted_communalities, initial_communalities);
2220 iterate_factor_matrix (idata->analysis_matrix, extracted_communalities, factor_matrix, fmw);
2223 extracted_eigenvalues = idata->eval;
2227 show_aic (factor, idata);
2228 show_communalities (factor, initial_communalities, extracted_communalities);
2230 if (factor->rotation != ROT_NONE)
2232 rotated_factors = gsl_matrix_calloc (factor_matrix->size1, factor_matrix->size2);
2233 rotated_loadings = gsl_vector_calloc (factor_matrix->size2);
2234 if (factor->rotation == ROT_PROMAX)
2236 pattern_matrix = gsl_matrix_calloc (factor_matrix->size1, factor_matrix->size2);
2237 fcm = gsl_matrix_calloc (factor_matrix->size2, factor_matrix->size2);
2241 rotate (factor, factor_matrix, extracted_communalities, rotated_factors, rotated_loadings, pattern_matrix, fcm);
2244 show_explained_variance (factor, idata, idata->eval, extracted_eigenvalues, rotated_loadings);
2246 factor_matrix_workspace_free (fmw);
2248 show_scree (factor, idata);
2250 show_factor_matrix (factor, idata,
2251 (factor->extraction == EXTRACTION_PC
2252 ? N_("Component Matrix") : N_("Factor Matrix")),
2255 if (factor->rotation == ROT_PROMAX)
2257 show_factor_matrix (factor, idata, N_("Pattern Matrix"),
2259 gsl_matrix_free (pattern_matrix);
2262 if (factor->rotation != ROT_NONE)
2264 show_factor_matrix (factor, idata,
2265 (factor->rotation == ROT_PROMAX
2266 ? N_("Structure Matrix")
2267 : factor->extraction == EXTRACTION_PC
2268 ? N_("Rotated Component Matrix")
2269 : N_("Rotated Factor Matrix")),
2272 gsl_matrix_free (rotated_factors);
2275 if (factor->rotation == ROT_PROMAX)
2277 show_factor_correlation (factor, fcm);
2278 gsl_matrix_free (fcm);
2281 gsl_matrix_free (factor_matrix);
2282 gsl_vector_free (rotated_loadings);
2283 gsl_vector_free (initial_communalities);
2284 gsl_vector_free (extracted_communalities);