openModeller  Version 1.4.0
csm.hh
Go to the documentation of this file.
00001 
00035 #ifndef CSM_H
00036 #define CSM_H
00037 
00038 #include <openmodeller/om.hh>
00039 #include <gsl/gsl_matrix.h>
00040 
00049 
00051 
00221 
00223 
00269 class Csm : public AlgorithmImpl
00270 {
00271     public:
00273         Csm(AlgMetadata const* metadata);
00275         ~Csm();
00276 
00277         
00278 
00279         //
00280         // Methods used to build the model
00281         //
00282 
00283 
00289         virtual int initialize();
00290 
00295         int iterate();
00296 
00303         int done() const;
00304 
00305         //
00306         // Methods used to project the model
00307         //
00308 
00309 
00316         Scalar getValue( const Sample& x ) const;
00317 
00324         int getConvergence( Scalar * const val ) const;
00325 
00326 
00327 
00328     private:
00329 
00330     protected:  
00334         int SamplerToMatrix();
00335 
00336 
00339         bool csm1();
00340 
00349         int calculateMeanAndSd( gsl_matrix * theMatrix, 
00350                 gsl_vector * theMeanVector,
00351                 gsl_vector * theStdDevVector);
00352 
00359         int center();
00360 
00368         virtual int discardComponents()=0;
00369 
00370 
00371 
00377         void displayVector(const gsl_vector * v, const char * name,const  bool roundFlag=true) const;
00378 
00384         void displayMatrix(const gsl_matrix * m, const char * name, const bool roundFlag=true) const;
00389         gsl_matrix * transpose (gsl_matrix * m);
00390 
00396         double product (gsl_vector * va, gsl_vector * vb) const;
00397 
00403         gsl_matrix * product (gsl_matrix * a, gsl_matrix * b) const;
00404 
00409         gsl_matrix * autoCovariance(gsl_matrix * m);
00410 
00414         virtual void _getConfiguration( ConfigurationPtr& config ) const;
00415 
00419         virtual void _setConfiguration( const ConstConfigurationPtr& config );
00420 
00421 
00423         int _initialized;
00426         int _done;
00430         gsl_matrix * _gsl_environment_matrix;
00433         gsl_matrix * _gsl_covariance_matrix;
00436         gsl_vector * _gsl_avg_vector;
00438         gsl_vector * _gsl_stddev_vector;
00440         gsl_vector * _gsl_eigenvalue_vector;
00442         gsl_matrix * _gsl_eigenvector_matrix;
00444         int _layer_count; 
00446         int _retained_components_count;
00448         int _localityCount; 
00449         
00451         int minComponentsInt;
00453         bool verboseDebuggingBool;
00454 
00455 };
00456 
00457 #endif