openModeller  Version 1.5.0
enfa.hh
Go to the documentation of this file.
1 
35 #ifndef ENFA_H
36 #define ENFA_H
37 
38 #include <openmodeller/om.hh>
40 #include <gsl/gsl_matrix.h>
41 
49 
51 
123 
125 
215 class Enfa : public AlgorithmImpl
216 {
217 public:
219  Enfa();
221  ~Enfa();
222 
223  //
224  // Methods used to build the model
225  //
226 
232  virtual int initialize();
233 
238  int iterate();
239 
246  int done() const;
247 
248  //
249  // Methods used to project the model
250  //
251 
252 
259  Scalar getValue( const Sample& x ) const;
260 
267  int getConvergence( Scalar * const val ) const;
268 
269 
270 
271 private:
272 
273 protected:
277  int SamplerToMatrix();
278 
282  int BackgroundToMatrix();
283 
286  bool enfa1();
287 
292  int calculateMeanAndSd( gsl_matrix * theMatrix,
293  gsl_vector * theMeanVector,
294  gsl_vector * theStdDevVector);
295 
302  int center(gsl_matrix * theMatrix,
303  int spCount);
304 
307  gsl_matrix* sqrtm(gsl_matrix* original_matrix) const;
308 
311  gsl_matrix* inverse(gsl_matrix* _m) const;
312 
316  double getGeomean(gsl_vector* v) const;
317 
326  int discardComponents() const;
327 
333  void displayVector(const gsl_vector * v, const char * name,const bool roundFlag=true) const;
334 
340  void displayMatrix(const gsl_matrix * m, const char * name, const bool roundFlag=true) const;
341 
347  void displayLoadings(const gsl_matrix * m, const int f) const;
348 
353  gsl_matrix * autoCovariance(gsl_matrix * m);
354 
358  virtual void _getConfiguration( ConfigurationPtr& ) const;
359 
363  virtual void _setConfiguration( const ConstConfigurationPtr& );
364 
365 
367  int _initialized;
370  int _done;
374  gsl_matrix * _gsl_environment_matrix;
375  /* Pointer to gsl matrix with environmental data for each locality
376  converted to factors using the enfa score matrix */
377  gsl_matrix * _gsl_environment_factor_matrix;
378 
382  gsl_matrix * _gsl_background_matrix;
383 
386  gsl_matrix * _gsl_covariance_matrix;
387 
390  gsl_matrix * _gsl_covariance_background_matrix;
391 
394  gsl_vector * _gsl_avg_vector;
396  gsl_vector * _gsl_stddev_vector;
397 
400  gsl_vector * _gsl_avg_background_vector;
401 
404  gsl_vector * _gsl_stddev_background_vector;
405 
407  gsl_vector * _gsl_eigenvalue_vector;
409  gsl_matrix * _gsl_eigenvector_matrix;
412  gsl_matrix * _gsl_score_matrix ;
415  gsl_matrix * _gsl_covariance_matrix_root_inverse;
417  gsl_matrix * _gsl_workspace_H;
419  gsl_matrix * _gsl_workspace_W;
421  gsl_matrix * _gsl_workspace_y;
423  gsl_vector * _gsl_workspace_z;
425  gsl_vector* _gsl_factor_weights_all_components;
427  gsl_vector* _gsl_factor_weights;
430  gsl_vector* _gsl_geomean_vector;
431 
432  int _layer_count;
434  int _retained_components_count;
436  int _localityCount;
437 
439  int _backgroundCount;
440 
442  int _backgroundProvided;
443 
445  int _discardMethod;
446 
448  int _retainComponents;
449 
451  double _retainVariation;
452 
454  int minComponentsInt;
455 
457  int _verboseDebug;
458 
460  int _displayLoadings;
461 
463  int _gsl_vector_min;
464 
466  double _marginality;
467 
469  double _specialisation;
470 
472  int _numRetries;
473 
475  int _retryCount;
476 
477 
478 };
479 
480 /* catch the case that generating an inverse fails
481  gsl doesnt do this well, and we want to catch,
482  ignore and retry when this happens */
483 class InverseFailedException : public OmException {
484 public:
485  InverseFailedException( const std::string& msg ) :
486  OmException( msg )
487  {}
488 };
489 
490 #endif
virtual int initialize()=0
double Scalar
Type of map values.
Definition: om_defs.hh:39
virtual int done() const
Definition: Algorithm.hh:140
virtual int iterate()
Definition: Algorithm.hh:132
virtual void _setConfiguration(const ConstConfigurationPtr &)
Definition: Algorithm.hh:200
virtual void _getConfiguration(ConfigurationPtr &) const
Definition: Algorithm.hh:199
virtual int getConvergence(Scalar *const val) const
Definition: Algorithm.hh:143
virtual Scalar getValue(const Sample &x) const =0
Definition: Sample.hh:25