c3a7e8ee2418b889499080e283d31c3dbe1c5cd3
[pspp] / src / language / stats / factor.c
1 /* PSPP - a program for statistical analysis.
2    Copyright (C) 2009, 2010, 2011, 2012 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 1
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
910   /* reflect negative sums and populate the rotated loadings vector*/
911   for (i = 0 ; i < result->size2; ++i)
912     {
913       double ssq = 0.0;
914       double sum = 0.0;
915       for (j = 0 ; j < result->size1; ++j)
916         {
917           double s = gsl_matrix_get (result, j, i);
918           ssq += s * s;
919           sum += s;
920         }
921
922       gsl_vector_set (rotated_loadings, i, ssq);
923
924       if ( sum < 0 )
925         for (j = 0 ; j < result->size1; ++j)
926           {
927             double *lambda = gsl_matrix_ptr (result, j, i);
928             *lambda = - *lambda;
929           }
930     }
931 }
932
933
934 /*
935   Get an approximation for the factor matrix into FACTORS, and the communalities into COMMUNALITIES.
936   R is the matrix to be analysed.
937   WS is a pointer to a structure which must have been initialised with factor_matrix_workspace_init.
938  */
939 static void
940 iterate_factor_matrix (const gsl_matrix *r, gsl_vector *communalities, gsl_matrix *factors, 
941                        struct factor_matrix_workspace *ws)
942 {
943   size_t i;
944   gsl_matrix_view mv ;
945
946   assert (r->size1 == r->size2);
947   assert (r->size1 == communalities->size);
948
949   assert (factors->size1 == r->size1);
950   assert (factors->size2 == ws->n_factors);
951
952   gsl_matrix_memcpy (ws->r, r);
953
954   /* Apply Communalities to diagonal of correlation matrix */
955   for (i = 0 ; i < communalities->size ; ++i)
956     {
957       double *x = gsl_matrix_ptr (ws->r, i, i);
958       *x = gsl_vector_get (communalities, i);
959     }
960
961   gsl_eigen_symmv (ws->r, ws->eval, ws->evec, ws->eigen_ws);
962
963   mv = gsl_matrix_submatrix (ws->evec, 0, 0, ws->evec->size1, ws->n_factors);
964
965   /* Gamma is the diagonal matrix containing the absolute values of the eigenvalues */
966   for (i = 0 ; i < ws->n_factors ; ++i)
967     {
968       double *ptr = gsl_matrix_ptr (ws->gamma, i, i);
969       *ptr = fabs (gsl_vector_get (ws->eval, i));
970     }
971
972   /* Take the square root of gamma */
973   gsl_linalg_cholesky_decomp (ws->gamma);
974
975   gsl_blas_dgemm (CblasNoTrans,  CblasNoTrans, 1.0, &mv.matrix, ws->gamma, 0.0, factors);
976
977   for (i = 0 ; i < r->size1 ; ++i)
978     {
979       double h = the_communality (ws->evec, ws->eval, i, ws->n_factors);
980       gsl_vector_set (communalities, i, h);
981     }
982 }
983
984
985
986 static bool run_factor (struct dataset *ds, const struct cmd_factor *factor);
987
988
989 int
990 cmd_factor (struct lexer *lexer, struct dataset *ds)
991 {
992   const struct dictionary *dict = dataset_dict (ds);
993   int n_iterations = 25;
994   struct cmd_factor factor;
995   factor.n_vars = 0;
996   factor.vars = NULL;
997   factor.method = METHOD_CORR;
998   factor.missing_type = MISS_LISTWISE;
999   factor.exclude = MV_ANY;
1000   factor.print = PRINT_INITIAL | PRINT_EXTRACTION | PRINT_ROTATION;
1001   factor.extraction = EXTRACTION_PC;
1002   factor.n_factors = 0;
1003   factor.min_eigen = SYSMIS;
1004   factor.extraction_iterations = 25;
1005   factor.rotation_iterations = 25;
1006   factor.econverge = 0.001;
1007
1008   factor.blank = 0;
1009   factor.sort = false;
1010   factor.plot = 0;
1011   factor.rotation = ROT_VARIMAX;
1012
1013   factor.rconverge = 0.0001;
1014
1015   factor.wv = dict_get_weight (dict);
1016
1017   lex_match (lexer, T_SLASH);
1018
1019   if (!lex_force_match_id (lexer, "VARIABLES"))
1020     {
1021       goto error;
1022     }
1023
1024   lex_match (lexer, T_EQUALS);
1025
1026   if (!parse_variables_const (lexer, dict, &factor.vars, &factor.n_vars,
1027                               PV_NO_DUPLICATE | PV_NUMERIC))
1028     goto error;
1029
1030   if (factor.n_vars < 2)
1031     msg (MW, _("Factor analysis on a single variable is not useful."));
1032
1033   while (lex_token (lexer) != T_ENDCMD)
1034     {
1035       lex_match (lexer, T_SLASH);
1036
1037       if (lex_match_id (lexer, "PLOT"))
1038         {
1039           lex_match (lexer, T_EQUALS);
1040           while (lex_token (lexer) != T_ENDCMD && lex_token (lexer) != T_SLASH)
1041             {
1042               if (lex_match_id (lexer, "EIGEN"))
1043                 {
1044                   factor.plot |= PLOT_SCREE;
1045                 }
1046 #if FACTOR_FULLY_IMPLEMENTED
1047               else if (lex_match_id (lexer, "ROTATION"))
1048                 {
1049                 }
1050 #endif
1051               else
1052                 {
1053                   lex_error (lexer, NULL);
1054                   goto error;
1055                 }
1056             }
1057         }
1058       else if (lex_match_id (lexer, "METHOD"))
1059         {
1060           lex_match (lexer, T_EQUALS);
1061           while (lex_token (lexer) != T_ENDCMD && lex_token (lexer) != T_SLASH)
1062             {
1063               if (lex_match_id (lexer, "COVARIANCE"))
1064                 {
1065                   factor.method = METHOD_COV;
1066                 }
1067               else if (lex_match_id (lexer, "CORRELATION"))
1068                 {
1069                   factor.method = METHOD_CORR;
1070                 }
1071               else
1072                 {
1073                   lex_error (lexer, NULL);
1074                   goto error;
1075                 }
1076             }
1077         }
1078       else if (lex_match_id (lexer, "ROTATION"))
1079         {
1080           lex_match (lexer, T_EQUALS);
1081           while (lex_token (lexer) != T_ENDCMD && lex_token (lexer) != T_SLASH)
1082             {
1083               /* VARIMAX and DEFAULT are defaults */
1084               if (lex_match_id (lexer, "VARIMAX") || lex_match_id (lexer, "DEFAULT"))
1085                 {
1086                   factor.rotation = ROT_VARIMAX;
1087                 }
1088               else if (lex_match_id (lexer, "EQUAMAX"))
1089                 {
1090                   factor.rotation = ROT_EQUAMAX;
1091                 }
1092               else if (lex_match_id (lexer, "QUARTIMAX"))
1093                 {
1094                   factor.rotation = ROT_QUARTIMAX;
1095                 }
1096               else if (lex_match_id (lexer, "PROMAX"))
1097                 {
1098                   factor.promax_power = 5;
1099                   if (lex_match (lexer, T_LPAREN))
1100                     {
1101                       lex_force_int (lexer);
1102                       factor.promax_power = lex_integer (lexer);
1103                       lex_get (lexer);
1104                       lex_force_match (lexer, T_RPAREN);
1105                     }
1106                   factor.rotation = ROT_PROMAX;
1107                 }
1108               else if (lex_match_id (lexer, "NOROTATE"))
1109                 {
1110                   factor.rotation = ROT_NONE;
1111                 }
1112               else
1113                 {
1114                   lex_error (lexer, NULL);
1115                   goto error;
1116                 }
1117             }
1118           factor.rotation_iterations = n_iterations;
1119         }
1120       else if (lex_match_id (lexer, "CRITERIA"))
1121         {
1122           lex_match (lexer, T_EQUALS);
1123           while (lex_token (lexer) != T_ENDCMD && lex_token (lexer) != T_SLASH)
1124             {
1125               if (lex_match_id (lexer, "FACTORS"))
1126                 {
1127                   if ( lex_force_match (lexer, T_LPAREN))
1128                     {
1129                       lex_force_int (lexer);
1130                       factor.n_factors = lex_integer (lexer);
1131                       lex_get (lexer);
1132                       lex_force_match (lexer, T_RPAREN);
1133                     }
1134                 }
1135               else if (lex_match_id (lexer, "MINEIGEN"))
1136                 {
1137                   if ( lex_force_match (lexer, T_LPAREN))
1138                     {
1139                       lex_force_num (lexer);
1140                       factor.min_eigen = lex_number (lexer);
1141                       lex_get (lexer);
1142                       lex_force_match (lexer, T_RPAREN);
1143                     }
1144                 }
1145               else if (lex_match_id (lexer, "ECONVERGE"))
1146                 {
1147                   if ( lex_force_match (lexer, T_LPAREN))
1148                     {
1149                       lex_force_num (lexer);
1150                       factor.econverge = lex_number (lexer);
1151                       lex_get (lexer);
1152                       lex_force_match (lexer, T_RPAREN);
1153                     }
1154                 }
1155               else if (lex_match_id (lexer, "RCONVERGE"))
1156                 {
1157                   if ( lex_force_match (lexer, T_LPAREN))
1158                     {
1159                       lex_force_num (lexer);
1160                       factor.rconverge = lex_number (lexer);
1161                       lex_get (lexer);
1162                       lex_force_match (lexer, T_RPAREN);
1163                     }
1164                 }
1165               else if (lex_match_id (lexer, "ITERATE"))
1166                 {
1167                   if ( lex_force_match (lexer, T_LPAREN))
1168                     {
1169                       lex_force_int (lexer);
1170                       n_iterations = lex_integer (lexer);
1171                       lex_get (lexer);
1172                       lex_force_match (lexer, T_RPAREN);
1173                     }
1174                 }
1175               else if (lex_match_id (lexer, "DEFAULT"))
1176                 {
1177                   factor.n_factors = 0;
1178                   factor.min_eigen = 1;
1179                   n_iterations = 25;
1180                 }
1181               else
1182                 {
1183                   lex_error (lexer, NULL);
1184                   goto error;
1185                 }
1186             }
1187         }
1188       else if (lex_match_id (lexer, "EXTRACTION"))
1189         {
1190           lex_match (lexer, T_EQUALS);
1191           while (lex_token (lexer) != T_ENDCMD && lex_token (lexer) != T_SLASH)
1192             {
1193               if (lex_match_id (lexer, "PAF"))
1194                 {
1195                   factor.extraction = EXTRACTION_PAF;
1196                 }
1197               else if (lex_match_id (lexer, "PC"))
1198                 {
1199                   factor.extraction = EXTRACTION_PC;
1200                 }
1201               else if (lex_match_id (lexer, "PA1"))
1202                 {
1203                   factor.extraction = EXTRACTION_PC;
1204                 }
1205               else if (lex_match_id (lexer, "DEFAULT"))
1206                 {
1207                   factor.extraction = EXTRACTION_PC;
1208                 }
1209               else
1210                 {
1211                   lex_error (lexer, NULL);
1212                   goto error;
1213                 }
1214             }
1215           factor.extraction_iterations = n_iterations;
1216         }
1217       else if (lex_match_id (lexer, "FORMAT"))
1218         {
1219           lex_match (lexer, T_EQUALS);
1220           while (lex_token (lexer) != T_ENDCMD && lex_token (lexer) != T_SLASH)
1221             {
1222               if (lex_match_id (lexer, "SORT"))
1223                 {
1224                   factor.sort = true;
1225                 }
1226               else if (lex_match_id (lexer, "BLANK"))
1227                 {
1228                   if ( lex_force_match (lexer, T_LPAREN))
1229                     {
1230                       lex_force_num (lexer);
1231                       factor.blank = lex_number (lexer);
1232                       lex_get (lexer);
1233                       lex_force_match (lexer, T_RPAREN);
1234                     }
1235                 }
1236               else if (lex_match_id (lexer, "DEFAULT"))
1237                 {
1238                   factor.blank = 0;
1239                   factor.sort = false;
1240                 }
1241               else
1242                 {
1243                   lex_error (lexer, NULL);
1244                   goto error;
1245                 }
1246             }
1247         }
1248       else if (lex_match_id (lexer, "PRINT"))
1249         {
1250           factor.print = 0;
1251           lex_match (lexer, T_EQUALS);
1252           while (lex_token (lexer) != T_ENDCMD && lex_token (lexer) != T_SLASH)
1253             {
1254               if (lex_match_id (lexer, "UNIVARIATE"))
1255                 {
1256                   factor.print |= PRINT_UNIVARIATE;
1257                 }
1258               else if (lex_match_id (lexer, "DET"))
1259                 {
1260                   factor.print |= PRINT_DETERMINANT;
1261                 }
1262 #if FACTOR_FULLY_IMPLEMENTED
1263               else if (lex_match_id (lexer, "INV"))
1264                 {
1265                 }
1266               else if (lex_match_id (lexer, "AIC"))
1267                 {
1268                 }
1269 #endif
1270               else if (lex_match_id (lexer, "SIG"))
1271                 {
1272                   factor.print |= PRINT_SIG;
1273                 }
1274               else if (lex_match_id (lexer, "CORRELATION"))
1275                 {
1276                   factor.print |= PRINT_CORRELATION;
1277                 }
1278 #if FACTOR_FULLY_IMPLEMENTED
1279               else if (lex_match_id (lexer, "COVARIANCE"))
1280                 {
1281                 }
1282 #endif
1283               else if (lex_match_id (lexer, "ROTATION"))
1284                 {
1285                   factor.print |= PRINT_ROTATION;
1286                 }
1287               else if (lex_match_id (lexer, "EXTRACTION"))
1288                 {
1289                   factor.print |= PRINT_EXTRACTION;
1290                 }
1291               else if (lex_match_id (lexer, "INITIAL"))
1292                 {
1293                   factor.print |= PRINT_INITIAL;
1294                 }
1295               else if (lex_match_id (lexer, "KMO"))
1296                 {
1297                   factor.print |= PRINT_KMO;
1298                 }
1299 #if FACTOR_FULLY_IMPLEMENTED
1300               else if (lex_match_id (lexer, "REPR"))
1301                 {
1302                 }
1303               else if (lex_match_id (lexer, "FSCORE"))
1304                 {
1305                 }
1306 #endif
1307               else if (lex_match (lexer, T_ALL))
1308                 {
1309                   factor.print = 0xFFFF;
1310                 }
1311               else if (lex_match_id (lexer, "DEFAULT"))
1312                 {
1313                   factor.print |= PRINT_INITIAL ;
1314                   factor.print |= PRINT_EXTRACTION ;
1315                   factor.print |= PRINT_ROTATION ;
1316                 }
1317               else
1318                 {
1319                   lex_error (lexer, NULL);
1320                   goto error;
1321                 }
1322             }
1323         }
1324       else if (lex_match_id (lexer, "MISSING"))
1325         {
1326           lex_match (lexer, T_EQUALS);
1327           while (lex_token (lexer) != T_ENDCMD && lex_token (lexer) != T_SLASH)
1328             {
1329               if (lex_match_id (lexer, "INCLUDE"))
1330                 {
1331                   factor.exclude = MV_SYSTEM;
1332                 }
1333               else if (lex_match_id (lexer, "EXCLUDE"))
1334                 {
1335                   factor.exclude = MV_ANY;
1336                 }
1337               else if (lex_match_id (lexer, "LISTWISE"))
1338                 {
1339                   factor.missing_type = MISS_LISTWISE;
1340                 }
1341               else if (lex_match_id (lexer, "PAIRWISE"))
1342                 {
1343                   factor.missing_type = MISS_PAIRWISE;
1344                 }
1345               else if (lex_match_id (lexer, "MEANSUB"))
1346                 {
1347                   factor.missing_type = MISS_MEANSUB;
1348                 }
1349               else
1350                 {
1351                   lex_error (lexer, NULL);
1352                   goto error;
1353                 }
1354             }
1355         }
1356       else
1357         {
1358           lex_error (lexer, NULL);
1359           goto error;
1360         }
1361     }
1362
1363   if ( factor.rotation == ROT_NONE )
1364     factor.print &= ~PRINT_ROTATION;
1365
1366   if ( ! run_factor (ds, &factor)) 
1367     goto error;
1368
1369   free (factor.vars);
1370   return CMD_SUCCESS;
1371
1372  error:
1373   free (factor.vars);
1374   return CMD_FAILURE;
1375 }
1376
1377 static void do_factor (const struct cmd_factor *factor, struct casereader *group);
1378
1379
1380 static bool
1381 run_factor (struct dataset *ds, const struct cmd_factor *factor)
1382 {
1383   struct dictionary *dict = dataset_dict (ds);
1384   bool ok;
1385   struct casereader *group;
1386
1387   struct casegrouper *grouper = casegrouper_create_splits (proc_open (ds), dict);
1388
1389   while (casegrouper_get_next_group (grouper, &group))
1390     {
1391       if ( factor->missing_type == MISS_LISTWISE )
1392         group  = casereader_create_filter_missing (group, factor->vars, factor->n_vars,
1393                                                    factor->exclude,
1394                                                    NULL,  NULL);
1395       do_factor (factor, group);
1396     }
1397
1398   ok = casegrouper_destroy (grouper);
1399   ok = proc_commit (ds) && ok;
1400
1401   return ok;
1402 }
1403
1404
1405 /* Return the communality of variable N, calculated to N_FACTORS */
1406 static double
1407 the_communality (const gsl_matrix *evec, const gsl_vector *eval, int n, int n_factors)
1408 {
1409   size_t i;
1410
1411   double comm = 0;
1412
1413   assert (n >= 0);
1414   assert (n < eval->size);
1415   assert (n < evec->size1);
1416   assert (n_factors <= eval->size);
1417
1418   for (i = 0 ; i < n_factors; ++i)
1419     {
1420       double evali = fabs (gsl_vector_get (eval, i));
1421
1422       double eveci = gsl_matrix_get (evec, n, i);
1423
1424       comm += pow2 (eveci) * evali;
1425     }
1426
1427   return comm;
1428 }
1429
1430 /* Return the communality of variable N, calculated to N_FACTORS */
1431 static double
1432 communality (struct idata *idata, int n, int n_factors)
1433 {
1434   return the_communality (idata->evec, idata->eval, n, n_factors);
1435 }
1436
1437
1438 static void
1439 show_scree (const struct cmd_factor *f, struct idata *idata)
1440 {
1441   struct scree *s;
1442   const char *label ;
1443
1444   if ( !(f->plot & PLOT_SCREE) )
1445     return;
1446
1447
1448   label = f->extraction == EXTRACTION_PC ? _("Component Number") : _("Factor Number");
1449
1450   s = scree_create (idata->eval, label);
1451
1452   scree_submit (s);
1453 }
1454
1455 static void
1456 show_communalities (const struct cmd_factor * factor,
1457                     const gsl_vector *initial, const gsl_vector *extracted)
1458 {
1459   int i;
1460   int c = 0;
1461   const int heading_columns = 1;
1462   int nc = heading_columns;
1463   const int heading_rows = 1;
1464   const int nr = heading_rows + factor->n_vars;
1465   struct tab_table *t;
1466
1467   if (factor->print & PRINT_EXTRACTION)
1468     nc++;
1469
1470   if (factor->print & PRINT_INITIAL)
1471     nc++;
1472
1473   /* No point having a table with only headings */
1474   if (nc <= 1)
1475     return;
1476
1477   t = tab_create (nc, nr);
1478
1479   tab_title (t, _("Communalities"));
1480
1481   tab_headers (t, heading_columns, 0, heading_rows, 0);
1482
1483   c = 1;
1484   if (factor->print & PRINT_INITIAL)
1485     tab_text (t, c++, 0, TAB_CENTER | TAT_TITLE, _("Initial"));
1486
1487   if (factor->print & PRINT_EXTRACTION)
1488     tab_text (t, c++, 0, TAB_CENTER | TAT_TITLE, _("Extraction"));
1489
1490   /* Outline the box */
1491   tab_box (t,
1492            TAL_2, TAL_2,
1493            -1, -1,
1494            0, 0,
1495            nc - 1, nr - 1);
1496
1497   /* Vertical lines */
1498   tab_box (t,
1499            -1, -1,
1500            -1, TAL_1,
1501            heading_columns, 0,
1502            nc - 1, nr - 1);
1503
1504   tab_hline (t, TAL_1, 0, nc - 1, heading_rows);
1505   tab_vline (t, TAL_2, heading_columns, 0, nr - 1);
1506
1507   for (i = 0 ; i < factor->n_vars; ++i)
1508     {
1509       c = 0;
1510       tab_text (t, c++, i + heading_rows, TAT_TITLE, var_to_string (factor->vars[i]));
1511
1512       if (factor->print & PRINT_INITIAL)
1513         tab_double (t, c++, i + heading_rows, 0, gsl_vector_get (initial, i), NULL, RC_OTHER);
1514
1515       if (factor->print & PRINT_EXTRACTION)
1516         tab_double (t, c++, i + heading_rows, 0, gsl_vector_get (extracted, i), NULL, RC_OTHER);
1517     }
1518
1519   tab_submit (t);
1520 }
1521
1522
1523 static void
1524 show_factor_matrix (const struct cmd_factor *factor, struct idata *idata, const char *title, const gsl_matrix *fm)
1525 {
1526   int i;
1527
1528   const int n_factors = idata->n_extractions;
1529
1530   const int heading_columns = 1;
1531   const int heading_rows = 2;
1532   const int nr = heading_rows + factor->n_vars;
1533   const int nc = heading_columns + n_factors;
1534   gsl_permutation *perm;
1535
1536   struct tab_table *t = tab_create (nc, nr);
1537
1538   /* 
1539   if ( factor->extraction == EXTRACTION_PC )
1540     tab_title (t, _("Component Matrix"));
1541   else 
1542     tab_title (t, _("Factor Matrix"));
1543   */
1544
1545   tab_title (t, "%s", title);
1546
1547   tab_headers (t, heading_columns, 0, heading_rows, 0);
1548
1549   if ( factor->extraction == EXTRACTION_PC )
1550     tab_joint_text (t,
1551                     1, 0,
1552                     nc - 1, 0,
1553                     TAB_CENTER | TAT_TITLE, _("Component"));
1554   else
1555     tab_joint_text (t,
1556                     1, 0,
1557                     nc - 1, 0,
1558                     TAB_CENTER | TAT_TITLE, _("Factor"));
1559
1560
1561   tab_hline (t, TAL_1, heading_columns, nc - 1, 1);
1562
1563
1564   /* Outline the box */
1565   tab_box (t,
1566            TAL_2, TAL_2,
1567            -1, -1,
1568            0, 0,
1569            nc - 1, nr - 1);
1570
1571   /* Vertical lines */
1572   tab_box (t,
1573            -1, -1,
1574            -1, TAL_1,
1575            heading_columns, 1,
1576            nc - 1, nr - 1);
1577
1578   tab_hline (t, TAL_1, 0, nc - 1, heading_rows);
1579   tab_vline (t, TAL_2, heading_columns, 0, nr - 1);
1580
1581
1582   /* Initialise to the identity permutation */
1583   perm = gsl_permutation_calloc (factor->n_vars);
1584
1585   if ( factor->sort)
1586     sort_matrix_indirect (fm, perm);
1587
1588   for (i = 0 ; i < n_factors; ++i)
1589     {
1590       tab_text_format (t, heading_columns + i, 1, TAB_CENTER | TAT_TITLE, _("%d"), i + 1);
1591     }
1592
1593   for (i = 0 ; i < factor->n_vars; ++i)
1594     {
1595       int j;
1596       const int matrix_row = perm->data[i];
1597       tab_text (t, 0, i + heading_rows, TAT_TITLE, var_to_string (factor->vars[matrix_row]));
1598
1599       for (j = 0 ; j < n_factors; ++j)
1600         {
1601           double x = gsl_matrix_get (fm, matrix_row, j);
1602
1603           if ( fabs (x) < factor->blank)
1604             continue;
1605
1606           tab_double (t, heading_columns + j, heading_rows + i, 0, x, NULL, RC_OTHER);
1607         }
1608     }
1609
1610   gsl_permutation_free (perm);
1611
1612   tab_submit (t);
1613 }
1614
1615
1616 static void
1617 show_explained_variance (const struct cmd_factor * factor, struct idata *idata,
1618                          const gsl_vector *initial_eigenvalues,
1619                          const gsl_vector *extracted_eigenvalues,
1620                          const gsl_vector *rotated_loadings)
1621 {
1622   size_t i;
1623   int c = 0;
1624   const int heading_columns = 1;
1625   const int heading_rows = 2;
1626   const int nr = heading_rows + factor->n_vars;
1627
1628   struct tab_table *t ;
1629
1630   double i_total = 0.0;
1631   double i_cum = 0.0;
1632
1633   double e_total = 0.0;
1634   double e_cum = 0.0;
1635
1636   double r_cum = 0.0;
1637
1638   int nc = heading_columns;
1639
1640   if (factor->print & PRINT_EXTRACTION)
1641     nc += 3;
1642
1643   if (factor->print & PRINT_INITIAL)
1644     nc += 3;
1645
1646   if (factor->print & PRINT_ROTATION)
1647     {
1648       nc += factor->rotation == ROT_PROMAX ? 1 : 3;
1649     }
1650
1651   /* No point having a table with only headings */
1652   if ( nc <= heading_columns)
1653     return;
1654
1655   t = tab_create (nc, nr);
1656
1657   tab_title (t, _("Total Variance Explained"));
1658
1659   tab_headers (t, heading_columns, 0, heading_rows, 0);
1660
1661   /* Outline the box */
1662   tab_box (t,
1663            TAL_2, TAL_2,
1664            -1, -1,
1665            0, 0,
1666            nc - 1, nr - 1);
1667
1668   /* Vertical lines */
1669   tab_box (t,
1670            -1, -1,
1671            -1, TAL_1,
1672            heading_columns, 0,
1673            nc - 1, nr - 1);
1674
1675   tab_hline (t, TAL_1, 0, nc - 1, heading_rows);
1676   tab_hline (t, TAL_1, 1, nc - 1, 1);
1677
1678   tab_vline (t, TAL_2, heading_columns, 0, nr - 1);
1679
1680
1681   if ( factor->extraction == EXTRACTION_PC)
1682     tab_text (t, 0, 1, TAB_LEFT | TAT_TITLE, _("Component"));
1683   else
1684     tab_text (t, 0, 1, TAB_LEFT | TAT_TITLE, _("Factor"));
1685
1686   c = 1;
1687   if (factor->print & PRINT_INITIAL)
1688     {
1689       tab_joint_text (t, c, 0, c + 2, 0, TAB_CENTER | TAT_TITLE, _("Initial Eigenvalues"));
1690       c += 3;
1691     }
1692
1693   if (factor->print & PRINT_EXTRACTION)
1694     {
1695       tab_joint_text (t, c, 0, c + 2, 0, TAB_CENTER | TAT_TITLE, _("Extraction Sums of Squared Loadings"));
1696       c += 3;
1697     }
1698
1699   if (factor->print & PRINT_ROTATION)
1700     {
1701       const int width = factor->rotation == ROT_PROMAX ? 0 : 2;
1702       tab_joint_text (t, c, 0, c + width, 0, TAB_CENTER | TAT_TITLE, _("Rotation Sums of Squared Loadings"));
1703       c += width + 1;
1704     }
1705
1706   for (i = 0; i < (nc - heading_columns + 2) / 3 ; ++i)
1707     {
1708       tab_text (t, i * 3 + 1, 1, TAB_CENTER | TAT_TITLE, _("Total"));
1709
1710       tab_vline (t, TAL_2, heading_columns + i * 3, 0, nr - 1);
1711
1712       if (i == 2 && factor->rotation == ROT_PROMAX)
1713         continue;
1714
1715       /* xgettext:no-c-format */
1716       tab_text (t, i * 3 + 2, 1, TAB_CENTER | TAT_TITLE, _("% of Variance"));
1717       tab_text (t, i * 3 + 3, 1, TAB_CENTER | TAT_TITLE, _("Cumulative %"));
1718     }
1719
1720   for (i = 0 ; i < initial_eigenvalues->size; ++i)
1721     i_total += gsl_vector_get (initial_eigenvalues, i);
1722
1723   if ( factor->extraction == EXTRACTION_PAF)
1724     {
1725       e_total = factor->n_vars;
1726     }
1727   else
1728     {
1729       e_total = i_total;
1730     }
1731
1732   for (i = 0 ; i < factor->n_vars; ++i)
1733     {
1734       const double i_lambda = gsl_vector_get (initial_eigenvalues, i);
1735       double i_percent = 100.0 * i_lambda / i_total ;
1736
1737       const double e_lambda = gsl_vector_get (extracted_eigenvalues, i);
1738       double e_percent = 100.0 * e_lambda / e_total ;
1739
1740       c = 0;
1741
1742       tab_text_format (t, c++, i + heading_rows, TAB_LEFT | TAT_TITLE, _("%zu"), i + 1);
1743
1744       i_cum += i_percent;
1745       e_cum += e_percent;
1746
1747       /* Initial Eigenvalues */
1748       if (factor->print & PRINT_INITIAL)
1749       {
1750         tab_double (t, c++, i + heading_rows, 0, i_lambda, NULL, RC_OTHER);
1751         tab_double (t, c++, i + heading_rows, 0, i_percent, NULL, RC_OTHER);
1752         tab_double (t, c++, i + heading_rows, 0, i_cum, NULL, RC_OTHER);
1753       }
1754
1755
1756       if (factor->print & PRINT_EXTRACTION)
1757         {
1758           if (i < idata->n_extractions)
1759             {
1760               /* Sums of squared loadings */
1761               tab_double (t, c++, i + heading_rows, 0, e_lambda, NULL, RC_OTHER);
1762               tab_double (t, c++, i + heading_rows, 0, e_percent, NULL, RC_OTHER);
1763               tab_double (t, c++, i + heading_rows, 0, e_cum, NULL, RC_OTHER);
1764             }
1765         }
1766
1767       if (rotated_loadings != NULL)
1768         {
1769           const double r_lambda = gsl_vector_get (rotated_loadings, i);
1770           double r_percent = 100.0 * r_lambda / e_total ;
1771
1772           if (factor->print & PRINT_ROTATION)
1773             {
1774               if (i < idata->n_extractions)
1775                 {
1776                   r_cum += r_percent;
1777                   tab_double (t, c++, i + heading_rows, 0, r_lambda, NULL, RC_OTHER);
1778                   if (factor->rotation != ROT_PROMAX)
1779                     {
1780                       tab_double (t, c++, i + heading_rows, 0, r_percent, NULL, RC_OTHER);
1781                       tab_double (t, c++, i + heading_rows, 0, r_cum, NULL, RC_OTHER);
1782                     }
1783                 }
1784             }
1785         }
1786     }
1787
1788   tab_submit (t);
1789 }
1790
1791
1792 static void
1793 show_factor_correlation (const struct cmd_factor * factor, const gsl_matrix *fcm)
1794 {
1795   size_t i, j;
1796   const int heading_columns = 1;
1797   const int heading_rows = 1;
1798   const int nr = heading_rows + fcm->size2;
1799   const int nc = heading_columns + fcm->size1;
1800   struct tab_table *t = tab_create (nc, nr);
1801
1802   tab_title (t, _("Factor Correlation Matrix"));
1803
1804   tab_headers (t, heading_columns, 0, heading_rows, 0);
1805
1806   /* Outline the box */
1807   tab_box (t,
1808            TAL_2, TAL_2,
1809            -1, -1,
1810            0, 0,
1811            nc - 1, nr - 1);
1812
1813   /* Vertical lines */
1814   tab_box (t,
1815            -1, -1,
1816            -1, TAL_1,
1817            heading_columns, 0,
1818            nc - 1, nr - 1);
1819
1820   tab_hline (t, TAL_1, 0, nc - 1, heading_rows);
1821   tab_hline (t, TAL_1, 1, nc - 1, 1);
1822
1823   tab_vline (t, TAL_2, heading_columns, 0, nr - 1);
1824
1825
1826   if ( factor->extraction == EXTRACTION_PC)
1827     tab_text (t, 0, 0, TAB_LEFT | TAT_TITLE, _("Component"));
1828   else
1829     tab_text (t, 0, 0, TAB_LEFT | TAT_TITLE, _("Factor"));
1830
1831   for (i = 0 ; i < fcm->size1; ++i)
1832     {
1833       tab_text_format (t, heading_columns + i, 0, TAB_CENTER | TAT_TITLE, _("%d"), i + 1);
1834     }
1835
1836   for (i = 0 ; i < fcm->size2; ++i)
1837     {
1838       tab_text_format (t, 0, heading_rows + i, TAB_CENTER | TAT_TITLE, _("%d"), i + 1);
1839     }
1840
1841
1842   for (i = 0 ; i < fcm->size1; ++i)
1843     {
1844       for (j = 0 ; j < fcm->size2; ++j)
1845         tab_double (t, heading_columns + i,  heading_rows +j, 0, 
1846                     gsl_matrix_get (fcm, i, j), NULL, RC_OTHER);
1847     }
1848
1849   tab_submit (t);
1850 }
1851
1852
1853 static void
1854 show_correlation_matrix (const struct cmd_factor *factor, const struct idata *idata)
1855 {
1856   struct tab_table *t ;
1857   size_t i, j;
1858   int y_pos_corr = -1;
1859   int y_pos_sig = -1;
1860   int suffix_rows = 0;
1861
1862   const int heading_rows = 1;
1863   const int heading_columns = 2;
1864
1865   int nc = heading_columns ;
1866   int nr = heading_rows ;
1867   int n_data_sets = 0;
1868
1869   if (factor->print & PRINT_CORRELATION)
1870     {
1871       y_pos_corr = n_data_sets;
1872       n_data_sets++;
1873       nc = heading_columns + factor->n_vars;
1874     }
1875
1876   if (factor->print & PRINT_SIG)
1877     {
1878       y_pos_sig = n_data_sets;
1879       n_data_sets++;
1880       nc = heading_columns + factor->n_vars;
1881     }
1882
1883   nr += n_data_sets * factor->n_vars;
1884
1885   if (factor->print & PRINT_DETERMINANT)
1886     suffix_rows = 1;
1887
1888   /* If the table would contain only headings, don't bother rendering it */
1889   if (nr <= heading_rows && suffix_rows == 0)
1890     return;
1891
1892   t = tab_create (nc, nr + suffix_rows);
1893
1894   tab_title (t, _("Correlation Matrix"));
1895
1896   tab_hline (t, TAL_1, 0, nc - 1, heading_rows);
1897
1898   if (nr > heading_rows)
1899     {
1900       tab_headers (t, heading_columns, 0, heading_rows, 0);
1901
1902       tab_vline (t, TAL_2, 2, 0, nr - 1);
1903
1904       /* Outline the box */
1905       tab_box (t,
1906                TAL_2, TAL_2,
1907                -1, -1,
1908                0, 0,
1909                nc - 1, nr - 1);
1910
1911       /* Vertical lines */
1912       tab_box (t,
1913                -1, -1,
1914                -1, TAL_1,
1915                heading_columns, 0,
1916                nc - 1, nr - 1);
1917
1918
1919       for (i = 0; i < factor->n_vars; ++i)
1920         tab_text (t, heading_columns + i, 0, TAT_TITLE, var_to_string (factor->vars[i]));
1921
1922
1923       for (i = 0 ; i < n_data_sets; ++i)
1924         {
1925           int y = heading_rows + i * factor->n_vars;
1926           size_t v;
1927           for (v = 0; v < factor->n_vars; ++v)
1928             tab_text (t, 1, y + v, TAT_TITLE, var_to_string (factor->vars[v]));
1929
1930           tab_hline (t, TAL_1, 0, nc - 1, y);
1931         }
1932
1933       if (factor->print & PRINT_CORRELATION)
1934         {
1935           const double y = heading_rows + y_pos_corr;
1936           tab_text (t, 0, y, TAT_TITLE, _("Correlations"));
1937
1938           for (i = 0; i < factor->n_vars; ++i)
1939             {
1940               for (j = 0; j < factor->n_vars; ++j)
1941                 tab_double (t, heading_columns + i,  y + j, 0, gsl_matrix_get (idata->corr, i, j), NULL, RC_OTHER);
1942             }
1943         }
1944
1945       if (factor->print & PRINT_SIG)
1946         {
1947           const double y = heading_rows + y_pos_sig * factor->n_vars;
1948           tab_text (t, 0, y, TAT_TITLE, _("Sig. (1-tailed)"));
1949
1950           for (i = 0; i < factor->n_vars; ++i)
1951             {
1952               for (j = 0; j < factor->n_vars; ++j)
1953                 {
1954                   double rho = gsl_matrix_get (idata->corr, i, j);
1955                   double w = gsl_matrix_get (idata->n, i, j);
1956
1957                   if (i == j)
1958                     continue;
1959
1960                   tab_double (t, heading_columns + i,  y + j, 0, significance_of_correlation (rho, w), NULL, RC_PVALUE);
1961                 }
1962             }
1963         }
1964     }
1965
1966   if (factor->print & PRINT_DETERMINANT)
1967     {
1968       tab_text (t, 0, nr, TAB_LEFT | TAT_TITLE, _("Determinant"));
1969
1970       tab_double (t, 1, nr, 0, idata->detR, NULL, RC_OTHER);
1971     }
1972
1973   tab_submit (t);
1974 }
1975
1976
1977
1978 static void
1979 do_factor (const struct cmd_factor *factor, struct casereader *r)
1980 {
1981   struct ccase *c;
1982   const gsl_matrix *var_matrix;
1983   const gsl_matrix *mean_matrix;
1984
1985   const gsl_matrix *analysis_matrix;
1986   struct idata *idata = idata_alloc (factor->n_vars);
1987
1988   struct covariance *cov = covariance_1pass_create (factor->n_vars, factor->vars,
1989                                               factor->wv, factor->exclude);
1990
1991   for ( ; (c = casereader_read (r) ); case_unref (c))
1992     {
1993       covariance_accumulate (cov, c);
1994     }
1995
1996   idata->cov = covariance_calculate (cov);
1997
1998   if (idata->cov == NULL)
1999     {
2000       msg (MW, _("The dataset contains no complete observations. No analysis will be performed."));
2001       covariance_destroy (cov);
2002       goto finish;
2003     }
2004
2005   var_matrix = covariance_moments (cov, MOMENT_VARIANCE);
2006   mean_matrix = covariance_moments (cov, MOMENT_MEAN);
2007   idata->n = covariance_moments (cov, MOMENT_NONE);
2008   
2009
2010   if ( factor->method == METHOD_CORR)
2011     {
2012       idata->corr = correlation_from_covariance (idata->cov, var_matrix);
2013       
2014       analysis_matrix = idata->corr;
2015     }
2016   else
2017     analysis_matrix = idata->cov;
2018
2019
2020   if (factor->print & PRINT_DETERMINANT
2021       || factor->print & PRINT_KMO)
2022     {
2023       int sign = 0;
2024
2025       const int size = idata->corr->size1;
2026       gsl_permutation *p = gsl_permutation_calloc (size);
2027       gsl_matrix *tmp = gsl_matrix_calloc (size, size);
2028       gsl_matrix_memcpy (tmp, idata->corr);
2029
2030       gsl_linalg_LU_decomp (tmp, p, &sign);
2031       idata->detR = gsl_linalg_LU_det (tmp, sign);
2032       gsl_permutation_free (p);
2033       gsl_matrix_free (tmp);
2034     }
2035
2036   if ( factor->print & PRINT_UNIVARIATE)
2037     {
2038       const struct fmt_spec *wfmt = factor->wv ? var_get_print_format (factor->wv) : & F_8_0;
2039       const int nc = 4;
2040       int i;
2041
2042       const int heading_columns = 1;
2043       const int heading_rows = 1;
2044
2045       const int nr = heading_rows + factor->n_vars;
2046
2047       struct tab_table *t = tab_create (nc, nr);
2048       tab_set_format (t, RC_WEIGHT, wfmt);
2049       tab_title (t, _("Descriptive Statistics"));
2050
2051       tab_headers (t, heading_columns, 0, heading_rows, 0);
2052
2053       /* Outline the box */
2054       tab_box (t,
2055                TAL_2, TAL_2,
2056                -1, -1,
2057                0, 0,
2058                nc - 1, nr - 1);
2059
2060       /* Vertical lines */
2061       tab_box (t,
2062                -1, -1,
2063                -1, TAL_1,
2064                heading_columns, 0,
2065                nc - 1, nr - 1);
2066
2067       tab_hline (t, TAL_1, 0, nc - 1, heading_rows);
2068       tab_vline (t, TAL_2, heading_columns, 0, nr - 1);
2069
2070       tab_text (t, 1, 0, TAB_CENTER | TAT_TITLE, _("Mean"));
2071       tab_text (t, 2, 0, TAB_CENTER | TAT_TITLE, _("Std. Deviation"));
2072       tab_text (t, 3, 0, TAB_CENTER | TAT_TITLE, _("Analysis N"));
2073
2074       for (i = 0 ; i < factor->n_vars; ++i)
2075         {
2076           const struct variable *v = factor->vars[i];
2077           tab_text (t, 0, i + heading_rows, TAB_LEFT | TAT_TITLE, var_to_string (v));
2078
2079           tab_double (t, 1, i + heading_rows, 0, gsl_matrix_get (mean_matrix, i, i), NULL, RC_OTHER);
2080           tab_double (t, 2, i + heading_rows, 0, sqrt (gsl_matrix_get (var_matrix, i, i)), NULL, RC_OTHER);
2081           tab_double (t, 3, i + heading_rows, 0, gsl_matrix_get (idata->n, i, i), NULL, RC_WEIGHT);
2082         }
2083
2084       tab_submit (t);
2085     }
2086
2087   if (factor->print & PRINT_KMO)
2088     {
2089       int i;
2090       double sum_ssq_r = 0;
2091       double sum_ssq_a = 0;
2092
2093       double df = factor->n_vars * ( factor->n_vars - 1) / 2;
2094
2095       double w = 0;
2096
2097
2098       double xsq;
2099
2100       const int heading_columns = 2;
2101       const int heading_rows = 0;
2102
2103       const int nr = heading_rows + 4;
2104       const int nc = heading_columns + 1;
2105
2106       gsl_matrix *a, *x;
2107
2108       struct tab_table *t = tab_create (nc, nr);
2109       tab_title (t, _("KMO and Bartlett's Test"));
2110
2111       x  = clone_matrix (idata->corr);
2112       gsl_linalg_cholesky_decomp (x);
2113       gsl_linalg_cholesky_invert (x);
2114
2115       a = anti_image (x);
2116
2117       for (i = 0; i < x->size1; ++i)
2118         {
2119           sum_ssq_r += ssq_od_n (x, i);
2120           sum_ssq_a += ssq_od_n (a, i);
2121         }
2122
2123       gsl_matrix_free (a);
2124       gsl_matrix_free (x);
2125
2126       tab_headers (t, heading_columns, 0, heading_rows, 0);
2127
2128       /* Outline the box */
2129       tab_box (t,
2130                TAL_2, TAL_2,
2131                -1, -1,
2132                0, 0,
2133                nc - 1, nr - 1);
2134
2135       tab_vline (t, TAL_2, heading_columns, 0, nr - 1);
2136
2137       tab_text (t, 0, 0, TAT_TITLE | TAB_LEFT, _("Kaiser-Meyer-Olkin Measure of Sampling Adequacy"));
2138
2139       tab_double (t, 2, 0, 0, sum_ssq_r /  (sum_ssq_r + sum_ssq_a), NULL, RC_OTHER);
2140
2141       tab_text (t, 0, 1, TAT_TITLE | TAB_LEFT, _("Bartlett's Test of Sphericity"));
2142
2143       tab_text (t, 1, 1, TAT_TITLE, _("Approx. Chi-Square"));
2144       tab_text (t, 1, 2, TAT_TITLE, _("df"));
2145       tab_text (t, 1, 3, TAT_TITLE, _("Sig."));
2146
2147
2148       /* The literature doesn't say what to do for the value of W when 
2149          missing values are involved.  The best thing I can think of
2150          is to take the mean average. */
2151       w = 0;
2152       for (i = 0; i < idata->n->size1; ++i)
2153         w += gsl_matrix_get (idata->n, i, i);
2154       w /= idata->n->size1;
2155
2156       xsq = w - 1 - (2 * factor->n_vars + 5) / 6.0;
2157       xsq *= -log (idata->detR);
2158
2159       tab_double (t, 2, 1, 0, xsq, NULL, RC_OTHER);
2160       tab_double (t, 2, 2, 0, df, NULL, RC_INTEGER);
2161       tab_double (t, 2, 3, 0, gsl_cdf_chisq_Q (xsq, df), NULL, RC_PVALUE);
2162       
2163
2164       tab_submit (t);
2165     }
2166
2167   show_correlation_matrix (factor, idata);
2168   covariance_destroy (cov);
2169
2170   {
2171     gsl_matrix *am = matrix_dup (analysis_matrix);
2172     gsl_eigen_symmv_workspace *workspace = gsl_eigen_symmv_alloc (factor->n_vars);
2173     
2174     gsl_eigen_symmv (am, idata->eval, idata->evec, workspace);
2175
2176     gsl_eigen_symmv_free (workspace);
2177     gsl_matrix_free (am);
2178   }
2179
2180   gsl_eigen_symmv_sort (idata->eval, idata->evec, GSL_EIGEN_SORT_ABS_DESC);
2181
2182   idata->n_extractions = n_extracted_factors (factor, idata);
2183
2184   if (idata->n_extractions == 0)
2185     {
2186       msg (MW, _("The %s criteria result in zero factors extracted. Therefore no analysis will be performed."), "FACTOR");
2187       goto finish;
2188     }
2189
2190   if (idata->n_extractions > factor->n_vars)
2191     {
2192       msg (MW, 
2193            _("The %s criteria result in more factors than variables, which is not meaningful. No analysis will be performed."), 
2194            "FACTOR");
2195       goto finish;
2196     }
2197     
2198   {
2199     gsl_matrix *rotated_factors = NULL;
2200     gsl_matrix *pattern_matrix = NULL;
2201     gsl_matrix *fcm = NULL;
2202     gsl_vector *rotated_loadings = NULL;
2203
2204     const gsl_vector *extracted_eigenvalues = NULL;
2205     gsl_vector *initial_communalities = gsl_vector_alloc (factor->n_vars);
2206     gsl_vector *extracted_communalities = gsl_vector_alloc (factor->n_vars);
2207     size_t i;
2208     struct factor_matrix_workspace *fmw = factor_matrix_workspace_alloc (idata->msr->size, idata->n_extractions);
2209     gsl_matrix *factor_matrix = gsl_matrix_calloc (factor->n_vars, fmw->n_factors);
2210
2211     if ( factor->extraction == EXTRACTION_PAF)
2212       {
2213         gsl_vector *diff = gsl_vector_alloc (idata->msr->size);
2214         struct smr_workspace *ws = ws_create (analysis_matrix);
2215
2216         for (i = 0 ; i < factor->n_vars ; ++i)
2217           {
2218             double r2 = squared_multiple_correlation (analysis_matrix, i, ws);
2219
2220             gsl_vector_set (idata->msr, i, r2);
2221           }
2222         ws_destroy (ws);
2223
2224         gsl_vector_memcpy (initial_communalities, idata->msr);
2225
2226         for (i = 0; i < factor->extraction_iterations; ++i)
2227           {
2228             double min, max;
2229             gsl_vector_memcpy (diff, idata->msr);
2230
2231             iterate_factor_matrix (analysis_matrix, idata->msr, factor_matrix, fmw);
2232       
2233             gsl_vector_sub (diff, idata->msr);
2234
2235             gsl_vector_minmax (diff, &min, &max);
2236       
2237             if ( fabs (min) < factor->econverge && fabs (max) < factor->econverge)
2238               break;
2239           }
2240         gsl_vector_free (diff);
2241
2242
2243
2244         gsl_vector_memcpy (extracted_communalities, idata->msr);
2245         extracted_eigenvalues = fmw->eval;
2246       }
2247     else if (factor->extraction == EXTRACTION_PC)
2248       {
2249         for (i = 0; i < factor->n_vars; ++i)
2250           gsl_vector_set (initial_communalities, i, communality (idata, i, factor->n_vars));
2251
2252         gsl_vector_memcpy (extracted_communalities, initial_communalities);
2253
2254         iterate_factor_matrix (analysis_matrix, extracted_communalities, factor_matrix, fmw);
2255
2256
2257         extracted_eigenvalues = idata->eval;
2258       }
2259
2260
2261     show_communalities (factor, initial_communalities, extracted_communalities);
2262
2263
2264     if ( factor->rotation != ROT_NONE)
2265       {
2266         rotated_factors = gsl_matrix_calloc (factor_matrix->size1, factor_matrix->size2);
2267         rotated_loadings = gsl_vector_calloc (factor_matrix->size2);
2268         if (factor->rotation == ROT_PROMAX)
2269           {
2270             pattern_matrix = gsl_matrix_calloc (factor_matrix->size1, factor_matrix->size2);
2271             fcm = gsl_matrix_calloc (factor_matrix->size2, factor_matrix->size2);
2272           }
2273           
2274
2275         rotate (factor, factor_matrix, extracted_communalities, rotated_factors, rotated_loadings, pattern_matrix, fcm);
2276       }
2277     
2278     show_explained_variance (factor, idata, idata->eval, extracted_eigenvalues, rotated_loadings);
2279
2280     factor_matrix_workspace_free (fmw);
2281
2282     show_scree (factor, idata);
2283
2284     show_factor_matrix (factor, idata,
2285                         factor->extraction == EXTRACTION_PC ? _("Component Matrix") : _("Factor Matrix"),
2286                         factor_matrix);
2287
2288     if ( factor->rotation == ROT_PROMAX)
2289       {
2290         show_factor_matrix (factor, idata, _("Pattern Matrix"),  pattern_matrix);
2291         gsl_matrix_free (pattern_matrix);
2292       }
2293
2294     if ( factor->rotation != ROT_NONE)
2295       {
2296         show_factor_matrix (factor, idata,
2297                             (factor->rotation == ROT_PROMAX) ? _("Structure Matrix") :
2298                             (factor->extraction == EXTRACTION_PC ? _("Rotated Component Matrix") : _("Rotated Factor Matrix")),
2299                             rotated_factors);
2300
2301         gsl_matrix_free (rotated_factors);
2302       }
2303
2304     if ( factor->rotation == ROT_PROMAX)
2305       {
2306         show_factor_correlation (factor, fcm);
2307         gsl_matrix_free (fcm);
2308       }
2309
2310     gsl_matrix_free (factor_matrix);
2311     gsl_vector_free (rotated_loadings);
2312     gsl_vector_free (initial_communalities);
2313     gsl_vector_free (extracted_communalities);
2314   }
2315
2316  finish:
2317
2318   idata_free (idata);
2319
2320   casereader_destroy (r);
2321 }
2322
2323