openModeller  Version 1.5.0
csm.hh
Go to the documentation of this file.
1 
35 #ifndef CSM_H
36 #define CSM_H
37 
38 #include <openmodeller/om.hh>
39 #include <gsl/gsl_matrix.h>
40 
49 
51 
221 
223 
269 class Csm : public AlgorithmImpl
270 {
271  public:
273  Csm(AlgMetadata const* metadata);
275  ~Csm();
276 
277 
278 
279  //
280  // Methods used to build the model
281  //
282 
283 
289  virtual int initialize();
290 
295  int iterate();
296 
303  int done() const;
304 
305  //
306  // Methods used to project the model
307  //
308 
309 
316  Scalar getValue( const Sample& x ) const;
317 
324  int getConvergence( Scalar * const val ) const;
325 
326 
327 
328  private:
329 
330  protected:
334  int SamplerToMatrix();
335 
336 
339  bool csm1();
340 
349  int calculateMeanAndSd( gsl_matrix * theMatrix,
350  gsl_vector * theMeanVector,
351  gsl_vector * theStdDevVector);
352 
359  int center();
360 
368  virtual int discardComponents()=0;
369 
370 
371 
377  void displayVector(const gsl_vector * v, const char * name,const bool roundFlag=true) const;
378 
384  void displayMatrix(const gsl_matrix * m, const char * name, const bool roundFlag=true) const;
389  gsl_matrix * transpose (gsl_matrix * m);
390 
396  double product (gsl_vector * va, gsl_vector * vb) const;
397 
403  gsl_matrix * product (gsl_matrix * a, gsl_matrix * b) const;
404 
409  gsl_matrix * autoCovariance(gsl_matrix * m);
410 
414  virtual void _getConfiguration( ConfigurationPtr& config ) const;
415 
419  virtual void _setConfiguration( const ConstConfigurationPtr& config );
420 
421 
426  int _done;
436  gsl_vector * _gsl_avg_vector;
438  gsl_vector * _gsl_stddev_vector;
449 
454 
455 };
456 
457 #endif
Scalar getValue(const Sample &x) const
Definition: csm.cpp:283
gsl_matrix * _gsl_eigenvector_matrix
Definition: csm.hh:442
bool csm1()
Definition: csm.cpp:648
double Scalar
Type of map values.
Definition: om_defs.hh:39
int done() const
Definition: csm.cpp:268
virtual int initialize()
Definition: csm.cpp:98
int center()
Definition: csm.cpp:217
int getConvergence(Scalar *const val) const
Definition: csm.cpp:403
virtual void _getConfiguration(ConfigurationPtr &config) const
Definition: csm.cpp:723
int calculateMeanAndSd(gsl_matrix *theMatrix, gsl_vector *theMeanVector, gsl_vector *theStdDevVector)
Definition: csm.cpp:175
int _localityCount
Definition: csm.hh:448
virtual void _setConfiguration(const ConstConfigurationPtr &config)
Definition: csm.cpp:770
bool verboseDebuggingBool
Definition: csm.hh:453
int _retained_components_count
Definition: csm.hh:446
virtual int discardComponents()=0
Csm(AlgMetadata const *metadata)
Definition: csm.cpp:57
gsl_vector * _gsl_eigenvalue_vector
Definition: csm.hh:440
~Csm()
Definition: csm.cpp:70
gsl_vector * _gsl_avg_vector
Definition: csm.hh:436
gsl_vector * _gsl_stddev_vector
Definition: csm.hh:438
void displayMatrix(const gsl_matrix *m, const char *name, const bool roundFlag=true) const
Definition: csm.cpp:451
int _layer_count
Definition: csm.hh:444
Definition: csm.hh:269
gsl_matrix * autoCovariance(gsl_matrix *m)
Definition: csm.cpp:578
void displayVector(const gsl_vector *v, const char *name, const bool roundFlag=true) const
Definition: csm.cpp:413
int SamplerToMatrix()
Definition: csm.cpp:142
int _initialized
Definition: csm.hh:423
gsl_matrix * _gsl_covariance_matrix
Definition: csm.hh:433
gsl_matrix * transpose(gsl_matrix *m)
Definition: csm.cpp:493
double product(gsl_vector *va, gsl_vector *vb) const
Definition: csm.cpp:510
int minComponentsInt
Definition: csm.hh:451
gsl_matrix * _gsl_environment_matrix
Definition: csm.hh:430
int _done
Definition: csm.hh:426
int iterate()
Definition: csm.cpp:255
AlgMetadata metadata
Definition: garp.cpp:134
Definition: Sample.hh:25