EMAN::SVDAnalyzer Class Reference

Singular Value Decomposition from GSL. More...

#include <analyzer.h>

Inheritance diagram for EMAN::SVDAnalyzer:

Inheritance graph
[legend]
Collaboration diagram for EMAN::SVDAnalyzer:

Collaboration graph
[legend]

List of all members.

Public Member Functions

 SVDAnalyzer ()
virtual int insert_image (EMData *image)
 insert a image to the list of input images
virtual int insert_images_list (vector< EMData * > image_list)
 insert a list of images to the list of input images
virtual vector< EMData * > analyze ()
 main function for Analyzer, analyze input images and create output images
string get_name () const
 Get the Analyzer's name.
string get_desc () const
 Get the Analyzer's description.
void set_params (const Dict &new_params)
 Set the Analyzer parameters using a key/value dictionary.
TypeDict get_param_types () const
 Get Analyzer parameter information in a dictionary.

Static Public Member Functions

static AnalyzerNEW ()

Protected Attributes

EMDatamask
int nvec
int pixels
int nimg

Private Attributes

int nsofar
gsl_matrix * A


Detailed Description

Singular Value Decomposition from GSL.

Comparable to pca

Parameters:
mask mask image
nvec number of desired basis vectors
nimg total number of input images, required even with insert_image()

Definition at line 202 of file analyzer.h.


Constructor & Destructor Documentation

EMAN::SVDAnalyzer::SVDAnalyzer (  )  [inline]

Definition at line 205 of file analyzer.h.

Referenced by NEW().

00205 : mask(0), nvec(0), nimg(0), A(NULL) {}


Member Function Documentation

int SVDAnalyzer::insert_image ( EMData image  )  [virtual]

insert a image to the list of input images

Parameters:
image 
Returns:
int 0 for success, <0 for fail

Implements EMAN::Analyzer.

Definition at line 728 of file analyzer.cpp.

References A, EMAN::EMData::get_data(), EMAN::EMData::get_xsize(), EMAN::EMData::get_ysize(), EMAN::EMData::get_zsize(), mask, nsofar, and NullPointerException.

00729 {
00730         if (mask==0)
00731                 throw NullPointerException("Null mask image pointer, set_params() first");
00732 
00733         // count pixels under mask
00734         size_t totpix=mask->get_xsize()*mask->get_ysize()*mask->get_zsize();
00735         float  *d=image->get_data();
00736         float *md=mask ->get_data();
00737         for (size_t i=0,j=0; i<totpix; ++i) {
00738                 if (md[i]) {
00739                         gsl_matrix_set(A,j,nsofar,d[i]);
00740                         j++;
00741                 }
00742         }
00743         nsofar++;
00744 
00745    return 0;
00746 }

virtual int EMAN::SVDAnalyzer::insert_images_list ( vector< EMData * >  image_list  )  [inline, virtual]

insert a list of images to the list of input images

Parameters:
image_list 
Returns:
int 0 for success, <0 for fail

Reimplemented from EMAN::Analyzer.

Definition at line 209 of file analyzer.h.

References EMAN::Analyzer::images.

00209                                                                             {
00210                         vector<EMData*>::const_iterator iter;
00211                         for(iter=image_list.begin(); iter!=image_list.end(); ++iter) {
00212                                 images.push_back(*iter);
00213                         }
00214                         return 0;
00215                 }

vector< EMData * > SVDAnalyzer::analyze (  )  [virtual]

main function for Analyzer, analyze input images and create output images

Returns:
vector<EMData *> result os images analysis

Implements EMAN::Analyzer.

Definition at line 750 of file analyzer.cpp.

References A, EMAN::EMData::get_data(), EMAN::EMData::get_xsize(), EMAN::EMData::get_ysize(), EMAN::EMData::get_zsize(), mask, nimg, nvec, EMAN::EMData::set_attr(), EMAN::EMData::set_size(), V, and X.

