c@244: /*********************************/ c@244: /* Principal Components Analysis */ c@244: /*********************************/ c@244: c@244: /*********************************************************************/ c@244: /* Principal Components Analysis or the Karhunen-Loeve expansion is a c@244: classical method for dimensionality reduction or exploratory data c@244: analysis. One reference among many is: F. Murtagh and A. Heck, c@244: Multivariate Data Analysis, Kluwer Academic, Dordrecht, 1987. c@244: c@244: Author: c@244: F. Murtagh c@244: Phone: + 49 89 32006298 (work) c@244: + 49 89 965307 (home) c@244: Earn/Bitnet: fionn@dgaeso51, fim@dgaipp1s, murtagh@stsci c@244: Span: esomc1::fionn c@244: Internet: murtagh@scivax.stsci.edu c@244: c@244: F. Murtagh, Munich, 6 June 1989 */ c@244: /*********************************************************************/ c@244: c@244: #include c@244: #include c@244: #include c@244: c@244: #include "pca.h" c@244: c@244: #define SIGN(a, b) ( (b) < 0 ? -fabs(a) : fabs(a) ) c@244: c@244: /** Variance-covariance matrix: creation *****************************/ c@244: c@244: /* Create m * m covariance matrix from given n * m data matrix. */ c@244: void covcol(double** data, int n, int m, double** symmat) c@244: { c@244: double *mean; c@244: int i, j, j1, j2; c@244: c@244: /* Allocate storage for mean vector */ c@244: c@244: mean = (double*) malloc(m*sizeof(double)); c@244: c@244: /* Determine mean of column vectors of input data matrix */ c@244: c@244: for (j = 0; j < m; j++) c@244: { c@244: mean[j] = 0.0; c@244: for (i = 0; i < n; i++) c@244: { c@244: mean[j] += data[i][j]; c@244: } c@244: mean[j] /= (double)n; c@244: } c@244: c@244: /* c@244: printf("\nMeans of column vectors:\n"); c@244: for (j = 0; j < m; j++) { c@244: printf("%12.1f",mean[j]); } printf("\n"); c@244: */ c@244: c@244: /* Center the column vectors. */ c@244: c@244: for (i = 0; i < n; i++) c@244: { c@244: for (j = 0; j < m; j++) c@244: { c@244: data[i][j] -= mean[j]; c@244: } c@244: } c@244: c@244: /* Calculate the m * m covariance matrix. */ c@244: for (j1 = 0; j1 < m; j1++) c@244: { c@244: for (j2 = j1; j2 < m; j2++) c@244: { c@244: symmat[j1][j2] = 0.0; c@244: for (i = 0; i < n; i++) c@244: { c@244: symmat[j1][j2] += data[i][j1] * data[i][j2]; c@244: } c@244: symmat[j2][j1] = symmat[j1][j2]; c@244: } c@244: } c@244: c@244: free(mean); c@244: c@244: return; c@244: c@244: } c@244: c@244: /** Error handler **************************************************/ c@244: c@244: void erhand(char* err_msg) c@244: { c@244: fprintf(stderr,"Run-time error:\n"); c@244: fprintf(stderr,"%s\n", err_msg); c@244: fprintf(stderr,"Exiting to system.\n"); c@244: exit(1); c@244: } c@244: c@244: c@244: /** Reduce a real, symmetric matrix to a symmetric, tridiag. matrix. */ c@244: c@244: /* Householder reduction of matrix a to tridiagonal form. c@244: Algorithm: Martin et al., Num. Math. 11, 181-195, 1968. c@244: Ref: Smith et al., Matrix Eigensystem Routines -- EISPACK Guide c@244: Springer-Verlag, 1976, pp. 489-494. c@244: W H Press et al., Numerical Recipes in C, Cambridge U P, c@244: 1988, pp. 373-374. */ c@244: void tred2(double** a, int n, double* d, double* e) c@244: { c@244: int l, k, j, i; c@244: double scale, hh, h, g, f; c@244: c@244: for (i = n-1; i >= 1; i--) c@244: { c@244: l = i - 1; c@244: h = scale = 0.0; c@244: if (l > 0) c@244: { c@244: for (k = 0; k <= l; k++) c@244: scale += fabs(a[i][k]); c@244: if (scale == 0.0) c@244: e[i] = a[i][l]; c@244: else c@244: { c@244: for (k = 0; k <= l; k++) c@244: { c@244: a[i][k] /= scale; c@244: h += a[i][k] * a[i][k]; c@244: } c@244: f = a[i][l]; c@244: g = f>0 ? -sqrt(h) : sqrt(h); c@244: e[i] = scale * g; c@244: h -= f * g; c@244: a[i][l] = f - g; c@244: f = 0.0; c@244: for (j = 0; j <= l; j++) c@244: { c@244: a[j][i] = a[i][j]/h; c@244: g = 0.0; c@244: for (k = 0; k <= j; k++) c@244: g += a[j][k] * a[i][k]; c@244: for (k = j+1; k <= l; k++) c@244: g += a[k][j] * a[i][k]; c@244: e[j] = g / h; c@244: f += e[j] * a[i][j]; c@244: } c@244: hh = f / (h + h); c@244: for (j = 0; j <= l; j++) c@244: { c@244: f = a[i][j]; c@244: e[j] = g = e[j] - hh * f; c@244: for (k = 0; k <= j; k++) c@244: a[j][k] -= (f * e[k] + g * a[i][k]); c@244: } c@244: } c@244: } c@244: else c@244: e[i] = a[i][l]; c@244: d[i] = h; c@244: } c@244: d[0] = 0.0; c@244: e[0] = 0.0; c@244: for (i = 0; i < n; i++) c@244: { c@244: l = i - 1; c@244: if (d[i]) c@244: { c@244: for (j = 0; j <= l; j++) c@244: { c@244: g = 0.0; c@244: for (k = 0; k <= l; k++) c@244: g += a[i][k] * a[k][j]; c@244: for (k = 0; k <= l; k++) c@244: a[k][j] -= g * a[k][i]; c@244: } c@244: } c@244: d[i] = a[i][i]; c@244: a[i][i] = 1.0; c@244: for (j = 0; j <= l; j++) c@244: a[j][i] = a[i][j] = 0.0; c@244: } c@244: } c@244: c@244: /** Tridiagonal QL algorithm -- Implicit **********************/ c@244: c@244: void tqli(double* d, double* e, int n, double** z) c@244: { c@244: int m, l, iter, i, k; c@244: double s, r, p, g, f, dd, c, b; c@244: c@244: for (i = 1; i < n; i++) c@244: e[i-1] = e[i]; c@244: e[n-1] = 0.0; c@244: for (l = 0; l < n; l++) c@244: { c@244: iter = 0; c@244: do c@244: { c@244: for (m = l; m < n-1; m++) c@244: { c@244: dd = fabs(d[m]) + fabs(d[m+1]); c@244: if (fabs(e[m]) + dd == dd) break; c@244: } c@244: if (m != l) c@244: { c@244: if (iter++ == 30) erhand("No convergence in TLQI."); c@244: g = (d[l+1] - d[l]) / (2.0 * e[l]); c@244: r = sqrt((g * g) + 1.0); c@244: g = d[m] - d[l] + e[l] / (g + SIGN(r, g)); c@244: s = c = 1.0; c@244: p = 0.0; c@244: for (i = m-1; i >= l; i--) c@244: { c@244: f = s * e[i]; c@244: b = c * e[i]; c@244: if (fabs(f) >= fabs(g)) c@244: { c@244: c = g / f; c@244: r = sqrt((c * c) + 1.0); c@244: e[i+1] = f * r; c@244: c *= (s = 1.0/r); c@244: } c@244: else c@244: { c@244: s = f / g; c@244: r = sqrt((s * s) + 1.0); c@244: e[i+1] = g * r; c@244: s *= (c = 1.0/r); c@244: } c@244: g = d[i+1] - p; c@244: r = (d[i] - g) * s + 2.0 * c * b; c@244: p = s * r; c@244: d[i+1] = g + p; c@244: g = c * r - b; c@244: for (k = 0; k < n; k++) c@244: { c@244: f = z[k][i+1]; c@244: z[k][i+1] = s * z[k][i] + c * f; c@244: z[k][i] = c * z[k][i] - s * f; c@244: } c@244: } c@244: d[l] = d[l] - p; c@244: e[l] = g; c@244: e[m] = 0.0; c@244: } c@244: } while (m != l); c@244: } c@244: } c@244: c@244: /* In place projection onto basis vectors */ c@244: void pca_project(double** data, int n, int m, int ncomponents) c@244: { c@244: int i, j, k, k2; c@414: double **symmat, /* **symmat2, */ *evals, *interm; c@244: c@244: //TODO: assert ncomponents < m c@244: c@244: symmat = (double**) malloc(m*sizeof(double*)); c@244: for (i = 0; i < m; i++) c@244: symmat[i] = (double*) malloc(m*sizeof(double)); c@244: c@244: covcol(data, n, m, symmat); c@244: c@244: /********************************************************************* c@244: Eigen-reduction c@244: **********************************************************************/ c@244: c@244: /* Allocate storage for dummy and new vectors. */ c@244: evals = (double*) malloc(m*sizeof(double)); /* Storage alloc. for vector of eigenvalues */ c@244: interm = (double*) malloc(m*sizeof(double)); /* Storage alloc. for 'intermediate' vector */ c@244: //MALLOC_ARRAY(symmat2,m,m,double); c@244: //for (i = 0; i < m; i++) { c@244: // for (j = 0; j < m; j++) { c@244: // symmat2[i][j] = symmat[i][j]; /* Needed below for col. projections */ c@244: // } c@244: //} c@244: tred2(symmat, m, evals, interm); /* Triangular decomposition */ c@244: tqli(evals, interm, m, symmat); /* Reduction of sym. trid. matrix */ c@244: /* evals now contains the eigenvalues, c@244: columns of symmat now contain the associated eigenvectors. */ c@244: c@244: /* c@244: printf("\nEigenvalues:\n"); c@244: for (j = m-1; j >= 0; j--) { c@244: printf("%18.5f\n", evals[j]); } c@244: printf("\n(Eigenvalues should be strictly positive; limited\n"); c@244: printf("precision machine arithmetic may affect this.\n"); c@244: printf("Eigenvalues are often expressed as cumulative\n"); c@244: printf("percentages, representing the 'percentage variance\n"); c@244: printf("explained' by the associated axis or principal component.)\n"); c@244: c@244: printf("\nEigenvectors:\n"); c@244: printf("(First three; their definition in terms of original vbes.)\n"); c@244: for (j = 0; j < m; j++) { c@244: for (i = 1; i <= 3; i++) { c@244: printf("%12.4f", symmat[j][m-i]); } c@244: printf("\n"); } c@244: */ c@244: c@244: /* Form projections of row-points on prin. components. */ c@244: /* Store in 'data', overwriting original data. */ c@244: for (i = 0; i < n; i++) { c@244: for (j = 0; j < m; j++) { c@244: interm[j] = data[i][j]; } /* data[i][j] will be overwritten */ c@244: for (k = 0; k < ncomponents; k++) { c@244: data[i][k] = 0.0; c@244: for (k2 = 0; k2 < m; k2++) { c@244: data[i][k] += interm[k2] * symmat[k2][m-k-1]; } c@244: } c@244: } c@244: c@244: /* c@244: printf("\nProjections of row-points on first 3 prin. comps.:\n"); c@244: for (i = 0; i < n; i++) { c@244: for (j = 0; j < 3; j++) { c@244: printf("%12.4f", data[i][j]); } c@244: printf("\n"); } c@244: */ c@244: c@244: /* Form projections of col.-points on first three prin. components. */ c@244: /* Store in 'symmat2', overwriting what was stored in this. */ c@244: //for (j = 0; j < m; j++) { c@244: // for (k = 0; k < m; k++) { c@244: // interm[k] = symmat2[j][k]; } /*symmat2[j][k] will be overwritten*/ c@244: // for (i = 0; i < 3; i++) { c@244: // symmat2[j][i] = 0.0; c@244: // for (k2 = 0; k2 < m; k2++) { c@244: // symmat2[j][i] += interm[k2] * symmat[k2][m-i-1]; } c@244: // if (evals[m-i-1] > 0.0005) /* Guard against zero eigenvalue */ c@244: // symmat2[j][i] /= sqrt(evals[m-i-1]); /* Rescale */ c@244: // else c@244: // symmat2[j][i] = 0.0; /* Standard kludge */ c@244: // } c@244: // } c@244: c@244: /* c@244: printf("\nProjections of column-points on first 3 prin. comps.:\n"); c@244: for (j = 0; j < m; j++) { c@244: for (k = 0; k < 3; k++) { c@244: printf("%12.4f", symmat2[j][k]); } c@244: printf("\n"); } c@244: */ c@244: c@244: c@244: for (i = 0; i < m; i++) c@244: free(symmat[i]); c@244: free(symmat); c@244: //FREE_ARRAY(symmat2,m); c@244: free(evals); c@244: free(interm); c@244: c@244: } c@244: c@244: c@244: