Cope with new location of upstream files
[pspp] / src / language / stats / factor.c
1 /* PSPP - a program for statistical analysis.
2    Copyright (C) 2009, 2010, 2011, 2012, 2014, 2015, 2016 Free Software Foundation, Inc.
3
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.
8
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.
13
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/>. */
16
17 #include <config.h>
18
19 #include <gsl/gsl_vector.h>
20 #include <gsl/gsl_linalg.h>
21 #include <gsl/gsl_matrix.h>
22 #include <gsl/gsl_eigen.h>
23 #include <gsl/gsl_blas.h>
24 #include <gsl/gsl_sort_vector.h>
25 #include <gsl/gsl_cdf.h>
26
27 #include "data/casegrouper.h"
28 #include "data/casereader.h"
29 #include "data/casewriter.h"
30 #include "data/dataset.h"
31 #include "data/dictionary.h"
32 #include "data/format.h"
33 #include "data/subcase.h"
34 #include "language/command.h"
35 #include "language/lexer/lexer.h"
36 #include "language/lexer/value-parser.h"
37 #include "language/lexer/variable-parser.h"
38 #include "libpspp/cast.h"
39 #include "libpspp/message.h"
40 #include "libpspp/misc.h"
41 #include "math/correlation.h"
42 #include "math/covariance.h"
43 #include "math/moments.h"
44 #include "output/chart-item.h"
45 #include "output/charts/scree.h"
46 #include "output/tab.h"
47
48 #include "gettext.h"
49 #define _(msgid) gettext (msgid)
50 #define N_(msgid) msgid
51
52 enum method
53   {
54     METHOD_CORR,
55     METHOD_COV
56   };
57
58 enum missing_type
59   {
60     MISS_LISTWISE,
61     MISS_PAIRWISE,
62     MISS_MEANSUB,
63   };
64
65 enum extraction_method
66   {
67     EXTRACTION_PC,
68     EXTRACTION_PAF,
69   };
70
71 enum plot_opts
72   {
73     PLOT_SCREE = 0x0001,
74     PLOT_ROTATION = 0x0002
75   };
76
77 enum print_opts
78   {
79     PRINT_UNIVARIATE  = 0x0001,
80     PRINT_DETERMINANT = 0x0002,
81     PRINT_INV         = 0x0004,
82     PRINT_AIC         = 0x0008,
83     PRINT_SIG         = 0x0010,
84     PRINT_COVARIANCE  = 0x0020,
85     PRINT_CORRELATION = 0x0040,
86     PRINT_ROTATION    = 0x0080,
87     PRINT_EXTRACTION  = 0x0100,
88     PRINT_INITIAL     = 0x0200,
89     PRINT_KMO         = 0x0400,
90     PRINT_REPR        = 0x0800,
91     PRINT_FSCORE      = 0x1000
92   };
93
94 enum rotation_type
95   {
96     ROT_VARIMAX = 0,
97     ROT_EQUAMAX,
98     ROT_QUARTIMAX,
99     ROT_PROMAX,
100     ROT_NONE
101   };
102
103 typedef void (*rotation_coefficients) (double *x, double *y,
104                                     double a, double b, double c, double d,
105                                     const gsl_matrix *loadings );
106
107
108 static void
109 varimax_coefficients (double *x, double *y,
110                       double a, double b, double c, double d,
111                       const gsl_matrix *loadings )
112 {
113   *x = d - 2 * a * b / loadings->size1;
114   *y = c - (a * a - b * b) / loadings->size1;
115 }
116
117 static void
118 equamax_coefficients (double *x, double *y,
119                       double a, double b, double c, double d,
120                       const gsl_matrix *loadings )
121 {
122   *x = d - loadings->size2 * a * b / loadings->size1;
123   *y = c - loadings->size2 * (a * a - b * b) / (2 * loadings->size1);
124 }
125
126 static void
127 quartimax_coefficients (double *x, double *y,
128                       double a UNUSED, double b UNUSED, double c, double d,
129                       const gsl_matrix *loadings UNUSED)
130 {
131   *x = d ;
132   *y = c ;
133 }
134
135 static const rotation_coefficients rotation_coeff[] = {
136   varimax_coefficients,
137   equamax_coefficients,
138   quartimax_coefficients,
139   varimax_coefficients  /* PROMAX is identical to VARIMAX */
140 };
141
142
143 /* return diag (C'C) ^ {-0.5} */
144 static gsl_matrix *
145 diag_rcp_sqrt (const gsl_matrix *C)
146 {
147   int j;
148   gsl_matrix *d =  gsl_matrix_calloc (C->size1, C->size2);
149   gsl_matrix *r =  gsl_matrix_calloc (C->size1, C->size2);
150
151   assert (C->size1 == C->size2);
152
153   gsl_linalg_matmult_mod (C,  GSL_LINALG_MOD_TRANSPOSE,
154                           C,  GSL_LINALG_MOD_NONE,
155                           d);
156
157   for (j = 0 ; j < d->size2; ++j)
158     {
159       double e = gsl_matrix_get (d, j, j);
160       e = 1.0 / sqrt (e);
161       gsl_matrix_set (r, j, j, e);
162     }
163
164   gsl_matrix_free (d);
165
166   return r;
167 }
168
169
170
171 /* return diag ((C'C)^-1) ^ {-0.5} */
172 static gsl_matrix *
173 diag_rcp_inv_sqrt (const gsl_matrix *CCinv)
174 {
175   int j;
176   gsl_matrix *r =  gsl_matrix_calloc (CCinv->size1, CCinv->size2);
177
178   assert (CCinv->size1 == CCinv->size2);
179
180   for (j = 0 ; j < CCinv->size2; ++j)
181     {
182       double e = gsl_matrix_get (CCinv, j, j);
183       e = 1.0 / sqrt (e);
184       gsl_matrix_set (r, j, j, e);
185     }
186
187   return r;
188 }
189
190
191
192
193
194 struct cmd_factor
195 {
196   size_t n_vars;
197   const struct variable **vars;
198
199   const struct variable *wv;
200
201   enum method method;
202   enum missing_type missing_type;
203   enum mv_class exclude;
204   enum print_opts print;
205   enum extraction_method extraction;
206   enum plot_opts plot;
207   enum rotation_type rotation;
208   int rotation_iterations;
209   int promax_power;
210
211   /* Extraction Criteria */
212   int n_factors;
213   double min_eigen;
214   double econverge;
215   int extraction_iterations;
216
217   double rconverge;
218
219   /* Format */
220   double blank;
221   bool sort;
222 };
223
224 struct idata
225 {
226   /* Intermediate values used in calculation */
227
228   const gsl_matrix *corr ;  /* The correlation matrix */
229   gsl_matrix *cov ;         /* The covariance matrix */
230   const gsl_matrix *n ;     /* Matrix of number of samples */
231
232   gsl_vector *eval ;  /* The eigenvalues */
233   gsl_matrix *evec ;  /* The eigenvectors */
234
235   int n_extractions;
236
237   gsl_vector *msr ;  /* Multiple Squared Regressions */
238
239   double detR;  /* The determinant of the correlation matrix */
240 };
241
242 static struct idata *
243 idata_alloc (size_t n_vars)
244 {
245   struct idata *id = xzalloc (sizeof (*id));
246
247   id->n_extractions = 0;
248   id->msr = gsl_vector_alloc (n_vars);
249
250   id->eval = gsl_vector_alloc (n_vars);
251   id->evec = gsl_matrix_alloc (n_vars, n_vars);
252
253   return id;
254 }
255
256 static void
257 idata_free (struct idata *id)
258 {
259   gsl_vector_free (id->msr);
260   gsl_vector_free (id->eval);
261   gsl_matrix_free (id->evec);
262   if (id->cov != NULL)
263     gsl_matrix_free (id->cov);
264   if (id->corr != NULL)
265     gsl_matrix_free (CONST_CAST (gsl_matrix *, id->corr));
266
267   free (id);
268 }
269
270
271 static gsl_matrix *
272 anti_image (const gsl_matrix *m)
273 {
274   int i, j;
275   gsl_matrix *a;
276   assert (m->size1 == m->size2);
277
278   a = gsl_matrix_alloc (m->size1, m->size2);
279
280   for (i = 0; i < m->size1; ++i)
281     {
282       for (j = 0; j < m->size2; ++j)
283         {
284           double *p = gsl_matrix_ptr (a, i, j);
285           *p = gsl_matrix_get (m, i, j);
286           *p /= gsl_matrix_get (m, i, i);
287           *p /= gsl_matrix_get (m, j, j);
288         }
289     }
290
291   return a;
292 }
293
294
295 /* Return the sum of all the elements excluding row N */
296 static double
297 ssq_od_n (const gsl_matrix *m, int n)
298 {
299   int i, j;
300   double ss = 0;
301   assert (m->size1 == m->size2);
302
303   assert (n < m->size1);
304
305   for (i = 0; i < m->size1; ++i)
306     {
307       if (i == n ) continue;
308       for (j = 0; j < m->size2; ++j)
309         {
310           ss += pow2 (gsl_matrix_get (m, i, j));
311         }
312     }
313
314   return ss;
315 }
316
317
318
319 #if 0
320 static void
321 dump_matrix (const gsl_matrix *m)
322 {
323   size_t i, j;
324
325   for (i = 0 ; i < m->size1; ++i)
326     {
327       for (j = 0 ; j < m->size2; ++j)
328         printf ("%02f ", gsl_matrix_get (m, i, j));
329       printf ("\n");
330     }
331 }
332
333 static void
334 dump_matrix_permute (const gsl_matrix *m, const gsl_permutation *p)
335 {
336   size_t i, j;
337
338   for (i = 0 ; i < m->size1; ++i)
339     {
340       for (j = 0 ; j < m->size2; ++j)
341         printf ("%02f ", gsl_matrix_get (m, gsl_permutation_get (p, i), j));
342       printf ("\n");
343     }
344 }
345
346
347 static void
348 dump_vector (const gsl_vector *v)
349 {
350   size_t i;
351   for (i = 0 ; i < v->size; ++i)
352     {
353       printf ("%02f\n", gsl_vector_get (v, i));
354     }
355   printf ("\n");
356 }
357 #endif
358
359
360 static int
361 n_extracted_factors (const struct cmd_factor *factor, struct idata *idata)
362 {
363   int i;
364
365   /* If there is a cached value, then return that. */
366   if ( idata->n_extractions != 0)
367     return idata->n_extractions;
368
369   /* Otherwise, if the number of factors has been explicitly requested,
370      use that. */
371   if (factor->n_factors > 0)
372     {
373       idata->n_extractions = factor->n_factors;
374       goto finish;
375     }
376
377   /* Use the MIN_EIGEN setting. */
378   for (i = 0 ; i < idata->eval->size; ++i)
379     {
380       double evali = fabs (gsl_vector_get (idata->eval, i));
381
382       idata->n_extractions = i;
383
384       if (evali < factor->min_eigen)
385         goto finish;
386     }
387
388  finish:
389   return idata->n_extractions;
390 }
391
392
393 /* Returns a newly allocated matrix identical to M.
394    It it the callers responsibility to free the returned value.
395 */
396 static gsl_matrix *
397 matrix_dup (const gsl_matrix *m)
398 {
399   gsl_matrix *n =  gsl_matrix_alloc (m->size1, m->size2);
400
401   gsl_matrix_memcpy (n, m);
402
403   return n;
404 }
405
406
407 struct smr_workspace
408 {
409   /* Copy of the subject */
410   gsl_matrix *m;
411
412   gsl_matrix *inverse;
413
414   gsl_permutation *perm;
415
416   gsl_matrix *result1;
417   gsl_matrix *result2;
418 };
419
420
421 static struct smr_workspace *ws_create (const gsl_matrix *input)
422 {
423   struct smr_workspace *ws = xmalloc (sizeof (*ws));
424
425   ws->m = gsl_matrix_alloc (input->size1, input->size2);
426   ws->inverse = gsl_matrix_calloc (input->size1 - 1, input->size2 - 1);
427   ws->perm = gsl_permutation_alloc (input->size1 - 1);
428   ws->result1 = gsl_matrix_calloc (input->size1 - 1, 1);
429   ws->result2 = gsl_matrix_calloc (1, 1);
430
431   return ws;
432 }
433
434 static void
435 ws_destroy (struct smr_workspace *ws)
436 {
437   gsl_matrix_free (ws->result2);
438   gsl_matrix_free (ws->result1);
439   gsl_permutation_free (ws->perm);
440   gsl_matrix_free (ws->inverse);
441   gsl_matrix_free (ws->m);
442
443   free (ws);
444 }
445
446
447 /*
448    Return the square of the regression coefficient for VAR regressed against all other variables.
449  */
450 static double
451 squared_multiple_correlation (const gsl_matrix *corr, int var, struct smr_workspace *ws)
452 {
453   /* For an explanation of what this is doing, see
454      http://www.visualstatistics.net/Visual%20Statistics%20Multimedia/multiple_regression_analysis.htm
455   */
456
457   int signum = 0;
458   gsl_matrix_view rxx;
459
460   gsl_matrix_memcpy (ws->m, corr);
461
462   gsl_matrix_swap_rows (ws->m, 0, var);
463   gsl_matrix_swap_columns (ws->m, 0, var);
464
465   rxx = gsl_matrix_submatrix (ws->m, 1, 1, ws->m->size1 - 1, ws->m->size1 - 1);
466
467   gsl_linalg_LU_decomp (&rxx.matrix, ws->perm, &signum);
468
469   gsl_linalg_LU_invert (&rxx.matrix, ws->perm, ws->inverse);
470
471   {
472     gsl_matrix_const_view rxy = gsl_matrix_const_submatrix (ws->m, 1, 0, ws->m->size1 - 1, 1);
473     gsl_matrix_const_view ryx = gsl_matrix_const_submatrix (ws->m, 0, 1, 1, ws->m->size1 - 1);
474
475     gsl_blas_dgemm (CblasNoTrans,  CblasNoTrans,
476                     1.0, ws->inverse, &rxy.matrix, 0.0, ws->result1);
477
478     gsl_blas_dgemm (CblasNoTrans,  CblasNoTrans,
479                     1.0, &ryx.matrix, ws->result1, 0.0, ws->result2);
480   }
481
482   return gsl_matrix_get (ws->result2, 0, 0);
483 }
484
485
486
487 static double the_communality (const gsl_matrix *evec, const gsl_vector *eval, int n, int n_factors);
488
489
490 struct factor_matrix_workspace
491 {
492   size_t n_factors;
493   gsl_eigen_symmv_workspace *eigen_ws;
494
495   gsl_vector *eval ;
496   gsl_matrix *evec ;
497
498   gsl_matrix *gamma ;
499
500   gsl_matrix *r;
501 };
502
503 static struct factor_matrix_workspace *
504 factor_matrix_workspace_alloc (size_t n, size_t nf)
505 {
506   struct factor_matrix_workspace *ws = xmalloc (sizeof (*ws));
507
508   ws->n_factors = nf;
509   ws->gamma = gsl_matrix_calloc (nf, nf);
510   ws->eigen_ws = gsl_eigen_symmv_alloc (n);
511   ws->eval = gsl_vector_alloc (n);
512   ws->evec = gsl_matrix_alloc (n, n);
513   ws->r  = gsl_matrix_alloc (n, n);
514
515   return ws;
516 }
517
518 static void
519 factor_matrix_workspace_free (struct factor_matrix_workspace *ws)
520 {
521   gsl_eigen_symmv_free (ws->eigen_ws);
522   gsl_vector_free (ws->eval);
523   gsl_matrix_free (ws->evec);
524   gsl_matrix_free (ws->gamma);
525   gsl_matrix_free (ws->r);
526   free (ws);
527 }
528
529 /*
530   Shift P left by OFFSET places, and overwrite TARGET
531   with the shifted result.
532   Positions in TARGET less than OFFSET are unchanged.
533 */
534 static void
535 perm_shift_apply (gsl_permutation *target, const gsl_permutation *p,
536                   size_t offset)
537 {
538   size_t i;
539   assert (target->size == p->size);
540   assert (offset <= target->size);
541
542   for (i = 0; i < target->size - offset; ++i)
543     {
544       target->data[i] = p->data [i + offset];
545     }
546 }
547
548
549 /*
550    Indirectly sort the rows of matrix INPUT, storing the sort order in PERM.
551    The sort criteria are as follows:
552
553    Rows are sorted on the first column, until the absolute value of an
554    element in a subsequent column  is greater than that of the first
555    column.  Thereafter, rows will be sorted on the second column,
556    until the absolute value of an element in a subsequent column
557    exceeds that of the second column ...
558 */
559 static void
560 sort_matrix_indirect (const gsl_matrix *input, gsl_permutation *perm)
561 {
562   const size_t n = perm->size;
563   const size_t m = input->size2;
564   int i, j;
565   gsl_matrix *mat ;
566   int column_n = 0;
567   int row_n = 0;
568   gsl_permutation *p;
569
570   assert (perm->size == input->size1);
571
572   p = gsl_permutation_alloc (n);
573
574   /* Copy INPUT into MAT, discarding the sign */
575   mat = gsl_matrix_alloc (n, m);
576   for (i = 0 ; i < mat->size1; ++i)
577     {
578       for (j = 0 ; j < mat->size2; ++j)
579         {
580           double x = gsl_matrix_get (input, i, j);
581           gsl_matrix_set (mat, i, j, fabs (x));
582         }
583     }
584
585   while (column_n < m && row_n < n)
586     {
587       gsl_vector_const_view columni = gsl_matrix_const_column (mat, column_n);
588       gsl_sort_vector_index (p, &columni.vector);
589
590       for (i = 0 ; i < n; ++i)
591         {
592           gsl_vector_view row = gsl_matrix_row (mat, p->data[n - 1 - i]);
593           size_t maxindex = gsl_vector_max_index (&row.vector);
594
595           if ( maxindex > column_n )
596             break;
597
598           /* All subsequent elements of this row, are of no interest.
599              So set them all to a highly negative value */
600           for (j = column_n + 1; j < row.vector.size ; ++j)
601             gsl_vector_set (&row.vector, j, -DBL_MAX);
602         }
603
604       perm_shift_apply (perm, p, row_n);
605       row_n += i;
606
607       column_n++;
608     }
609
610   gsl_permutation_free (p);
611   gsl_matrix_free (mat);
612
613   assert ( 0 == gsl_permutation_valid (perm));
614
615   /* We want the biggest value to be first */
616   gsl_permutation_reverse (perm);
617 }
618
619
620 static void
621 drot_go (double phi, double *l0, double *l1)
622 {
623   double r0 = cos (phi) * *l0 + sin (phi) * *l1;
624   double r1 = - sin (phi) * *l0 + cos (phi) * *l1;
625
626   *l0 = r0;
627   *l1 = r1;
628 }
629
630
631 static gsl_matrix *
632 clone_matrix (const gsl_matrix *m)
633 {
634   int j, k;
635   gsl_matrix *c = gsl_matrix_calloc (m->size1, m->size2);
636
637   for (j = 0 ; j < c->size1; ++j)
638     {
639       for (k = 0 ; k < c->size2; ++k)
640         {
641           const double *v = gsl_matrix_const_ptr (m, j, k);
642           gsl_matrix_set (c, j, k, *v);
643         }
644     }
645
646   return c;
647 }
648
649
650 static double
651 initial_sv (const gsl_matrix *fm)
652 {
653   int j, k;
654
655   double sv = 0.0;
656   for (j = 0 ; j < fm->size2; ++j)
657     {
658       double l4s = 0;
659       double l2s = 0;
660
661       for (k = j + 1 ; k < fm->size2; ++k)
662         {
663           double lambda = gsl_matrix_get (fm, k, j);
664           double lambda_sq = lambda * lambda;
665           double lambda_4 = lambda_sq * lambda_sq;
666
667           l4s += lambda_4;
668           l2s += lambda_sq;
669         }
670       sv += ( fm->size1 * l4s - (l2s * l2s) ) / (fm->size1 * fm->size1 );
671     }
672   return sv;
673 }
674
675 static void
676 rotate (const struct cmd_factor *cf, const gsl_matrix *unrot,
677         const gsl_vector *communalities,
678         gsl_matrix *result,
679         gsl_vector *rotated_loadings,
680         gsl_matrix *pattern_matrix,
681         gsl_matrix *factor_correlation_matrix
682         )
683 {
684   int j, k;
685   int i;
686   double prev_sv;
687
688   /* First get a normalised version of UNROT */
689   gsl_matrix *normalised = gsl_matrix_calloc (unrot->size1, unrot->size2);
690   gsl_matrix *h_sqrt = gsl_matrix_calloc (communalities->size, communalities->size);
691   gsl_matrix *h_sqrt_inv ;
692
693   /* H is the diagonal matrix containing the absolute values of the communalities */
694   for (i = 0 ; i < communalities->size ; ++i)
695     {
696       double *ptr = gsl_matrix_ptr (h_sqrt, i, i);
697       *ptr = fabs (gsl_vector_get (communalities, i));
698     }
699
700   /* Take the square root of the communalities */
701   gsl_linalg_cholesky_decomp (h_sqrt);
702
703
704   /* Save a copy of h_sqrt and invert it */
705   h_sqrt_inv = clone_matrix (h_sqrt);
706   gsl_linalg_cholesky_decomp (h_sqrt_inv);
707   gsl_linalg_cholesky_invert (h_sqrt_inv);
708
709   /* normalised vertion is H^{1/2} x UNROT */
710   gsl_blas_dgemm (CblasNoTrans,  CblasNoTrans, 1.0, h_sqrt_inv, unrot, 0.0, normalised);
711
712   gsl_matrix_free (h_sqrt_inv);
713
714
715   /* Now perform the rotation iterations */
716
717   prev_sv = initial_sv (normalised);
718   for (i = 0 ; i < cf->rotation_iterations ; ++i)
719     {
720       double sv = 0.0;
721       for (j = 0 ; j < normalised->size2; ++j)
722         {
723           /* These variables relate to the convergence criterium */
724           double l4s = 0;
725           double l2s = 0;
726
727           for (k = j + 1 ; k < normalised->size2; ++k)
728             {
729               int p;
730               double a = 0.0;
731               double b = 0.0;
732               double c = 0.0;
733               double d = 0.0;
734               double x, y;
735               double phi;
736
737               for (p = 0; p < normalised->size1; ++p)
738                 {
739                   double jv = gsl_matrix_get (normalised, p, j);
740                   double kv = gsl_matrix_get (normalised, p, k);
741
742                   double u = jv * jv - kv * kv;
743                   double v = 2 * jv * kv;
744                   a += u;
745                   b += v;
746                   c +=  u * u - v * v;
747                   d += 2 * u * v;
748                 }
749
750               rotation_coeff [cf->rotation] (&x, &y, a, b, c, d, normalised);
751
752               phi = atan2 (x,  y) / 4.0 ;
753
754               /* Don't bother rotating if the angle is small */
755               if ( fabs (sin (phi) ) <= pow (10.0, -15.0))
756                   continue;
757
758               for (p = 0; p < normalised->size1; ++p)
759                 {
760                   double *lambda0 = gsl_matrix_ptr (normalised, p, j);
761                   double *lambda1 = gsl_matrix_ptr (normalised, p, k);
762                   drot_go (phi, lambda0, lambda1);
763                 }
764
765               /* Calculate the convergence criterium */
766               {
767                 double lambda = gsl_matrix_get (normalised, k, j);
768                 double lambda_sq = lambda * lambda;
769                 double lambda_4 = lambda_sq * lambda_sq;
770
771                 l4s += lambda_4;
772                 l2s += lambda_sq;
773               }
774             }
775           sv += ( normalised->size1 * l4s - (l2s * l2s) ) / (normalised->size1 * normalised->size1 );
776         }
777
778       if ( fabs (sv - prev_sv) <= cf->rconverge)
779         break;
780
781       prev_sv = sv;
782     }
783
784   gsl_blas_dgemm (CblasNoTrans,  CblasNoTrans, 1.0,
785                   h_sqrt, normalised,  0.0,   result);
786
787   gsl_matrix_free (h_sqrt);
788   gsl_matrix_free (normalised);
789
790   if (cf->rotation == ROT_PROMAX)
791     {
792       /* general purpose m by m matrix, where m is the number of factors */
793       gsl_matrix *mm1 =  gsl_matrix_calloc (unrot->size2, unrot->size2);
794       gsl_matrix *mm2 =  gsl_matrix_calloc (unrot->size2, unrot->size2);
795
796       /* general purpose m by p matrix, where p is the number of variables */
797       gsl_matrix *mp1 =  gsl_matrix_calloc (unrot->size2, unrot->size1);
798
799       gsl_matrix *pm1 =  gsl_matrix_calloc (unrot->size1, unrot->size2);
800
801       gsl_permutation *perm = gsl_permutation_alloc (unrot->size2);
802
803       int signum;
804
805       int i, j;
806
807       /* The following variables follow the notation by SPSS Statistical Algorithms
808          page 342 */
809       gsl_matrix *L =  gsl_matrix_calloc (unrot->size2, unrot->size2);
810       gsl_matrix *P = clone_matrix (result);
811       gsl_matrix *D ;
812       gsl_matrix *Q ;
813
814
815       /* Vector of length p containing (indexed by i)
816          \Sum^m_j {\lambda^2_{ij}} */
817       gsl_vector *rssq = gsl_vector_calloc (unrot->size1);
818
819       for (i = 0; i < P->size1; ++i)
820         {
821           double sum = 0;
822           for (j = 0; j < P->size2; ++j)
823             {
824               sum += gsl_matrix_get (result, i, j)
825                 * gsl_matrix_get (result, i, j);
826
827             }
828
829           gsl_vector_set (rssq, i, sqrt (sum));
830         }
831
832       for (i = 0; i < P->size1; ++i)
833         {
834           for (j = 0; j < P->size2; ++j)
835             {
836               double l = gsl_matrix_get (result, i, j);
837               double r = gsl_vector_get (rssq, i);
838               gsl_matrix_set (P, i, j, pow (fabs (l / r), cf->promax_power + 1) * r / l);
839             }
840         }
841
842       gsl_vector_free (rssq);
843
844       gsl_linalg_matmult_mod (result,
845                               GSL_LINALG_MOD_TRANSPOSE,
846                               result,
847                               GSL_LINALG_MOD_NONE,
848                               mm1);
849
850       gsl_linalg_LU_decomp (mm1, perm, &signum);
851       gsl_linalg_LU_invert (mm1, perm, mm2);
852
853       gsl_linalg_matmult_mod (mm2,   GSL_LINALG_MOD_NONE,
854                               result,  GSL_LINALG_MOD_TRANSPOSE,
855                               mp1);
856
857       gsl_linalg_matmult_mod (mp1, GSL_LINALG_MOD_NONE,
858                               P,   GSL_LINALG_MOD_NONE,
859                               L);
860
861       D = diag_rcp_sqrt (L);
862       Q = gsl_matrix_calloc (unrot->size2, unrot->size2);
863
864       gsl_linalg_matmult_mod (L, GSL_LINALG_MOD_NONE,
865                               D, GSL_LINALG_MOD_NONE,
866                               Q);
867
868       gsl_matrix *QQinv = gsl_matrix_calloc (unrot->size2, unrot->size2);
869
870       gsl_linalg_matmult_mod (Q, GSL_LINALG_MOD_TRANSPOSE,
871                               Q,  GSL_LINALG_MOD_NONE,
872                               QQinv);
873
874       gsl_linalg_cholesky_decomp (QQinv);
875       gsl_linalg_cholesky_invert (QQinv);
876
877
878       gsl_matrix *C = diag_rcp_inv_sqrt (QQinv);
879       gsl_matrix *Cinv =  clone_matrix (C);
880
881       gsl_linalg_cholesky_decomp (Cinv);
882       gsl_linalg_cholesky_invert (Cinv);
883
884
885       gsl_linalg_matmult_mod (result, GSL_LINALG_MOD_NONE,
886                               Q,      GSL_LINALG_MOD_NONE,
887                               pm1);
888
889       gsl_linalg_matmult_mod (pm1,      GSL_LINALG_MOD_NONE,
890                               Cinv,         GSL_LINALG_MOD_NONE,
891                               pattern_matrix);
892
893
894       gsl_linalg_matmult_mod (C,      GSL_LINALG_MOD_NONE,
895                               QQinv,  GSL_LINALG_MOD_NONE,
896                               mm1);
897
898       gsl_linalg_matmult_mod (mm1,      GSL_LINALG_MOD_NONE,
899                               C,  GSL_LINALG_MOD_TRANSPOSE,
900                               factor_correlation_matrix);
901
902       gsl_linalg_matmult_mod (pattern_matrix,      GSL_LINALG_MOD_NONE,
903                               factor_correlation_matrix,  GSL_LINALG_MOD_NONE,
904                               pm1);
905
906       gsl_matrix_memcpy (result, pm1);
907
908
909       gsl_matrix_free (QQinv);
910       gsl_matrix_free (C);
911       gsl_matrix_free (Cinv);
912
913       gsl_matrix_free (D);
914       gsl_matrix_free (Q);
915       gsl_matrix_free (L);
916       gsl_matrix_free (P);
917
918       gsl_permutation_free (perm);
919
920       gsl_matrix_free (mm1);
921       gsl_matrix_free (mm2);
922       gsl_matrix_free (mp1);
923       gsl_matrix_free (pm1);
924     }
925
926
927   /* reflect negative sums and populate the rotated loadings vector*/
928   for (i = 0 ; i < result->size2; ++i)
929     {
930       double ssq = 0.0;
931       double sum = 0.0;
932       for (j = 0 ; j < result->size1; ++j)
933         {
934           double s = gsl_matrix_get (result, j, i);
935           ssq += s * s;
936           sum += s;
937         }
938
939       gsl_vector_set (rotated_loadings, i, ssq);
940
941       if ( sum < 0 )
942         for (j = 0 ; j < result->size1; ++j)
943           {
944             double *lambda = gsl_matrix_ptr (result, j, i);
945             *lambda = - *lambda;
946           }
947     }
948 }
949
950
951 /*
952   Get an approximation for the factor matrix into FACTORS, and the communalities into COMMUNALITIES.
953   R is the matrix to be analysed.
954   WS is a pointer to a structure which must have been initialised with factor_matrix_workspace_init.
955  */
956 static void
957 iterate_factor_matrix (const gsl_matrix *r, gsl_vector *communalities, gsl_matrix *factors,
958                        struct factor_matrix_workspace *ws)
959 {
960   size_t i;
961   gsl_matrix_view mv ;
962
963   assert (r->size1 == r->size2);
964   assert (r->size1 == communalities->size);
965
966   assert (factors->size1 == r->size1);
967   assert (factors->size2 == ws->n_factors);
968
969   gsl_matrix_memcpy (ws->r, r);
970
971   /* Apply Communalities to diagonal of correlation matrix */
972   for (i = 0 ; i < communalities->size ; ++i)
973     {
974       double *x = gsl_matrix_ptr (ws->r, i, i);
975       *x = gsl_vector_get (communalities, i);
976     }
977
978   gsl_eigen_symmv (ws->r, ws->eval, ws->evec, ws->eigen_ws);
979
980   mv = gsl_matrix_submatrix (ws->evec, 0, 0, ws->evec->size1, ws->n_factors);
981
982   /* Gamma is the diagonal matrix containing the absolute values of the eigenvalues */
983   for (i = 0 ; i < ws->n_factors ; ++i)
984     {
985       double *ptr = gsl_matrix_ptr (ws->gamma, i, i);
986       *ptr = fabs (gsl_vector_get (ws->eval, i));
987     }
988
989   /* Take the square root of gamma */
990   gsl_linalg_cholesky_decomp (ws->gamma);
991
992   gsl_blas_dgemm (CblasNoTrans,  CblasNoTrans, 1.0, &mv.matrix, ws->gamma, 0.0, factors);
993
994   for (i = 0 ; i < r->size1 ; ++i)
995     {
996       double h = the_communality (ws->evec, ws->eval, i, ws->n_factors);
997       gsl_vector_set (communalities, i, h);
998     }
999 }
1000
1001
1002
1003 static bool run_factor (struct dataset *ds, const struct cmd_factor *factor);
1004
1005
1006 int
1007 cmd_factor (struct lexer *lexer, struct dataset *ds)
1008 {
1009   const struct dictionary *dict = dataset_dict (ds);
1010   int n_iterations = 25;
1011   struct cmd_factor factor;
1012   factor.n_vars = 0;
1013   factor.vars = NULL;
1014   factor.method = METHOD_CORR;
1015   factor.missing_type = MISS_LISTWISE;
1016   factor.exclude = MV_ANY;
1017   factor.print = PRINT_INITIAL | PRINT_EXTRACTION | PRINT_ROTATION;
1018   factor.extraction = EXTRACTION_PC;
1019   factor.n_factors = 0;
1020   factor.min_eigen = SYSMIS;
1021   factor.extraction_iterations = 25;
1022   factor.rotation_iterations = 25;
1023   factor.econverge = 0.001;
1024
1025   factor.blank = 0;
1026   factor.sort = false;
1027   factor.plot = 0;
1028   factor.rotation = ROT_VARIMAX;
1029
1030   factor.rconverge = 0.0001;
1031
1032   factor.wv = dict_get_weight (dict);
1033
1034   lex_match (lexer, T_SLASH);
1035
1036   if (!lex_force_match_id (lexer, "VARIABLES"))
1037     {
1038       goto error;
1039     }
1040
1041   lex_match (lexer, T_EQUALS);
1042
1043   if (!parse_variables_const (lexer, dict, &factor.vars, &factor.n_vars,
1044                               PV_NO_DUPLICATE | PV_NUMERIC))
1045     goto error;
1046
1047   if (factor.n_vars < 2)
1048     msg (MW, _("Factor analysis on a single variable is not useful."));
1049
1050   while (lex_token (lexer) != T_ENDCMD)
1051     {
1052       lex_match (lexer, T_SLASH);
1053
1054       if (lex_match_id (lexer, "ANALYSIS"))
1055         {
1056           struct const_var_set *vs;
1057           const struct variable **vars;
1058           size_t n_vars;
1059           bool ok;
1060
1061           lex_match (lexer, T_EQUALS);
1062
1063           vs = const_var_set_create_from_array (factor.vars, factor.n_vars);
1064           ok = parse_const_var_set_vars (lexer, vs, &vars, &n_vars,
1065                                          PV_NO_DUPLICATE | PV_NUMERIC);
1066           const_var_set_destroy (vs);
1067
1068           if (!ok)
1069             goto error;
1070
1071           free (factor.vars);
1072           factor.vars = vars;
1073           factor.n_vars = n_vars;
1074         }
1075       else if (lex_match_id (lexer, "PLOT"))
1076         {
1077           lex_match (lexer, T_EQUALS);
1078           while (lex_token (lexer) != T_ENDCMD && lex_token (lexer) != T_SLASH)
1079             {
1080               if (lex_match_id (lexer, "EIGEN"))
1081                 {
1082                   factor.plot |= PLOT_SCREE;
1083                 }
1084 #if FACTOR_FULLY_IMPLEMENTED
1085               else if (lex_match_id (lexer, "ROTATION"))
1086                 {
1087                 }
1088 #endif
1089               else
1090                 {
1091                   lex_error (lexer, NULL);
1092                   goto error;
1093                 }
1094             }
1095         }
1096       else if (lex_match_id (lexer, "METHOD"))
1097         {
1098           lex_match (lexer, T_EQUALS);
1099           while (lex_token (lexer) != T_ENDCMD && lex_token (lexer) != T_SLASH)
1100             {
1101               if (lex_match_id (lexer, "COVARIANCE"))
1102                 {
1103                   factor.method = METHOD_COV;
1104                 }
1105               else if (lex_match_id (lexer, "CORRELATION"))
1106                 {
1107                   factor.method = METHOD_CORR;
1108                 }
1109               else
1110                 {
1111                   lex_error (lexer, NULL);
1112                   goto error;
1113                 }
1114             }
1115         }
1116       else if (lex_match_id (lexer, "ROTATION"))
1117         {
1118           lex_match (lexer, T_EQUALS);
1119           while (lex_token (lexer) != T_ENDCMD && lex_token (lexer) != T_SLASH)
1120             {
1121               /* VARIMAX and DEFAULT are defaults */
1122               if (lex_match_id (lexer, "VARIMAX") || lex_match_id (lexer, "DEFAULT"))
1123                 {
1124                   factor.rotation = ROT_VARIMAX;
1125                 }
1126               else if (lex_match_id (lexer, "EQUAMAX"))
1127                 {
1128                   factor.rotation = ROT_EQUAMAX;
1129                 }
1130               else if (lex_match_id (lexer, "QUARTIMAX"))
1131                 {
1132                   factor.rotation = ROT_QUARTIMAX;
1133                 }
1134               else if (lex_match_id (lexer, "PROMAX"))
1135                 {
1136                   factor.promax_power = 5;
1137                   if (lex_match (lexer, T_LPAREN)
1138                       && lex_force_int (lexer))
1139                     {
1140                       factor.promax_power = lex_integer (lexer);
1141                       lex_get (lexer);
1142                       if (! lex_force_match (lexer, T_RPAREN))
1143                         goto error;
1144                     }
1145                   factor.rotation = ROT_PROMAX;
1146                 }
1147               else if (lex_match_id (lexer, "NOROTATE"))
1148                 {
1149                   factor.rotation = ROT_NONE;
1150                 }
1151               else
1152                 {
1153                   lex_error (lexer, NULL);
1154                   goto error;
1155                 }
1156             }
1157           factor.rotation_iterations = n_iterations;
1158         }
1159       else if (lex_match_id (lexer, "CRITERIA"))
1160         {
1161           lex_match (lexer, T_EQUALS);
1162           while (lex_token (lexer) != T_ENDCMD && lex_token (lexer) != T_SLASH)
1163             {
1164               if (lex_match_id (lexer, "FACTORS"))
1165                 {
1166                   if ( lex_force_match (lexer, T_LPAREN)
1167                        && lex_force_int (lexer))
1168                     {
1169                       factor.n_factors = lex_integer (lexer);
1170                       lex_get (lexer);
1171                       if (! lex_force_match (lexer, T_RPAREN))
1172                         goto error;
1173                     }
1174                 }
1175               else if (lex_match_id (lexer, "MINEIGEN"))
1176                 {
1177                   if ( lex_force_match (lexer, T_LPAREN)
1178                        && lex_force_num (lexer))
1179                     {
1180                       factor.min_eigen = lex_number (lexer);
1181                       lex_get (lexer);
1182                       if (! lex_force_match (lexer, T_RPAREN))
1183                         goto error;
1184                     }
1185                 }
1186               else if (lex_match_id (lexer, "ECONVERGE"))
1187                 {
1188                   if ( lex_force_match (lexer, T_LPAREN)
1189                        && lex_force_num (lexer))
1190                     {
1191                       factor.econverge = lex_number (lexer);
1192                       lex_get (lexer);
1193                       if (! lex_force_match (lexer, T_RPAREN))
1194                         goto error;
1195                     }
1196                 }
1197               else if (lex_match_id (lexer, "RCONVERGE"))
1198                 {
1199                   if (lex_force_match (lexer, T_LPAREN)
1200                       && lex_force_num (lexer))
1201                     {
1202                       factor.rconverge = lex_number (lexer);
1203                       lex_get (lexer);
1204                       if (! lex_force_match (lexer, T_RPAREN))
1205                         goto error;
1206                     }
1207                 }
1208               else if (lex_match_id (lexer, "ITERATE"))
1209                 {
1210                   if ( lex_force_match (lexer, T_LPAREN)
1211                        && lex_force_int (lexer))
1212                     {
1213                       n_iterations = lex_integer (lexer);
1214                       lex_get (lexer);
1215                       if (! lex_force_match (lexer, T_RPAREN))
1216                         goto error;
1217                     }
1218                 }
1219               else if (lex_match_id (lexer, "DEFAULT"))
1220                 {
1221                   factor.n_factors = 0;
1222                   factor.min_eigen = 1;
1223                   n_iterations = 25;
1224                 }
1225               else
1226                 {
1227                   lex_error (lexer, NULL);
1228                   goto error;
1229                 }
1230             }
1231         }
1232       else if (lex_match_id (lexer, "EXTRACTION"))
1233         {
1234           lex_match (lexer, T_EQUALS);
1235           while (lex_token (lexer) != T_ENDCMD && lex_token (lexer) != T_SLASH)
1236             {
1237               if (lex_match_id (lexer, "PAF"))
1238                 {
1239                   factor.extraction = EXTRACTION_PAF;
1240                 }
1241               else if (lex_match_id (lexer, "PC"))
1242                 {
1243                   factor.extraction = EXTRACTION_PC;
1244                 }
1245               else if (lex_match_id (lexer, "PA1"))
1246                 {
1247                   factor.extraction = EXTRACTION_PC;
1248                 }
1249               else if (lex_match_id (lexer, "DEFAULT"))
1250                 {
1251                   factor.extraction = EXTRACTION_PC;
1252                 }
1253               else
1254                 {
1255                   lex_error (lexer, NULL);
1256                   goto error;
1257                 }
1258             }
1259           factor.extraction_iterations = n_iterations;
1260         }
1261       else if (lex_match_id (lexer, "FORMAT"))
1262         {
1263           lex_match (lexer, T_EQUALS);
1264           while (lex_token (lexer) != T_ENDCMD && lex_token (lexer) != T_SLASH)
1265             {
1266               if (lex_match_id (lexer, "SORT"))
1267                 {
1268                   factor.sort = true;
1269                 }
1270               else if (lex_match_id (lexer, "BLANK"))
1271                 {
1272                   if ( lex_force_match (lexer, T_LPAREN)
1273                        && lex_force_num (lexer))
1274                     {
1275                       factor.blank = lex_number (lexer);
1276                       lex_get (lexer);
1277                       if (! lex_force_match (lexer, T_RPAREN))
1278                         goto error;
1279                     }
1280                 }
1281               else if (lex_match_id (lexer, "DEFAULT"))
1282                 {
1283                   factor.blank = 0;
1284                   factor.sort = false;
1285                 }
1286               else
1287                 {
1288                   lex_error (lexer, NULL);
1289                   goto error;
1290                 }
1291             }
1292         }
1293       else if (lex_match_id (lexer, "PRINT"))
1294         {
1295           factor.print = 0;
1296           lex_match (lexer, T_EQUALS);
1297           while (lex_token (lexer) != T_ENDCMD && lex_token (lexer) != T_SLASH)
1298             {
1299               if (lex_match_id (lexer, "UNIVARIATE"))
1300                 {
1301                   factor.print |= PRINT_UNIVARIATE;
1302                 }
1303               else if (lex_match_id (lexer, "DET"))
1304                 {
1305                   factor.print |= PRINT_DETERMINANT;
1306                 }
1307 #if FACTOR_FULLY_IMPLEMENTED
1308               else if (lex_match_id (lexer, "INV"))
1309                 {
1310                 }
1311               else if (lex_match_id (lexer, "AIC"))
1312                 {
1313                 }
1314 #endif
1315               else if (lex_match_id (lexer, "SIG"))
1316                 {
1317                   factor.print |= PRINT_SIG;
1318                 }
1319               else if (lex_match_id (lexer, "CORRELATION"))
1320                 {
1321                   factor.print |= PRINT_CORRELATION;
1322                 }
1323 #if FACTOR_FULLY_IMPLEMENTED
1324               else if (lex_match_id (lexer, "COVARIANCE"))
1325                 {
1326                 }
1327 #endif
1328               else if (lex_match_id (lexer, "ROTATION"))
1329                 {
1330                   factor.print |= PRINT_ROTATION;
1331                 }
1332               else if (lex_match_id (lexer, "EXTRACTION"))
1333                 {
1334                   factor.print |= PRINT_EXTRACTION;
1335                 }
1336               else if (lex_match_id (lexer, "INITIAL"))
1337                 {
1338                   factor.print |= PRINT_INITIAL;
1339                 }
1340               else if (lex_match_id (lexer, "KMO"))
1341                 {
1342                   factor.print |= PRINT_KMO;
1343                 }
1344 #if FACTOR_FULLY_IMPLEMENTED
1345               else if (lex_match_id (lexer, "REPR"))
1346                 {
1347                 }
1348               else if (lex_match_id (lexer, "FSCORE"))
1349                 {
1350                 }
1351 #endif
1352               else if (lex_match (lexer, T_ALL))
1353                 {
1354                   factor.print = 0xFFFF;
1355                 }
1356               else if (lex_match_id (lexer, "DEFAULT"))
1357                 {
1358                   factor.print |= PRINT_INITIAL ;
1359                   factor.print |= PRINT_EXTRACTION ;
1360                   factor.print |= PRINT_ROTATION ;
1361                 }
1362               else
1363                 {
1364                   lex_error (lexer, NULL);
1365                   goto error;
1366                 }
1367             }
1368         }
1369       else if (lex_match_id (lexer, "MISSING"))
1370         {
1371           lex_match (lexer, T_EQUALS);
1372           while (lex_token (lexer) != T_ENDCMD && lex_token (lexer) != T_SLASH)
1373             {
1374               if (lex_match_id (lexer, "INCLUDE"))
1375                 {
1376                   factor.exclude = MV_SYSTEM;
1377                 }
1378               else if (lex_match_id (lexer, "EXCLUDE"))
1379                 {
1380                   factor.exclude = MV_ANY;
1381                 }
1382               else if (lex_match_id (lexer, "LISTWISE"))
1383                 {
1384                   factor.missing_type = MISS_LISTWISE;
1385                 }
1386               else if (lex_match_id (lexer, "PAIRWISE"))
1387                 {
1388                   factor.missing_type = MISS_PAIRWISE;
1389                 }
1390               else if (lex_match_id (lexer, "MEANSUB"))
1391                 {
1392                   factor.missing_type = MISS_MEANSUB;
1393                 }
1394               else
1395                 {
1396                   lex_error (lexer, NULL);
1397                   goto error;
1398                 }
1399             }
1400         }
1401       else
1402         {
1403           lex_error (lexer, NULL);
1404           goto error;
1405         }
1406     }
1407
1408   if ( factor.rotation == ROT_NONE )
1409     factor.print &= ~PRINT_ROTATION;
1410
1411   if ( ! run_factor (ds, &factor))
1412     goto error;
1413
1414   free (factor.vars);
1415   return CMD_SUCCESS;
1416
1417  error:
1418   free (factor.vars);
1419   return CMD_FAILURE;
1420 }
1421
1422 static void do_factor (const struct cmd_factor *factor, struct casereader *group);
1423
1424
1425 static bool
1426 run_factor (struct dataset *ds, const struct cmd_factor *factor)
1427 {
1428   struct dictionary *dict = dataset_dict (ds);
1429   bool ok;
1430   struct casereader *group;
1431
1432   struct casegrouper *grouper = casegrouper_create_splits (proc_open (ds), dict);
1433
1434   while (casegrouper_get_next_group (grouper, &group))
1435     {
1436       if ( factor->missing_type == MISS_LISTWISE )
1437         group  = casereader_create_filter_missing (group, factor->vars, factor->n_vars,
1438                                                    factor->exclude,
1439                                                    NULL,  NULL);
1440       do_factor (factor, group);
1441     }
1442
1443   ok = casegrouper_destroy (grouper);
1444   ok = proc_commit (ds) && ok;
1445
1446   return ok;
1447 }
1448
1449
1450 /* Return the communality of variable N, calculated to N_FACTORS */
1451 static double
1452 the_communality (const gsl_matrix *evec, const gsl_vector *eval, int n, int n_factors)
1453 {
1454   size_t i;
1455
1456   double comm = 0;
1457
1458   assert (n >= 0);
1459   assert (n < eval->size);
1460   assert (n < evec->size1);
1461   assert (n_factors <= eval->size);
1462
1463   for (i = 0 ; i < n_factors; ++i)
1464     {
1465       double evali = fabs (gsl_vector_get (eval, i));
1466
1467       double eveci = gsl_matrix_get (evec, n, i);
1468
1469       comm += pow2 (eveci) * evali;
1470     }
1471
1472   return comm;
1473 }
1474
1475 /* Return the communality of variable N, calculated to N_FACTORS */
1476 static double
1477 communality (struct idata *idata, int n, int n_factors)
1478 {
1479   return the_communality (idata->evec, idata->eval, n, n_factors);
1480 }
1481
1482
1483 static void
1484 show_scree (const struct cmd_factor *f, struct idata *idata)
1485 {
1486   struct scree *s;
1487   const char *label ;
1488
1489   if ( !(f->plot & PLOT_SCREE) )
1490     return;
1491
1492
1493   label = f->extraction == EXTRACTION_PC ? _("Component Number") : _("Factor Number");
1494
1495   s = scree_create (idata->eval, label);
1496
1497   scree_submit (s);
1498 }
1499
1500 static void
1501 show_communalities (const struct cmd_factor * factor,
1502                     const gsl_vector *initial, const gsl_vector *extracted)
1503 {
1504   int i;
1505   int c = 0;
1506   const int heading_columns = 1;
1507   int nc = heading_columns;
1508   const int heading_rows = 1;
1509   const int nr = heading_rows + factor->n_vars;
1510   struct tab_table *t;
1511
1512   if (factor->print & PRINT_EXTRACTION)
1513     nc++;
1514
1515   if (factor->print & PRINT_INITIAL)
1516     nc++;
1517
1518   /* No point having a table with only headings */
1519   if (nc <= 1)
1520     return;
1521
1522   t = tab_create (nc, nr);
1523
1524   tab_title (t, _("Communalities"));
1525
1526   tab_headers (t, heading_columns, 0, heading_rows, 0);
1527
1528   c = 1;
1529   if (factor->print & PRINT_INITIAL)
1530     tab_text (t, c++, 0, TAB_CENTER | TAT_TITLE, _("Initial"));
1531
1532   if (factor->print & PRINT_EXTRACTION)
1533     tab_text (t, c++, 0, TAB_CENTER | TAT_TITLE, _("Extraction"));
1534
1535   /* Outline the box */
1536   tab_box (t,
1537            TAL_2, TAL_2,
1538            -1, -1,
1539            0, 0,
1540            nc - 1, nr - 1);
1541
1542   /* Vertical lines */
1543   tab_box (t,
1544            -1, -1,
1545            -1, TAL_1,
1546            heading_columns, 0,
1547            nc - 1, nr - 1);
1548
1549   tab_hline (t, TAL_1, 0, nc - 1, heading_rows);
1550   tab_vline (t, TAL_2, heading_columns, 0, nr - 1);
1551
1552   for (i = 0 ; i < factor->n_vars; ++i)
1553     {
1554       c = 0;
1555       tab_text (t, c++, i + heading_rows, TAT_TITLE, var_to_string (factor->vars[i]));
1556
1557       if (factor->print & PRINT_INITIAL)
1558         tab_double (t, c++, i + heading_rows, 0, gsl_vector_get (initial, i), NULL, RC_OTHER);
1559
1560       if (factor->print & PRINT_EXTRACTION)
1561         tab_double (t, c++, i + heading_rows, 0, gsl_vector_get (extracted, i), NULL, RC_OTHER);
1562     }
1563
1564   tab_submit (t);
1565 }
1566
1567
1568 static void
1569 show_factor_matrix (const struct cmd_factor *factor, struct idata *idata, const char *title, const gsl_matrix *fm)
1570 {
1571   int i;
1572
1573   const int n_factors = idata->n_extractions;
1574
1575   const int heading_columns = 1;
1576   const int heading_rows = 2;
1577   const int nr = heading_rows + factor->n_vars;
1578   const int nc = heading_columns + n_factors;
1579   gsl_permutation *perm;
1580
1581   struct tab_table *t = tab_create (nc, nr);
1582
1583   /*
1584   if ( factor->extraction == EXTRACTION_PC )
1585     tab_title (t, _("Component Matrix"));
1586   else
1587     tab_title (t, _("Factor Matrix"));
1588   */
1589
1590   tab_title (t, "%s", title);
1591
1592   tab_headers (t, heading_columns, 0, heading_rows, 0);
1593
1594   if ( factor->extraction == EXTRACTION_PC )
1595     tab_joint_text (t,
1596                     1, 0,
1597                     nc - 1, 0,
1598                     TAB_CENTER | TAT_TITLE, _("Component"));
1599   else
1600     tab_joint_text (t,
1601                     1, 0,
1602                     nc - 1, 0,
1603                     TAB_CENTER | TAT_TITLE, _("Factor"));
1604
1605
1606   tab_hline (t, TAL_1, heading_columns, nc - 1, 1);
1607
1608
1609   /* Outline the box */
1610   tab_box (t,
1611            TAL_2, TAL_2,
1612            -1, -1,
1613            0, 0,
1614            nc - 1, nr - 1);
1615
1616   /* Vertical lines */
1617   tab_box (t,
1618            -1, -1,
1619            -1, TAL_1,
1620            heading_columns, 1,
1621            nc - 1, nr - 1);
1622
1623   tab_hline (t, TAL_1, 0, nc - 1, heading_rows);
1624   tab_vline (t, TAL_2, heading_columns, 0, nr - 1);
1625
1626
1627   /* Initialise to the identity permutation */
1628   perm = gsl_permutation_calloc (factor->n_vars);
1629
1630   if ( factor->sort)
1631     sort_matrix_indirect (fm, perm);
1632
1633   for (i = 0 ; i < n_factors; ++i)
1634     {
1635       tab_text_format (t, heading_columns + i, 1, TAB_CENTER | TAT_TITLE, _("%d"), i + 1);
1636     }
1637
1638   for (i = 0 ; i < factor->n_vars; ++i)
1639     {
1640       int j;
1641       const int matrix_row = perm->data[i];
1642       tab_text (t, 0, i + heading_rows, TAT_TITLE, var_to_string (factor->vars[matrix_row]));
1643
1644       for (j = 0 ; j < n_factors; ++j)
1645         {
1646           double x = gsl_matrix_get (fm, matrix_row, j);
1647
1648           if ( fabs (x) < factor->blank)
1649             continue;
1650
1651           tab_double (t, heading_columns + j, heading_rows + i, 0, x, NULL, RC_OTHER);
1652         }
1653     }
1654
1655   gsl_permutation_free (perm);
1656
1657   tab_submit (t);
1658 }
1659
1660
1661 static void
1662 show_explained_variance (const struct cmd_factor * factor, struct idata *idata,
1663                          const gsl_vector *initial_eigenvalues,
1664                          const gsl_vector *extracted_eigenvalues,
1665                          const gsl_vector *rotated_loadings)
1666 {
1667   size_t i;
1668   int c = 0;
1669   const int heading_columns = 1;
1670   const int heading_rows = 2;
1671   const int nr = heading_rows + factor->n_vars;
1672
1673   struct tab_table *t ;
1674
1675   double i_total = 0.0;
1676   double i_cum = 0.0;
1677
1678   double e_total = 0.0;
1679   double e_cum = 0.0;
1680
1681   double r_cum = 0.0;
1682
1683   int nc = heading_columns;
1684
1685   if (factor->print & PRINT_EXTRACTION)
1686     nc += 3;
1687
1688   if (factor->print & PRINT_INITIAL)
1689     nc += 3;
1690
1691   if (factor->print & PRINT_ROTATION)
1692     {
1693       nc += factor->rotation == ROT_PROMAX ? 1 : 3;
1694     }
1695
1696   /* No point having a table with only headings */
1697   if ( nc <= heading_columns)
1698     return;
1699
1700   t = tab_create (nc, nr);
1701
1702   tab_title (t, _("Total Variance Explained"));
1703
1704   tab_headers (t, heading_columns, 0, heading_rows, 0);
1705
1706   /* Outline the box */
1707   tab_box (t,
1708            TAL_2, TAL_2,
1709            -1, -1,
1710            0, 0,
1711            nc - 1, nr - 1);
1712
1713   /* Vertical lines */
1714   tab_box (t,
1715            -1, -1,
1716            -1, TAL_1,
1717            heading_columns, 0,
1718            nc - 1, nr - 1);
1719
1720   tab_hline (t, TAL_1, 0, nc - 1, heading_rows);
1721   tab_hline (t, TAL_1, 1, nc - 1, 1);
1722
1723   tab_vline (t, TAL_2, heading_columns, 0, nr - 1);
1724
1725
1726   if ( factor->extraction == EXTRACTION_PC)
1727     tab_text (t, 0, 1, TAB_LEFT | TAT_TITLE, _("Component"));
1728   else
1729     tab_text (t, 0, 1, TAB_LEFT | TAT_TITLE, _("Factor"));
1730
1731   c = 1;
1732   if (factor->print & PRINT_INITIAL)
1733     {
1734       tab_joint_text (t, c, 0, c + 2, 0, TAB_CENTER | TAT_TITLE, _("Initial Eigenvalues"));
1735       c += 3;
1736     }
1737
1738   if (factor->print & PRINT_EXTRACTION)
1739     {
1740       tab_joint_text (t, c, 0, c + 2, 0, TAB_CENTER | TAT_TITLE, _("Extraction Sums of Squared Loadings"));
1741       c += 3;
1742     }
1743
1744   if (factor->print & PRINT_ROTATION)
1745     {
1746       const int width = factor->rotation == ROT_PROMAX ? 0 : 2;
1747       tab_joint_text (t, c, 0, c + width, 0, TAB_CENTER | TAT_TITLE, _("Rotation Sums of Squared Loadings"));
1748       c += width + 1;
1749     }
1750
1751   for (i = 0; i < (nc - heading_columns + 2) / 3 ; ++i)
1752     {
1753       tab_text (t, i * 3 + 1, 1, TAB_CENTER | TAT_TITLE, _("Total"));
1754
1755       tab_vline (t, TAL_2, heading_columns + i * 3, 0, nr - 1);
1756
1757       if (i == 2 && factor->rotation == ROT_PROMAX)
1758         continue;
1759
1760       /* xgettext:no-c-format */
1761       tab_text (t, i * 3 + 2, 1, TAB_CENTER | TAT_TITLE, _("% of Variance"));
1762       tab_text (t, i * 3 + 3, 1, TAB_CENTER | TAT_TITLE, _("Cumulative %"));
1763     }
1764
1765   for (i = 0 ; i < initial_eigenvalues->size; ++i)
1766     i_total += gsl_vector_get (initial_eigenvalues, i);
1767
1768   if ( factor->extraction == EXTRACTION_PAF)
1769     {
1770       e_total = factor->n_vars;
1771     }
1772   else
1773     {
1774       e_total = i_total;
1775     }
1776
1777   for (i = 0 ; i < factor->n_vars; ++i)
1778     {
1779       const double i_lambda = gsl_vector_get (initial_eigenvalues, i);
1780       double i_percent = 100.0 * i_lambda / i_total ;
1781
1782       const double e_lambda = gsl_vector_get (extracted_eigenvalues, i);
1783       double e_percent = 100.0 * e_lambda / e_total ;
1784
1785       c = 0;
1786
1787       tab_text_format (t, c++, i + heading_rows, TAB_LEFT | TAT_TITLE, _("%zu"), i + 1);
1788
1789       i_cum += i_percent;
1790       e_cum += e_percent;
1791
1792       /* Initial Eigenvalues */
1793       if (factor->print & PRINT_INITIAL)
1794       {
1795         tab_double (t, c++, i + heading_rows, 0, i_lambda, NULL, RC_OTHER);
1796         tab_double (t, c++, i + heading_rows, 0, i_percent, NULL, RC_OTHER);
1797         tab_double (t, c++, i + heading_rows, 0, i_cum, NULL, RC_OTHER);
1798       }
1799
1800
1801       if (factor->print & PRINT_EXTRACTION)
1802         {
1803           if (i < idata->n_extractions)
1804             {
1805               /* Sums of squared loadings */
1806               tab_double (t, c++, i + heading_rows, 0, e_lambda, NULL, RC_OTHER);
1807               tab_double (t, c++, i + heading_rows, 0, e_percent, NULL, RC_OTHER);
1808               tab_double (t, c++, i + heading_rows, 0, e_cum, NULL, RC_OTHER);
1809             }
1810         }
1811
1812       if (rotated_loadings != NULL)
1813         {
1814           const double r_lambda = gsl_vector_get (rotated_loadings, i);
1815           double r_percent = 100.0 * r_lambda / e_total ;
1816
1817           if (factor->print & PRINT_ROTATION)
1818             {
1819               if (i < idata->n_extractions)
1820                 {
1821                   r_cum += r_percent;
1822                   tab_double (t, c++, i + heading_rows, 0, r_lambda, NULL, RC_OTHER);
1823                   if (factor->rotation != ROT_PROMAX)
1824                     {
1825                       tab_double (t, c++, i + heading_rows, 0, r_percent, NULL, RC_OTHER);
1826                       tab_double (t, c++, i + heading_rows, 0, r_cum, NULL, RC_OTHER);
1827                     }
1828                 }
1829             }
1830         }
1831     }
1832
1833   tab_submit (t);
1834 }
1835
1836
1837 static void
1838 show_factor_correlation (const struct cmd_factor * factor, const gsl_matrix *fcm)
1839 {
1840   size_t i, j;
1841   const int heading_columns = 1;
1842   const int heading_rows = 1;
1843   const int nr = heading_rows + fcm->size2;
1844   const int nc = heading_columns + fcm->size1;
1845   struct tab_table *t = tab_create (nc, nr);
1846
1847   tab_title (t, _("Factor Correlation Matrix"));
1848
1849   tab_headers (t, heading_columns, 0, heading_rows, 0);
1850
1851   /* Outline the box */
1852   tab_box (t,
1853            TAL_2, TAL_2,
1854            -1, -1,
1855            0, 0,
1856            nc - 1, nr - 1);
1857
1858   /* Vertical lines */
1859   tab_box (t,
1860            -1, -1,
1861            -1, TAL_1,
1862            heading_columns, 0,
1863            nc - 1, nr - 1);
1864
1865   tab_hline (t, TAL_1, 0, nc - 1, heading_rows);
1866   tab_hline (t, TAL_1, 1, nc - 1, 1);
1867
1868   tab_vline (t, TAL_2, heading_columns, 0, nr - 1);
1869
1870
1871   if ( factor->extraction == EXTRACTION_PC)
1872     tab_text (t, 0, 0, TAB_LEFT | TAT_TITLE, _("Component"));
1873   else
1874     tab_text (t, 0, 0, TAB_LEFT | TAT_TITLE, _("Factor"));
1875
1876   for (i = 0 ; i < fcm->size1; ++i)
1877     {
1878       tab_text_format (t, heading_columns + i, 0, TAB_CENTER | TAT_TITLE, _("%zu"), i + 1);
1879     }
1880
1881   for (i = 0 ; i < fcm->size2; ++i)
1882     {
1883       tab_text_format (t, 0, heading_rows + i, TAB_CENTER | TAT_TITLE, _("%zu"), i + 1);
1884     }
1885
1886
1887   for (i = 0 ; i < fcm->size1; ++i)
1888     {
1889       for (j = 0 ; j < fcm->size2; ++j)
1890         tab_double (t, heading_columns + i,  heading_rows +j, 0,
1891                     gsl_matrix_get (fcm, i, j), NULL, RC_OTHER);
1892     }
1893
1894   tab_submit (t);
1895 }
1896
1897
1898 static void
1899 show_correlation_matrix (const struct cmd_factor *factor, const struct idata *idata)
1900 {
1901   struct tab_table *t ;
1902   size_t i, j;
1903   int y_pos_corr = -1;
1904   int y_pos_sig = -1;
1905   int suffix_rows = 0;
1906
1907   const int heading_rows = 1;
1908   const int heading_columns = 2;
1909
1910   int nc = heading_columns ;
1911   int nr = heading_rows ;
1912   int n_data_sets = 0;
1913
1914   if (factor->print & PRINT_CORRELATION)
1915     {
1916       y_pos_corr = n_data_sets;
1917       n_data_sets++;
1918       nc = heading_columns + factor->n_vars;
1919     }
1920
1921   if (factor->print & PRINT_SIG)
1922     {
1923       y_pos_sig = n_data_sets;
1924       n_data_sets++;
1925       nc = heading_columns + factor->n_vars;
1926     }
1927
1928   nr += n_data_sets * factor->n_vars;
1929
1930   if (factor->print & PRINT_DETERMINANT)
1931     suffix_rows = 1;
1932
1933   /* If the table would contain only headings, don't bother rendering it */
1934   if (nr <= heading_rows && suffix_rows == 0)
1935     return;
1936
1937   t = tab_create (nc, nr + suffix_rows);
1938
1939   tab_title (t, _("Correlation Matrix"));
1940
1941   tab_hline (t, TAL_1, 0, nc - 1, heading_rows);
1942
1943   if (nr > heading_rows)
1944     {
1945       tab_headers (t, heading_columns, 0, heading_rows, 0);
1946
1947       tab_vline (t, TAL_2, 2, 0, nr - 1);
1948
1949       /* Outline the box */
1950       tab_box (t,
1951                TAL_2, TAL_2,
1952                -1, -1,
1953                0, 0,
1954                nc - 1, nr - 1);
1955
1956       /* Vertical lines */
1957       tab_box (t,
1958                -1, -1,
1959                -1, TAL_1,
1960                heading_columns, 0,
1961                nc - 1, nr - 1);
1962
1963
1964       for (i = 0; i < factor->n_vars; ++i)
1965         tab_text (t, heading_columns + i, 0, TAT_TITLE, var_to_string (factor->vars[i]));
1966
1967
1968       for (i = 0 ; i < n_data_sets; ++i)
1969         {
1970           int y = heading_rows + i * factor->n_vars;
1971           size_t v;
1972           for (v = 0; v < factor->n_vars; ++v)
1973             tab_text (t, 1, y + v, TAT_TITLE, var_to_string (factor->vars[v]));
1974
1975           tab_hline (t, TAL_1, 0, nc - 1, y);
1976         }
1977
1978       if (factor->print & PRINT_CORRELATION)
1979         {
1980           const double y = heading_rows + y_pos_corr;
1981           tab_text (t, 0, y, TAT_TITLE, _("Correlations"));
1982
1983           for (i = 0; i < factor->n_vars; ++i)
1984             {
1985               for (j = 0; j < factor->n_vars; ++j)
1986                 tab_double (t, heading_columns + i,  y + j, 0, gsl_matrix_get (idata->corr, i, j), NULL, RC_OTHER);
1987             }
1988         }
1989
1990       if (factor->print & PRINT_SIG)
1991         {
1992           const double y = heading_rows + y_pos_sig * factor->n_vars;
1993           tab_text (t, 0, y, TAT_TITLE, _("Sig. (1-tailed)"));
1994
1995           for (i = 0; i < factor->n_vars; ++i)
1996             {
1997               for (j = 0; j < factor->n_vars; ++j)
1998                 {
1999                   double rho = gsl_matrix_get (idata->corr, i, j);
2000                   double w = gsl_matrix_get (idata->n, i, j);
2001
2002                   if (i == j)
2003                     continue;
2004
2005                   tab_double (t, heading_columns + i,  y + j, 0, significance_of_correlation (rho, w), NULL, RC_PVALUE);
2006                 }
2007             }
2008         }
2009     }
2010
2011   if (factor->print & PRINT_DETERMINANT)
2012     {
2013       tab_text (t, 0, nr, TAB_LEFT | TAT_TITLE, _("Determinant"));
2014
2015       tab_double (t, 1, nr, 0, idata->detR, NULL, RC_OTHER);
2016     }
2017
2018   tab_submit (t);
2019 }
2020
2021
2022
2023 static void
2024 do_factor (const struct cmd_factor *factor, struct casereader *r)
2025 {
2026   struct ccase *c;
2027   const gsl_matrix *var_matrix;
2028   const gsl_matrix *mean_matrix;
2029
2030   const gsl_matrix *analysis_matrix;
2031   struct idata *idata = idata_alloc (factor->n_vars);
2032
2033   struct covariance *cov = covariance_1pass_create (factor->n_vars, factor->vars,
2034                                               factor->wv, factor->exclude);
2035
2036   for ( ; (c = casereader_read (r) ); case_unref (c))
2037     {
2038       covariance_accumulate (cov, c);
2039     }
2040
2041   idata->cov = covariance_calculate (cov);
2042
2043   if (idata->cov == NULL)
2044     {
2045       msg (MW, _("The dataset contains no complete observations. No analysis will be performed."));
2046       covariance_destroy (cov);
2047       goto finish;
2048     }
2049
2050   var_matrix = covariance_moments (cov, MOMENT_VARIANCE);
2051   mean_matrix = covariance_moments (cov, MOMENT_MEAN);
2052   idata->n = covariance_moments (cov, MOMENT_NONE);
2053
2054
2055   if ( factor->method == METHOD_CORR)
2056     {
2057       idata->corr = correlation_from_covariance (idata->cov, var_matrix);
2058
2059       analysis_matrix = idata->corr;
2060     }
2061   else
2062     analysis_matrix = idata->cov;
2063
2064
2065   if (factor->print & PRINT_DETERMINANT
2066       || factor->print & PRINT_KMO)
2067     {
2068       int sign = 0;
2069
2070       const int size = idata->corr->size1;
2071       gsl_permutation *p = gsl_permutation_calloc (size);
2072       gsl_matrix *tmp = gsl_matrix_calloc (size, size);
2073       gsl_matrix_memcpy (tmp, idata->corr);
2074
2075       gsl_linalg_LU_decomp (tmp, p, &sign);
2076       idata->detR = gsl_linalg_LU_det (tmp, sign);
2077       gsl_permutation_free (p);
2078       gsl_matrix_free (tmp);
2079     }
2080
2081   if ( factor->print & PRINT_UNIVARIATE)
2082     {
2083       const struct fmt_spec *wfmt = factor->wv ? var_get_print_format (factor->wv) : & F_8_0;
2084       const int nc = 4;
2085       int i;
2086
2087       const int heading_columns = 1;
2088       const int heading_rows = 1;
2089
2090       const int nr = heading_rows + factor->n_vars;
2091
2092       struct tab_table *t = tab_create (nc, nr);
2093       tab_set_format (t, RC_WEIGHT, wfmt);
2094       tab_title (t, _("Descriptive Statistics"));
2095
2096       tab_headers (t, heading_columns, 0, heading_rows, 0);
2097
2098       /* Outline the box */
2099       tab_box (t,
2100                TAL_2, TAL_2,
2101                -1, -1,
2102                0, 0,
2103                nc - 1, nr - 1);
2104
2105       /* Vertical lines */
2106       tab_box (t,
2107                -1, -1,
2108                -1, TAL_1,
2109                heading_columns, 0,
2110                nc - 1, nr - 1);
2111
2112       tab_hline (t, TAL_1, 0, nc - 1, heading_rows);
2113       tab_vline (t, TAL_2, heading_columns, 0, nr - 1);
2114
2115       tab_text (t, 1, 0, TAB_CENTER | TAT_TITLE, _("Mean"));
2116       tab_text (t, 2, 0, TAB_CENTER | TAT_TITLE, _("Std. Deviation"));
2117       tab_text (t, 3, 0, TAB_CENTER | TAT_TITLE, _("Analysis N"));
2118
2119       for (i = 0 ; i < factor->n_vars; ++i)
2120         {
2121           const struct variable *v = factor->vars[i];
2122           tab_text (t, 0, i + heading_rows, TAB_LEFT | TAT_TITLE, var_to_string (v));
2123
2124           tab_double (t, 1, i + heading_rows, 0, gsl_matrix_get (mean_matrix, i, i), NULL, RC_OTHER);
2125           tab_double (t, 2, i + heading_rows, 0, sqrt (gsl_matrix_get (var_matrix, i, i)), NULL, RC_OTHER);
2126           tab_double (t, 3, i + heading_rows, 0, gsl_matrix_get (idata->n, i, i), NULL, RC_WEIGHT);
2127         }
2128
2129       tab_submit (t);
2130     }
2131
2132   if (factor->print & PRINT_KMO)
2133     {
2134       int i;
2135       double sum_ssq_r = 0;
2136       double sum_ssq_a = 0;
2137
2138       double df = factor->n_vars * ( factor->n_vars - 1) / 2;
2139
2140       double w = 0;
2141
2142
2143       double xsq;
2144
2145       const int heading_columns = 2;
2146       const int heading_rows = 0;
2147
2148       const int nr = heading_rows + 4;
2149       const int nc = heading_columns + 1;
2150
2151       gsl_matrix *a, *x;
2152
2153       struct tab_table *t = tab_create (nc, nr);
2154       tab_title (t, _("KMO and Bartlett's Test"));
2155
2156       x  = clone_matrix (idata->corr);
2157       gsl_linalg_cholesky_decomp (x);
2158       gsl_linalg_cholesky_invert (x);
2159
2160       a = anti_image (x);
2161
2162       for (i = 0; i < x->size1; ++i)
2163         {
2164           sum_ssq_r += ssq_od_n (x, i);
2165           sum_ssq_a += ssq_od_n (a, i);
2166         }
2167
2168       gsl_matrix_free (a);
2169       gsl_matrix_free (x);
2170
2171       tab_headers (t, heading_columns, 0, heading_rows, 0);
2172
2173       /* Outline the box */
2174       tab_box (t,
2175                TAL_2, TAL_2,
2176                -1, -1,
2177                0, 0,
2178                nc - 1, nr - 1);
2179
2180       tab_vline (t, TAL_2, heading_columns, 0, nr - 1);
2181
2182       tab_text (t, 0, 0, TAT_TITLE | TAB_LEFT, _("Kaiser-Meyer-Olkin Measure of Sampling Adequacy"));
2183
2184       tab_double (t, 2, 0, 0, sum_ssq_r /  (sum_ssq_r + sum_ssq_a), NULL, RC_OTHER);
2185
2186       tab_text (t, 0, 1, TAT_TITLE | TAB_LEFT, _("Bartlett's Test of Sphericity"));
2187
2188       tab_text (t, 1, 1, TAT_TITLE, _("Approx. Chi-Square"));
2189       tab_text (t, 1, 2, TAT_TITLE, _("df"));
2190       tab_text (t, 1, 3, TAT_TITLE, _("Sig."));
2191
2192
2193       /* The literature doesn't say what to do for the value of W when
2194          missing values are involved.  The best thing I can think of
2195          is to take the mean average. */
2196       w = 0;
2197       for (i = 0; i < idata->n->size1; ++i)
2198         w += gsl_matrix_get (idata->n, i, i);
2199       w /= idata->n->size1;
2200
2201       xsq = w - 1 - (2 * factor->n_vars + 5) / 6.0;
2202       xsq *= -log (idata->detR);
2203
2204       tab_double (t, 2, 1, 0, xsq, NULL, RC_OTHER);
2205       tab_double (t, 2, 2, 0, df, NULL, RC_INTEGER);
2206       tab_double (t, 2, 3, 0, gsl_cdf_chisq_Q (xsq, df), NULL, RC_PVALUE);
2207
2208
2209       tab_submit (t);
2210     }
2211
2212   show_correlation_matrix (factor, idata);
2213   covariance_destroy (cov);
2214
2215   {
2216     gsl_matrix *am = matrix_dup (analysis_matrix);
2217     gsl_eigen_symmv_workspace *workspace = gsl_eigen_symmv_alloc (factor->n_vars);
2218
2219     gsl_eigen_symmv (am, idata->eval, idata->evec, workspace);
2220
2221     gsl_eigen_symmv_free (workspace);
2222     gsl_matrix_free (am);
2223   }
2224
2225   gsl_eigen_symmv_sort (idata->eval, idata->evec, GSL_EIGEN_SORT_ABS_DESC);
2226
2227   idata->n_extractions = n_extracted_factors (factor, idata);
2228
2229   if (idata->n_extractions == 0)
2230     {
2231       msg (MW, _("The %s criteria result in zero factors extracted. Therefore no analysis will be performed."), "FACTOR");
2232       goto finish;
2233     }
2234
2235   if (idata->n_extractions > factor->n_vars)
2236     {
2237       msg (MW,
2238            _("The %s criteria result in more factors than variables, which is not meaningful. No analysis will be performed."),
2239            "FACTOR");
2240       goto finish;
2241     }
2242
2243   {
2244     gsl_matrix *rotated_factors = NULL;
2245     gsl_matrix *pattern_matrix = NULL;
2246     gsl_matrix *fcm = NULL;
2247     gsl_vector *rotated_loadings = NULL;
2248
2249     const gsl_vector *extracted_eigenvalues = NULL;
2250     gsl_vector *initial_communalities = gsl_vector_alloc (factor->n_vars);
2251     gsl_vector *extracted_communalities = gsl_vector_alloc (factor->n_vars);
2252     size_t i;
2253     struct factor_matrix_workspace *fmw = factor_matrix_workspace_alloc (idata->msr->size, idata->n_extractions);
2254     gsl_matrix *factor_matrix = gsl_matrix_calloc (factor->n_vars, fmw->n_factors);
2255
2256     if ( factor->extraction == EXTRACTION_PAF)
2257       {
2258         gsl_vector *diff = gsl_vector_alloc (idata->msr->size);
2259         struct smr_workspace *ws = ws_create (analysis_matrix);
2260
2261         for (i = 0 ; i < factor->n_vars ; ++i)
2262           {
2263             double r2 = squared_multiple_correlation (analysis_matrix, i, ws);
2264
2265             gsl_vector_set (idata->msr, i, r2);
2266           }
2267         ws_destroy (ws);
2268
2269         gsl_vector_memcpy (initial_communalities, idata->msr);
2270
2271         for (i = 0; i < factor->extraction_iterations; ++i)
2272           {
2273             double min, max;
2274             gsl_vector_memcpy (diff, idata->msr);
2275
2276             iterate_factor_matrix (analysis_matrix, idata->msr, factor_matrix, fmw);
2277
2278             gsl_vector_sub (diff, idata->msr);
2279
2280             gsl_vector_minmax (diff, &min, &max);
2281
2282             if ( fabs (min) < factor->econverge && fabs (max) < factor->econverge)
2283               break;
2284           }
2285         gsl_vector_free (diff);
2286
2287
2288
2289         gsl_vector_memcpy (extracted_communalities, idata->msr);
2290         extracted_eigenvalues = fmw->eval;
2291       }
2292     else if (factor->extraction == EXTRACTION_PC)
2293       {
2294         for (i = 0; i < factor->n_vars; ++i)
2295           gsl_vector_set (initial_communalities, i, communality (idata, i, factor->n_vars));
2296
2297         gsl_vector_memcpy (extracted_communalities, initial_communalities);
2298
2299         iterate_factor_matrix (analysis_matrix, extracted_communalities, factor_matrix, fmw);
2300
2301
2302         extracted_eigenvalues = idata->eval;
2303       }
2304
2305
2306     show_communalities (factor, initial_communalities, extracted_communalities);
2307
2308
2309     if ( factor->rotation != ROT_NONE)
2310       {
2311         rotated_factors = gsl_matrix_calloc (factor_matrix->size1, factor_matrix->size2);
2312         rotated_loadings = gsl_vector_calloc (factor_matrix->size2);
2313         if (factor->rotation == ROT_PROMAX)
2314           {
2315             pattern_matrix = gsl_matrix_calloc (factor_matrix->size1, factor_matrix->size2);
2316             fcm = gsl_matrix_calloc (factor_matrix->size2, factor_matrix->size2);
2317           }
2318
2319
2320         rotate (factor, factor_matrix, extracted_communalities, rotated_factors, rotated_loadings, pattern_matrix, fcm);
2321       }
2322
2323     show_explained_variance (factor, idata, idata->eval, extracted_eigenvalues, rotated_loadings);
2324
2325     factor_matrix_workspace_free (fmw);
2326
2327     show_scree (factor, idata);
2328
2329     show_factor_matrix (factor, idata,
2330                         factor->extraction == EXTRACTION_PC ? _("Component Matrix") : _("Factor Matrix"),
2331                         factor_matrix);
2332
2333     if ( factor->rotation == ROT_PROMAX)
2334       {
2335         show_factor_matrix (factor, idata, _("Pattern Matrix"),  pattern_matrix);
2336         gsl_matrix_free (pattern_matrix);
2337       }
2338
2339     if ( factor->rotation != ROT_NONE)
2340       {
2341         show_factor_matrix (factor, idata,
2342                             (factor->rotation == ROT_PROMAX) ? _("Structure Matrix") :
2343                             (factor->extraction == EXTRACTION_PC ? _("Rotated Component Matrix") : _("Rotated Factor Matrix")),
2344                             rotated_factors);
2345
2346         gsl_matrix_free (rotated_factors);
2347       }
2348
2349     if ( factor->rotation == ROT_PROMAX)
2350       {
2351         show_factor_correlation (factor, fcm);
2352         gsl_matrix_free (fcm);
2353       }
2354
2355     gsl_matrix_free (factor_matrix);
2356     gsl_vector_free (rotated_loadings);
2357     gsl_vector_free (initial_communalities);
2358     gsl_vector_free (extracted_communalities);
2359   }
2360
2361  finish:
2362
2363   idata_free (idata);
2364
2365   casereader_destroy (r);
2366 }
2367
2368