Added symmetric poweriterations to sandbox
This commit is contained in:
27
sandbox/powerit/Makefile
Normal file
27
sandbox/powerit/Makefile
Normal file
@@ -0,0 +1,27 @@
|
||||
CC=gcc
|
||||
PYTHON_INCLUDE=/usr/include/python2.5
|
||||
PYTHON_LIBRARY=/usr/lib/python2.5
|
||||
NUMPY_INCLUDE=$(PYTHON_LIBRARY)/site-packages/numpy/core/include/numpy/
|
||||
GOTO_LIBRARY=/home/flatberg/GotoBLAS
|
||||
BLAS_LIBRARY=/usr/lib
|
||||
BLAS=blas
|
||||
GOTO=goto
|
||||
|
||||
|
||||
all: numpy_module.so
|
||||
|
||||
numpy_module.so: numpy_module.o sympowerit.o
|
||||
$(CC) -Wall -shared numpy_module.o sympowerit.o -o _numpy_module.so -L$(PYTHON_LIBRARY) -lpython2.5 -L$(BLAS_LIBRARY) -l$(BLAS) -llapack -lg2c -lm
|
||||
|
||||
numpy_module.o: numpy_module.c numpy_module.h
|
||||
$(CC) -fPIC -Wall -O2 -g -c -I$(PYTHON_INCLUDE) -I$(NUMPY_INCLUDE) numpy_module.c
|
||||
|
||||
c_egines.o: sympowerit.c
|
||||
$(CC) -Wall -O2 -g -c sympowerit.c
|
||||
|
||||
clean:
|
||||
-rm numpy_module.o sympowerit.o
|
||||
-rm -rf _numpy_module.so
|
||||
|
||||
|
||||
|
70
sandbox/powerit/numpy_module.c
Executable file
70
sandbox/powerit/numpy_module.c
Executable file
@@ -0,0 +1,70 @@
|
||||
/* A file to test imorting C modules for handling arrays to Python */
|
||||
|
||||
/* Python.h includes <stdio.h>, <string.h>, <errno.h>, <limits.h>, and <stdlib.h> (if available)*/
|
||||
#include "Python.h"
|
||||
#include "arrayobject.h"
|
||||
#include "numpy_module.h"
|
||||
#include "sympowerit.h"
|
||||
#include <cblas.h>
|
||||
|
||||
|
||||
/* ==== Set up the methods table ====================== */
|
||||
static PyMethodDef _numpy_moduleMethods[] = {
|
||||
{"sym_powerit", sym_powerit, METH_VARARGS},
|
||||
{NULL, NULL} /* Sentinel - marks the end of this structure */
|
||||
};
|
||||
|
||||
|
||||
/* ==== Initialize the numpy_module functions =============== */
|
||||
// Module name must be _numpy_module in compile and linked
|
||||
void init_numpy_module() {
|
||||
(void) Py_InitModule("_numpy_module", _numpy_moduleMethods);
|
||||
import_array(); // Must be present for NumPy. Called first after above line.
|
||||
}
|
||||
|
||||
|
||||
/* =========== Power iteration on symmetric matrix ========== */
|
||||
static PyObject *sym_powerit(PyObject *self, PyObject *args)
|
||||
{
|
||||
PyArrayObject *X=NULL, *T=NULL, *E=NULL,*U=NULL;
|
||||
double *tin, *ein, tol;
|
||||
int amax, n, info, maxiter, verbose;
|
||||
int dims[2];
|
||||
|
||||
/* Parse tuple of input arguments*/
|
||||
if (!PyArg_ParseTuple(args, "O!O!iidii",
|
||||
&PyArray_Type, &X,
|
||||
&PyArray_Type, &T,
|
||||
&n,
|
||||
&amax,
|
||||
&tol,
|
||||
&maxiter,
|
||||
&verbose)
|
||||
)
|
||||
return NULL;
|
||||
if (NULL == X) return NULL;
|
||||
|
||||
/* Get the dimensions of the input */
|
||||
n = X->dimensions[0];
|
||||
|
||||
/* Create output/ work arrays, no inplace calculations */
|
||||
dims[0] = n;
|
||||
dims[1] = amax;
|
||||
U = (PyArrayObject*) PyArray_NewCopy(T, NPY_CORDER);
|
||||
dims[1] = n;
|
||||
E = (PyArrayObject*) PyArray_NewCopy(X, NPY_CORDER);
|
||||
|
||||
/* Get pointers to contigous data (row-major)*/
|
||||
ein = (double *)E->data;
|
||||
tin = (double *)U->data;
|
||||
|
||||
/* Call sympower method */
|
||||
info = sympowerit (ein, n, tin, amax, tol, maxiter, verbose);
|
||||
|
||||
Py_DECREF(E);
|
||||
|
||||
return Py_BuildValue("Ni", U, info);
|
||||
|
||||
}
|
||||
|
||||
|
1
sandbox/powerit/numpy_module.h
Executable file
1
sandbox/powerit/numpy_module.h
Executable file
@@ -0,0 +1 @@
|
||||
/* Header to test of C modules for arrays for Python: sympowerit.c */
|
116
sandbox/powerit/sympowerit.c
Normal file
116
sandbox/powerit/sympowerit.c
Normal file
@@ -0,0 +1,116 @@
|
||||
#include <cblas.h>
|
||||
#include <math.h>
|
||||
#include <stdio.h>
|
||||
#include <stdlib.h>
|
||||
|
||||
/* =============== Main functions body ====================*/
|
||||
|
||||
int sympowerit(double *E, int n, double *T, int amax, double tol,
|
||||
int maxiter, int verbose)
|
||||
{
|
||||
/*
|
||||
PURPOSE: Estimate eigenvectos of a symmetric matrix using the power method.
|
||||
|
||||
CALLING SEQUENCE:
|
||||
info = sympowerit(*E, n, *T, amax, tol, maxiter);
|
||||
|
||||
INPUTS:
|
||||
E symmetric matrix from XX^T/X^TX of centered data matrix
|
||||
type: *double
|
||||
n elements in X first dimension
|
||||
type: int
|
||||
type: *double
|
||||
T workspace for scores (m, amax)
|
||||
type: *double
|
||||
amax maximum number of principal components
|
||||
type: int
|
||||
tol convergence limit
|
||||
type: double
|
||||
maxiter maximum number of iteraitons
|
||||
type: int
|
||||
verbose the amount of debug output
|
||||
*/
|
||||
|
||||
int iter, a, j;
|
||||
int info = 0;
|
||||
double sumsq, l2norm, *y,*x, sgn, diff, alpha;
|
||||
double lambda = 0.0;
|
||||
|
||||
if (verbose>1) verbose=1;
|
||||
|
||||
/* Allocate space for working vector and eigenvector */
|
||||
y = (double *) malloc(n*sizeof(double));
|
||||
x = (double *) malloc(n*sizeof(double));
|
||||
|
||||
/* Loop over all components to be estimated*/
|
||||
for (a=0; a<=amax-1; a++){
|
||||
/* Reset work-/eigen-vector*/
|
||||
for (j=0; j<n; j++) {
|
||||
y[j] = 0.0;
|
||||
x[j] = T[(j*amax)+a];
|
||||
}
|
||||
/*x[0] = 1.0;*/
|
||||
|
||||
/* Main power-iteration loop */
|
||||
for ( iter = 0; iter <= maxiter; iter++ ) {
|
||||
/* Matrix-vector product */
|
||||
/*cblas_dsymv (CblasRowMajor, CblasUpper, n, 1.0, E, n,
|
||||
x, 1, 0.0, y, 1);*/
|
||||
cblas_dgemv (CblasRowMajor, CblasNoTrans, n, n, 1.0, E, n,
|
||||
x, 1, 0.0, y, 1);
|
||||
|
||||
/* Normalise y */
|
||||
sumsq = y[0] * y[0];
|
||||
lambda = x[0] * y[0];
|
||||
for ( j = 1; j < n; j++ ) {
|
||||
sumsq += y[j] * y[j];
|
||||
lambda += x[j] * y[j];
|
||||
}
|
||||
l2norm = sqrt ( sumsq );
|
||||
for ( j = 0; j < n; j++ ) y[j] /= l2norm;
|
||||
|
||||
/*Check for convergence */
|
||||
sgn = ( lambda < 0 ? -1.0 : 1.0 );
|
||||
diff = x[0] - sgn * y[0];
|
||||
sumsq = diff * diff;
|
||||
x[0] = y[0];
|
||||
for ( j = 0; j < n; j++ ) {
|
||||
diff = x[j] - sgn * y[j];
|
||||
sumsq += diff * diff;
|
||||
x[j] = y[j];
|
||||
}
|
||||
if ( sqrt ( sumsq ) < tol ) {
|
||||
if (verbose == 1){
|
||||
printf("\nComp: %d\n", a);
|
||||
printf("Converged in %d iterations\n", iter);
|
||||
}
|
||||
break;
|
||||
}
|
||||
if (iter >= maxiter){
|
||||
if (verbose == 1){
|
||||
printf("\nComp: %d\n", a);
|
||||
printf("Max iter reached.\n");
|
||||
printf("Error: %.2E\n", sumsq);
|
||||
}
|
||||
info = 1;
|
||||
break;
|
||||
}
|
||||
}
|
||||
/* Calculate T */
|
||||
for (j=0; j<n; j++){
|
||||
y[j] = sgn*sqrt(sgn*lambda)*x[j];
|
||||
T[(j*amax)+a] = y[j];
|
||||
}
|
||||
|
||||
/* rank one deflation of residual matrix */
|
||||
/*cblas_dsyr (CblasRowMajor, CblasUpper, n, -1.0, x, 1, E, n);*/
|
||||
alpha = -1.0*lambda;
|
||||
cblas_dger (CblasRowMajor, n, n, alpha, x, 1, x, 1, E, n);
|
||||
}
|
||||
|
||||
/* Free used space */
|
||||
free(x);
|
||||
free(y);
|
||||
|
||||
return (info);
|
||||
}
|
2
sandbox/powerit/sympowerit.h
Normal file
2
sandbox/powerit/sympowerit.h
Normal file
@@ -0,0 +1,2 @@
|
||||
int sympowerit (double *X, int n, double *T, int amax, double tol, int maxiter, int verbose);
|
||||
/* Simple symmetric power iterations */
|
100
sandbox/powerit/sympowerit.py
Normal file
100
sandbox/powerit/sympowerit.py
Normal file
@@ -0,0 +1,100 @@
|
||||
from numpy import *
|
||||
|
||||
HAS_SYMPOWER=True
|
||||
try:
|
||||
from _numpy_module import sym_powerit
|
||||
except:
|
||||
raise ImportError("Sym_powerit module not present")
|
||||
HAS_SYMPOWER = False
|
||||
|
||||
class SymPowerException(Exception):
|
||||
pass
|
||||
|
||||
_ERRORCODES = {1: "Some eigenvectors did not converge, try to increase \nthe number of iterations or lower the tolerance level",
|
||||
0: ""}
|
||||
|
||||
|
||||
def sympowerit(xx, T0=None, mn_center=False, a_max=10, tol=1e-7, maxiter=100,
|
||||
verbose=0):
|
||||
"""Estimate eigenvectos of a symmetric matrix using the power method.
|
||||
|
||||
*Parameters*:
|
||||
|
||||
xx : {array}
|
||||
Symmetric square array (m, m)
|
||||
T0 : {array}
|
||||
Initial solution (m, a_max), optional
|
||||
mn_center : {boolean}, optional
|
||||
Mean centering
|
||||
a_max : {integer}, optional
|
||||
Number of components to extract
|
||||
tol : {float}, optional
|
||||
Tolerance level of eigenvector solver
|
||||
maxiter : {integer}
|
||||
Maximum number of poweriterations to use
|
||||
verbose : {integer}
|
||||
Debug output (==1)
|
||||
|
||||
*Returns*:
|
||||
v : {array}
|
||||
Eigenvectors of xx, (m , a_max)
|
||||
"""
|
||||
|
||||
valid_types = ['D','d','F','f']
|
||||
dtype = xx.dtype.char
|
||||
n, m = xx.shape
|
||||
if not(dtype in valid_types):
|
||||
msg = "Array type: (%s) needs to be a float or double" %dtype
|
||||
raise SymPowerException(msg)
|
||||
if not (m==n):
|
||||
msg = "Input array needs to be square, input: (%d,%d)" %(m,n)
|
||||
raise SymPowerException(msg)
|
||||
# small test of symmetry
|
||||
N = 5
|
||||
num = random.randint(0,n,N)
|
||||
for i in range(5):
|
||||
j = N-5
|
||||
if abs(xx[num[i],num[j]] - xx[num[j],num[i]])>1e-15:
|
||||
msg = "Array needs to be symmetric"
|
||||
raise SymPowerException(msg)
|
||||
|
||||
if not a_max:
|
||||
a_max = 10
|
||||
|
||||
if T0 !=None:
|
||||
tn, tm = T0.shape
|
||||
if not (tn==n):
|
||||
msg = "Start eigenvectors need to match input array ()"
|
||||
raise SymPowerException(msg)
|
||||
if not (tm==a_max):
|
||||
msg = "Start eigenvectors need to match input a_max ()"
|
||||
raise SymPowerException(msg)
|
||||
else:
|
||||
T0 = zeros((n, a_max), 'd')
|
||||
T0[0,:] = ones((a_max,),'d')
|
||||
|
||||
if mn_center:
|
||||
xx = _center(xx)
|
||||
|
||||
# call c-function
|
||||
T, info = sym_powerit(xx, T0, n, a_max, tol, maxiter, verbose)
|
||||
|
||||
if info != 0:
|
||||
if verbose:
|
||||
print _ERRORCODES.get(info, "Dont know this error")
|
||||
return T
|
||||
|
||||
|
||||
def _center(xx, ret_mn=False):
|
||||
"""Returns mean centered symmetric kernel matrix.
|
||||
"""
|
||||
n = xx.shape[0]
|
||||
h = xx.sum(0)[:,newaxis]
|
||||
h = (h - mean(h)/2)/n
|
||||
mn_a = h + h.T
|
||||
xxc = xx - mn_a
|
||||
if ret_mn:
|
||||
return xxc, mn_a
|
||||
return xxc
|
||||
|
||||
|
Reference in New Issue
Block a user