00751 {
00752 // Allocate the working space
00753 gsl_vector *work=gsl_vector_alloc(nimg);
00754 gsl_vector *S=gsl_vector_alloc(nimg);
00755 gsl_matrix *V=gsl_matrix_alloc(nimg,nimg);
00756 gsl_matrix *X=gsl_matrix_alloc(nimg,nimg);
00757 
00758 
00759 // Do the decomposition. All the real work is here
00760 gsl_linalg_SV_decomp_mod (A,X, V, S, work);
00761 //else gsl_linalg_SV_decomp_jacobi(A,V,S);
00762 
00763 vector<EMData*> ret;
00764 //unpack the results and write the output file
00765 float *md=mask->get_data();
00766 size_t totpix=mask->get_xsize()*mask->get_ysize()*mask->get_zsize();
00767 for (int k=0; k<nvec; k++) {
00768         EMData *img = new EMData;
00769         img->set_size(mask->get_xsize(),mask->get_ysize(),mask->get_zsize());
00770 
00771         float  *d=img->get_data();
00772         for (size_t i=0,j=0; i<totpix; ++i) {
00773                 if (md[i]) {
00774                         d[i]=(float)gsl_matrix_get(A,j,k);
00775                         j++;
00776                 }
00777         }
00778         img->set_attr( "eigval", gsl_vector_get(S,k));
00779         ret.push_back(img);
00780 }
00781 
00782 gsl_vector_free(work);
00783 gsl_vector_free(S);
00784 gsl_matrix_free(V);
00785 gsl_matrix_free(X);
00786 
00787 gsl_matrix_free(A);
00788 A=NULL;
00789 mask=NULL;
00790 
00791 return ret;
00792 }

string EMAN::SVDAnalyzer::get_name (  )  const [inline, virtual]

Get the Analyzer's name.

Each Analyzer is identified by a unique name.

Returns:
The Analyzer's name.

Implements EMAN::Analyzer.

Definition at line 219 of file analyzer.h.

00220                 {
00221                         return "svd_gsl";
00222                 }

string EMAN::SVDAnalyzer::get_desc (  )  const [inline, virtual]

Get the Analyzer's description.

Returns:
The Analyzer's description.

Implements EMAN::Analyzer.

Definition at line 224 of file analyzer.h.

00225                 {
00226                         return "Singular Value Decomposition from GSL. Comparable to pca";
00227                 }

static Analyzer* EMAN::SVDAnalyzer::NEW (  )  [inline, static]

Definition at line 229 of file analyzer.h.

References SVDAnalyzer().

Referenced by EMAN::Factory< T >::Factory().

00230                 {
00231                         return new SVDAnalyzer();
00232                 }

void SVDAnalyzer::set_params ( const Dict new_params  )  [virtual]

Set the Analyzer parameters using a key/value dictionary.

Parameters:
new_params A dictionary containing the new parameters.

Reimplemented from EMAN::Analyzer.

Definition at line 794 of file analyzer.cpp.

References A, EMAN::EMData::get_data(), EMAN::EMData::get_xsize(), EMAN::EMData::get_ysize(), EMAN::EMData::get_zsize(), mask, nimg, nsofar, nvec, EMAN::Analyzer::params, and pixels.

00795 {
00796         params = new_params;
00797         mask = params["mask"];
00798         nvec = params["nvec"];
00799         nimg = params["nimg"];
00800 
00801         // count pixels under mask
00802         pixels=0;
00803         size_t totpix=mask->get_xsize()*mask->get_ysize()*mask->get_zsize();
00804         float *d=mask->get_data();
00805         for (size_t i=0; i<totpix; ++i) if (d[i]) ++pixels;
00806 
00807         printf("%d,%d\n",pixels,nimg);
00808         A=gsl_matrix_alloc(pixels,nimg);
00809         nsofar=0;
00810 }

TypeDict EMAN::SVDAnalyzer::get_param_types (  )  const [inline, virtual]

Get Analyzer parameter information in a dictionary.

Each parameter has one record in the dictionary. Each record contains its name, data-type, and description.

Returns:
A dictionary containing the parameter info.

Implements EMAN::Analyzer.

Definition at line 236 of file analyzer.h.

References EMAN::EMObject::EMDATA, EMAN::EMObject::INT, and EMAN::TypeDict::put().

00237                 {
00238                         TypeDict d;
00239                         d.put("mask", EMObject::EMDATA, "mask image");
00240                         d.put("nvec", EMObject::INT, "number of desired basis vectors");
00241                         d.put("nimg", EMObject::INT, "total number of input images, required even with insert_image()");
00242                         return d;
00243                 }


Member Data Documentation

Definition at line 246 of file analyzer.h.

Referenced by analyze(), insert_image(), and set_params().

int EMAN::SVDAnalyzer::nvec [protected]

Definition at line 247 of file analyzer.h.

Referenced by analyze(), and set_params().

int EMAN::SVDAnalyzer::pixels [protected]

Definition at line 248 of file analyzer.h.

Referenced by set_params().

int EMAN::SVDAnalyzer::nimg [protected]

Definition at line 249 of file analyzer.h.

Referenced by analyze(), and set_params().

Definition at line 252 of file analyzer.h.

Referenced by insert_image(), and set_params().

gsl_matrix* EMAN::SVDAnalyzer::A [private]

Definition at line 253 of file analyzer.h.

Referenced by analyze(), insert_image(), and set_params().


The documentation for this class was generated from the following files:

Generated on Sat Nov 21 02:20:15 2009 for EMAN2 by  doxygen 1.5.6